chgres_cube  1.4.0
 All Data Structures Files Functions Variables
static_data.F90
Go to the documentation of this file.
1 
4 
13  module static_data
14 
15  use esmf
16 
17  implicit none
18 
19  private
20 
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
35 
36  public :: get_static_fields
37  public :: cleanup_static_fields
38 
39  contains
40 
46  subroutine get_static_fields(localpet)
47 
48  use model_grid, only : target_grid, &
49  num_tiles_target_grid, &
50  i_target, j_target
51 
52  implicit none
53 
54  integer, intent(in) :: localpet
55 
56  integer :: error, tile, i, j
57 
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(:,:)
61 
62  if (localpet==0) then
63  allocate(data_one_tile(i_target,j_target))
64  else
65  allocate(data_one_tile(0,0))
66  endif
67 
68 !------------------------------------------------------------------------------
69 ! Slope type
70 !------------------------------------------------------------------------------
71 
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__)) &
77  call error_handler("IN FieldCreate", error)
78 
79  do tile = 1, num_tiles_target_grid
80  if (localpet == 0) then
81  call read_static_file('slope_type', i_target, j_target, tile, data_one_tile)
82  endif
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__)) &
86  call error_handler("IN FieldScatter", error)
87  enddo
88 
89 !------------------------------------------------------------------------------
90 ! Maximum snow albedo.
91 !------------------------------------------------------------------------------
92 
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__)) &
98  call error_handler("IN FieldCreate", error)
99 
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)
103  endif
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__)) &
107  call error_handler("IN FieldScatter", error)
108  enddo
109 
110 !------------------------------------------------------------------------------
111 ! Soil type
112 !------------------------------------------------------------------------------
113 
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__)) &
119  call error_handler("IN FieldCreate", error)
120 
121  do tile = 1, num_tiles_target_grid
122  if (localpet == 0) then
123  call read_static_file('soil_type', i_target, j_target, tile, data_one_tile)
124  endif
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__)) &
128  call error_handler("IN FieldScatter", error)
129  enddo
130 
131 !------------------------------------------------------------------------------
132 ! Vegetation type
133 !------------------------------------------------------------------------------
134 
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__)) &
140  call error_handler("IN FieldCreate", error)
141 
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)
145  endif
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__)) &
149  call error_handler("IN FieldScatter", error)
150  enddo
151 
152 !------------------------------------------------------------------------------
153 ! Vegetation greenness
154 !------------------------------------------------------------------------------
155 
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__)) &
161  call error_handler("IN FieldCreate", error)
162 
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__)) &
168  call error_handler("IN FieldCreate", error)
169 
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__)) &
175  call error_handler("IN FieldCreate", error)
176 
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))
180  else
181  allocate(max_data_one_tile(0,0))
182  allocate(min_data_one_tile(0,0))
183  endif
184 
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)
189  endif
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__)) &
193  call error_handler("IN FieldScatter", error)
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__)) &
197  call error_handler("IN FieldScatter", error)
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__)) &
201  call error_handler("IN FieldScatter", error)
202  enddo
203 
204  deallocate(max_data_one_tile, min_data_one_tile)
205 
206 !------------------------------------------------------------------------------
207 ! Soil substrate temperature
208 !------------------------------------------------------------------------------
209 
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__)) &
215  call error_handler("IN FieldCreate", error)
216 
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)
220  endif
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__)) &
224  call error_handler("IN FieldScatter", error)
225  enddo
226 
227 !------------------------------------------------------------------------------
228 ! Four-component albedo.
229 !------------------------------------------------------------------------------
230 
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__)) &
236  call error_handler("IN FieldCreate", error)
237 
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)
241  endif
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__)) &
245  call error_handler("IN FieldScatter", error)
246  enddo
247 
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__)) &
253  call error_handler("IN FieldCreate", error)
254 
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)
258  endif
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__)) &
262  call error_handler("IN FieldScatter", error)
263  enddo
264 
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__)) &
270  call error_handler("IN FieldCreate", error)
271 
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)
275  endif
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__)) &
279  call error_handler("IN FieldScatter", error)
280  enddo
281 
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__)) &
287  call error_handler("IN FieldCreate", error)
288 
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)
292  endif
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__)) &
296  call error_handler("IN FieldScatter", error)
297  enddo
298 
299 !------------------------------------------------------------------------------
300 ! facsf and facwf
301 !------------------------------------------------------------------------------
302 
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__)) &
308  call error_handler("IN FieldCreate", error)
309 
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__)) &
315  call error_handler("IN FieldCreate", error)
316 
317  do tile = 1, num_tiles_target_grid
318  if (localpet == 0) then
319  call read_static_file('facsf', i_target, j_target, tile, data_one_tile)
320  endif
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__)) &
324  call error_handler("IN FieldScatter", error)
325  if (localpet == 0) then
326  do j = 1, j_target
327  do i = 1, i_target
328  if (data_one_tile(i,j) >= 0.0) then
329  data_one_tile(i,j) = 1.0 - data_one_tile(i,j)
330  endif
331  enddo
332  enddo
333  endif
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__)) &
336  call error_handler("IN FieldScatter", error)
337  enddo
338 
339  deallocate(data_one_tile)
340 
341  end subroutine get_static_fields
342 
355  subroutine read_static_file(field, i_target, j_target, tile, &
356  data_one_tile, max_data_one_tile, &
357  min_data_one_tile)
358 
359  use netcdf
360  use model_grid, only : tiles_target_grid
361  use program_setup, only : fix_dir_target_grid, cres_target_grid, &
362  cycle_mon, cycle_day, cycle_hour
363 
364  implicit none
365 
366  character(len=*), intent(in) :: field
367  character(len=100) :: filename
368  character(len=500) :: the_file
369 
370  integer, intent(in) :: i_target, j_target, tile
371 
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)
375 
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(:)
381 
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
385 
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"
397 
398  the_file = trim(fix_dir_target_grid) // trim(filename)
399 
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) )
403 
404  error=nf90_inq_dimid(ncid, 'time', id_time)
405  call netcdf_err(error, 'INQ TIME DIMENSION')
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.'
409 
410  allocate(dummy(i_target,j_target,num_times))
411  error=nf90_inq_varid(ncid, field, id_var)
412  call netcdf_err(error, 'READING FIELD ID' )
413  error=nf90_get_var(ncid, id_var, dummy)
414  call netcdf_err(error, 'READING FIELD' )
415 
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
421  idat = 0
422  idat(1) = 2015
423  idat(2) = 1
424  idat(3) = 1
425  idat(5) = 0
426  jdat = 0
427  jdat(1) = 2015
428  jdat(2) = cycle_mon
429  jdat(3) = cycle_day
430  jdat(5) = cycle_hour
431  call w3difdat(jdat,idat,1,rinc)
432  do n = 1, num_times
433  if (rinc(1) <= days_since(n)) exit
434  enddo
435  bound2 = n
436  bound1 = n - 1
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
446  else
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))
450  else
451  num_days_rec1 = (365.0 - float(days_since(bound1))) + rinc(1)
452  endif
453  weight_rec2 = num_days_rec1 / num_days
454  weight_rec1 = 1.0 - weight_rec2
455  print*,"- BOUNDING WEIGHTS ", weight_rec1, weight_rec2
456  endif
457 
458  do j = 1, j_target
459  do i = 1, i_target
460  data_one_tile(i,j) = (weight_rec1*dummy(i,j,bound1)) + (weight_rec2*dummy(i,j,bound2))
461  enddo
462  enddo
463 
464  deallocate(days_since)
465 
466  else ! file contains only one time record
467 
468  data_one_tile = dummy(:,:,1)
469 
470  endif
471 
472  if (trim(field) == 'vegetation_greenness') then
473 
474  do j = 1, j_target
475  do i = 1, i_target
476  max_data_one_tile(i,j) = maxval(dummy(i,j,:))
477  min_data_one_tile(i,j) = minval(dummy(i,j,:))
478  enddo
479  enddo
480 
481  endif
482 
483  deallocate(dummy)
484 
485  error = nf90_close(ncid)
486 
487  end subroutine read_static_file
488 
493 
494  implicit none
495 
496  integer :: rc
497 
498  print*,"- DESTROY STATIC FIELDS."
499 
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)
514 
515  end subroutine cleanup_static_fields
516 
517  end module static_data
subroutine, public get_static_fields(localpet)
Driver routine to read/time interpolate static/climo fields on the fv3 target grid.
Definition: static_data.F90:46
Sets up the ESMF grid objects for the input data grid and target FV3 grid.
Definition: model_grid.F90:9
subroutine, public cleanup_static_fields
Free up memory for fields in this module.
subroutine netcdf_err(err, string)
Error handler for netcdf.
Definition: utils.F90:31
subroutine error_handler(string, rc)
General error handler.
Definition: utils.F90:9
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...
Definition: static_data.F90:13