chgres_cube 1.14.0
Loading...
Searching...
No Matches
static_data.F90
Go to the documentation of this file.
1
4
14
15 use esmf
16
17 use utilities, only : error_handler, netcdf_err
18
19implicit none
20
21 private
22
23 type(esmf_field), public :: alvsf_target_grid
24 type(esmf_field), public :: alvwf_target_grid
25 type(esmf_field), public :: alnsf_target_grid
26 type(esmf_field), public :: alnwf_target_grid
27 type(esmf_field), public :: facsf_target_grid
28 type(esmf_field), public :: facwf_target_grid
29 type(esmf_field), public :: max_veg_greenness_target_grid
30 type(esmf_field), public :: min_veg_greenness_target_grid
31 type(esmf_field), public :: mxsno_albedo_target_grid
32 type(esmf_field), public :: slope_type_target_grid
33 type(esmf_field), public :: soil_type_target_grid
34 type(esmf_field), public :: substrate_temp_target_grid
35 type(esmf_field), public :: veg_greenness_target_grid
36 type(esmf_field), public :: veg_type_target_grid
37
38 public :: get_static_fields
39 public :: create_static_fields
40 public :: cleanup_static_fields
41
42 contains
43
49 subroutine get_static_fields(localpet)
50
54
55 implicit none
56
57 integer, intent(in) :: localpet
58
59 integer :: error, tile, i, j
60
61 real(esmf_kind_r8), allocatable :: data_one_tile(:,:)
62 real(esmf_kind_r8), allocatable :: max_data_one_tile(:,:)
63 real(esmf_kind_r8), allocatable :: min_data_one_tile(:,:)
64 real(esmf_kind_r8), allocatable :: land_frac_target_tile(:,:)
65 if (localpet==0) then
66 allocate(data_one_tile(i_target,j_target))
67 allocate(land_frac_target_tile(i_target,j_target))
68 else
69 allocate(data_one_tile(0,0))
70 allocate(land_frac_target_tile(0,0))
71 endif
72
74
75!------------------------------------------------------------------------------
76! Slope type
77!------------------------------------------------------------------------------
78
79 do tile = 1, num_tiles_target_grid
80 call esmf_fieldgather(land_frac_target_grid, land_frac_target_tile, rootpet=0,tile=tile, rc=error)
81 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__))&
82 call error_handler("IN FieldGather", error)
83 if (localpet == 0) then
84 call read_static_file('slope_type', i_target, j_target, tile, data_one_tile, land_frac_target_tile)
85 endif
86 print*,"- CALL FieldScatter FOR TARGET GRID SLOPE TYPE."
87 call esmf_fieldscatter(slope_type_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
88 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
89 call error_handler("IN FieldScatter", error)
90 enddo
91
92!------------------------------------------------------------------------------
93! Maximum snow albedo.
94!------------------------------------------------------------------------------
95
96 do tile = 1, num_tiles_target_grid
97 call esmf_fieldgather(land_frac_target_grid, land_frac_target_tile, rootpet=0,tile=tile, rc=error)
98 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__))&
99 call error_handler("IN FieldGather", error)
100 if (localpet == 0) then
101 call read_static_file('maximum_snow_albedo', i_target, j_target, tile, data_one_tile, land_frac_target_tile)
102 endif
103 print*,"- CALL FieldScatter FOR TARGET GRID MAXIMUM SNOW ALBEDO."
104 call esmf_fieldscatter(mxsno_albedo_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
105 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
106 call error_handler("IN FieldScatter", error)
107 enddo
108
109!------------------------------------------------------------------------------
110! Soil type
111!------------------------------------------------------------------------------
112
113 do tile = 1, num_tiles_target_grid
114 call esmf_fieldgather(land_frac_target_grid, land_frac_target_tile, rootpet=0,tile=tile, rc=error)
115 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__))&
116 call error_handler("IN FieldGather", error)
117 if (localpet == 0) then
118 call read_static_file('soil_type', i_target, j_target, tile, data_one_tile, land_frac_target_tile)
119 endif
120 print*,"- CALL FieldScatter FOR TARGET GRID SOIL TYPE."
121 call esmf_fieldscatter(soil_type_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
122 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
123 call error_handler("IN FieldScatter", error)
124 enddo
125
126!------------------------------------------------------------------------------
127! Vegetation type
128!------------------------------------------------------------------------------
129
130 do tile = 1, num_tiles_target_grid
131 call esmf_fieldgather(land_frac_target_grid, land_frac_target_tile, rootpet=0,tile=tile, rc=error)
132 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__))&
133 call error_handler("IN FieldGather", error)
134 if (localpet == 0) then
135 call read_static_file('vegetation_type', i_target, j_target, tile, data_one_tile, land_frac_target_tile)
136 endif
137 print*,"- CALL FieldScatter FOR TARGET GRID VEGETATION TYPE."
138 call esmf_fieldscatter(veg_type_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
139 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
140 call error_handler("IN FieldScatter", error)
141 enddo
142
143!------------------------------------------------------------------------------
144! Vegetation greenness
145!------------------------------------------------------------------------------
146
147 if (localpet == 0) then
148 allocate(max_data_one_tile(i_target,j_target))
149 allocate(min_data_one_tile(i_target,j_target))
150 else
151 allocate(max_data_one_tile(0,0))
152 allocate(min_data_one_tile(0,0))
153 endif
154
155 do tile = 1, num_tiles_target_grid
156 call esmf_fieldgather(land_frac_target_grid, land_frac_target_tile, rootpet=0,tile=tile, rc=error)
157 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__))&
158 call error_handler("IN FieldGather", error)
159 if (localpet == 0) then
160 call read_static_file('vegetation_greenness', i_target, j_target, tile, data_one_tile, &
161 land_frac_target_tile, max_data_one_tile, min_data_one_tile)
162 endif
163 print*,"- CALL FieldScatter FOR TARGET GRID VEGETATION GREENNESS."
164 call esmf_fieldscatter(veg_greenness_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
165 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
166 call error_handler("IN FieldScatter", error)
167 print*,"- CALL FieldScatter FOR TARGET GRID MAXIMUM VEGETATION GREENNESS."
168 call esmf_fieldscatter(max_veg_greenness_target_grid, max_data_one_tile, rootpet=0, tile=tile, rc=error)
169 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
170 call error_handler("IN FieldScatter", error)
171 print*,"- CALL FieldScatter FOR TARGET GRID MINIMUM VEGETATION GREENNESS."
172 call esmf_fieldscatter(min_veg_greenness_target_grid, min_data_one_tile, rootpet=0, tile=tile, rc=error)
173 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
174 call error_handler("IN FieldScatter", error)
175 enddo
176
177 deallocate(max_data_one_tile, min_data_one_tile)
178
179!------------------------------------------------------------------------------
180! Soil substrate temperature
181!------------------------------------------------------------------------------
182
183 do tile = 1, num_tiles_target_grid
184 call esmf_fieldgather(land_frac_target_grid, land_frac_target_tile, rootpet=0,tile=tile, rc=error)
185 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__))&
186 call error_handler("IN FieldGather", error)
187 if (localpet == 0) then
188 call read_static_file('substrate_temperature', i_target, j_target, tile, data_one_tile, land_frac_target_tile)
189 endif
190 print*,"- CALL FieldScatter FOR TARGET GRID SUBSTRATE TEMPERATURE."
191 call esmf_fieldscatter(substrate_temp_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__)) &
193 call error_handler("IN FieldScatter", error)
194 enddo
195
196!------------------------------------------------------------------------------
197! Four-component albedo.
198!------------------------------------------------------------------------------
199
200 do tile = 1, num_tiles_target_grid
201 call esmf_fieldgather(land_frac_target_grid, land_frac_target_tile, rootpet=0,tile=tile, rc=error)
202 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__))&
203 call error_handler("IN FieldGather", error)
204 if (localpet == 0) then
205 call read_static_file('visible_black_sky_albedo', i_target, j_target, tile, data_one_tile, land_frac_target_tile)
206 endif
207 print*,"- CALL FieldScatter FOR TARGET GRID ALVSF."
208 call esmf_fieldscatter(alvsf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
209 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
210 call error_handler("IN FieldScatter", error)
211 enddo
212
213 do tile = 1, num_tiles_target_grid
214 call esmf_fieldgather(land_frac_target_grid, land_frac_target_tile, rootpet=0,tile=tile, rc=error)
215 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__))&
216 call error_handler("IN FieldGather", error)
217 if (localpet == 0) then
218 call read_static_file('visible_white_sky_albedo', i_target, j_target, tile, data_one_tile, land_frac_target_tile)
219 endif
220 print*,"- CALL FieldScatter FOR TARGET GRID ALVWF."
221 call esmf_fieldscatter(alvwf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
222 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
223 call error_handler("IN FieldScatter", error)
224 enddo
225
226 do tile = 1, num_tiles_target_grid
227 call esmf_fieldgather(land_frac_target_grid, land_frac_target_tile, rootpet=0,tile=tile, rc=error)
228 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__))&
229 call error_handler("IN FieldGather", error)
230 if (localpet == 0) then
231 call read_static_file('near_IR_black_sky_albedo', i_target, j_target, tile, data_one_tile, land_frac_target_tile)
232 endif
233 print*,"- CALL FieldScatter FOR TARGET GRID ALNSF."
234 call esmf_fieldscatter(alnsf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
235 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
236 call error_handler("IN FieldScatter", error)
237 enddo
238
239 do tile = 1, num_tiles_target_grid
240 call esmf_fieldgather(land_frac_target_grid, land_frac_target_tile, rootpet=0,tile=tile, rc=error)
241 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__))&
242 call error_handler("IN FieldGather", error)
243 if (localpet == 0) then
244 call read_static_file('near_IR_white_sky_albedo', i_target, j_target, tile, data_one_tile, land_frac_target_tile)
245 endif
246 print*,"- CALL FieldScatter FOR TARGET GRID ALNWF."
247 call esmf_fieldscatter(alnwf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
248 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
249 call error_handler("IN FieldScatter", error)
250 enddo
251
252!------------------------------------------------------------------------------
253! facsf and facwf
254!------------------------------------------------------------------------------
255
256 do tile = 1, num_tiles_target_grid
257 call esmf_fieldgather(land_frac_target_grid, land_frac_target_tile, rootpet=0,tile=tile, rc=error)
258 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__))&
259 call error_handler("IN FieldGather", error)
260 if (localpet == 0) then
261 call read_static_file('facsf', i_target, j_target, tile, data_one_tile, land_frac_target_tile)
262 endif
263 print*,"- CALL FieldScatter FOR TARGET GRID FACSF."
264 call esmf_fieldscatter(facsf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
265 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
266 call error_handler("IN FieldScatter", error)
267 if (localpet == 0) then
268 do j = 1, j_target
269 do i = 1, i_target
270 if (data_one_tile(i,j) >= 0.0) then
271 data_one_tile(i,j) = 1.0 - data_one_tile(i,j)
272 endif
273 enddo
274 enddo
275 endif
276 call esmf_fieldscatter(facwf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
277 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
278 call error_handler("IN FieldScatter", error)
279 enddo
280
281 deallocate(data_one_tile)
282
283 end subroutine get_static_fields
284
298 subroutine read_static_file(field, i_target, j_target, tile, &
299 data_one_tile, land_frac, max_data_one_tile, &
300 min_data_one_tile)
301
302 use netcdf
303 use model_grid, only : tiles_target_grid
306
307 implicit none
308
309 character(len=*), intent(in) :: field
310 character(len=100) :: filename
311 character(len=500) :: the_file, err_msg
312
313 integer, intent(in) :: i_target, j_target, tile
314
315 real(esmf_kind_r8), intent(in) :: land_frac(i_target,j_target)
316 real(esmf_kind_r8), intent(out) :: data_one_tile(i_target,j_target)
317 real(esmf_kind_r8), intent(out), optional :: max_data_one_tile(i_target,j_target)
318 real(esmf_kind_r8), intent(out), optional :: min_data_one_tile(i_target,j_target)
319
320 integer :: bound1, bound2
321 integer :: error, ncid, id_var, n
322 integer :: i, j, id_time, num_times
323 integer :: idat(8), jdat(8)
324 integer, allocatable :: days_since(:)
325
326 real(kind=4), allocatable :: dummy(:,:,:)
327 real(esmf_kind_r8) :: num_days, num_days_rec1, rinc(5)
328 real(esmf_kind_r8) :: weight_rec1, weight_rec2
329
330 if (trim(field) == 'facsf') filename = "/" // trim(cres_target_grid) // ".facsf." // trim(tiles_target_grid(tile)) // ".nc"
331 if (trim(field) == 'maximum_snow_albedo') filename = "/" // trim(cres_target_grid) // ".maximum_snow_albedo." // trim(tiles_target_grid(tile)) // ".nc"
332 if (trim(field) == 'slope_type') filename = "/" // trim(cres_target_grid) // ".slope_type." // trim(tiles_target_grid(tile)) // ".nc"
333 if (trim(field) == 'soil_type') filename = "/" // trim(cres_target_grid) // ".soil_type." // trim(tiles_target_grid(tile)) // ".nc"
334 if (trim(field) == 'substrate_temperature') filename = "/" // trim(cres_target_grid) // ".substrate_temperature." // trim(tiles_target_grid(tile)) // ".nc"
335 if (trim(field) == 'vegetation_greenness') filename = "/" // trim(cres_target_grid) // ".vegetation_greenness." // trim(tiles_target_grid(tile)) // ".nc"
336 if (trim(field) == 'vegetation_type') filename = "/" // trim(cres_target_grid) // ".vegetation_type." // trim(tiles_target_grid(tile)) // ".nc"
337 if (trim(field) == 'visible_black_sky_albedo') filename = "/" // trim(cres_target_grid) // ".snowfree_albedo." // trim(tiles_target_grid(tile)) // ".nc"
338 if (trim(field) == 'visible_white_sky_albedo') filename = "/" // trim(cres_target_grid) // ".snowfree_albedo." // trim(tiles_target_grid(tile)) // ".nc"
339 if (trim(field) == 'near_IR_black_sky_albedo') filename = "/" // trim(cres_target_grid) // ".snowfree_albedo." // trim(tiles_target_grid(tile)) // ".nc"
340 if (trim(field) == 'near_IR_white_sky_albedo') filename = "/" // trim(cres_target_grid) // ".snowfree_albedo." // trim(tiles_target_grid(tile)) // ".nc"
341
342 the_file = trim(fix_dir_target_grid) // trim(filename)
343
344 print*,'- OPEN FILE ',trim(the_file)
345 error=nf90_open(trim(the_file),nf90_nowrite,ncid)
346 call netcdf_err(error, 'OPENING: '//trim(the_file) )
347
348 error=nf90_inq_dimid(ncid, 'time', id_time)
349 call netcdf_err(error, 'INQ TIME DIMENSION')
350 error=nf90_inquire_dimension(ncid, id_time, len=num_times)
351 call netcdf_err(error, 'READING TIME DIMENSION')
352 print*,'- FILE CONTAINS ', num_times, ' TIME RECORDS.'
353
354 allocate(dummy(i_target,j_target,num_times))
355 error=nf90_inq_varid(ncid, field, id_var)
356 call netcdf_err(error, 'READING FIELD ID' )
357 error=nf90_get_var(ncid, id_var, dummy)
358 call netcdf_err(error, 'READING FIELD' )
359
360 if (num_times > 1) then
361 allocate (days_since(num_times))
362 error=nf90_inq_varid(ncid, 'time', id_time)
363 error=nf90_get_var(ncid, id_time, days_since)
364 print*,'- TIME RECORDS (DAYS SINCE): ', days_since
365 idat = 0
366 idat(1) = 2015
367 idat(2) = 1
368 idat(3) = 1
369 idat(5) = 0
370 jdat = 0
371 jdat(1) = 2015
372 jdat(2) = cycle_mon
373 jdat(3) = cycle_day
374 jdat(5) = cycle_hour
375 call w3difdat(jdat,idat,1,rinc)
376 do n = 1, num_times
377 if (rinc(1) <= days_since(n)) exit
378 enddo
379 bound2 = n
380 bound1 = n - 1
381 if (bound1 == 0) bound1 = num_times
382 if (bound2 == num_times+1) bound2 = 1
383 print*,"- BOUNDING TIME RECORDS: ", bound1, bound2
384 if (bound2 /= 1) then
385 num_days = float(days_since(bound2)) - float(days_since(bound1))
386 num_days_rec1 = rinc(1) - float(days_since(bound1))
387 weight_rec2 = num_days_rec1 / num_days
388 weight_rec1 = 1.0 - weight_rec2
389 print*,"- BOUNDING WEIGHTS ", weight_rec1, weight_rec2
390 else
391 num_days = (float(days_since(bound2)) + 1.0) + (365.0 - float(days_since(bound1)) - 1.0)
392 if (rinc(1) >= days_since(bound1)) then
393 num_days_rec1 = rinc(1) - float(days_since(bound1))
394 else
395 num_days_rec1 = (365.0 - float(days_since(bound1))) + rinc(1)
396 endif
397 weight_rec2 = num_days_rec1 / num_days
398 weight_rec1 = 1.0 - weight_rec2
399 print*,"- BOUNDING WEIGHTS ", weight_rec1, weight_rec2
400 endif
401
402 do j = 1, j_target
403 do i = 1, i_target
404 data_one_tile(i,j) = (weight_rec1*dummy(i,j,bound1)) + (weight_rec2*dummy(i,j,bound2))
405 enddo
406 enddo
407
408 deallocate(days_since)
409
410 else ! file contains only one time record
411
412 do j = 1, j_target
413 do i = 1, i_target
414 if(land_frac(i,j) > 0.0 .and. dummy(i,j,1) == -999) then
415 err_msg = "Detected missing data point in " // trim(filename) // ". Static data may be outdated. Please use static data created with the latest UFS_UTILS release that supports fractional land coverage"
416 call error_handler(err_msg,-1)
417 endif
418 data_one_tile(i,j) = dummy(i,j,1)
419 enddo
420 enddo
421 endif
422
423 if (trim(field) == 'vegetation_greenness') then
424
425 do j = 1, j_target
426 do i = 1, i_target
427 max_data_one_tile(i,j) = maxval(dummy(i,j,:))
428 min_data_one_tile(i,j) = minval(dummy(i,j,:))
429 enddo
430 enddo
431
432 endif
433
434 deallocate(dummy)
435
436 error = nf90_close(ncid)
437
438 end subroutine read_static_file
439
444
445 use model_grid, only : target_grid
446
447 implicit none
448
449 integer :: error
450
451 print*,"- CALL FieldCreate FOR TARGET GRID SLOPE TYPE."
452 slope_type_target_grid = esmf_fieldcreate(target_grid, &
453 typekind=esmf_typekind_r8, &
454 staggerloc=esmf_staggerloc_center, rc=error)
455 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
456 call error_handler("IN FieldCreate", error)
457
458 print*,"- CALL FieldCreate FOR TARGET GRID MAXIMUM SNOW ALBEDO."
459 mxsno_albedo_target_grid = esmf_fieldcreate(target_grid, &
460 typekind=esmf_typekind_r8, &
461 staggerloc=esmf_staggerloc_center, rc=error)
462 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
463 call error_handler("IN FieldCreate", error)
464
465 print*,"- CALL FieldCreate FOR TARGET GRID SOIL TYPE."
466 soil_type_target_grid = esmf_fieldcreate(target_grid, &
467 typekind=esmf_typekind_r8, &
468 staggerloc=esmf_staggerloc_center, rc=error)
469 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
470 call error_handler("IN FieldCreate", error)
471
472 print*,"- CALL FieldCreate FOR TARGET GRID VEGETATION TYPE."
473 veg_type_target_grid = esmf_fieldcreate(target_grid, &
474 typekind=esmf_typekind_r8, &
475 staggerloc=esmf_staggerloc_center, rc=error)
476 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
477 call error_handler("IN FieldCreate", error)
478
479 print*,"- CALL FieldCreate FOR TARGET GRID VEGETATION GREENNESS."
480 veg_greenness_target_grid = esmf_fieldcreate(target_grid, &
481 typekind=esmf_typekind_r8, &
482 staggerloc=esmf_staggerloc_center, rc=error)
483 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
484 call error_handler("IN FieldCreate", error)
485
486 print*,"- CALL FieldCreate FOR TARGET GRID MAXIMUM VEGETATION GREENNESS."
487 max_veg_greenness_target_grid = esmf_fieldcreate(target_grid, &
488 typekind=esmf_typekind_r8, &
489 staggerloc=esmf_staggerloc_center, rc=error)
490 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
491 call error_handler("IN FieldCreate", error)
492
493 print*,"- CALL FieldCreate FOR TARGET GRID MINIMUM VEGETATION GREENNESS."
494 min_veg_greenness_target_grid = esmf_fieldcreate(target_grid, &
495 typekind=esmf_typekind_r8, &
496 staggerloc=esmf_staggerloc_center, rc=error)
497 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
498 call error_handler("IN FieldCreate", error)
499
500 print*,"- CALL FieldCreate FOR TARGET GRID SUBSTRATE TEMPERATURE."
501 substrate_temp_target_grid = esmf_fieldcreate(target_grid, &
502 typekind=esmf_typekind_r8, &
503 staggerloc=esmf_staggerloc_center, rc=error)
504 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
505 call error_handler("IN FieldCreate", error)
506
507 print*,"- CALL FieldCreate FOR ALVSF."
508 alvsf_target_grid = esmf_fieldcreate(target_grid, &
509 typekind=esmf_typekind_r8, &
510 staggerloc=esmf_staggerloc_center, rc=error)
511 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
512 call error_handler("IN FieldCreate", error)
513
514 print*,"- CALL FieldCreate FOR ALVWF."
515 alvwf_target_grid = esmf_fieldcreate(target_grid, &
516 typekind=esmf_typekind_r8, &
517 staggerloc=esmf_staggerloc_center, rc=error)
518 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
519 call error_handler("IN FieldCreate", error)
520
521 print*,"- CALL FieldCreate FOR ALNSF."
522 alnsf_target_grid = esmf_fieldcreate(target_grid, &
523 typekind=esmf_typekind_r8, &
524 staggerloc=esmf_staggerloc_center, rc=error)
525 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
526 call error_handler("IN FieldCreate", error)
527
528 print*,"- CALL FieldCreate FOR ALNWF."
529 alnwf_target_grid = esmf_fieldcreate(target_grid, &
530 typekind=esmf_typekind_r8, &
531 staggerloc=esmf_staggerloc_center, rc=error)
532 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
533 call error_handler("IN FieldCreate", error)
534
535 print*,"- CALL FieldCreate FOR TARGET GRID FACSF."
536 facsf_target_grid = esmf_fieldcreate(target_grid, &
537 typekind=esmf_typekind_r8, &
538 staggerloc=esmf_staggerloc_center, rc=error)
539 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
540 call error_handler("IN FieldCreate", error)
541
542 print*,"- CALL FieldCreate FOR TARGET GRID FACWF."
543 facwf_target_grid = esmf_fieldcreate(target_grid, &
544 typekind=esmf_typekind_r8, &
545 staggerloc=esmf_staggerloc_center, rc=error)
546 if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
547 call error_handler("IN FieldCreate", error)
548
549 end subroutine create_static_fields
550
555
556 implicit none
557
558 integer :: rc
559
560 print*,"- DESTROY STATIC FIELDS."
561
562 call esmf_fielddestroy(alvsf_target_grid, rc=rc)
563 call esmf_fielddestroy(alvwf_target_grid, rc=rc)
564 call esmf_fielddestroy(alnsf_target_grid, rc=rc)
565 call esmf_fielddestroy(alnwf_target_grid, rc=rc)
566 call esmf_fielddestroy(facsf_target_grid, rc=rc)
567 call esmf_fielddestroy(facwf_target_grid, rc=rc)
568 call esmf_fielddestroy(max_veg_greenness_target_grid, rc=rc)
569 call esmf_fielddestroy(min_veg_greenness_target_grid, rc=rc)
570 call esmf_fielddestroy(mxsno_albedo_target_grid, rc=rc)
571 call esmf_fielddestroy(slope_type_target_grid, rc=rc)
572 call esmf_fielddestroy(soil_type_target_grid, rc=rc)
573 call esmf_fielddestroy(substrate_temp_target_grid, rc=rc)
574 call esmf_fielddestroy(veg_greenness_target_grid, rc=rc)
575 call esmf_fielddestroy(veg_type_target_grid, rc=rc)
576
577 end subroutine cleanup_static_fields
578
579 end module static_data
Sets up the ESMF grid objects for the input data grid and target FV3 grid.
Definition model_grid.F90:9
integer, public i_target
i dimension of each global tile, or of a nest, target grid.
type(esmf_field), public land_frac_target_grid
land fraction, target grid
integer, public num_tiles_target_grid
Number of tiles, target grid.
type(esmf_grid), public target_grid
target grid esmf grid object.
character(len=5), dimension(:), allocatable, public tiles_target_grid
Tile names of target grid.
integer, public j_target
j dimension of each global tile, or of a nest, target grid.
This module contains code to read the setup namelist file, handle the varmap file for GRIB2 data,...
integer, public cycle_mon
Cycle month.
integer, public cycle_day
Cycle day.
character(len=15), public cres_target_grid
Target grid resolution, i.e., C768.
integer, public cycle_hour
Cycle hour.
character(len=500), public fix_dir_target_grid
Directory containing target grid pre-computed fixed data (ex: soil type).
Reads static surface climatological data for the target FV3 grid (such as soil type and vegetation ty...
type(esmf_field), public veg_greenness_target_grid
vegetation greenness fraction
type(esmf_field), public veg_type_target_grid
vegetation type
subroutine read_static_file(field, i_target, j_target, tile, data_one_tile, land_frac, max_data_one_tile, min_data_one_tile)
Read static climatological data file.
type(esmf_field), public alvwf_target_grid
visible white sky albedo
subroutine, public cleanup_static_fields
Free up memory for fields in this module.
type(esmf_field), public substrate_temp_target_grid
soil subtrate temperature
type(esmf_field), public slope_type_target_grid
slope type
type(esmf_field), public alvsf_target_grid
visible black sky albedo
type(esmf_field), public alnsf_target_grid
near ir black sky albedo
type(esmf_field), public facsf_target_grid
fractional coverage for strong zenith angle dependent albedo
type(esmf_field), public facwf_target_grid
fractional coverage for weak zenith angle dependent albedo
subroutine, public get_static_fields(localpet)
Driver routine to read/time interpolate static/climo fields on the fv3 target grid.
type(esmf_field), public min_veg_greenness_target_grid
minimum annual greenness fraction
type(esmf_field), public alnwf_target_grid
near ir white sky albedo
type(esmf_field), public max_veg_greenness_target_grid
maximum annual greenness fraction
type(esmf_field), public mxsno_albedo_target_grid
maximum snow albedo
subroutine, public create_static_fields
Create ESMF fields for static target grid data.
type(esmf_field), public soil_type_target_grid
soil type