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