21 character(len=50),
allocatable,
public :: field_names(:)
22 character(len=75),
public :: source
24 integer,
public :: i_src
25 integer,
public :: j_src
26 integer,
public :: num_records
27 integer,
public :: num_time_recs
28 integer,
public :: num_fields
31 integer,
allocatable,
public :: day_of_rec(:)
34 type(esmf_grid
),
public :: grid_src
58 character(len=*),
intent(in) :: input_file
60 integer,
intent(in) :: localpet, npets
62 character(len=50) :: vname
64 integer :: dimid, dims(1), ncid, status
65 integer :: count, num_vars, n
66 integer :: rc, varid, i, j, i_src_corner
67 integer(esmf_kind_i4),
pointer :: mask_ptr(:,:)
68 integer :: clb(2), cub(2)
69 integer :: clb_corner(2), cub_corner(2)
71 real(esmf_kind_r4),
allocatable :: mask_global(:,:)
72 real(esmf_kind_r8),
allocatable :: lat_global(:)
73 real(esmf_kind_r8),
allocatable :: lon_global(:)
74 real(esmf_kind_r8),
allocatable :: lat_corner_global(:)
75 real(esmf_kind_r8),
allocatable :: lon_corner_global(:)
76 real(esmf_kind_r4),
pointer :: mask_field_ptr(:,:)
77 real(esmf_kind_r8),
pointer :: lat_ptr(:,:)
78 real(esmf_kind_r8),
pointer :: lon_ptr(:,:)
79 real(esmf_kind_r8),
pointer :: lat_corner_ptr(:,:)
80 real(esmf_kind_r8),
pointer :: lon_corner_ptr(:,:)
83 type(esmf_field
) :: mask_field
84 type(esmf_polekind_flag
) :: polekindflag(2)
86 print*,
"- OPEN FILE: ", trim(input_file)
87 status = nf90_open(input_file, nf90_nowrite, ncid)
88 call
netcdf_err(status,
"OPENING INPUT SOURCE FILE")
90 status = nf90_get_att(ncid, nf90_global,
'source', source)
91 if (status /= nf90_noerr) source=
"unknown"
93 if(localpet == 0) print*,
'- SOURCE OF DATA IS: ', trim(source)
95 status = nf90_inq_dimid(ncid,
'idim', dimid)
96 call
netcdf_err(status,
"READING I DIMENSION ID OF SOURCE DATA")
97 status = nf90_inquire_dimension(ncid, dimid, len=i_src)
98 call
netcdf_err(status,
"READING I DIMENSION OF SOURCE DATA")
100 status = nf90_inq_dimid(ncid,
'jdim', dimid)
101 call
netcdf_err(status,
"READING J DIMENSION ID OF SOURCE DATA")
102 status = nf90_inquire_dimension(ncid, dimid, len=j_src)
103 call
netcdf_err(status,
"READING J DIMENSION OF SOURCE DATA")
105 if(localpet == 0) print*,
'- I/J DIMENSION OF SOURCE DATA: ', i_src, j_src
107 allocate(lat_global(j_src))
108 status = nf90_inq_varid(ncid,
'lat', varid)
109 call
netcdf_err(status,
"READING SOURCE DATA LATITUDE ID")
110 status = nf90_get_var(ncid, varid, lat_global)
111 call
netcdf_err(status,
"READING SOURCE DATA LATITUDES")
113 allocate(lon_global(i_src))
114 status = nf90_inq_varid(ncid,
'lon', varid)
115 call
netcdf_err(status,
"READING SOURCE DATA LONGITUDE ID")
116 status = nf90_get_var(ncid, varid, lon_global)
117 call
netcdf_err(status,
"READING SOURCE DATA LONGITUDE")
119 allocate(lat_corner_global(j_src+1))
120 status = nf90_inq_varid(ncid,
'lat_corner', varid)
121 call
netcdf_err(status,
"READING SOURCE DATA CORNER LATITUDE ID")
122 status = nf90_get_var(ncid, varid, lat_corner_global)
123 call
netcdf_err(status,
"READING SOURCE DATA CORNER LATITUDE")
130 status = nf90_inq_varid(ncid,
'lon_corner', varid)
131 call
netcdf_err(status,
"READING SOURCE DATA CORNER LONGITUDE ID")
132 status = nf90_inquire_variable(ncid, varid, dimids=dims)
133 call
netcdf_err(status,
"READING SOURCE DATA CORNER LONGITUDE DIMIDS")
134 status = nf90_inquire_dimension(ncid, dims(1), len=i_src_corner)
135 call
netcdf_err(status,
"READING SOURCE DATA CORNER LONGITUDE DIMS")
136 allocate(lon_corner_global(i_src_corner))
137 status = nf90_get_var(ncid, varid, lon_corner_global)
138 call
netcdf_err(status,
"READING SOURCE DATA CORNER LONGITUDE")
140 status = nf90_inq_dimid(ncid,
'time', dimid)
141 call
netcdf_err(status,
"READING SOURCE DATA TIME ID")
142 status = nf90_inquire_dimension(ncid, dimid, len=num_time_recs)
143 call
netcdf_err(status,
"READING SOURCE DATA NUM TIME PERIODS")
145 allocate(day_of_rec(num_time_recs))
146 status = nf90_inq_varid(ncid,
'time', varid)
147 call
netcdf_err(status,
"READING SOURCE DATA TIME VARID")
148 status = nf90_get_var(ncid, varid, day_of_rec)
149 call
netcdf_err(status,
"READING SOURCE DATA DAY OF RECORD")
151 print*,
'- SOURCE DATA DAYS OF RECORD: ', day_of_rec
153 status = nf90_inquire(ncid, nvariables=num_vars)
154 call
netcdf_err(status,
"READING NUMBER OF VARIABLES SOURCE DATA")
162 num_fields = num_vars - 5
163 num_records = num_vars * num_time_recs
165 allocate(field_names(num_fields))
170 status = nf90_inquire_variable(ncid, n, name=vname)
171 call
netcdf_err(status,
"READING FIELD NAMES")
173 if (trim(vname) ==
'time') cycle
174 if (trim(vname) ==
'lon') cycle
175 if (trim(vname) ==
'lon_corner') cycle
176 if (trim(vname) ==
'lat') cycle
177 if (trim(vname) ==
'lat_corner') cycle
180 field_names(count) = vname
184 if(localpet==0) print*,
'- FIELDS TO BE PROCESSED: ', field_names
186 if (localpet == 0)
then
187 allocate(mask_global(i_src,j_src))
188 status = nf90_inq_varid(ncid, field_names(1), varid)
190 status = nf90_get_var(ncid, varid, mask_global)
193 allocate(mask_global(0,0))
196 status = nf90_close(ncid)
209 lon_extent = mod((lon_global(i_src)-lon_global(1))-1+3600,360.)+1.0
211 if (lon_extent < 350.0)
then
213 print*,
"- CALL GridCreateNoPeriDim FOR SOURCE GRID"
214 grid_src = esmf_gridcreatenoperidim(minindex=(/1,1/), &
215 maxindex=(/i_src,j_src/), &
216 coordsys=esmf_coordsys_sph_deg, &
217 regdecomp=(/1,npets/), &
218 name=
"source_grid", &
219 indexflag=esmf_index_global, rc=rc)
220 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
225 polekindflag(1:2) = esmf_polekind_monopole
227 print*,
"- CALL GridCreate1PeriDim FOR SOURCE GRID"
228 grid_src = esmf_gridcreate1peridim(minindex=(/1,1/), &
229 maxindex=(/i_src,j_src/), &
230 polekindflag=polekindflag, &
233 coordsys=esmf_coordsys_sph_deg, &
234 regdecomp=(/1,npets/), &
235 name=
"source_grid", &
236 indexflag=esmf_index_global, rc=rc)
237 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
242 print*,
"- CALL FieldCreate FOR SOURCE GRID LANDMASK."
243 mask_field = esmf_fieldcreate(grid_src, &
244 typekind=esmf_typekind_r4, &
245 staggerloc=esmf_staggerloc_center, &
246 name=
"source grid land mask", &
248 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
251 print*,
"- CALL FieldScatter FOR SOURCE GRID MASK."
252 call esmf_fieldscatter(mask_field, mask_global, rootpet=0, rc=rc)
253 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
256 print*,
"- CALL GridAddItem FOR SOURCE GRID MASK."
257 call esmf_gridadditem(grid_src, &
258 itemflag=esmf_griditem_mask, &
259 staggerloc=esmf_staggerloc_center, &
261 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
264 print*,
"- CALL GridGetItem FOR SOURCE GRID MASK."
266 call esmf_gridgetitem(grid_src, &
267 itemflag=esmf_griditem_mask, &
268 farrayptr=mask_ptr, &
272 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
275 print*,
"- CALL FieldGet FOR SOURCE GRID LANDMASK."
276 nullify(mask_field_ptr)
277 call esmf_fieldget(mask_field, &
278 farrayptr=mask_field_ptr, &
280 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
283 do j = clb(2), cub(2)
284 do i = clb(1), cub(1)
285 if (mask_field_ptr(i,j) < -1.0)
then
293 deallocate(mask_global)
295 print*,
"- CALL FieldDestroy FOR SOURCE GRID LAND MASK."
296 call esmf_fielddestroy(mask_field,rc=rc)
297 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
302 print*,
"- CALL GridAddCoord FOR SOURCE GRID CENTER LOCATION."
303 call esmf_gridaddcoord(grid_src, &
304 staggerloc=esmf_staggerloc_center, rc=rc)
305 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
308 print*,
"- CALL GridGetCoord FOR SOURCE GRID CENTER LONGITUDE."
310 call esmf_gridgetcoord(grid_src, &
311 staggerloc=esmf_staggerloc_center, &
313 farrayptr=lon_ptr, rc=rc)
314 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
317 print*,
"- CALL GridGetCoord FOR SOURCE GRID CENTER LATITUDE."
319 call esmf_gridgetcoord(grid_src, &
320 staggerloc=esmf_staggerloc_center, &
322 farrayptr=lat_ptr, rc=rc)
323 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
326 do j = clb(2), cub(2)
327 lat_ptr(:,j) = lat_global(j)
330 do i = clb(1), cub(1)
331 lon_ptr(i,:) = lon_global(i)
334 print*,
"- CALL GridAddCoord FOR SOURCE GRID CORNER LOCATION."
335 call esmf_gridaddcoord(grid_src, &
336 staggerloc=esmf_staggerloc_corner, &
338 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
341 print*,
"- CALL GridGetCoord FOR SOURCE GRID CORNER LONGITUDE."
342 nullify(lon_corner_ptr)
343 call esmf_gridgetcoord(grid_src, &
344 staggerloc=esmf_staggerloc_corner, &
346 farrayptr=lon_corner_ptr, &
348 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
351 print*,
"- CALL GridGetCoord FOR SOURCE GRID CORNER LATITUDE."
352 nullify(lat_corner_ptr)
353 call esmf_gridgetcoord(grid_src, &
354 staggerloc=esmf_staggerloc_corner, &
356 computationallbound=clb_corner, &
357 computationalubound=cub_corner, &
358 farrayptr=lat_corner_ptr, &
360 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
363 do j = clb_corner(2), cub_corner(2)
364 lat_corner_ptr(:,j) = lat_corner_global(j)
367 do i = clb_corner(1), cub_corner(1)
368 lon_corner_ptr(i,:) = lon_corner_global(i)
371 deallocate(lat_global)
372 deallocate(lon_global)
373 deallocate(lat_corner_global)
374 deallocate(lon_corner_global)
387 print*,
"- CALL GridDestroy FOR SOURCE GRID."
388 call esmf_griddestroy(grid_src,rc=rc)
390 deallocate(day_of_rec)
391 deallocate(field_names)
subroutine, public source_grid_cleanup
Free up memory associated with this module.
subroutine, public define_source_grid(localpet, npets, input_file)
Defines esmf grid object for source grid.
subroutine, public netcdf_err(err, string)
Handle netCDF error codes.
subroutine, public error_handler(string, rc)
Handle errors.
Read grid specs, date information and land/sea mask for the source data that will be interpolated to ...