sfc_climo_gen 1.14.0
Loading...
Searching...
No Matches
driver.F90
Go to the documentation of this file.
1
5
24 program driver
25
26 use model_grid
27 use source_grid
29 use utils
30 use esmf
31
32 implicit none
33
34 integer :: rc
35 integer :: localpet
36 integer :: npets
37
38 type(esmf_vm) :: vm
39 type(esmf_regridmethod_flag) :: method
40
41!-------------------------------------------------------------------------
42! Initialize MPI and ESMF environments.
43!-------------------------------------------------------------------------
44
45 call mpi_init(rc)
46
47 print*,"- INITIALIZE ESMF"
48 call esmf_initialize(rc=rc)
49 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
50 call error_handler("INITIALIZING ESMF", rc)
51
52 print*,"- CALL VMGetGlobal"
53 call esmf_vmgetglobal(vm, rc=rc)
54 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
55 call error_handler("IN VMGetGlobal.", rc)
56
57 print*,"- CALL VMGet"
58 call esmf_vmget(vm, localpet=localpet, petcount=npets, rc=rc)
59 if(esmf_logfounderror(rctocheck=rc,msg=esmf_logerr_passthru,line=__line__,file=__file__)) &
60 call error_handler("IN VMGet.", rc)
61
62 print*,'- NPETS IS ',npets
63 print*,'- LOCAL PET ',localpet
64
65!-------------------------------------------------------------------------
66! Program set up step.
67!-------------------------------------------------------------------------
68
69 call read_setup_namelist(localpet)
70
71!-------------------------------------------------------------------------
72! Define fv3 model grid.
73!-------------------------------------------------------------------------
74
75 call define_model_grid(localpet, npets)
76
77!-------------------------------------------------------------------------
78! Interpolate fields. Vegetation type must be done first as
79! it defines which points are permanent land ice.
80!-------------------------------------------------------------------------
81
83 if (fract_vegsoil_type) then
84 print*,'- WILL OUTPUT VEGETATION TYPE FRACTION AND DOMINANT CATEGORY.'
86 else
87 print*,'- WILL OUTPUT DOMINANT VEGETATION TYPE.'
88 method=esmf_regridmethod_nearest_stod
89 call interp(localpet, method, input_vegetation_type_file)
90 endif
92
93! Snow free albedo
94
95 if (trim(input_snowfree_albedo_file) /= "NULL") then
97 method=esmf_regridmethod_bilinear
98 if (trim(snowfree_albedo_method)=="conserve") method=esmf_regridmethod_conserve
99 call interp(localpet, method, input_snowfree_albedo_file)
101 endif
102
103! Maximum snow albedo
104
105 if (trim(input_maximum_snow_albedo_file) /= "NULL") then
107 method=esmf_regridmethod_bilinear
108 if (trim(maximum_snow_albedo_method)=="conserve") method=esmf_regridmethod_conserve
109 call interp(localpet, method, input_maximum_snow_albedo_file)
111 endif
112
113! FACSF - fractional coverage for strong zenith angle
114! dependant albedo.
115
116 if (trim(input_facsf_file) /= "NULL") then
117 call define_source_grid(localpet, npets, input_facsf_file)
118 method=esmf_regridmethod_bilinear
119 call interp(localpet, method, input_facsf_file)
121 endif
122
123! Soil substrate temperature
124
125 if (trim(input_substrate_temperature_file) /= "NULL") then
127 method=esmf_regridmethod_bilinear
128 call interp(localpet, method, input_substrate_temperature_file)
130 endif
131
132! Slope type
133
134 if (trim(input_slope_type_file) /= "NULL") then
135 call define_source_grid(localpet, npets, input_slope_type_file)
136 method=esmf_regridmethod_nearest_stod
137 call interp(localpet, method, input_slope_type_file)
139 endif
140
141! Soil type
142
143 if (trim(input_soil_type_file) /= "NULL") then
144 call define_source_grid(localpet, npets, input_soil_type_file)
145 if (fract_vegsoil_type) then
146 print*,'- WILL OUTPUT SOIL TYPE FRACTION AND DOMINANT CATEGORY.'
148 else
149 print*,'- WILL OUTPUT DOMINANT SOIL TYPE.'
150 method=esmf_regridmethod_nearest_stod
151 call interp(localpet, method, input_soil_type_file)
152 endif
154 endif
155
156! Soil color
157
158 if (trim(input_soil_color_file) /= "NULL") then
159 call define_source_grid(localpet, npets, input_soil_color_file)
160 method=esmf_regridmethod_nearest_stod
161 call interp(localpet, method, input_soil_color_file)
163 endif
164
165
166! Vegetation greenness
167
168 if (trim(input_vegetation_greenness_file) /= "NULL") then
170 method=esmf_regridmethod_bilinear
171 if (trim(vegetation_greenness_method)=="conserve") method=esmf_regridmethod_conserve
172 call interp(localpet, method, input_vegetation_greenness_file)
174 endif
175
176! Leaf Area Index
177
178 if (trim(input_leaf_area_index_file) /= "NULL") then
180 method=esmf_regridmethod_bilinear
181 if (trim(leaf_area_index_method)=="conserve") method=esmf_regridmethod_conserve
182 call interp(localpet, method, input_leaf_area_index_file)
184 endif
185
187
188 print*,"- CALL ESMF_finalize"
189 call esmf_finalize(endflag=esmf_end_keepmpi, rc=rc)
190
191 call mpi_finalize(rc)
192
193 print*,'- DONE.'
194 stop
195
196 end program driver
program driver
Reads static surface data on a global lat/lon grid, interpolates the data to the fv3 model grid,...
Definition driver.F90:24
subroutine interp(localpet, method, input_file)
Read the input source data and interpolate it to the model grid.
Definition interp.F90:14
subroutine interp_frac_cats(localpet, input_file)
Read the input source data and interpolate it to the model grid.
This module defines the model grid.
subroutine, public define_model_grid(localpet, npets)
Define model grid.
subroutine, public model_grid_cleanup
Model grid cleanup.
Set up program execution.
character(len=500), public input_substrate_temperature_file
File containing input soil substrate temperature data.
character(len=500), public input_leaf_area_index_file
File containing input leaf area index data.
logical, public fract_vegsoil_type
When true, output the percentage of each soil and vegetation type category, and the dominant category...
character(len=500), public input_maximum_snow_albedo_file
File containing input maximum snow albedo data.
character(len=50), public maximum_snow_albedo_method
Interpolation method for max snow albedo.
character(len=500), public input_vegetation_greenness_file
File containing input vegetation greenness data.
character(len=500), public input_snowfree_albedo_file
File containing input snow-free albedo data.
character(len=50), public leaf_area_index_method
Interpolation method for leaf area index.
character(len=500), public input_slope_type_file
File containing input slope type data.
character(len=500), public input_facsf_file
File containing input fractional coverage data for strong zenith angle dependent albedo.
character(len=500), public input_soil_type_file
File containing input soil type data.
subroutine, public read_setup_namelist(localpet)
Read program setup namelist.
character(len=500), public input_vegetation_type_file
File containing input vegetation type data.
character(len=50), public snowfree_albedo_method
Interpolation method for snowfree albedo.
character(len=500), public input_soil_color_file
File containing input soil color data.
character(len=50), public vegetation_greenness_method
Interpolation method for vegetation greenness.
Read grid specs, date information and land/sea mask for the source data that will be interpolated to ...
subroutine, public source_grid_cleanup
Free up memory associated with this module.
subroutine, public define_source_grid(localpet, npets, input_file)
Defines esmf grid object for source grid.
Utilities.
Definition utils.f90:8
subroutine, public error_handler(string, rc)
Handle errors.
Definition utils.f90:49