58 character(len=*),
intent(in) :: input_file
60 integer,
intent(in) :: localpet, npets
62 character(len=50) :: field_names_save(100)
63 character(len=50) :: vname
65 integer :: dimid, dims(1), ncid, status
66 integer :: count, num_vars, n
67 integer :: rc, varid, i, j, i_src_corner
68 integer(esmf_kind_i4),
pointer :: mask_ptr(:,:)
69 integer :: clb(2), cub(2)
70 integer :: clb_corner(2), cub_corner(2)
72 real(esmf_kind_r4),
allocatable :: mask_global(:,:)
73 real(esmf_kind_r8),
allocatable :: lat_global(:)
74 real(esmf_kind_r8),
allocatable :: lon_global(:)
75 real(esmf_kind_r8),
allocatable :: lat_corner_global(:)
76 real(esmf_kind_r8),
allocatable :: lon_corner_global(:)
77 real(esmf_kind_r4),
pointer :: mask_field_ptr(:,:)
78 real(esmf_kind_r8),
pointer :: lat_ptr(:,:)
79 real(esmf_kind_r8),
pointer :: lon_ptr(:,:)
80 real(esmf_kind_r8),
pointer :: lat_corner_ptr(:,:)
81 real(esmf_kind_r8),
pointer :: lon_corner_ptr(:,:)
83 real(esmf_kind_r4) :: missing
85 type(esmf_field) :: mask_field
86 type(esmf_polekind_flag) :: polekindflag(2)
88 print*,
"- OPEN FILE: ", trim(input_file)
89 status = nf90_open(input_file, nf90_nowrite, ncid)
90 call netcdf_err(status,
"OPENING INPUT SOURCE FILE")
92 status = nf90_get_att(ncid, nf90_global,
'source',
source)
93 if (status /= nf90_noerr)
source=
"unknown" 95 if(localpet == 0) print*,
'- SOURCE OF DATA IS: ', trim(
source)
97 status = nf90_inq_dimid(ncid,
'idim', dimid)
98 call netcdf_err(status,
"READING I DIMENSION ID OF SOURCE DATA")
99 status = nf90_inquire_dimension(ncid, dimid, len=
i_src)
100 call netcdf_err(status,
"READING I DIMENSION OF SOURCE DATA")
102 status = nf90_inq_dimid(ncid,
'jdim', dimid)
103 call netcdf_err(status,
"READING J DIMENSION ID OF SOURCE DATA")
104 status = nf90_inquire_dimension(ncid, dimid, len=
j_src)
105 call netcdf_err(status,
"READING J DIMENSION OF SOURCE DATA")
107 if(localpet == 0) print*,
'- I/J DIMENSION OF SOURCE DATA: ',
i_src,
j_src 109 allocate(lat_global(
j_src))
110 status = nf90_inq_varid(ncid,
'lat', varid)
111 call netcdf_err(status,
"READING SOURCE DATA LATITUDE ID")
112 status = nf90_get_var(ncid, varid, lat_global)
113 call netcdf_err(status,
"READING SOURCE DATA LATITUDES")
115 allocate(lon_global(
i_src))
116 status = nf90_inq_varid(ncid,
'lon', varid)
117 call netcdf_err(status,
"READING SOURCE DATA LONGITUDE ID")
118 status = nf90_get_var(ncid, varid, lon_global)
119 call netcdf_err(status,
"READING SOURCE DATA LONGITUDE")
121 allocate(lat_corner_global(
j_src+1))
122 status = nf90_inq_varid(ncid,
'lat_corner', varid)
123 call netcdf_err(status,
"READING SOURCE DATA CORNER LATITUDE ID")
124 status = nf90_get_var(ncid, varid, lat_corner_global)
125 call netcdf_err(status,
"READING SOURCE DATA CORNER LATITUDE")
132 status = nf90_inq_varid(ncid,
'lon_corner', varid)
133 call netcdf_err(status,
"READING SOURCE DATA CORNER LONGITUDE ID")
134 status = nf90_inquire_variable(ncid, varid, dimids=dims)
135 call netcdf_err(status,
"READING SOURCE DATA CORNER LONGITUDE DIMIDS")
136 status = nf90_inquire_dimension(ncid, dims(1), len=i_src_corner)
137 call netcdf_err(status,
"READING SOURCE DATA CORNER LONGITUDE DIMS")
138 allocate(lon_corner_global(i_src_corner))
139 status = nf90_get_var(ncid, varid, lon_corner_global)
140 call netcdf_err(status,
"READING SOURCE DATA CORNER LONGITUDE")
142 status = nf90_inq_dimid(ncid,
'time', dimid)
143 call netcdf_err(status,
"READING SOURCE DATA TIME ID")
144 status = nf90_inquire_dimension(ncid, dimid, len=
num_time_recs)
145 call netcdf_err(status,
"READING SOURCE DATA NUM TIME PERIODS")
148 status = nf90_inq_varid(ncid,
'time', varid)
149 call netcdf_err(status,
"READING SOURCE DATA TIME VARID")
150 status = nf90_get_var(ncid, varid,
day_of_rec)
151 call netcdf_err(status,
"READING SOURCE DATA DAY OF RECORD")
153 print*,
'- SOURCE DATA DAYS OF RECORD: ',
day_of_rec 155 status = nf90_inquire(ncid, nvariables=num_vars)
156 call netcdf_err(status,
"READING NUMBER OF VARIABLES SOURCE DATA")
171 status = nf90_inquire_variable(ncid, n, name=vname)
172 call netcdf_err(status,
"READING FIELD NAMES")
174 if (trim(vname) ==
'time') cycle
175 if (trim(vname) ==
'lon') cycle
176 if (trim(vname) ==
'lon_corner') cycle
177 if (trim(vname) ==
'lat') cycle
178 if (trim(vname) ==
'lat_corner') cycle
179 if (trim(vname) ==
'clay_lev1') cycle
180 if (trim(vname) ==
'clay_top') cycle
181 if (trim(vname) ==
'sand_lev1') cycle
182 if (trim(vname) ==
'sand_top') cycle
185 field_names_save(count) = vname
196 if(localpet==0) print*,
'- FIELDS TO BE PROCESSED: ',
field_names 198 if (localpet == 0)
then 200 status = nf90_inq_varid(ncid,
field_names(1), varid)
202 status = nf90_get_var(ncid, varid, mask_global)
205 allocate(mask_global(0,0))
213 status = nf90_inq_varid(ncid,
field_names(1), varid)
215 status=nf90_get_att(ncid, varid,
'missing_value', missing)
216 call netcdf_err(status,
"READING MISSING VALUE")
218 status = nf90_close(ncid)
231 lon_extent = mod((lon_global(
i_src)-lon_global(1))-1+3600,360.)+1.0
233 if (lon_extent < 350.0)
then 235 print*,
"- CALL GridCreateNoPeriDim FOR SOURCE GRID" 236 grid_src = esmf_gridcreatenoperidim(minindex=(/1,1/), &
238 coordsys=esmf_coordsys_sph_deg, &
239 regdecomp=(/1,npets/), &
240 name=
"source_grid", &
241 indexflag=esmf_index_global, rc=rc)
242 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
247 polekindflag(1:2) = esmf_polekind_monopole
249 print*,
"- CALL GridCreate1PeriDim FOR SOURCE GRID" 250 grid_src = esmf_gridcreate1peridim(minindex=(/1,1/), &
252 polekindflag=polekindflag, &
255 coordsys=esmf_coordsys_sph_deg, &
256 regdecomp=(/1,npets/), &
257 name=
"source_grid", &
258 indexflag=esmf_index_global, rc=rc)
259 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
264 print*,
"- CALL FieldCreate FOR SOURCE GRID LANDMASK." 265 mask_field = esmf_fieldcreate(
grid_src, &
266 typekind=esmf_typekind_r4, &
267 staggerloc=esmf_staggerloc_center, &
268 name=
"source grid land mask", &
270 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
273 print*,
"- CALL FieldScatter FOR SOURCE GRID MASK." 274 call esmf_fieldscatter(mask_field, mask_global, rootpet=0, rc=rc)
275 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
278 print*,
"- CALL GridAddItem FOR SOURCE GRID MASK." 280 itemflag=esmf_griditem_mask, &
281 staggerloc=esmf_staggerloc_center, &
283 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
286 print*,
"- CALL GridGetItem FOR SOURCE GRID MASK." 289 itemflag=esmf_griditem_mask, &
290 farrayptr=mask_ptr, &
294 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
297 print*,
"- CALL FieldGet FOR SOURCE GRID LANDMASK." 298 nullify(mask_field_ptr)
299 call esmf_fieldget(mask_field, &
300 farrayptr=mask_field_ptr, &
302 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
305 do j = clb(2), cub(2)
306 do i = clb(1), cub(1)
307 if ( abs(mask_field_ptr(i,j)-missing) < 0.001)
then 315 deallocate(mask_global)
317 print*,
"- CALL FieldDestroy FOR SOURCE GRID LAND MASK." 318 call esmf_fielddestroy(mask_field,rc=rc)
319 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
324 print*,
"- CALL GridAddCoord FOR SOURCE GRID CENTER LOCATION." 326 staggerloc=esmf_staggerloc_center, rc=rc)
327 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
330 print*,
"- CALL GridGetCoord FOR SOURCE GRID CENTER LONGITUDE." 333 staggerloc=esmf_staggerloc_center, &
335 farrayptr=lon_ptr, rc=rc)
336 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
339 print*,
"- CALL GridGetCoord FOR SOURCE GRID CENTER LATITUDE." 342 staggerloc=esmf_staggerloc_center, &
344 farrayptr=lat_ptr, rc=rc)
345 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
348 do j = clb(2), cub(2)
349 lat_ptr(:,j) = lat_global(j)
352 do i = clb(1), cub(1)
353 lon_ptr(i,:) = lon_global(i)
356 print*,
"- CALL GridAddCoord FOR SOURCE GRID CORNER LOCATION." 358 staggerloc=esmf_staggerloc_corner, &
360 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
363 print*,
"- CALL GridGetCoord FOR SOURCE GRID CORNER LONGITUDE." 364 nullify(lon_corner_ptr)
366 staggerloc=esmf_staggerloc_corner, &
368 farrayptr=lon_corner_ptr, &
370 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
373 print*,
"- CALL GridGetCoord FOR SOURCE GRID CORNER LATITUDE." 374 nullify(lat_corner_ptr)
376 staggerloc=esmf_staggerloc_corner, &
378 computationallbound=clb_corner, &
379 computationalubound=cub_corner, &
380 farrayptr=lat_corner_ptr, &
382 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
385 do j = clb_corner(2), cub_corner(2)
386 lat_corner_ptr(:,j) = lat_corner_global(j)
389 do i = clb_corner(1), cub_corner(1)
390 lon_corner_ptr(i,:) = lon_corner_global(i)
393 deallocate(lat_global)
394 deallocate(lon_global)
395 deallocate(lat_corner_global)
396 deallocate(lon_corner_global)
409 print*,
"- CALL GridDestroy FOR SOURCE GRID." 410 call esmf_griddestroy(
grid_src,rc=rc)
integer, public j_src
j dimension of the source grid.
character(len=50), dimension(:), allocatable, public field_names
Names of fields to be processed.
Read grid specs, date information and land/sea mask for the source data that will be interpolated to ...
integer, public num_time_recs
Number of time records.
subroutine, public define_source_grid(localpet, npets, input_file)
Defines esmf grid object for source grid.
subroutine, public source_grid_cleanup
Free up memory associated with this module.
integer, dimension(:), allocatable, public day_of_rec
Day of each time record with respect to Jan 1.
integer, public i_src
i dimension of the source grid.
integer, public num_records
Number of fields times time records.
subroutine, public error_handler(string, rc)
Handle errors.
integer, public num_fields
Number of fields in the file.
subroutine, public netcdf_err(err, string)
Handle netCDF error codes.
type(esmf_grid), public grid_src
ESMF grid object for the source grid.
character(len=75), public source
Original source of the data.