OB.DAAC Logo
NASA Logo
Ocean Color Science Software

ocssw V2022
interp_hycom.f
Go to the documentation of this file.
1  program interp_hycom
2 
3  use netcdf
4  USE hdf5
5 c
6 c This program interpolate a global 1/12 degree hycom temperature
7 c or salinity field at selected depth into a regular grid of
8 c a quarter degree (1440x720).
9 c
10 c This program is based on isubaregion.f written by Alan J. Wallcraft
11 c It takes two pre-generated NetCDF files: hycom_weight.nc
12 c that stores the interpolation weight table from the hycom grid
13 c to the regular grid based on a bilinear interpolation
14 c The weight table was generated using a modified program based
15 c on isuba_gmapi.f also written by Wallcraft.
16 c
17 c The second file is mask1.nc. It provides the landmask for the
18 c output grid. This two files have to reside in the directory
19 c where the program is executed.
20 c
21 c The input file is a NetCDF file extracted from the HYCOM output
22 c file via their OpenDAP server. The variable field contains only
23 c one depth only. If you change the input file format, you will
24 c have to modify the code to read it properly.
25 c
26 c The output file is a NetCDF file containing 1D arrays for latitude
27 c and longitude, two 3D arrays for temp and salt. There are
28 c 9 pre-defined depth at 0, 10, 20, 30, 50, 100, 250, 500, 1000 meter
29 c below the sea surface.
30 c
31 c Original Author:
32 c Alan J. Wallcraft, NRL, January 2005 and January 2009.c
33 c Modified by
34 c Peggy Li, JPL, August 31, 2010
35 c Modified by
36 c Joel Gales, GSFC, September 5, 2010
37 c Add HDF5 and binary output
38 c
39 c
40  character*40 :: field, timetag, depth, varname
41  character*128 :: dir_in, dir_out
42  character*256 :: file_in, file_out
43  integer :: nargs,iargc
44  integer :: level
45  integer :: idm, jdm,idm_out,jdm_out,ddim
46  integer :: i,ii,ip,j,jj,jp
47  integer :: if_sm,il_sm,jf_sm,jl_sm
48  integer :: status,ncid,ncid2,vid
49  integer :: latid, lonid, depthid
50  integer :: varid, xoutid, youtid
51  real, allocatable :: lats(:), lons(:), dep(:)
52  character*80 :: attstr
53  integer :: start_3d(4),count_3d(4)
54 ! logical :: lmask
55  integer, allocatable :: m_sm(:,:), iv_sm(:,:)
56  integer, allocatable :: m_out(:,:), m_osm(:,:)
57  real, allocatable :: a_in(:,:)
58  real, allocatable, target :: a_out(:,:)
59 c
60  integer, allocatable :: i_out(:,:),j_out(:,:)
61  real, allocatable :: x_out(:,:),y_out(:,:)
62 c
63  real, parameter :: spval=2.0**100 ! spval
64  real, parameter :: hspval=0.5*2.0**100 ! half spval
65 ! DOUBLE PRECISION :: time
66  logical :: file_exists
67  integer :: count, count1
68 
69  INTEGER(HID_T) :: file_id
70  INTEGER(HID_T) :: dset_id
71  INTEGER :: rank
72  INTEGER(HSIZE_T), DIMENSION(2) :: dims_sal
73  INTEGER(HSIZE_T), DIMENSION(3) :: offset=(/0,0,0/)
74  INTEGER(HID_T) :: filespace
75  INTEGER(HID_T) :: dataspace
76  real(4), POINTER :: pt
77 
78  integer slash_pos
79  CHARACTER datadir*255
80 
81  include "aquarius.fin"
82 c
83 c define input/output dimension
84 c parse input arguments
85  idm=4500;jdm=3298
86  idm_out=1440;jdm_out=720
87  nargs=iargc()
88  if (nargs.eq.2) then
89 ! call getarg(1,field)
90 ! call getarg(2,timetag)
91 ! call getarg(2,depth)
92  call getarg(1,file_in)
93  call getarg(2,dir_out)
94  else
95  print *,
96 ! & "Command:get_hycom [temp|salt] yyyymmdd depth infile outdir"
97  & "Command:get_hycom [temp|salt] infile outdir"
98 ! print *," depth: 0,10,20,30,50,100,250,500,and 1000"
99 ! stop
100  call exit(1)
101  endif
102  field = 'salt'
103  depth = '0'
104 
105  slash_pos = index(file_in,'/',back=.true.)
106  timetag = file_in(slash_pos+7:slash_pos+10)//file_in(slash_pos+12:slash_pos+14)
107 
108  if (field .eq. "temp") then
109  varname="temperature"
110  elseif (field .eq. "salt") then
111  varname="salinity"
112  else
113  print *, "field has to be either temp or salt"
114 ! stop
115  call exit(1)
116  endif
117  level = -1
118  if (depth .eq. '0') level=1
119  if (depth .eq. '10') level=2
120  if (depth .eq. '20') level=3
121  if (depth .eq. '30') level=4
122  if (depth .eq. '50') level=5
123  if (depth .eq. '100') level=6
124  if (depth .eq. '250') level=7
125  if (depth .eq. '500') level=8
126  if (depth .eq. '1000') level=9
127  if (level < 0) then
128  print *, "depth value is invalid"
129 ! stop
130  call exit(1)
131  endif
132 
133  allocate( iv_sm(jdm,2) )
134 c
135  allocate( m_sm(idm,jdm), m_osm(idm_out,jdm_out) )
136  allocate( m_out(idm_out,jdm_out) )
137  allocate( a_in(idm,jdm), a_out(idm_out,jdm_out) )
138 c
139  allocate( i_out(idm_out,jdm_out), j_out(idm_out,jdm_out) )
140  allocate( x_out(idm_out,jdm_out), y_out(idm_out,jdm_out) )
141 
142  a_out(:,:)=spval
143 ! write(*,*) 'interpolating ',trim(file_in),' to ',
144 ! & trim(file_out)
145 
146  call getenv('OCDATAROOT', datadir)
147 
148 c readin the land mask for the output grid
149  status=nf90_open(trim(datadir)//'/aquarius/static/hycom_landmask.nc', 0, ncid)
150  if (status /= nf90_noerr) then
151  write(*,*) 'Land mask "hycom_landmask.nc" cannot be opened.'
152  call handle_err(status)
153  endif
154 
155  status=nf90_inq_varid(ncid,"mask",vid)
156  if (status /= nf90_noerr) then
157  call handle_err(status)
158  endif
159 
160  status=nf90_get_var(ncid,vid,m_out)
161  if (status /= nf90_noerr) then
162  call handle_err(status)
163  endif
164 
165  status=nf90_close(ncid)
166 
167 
168 c stop if input file does not exit
169 c and create the output file if not exist
170  inquire(file=file_in, exist=file_exists)
171  if (.not. file_exists) then
172  print *, file_in, " does not exist"
173 ! stop
174  call exit(1)
175  endif
176 ! inquire(FILE=file_out, EXIST=file_exists)
177 ! if (.not. file_exists) then
178 
179  file_out=trim(dir_out)//'/N'//trim(timetag)//'00_SALINITY_HYCOM_24h.nc'
180 
181 ! print *, "Create ",trim(file_out)
182  ddim = 9
183 
184  status=nf90_create(file_out,0, ncid)
185  if (status /= nf90_noerr) then
186  write(*,*) 'NetCDF output file cannot be created.'
187  call handle_err(status)
188  endif
189 
190  status=nf90_def_dim(ncid,'lat', jdm_out, latid)
191  if (status /= nf90_noerr) call handle_err(status)
192  status=nf90_def_dim(ncid,'lon', idm_out, lonid)
193  if (status /= nf90_noerr) call handle_err(status)
194 
195  status=nf90_def_var(ncid, 'salt', nf90_float,
196  & (/lonid,latid/), varid)
197  if (status /= nf90_noerr) call handle_err(status)
198  attstr="Salinity"
199  status=nf90_put_att(ncid,varid,'long_name',attstr)
200  attstr="PS"
201  status=nf90_put_att(ncid,varid,'units',attstr)
202  if (status /= nf90_noerr) call handle_err(status)
203  status=nf90_put_att(ncid,varid,'_FillValue',spval)
204  if (status /= nf90_noerr) call handle_err(status)
205  status=nf90_enddef(ncid)
206  if (status /= nf90_noerr) call handle_err(status)
207 
208  allocate(lats(jdm_out), lons(idm_out), dep(ddim))
209 
210  do i=1,idm_out
211  lons(i)=-180.125+i*0.25
212  enddo
213  do i=1,jdm_out
214  lats(i)=-90.125+i*0.25
215  enddo
216  dep(1)=0; dep(2)=10; dep(3)=20; dep(4)=30
217  dep(5)=50; dep(6)=100; dep(7)=250; dep(8)=500
218  dep(9)=1000
219 
220  deallocate(lats,lons,dep)
221 
222 
223 c Open input file and read in the field
224  status=nf90_open(file_in, nf90_nowrite, ncid)
225  if (status /= nf90_noerr) call handle_err(status)
226  status=nf90_inq_varid(ncid,varname,vid)
227  status=nf90_get_var(ncid,vid,a_in)
228  if (status /= nf90_noerr) call handle_err(status)
229  status=nf90_close(ncid)
230  if (status /= nf90_noerr) call handle_err(status)
231 
232  count = 0
233  do j= 1,jdm
234  do i= 1,idm
235  if (a_in(i,j).gt.hspval) then
236  count = count+1
237  m_sm(i,j) = 0
238  else
239  m_sm(i,j) = 1
240  endif
241  enddo
242  enddo
243  print *, "hspval", hspval
244  print *, "Input data ", count, " land point"
245 
246 c read in the weight file and output mask
247 
248  status=nf90_open(trim(datadir)//'/aquarius/static/hycom_weight.nc',0, ncid)
249 
250  if (status /= nf90_noerr) then
251  write(*,*) 'HYCOM weight file "hycom_weight.nc" cannot be opened.'
252  call handle_err(status)
253  endif
254 
255  status=nf90_inq_varid(ncid, 'x_out', xoutid)
256  if (status /= nf90_noerr) call handle_err(status)
257  status=nf90_inq_varid(ncid, 'y_out', youtid)
258  if (status /= nf90_noerr) call handle_err(status)
259  status=nf90_get_var(ncid,xoutid,x_out)
260  if (status /= nf90_noerr) call handle_err(status)
261  status=nf90_get_var(ncid,youtid,y_out)
262  if (status /= nf90_noerr) call handle_err(status)
263  status=nf90_close(ncid)
264  if (status /= nf90_noerr) call handle_err(status)
265 
266 c m_out(:,:)=1
267  do jj= 1,jdm_out
268  do ii= 1,idm_out
269  if (x_out(ii,jj).lt.hspval) then
270  i_out(ii,jj) = int(x_out(ii,jj))
271  x_out(ii,jj) = x_out(ii,jj) - i_out(ii,jj)
272  j_out(ii,jj) = int(y_out(ii,jj))
273  y_out(ii,jj) = y_out(ii,jj) - j_out(ii,jj)
274  else !update output mask???
275  m_out(ii,jj) = 0
276  endif
277  enddo !ii
278  enddo !jj
279 c
280 c --- form the p-grid smoother mask,
281 c --- set to 2 if a land point is needed for interpolation.
282 c --- we are assuming that "2" is never needed outside the
283 c --- target subregion, which will be the case unless the
284 c --- subregion rectangle is poorly chosen.
285 c
286 
287  count = 0
288  do jj= 1,jdm_out
289  do ii= 1,idm_out
290  count = 0
291  if (m_out(ii,jj).eq.1) then
292  j = j_out(ii,jj)
293  i = i_out(ii,jj)
294  jp = min(j+1,jdm)
295  if (i.ne.idm) then
296  ip = i+1
297  else !lperiod
298  ip = 1
299  endif
300  if (m_sm(i,j).eq.0) then
301  count = count+1
302  m_sm(i, j ) = 2
303  endif
304  if (m_sm(i,jp).eq.0) then
305  count = count+1
306  m_sm(i, jp) = 2
307  endif
308  if (m_sm(ip,j).eq.0) then
309  count = count+1
310  m_sm(ip,j) = 2
311  endif
312  if (m_sm(ip,jp).eq.0) then
313  count = count+1
314  m_sm(ip,jp) = 2
315  endif
316  if (count > 0) count1=count1+1
317  endif !sea-point
318  enddo !ii
319  enddo !jj
320 c
321  print *, "Total ", count1, "ocean point hits land"
322 
323  do j= 1,jdm
324  iv_sm(j,1) = idm
325  do i= 1,idm
326  if (m_sm(i,j).eq.2) then
327  iv_sm(j,1) = i
328  exit
329  endif
330  enddo
331  iv_sm(j,2) = 1
332  do i= idm,iv_sm(j,1),-1
333  if (m_sm(i,j).eq.2) then
334  iv_sm(j,2) = i
335  exit
336  endif
337  enddo
338  enddo
339  jf_sm = jdm
340  do j= 1,jdm
341  if (iv_sm(j,1).le.iv_sm(j,2)) then
342  jf_sm = j
343  exit
344  endif
345  enddo
346  jl_sm = 1
347  do j= jdm,jf_sm,-1
348  if (iv_sm(j,1).le.iv_sm(j,2)) then
349  jl_sm = j
350  exit
351  endif
352  enddo
353  if_sm = minval(iv_sm(:,1))
354  il_sm = maxval(iv_sm(:,2))
355 
356  call landfill( a_in,m_sm,idm, jdm,
357  & iv_sm,if_sm,il_sm,jf_sm,jl_sm)
358  call bilinear_p(a_in, idm, jdm,
359  & a_out, idm_out,jdm_out,
360  & m_sm, m_out,i_out,j_out,x_out,y_out)
361 
362 
363 c output the interpolated field
364  status = nf90_open(file_out,nf90_write,ncid2)
365  if (status /= nf90_noerr) call handle_err(status)
366  start_3d(1)=1
367  start_3d(2)=1
368  start_3d(3)=level
369  count_3d(1)=idm_out
370  count_3d(2)=jdm_out
371  count_3d(3)=1
372  status=nf90_inq_varid(ncid2,field,vid)
373  if (status /= nf90_noerr) call handle_err(status)
374  status=nf90_put_var(ncid2, vid,a_out,
375  & start=start_3d,count=count_3d)
376  if (status /= nf90_noerr) call handle_err(status)
377  status=nf90_close(ncid2)
378  if (status /= nf90_noerr) call handle_err(status)
379 
380 
381 c write out a land mask based on the output
382  do j=1,jdm_out
383  do i=1,idm_out
384  if (a_out(i,j) .eq. spval) then
385  m_out(i,j)=0
386  else
387  m_out(i,j)=1
388  endif
389  enddo
390  enddo
391 
392 C Write out HD5 file
393  file_out=trim(dir_out)//'/N'//trim(timetag)//'00_SALINITY_HYCOM_24h.h5'
394 
395  call create_hdf5(file_out, file_id)
396 
397  rank = 2
398  dims_sal(1) = 360*4
399  dims_sal(2) = 180*4
400 
401  call create_hdf5_real(file_id, 'salinity', rank, dims_sal, dset_id)
402  call set_space_hdf5(dset_id, rank, dims_sal, filespace, dataspace)
403  pt => a_out(1,1)
404  call write_hdf5_real(dset_id, filespace, dataspace, offset, dims_sal, pt)
405  call close_hdf5_ds(dset_id, dataspace, filespace)
406  call close_hdf5_df(file_id)
407 
408 
409 C Write binary file
410  file_out=trim(dir_out)//'/N'//trim(timetag)//'00_SALINITY_HYCOM_24h.dat'
411  open(unit=8,file=file_out, form='UNFORMATTED')
412  write(8) a_out
413  close(8)
414 
415  deallocate(iv_sm, m_sm, m_osm, m_out, a_in, a_out)
416  deallocate(i_out, j_out, x_out, y_out)
417 
418  write(*,*) ' '
419  write(*,*) 'NETCDF=N'//trim(timetag)//'00_SALINITY_HYCOM_24h.nc'
420  write(*,*) 'HDF5 =N'//trim(timetag)//'00_SALINITY_HYCOM_24h.h5'
421  write(*,*) 'BINARY=N'//trim(timetag)//'00_SALINITY_HYCOM_24h.dat'
422 
423  end program interp_hycom
424 
425  subroutine bilinear_p(a_in, idm_in, jdm_in,
426  & a_out,idm_out,jdm_out,
427  & m_in, m_out,i_out,j_out,x_out,y_out)
428  implicit none
429 c
430  integer idm_in, jdm_in,
431  & idm_out,jdm_out
432  integer m_out(idm_out,jdm_out),
433  & m_in(idm_in,jdm_in),
434  & i_out(idm_out,jdm_out),
435  & j_out(idm_out,jdm_out)
436  real a_in( idm_in, jdm_in ),
437  & a_out(idm_out,jdm_out),
438  & x_out(idm_out,jdm_out),
439  & y_out(idm_out,jdm_out)
440 c
441 c --- interpolate from a_in to a_out.
442 c
443  integer i,ii,ip,j,jj,jp
444  integer ki,kj,iii,jjj
445  real sa, ss
446  real sx,sy
447  integer count
448  real s(-1:1, -1:1)
449  data s / 1.0, 2.0, 1.0,
450  & 2.0, 4.0, 2.0,
451  & 1.0, 2.0, 1.0 /
452  real, parameter :: spval=2.0**100 ! spval
453 
454 c
455  count = 0
456  do jj= 1,jdm_out
457  do ii= 1,idm_out
458  if (m_out(ii,jj).eq.1) then
459  count = count+1
460  sx = x_out(ii,jj)
461  sy = y_out(ii,jj)
462  i = i_out(ii,jj)
463  if (i.ne.idm_in) then
464  ip = i+1
465  else
466  ip = 1
467  endif
468  j = j_out(ii,jj)
469  jp = j+1
470 c
471 c doing a 3x3 interpolation for jj<552
472  if (jj<552) then
473  sa = 0.0
474  ss = 0.0
475  do kj= -1,1
476  jjj = j+kj
477  if (jjj==0) jjj=1
478  if (jjj>jdm_in) jjj=jdm_in
479  do ki= -1,1
480  iii = i+ki
481  if (iii==0) iii=idm_in
482  if (iii>idm_in) iii=1
483  if (m_in(iii,jjj).eq.1) then
484  sa = sa + s(ki,kj)*a_in(iii,jjj)
485  ss = ss + s(ki,kj)
486  endif
487  enddo
488  enddo
489  if (ss.ne.0.0) then
490 c
491 c at least one ocean point within stencil.
492 c
493  a_out( ii,jj) = sa/ss
494  else
495  a_out(ii,jj) = spval
496  endif
497  else
498  a_out(ii,jj) = (1.0-sx)*(1.0-sy)*a_in(i, j ) +
499  & (1.0-sx)* sy *a_in(i, jp) +
500  & sx *(1.0-sy)*a_in(ip,j ) +
501  & sx * sy *a_in(ip,jp)
502 * if (a_out(ii,jj).lt.0.0) then
503 * write(6,'(a,6i5,2f7.3,5f9.2)')
504 * & 'ii,jj,i,ip,j,jp,sx,sy,a_out,a_in',
505 * & ii,jj,i,ip,j,jp,
506 * & sx,sy,a_out(ii,jj),
507 * & a_in(i, j ),
508 * & a_in(i, jp),
509 * & a_in(ip,j ),
510 * & a_in(ip,jp)
511 * endif
512  endif
513  endif
514  enddo
515  enddo
516  print *, "Output ocean point:", count
517  return
518  end
519 
520  subroutine landfill(a,mask,m,n, iv,if,il,jf,jl)
521  implicit none
522 c
523  integer m,n,mask(m,n), iv(n,2),if,il,jf,jl
524  real a(m,n)
525 c
526 c --- extrapolate a 1-grid cell into the land mask,
527 c --- mask == 0 for land.
528 c --- mask == 1 for ocean.
529 c --- mask == 2 for land to be extrapolated to ocean.
530 c
531  integer, allocatable :: mm(:,:,:)
532 c
533  integer i,ii,ip0,ip1,ipass,j,jj,ki,kj,nleft,nup
534  real sa,ss
535 c
536  logical lfirst
537 c real s(-1:1,-1:1)
538  real s(-2:2,-2:2)
539  save lfirst,s
540 c
541  data lfirst / .true. /
542 c data s / 1.0, 2.0, 1.0,
543 c & 2.0, 4.0, 2.0,
544 c & 1.0, 2.0, 1.0 /
545  data s /1.0, 1.5, 2.0, 1.5, 1.0,
546  & 1.5, 2.0, 3.0, 2.0, 1.5,
547  & 2.0, 3.0, 4.0, 3.0, 2.0,
548  & 1.5, 2.0, 3.0, 2.0, 1.5,
549  & 1.0, 1.5, 2.0, 1.5, 1.0 /
550 c
551 c adding a halo to mm simplifies ocean selection logic.
552 c
553  allocate( mm(0:m+1,0:n+1,0:1) )
554 c
555  mm( : , : ,0) = 0
556  mm(1:m,1:n,0) = mask
557  mm( : , : ,1) = mm(:,:,0)
558 c
559 c --- repeated passes of 9-point "smoother" to
560 c --- convert all mask==2 points to mask==1.
561 c --- double-buffering mm allows in-place use of a.
562 c
563  if (lfirst) then
564  write(6,'(/a,6i5/)')
565  & 'landfill - m,n,if,il,jf,jl =',m,n,if,il,jf,jl
566  endif
567  do ipass= 1,n+m
568  ip0 = mod(ipass+1,2)
569  ip1 = mod(ipass, 2)
570  nup = 0
571  nleft = 0
572  do j= jf,jl
573  do i= iv(j,1),iv(j,2)
574  if (mm(i,j,ip0).eq.2) then
575  sa = 0.0
576  ss = 0.0
577  do kj= -2,2
578  jj = j+kj
579  do ki= -2,2
580  ii = i+ki
581  if ((ii>0) .and. (ii.le.m) .and. (jj>0)
582  & .and. (jj.le.n)) then
583  if (mm(ii,jj,ip0).eq.1) then
584  sa = sa + s(ki,kj)*a(ii,jj)
585  ss = ss + s(ki,kj)
586  endif
587  endif
588  enddo
589  enddo
590  if (ss.ne.0.0) then
591 c
592 c at least one ocean point within stencil.
593 c
594  a( i,j) = sa/ss
595  mm(i,j,ip1) = 1
596  nup = nup + 1
597 * if (mask(i,j).eq.1) then
598 * write(6,*) 'error - i,j,ip0,ip1,mask,mm = ',
599 * & i,j,ip0,ip1,mask(i,j),mm(i,j,ip0)
600 * stop
601 * endif
602 * if (mod(nup,1000).eq.1) then
603 * write(6,'(a,2i5,f5.1,f10.3)')
604 * & ' i,j,ss,a = ',i,j,ss,a(i,j)
605 * endif
606  else
607  nleft = nleft + 1
608  endif
609  endif
610  enddo
611  enddo
612  if (lfirst) then
613  write(6,'(a,i4,a,i6,a,i6,a)')
614  & 'landfill: pass',ipass,
615  & ' filled in',nup,
616  & ' points, with',nleft,' still to fill'
617  call flush(6)
618  endif
619  if (nup.eq.0) then
620  exit
621  endif
622  mm(if:il,jf:jl,ip0) = mm(if:il,jf:jl,ip1)
623  enddo ! ipass=1,...
624  if (lfirst) then
625  write(6,*)
626  lfirst = .false.
627  endif
628  if (nleft.ne.0) then
629  write(6,'(/a,i6,a/a/)')
630  & 'error in landfill - ',
631  & nleft,' "mask==2" values are not fillable',
632  & 'probably a mismatch between input and output land masks'
633  call flush(6)
634 c do j= jf,jl
635 c do i= iv(j,1),iv(j,2)
636 c if (mm(i,j,ip1).eq.2) then
637 c write(6,'(a,2i5)') 'mask==2 at (input) i,j = ',i,j
638 c endif
639 c enddo
640 c enddo
641 c write(6,*)
642 c call flush(6)
643 c stop
644  endif
645 c
646  mask=mm(1:m,1:n,ip1)
647  deallocate( mm )
648 c
649  return
650  end subroutine landfill
651 
652  subroutine psmooth(a,mask,amn,amx,m,n)
653  implicit none
654 c
655  integer m,n,mask(m,n)
656  real a(m,n),amn(m,n),amx(m,n)
657 c
658 c --- smooth under mask and within amn,amx range.
659 c
660  integer, allocatable :: mm(:,:)
661  real, allocatable :: aa(:,:)
662 c
663  integer i,ii,j,jj,ki,kj
664  real rss,sa
665 c
666  real s(-1:1,-1:1)
667  save s
668  data s / 1.0, 2.0, 1.0,
669  & 2.0, 4.0, 2.0,
670  & 1.0, 2.0, 1.0 /
671 c
672  rss = 1.0/sum(s(:,:))
673 c
674 c local copy of a.
675 c
676  allocate( aa(m,n) )
677  aa = a
678 c
679 c adding a halo to mm simplifies ocean selection logic.
680 c
681  allocate( mm(0:m+1,0:n+1) )
682  mm( 0, : ) = 0
683  mm(m+1, : ) = 0
684  mm( : , 0) = 0
685  mm( : ,n+1) = 0
686  mm(1:m,1:n) = mask
687 c
688  do j= 1,n
689  do i= 1,m
690  if (mm(i,j).eq.1) then
691  sa = 0.0
692  do kj= -1,1
693  jj = j+kj
694  do ki= -1,1
695  ii = i+ki
696  if (mm(ii,jj).eq.1) then
697  sa = sa + s(ki,kj)*aa(ii,jj) ! must use local copy of a
698  else
699  sa = sa + s(ki,kj)*aa(i ,j )
700  endif
701  enddo
702  enddo
703  a(i,j) = max( amn(i,j),
704  & min( amx(i,j), sa*rss ) )
705  endif
706  enddo
707  enddo
708 c
709  deallocate( aa, mm )
710 c
711  return
712  end subroutine psmooth
713 
714 c
715  subroutine handle_err(status)
716  integer :: status
717  print *, "NetCDF return error ", status
718 ! stop "Stopped"
719  call exit(1)
720  end subroutine handle_err
721 
722 
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 the first two or last two tables respectively will be used for the interpolation If there is only one LUT in the set of it will be treated as a constant LUT The manner in which Earth View data is checked for saturation was changed Previously the raw Earth View DNs and Space View DNs were checked against the lookup table values contained in the table dn_sat The change made is to check the raw Earth and Space View DNs to be sure they are less than the maximum saturation value and to check the Space View subtracted Earth View dns against a set of values contained in the new lookup table dn_sat_ev The metadata configuration and ASSOCIATEDINSTRUMENTSHORTNAME from the MOD02HKM product The same metatdata with extensions and were removed from the MOD021KM and MOD02OBC products ASSOCIATEDSENSORSHORTNAME was set to MODIS in all products These changes are reflected in new File Specification which users may consult for exact the pow functions were eliminated in Emissive_Cal and Emissive bands replaced by more efficient code Other calculations throughout the code were also made more efficient Aside from a few round off there was no difference to the product The CPU time decreased by about for a day case and for a night case A minor bug in calculating the uncertainty index for emissive bands was corrected The frame index(0-based) was previously being used the frame number(1-based) should have been used. There were only a few minor changes to the uncertainty index(maximum of 1 digit). 3. Some inefficient arrays(Sigma_RVS_norm_sq) were eliminated and some code lines in Preprocess_L1A_Data were moved into Process_OBCEng_Emiss. There were no changes to the product. Required RAM was reduced by 20 MB. Now
#define max(A, B)
Definition: main_biosmap.c:61
program interp_hycom
Definition: interp_hycom.f:1
subroutine handle_err(status)
Definition: interp_hycom.f:716
string & trim(string &s, const string &delimiters)
Definition: EnvsatUtil.cpp:29
#define min(A, B)
Definition: main_biosmap.c:62
subroutine landfill(a, mask, m, n, iv, if, il, jf, jl)
Definition: interp_hycom.f:521
subroutine bilinear_p(a_in, idm_in, jdm_in, a_out, idm_out, jdm_out, m_in, m_out, i_out, j_out, x_out, y_out)
Definition: interp_hycom.f:428
subroutine psmooth(a, mask, amn, amx, m, n)
Definition: interp_hycom.f:653