OB.DAAC Logo
NASA Logo
Ocean Color Science Software

ocssw V2022
viirs_ler_luts_nc4.f95
Go to the documentation of this file.
1 module viirs_ler_luts
2 
3 !private
4 
5 public :: load_viirs_ler_luts
6 !public :: unload_viirs_ler_luts
7 
8 real, dimension(20800) :: logi0
9 real, dimension(20800) :: z1i0
10 real, dimension(20800) :: z2i0
11 real, dimension(20800) :: ti0
12 real, dimension(260) :: sb
13 real, dimension(160) :: li0r
14 real, dimension(160) :: z1i0r
15 real, dimension(160) :: z2i0r
16 real, dimension(160) :: ti0r
17 real, dimension(2) :: sbr
18 
19 real, dimension(10) :: xxzlog = (/0.0000000, 0.00977964, 0.0395086, 0.0904221, &
20 & 0.164818, 0.266515, 0.401776, 0.581261, &
21 & 0.824689, 1.17436/)
22 ! -- these numbers correspond to satellite zenith angle
23 ! -- node points of 0.,16.,30.,40.,48.,54.,58.,60. degrees
24 real, dimension(8) :: xxlog = (/0.000000,0.0395086,0.143841,0.266515,0.401776, &
25 & 0.531394,0.635031,0.693147/)
26 
27 real, dimension(10) :: xzlog
28 real, dimension(8) :: xlog
29 real, dimension(4,7) :: densol
30 real, dimension(4,5) :: denscn
31 real, dimension(4) :: cthet0
32 real, dimension(4) :: ctheta
33 real, dimension(16) :: cofs
34 integer :: indsol
35 integer :: indscn
36 integer :: iofset
37 real :: p1
38 real :: pr
39 real, dimension(10) :: dum10
40 
41 common /lpoly/ xzlog, xlog, densol, denscn, cthet0, ctheta, cofs, &
42 & indsol, indscn, iofset, p1, pr, dum10
43 
44 real :: xlat
45 real :: xlong
46 real :: sza
47 real :: xthet
48 real :: xphi
49 real :: cphi
50 real :: c2phi
51 real :: pteran
52 real :: xnvalm
53 integer :: isnow
54 integer :: partial
55 real :: so2ind
56 real :: resn
57 real :: sens
58 real :: rsens
59 real :: ozbst
60 real :: ref
61 real :: estozn
62 real :: ozcld
63 real :: pcloud
64 real :: prfrac
65 real :: clfrac
66 real :: rayval
67 real :: r412
68 real :: r470
69 real :: sfref412
70 real :: sfref470
71 real :: sfref650
72 real :: qdif412
73 real :: qdif470
74 real :: qdif650
75 real :: stdv
76 
77 common /sample/ xlat,xlong,sza,xthet,xphi,cphi,c2phi,pteran, &
78 & xnvalm(6),isnow,so2ind,resn(5),sens(5), &
79 & rsens,ozbst,ref,estozn,ozcld,pcloud,prfrac,clfrac,partial, &
80 & rayval(5),r412,r470,sfref412,sfref470,sfref650,qdif412, &
81 & qdif470,qdif650, stdv
82 
83 real, dimension(26) :: realbuf
84 common /bufout/ realbuf
85 
86 common /table/ logi0, z1i0, z2i0, ti0, sb, li0r, z1i0r, z2i0r, ti0r, sbr
87 
88 contains
89 
90 subroutine load_viirs_ler_luts(lut_ler_file, status)
91  use netcdf
92 
93  implicit none
94 
95  character(len=255), intent(in) :: lut_ler_file
96 
97  real :: xdenom
98  integer :: i, j, k
99  integer :: status
100 
101  character(len=255) :: sds_name
102  integer :: sd_id, sds_id, sds_index
103  integer, dimension(1) :: start1, stride1, dims1
104  integer, dimension(32):: dimids
105  integer :: rank, dtype, nattrs
106 
107  integer :: nc_id
108  integer :: dim_id
109  integer :: dset_id
110  integer :: grp_id
111 
112  character(len=255) :: dset_name
113  character(len=255) :: attr_name
114  character(len=255) :: group_name
115  character(len=255) :: err_msg
116 
117  status = nf90_open(trim(lut_ler_file), nf90_nowrite, nc_id)
118  if (status /= nf90_noerr) then
119  print *, "ERROR: Failed to open deepblue lut_nc4 file: ", status
120  return
121  end if
122 
123  group_name = 'LER_TABLES'
124  status = nf90_inq_ncid(nc_id, group_name, grp_id)
125  if (status /= nf90_noerr) then
126  print *, "ERROR: Failed to get ID of group "//trim(group_name)//": ", status
127  return
128  end if
129 
130  start1=(/1/)
131  stride1=(/1/)
132 
133  dset_name = 'LOGI0'
134  status = nf90_inq_varid(grp_id, dset_name, dset_id)
135  if (status /= nf90_noerr) then
136  print *, "ERROR: Failed to get ID of dataset "//trim(dset_name)//": ", status
137  return
138  end if
139  status = nf90_inquire_variable(grp_id, dset_id, dimids=dimids)
140  status = nf90_inquire_dimension(grp_id, dimids(1), len = dims1(1))
141  status = nf90_get_var(grp_id, dset_id, logi0, start=start1, &
142  stride=stride1, count=dims1)
143  err_msg = nf90_strerror(status)
144  if (status /= nf90_noerr) then
145  print *, "ERROR: Failed to read dataset "//trim(dset_name)//": ", status
146  return
147  end if
148 
149  dset_name = 'Z1I0'
150  status = nf90_inq_varid(grp_id, dset_name, dset_id)
151  if (status /= nf90_noerr) then
152  print *, "ERROR: Failed to get ID of dataset "//trim(dset_name)//": ", status
153  return
154  end if
155  status = nf90_inquire_variable(grp_id, dset_id, dimids=dimids)
156  status = nf90_inquire_dimension(grp_id, dimids(1), len = dims1(1))
157  status = nf90_get_var(grp_id, dset_id, z1i0, start=start1, &
158  stride=stride1, count=dims1)
159  err_msg = nf90_strerror(status)
160  if (status /= nf90_noerr) then
161  print *, "ERROR: Failed to read dataset "//trim(dset_name)//": ", status
162  return
163  end if
164 
165  dset_name = 'Z2I0'
166  status = nf90_inq_varid(grp_id, dset_name, dset_id)
167  if (status /= nf90_noerr) then
168  print *, "ERROR: Failed to get ID of dataset "//trim(dset_name)//": ", status
169  return
170  end if
171  status = nf90_inquire_variable(grp_id, dset_id, dimids=dimids)
172  status = nf90_inquire_dimension(grp_id, dimids(1), len = dims1(1))
173  status = nf90_get_var(grp_id, dset_id, z2i0, start=start1, &
174  stride=stride1, count=dims1)
175  err_msg = nf90_strerror(status)
176  if (status /= nf90_noerr) then
177  print *, "ERROR: Failed to read dataset "//trim(dset_name)//": ", status
178  return
179  end if
180 
181  dset_name = 'TI0'
182  status = nf90_inq_varid(grp_id, dset_name, dset_id)
183  if (status /= nf90_noerr) then
184  print *, "ERROR: Failed to get ID of dataset "//trim(dset_name)//": ", status
185  return
186  end if
187  status = nf90_inquire_variable(grp_id, dset_id, dimids=dimids)
188  status = nf90_inquire_dimension(grp_id, dimids(1), len = dims1(1))
189  status = nf90_get_var(grp_id, dset_id, ti0, start=start1, &
190  stride=stride1, count=dims1)
191  err_msg = nf90_strerror(status)
192  if (status /= nf90_noerr) then
193  print *, "ERROR: Failed to read dataset "//trim(dset_name)//": ", status
194  return
195  end if
196 
197  dset_name = 'SB'
198  status = nf90_inq_varid(grp_id, dset_name, dset_id)
199  if (status /= nf90_noerr) then
200  print *, "ERROR: Failed to get ID of dataset "//trim(dset_name)//": ", status
201  return
202  end if
203  status = nf90_inquire_variable(grp_id, dset_id, dimids=dimids)
204  status = nf90_inquire_dimension(grp_id, dimids(1), len = dims1(1))
205  status = nf90_get_var(grp_id, dset_id, sb, start=start1, &
206  stride=stride1, count=dims1)
207  err_msg = nf90_strerror(status)
208  if (status /= nf90_noerr) then
209  print *, "ERROR: Failed to read dataset "//trim(dset_name)//": ", status
210  return
211  end if
212 
213  dset_name = 'LOGI0R'
214  status = nf90_inq_varid(grp_id, dset_name, dset_id)
215  if (status /= nf90_noerr) then
216  print *, "ERROR: Failed to get ID of dataset "//trim(dset_name)//": ", status
217  return
218  end if
219  status = nf90_inquire_variable(grp_id, dset_id, dimids=dimids)
220  status = nf90_inquire_dimension(grp_id, dimids(1), len = dims1(1))
221  status = nf90_get_var(grp_id, dset_id, li0r, start=start1, &
222  stride=stride1, count=dims1)
223  err_msg = nf90_strerror(status)
224  if (status /= nf90_noerr) then
225  print *, "ERROR: Failed to read dataset "//trim(dset_name)//": ", status
226  return
227  end if
228 
229  dset_name = 'Z1I0R'
230  status = nf90_inq_varid(grp_id, dset_name, dset_id)
231  if (status /= nf90_noerr) then
232  print *, "ERROR: Failed to get ID of dataset "//trim(dset_name)//": ", status
233  return
234  end if
235  status = nf90_inquire_variable(grp_id, dset_id, dimids=dimids)
236  status = nf90_inquire_dimension(grp_id, dimids(1), len = dims1(1))
237  status = nf90_get_var(grp_id, dset_id, z1i0r, start=start1, &
238  stride=stride1, count=dims1)
239  err_msg = nf90_strerror(status)
240  if (status /= nf90_noerr) then
241  print *, "ERROR: Failed to read dataset "//trim(dset_name)//": ", status
242  return
243  end if
244 
245  dset_name = 'Z2I0R'
246  status = nf90_inq_varid(grp_id, dset_name, dset_id)
247  if (status /= nf90_noerr) then
248  print *, "ERROR: Failed to get ID of dataset "//trim(dset_name)//": ", status
249  return
250  end if
251  status = nf90_inquire_variable(grp_id, dset_id, dimids=dimids)
252  status = nf90_inquire_dimension(grp_id, dimids(1), len = dims1(1))
253  status = nf90_get_var(grp_id, dset_id, z2i0r, start=start1, &
254  stride=stride1, count=dims1)
255  err_msg = nf90_strerror(status)
256  if (status /= nf90_noerr) then
257  print *, "ERROR: Failed to read dataset "//trim(dset_name)//": ", status
258  return
259  end if
260 
261  dset_name = 'TI0R'
262  status = nf90_inq_varid(grp_id, dset_name, dset_id)
263  if (status /= nf90_noerr) then
264  print *, "ERROR: Failed to get ID of dataset "//trim(dset_name)//": ", status
265  return
266  end if
267  status = nf90_inquire_variable(grp_id, dset_id, dimids=dimids)
268  status = nf90_inquire_dimension(grp_id, dimids(1), len = dims1(1))
269  status = nf90_get_var(grp_id, dset_id, ti0r, start=start1, &
270  stride=stride1, count=dims1)
271  err_msg = nf90_strerror(status)
272  if (status /= nf90_noerr) then
273  print *, "ERROR: Failed to read dataset "//trim(dset_name)//": ", status
274  return
275  end if
276 
277  dset_name = 'SBR'
278  status = nf90_inq_varid(grp_id, dset_name, dset_id)
279  if (status /= nf90_noerr) then
280  print *, "ERROR: Failed to get ID of dataset "//trim(dset_name)//": ", status
281  return
282  end if
283  status = nf90_inquire_variable(grp_id, dset_id, dimids=dimids)
284  status = nf90_inquire_dimension(grp_id, dimids(1), len = dims1(1))
285  status = nf90_get_var(grp_id, dset_id, sbr, start=start1, &
286  stride=stride1, count=dims1)
287  err_msg = nf90_strerror(status)
288  if (status /= nf90_noerr) then
289  print *, "ERROR: Failed to read dataset "//trim(dset_name)//": ", status
290  return
291  end if
292 
293  status = nf90_close(nc_id)
294  if (status /= nf90_noerr) then
295  print *, "ERROR: Failed to close lut_nc4 file: ", status
296  return
297  end if
298 
299  xzlog(:) = xxzlog(:)
300  xlog(:) = xxlog(:)
301 
302 ! -- calculate values needed in table interpolation
303  do j = 1, 7
304  do k = j, j + 3
305  xdenom = 1.0
306  do i = j, j + 3
307  if (i .ne. k) xdenom = xdenom * (xzlog(k) - xzlog(i))
308  end do
309  densol(k - j + 1, j) = xdenom
310  end do
311  end do
312  do j = 1, 5
313  do k = j, j + 3
314  xdenom = 1.0
315  do i = j, j + 3
316  if (i .ne. k) xdenom = xdenom * (xlog(k) - xlog(i))
317  end do
318  denscn(k - j + 1, j) = xdenom
319  end do
320  end do
321 
322  return
323 end subroutine load_viirs_ler_luts
324 
325 integer function wrap_total() result(status)
326  implicit none
327 
328  status = 0
329  return
330 
331 end function wrap_total
332 
333 end module viirs_ler_luts
an array had not been initialized Several spelling and grammar corrections were which is read from the appropriate MCF the above metadata values were hard coded A problem calculating the average background DN for SWIR bands when the moon is in the space view port was corrected The new algorithm used to calculate the average background DN for all reflective bands when the moon is in the space view port is now the same as the algorithm employed by the thermal bands For non SWIR changes in the averages are typically less than Also for non SWIR the black body DNs remain a backup in case the SV DNs are not available For SWIR the changes in computed averages were larger because the old which used the black body suffered from contamination by the micron leak As a consequence of the if SV DNs are not available for the SWIR the EV pixels will not be the granule time is used to identify the appropriate tables within the set given for one LUT table(e.g. m1) to be used for interpolation. The table values will be linearly interpolated using the tables corresponding to the node times bracketing the granule time. If the granule time falls before the time of the first node or after the time of the last node
string & trim(string &s, const string &delimiters)
Definition: EnvsatUtil.cpp:29
subroutine, public load_viirs_ler_luts(lut_ler_file, status)
integer function wrap_total()