10 use fckit_configuration_module,
only: fckit_configuration
11 use fckit_log_module,
only: fckit_log
14 use kinds,
only : kind_real
17 use gnssro_mod_transform,
only: geometric2geop
21 use mpas_derived_types
22 use mpas_field_routines
23 use mpas_kind_types,
only: rkind
24 use mpas_pool_routines
58 type(fckit_configuration),
intent(in) :: conf
80 type(mpas_pool_type),
pointer :: mFields
81 type(mpas_pool_data_type),
pointer :: mdata, gdata
84 type(field2dreal),
pointer :: fieldr2_a, fieldr2_b
87 real(kind=kind_real),
dimension(:),
pointer :: ptrr1_a, ptrr1_b
88 real(kind=kind_real),
dimension(:,:),
pointer :: ptrr2_a, ptrr2_b
89 real(kind=kind_real),
dimension(:,:),
allocatable :: r2_a, r2_b
92 character(len=MAXVARLEN) :: geovar
93 integer :: nCells, nVertLevels, nVertLevelsP1
94 integer :: iVar, iCell, iLevel
95 real (kind=kind_real) :: lat
98 character(len=MAXVARLEN),
parameter :: &
99 MPASSfcClassifyNames(5) = &
100 [ character(len=maxvarlen) ::
'ivgtyp',
'isltyp',
'landmask',
'xice',
'snowc']
101 character(len=MAXVARLEN),
parameter :: &
102 CRTMSfcClassifyNames(7) = &
103 [var_sfc_vegtyp, var_sfc_landtyp, var_sfc_soiltyp, &
104 var_sfc_wfrac, var_sfc_lfrac, var_sfc_ifrac, var_sfc_sfrac]
105 type(mpas_pool_type),
pointer :: CRTMSfcClassifyFields => null()
106 integer,
dimension(:),
pointer :: vegtyp, soiltyp
107 real(kind=kind_real),
dimension(:),
pointer :: wfrac, lfrac, ifrac, sfrac
110 real(kind=kind_real),
allocatable :: plevels(:,:)
113 character(len=StrKIND),
pointer :: &
114 config_microp_scheme, config_radt_cld_scheme
115 logical,
pointer :: config_microp_re
118 mfields => xm % subFields
119 ncells = geom%nCellsSolve
120 nvertlevels = geom%nVertLevels
121 nvertlevelsp1 = geom%nVertLevelsP1
125 if ( any(xg%has(crtmsfcclassifynames)) )
then
127 if (.not. all(xm%has(mpassfcclassifynames)))
then
128 call abor1_ftn(
'mpasjedi_vc_model2geovars::changevar: xm must include MPASSfcClassifyNames!')
131 call da_template_pool(geom, crtmsfcclassifyfields,
size(crtmsfcclassifynames), crtmsfcclassifynames)
135 call xm%copy_to(
'ivgtyp', crtmsfcclassifyfields, var_sfc_landtyp)
139 call mpas_pool_get_array(crtmsfcclassifyfields, var_sfc_vegtyp, vegtyp)
140 call xm%get(
'ivgtyp', mdata)
146 call mpas_pool_get_array(crtmsfcclassifyfields, var_sfc_soiltyp, soiltyp)
147 call xm%get(
'isltyp', mdata)
155 call mpas_pool_get_array(crtmsfcclassifyfields, var_sfc_lfrac, lfrac)
156 call xm%get(
'landmask', mdata)
157 lfrac(1:ncells) = real(mdata%i1%array(1:ncells), rkind)
160 call mpas_pool_get_array(crtmsfcclassifyfields, var_sfc_ifrac, ifrac)
161 call xm%get(
'xice', mdata)
162 ifrac(1:ncells) = mdata%r1%array(1:ncells)
167 call mpas_pool_get_array(crtmsfcclassifyfields, var_sfc_sfrac, sfrac)
168 call xm%get(
'snowc', mdata)
169 sfrac(1:ncells) = mdata%r1%array(1:ncells)
172 call mpas_pool_get_array(crtmsfcclassifyfields, var_sfc_wfrac, wfrac)
190 allocate(plevels(1:nvertlevelsp1,1:ncells))
191 call xm%get(
'pressure', ptrr2_a)
192 call xm%get(
'surface_pressure', ptrr1_a)
194 ncells, nvertlevels, plevels)
200 geovar = trim(xg % fldnames(ivar))
202 if (geom%has_identity(geovar))
then
205 if (xm%has(geom%identity(geovar)))
then
206 call xm%copy_to(geom%identity(geovar), xg, geovar)
208 call abor1_ftn(
'mpasjedi_vc_model2geovars::changevar: '&
209 &
'state missing identity field for geovar => '//trim(geovar))
213 call xg%get(geovar, gdata)
215 select case (trim(geovar))
218 call xm%get(
'temperature', ptrr2_a)
219 call xm%get(
'spechum', ptrr2_b)
220 allocate(r2_a(1:nvertlevels, 1:ncells))
221 call q_to_w( ptrr2_b(:,1:ncells), r2_a(:,1:ncells) )
222 call tw_to_tv( ptrr2_a(:,1:ncells), r2_a(:,1:ncells), gdata%r2%array(:,1:ncells) )
226 call xm%get(
'spechum', ptrr2_a)
227 call q_to_w(ptrr2_a(:,1:ncells), gdata%r2%array(:,1:ncells))
239 gdata%r2%array(1:nvertlevelsp1,1:ncells) = plevels(1:nvertlevelsp1,1:ncells)
250 call q_fields_forward(
'qc', mfields, gdata%r2, plevels, ncells, nvertlevels)
253 call q_fields_forward(
'qi', mfields, gdata%r2, plevels, ncells, nvertlevels)
256 call q_fields_forward(
'qr', mfields, gdata%r2, plevels, ncells, nvertlevels)
259 call q_fields_forward(
'qs', mfields, gdata%r2, plevels, ncells, nvertlevels)
262 call q_fields_forward(
'qg', mfields, gdata%r2, plevels, ncells, nvertlevels)
265 call q_fields_forward(
'qh', mfields, gdata%r2, plevels, ncells, nvertlevels)
268 call mpas_pool_get_config(geom % domain % blocklist % configs,
'config_microp_re', config_microp_re)
269 if (config_microp_re)
then
270 call xm%get(
're_cloud', ptrr2_a)
273 gdata%r2%array(:,1:ncells) = 10.0_kind_real
277 call mpas_pool_get_config(geom % domain % blocklist % configs,
'config_microp_re', config_microp_re)
278 if (config_microp_re)
then
279 call xm%get(
're_ice', ptrr2_a)
282 gdata%r2%array(:,1:ncells) = 30.0_kind_real
287 call xm%get(
'qr', fieldr2_a)
288 call mpas_pool_get_config(geom % domain % blocklist % configs,
'config_microp_scheme', config_microp_scheme)
289 call mpas_pool_get_config(geom % domain % blocklist % configs,
'config_microp_re', config_microp_re)
290 if (config_microp_re)
then
291 allocate(r2_a(1:nvertlevels, 1:ncells))
292 allocate(r2_b(1:nvertlevels, 1:ncells))
293 if (trim(config_microp_scheme) ==
'mp_thompson')
then
294 call xm%get(
'nr', ptrr2_a)
295 r2_b(:,1:ncells) = ptrr2_a(:,1:ncells)
299 call xm%get(
'rho', fieldr2_b)
302 r2_b(:,1:ncells), r2_a(:,1:ncells), config_microp_scheme, &
308 gdata%r2%array(:,1:ncells) = 999.0_kind_real
312 call mpas_pool_get_config(geom % domain % blocklist % configs,
'config_microp_re', config_microp_re)
313 if (config_microp_re)
then
314 call xm%get(
're_snow', mdata)
317 gdata%r2%array(:,1:ncells) = 600.0_kind_real
322 call xm%get(
'qg', fieldr2_a)
323 call mpas_pool_get_config(geom % domain % blocklist % configs,
'config_microp_scheme', config_microp_scheme)
324 call mpas_pool_get_config(geom % domain % blocklist % configs,
'config_microp_re', config_microp_re)
325 if (config_microp_re)
then
326 allocate(r2_a(1:nvertlevels, 1:ncells))
327 call xm%get(
'rho', fieldr2_b)
328 call effectrad_graupel(fieldr2_a%array(:,1:ncells), fieldr2_b%array(:,1:ncells), &
329 r2_a(:,1:ncells), config_microp_scheme, &
334 gdata%r2%array(:,1:ncells) = 600.0_kind_real
341 gdata%r2%array(:,1:ncells) = 600.0_kind_real
344 call mpas_pool_get_config(geom % domain % blocklist % configs, &
345 'config_radt_cld_scheme', config_radt_cld_scheme)
346 call mpas_pool_get_config(geom % domain % blocklist % configs, &
347 'config_microp_scheme', config_microp_scheme)
352 else if (xm%has(
'cldfrac'))
then
353 call xm%copy_to(
'cldfrac', xg, geovar)
355 call abor1_ftn(
'mpasjedi_vc_model2geovars::changevar: cldfrac must be added to the state &
356 & variables in order to populate the var_cldfrac geovar with the MPAS diagnostic cloud &
362 allocate(r2_a(1:nvertlevels,1:ncells))
364 nvertlevels, r2_a(:,1:ncells))
367 do ilevel = 1, nvertlevels
368 call geometric2geop(lat, r2_a(ilevel,icell), gdata%r2%array(ilevel,icell))
376 nvertlevels, gdata%r2%array(:,1:ncells))
382 call geometric2geop(lat, geom%zgrid(1,icell), gdata%r1%array(icell))
385 case ( var_sfc_geomz )
386 gdata%r1%array(1:ncells) = geom%zgrid(1,1:ncells)
388 case ( var_sfc_sdepth )
389 call xm%get(
'snowh', mdata)
392 case ( var_sfc_vegfrac )
393 call xm%get(
'vegfra', mdata)
394 gdata%r1%array(1:ncells) = mdata%r1%array(1:ncells) / 100.0_kind_real
396 case ( var_sfc_soilm )
399 call xm%get(
'smois', mdata)
400 gdata%r1%array(1:ncells) = mdata%r2%array(1,1:ncells)
402 case ( var_sfc_soilt )
404 call xm%get(
'tslb', mdata)
405 gdata%r1%array(1:ncells) = mdata%r2%array(1,1:ncells)
407 case ( var_sfc_landtyp, var_sfc_vegtyp, var_sfc_soiltyp, &
408 var_sfc_wfrac, var_sfc_lfrac, var_sfc_ifrac, var_sfc_sfrac )
412 call xg%copy_from(geovar, crtmsfcclassifyfields)
414 case ( var_sfc_wspeed )
415 call xm%get(
'u10', ptrr1_a)
416 call xm%get(
'v10', ptrr1_b)
417 gdata%r1%array(1:ncells)=sqrt( ptrr1_a(1:ncells)**2 + ptrr1_b(1:ncells)**2 )
422 call abor1_ftn(
'mpasjedi_vc_model2geovars::changevar: geovar not implemented => '//trim(geovar))
430 if (
associated(crtmsfcclassifyfields))
then
431 call mpas_pool_destroy_pool(crtmsfcclassifyfields)
elemental subroutine, public q_to_w(specific_humidity, mixing_ratio)
integer function, public convert_type_soil(type_in)
subroutine, public pressure_half_to_full(pressure, zgrid, surface_pressure, nC, nV, pressure_f)
elemental subroutine, public tw_to_tv(temperature, mixing_ratio, virtual_temperature)
subroutine, public geometricz_full_to_half(zgrid_f, nC, nV, zgrid)
subroutine, public effectrad_rainwater(qr, rho, nr, re_qr, mp_scheme, ngrid, nVertLevels)
integer function, public convert_type_veg(type_in)
subroutine, public effectrad_graupel(qg, rho, re_qg, mp_scheme, ngrid, nVertLevels)
subroutine, public q_fields_forward(mqName, modelFields, qGeo, plevels, nCells, nVertLevels)
subroutine, public da_template_pool(geom, templatePool, nf, fieldnames)
Subset a pool from fields described in geom.
real(kind=kind_real), parameter mpas_jedi_thousand_kr
real(kind=kind_real), parameter mpas_jedi_zero_kr
real(kind=kind_real), parameter mpas_jedi_rad2deg_kr
real(kind=kind_real), parameter mpas_jedi_lessone_kr
character(len=3), parameter mpas_jedi_off
real(kind=kind_real), parameter mpas_jedi_million_kr
real(kind=kind_real), parameter mpas_jedi_one_kr
subroutine create(self, geom, conf)
subroutine changevar(self, geom, xm, xg)
Fortran derived type to hold MPAS field.
Fortran derived type to hold geometry definition.