chgres_cube  1.13.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 
18 
19  implicit 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 
51  use model_grid, only : num_tiles_target_grid, &
52  i_target, j_target
53 
54  implicit none
55 
56  integer, intent(in) :: localpet
57 
58  integer :: error, tile, i, j
59 
60  real(esmf_kind_r8), allocatable :: data_one_tile(:,:)
61  real(esmf_kind_r8), allocatable :: max_data_one_tile(:,:)
62  real(esmf_kind_r8), allocatable :: min_data_one_tile(:,:)
63 
64  if (localpet==0) then
65  allocate(data_one_tile(i_target,j_target))
66  else
67  allocate(data_one_tile(0,0))
68  endif
69 
71 
72 !------------------------------------------------------------------------------
73 ! Slope type
74 !------------------------------------------------------------------------------
75 
76  do tile = 1, num_tiles_target_grid
77  if (localpet == 0) then
78  call read_static_file('slope_type', i_target, j_target, tile, data_one_tile)
79  endif
80  print*,"- CALL FieldScatter FOR TARGET GRID SLOPE TYPE."
81  call esmf_fieldscatter(slope_type_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
82  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
83  call error_handler("IN FieldScatter", error)
84  enddo
85 
86 !------------------------------------------------------------------------------
87 ! Maximum snow albedo.
88 !------------------------------------------------------------------------------
89 
90  do tile = 1, num_tiles_target_grid
91  if (localpet == 0) then
92  call read_static_file('maximum_snow_albedo', i_target, j_target, tile, data_one_tile)
93  endif
94  print*,"- CALL FieldScatter FOR TARGET GRID MAXIMUM SNOW ALBEDO."
95  call esmf_fieldscatter(mxsno_albedo_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
96  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
97  call error_handler("IN FieldScatter", error)
98  enddo
99 
100 !------------------------------------------------------------------------------
101 ! Soil type
102 !------------------------------------------------------------------------------
103 
104  do tile = 1, num_tiles_target_grid
105  if (localpet == 0) then
106  call read_static_file('soil_type', i_target, j_target, tile, data_one_tile)
107  endif
108  print*,"- CALL FieldScatter FOR TARGET GRID SOIL TYPE."
109  call esmf_fieldscatter(soil_type_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
110  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
111  call error_handler("IN FieldScatter", error)
112  enddo
113 
114 !------------------------------------------------------------------------------
115 ! Vegetation type
116 !------------------------------------------------------------------------------
117 
118  do tile = 1, num_tiles_target_grid
119  if (localpet == 0) then
120  call read_static_file('vegetation_type', i_target, j_target, tile, data_one_tile)
121  endif
122  print*,"- CALL FieldScatter FOR TARGET GRID VEGETATION TYPE."
123  call esmf_fieldscatter(veg_type_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
124  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
125  call error_handler("IN FieldScatter", error)
126  enddo
127 
128 !------------------------------------------------------------------------------
129 ! Vegetation greenness
130 !------------------------------------------------------------------------------
131 
132  if (localpet == 0) then
133  allocate(max_data_one_tile(i_target,j_target))
134  allocate(min_data_one_tile(i_target,j_target))
135  else
136  allocate(max_data_one_tile(0,0))
137  allocate(min_data_one_tile(0,0))
138  endif
139 
140  do tile = 1, num_tiles_target_grid
141  if (localpet == 0) then
142  call read_static_file('vegetation_greenness', i_target, j_target, tile, data_one_tile, &
143  max_data_one_tile, min_data_one_tile)
144  endif
145  print*,"- CALL FieldScatter FOR TARGET GRID VEGETATION GREENNESS."
146  call esmf_fieldscatter(veg_greenness_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
147  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
148  call error_handler("IN FieldScatter", error)
149  print*,"- CALL FieldScatter FOR TARGET GRID MAXIMUM VEGETATION GREENNESS."
150  call esmf_fieldscatter(max_veg_greenness_target_grid, max_data_one_tile, rootpet=0, tile=tile, rc=error)
151  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
152  call error_handler("IN FieldScatter", error)
153  print*,"- CALL FieldScatter FOR TARGET GRID MINIMUM VEGETATION GREENNESS."
154  call esmf_fieldscatter(min_veg_greenness_target_grid, min_data_one_tile, rootpet=0, tile=tile, rc=error)
155  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
156  call error_handler("IN FieldScatter", error)
157  enddo
158 
159  deallocate(max_data_one_tile, min_data_one_tile)
160 
161 !------------------------------------------------------------------------------
162 ! Soil substrate temperature
163 !------------------------------------------------------------------------------
164 
165  do tile = 1, num_tiles_target_grid
166  if (localpet == 0) then
167  call read_static_file('substrate_temperature', i_target, j_target, tile, data_one_tile)
168  endif
169  print*,"- CALL FieldScatter FOR TARGET GRID SUBSTRATE TEMPERATURE."
170  call esmf_fieldscatter(substrate_temp_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
171  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
172  call error_handler("IN FieldScatter", error)
173  enddo
174 
175 !------------------------------------------------------------------------------
176 ! Four-component albedo.
177 !------------------------------------------------------------------------------
178 
179  do tile = 1, num_tiles_target_grid
180  if (localpet == 0) then
181  call read_static_file('visible_black_sky_albedo', i_target, j_target, tile, data_one_tile)
182  endif
183  print*,"- CALL FieldScatter FOR TARGET GRID ALVSF."
184  call esmf_fieldscatter(alvsf_target_grid, data_one_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 FieldScatter", error)
187  enddo
188 
189  do tile = 1, num_tiles_target_grid
190  if (localpet == 0) then
191  call read_static_file('visible_white_sky_albedo', i_target, j_target, tile, data_one_tile)
192  endif
193  print*,"- CALL FieldScatter FOR TARGET GRID ALVWF."
194  call esmf_fieldscatter(alvwf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
195  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
196  call error_handler("IN FieldScatter", error)
197  enddo
198 
199  do tile = 1, num_tiles_target_grid
200  if (localpet == 0) then
201  call read_static_file('near_IR_black_sky_albedo', i_target, j_target, tile, data_one_tile)
202  endif
203  print*,"- CALL FieldScatter FOR TARGET GRID ALNSF."
204  call esmf_fieldscatter(alnsf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
205  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
206  call error_handler("IN FieldScatter", error)
207  enddo
208 
209  do tile = 1, num_tiles_target_grid
210  if (localpet == 0) then
211  call read_static_file('near_IR_white_sky_albedo', i_target, j_target, tile, data_one_tile)
212  endif
213  print*,"- CALL FieldScatter FOR TARGET GRID ALNWF."
214  call esmf_fieldscatter(alnwf_target_grid, data_one_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 FieldScatter", error)
217  enddo
218 
219 !------------------------------------------------------------------------------
220 ! facsf and facwf
221 !------------------------------------------------------------------------------
222 
223  do tile = 1, num_tiles_target_grid
224  if (localpet == 0) then
225  call read_static_file('facsf', i_target, j_target, tile, data_one_tile)
226  endif
227  print*,"- CALL FieldScatter FOR TARGET GRID FACSF."
228  call esmf_fieldscatter(facsf_target_grid, data_one_tile, rootpet=0, tile=tile, rc=error)
229  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
230  call error_handler("IN FieldScatter", error)
231  if (localpet == 0) then
232  do j = 1, j_target
233  do i = 1, i_target
234  if (data_one_tile(i,j) >= 0.0) then
235  data_one_tile(i,j) = 1.0 - data_one_tile(i,j)
236  endif
237  enddo
238  enddo
239  endif
240  call esmf_fieldscatter(facwf_target_grid, data_one_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 FieldScatter", error)
243  enddo
244 
245  deallocate(data_one_tile)
246 
247  end subroutine get_static_fields
248 
261  subroutine read_static_file(field, i_target, j_target, tile, &
262  data_one_tile, max_data_one_tile, &
263  min_data_one_tile)
264 
265  use netcdf
266  use model_grid, only : tiles_target_grid
267  use program_setup, only : fix_dir_target_grid, cres_target_grid, &
268  cycle_mon, cycle_day, cycle_hour
269 
270  implicit none
271 
272  character(len=*), intent(in) :: field
273  character(len=100) :: filename
274  character(len=500) :: the_file
275 
276  integer, intent(in) :: i_target, j_target, tile
277 
278  real(esmf_kind_r8), intent(out) :: data_one_tile(i_target,j_target)
279  real(esmf_kind_r8), intent(out), optional :: max_data_one_tile(i_target,j_target)
280  real(esmf_kind_r8), intent(out), optional :: min_data_one_tile(i_target,j_target)
281 
282  integer :: bound1, bound2
283  integer :: error, ncid, id_var, n
284  integer :: i, j, id_time, num_times
285  integer :: idat(8), jdat(8)
286  integer, allocatable :: days_since(:)
287 
288  real(kind=4), allocatable :: dummy(:,:,:)
289  real(esmf_kind_r8) :: num_days, num_days_rec1, rinc(5)
290  real(esmf_kind_r8) :: weight_rec1, weight_rec2
291 
292  if (trim(field) == 'facsf') filename = "/" // trim(cres_target_grid) // ".facsf." // trim(tiles_target_grid(tile)) // ".nc"
293  if (trim(field) == 'maximum_snow_albedo') filename = "/" // trim(cres_target_grid) // ".maximum_snow_albedo." // trim(tiles_target_grid(tile)) // ".nc"
294  if (trim(field) == 'slope_type') filename = "/" // trim(cres_target_grid) // ".slope_type." // trim(tiles_target_grid(tile)) // ".nc"
295  if (trim(field) == 'soil_type') filename = "/" // trim(cres_target_grid) // ".soil_type." // trim(tiles_target_grid(tile)) // ".nc"
296  if (trim(field) == 'substrate_temperature') filename = "/" // trim(cres_target_grid) // ".substrate_temperature." // trim(tiles_target_grid(tile)) // ".nc"
297  if (trim(field) == 'vegetation_greenness') filename = "/" // trim(cres_target_grid) // ".vegetation_greenness." // trim(tiles_target_grid(tile)) // ".nc"
298  if (trim(field) == 'vegetation_type') filename = "/" // trim(cres_target_grid) // ".vegetation_type." // trim(tiles_target_grid(tile)) // ".nc"
299  if (trim(field) == 'visible_black_sky_albedo') filename = "/" // trim(cres_target_grid) // ".snowfree_albedo." // trim(tiles_target_grid(tile)) // ".nc"
300  if (trim(field) == 'visible_white_sky_albedo') filename = "/" // trim(cres_target_grid) // ".snowfree_albedo." // trim(tiles_target_grid(tile)) // ".nc"
301  if (trim(field) == 'near_IR_black_sky_albedo') filename = "/" // trim(cres_target_grid) // ".snowfree_albedo." // trim(tiles_target_grid(tile)) // ".nc"
302  if (trim(field) == 'near_IR_white_sky_albedo') filename = "/" // trim(cres_target_grid) // ".snowfree_albedo." // trim(tiles_target_grid(tile)) // ".nc"
303 
304  the_file = trim(fix_dir_target_grid) // trim(filename)
305 
306  print*,'- OPEN FILE ',trim(the_file)
307  error=nf90_open(trim(the_file),nf90_nowrite,ncid)
308  call netcdf_err(error, 'OPENING: '//trim(the_file) )
309 
310  error=nf90_inq_dimid(ncid, 'time', id_time)
311  call netcdf_err(error, 'INQ TIME DIMENSION')
312  error=nf90_inquire_dimension(ncid, id_time, len=num_times)
313  call netcdf_err(error, 'READING TIME DIMENSION')
314  print*,'- FILE CONTAINS ', num_times, ' TIME RECORDS.'
315 
316  allocate(dummy(i_target,j_target,num_times))
317  error=nf90_inq_varid(ncid, field, id_var)
318  call netcdf_err(error, 'READING FIELD ID' )
319  error=nf90_get_var(ncid, id_var, dummy)
320  call netcdf_err(error, 'READING FIELD' )
321 
322  if (num_times > 1) then
323  allocate (days_since(num_times))
324  error=nf90_inq_varid(ncid, 'time', id_time)
325  error=nf90_get_var(ncid, id_time, days_since)
326  print*,'- TIME RECORDS (DAYS SINCE): ', days_since
327  idat = 0
328  idat(1) = 2015
329  idat(2) = 1
330  idat(3) = 1
331  idat(5) = 0
332  jdat = 0
333  jdat(1) = 2015
334  jdat(2) = cycle_mon
335  jdat(3) = cycle_day
336  jdat(5) = cycle_hour
337  call w3difdat(jdat,idat,1,rinc)
338  do n = 1, num_times
339  if (rinc(1) <= days_since(n)) exit
340  enddo
341  bound2 = n
342  bound1 = n - 1
343  if (bound1 == 0) bound1 = num_times
344  if (bound2 == num_times+1) bound2 = 1
345  print*,"- BOUNDING TIME RECORDS: ", bound1, bound2
346  if (bound2 /= 1) then
347  num_days = float(days_since(bound2)) - float(days_since(bound1))
348  num_days_rec1 = rinc(1) - float(days_since(bound1))
349  weight_rec2 = num_days_rec1 / num_days
350  weight_rec1 = 1.0 - weight_rec2
351  print*,"- BOUNDING WEIGHTS ", weight_rec1, weight_rec2
352  else
353  num_days = (float(days_since(bound2)) + 1.0) + (365.0 - float(days_since(bound1)) - 1.0)
354  if (rinc(1) >= days_since(bound1)) then
355  num_days_rec1 = rinc(1) - float(days_since(bound1))
356  else
357  num_days_rec1 = (365.0 - float(days_since(bound1))) + rinc(1)
358  endif
359  weight_rec2 = num_days_rec1 / num_days
360  weight_rec1 = 1.0 - weight_rec2
361  print*,"- BOUNDING WEIGHTS ", weight_rec1, weight_rec2
362  endif
363 
364  do j = 1, j_target
365  do i = 1, i_target
366  data_one_tile(i,j) = (weight_rec1*dummy(i,j,bound1)) + (weight_rec2*dummy(i,j,bound2))
367  enddo
368  enddo
369 
370  deallocate(days_since)
371 
372  else ! file contains only one time record
373 
374  data_one_tile = dummy(:,:,1)
375 
376  endif
377 
378  if (trim(field) == 'vegetation_greenness') then
379 
380  do j = 1, j_target
381  do i = 1, i_target
382  max_data_one_tile(i,j) = maxval(dummy(i,j,:))
383  min_data_one_tile(i,j) = minval(dummy(i,j,:))
384  enddo
385  enddo
386 
387  endif
388 
389  deallocate(dummy)
390 
391  error = nf90_close(ncid)
392 
393  end subroutine read_static_file
394 
399 
400  use model_grid, only : target_grid
401 
402  implicit none
403 
404  integer :: error
405 
406  print*,"- CALL FieldCreate FOR TARGET GRID SLOPE TYPE."
407  slope_type_target_grid = esmf_fieldcreate(target_grid, &
408  typekind=esmf_typekind_r8, &
409  staggerloc=esmf_staggerloc_center, rc=error)
410  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
411  call error_handler("IN FieldCreate", error)
412 
413  print*,"- CALL FieldCreate FOR TARGET GRID MAXIMUM SNOW ALBEDO."
414  mxsno_albedo_target_grid = esmf_fieldcreate(target_grid, &
415  typekind=esmf_typekind_r8, &
416  staggerloc=esmf_staggerloc_center, rc=error)
417  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
418  call error_handler("IN FieldCreate", error)
419 
420  print*,"- CALL FieldCreate FOR TARGET GRID SOIL TYPE."
421  soil_type_target_grid = esmf_fieldcreate(target_grid, &
422  typekind=esmf_typekind_r8, &
423  staggerloc=esmf_staggerloc_center, rc=error)
424  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
425  call error_handler("IN FieldCreate", error)
426 
427  print*,"- CALL FieldCreate FOR TARGET GRID VEGETATION TYPE."
428  veg_type_target_grid = esmf_fieldcreate(target_grid, &
429  typekind=esmf_typekind_r8, &
430  staggerloc=esmf_staggerloc_center, rc=error)
431  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
432  call error_handler("IN FieldCreate", error)
433 
434  print*,"- CALL FieldCreate FOR TARGET GRID VEGETATION GREENNESS."
435  veg_greenness_target_grid = esmf_fieldcreate(target_grid, &
436  typekind=esmf_typekind_r8, &
437  staggerloc=esmf_staggerloc_center, rc=error)
438  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
439  call error_handler("IN FieldCreate", error)
440 
441  print*,"- CALL FieldCreate FOR TARGET GRID MAXIMUM VEGETATION GREENNESS."
442  max_veg_greenness_target_grid = esmf_fieldcreate(target_grid, &
443  typekind=esmf_typekind_r8, &
444  staggerloc=esmf_staggerloc_center, rc=error)
445  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
446  call error_handler("IN FieldCreate", error)
447 
448  print*,"- CALL FieldCreate FOR TARGET GRID MINIMUM VEGETATION GREENNESS."
449  min_veg_greenness_target_grid = esmf_fieldcreate(target_grid, &
450  typekind=esmf_typekind_r8, &
451  staggerloc=esmf_staggerloc_center, rc=error)
452  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
453  call error_handler("IN FieldCreate", error)
454 
455  print*,"- CALL FieldCreate FOR TARGET GRID SUBSTRATE TEMPERATURE."
456  substrate_temp_target_grid = esmf_fieldcreate(target_grid, &
457  typekind=esmf_typekind_r8, &
458  staggerloc=esmf_staggerloc_center, rc=error)
459  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
460  call error_handler("IN FieldCreate", error)
461 
462  print*,"- CALL FieldCreate FOR ALVSF."
463  alvsf_target_grid = esmf_fieldcreate(target_grid, &
464  typekind=esmf_typekind_r8, &
465  staggerloc=esmf_staggerloc_center, rc=error)
466  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
467  call error_handler("IN FieldCreate", error)
468 
469  print*,"- CALL FieldCreate FOR ALVWF."
470  alvwf_target_grid = esmf_fieldcreate(target_grid, &
471  typekind=esmf_typekind_r8, &
472  staggerloc=esmf_staggerloc_center, rc=error)
473  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
474  call error_handler("IN FieldCreate", error)
475 
476  print*,"- CALL FieldCreate FOR ALNSF."
477  alnsf_target_grid = esmf_fieldcreate(target_grid, &
478  typekind=esmf_typekind_r8, &
479  staggerloc=esmf_staggerloc_center, rc=error)
480  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
481  call error_handler("IN FieldCreate", error)
482 
483  print*,"- CALL FieldCreate FOR ALNWF."
484  alnwf_target_grid = esmf_fieldcreate(target_grid, &
485  typekind=esmf_typekind_r8, &
486  staggerloc=esmf_staggerloc_center, rc=error)
487  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
488  call error_handler("IN FieldCreate", error)
489 
490  print*,"- CALL FieldCreate FOR TARGET GRID FACSF."
491  facsf_target_grid = esmf_fieldcreate(target_grid, &
492  typekind=esmf_typekind_r8, &
493  staggerloc=esmf_staggerloc_center, rc=error)
494  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
495  call error_handler("IN FieldCreate", error)
496 
497  print*,"- CALL FieldCreate FOR TARGET GRID FACWF."
498  facwf_target_grid = esmf_fieldcreate(target_grid, &
499  typekind=esmf_typekind_r8, &
500  staggerloc=esmf_staggerloc_center, rc=error)
501  if(esmf_logfounderror(rctocheck=error,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
502  call error_handler("IN FieldCreate", error)
503 
504  end subroutine create_static_fields
505 
510 
511  implicit none
512 
513  integer :: rc
514 
515  print*,"- DESTROY STATIC FIELDS."
516 
517  call esmf_fielddestroy(alvsf_target_grid, rc=rc)
518  call esmf_fielddestroy(alvwf_target_grid, rc=rc)
519  call esmf_fielddestroy(alnsf_target_grid, rc=rc)
520  call esmf_fielddestroy(alnwf_target_grid, rc=rc)
521  call esmf_fielddestroy(facsf_target_grid, rc=rc)
522  call esmf_fielddestroy(facwf_target_grid, rc=rc)
523  call esmf_fielddestroy(max_veg_greenness_target_grid, rc=rc)
524  call esmf_fielddestroy(min_veg_greenness_target_grid, rc=rc)
525  call esmf_fielddestroy(mxsno_albedo_target_grid, rc=rc)
526  call esmf_fielddestroy(slope_type_target_grid, rc=rc)
527  call esmf_fielddestroy(soil_type_target_grid, rc=rc)
528  call esmf_fielddestroy(substrate_temp_target_grid, rc=rc)
529  call esmf_fielddestroy(veg_greenness_target_grid, rc=rc)
530  call esmf_fielddestroy(veg_type_target_grid, rc=rc)
531 
532  end subroutine cleanup_static_fields
533 
534  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:49
subroutine, public create_static_fields
Create ESMF fields for static target grid data.
subroutine netcdf_err(err, string)
Error handler for netcdf.
Definition: utils.F90:34
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 error_handler(string, rc)
General error handler.
Definition: utils.F90:12
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