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