21 type(esmf_field
),
public :: alvsf_target_grid
22 type(esmf_field),
public :: alvwf_target_grid
23 type(esmf_field),
public :: alnsf_target_grid
24 type(esmf_field),
public :: alnwf_target_grid
25 type(esmf_field),
public :: facsf_target_grid
26 type(esmf_field),
public :: facwf_target_grid
27 type(esmf_field),
public :: max_veg_greenness_target_grid
28 type(esmf_field),
public :: min_veg_greenness_target_grid
29 type(esmf_field),
public :: mxsno_albedo_target_grid
30 type(esmf_field),
public :: slope_type_target_grid
31 type(esmf_field),
public :: soil_type_target_grid
32 type(esmf_field),
public :: substrate_temp_target_grid
33 type(esmf_field),
public :: veg_greenness_target_grid
34 type(esmf_field),
public :: veg_type_target_grid
49 num_tiles_target_grid, &
54 integer,
intent(in) :: localpet
56 integer :: error, tile, i, j
58 real(esmf_kind_r8),
allocatable :: data_one_tile(:,:)
59 real(esmf_kind_r8),
allocatable :: max_data_one_tile(:,:)
60 real(esmf_kind_r8),
allocatable :: min_data_one_tile(:,:)
63 allocate(data_one_tile(i_target,j_target))
65 allocate(data_one_tile(0,0))
72 print*,
"- CALL FieldCreate FOR TARGET GRID SLOPE TYPE."
73 slope_type_target_grid = esmf_fieldcreate(target_grid, &
74 typekind=esmf_typekind_r8, &
75 staggerloc=esmf_staggerloc_center, rc=error)
76 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
79 do tile = 1, num_tiles_target_grid
80 if (localpet == 0)
then
83 print*,
"- CALL FieldScatter FOR TARGET GRID SLOPE TYPE."
84 call esmf_fieldscatter(slope_type_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
85 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
93 print*,
"- CALL FieldCreate FOR TARGET GRID MAXIMUM SNOW ALBEDO."
94 mxsno_albedo_target_grid = esmf_fieldcreate(target_grid, &
95 typekind=esmf_typekind_r8, &
96 staggerloc=esmf_staggerloc_center, rc=error)
97 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
100 do tile = 1, num_tiles_target_grid
101 if (localpet == 0)
then
102 call
read_static_file(
'maximum_snow_albedo', i_target, j_target, tile, data_one_tile)
104 print*,
"- CALL FieldScatter FOR TARGET GRID MAXIMUM SNOW ALBEDO."
105 call esmf_fieldscatter(mxsno_albedo_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
106 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
114 print*,
"- CALL FieldCreate FOR TARGET GRID SOIL TYPE."
115 soil_type_target_grid = esmf_fieldcreate(target_grid, &
116 typekind=esmf_typekind_r8, &
117 staggerloc=esmf_staggerloc_center, rc=error)
118 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
121 do tile = 1, num_tiles_target_grid
122 if (localpet == 0)
then
125 print*,
"- CALL FieldScatter FOR TARGET GRID SOIL TYPE."
126 call esmf_fieldscatter(soil_type_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
127 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
135 print*,
"- CALL FieldCreate FOR TARGET GRID VEGETATION TYPE."
136 veg_type_target_grid = esmf_fieldcreate(target_grid, &
137 typekind=esmf_typekind_r8, &
138 staggerloc=esmf_staggerloc_center, rc=error)
139 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
142 do tile = 1, num_tiles_target_grid
143 if (localpet == 0)
then
144 call
read_static_file(
'vegetation_type', i_target, j_target, tile, data_one_tile)
146 print*,
"- CALL FieldScatter FOR TARGET GRID VEGETATION TYPE."
147 call esmf_fieldscatter(veg_type_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
148 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
156 print*,
"- CALL FieldCreate FOR TARGET GRID VEGETATION GREENNESS."
157 veg_greenness_target_grid = esmf_fieldcreate(target_grid, &
158 typekind=esmf_typekind_r8, &
159 staggerloc=esmf_staggerloc_center, rc=error)
160 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
163 print*,
"- CALL FieldCreate FOR TARGET GRID MAXIMUM VEGETATION GREENNESS."
164 max_veg_greenness_target_grid = esmf_fieldcreate(target_grid, &
165 typekind=esmf_typekind_r8, &
166 staggerloc=esmf_staggerloc_center, rc=error)
167 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
170 print*,
"- CALL FieldCreate FOR TARGET GRID MINIMUM VEGETATION GREENNESS."
171 min_veg_greenness_target_grid = esmf_fieldcreate(target_grid, &
172 typekind=esmf_typekind_r8, &
173 staggerloc=esmf_staggerloc_center, rc=error)
174 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
177 if (localpet == 0)
then
178 allocate(max_data_one_tile(i_target,j_target))
179 allocate(min_data_one_tile(i_target,j_target))
181 allocate(max_data_one_tile(0,0))
182 allocate(min_data_one_tile(0,0))
185 do tile = 1, num_tiles_target_grid
186 if (localpet == 0)
then
187 call
read_static_file(
'vegetation_greenness', i_target, j_target, tile, data_one_tile, &
188 max_data_one_tile, min_data_one_tile)
190 print*,
"- CALL FieldScatter FOR TARGET GRID VEGETATION GREENNESS."
191 call esmf_fieldscatter(veg_greenness_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
192 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
194 print*,
"- CALL FieldScatter FOR TARGET GRID MAXIMUM VEGETATION GREENNESS."
195 call esmf_fieldscatter(max_veg_greenness_target_grid, max_data_one_tile, rootpet=0, tile=tile, rc=error)
196 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
198 print*,
"- CALL FieldScatter FOR TARGET GRID MINIMUM VEGETATION GREENNESS."
199 call esmf_fieldscatter(min_veg_greenness_target_grid, min_data_one_tile, rootpet=0, tile=tile, rc=error)
200 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
204 deallocate(max_data_one_tile, min_data_one_tile)
210 print*,
"- CALL FieldCreate FOR TARGET GRID SUBSTRATE TEMPERATURE."
211 substrate_temp_target_grid = esmf_fieldcreate(target_grid, &
212 typekind=esmf_typekind_r8, &
213 staggerloc=esmf_staggerloc_center, rc=error)
214 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
217 do tile = 1, num_tiles_target_grid
218 if (localpet == 0)
then
219 call
read_static_file(
'substrate_temperature', i_target, j_target, tile, data_one_tile)
221 print*,
"- CALL FieldScatter FOR TARGET GRID SUBSTRATE TEMPERATURE."
222 call esmf_fieldscatter(substrate_temp_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
223 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
231 print*,
"- CALL FieldCreate FOR ALVSF."
232 alvsf_target_grid = esmf_fieldcreate(target_grid, &
233 typekind=esmf_typekind_r8, &
234 staggerloc=esmf_staggerloc_center, rc=error)
235 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
238 do tile = 1, num_tiles_target_grid
239 if (localpet == 0)
then
240 call
read_static_file(
'visible_black_sky_albedo', i_target, j_target, tile, data_one_tile)
242 print*,
"- CALL FieldScatter FOR TARGET GRID ALVSF."
243 call esmf_fieldscatter(alvsf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
244 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
248 print*,
"- CALL FieldCreate FOR ALVWF."
249 alvwf_target_grid = esmf_fieldcreate(target_grid, &
250 typekind=esmf_typekind_r8, &
251 staggerloc=esmf_staggerloc_center, rc=error)
252 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
255 do tile = 1, num_tiles_target_grid
256 if (localpet == 0)
then
257 call
read_static_file(
'visible_white_sky_albedo', i_target, j_target, tile, data_one_tile)
259 print*,
"- CALL FieldScatter FOR TARGET GRID ALVWF."
260 call esmf_fieldscatter(alvwf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
261 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
265 print*,
"- CALL FieldCreate FOR ALNSF."
266 alnsf_target_grid = esmf_fieldcreate(target_grid, &
267 typekind=esmf_typekind_r8, &
268 staggerloc=esmf_staggerloc_center, rc=error)
269 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
272 do tile = 1, num_tiles_target_grid
273 if (localpet == 0)
then
274 call
read_static_file(
'near_IR_black_sky_albedo', i_target, j_target, tile, data_one_tile)
276 print*,
"- CALL FieldScatter FOR TARGET GRID ALNSF."
277 call esmf_fieldscatter(alnsf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
278 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
282 print*,
"- CALL FieldCreate FOR ALNWF."
283 alnwf_target_grid = esmf_fieldcreate(target_grid, &
284 typekind=esmf_typekind_r8, &
285 staggerloc=esmf_staggerloc_center, rc=error)
286 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
289 do tile = 1, num_tiles_target_grid
290 if (localpet == 0)
then
291 call
read_static_file(
'near_IR_white_sky_albedo', i_target, j_target, tile, data_one_tile)
293 print*,
"- CALL FieldScatter FOR TARGET GRID ALNWF."
294 call esmf_fieldscatter(alnwf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
295 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
303 print*,
"- CALL FieldCreate FOR TARGET GRID FACSF."
304 facsf_target_grid = esmf_fieldcreate(target_grid, &
305 typekind=esmf_typekind_r8, &
306 staggerloc=esmf_staggerloc_center, rc=error)
307 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
310 print*,
"- CALL FieldCreate FOR TARGET GRID FACWF."
311 facwf_target_grid = esmf_fieldcreate(target_grid, &
312 typekind=esmf_typekind_r8, &
313 staggerloc=esmf_staggerloc_center, rc=error)
314 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
317 do tile = 1, num_tiles_target_grid
318 if (localpet == 0)
then
321 print*,
"- CALL FieldScatter FOR TARGET GRID FACSF."
322 call esmf_fieldscatter(facsf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
323 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
325 if (localpet == 0)
then
328 if (data_one_tile(i,j) >= 0.0)
then
329 data_one_tile(i,j) = 1.0 - data_one_tile(i,j)
334 call esmf_fieldscatter(facwf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
335 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
339 deallocate(data_one_tile)
356 data_one_tile, max_data_one_tile, &
361 use program_setup, only : fix_dir_target_grid, cres_target_grid, &
362 cycle_mon, cycle_day, cycle_hour
366 character(len=*),
intent(in) :: field
367 character(len=100) :: filename
368 character(len=500) :: the_file
370 integer,
intent(in) :: i_target, j_target, tile
372 real(esmf_kind_r8),
intent(out) :: data_one_tile(i_target,j_target)
373 real(esmf_kind_r8),
intent(out),
optional :: max_data_one_tile(i_target,j_target)
374 real(esmf_kind_r8),
intent(out),
optional :: min_data_one_tile(i_target,j_target)
376 integer :: bound1, bound2
377 integer :: error, ncid, id_var, n
378 integer :: i, j, id_time, num_times
379 integer :: idat(8), jdat(8)
380 integer,
allocatable :: days_since(:)
382 real(kind=4),
allocatable :: dummy(:,:,:)
383 real(esmf_kind_r8) :: num_days, num_days_rec1, rinc(5)
384 real(esmf_kind_r8) :: weight_rec1, weight_rec2
386 if (trim(field) ==
'facsf') filename =
"/" // trim(cres_target_grid) //
".facsf." // trim(tiles_target_grid(tile)) //
".nc"
387 if (trim(field) ==
'maximum_snow_albedo') filename =
"/" // trim(cres_target_grid) //
".maximum_snow_albedo." // trim(tiles_target_grid(tile)) //
".nc"
388 if (trim(field) ==
'slope_type') filename =
"/" // trim(cres_target_grid) //
".slope_type." // trim(tiles_target_grid(tile)) //
".nc"
389 if (trim(field) ==
'soil_type') filename =
"/" // trim(cres_target_grid) //
".soil_type." // trim(tiles_target_grid(tile)) //
".nc"
390 if (trim(field) ==
'substrate_temperature') filename =
"/" // trim(cres_target_grid) //
".substrate_temperature." // trim(tiles_target_grid(tile)) //
".nc"
391 if (trim(field) ==
'vegetation_greenness') filename =
"/" // trim(cres_target_grid) //
".vegetation_greenness." // trim(tiles_target_grid(tile)) //
".nc"
392 if (trim(field) ==
'vegetation_type') filename =
"/" // trim(cres_target_grid) //
".vegetation_type." // trim(tiles_target_grid(tile)) //
".nc"
393 if (trim(field) ==
'visible_black_sky_albedo') filename =
"/" // trim(cres_target_grid) //
".snowfree_albedo." // trim(tiles_target_grid(tile)) //
".nc"
394 if (trim(field) ==
'visible_white_sky_albedo') filename =
"/" // trim(cres_target_grid) //
".snowfree_albedo." // trim(tiles_target_grid(tile)) //
".nc"
395 if (trim(field) ==
'near_IR_black_sky_albedo') filename =
"/" // trim(cres_target_grid) //
".snowfree_albedo." // trim(tiles_target_grid(tile)) //
".nc"
396 if (trim(field) ==
'near_IR_white_sky_albedo') filename =
"/" // trim(cres_target_grid) //
".snowfree_albedo." // trim(tiles_target_grid(tile)) //
".nc"
398 the_file = trim(fix_dir_target_grid) // trim(filename)
400 print*,
'- OPEN FILE ',trim(the_file)
401 error=nf90_open(trim(the_file),nf90_nowrite,ncid)
402 call
netcdf_err(error,
'OPENING: '//trim(the_file) )
404 error=nf90_inq_dimid(ncid,
'time', id_time)
406 error=nf90_inquire_dimension(ncid, id_time, len=num_times)
407 call
netcdf_err(error,
'READING TIME DIMENSION')
408 print*,
'- FILE CONTAINS ', num_times,
' TIME RECORDS.'
410 allocate(dummy(i_target,j_target,num_times))
411 error=nf90_inq_varid(ncid, field, id_var)
413 error=nf90_get_var(ncid, id_var, dummy)
416 if (num_times > 1)
then
417 allocate (days_since(num_times))
418 error=nf90_inq_varid(ncid,
'time', id_time)
419 error=nf90_get_var(ncid, id_time, days_since)
420 print*,
'- TIME RECORDS (DAYS SINCE): ', days_since
431 call w3difdat(jdat,idat,1,rinc)
433 if (rinc(1) <= days_since(n))
exit
437 if (bound1 == 0) bound1 = num_times
438 if (bound2 == num_times+1) bound2 = 1
439 print*,
"- BOUNDING TIME RECORDS: ", bound1, bound2
440 if (bound2 /= 1)
then
441 num_days = float(days_since(bound2)) - float(days_since(bound1))
442 num_days_rec1 = rinc(1) - float(days_since(bound1))
443 weight_rec2 = num_days_rec1 / num_days
444 weight_rec1 = 1.0 - weight_rec2
445 print*,
"- BOUNDING WEIGHTS ", weight_rec1, weight_rec2
447 num_days = (float(days_since(bound2)) + 1.0) + (365.0 - float(days_since(bound1)) - 1.0)
448 if (rinc(1) >= days_since(bound1))
then
449 num_days_rec1 = rinc(1) - float(days_since(bound1))
451 num_days_rec1 = (365.0 - float(days_since(bound1))) + rinc(1)
453 weight_rec2 = num_days_rec1 / num_days
454 weight_rec1 = 1.0 - weight_rec2
455 print*,
"- BOUNDING WEIGHTS ", weight_rec1, weight_rec2
460 data_one_tile(i,j) = (weight_rec1*dummy(i,j,bound1)) + (weight_rec2*dummy(i,j,bound2))
464 deallocate(days_since)
468 data_one_tile = dummy(:,:,1)
472 if (trim(field) ==
'vegetation_greenness')
then
476 max_data_one_tile(i,j) = maxval(dummy(i,j,:))
477 min_data_one_tile(i,j) = minval(dummy(i,j,:))
485 error = nf90_close(ncid)
498 print*,
"- DESTROY STATIC FIELDS."
500 call esmf_fielddestroy(alvsf_target_grid, rc=rc)
501 call esmf_fielddestroy(alvwf_target_grid, rc=rc)
502 call esmf_fielddestroy(alnsf_target_grid, rc=rc)
503 call esmf_fielddestroy(alnwf_target_grid, rc=rc)
504 call esmf_fielddestroy(facsf_target_grid, rc=rc)
505 call esmf_fielddestroy(facwf_target_grid, rc=rc)
506 call esmf_fielddestroy(max_veg_greenness_target_grid, rc=rc)
507 call esmf_fielddestroy(min_veg_greenness_target_grid, rc=rc)
508 call esmf_fielddestroy(mxsno_albedo_target_grid, rc=rc)
509 call esmf_fielddestroy(slope_type_target_grid, rc=rc)
510 call esmf_fielddestroy(soil_type_target_grid, rc=rc)
511 call esmf_fielddestroy(substrate_temp_target_grid, rc=rc)
512 call esmf_fielddestroy(veg_greenness_target_grid, rc=rc)
513 call esmf_fielddestroy(veg_type_target_grid, rc=rc)
subroutine, public get_static_fields(localpet)
Driver routine to read/time interpolate static/climo fields on the fv3 target grid.
Sets up the ESMF grid objects for the input data grid and target FV3 grid.
subroutine, public cleanup_static_fields
Free up memory for fields in this module.
subroutine netcdf_err(err, string)
Error handler for netcdf.
subroutine error_handler(string, rc)
General error handler.
This module contains code to read the setup namelist file, handle the varmap file for GRIB2 data...
subroutine read_static_file(field, i_target, j_target, tile, data_one_tile, max_data_one_tile, min_data_one_tile)
Read static climatological data file.
Reads static surface climatological data for the target FV3 grid (such as soil type and vegetation ty...