Line data Source code
1 : module atm_import_export
2 :
3 : use NUOPC , only : NUOPC_CompAttributeGet, NUOPC_Advertise, NUOPC_IsConnected
4 : use NUOPC_Model , only : NUOPC_ModelGet
5 : use ESMF , only : ESMF_GridComp, ESMF_State, ESMF_Mesh, ESMF_StateGet, ESMF_Field
6 : use ESMF , only : ESMF_Clock
7 : use ESMF , only : ESMF_KIND_R8, ESMF_SUCCESS, ESMF_MAXSTR, ESMF_LOGMSG_INFO
8 : use ESMF , only : ESMF_LogWrite, ESMF_LOGMSG_ERROR, ESMF_LogFoundError
9 : use ESMF , only : ESMF_STATEITEM_NOTFOUND, ESMF_StateItem_Flag
10 : use ESMF , only : operator(/=), operator(==)
11 : use shr_kind_mod , only : r8=>shr_kind_r8, i8=>shr_kind_i8, cl=>shr_kind_cl, cs=>shr_kind_cs, cx=>shr_kind_cx
12 : use shr_sys_mod , only : shr_sys_abort
13 : use shr_mpi_mod , only : shr_mpi_min, shr_mpi_max
14 : use nuopc_shr_methods , only : chkerr
15 : use cam_logfile , only : iulog
16 : use spmd_utils , only : masterproc, mpicom
17 : use srf_field_check , only : set_active_Sl_ram1
18 : use srf_field_check , only : set_active_Sl_fv
19 : use srf_field_check , only : set_active_Sl_soilw
20 : use srf_field_check , only : set_active_Fall_flxdst1
21 : use srf_field_check , only : set_active_Fall_flxvoc
22 : use srf_field_check , only : set_active_Fall_flxfire
23 : use srf_field_check , only : set_active_Fall_fco2_lnd
24 : use srf_field_check , only : set_active_Faoo_fco2_ocn
25 : use atm_stream_ndep , only : stream_ndep_init, stream_ndep_interp, stream_ndep_is_initialized
26 : use atm_stream_ndep , only : ndep_stream_active
27 : use chemistry , only : chem_has_ndep_flx
28 : use cam_control_mod , only : aqua_planet, simple_phys
29 :
30 : implicit none
31 : private ! except
32 :
33 : public :: read_surface_fields_namelists
34 : public :: advertise_fields
35 : public :: realize_fields
36 : public :: import_fields
37 : public :: export_fields
38 :
39 : private :: fldlist_add
40 : private :: fldlist_realize
41 : private :: state_getfldptr
42 :
43 : type fldlist_type
44 : character(len=128) :: stdname
45 : integer :: ungridded_lbound = 0
46 : integer :: ungridded_ubound = 0
47 : end type fldlist_type
48 :
49 : integer , parameter :: fldsMax = 100
50 : integer , public, protected :: fldsToAtm_num = 0
51 : integer , public, protected :: fldsFrAtm_num = 0
52 : type (fldlist_type) , public, protected :: fldsToAtm(fldsMax)
53 : type (fldlist_type) , public, protected :: fldsFrAtm(fldsMax)
54 :
55 : ! area correction factors for fluxes send and received from mediator
56 : real(r8), allocatable :: mod2med_areacor(:)
57 : real(r8), allocatable :: med2mod_areacor(:)
58 :
59 : character(len=cx) :: carma_fields = ' ' ! list of CARMA fields from lnd->atm
60 : integer :: drydep_nflds = -huge(1) ! number of dry deposition velocity fields lnd-> atm
61 : integer :: megan_nflds = -huge(1) ! number of MEGAN voc fields from lnd-> atm
62 : integer :: emis_nflds = -huge(1) ! number of fire emission fields from lnd-> atm
63 : logical :: atm_provides_lightning = .false. ! cld to grnd lightning flash freq (min-1)
64 : character(*),parameter :: F01 = "('(cam_import_export) ',a,i8,2x,i8,2x,d21.14)"
65 : character(*),parameter :: F02 = "('(cam_import_export) ',a,i8,2x,i8,2x,i8,2x,d21.14)"
66 : character(*),parameter :: u_FILE_u = __FILE__
67 :
68 : !===============================================================================
69 : contains
70 : !===============================================================================
71 :
72 : !-----------------------------------------------------------
73 : ! read mediator fields namelist file
74 : !-----------------------------------------------------------
75 2304 : subroutine read_surface_fields_namelists()
76 :
77 : use shr_drydep_mod , only : shr_drydep_readnl
78 : use shr_megan_mod , only : shr_megan_readnl
79 : use shr_fire_emis_mod , only : shr_fire_emis_readnl
80 : use shr_carma_mod , only : shr_carma_readnl
81 : use shr_lightning_coupling_mod, only : shr_lightning_coupling_readnl
82 :
83 : character(len=*), parameter :: nl_file_name = 'drv_flds_in'
84 :
85 : ! read mediator fields options
86 2304 : call shr_drydep_readnl(nl_file_name, drydep_nflds)
87 2304 : call shr_megan_readnl(nl_file_name, megan_nflds)
88 2304 : call shr_fire_emis_readnl(nl_file_name, emis_nflds)
89 2304 : call shr_carma_readnl(nl_file_name, carma_fields)
90 2304 : call shr_lightning_coupling_readnl(nl_file_name, atm_provides_lightning)
91 :
92 2304 : end subroutine read_surface_fields_namelists
93 :
94 : !-----------------------------------------------------------
95 : ! advertise fields
96 : !-----------------------------------------------------------
97 2304 : subroutine advertise_fields(gcomp, flds_scalar_name, rc)
98 :
99 : ! input/output variables
100 : type(ESMF_GridComp) :: gcomp
101 : character(len=*) , intent(in) :: flds_scalar_name
102 : integer , intent(out) :: rc
103 :
104 : ! local variables
105 : type(ESMF_State) :: importState
106 : type(ESMF_State) :: exportState
107 : character(ESMF_MAXSTR) :: stdname
108 : character(ESMF_MAXSTR) :: cvalue
109 : integer :: n, num
110 : logical :: flds_co2a ! use case
111 : logical :: flds_co2b ! use case
112 : logical :: flds_co2c ! use case
113 : character(len=128) :: fldname
114 : character(len=*), parameter :: subname='(atm_import_export:advertise_fields)'
115 : !-------------------------------------------------------------------------------
116 :
117 2304 : rc = ESMF_SUCCESS
118 :
119 2304 : call NUOPC_ModelGet(gcomp, importState=importState, exportState=exportState, rc=rc)
120 2304 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
121 :
122 : !--------------------------------
123 : ! determine necessary toggles for below
124 : !--------------------------------
125 :
126 2304 : call NUOPC_CompAttributeGet(gcomp, name='flds_co2a', value=cvalue, rc=rc)
127 2304 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
128 2304 : read(cvalue,*) flds_co2a
129 2304 : if (masterproc) write(iulog,'(a)') trim(subname)//'flds_co2a = '// trim(cvalue)
130 :
131 2304 : call NUOPC_CompAttributeGet(gcomp, name='flds_co2b', value=cvalue, rc=rc)
132 2304 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
133 2304 : read(cvalue,*) flds_co2b
134 2304 : if (masterproc) write(iulog,'(a)') trim(subname)//'flds_co2b = '// trim(cvalue)
135 :
136 2304 : call NUOPC_CompAttributeGet(gcomp, name='flds_co2c', value=cvalue, rc=rc)
137 2304 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
138 2304 : read(cvalue,*) flds_co2c
139 2304 : if (masterproc) write(iulog,'(a)') trim(subname)//'flds_co2c = '// trim(cvalue)
140 :
141 : !--------------------------------
142 : ! Export fields
143 : !--------------------------------
144 :
145 2304 : if (masterproc) write(iulog,'(a)') trim(subname)//'export_fields '
146 :
147 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, trim(flds_scalar_name))
148 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Sa_topo' )
149 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Sa_z' )
150 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Sa_u' )
151 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Sa_v' )
152 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Sa_tbot' )
153 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Sa_ptem' )
154 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Sa_shum' )
155 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Sa_pbot' )
156 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Sa_dens' )
157 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Sa_pslv' )
158 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Sa_o3' )
159 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_rainc' )
160 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_rainl' )
161 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_snowc' )
162 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_snowl' )
163 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_lwdn' )
164 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_swndr' )
165 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_swvdr' )
166 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_swndf' )
167 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_swvdf' )
168 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_swnet' ) ! only diagnostic
169 :
170 : ! from atm - black carbon deposition fluxes (3)
171 : ! (1) => bcphidry, (2) => bcphodry, (3) => bcphiwet
172 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_bcph', ungridded_lbound=1, ungridded_ubound=3)
173 :
174 : ! from atm - organic carbon deposition fluxes (3)
175 : ! (1) => ocphidry, (2) => ocphodry, (3) => ocphiwet
176 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_ocph', ungridded_lbound=1, ungridded_ubound=3)
177 :
178 : ! from atm - wet dust deposition frluxes (4 sizes)
179 : ! (1) => dstwet1, (2) => dstwet2, (3) => dstwet3, (4) => dstwet4
180 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_dstwet', ungridded_lbound=1, ungridded_ubound=4)
181 :
182 : ! from atm - dry dust deposition frluxes (4 sizes)
183 : ! (1) => dstdry1, (2) => dstdry2, (3) => dstdry3, (4) => dstdry4
184 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_dstdry', ungridded_lbound=1, ungridded_ubound=4)
185 :
186 2304 : call ESMF_LogWrite(subname//' export fields co2', ESMF_LOGMSG_INFO)
187 :
188 : ! from atm co2 fields
189 2304 : if (flds_co2a .or. flds_co2b .or. flds_co2c) then
190 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Sa_co2prog' )
191 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Sa_co2diag' )
192 : end if
193 :
194 : ! Nitrogen deposition fluxes
195 : ! Assume that 2 fields are always sent as part of Faxa_ndep
196 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Faxa_ndep', ungridded_lbound=1, ungridded_ubound=2)
197 :
198 : ! lightning flash freq
199 2304 : if (atm_provides_lightning) then
200 2304 : call fldlist_add(fldsFrAtm_num, fldsFrAtm, 'Sa_lightning')
201 : end if
202 :
203 : ! Now advertise above export fields
204 2304 : if (masterproc) write(iulog,*) trim(subname)//' advertise export fields'
205 71424 : do n = 1,fldsFrAtm_num
206 69120 : call NUOPC_Advertise(exportState, standardName=fldsFrAtm(n)%stdname, &
207 69120 : TransferOfferGeomObject='will provide', rc=rc)
208 71424 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
209 : enddo
210 :
211 : !-----------------
212 : ! Import fields
213 : !-----------------
214 :
215 2304 : if (masterproc) write(iulog,'(a)') trim(subname)//' import fields '
216 :
217 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, trim(flds_scalar_name))
218 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sx_anidr' )
219 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sx_avsdf' )
220 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sx_anidf' )
221 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sx_avsdr' )
222 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sl_lfrac' )
223 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Si_ifrac' )
224 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'So_ofrac' )
225 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sx_tref' )
226 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sx_qref' )
227 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sx_t' )
228 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'So_t' )
229 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sl_fv' ); call set_active_Sl_fv(.true.)
230 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sl_ram1' ); call set_active_Sl_ram1(.true.)
231 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sl_snowh' )
232 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Si_snowh' )
233 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'So_ssq' )
234 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'So_re' )
235 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'So_ustar' )
236 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sx_u10' )
237 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'So_ugustOut')
238 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'So_u10withGust')
239 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Faxx_taux' )
240 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Faxx_tauy' )
241 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Faxx_lat' )
242 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Faxx_sen' )
243 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Faxx_lwup' )
244 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Faxx_evap' )
245 :
246 : ! dust fluxes from land (4 sizes)
247 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Fall_flxdst', ungridded_lbound=1, ungridded_ubound=4)
248 2304 : call set_active_Fall_flxdst1(.true.)
249 :
250 : ! co2 fields from land and ocean
251 2304 : if (flds_co2b .or. flds_co2c) then
252 0 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Fall_fco2_lnd')
253 0 : call set_active_Fall_fco2_lnd(.true.)
254 : end if
255 2304 : if (flds_co2c) then
256 0 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Faoo_fco2_ocn')
257 0 : call set_active_Faoo_fco2_ocn(.true.)
258 : end if
259 :
260 : ! dry deposition velocities from land - ALSO initialize drydep here
261 2304 : if (drydep_nflds > 0) then
262 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sl_ddvel', ungridded_lbound=1, ungridded_ubound=drydep_nflds)
263 : end if
264 :
265 : ! MEGAN VOC emissions fluxes from land
266 2304 : if (megan_nflds > 0) then
267 2304 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Fall_voc', ungridded_lbound=1, ungridded_ubound=megan_nflds)
268 2304 : call set_active_Fall_flxvoc(.true.)
269 : end if
270 :
271 : ! fire emissions fluxes from land
272 2304 : if (emis_nflds > 0) then
273 0 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Fall_fire', ungridded_lbound=1, ungridded_ubound=emis_nflds)
274 0 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sl_fztop')
275 0 : call set_active_Fall_flxfire(.true.)
276 : end if
277 :
278 : ! CARMA volumetric soil water from land
279 2304 : if (carma_fields /= ' ') then
280 0 : call fldlist_add(fldsToAtm_num, fldsToAtm, 'Sl_soilw') ! optional for carma
281 0 : call set_active_Sl_soilw(.true.) ! check for carma
282 : end if
283 :
284 : ! ------------------------------------------
285 : ! Now advertise above import fields
286 : ! ------------------------------------------
287 2304 : call ESMF_LogWrite(trim(subname)//' advertise import fields ', ESMF_LOGMSG_INFO)
288 73728 : do n = 1,fldsToAtm_num
289 71424 : call NUOPC_Advertise(importState, standardName=fldsToAtm(n)%stdname, &
290 71424 : TransferOfferGeomObject='will provide', rc=rc)
291 73728 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
292 : enddo
293 :
294 4608 : end subroutine advertise_fields
295 :
296 : !===============================================================================
297 :
298 2304 : subroutine realize_fields(gcomp, Emesh, flds_scalar_name, flds_scalar_num, single_column, rc)
299 :
300 : use ESMF , only : ESMF_MeshGet, ESMF_StateGet
301 : use ESMF , only : ESMF_FieldRegridGetArea,ESMF_FieldGet
302 : use ppgrid , only : pcols, begchunk, endchunk
303 : use phys_grid , only : get_area_all_p, get_ncols_p
304 :
305 : ! input/output variables
306 : type(ESMF_GridComp) , intent(inout) :: gcomp
307 : type(ESMF_Mesh) , intent(in) :: Emesh
308 : character(len=*) , intent(in) :: flds_scalar_name
309 : integer , intent(in) :: flds_scalar_num
310 : logical , intent(in) :: single_column
311 : integer , intent(out) :: rc
312 :
313 : ! local variables
314 : type(ESMF_State) :: importState
315 : type(ESMF_State) :: exportState
316 : type(ESMF_Field) :: lfield
317 : integer :: numOwnedElements
318 : integer :: c,i,n,ncols
319 2304 : real(r8), allocatable :: mesh_areas(:)
320 2304 : real(r8), allocatable :: model_areas(:)
321 2304 : real(r8), allocatable :: area(:)
322 2304 : real(r8), pointer :: dataptr(:)
323 : real(r8) :: max_mod2med_areacor
324 : real(r8) :: max_med2mod_areacor
325 : real(r8) :: min_mod2med_areacor
326 : real(r8) :: min_med2mod_areacor
327 : real(r8) :: max_mod2med_areacor_glob
328 : real(r8) :: max_med2mod_areacor_glob
329 : real(r8) :: min_mod2med_areacor_glob
330 : real(r8) :: min_med2mod_areacor_glob
331 : character(len=cl) :: cvalue
332 : character(len=cl) :: mesh_atm
333 : character(len=cl) :: mesh_lnd
334 : character(len=cl) :: mesh_ocn
335 : logical :: samegrid_atm_lnd_ocn
336 : character(len=*), parameter :: subname='(atm_import_export:realize_fields)'
337 : !---------------------------------------------------------------------------
338 :
339 2304 : rc = ESMF_SUCCESS
340 :
341 2304 : call ESMF_LogWrite(subname//' called', ESMF_LOGMSG_INFO)
342 :
343 2304 : call NUOPC_ModelGet(gcomp, importState=importState, exportState=exportState, rc=rc)
344 2304 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
345 :
346 : call fldlist_realize( &
347 : state=ExportState, &
348 : fldList=fldsFrAtm, &
349 : numflds=fldsFrAtm_num, &
350 : flds_scalar_name=flds_scalar_name, &
351 : flds_scalar_num=flds_scalar_num, &
352 : tag=subname//':camExport',&
353 2304 : mesh=Emesh, rc=rc)
354 2304 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
355 :
356 : call fldlist_realize( &
357 : state=importState, &
358 : fldList=fldsToAtm, &
359 : numflds=fldsToAtm_num, &
360 : flds_scalar_name=flds_scalar_name, &
361 : flds_scalar_num=flds_scalar_num, &
362 : tag=subname//':camImport',&
363 2304 : mesh=Emesh, rc=rc)
364 2304 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
365 :
366 : ! Determine if atm/lnd/ocn are on the same grid - if so set area correction factors to 1
367 2304 : call NUOPC_CompAttributeGet(gcomp, name='mesh_atm', value=mesh_atm, rc=rc)
368 2304 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
369 2304 : call NUOPC_CompAttributeGet(gcomp, name='mesh_lnd', value=mesh_lnd, rc=rc)
370 2304 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
371 2304 : call NUOPC_CompAttributeGet(gcomp, name='mesh_ocn', value=mesh_ocn, rc=rc)
372 2304 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
373 :
374 2304 : samegrid_atm_lnd_ocn = .false.
375 : if ( trim(mesh_lnd) /= 'UNSET' .and. trim(mesh_atm) == trim(mesh_lnd) .and. &
376 2304 : trim(mesh_ocn) /= 'UNSET' .and. trim(mesh_atm) == trim(mesh_ocn)) then
377 : samegrid_atm_lnd_ocn = .true.
378 0 : elseif ( trim(mesh_lnd) == 'UNSET' .and. trim(mesh_atm) == trim(mesh_ocn)) then
379 : samegrid_atm_lnd_ocn = .true.
380 0 : elseif ( trim(mesh_ocn) == 'UNSET' .and. trim(mesh_atm) == trim(mesh_lnd)) then
381 0 : samegrid_atm_lnd_ocn = .true.
382 : end if
383 :
384 : ! allocate area correction factors
385 2304 : call ESMF_MeshGet(Emesh, numOwnedElements=numOwnedElements, rc=rc)
386 2304 : if (chkerr(rc,__LINE__,u_FILE_u)) return
387 6912 : allocate (mod2med_areacor(numOwnedElements))
388 4608 : allocate (med2mod_areacor(numOwnedElements))
389 :
390 2304 : if (single_column .or. samegrid_atm_lnd_ocn) then
391 :
392 148104 : mod2med_areacor(:) = 1._r8
393 148104 : med2mod_areacor(:) = 1._r8
394 :
395 : else
396 :
397 : ! Determine areas for regridding
398 0 : call ESMF_StateGet(exportState, itemName=trim(fldsFrAtm(2)%stdname), field=lfield, rc=rc)
399 0 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
400 0 : call ESMF_FieldRegridGetArea(lfield, rc=rc)
401 0 : if (chkerr(rc,__LINE__,u_FILE_u)) return
402 0 : call ESMF_FieldGet(lfield, farrayPtr=dataptr, rc=rc)
403 0 : if (chkerr(rc,__LINE__,u_FILE_u)) return
404 0 : allocate(mesh_areas(numOwnedElements))
405 0 : mesh_areas(:) = dataptr(:)
406 :
407 : ! Determine model areas
408 0 : allocate(model_areas(numOwnedElements))
409 0 : allocate(area(numOwnedElements))
410 0 : n = 0
411 0 : do c = begchunk, endchunk
412 0 : ncols = get_ncols_p(c)
413 0 : call get_area_all_p(c, ncols, area)
414 0 : do i = 1,ncols
415 0 : n = n + 1
416 0 : model_areas(n) = area(i)
417 : end do
418 : end do
419 0 : deallocate(area)
420 :
421 : ! Determine flux correction factors (module variables)
422 0 : do n = 1,numOwnedElements
423 0 : mod2med_areacor(n) = model_areas(n) / mesh_areas(n)
424 0 : med2mod_areacor(n) = 1._r8 / mod2med_areacor(n)
425 : end do
426 0 : deallocate(model_areas)
427 0 : deallocate(mesh_areas)
428 :
429 : end if
430 :
431 150408 : min_mod2med_areacor = minval(mod2med_areacor)
432 150408 : max_mod2med_areacor = maxval(mod2med_areacor)
433 150408 : min_med2mod_areacor = minval(med2mod_areacor)
434 150408 : max_med2mod_areacor = maxval(med2mod_areacor)
435 2304 : call shr_mpi_max(max_mod2med_areacor, max_mod2med_areacor_glob, mpicom)
436 2304 : call shr_mpi_min(min_mod2med_areacor, min_mod2med_areacor_glob, mpicom)
437 2304 : call shr_mpi_max(max_med2mod_areacor, max_med2mod_areacor_glob, mpicom)
438 2304 : call shr_mpi_min(min_med2mod_areacor, min_med2mod_areacor_glob, mpicom)
439 :
440 2304 : if (masterproc) then
441 3 : write(iulog,'(2A,2g23.15,A )') trim(subname),' : min_mod2med_areacor, max_mod2med_areacor ',&
442 6 : min_mod2med_areacor_glob, max_mod2med_areacor_glob, 'CAM'
443 3 : write(iulog,'(2A,2g23.15,A )') trim(subname),' : min_med2mod_areacor, max_med2mod_areacor ',&
444 6 : min_med2mod_areacor_glob, max_med2mod_areacor_glob, 'CAM'
445 : end if
446 :
447 2304 : call ESMF_LogWrite(trim(subname)//' done', ESMF_LOGMSG_INFO)
448 :
449 9216 : end subroutine realize_fields
450 :
451 : !===============================================================================
452 :
453 46080 : subroutine import_fields( gcomp, cam_in, restart_init, rc)
454 :
455 : ! -----------------------------------------------------
456 : ! Set field pointers in import state and
457 : ! copy from field pointer to chunk array data structure
458 : ! -----------------------------------------------------
459 :
460 2304 : use camsrfexch , only : cam_in_t
461 : use phys_grid , only : get_ncols_p
462 : use ppgrid , only : begchunk, endchunk
463 : use shr_const_mod , only : shr_const_stebol
464 : use co2_cycle , only : c_i, co2_readFlux_ocn, co2_readFlux_fuel
465 : use co2_cycle , only : co2_transport, co2_time_interp_ocn, co2_time_interp_fuel
466 : use co2_cycle , only : data_flux_ocn, data_flux_fuel
467 : use physconst , only : mwco2
468 : use time_manager , only : is_first_step, get_nstep
469 :
470 : ! input/output variabes
471 : type(ESMF_GridComp) :: gcomp
472 : type(cam_in_t) , intent(inout) :: cam_in(begchunk:endchunk)
473 : logical, optional , intent(in) :: restart_init
474 : integer , intent(out) :: rc
475 :
476 : ! local variables
477 : type(ESMF_State) :: importState
478 : integer :: i,n,c,g, num ! indices
479 : integer :: ncols ! number of columns
480 : integer :: nstep
481 : logical :: overwrite_flds
482 : logical :: exists
483 : logical :: exists_fco2_ocn
484 : logical :: exists_fco2_lnd
485 : character(len=128) :: fldname
486 23040 : real(r8), pointer :: fldptr2d(:,:)
487 23040 : real(r8), pointer :: fldptr1d(:)
488 23040 : real(r8), pointer :: fldptr_lat(:)
489 23040 : real(r8), pointer :: fldptr_lwup(:)
490 23040 : real(r8), pointer :: fldptr_avsdr(:)
491 23040 : real(r8), pointer :: fldptr_anidr(:)
492 23040 : real(r8), pointer :: fldptr_avsdf(:)
493 23040 : real(r8), pointer :: fldptr_anidf(:)
494 23040 : real(r8), pointer :: fldptr_tsurf(:)
495 23040 : real(r8), pointer :: fldptr_tocn(:)
496 23040 : real(r8), pointer :: fldptr_tref(:)
497 23040 : real(r8), pointer :: fldptr_qref(:)
498 23040 : real(r8), pointer :: fldptr_u10(:)
499 23040 : real(r8), pointer :: fldptr_snowhland(:)
500 23040 : real(r8), pointer :: fldptr_snowhice(:)
501 23040 : real(r8), pointer :: fldptr_ifrac(:)
502 23040 : real(r8), pointer :: fldptr_ofrac(:)
503 23040 : real(r8), pointer :: fldptr_lfrac(:)
504 23040 : real(r8), pointer :: fldptr_taux(:)
505 23040 : real(r8), pointer :: fldptr_tauy(:)
506 23040 : real(r8), pointer :: fldptr_sen(:)
507 23040 : real(r8), pointer :: fldptr_evap(:)
508 : logical, save :: first_time = .true.
509 : character(len=*), parameter :: subname='(atm_import_export:import_fields)'
510 : !---------------------------------------------------------------------------
511 :
512 23040 : rc = ESMF_SUCCESS
513 :
514 : ! Get import state
515 23040 : call NUOPC_ModelGet(gcomp, importState=importState, rc=rc)
516 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
517 :
518 : ! don't overwrite fields if invoked during the initialization phase
519 : ! of a 'continue' or 'branch' run type with data from .rs file
520 23040 : overwrite_flds = .true.
521 23040 : if (present(restart_init)) overwrite_flds = .not. restart_init
522 :
523 : !--------------------------
524 : ! Required atmosphere input fields
525 : !--------------------------
526 :
527 768 : if (overwrite_flds) then
528 22272 : call state_getfldptr(importState, 'Faxx_taux', fldptr=fldptr_taux, rc=rc)
529 22272 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
530 22272 : call state_getfldptr(importState, 'Faxx_tauy', fldptr=fldptr_tauy, rc=rc)
531 22272 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
532 22272 : call state_getfldptr(importState, 'Faxx_sen' , fldptr=fldptr_sen, rc=rc)
533 22272 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
534 22272 : call state_getfldptr(importState, 'Faxx_evap', fldptr=fldptr_evap, rc=rc)
535 22272 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
536 22272 : g = 1
537 112056 : do c = begchunk,endchunk
538 1521456 : do i = 1,get_ncols_p(c)
539 1409400 : cam_in(c)%wsx(i) = -fldptr_taux(g) * med2mod_areacor(g)
540 1409400 : cam_in(c)%wsy(i) = -fldptr_tauy(g) * med2mod_areacor(g)
541 1409400 : cam_in(c)%shf(i) = -fldptr_sen(g) * med2mod_areacor(g)
542 1409400 : cam_in(c)%cflx(i,1) = -fldptr_evap(g) * med2mod_areacor(g)
543 1499184 : g = g + 1
544 : end do
545 : end do
546 : end if ! end of overwrite_flds
547 :
548 23040 : call state_getfldptr(importState, 'Faxx_lat', fldptr=fldptr_lat, rc=rc)
549 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
550 23040 : call state_getfldptr(importState, 'Faxx_lwup', fldptr=fldptr_lwup, rc=rc)
551 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
552 23040 : call state_getfldptr(importState, 'Sx_avsdr', fldptr=fldptr_avsdr, rc=rc)
553 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
554 23040 : call state_getfldptr(importState, 'Sx_anidr', fldptr=fldptr_anidr, rc=rc)
555 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
556 23040 : call state_getfldptr(importState, 'Sx_avsdf', fldptr=fldptr_avsdf, rc=rc)
557 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
558 23040 : call state_getfldptr(importState, 'Sx_anidf', fldptr=fldptr_anidf, rc=rc)
559 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
560 23040 : call state_getfldptr(importState, 'Sx_t', fldptr=fldptr_tsurf, rc=rc)
561 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
562 23040 : call state_getfldptr(importState, 'So_t', fldptr=fldptr_tocn, rc=rc)
563 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
564 23040 : call state_getfldptr(importState, 'Sl_snowh', fldptr=fldptr_snowhland, rc=rc)
565 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
566 23040 : call state_getfldptr(importState, 'Si_snowh', fldptr=fldptr_snowhice, rc=rc)
567 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
568 23040 : call state_getfldptr(importState, 'Sx_tref', fldptr=fldptr_tref, rc=rc)
569 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
570 23040 : call state_getfldptr(importState, 'Sx_qref', fldptr=fldptr_qref, rc=rc)
571 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
572 23040 : call state_getfldptr(importState, 'Sx_u10', fldptr=fldptr_u10, rc=rc)
573 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
574 23040 : call state_getfldptr(importState, 'Si_ifrac', fldptr=fldptr_ifrac, rc=rc)
575 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
576 23040 : call state_getfldptr(importState, 'So_ofrac', fldptr=fldptr_ofrac, rc=rc)
577 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
578 23040 : call state_getfldptr(importState, 'Sl_lfrac', fldptr=fldptr_lfrac, rc=rc)
579 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
580 :
581 : ! Only do area correction on fluxes
582 23040 : g = 1
583 115920 : do c = begchunk,endchunk
584 1573920 : do i = 1,get_ncols_p(c)
585 1458000 : cam_in(c)%lhf(i) = -fldptr_lat(g) * med2mod_areacor(g)
586 1458000 : cam_in(c)%lwup(i) = -fldptr_lwup(g) * med2mod_areacor(g)
587 1458000 : cam_in(c)%asdir(i) = fldptr_avsdr(g)
588 1458000 : cam_in(c)%aldir(i) = fldptr_anidr(g)
589 1458000 : cam_in(c)%asdif(i) = fldptr_avsdf(g)
590 1458000 : cam_in(c)%aldif(i) = fldptr_anidf(g)
591 1458000 : cam_in(c)%ts(i) = fldptr_tsurf(g)
592 1458000 : cam_in(c)%sst(i) = fldptr_tocn(g)
593 1458000 : cam_in(c)%tref(i) = fldptr_tref(g)
594 1458000 : cam_in(c)%qref(i) = fldptr_qref(g)
595 1458000 : cam_in(c)%u10(i) = fldptr_u10(g)
596 1458000 : cam_in(c)%snowhland(i) = fldptr_snowhland(g)
597 1458000 : cam_in(c)%snowhice(i) = fldptr_snowhice(g)
598 1458000 : cam_in(c)%icefrac(i) = fldptr_ifrac(g)
599 1458000 : cam_in(c)%ocnfrac(i) = fldptr_ofrac(g)
600 1458000 : cam_in(c)%landfrac(i) = fldptr_lfrac(g)
601 1550880 : g = g + 1
602 : end do
603 : end do
604 :
605 : ! Optional fields
606 :
607 23040 : call state_getfldptr(importState, 'Sl_ram1', fldptr=fldptr1d, exists=exists, rc=rc)
608 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
609 23040 : if (exists) then
610 23040 : g = 1
611 115920 : do c = begchunk,endchunk
612 115920 : if ( associated(cam_in(c)%ram1) ) then
613 1550880 : do i = 1, get_ncols_p(c)
614 1458000 : cam_in(c)%ram1(i) = fldptr1d(g)
615 1550880 : g = g + 1
616 : end do
617 : end if
618 : end do
619 : end if
620 :
621 23040 : call state_getfldptr(importState, 'Sl_fv', fldptr=fldptr1d, exists=exists, rc=rc)
622 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
623 23040 : if (exists) then
624 23040 : g = 1
625 115920 : do c = begchunk,endchunk
626 115920 : if ( associated(cam_in(c)%fv) ) then
627 1550880 : do i = 1,get_ncols_p(c)
628 1458000 : cam_in(c)%fv(i) = fldptr1d(g)
629 1550880 : g = g + 1
630 : end do
631 : end if
632 : end do
633 : end if
634 :
635 : ! For CARMA - soil water from land
636 23040 : call state_getfldptr(importState, 'Sl_soilw', fldptr=fldptr1d, exists=exists, rc=rc)
637 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
638 23040 : if (exists) then
639 0 : g = 1
640 0 : do c = begchunk,endchunk
641 0 : if ( associated(cam_in(c)%soilw)) then
642 0 : do i = 1,get_ncols_p(c)
643 0 : cam_in(c)%soilw(i) = fldptr1d(g)
644 0 : g = g+1
645 : end do
646 : end if
647 : end do
648 : end if
649 :
650 : ! dry deposition fluxes from land
651 23040 : call state_getfldptr(importState, 'Fall_flxdst', fldptr2d=fldptr2d, exists=exists, rc=rc)
652 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
653 23040 : if (exists) then
654 23040 : g = 1
655 115920 : do c = begchunk,endchunk
656 115920 : if ( associated(cam_in(c)%dstflx) ) then
657 1550880 : do i = 1,get_ncols_p(c)
658 7290000 : do n = 1, size(fldptr2d, dim=1)
659 7290000 : cam_in(c)%dstflx(i,n) = fldptr2d(n,g) * med2mod_areacor(g)
660 : end do
661 1550880 : g = g + 1
662 : end do
663 : end if
664 : end do
665 : end if
666 :
667 : ! MEGAN VOC emis fluxes from land
668 23040 : call state_getfldptr(importState, 'Fall_voc', fldptr2d=fldptr2d, exists=exists, rc=rc)
669 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
670 23040 : if (exists) then
671 23040 : g = 1
672 115920 : do c=begchunk,endchunk
673 115920 : if ( associated(cam_in(c)%meganflx) ) then
674 1550880 : do i = 1,get_ncols_p(c)
675 2916000 : do n = 1, size(fldptr2d, dim=1)
676 2916000 : cam_in(c)%meganflx(i,n) = fldptr2d(n,g) * med2mod_areacor(g)
677 : end do
678 1550880 : g = g + 1
679 : end do
680 : end if
681 : end do
682 : end if
683 :
684 : ! fire emission fluxes from land
685 23040 : call state_getfldptr(importState, 'Fall_fire', fldptr2d=fldptr2d, exists=exists, rc=rc)
686 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
687 23040 : if (exists) then
688 0 : g = 1
689 0 : do c = begchunk,endchunk
690 0 : if ( associated(cam_in(c)%fireflx) .and. associated(cam_in(c)%fireztop) ) then
691 0 : do i = 1,get_ncols_p(c)
692 0 : do n = 1, size(fldptr2d, dim=1)
693 0 : cam_in(c)%fireflx(i,n) = fldptr2d(n,g) * med2mod_areacor(g)
694 : end do
695 0 : g = g + 1
696 : end do
697 : end if
698 : end do
699 : end if
700 23040 : call state_getfldptr(importState, 'Sl_fztop', fldptr=fldptr1d, exists=exists, rc=rc)
701 23040 : if (exists) then
702 0 : g = 1
703 0 : do c = begchunk,endchunk
704 0 : do i = 1,get_ncols_p(c)
705 0 : cam_in(c)%fireztop(i) = fldptr1d(g)
706 0 : g = g + 1
707 : end do
708 : end do
709 : end if
710 :
711 : ! dry dep velocities
712 23040 : call state_getfldptr(importState, 'Sl_ddvel', fldptr2d=fldptr2d, exists=exists, rc=rc)
713 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
714 23040 : if (exists) then
715 23040 : g = 1
716 115920 : do c = begchunk,endchunk
717 1573920 : do i = 1,get_ncols_p(c)
718 8748000 : do n = 1, size(fldptr2d, dim=1)
719 8748000 : cam_in(c)%depvel(i,n) = fldptr2d(n,g)
720 : end do
721 1550880 : g = g + 1
722 : end do
723 : end do
724 : end if
725 :
726 : ! fields needed to calculate water isotopes to ocean evaporation processes
727 23040 : call state_getfldptr(importState, 'So_ustar', fldptr=fldptr1d, exists=exists, rc=rc)
728 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
729 23040 : if (exists) then
730 23040 : g = 1
731 115920 : do c = begchunk,endchunk
732 1573920 : do i = 1,get_ncols_p(c)
733 1458000 : cam_in(c)%ustar(i) = fldptr1d(g)
734 1550880 : g = g + 1
735 : end do
736 : end do
737 : end if
738 23040 : call state_getfldptr(importState, 'So_re', fldptr=fldptr1d, exists=exists, rc=rc)
739 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
740 23040 : if (exists) then
741 23040 : g = 1
742 115920 : do c = begchunk,endchunk
743 1573920 : do i = 1,get_ncols_p(c)
744 1458000 : cam_in(c)%re(i)= fldptr1d(g)
745 1550880 : g = g + 1
746 : end do
747 : end do
748 : end if
749 23040 : call state_getfldptr(importState, 'So_ssq', fldptr=fldptr1d, exists=exists, rc=rc)
750 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
751 23040 : if (exists) then
752 23040 : g = 1
753 115920 : do c = begchunk,endchunk
754 1573920 : do i = 1,get_ncols_p(c)
755 1458000 : cam_in(c)%ssq(i) = fldptr1d(g)
756 1550880 : g = g + 1
757 : end do
758 : end do
759 : end if
760 :
761 23040 : call state_getfldptr(importState, 'So_ugustOut', fldptr=fldptr1d, exists=exists, rc=rc)
762 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
763 23040 : if (exists) then
764 23040 : g = 1
765 115920 : do c = begchunk,endchunk
766 1573920 : do i = 1,get_ncols_p(c)
767 1458000 : cam_in(c)%ugustOut(i) = fldptr1d(g)
768 1550880 : g = g + 1
769 : end do
770 : end do
771 : end if
772 :
773 23040 : call state_getfldptr(importState, 'So_u10withGust', fldptr=fldptr1d, exists=exists, rc=rc)
774 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
775 23040 : if (exists) then
776 23040 : g = 1
777 115920 : do c = begchunk,endchunk
778 1573920 : do i = 1,get_ncols_p(c)
779 1458000 : cam_in(c)%u10withGusts(i) = fldptr1d(g)
780 1550880 : g = g + 1
781 : end do
782 : end do
783 : end if
784 :
785 : ! bgc scenarios
786 23040 : call state_getfldptr(importState, 'Fall_fco2_lnd', fldptr=fldptr1d, exists=exists_fco2_lnd, rc=rc)
787 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
788 23040 : if (exists_fco2_lnd) then
789 0 : g = 1
790 0 : do c = begchunk,endchunk
791 0 : do i = 1,get_ncols_p(c)
792 0 : cam_in(c)%fco2_lnd(i) = -fldptr1d(g) * med2mod_areacor(g)
793 0 : g = g + 1
794 : end do
795 : end do
796 : end if
797 23040 : call state_getfldptr(importState, 'Faoo_fco2_ocn', fldptr=fldptr1d, exists=exists_fco2_ocn, rc=rc)
798 23040 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
799 23040 : if (exists_fco2_ocn) then
800 0 : g = 1
801 0 : do c = begchunk,endchunk
802 0 : do i = 1,get_ncols_p(c)
803 0 : cam_in(c)%fco2_ocn(i) = -fldptr1d(g) * med2mod_areacor(g)
804 0 : g = g + 1
805 : end do
806 : end do
807 : else
808 : ! Consistency check
809 23040 : if (co2_readFlux_ocn) then
810 0 : call shr_sys_abort(subname // ':: co2_readFlux_ocn and x2a_Faoo_fco2_ocn cannot both be active')
811 : end if
812 : end if
813 23040 : call state_getfldptr(importState, 'Faoo_dms_ocn', fldptr=fldptr1d, exists=exists, rc=rc)
814 23040 : if (exists) then
815 0 : g = 1
816 0 : do c = begchunk,endchunk
817 0 : do i = 1,get_ncols_p(c)
818 0 : cam_in(c)%fdms(i) = -fldptr1d(g) * med2mod_areacor(g)
819 0 : g = g + 1
820 : end do
821 : end do
822 : end if
823 :
824 : ! -----------------------------------
825 : ! Get total co2 flux from components,
826 : ! -----------------------------------
827 :
828 : ! Note - co2_transport determines if cam_in(c)%cflx(i,c_i(1:4)) is allocated
829 :
830 23040 : if (co2_transport() .and. overwrite_flds) then
831 :
832 : ! Interpolate in time for flux data read in
833 0 : if (co2_readFlux_ocn) then
834 0 : call co2_time_interp_ocn
835 : end if
836 0 : if (co2_readFlux_fuel) then
837 0 : call co2_time_interp_fuel
838 : end if
839 :
840 : ! from ocn : data read in or from coupler or zero
841 : ! from fuel: data read in or zero
842 : ! from lnd : through coupler or zero
843 : ! all co2 fluxes in unit kgCO2/m2/s
844 :
845 0 : do c=begchunk,endchunk
846 0 : do i=1, get_ncols_p(c)
847 :
848 : ! co2 flux from ocn
849 0 : if (exists_fco2_ocn) then
850 0 : cam_in(c)%cflx(i,c_i(1)) = cam_in(c)%fco2_ocn(i)
851 0 : else if (co2_readFlux_ocn) then
852 : ! convert from molesCO2/m2/s to kgCO2/m2/s
853 0 : cam_in(c)%cflx(i,c_i(1)) = &
854 0 : -data_flux_ocn%co2flx(i,c)*(1._r8- cam_in(c)%landfrac(i))*mwco2*1.0e-3_r8
855 : else
856 0 : cam_in(c)%cflx(i,c_i(1)) = 0._r8
857 : end if
858 :
859 : ! co2 flux from fossil fuel
860 0 : if (co2_readFlux_fuel) then
861 0 : cam_in(c)%cflx(i,c_i(2)) = data_flux_fuel%co2flx(i,c)
862 : else
863 0 : cam_in(c)%cflx(i,c_i(2)) = 0._r8
864 : end if
865 :
866 : ! co2 flux from land (cpl already multiplies flux by land fraction)
867 0 : if (exists_fco2_lnd) then
868 0 : cam_in(c)%cflx(i,c_i(3)) = cam_in(c)%fco2_lnd(i)
869 : else
870 0 : cam_in(c)%cflx(i,c_i(3)) = 0._r8
871 : end if
872 :
873 : ! merged co2 flux
874 0 : cam_in(c)%cflx(i,c_i(4)) = cam_in(c)%cflx(i,c_i(1)) + cam_in(c)%cflx(i,c_i(2)) + cam_in(c)%cflx(i,c_i(3))
875 : end do
876 : end do
877 : end if
878 :
879 : ! if first step, determine longwave up flux from the surface temperature
880 23040 : if (first_time) then
881 2304 : if (is_first_step()) then
882 7728 : do c=begchunk, endchunk
883 104928 : do i=1, get_ncols_p(c)
884 103392 : cam_in(c)%lwup(i) = shr_const_stebol*(cam_in(c)%ts(i)**4)
885 : end do
886 : end do
887 : end if
888 2304 : first_time = .false.
889 : end if
890 :
891 46080 : end subroutine import_fields
892 :
893 : !===============================================================================
894 :
895 50688 : subroutine export_fields( gcomp, model_mesh, model_clock, cam_out, rc)
896 :
897 : ! -----------------------------------------------------
898 : ! Set field pointers in export set
899 : ! Copy from chunk array data structure into state fldptr
900 : ! -----------------------------------------------------
901 :
902 23040 : use camsrfexch , only : cam_out_t
903 : use phys_grid , only : get_ncols_p
904 : use ppgrid , only : begchunk, endchunk
905 : use time_manager , only : is_first_step, get_nstep
906 : use spmd_utils , only : masterproc
907 :
908 : !-------------------------------
909 : ! Pack the export state
910 : !-------------------------------
911 :
912 : ! input/output variables
913 : type(ESMF_GridComp) :: gcomp
914 : type(ESMF_Mesh) , intent(in) :: model_mesh
915 : type(ESMF_Clock), intent(in) :: model_clock
916 : type(cam_out_t) , intent(inout) :: cam_out(begchunk:endchunk)
917 : integer , intent(out) :: rc
918 :
919 : ! local variables
920 : type(ESMF_State) :: exportState
921 : type(ESMF_Clock) :: clock
922 : integer :: i,m,c,n,g ! indices
923 : integer :: ncols ! Number of columns
924 : integer :: nstep
925 : logical :: exists
926 : ! 2d pointers
927 25344 : real(r8), pointer :: fldptr_ndep(:,:)
928 25344 : real(r8), pointer :: fldptr_bcph(:,:) , fldptr_ocph(:,:)
929 25344 : real(r8), pointer :: fldptr_dstwet(:,:), fldptr_dstdry(:,:)
930 : ! 1d pointers
931 25344 : real(r8), pointer :: fldptr_soll(:) , fldptr_sols(:)
932 25344 : real(r8), pointer :: fldptr_solld(:) , fldptr_solsd(:)
933 25344 : real(r8), pointer :: fldptr_snowc(:) , fldptr_snowl(:)
934 25344 : real(r8), pointer :: fldptr_rainc(:) , fldptr_rainl(:)
935 25344 : real(r8), pointer :: fldptr_lwdn(:) , fldptr_swnet(:)
936 25344 : real(r8), pointer :: fldptr_topo(:) , fldptr_zbot(:)
937 25344 : real(r8), pointer :: fldptr_ubot(:) , fldptr_vbot(:)
938 25344 : real(r8), pointer :: fldptr_pbot(:) , fldptr_tbot(:)
939 25344 : real(r8), pointer :: fldptr_shum(:) , fldptr_dens(:)
940 25344 : real(r8), pointer :: fldptr_ptem(:) , fldptr_pslv(:)
941 25344 : real(r8), pointer :: fldptr_co2prog(:) , fldptr_co2diag(:)
942 25344 : real(r8), pointer :: fldptr_ozone(:)
943 25344 : real(r8), pointer :: fldptr_lght(:)
944 : character(len=*), parameter :: subname='(atm_import_export:export_fields)'
945 : !---------------------------------------------------------------------------
946 :
947 25344 : rc = ESMF_SUCCESS
948 :
949 : ! Get export state
950 25344 : call NUOPC_ModelGet(gcomp, exportState=exportState, rc=rc)
951 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
952 :
953 : ! required export state variables
954 25344 : call state_getfldptr(exportState, 'Sa_topo', fldptr=fldptr_topo, rc=rc)
955 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
956 25344 : call state_getfldptr(exportState, 'Sa_z' , fldptr=fldptr_zbot, rc=rc)
957 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
958 25344 : call state_getfldptr(exportState, 'Sa_u' , fldptr=fldptr_ubot, rc=rc)
959 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
960 25344 : call state_getfldptr(exportState, 'Sa_v' , fldptr=fldptr_vbot, rc=rc)
961 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
962 25344 : call state_getfldptr(exportState, 'Sa_tbot', fldptr=fldptr_tbot, rc=rc)
963 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
964 25344 : call state_getfldptr(exportState, 'Sa_pbot', fldptr=fldptr_pbot, rc=rc)
965 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
966 25344 : call state_getfldptr(exportState, 'Sa_shum', fldptr=fldptr_shum, rc=rc)
967 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
968 25344 : call state_getfldptr(exportState, 'Sa_dens', fldptr=fldptr_dens, rc=rc)
969 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
970 25344 : call state_getfldptr(exportState, 'Sa_ptem', fldptr=fldptr_ptem, rc=rc)
971 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
972 25344 : call state_getfldptr(exportState, 'Sa_pslv', fldptr=fldptr_pslv, rc=rc)
973 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
974 25344 : g = 1
975 127512 : do c = begchunk,endchunk
976 1731312 : do i = 1,get_ncols_p(c)
977 1603800 : fldptr_topo(g) = cam_out(c)%topo(i)
978 1603800 : fldptr_zbot(g) = cam_out(c)%zbot(i)
979 1603800 : fldptr_ubot(g) = cam_out(c)%ubot(i)
980 1603800 : fldptr_vbot(g) = cam_out(c)%vbot(i)
981 1603800 : fldptr_pbot(g) = cam_out(c)%pbot(i)
982 1603800 : fldptr_tbot(g) = cam_out(c)%tbot(i)
983 1603800 : fldptr_shum(g) = cam_out(c)%qbot(i,1)
984 1603800 : fldptr_dens(g) = cam_out(c)%rho(i)
985 1603800 : fldptr_ptem(g) = cam_out(c)%thbot(i)
986 1603800 : fldptr_pslv(g) = cam_out(c)%psl(i)
987 1705968 : g = g + 1
988 : end do
989 : end do
990 :
991 : ! required export flux variables
992 25344 : call state_getfldptr(exportState, 'Faxa_swnet', fldptr=fldptr_swnet, rc=rc)
993 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
994 25344 : call state_getfldptr(exportState, 'Faxa_lwdn' , fldptr=fldptr_lwdn , rc=rc)
995 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
996 25344 : call state_getfldptr(exportState, 'Faxa_rainc', fldptr=fldptr_rainc, rc=rc)
997 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
998 25344 : call state_getfldptr(exportState, 'Faxa_rainl', fldptr=fldptr_rainl, rc=rc)
999 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1000 25344 : call state_getfldptr(exportState, 'Faxa_snowc', fldptr=fldptr_snowc, rc=rc)
1001 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1002 25344 : call state_getfldptr(exportState, 'Faxa_snowl', fldptr=fldptr_snowl, rc=rc)
1003 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1004 25344 : call state_getfldptr(exportState, 'Faxa_swndr', fldptr=fldptr_soll, rc=rc)
1005 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1006 25344 : call state_getfldptr(exportState, 'Faxa_swvdr', fldptr=fldptr_sols, rc=rc)
1007 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1008 25344 : call state_getfldptr(exportState, 'Faxa_swndf', fldptr=fldptr_solld, rc=rc)
1009 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1010 25344 : call state_getfldptr(exportState, 'Faxa_swvdf', fldptr=fldptr_solsd, rc=rc)
1011 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1012 25344 : g = 1
1013 127512 : do c = begchunk,endchunk
1014 1731312 : do i = 1,get_ncols_p(c)
1015 1603800 : fldptr_lwdn(g) = cam_out(c)%flwds(i) * mod2med_areacor(g)
1016 1603800 : fldptr_swnet(g) = cam_out(c)%netsw(i) * mod2med_areacor(g)
1017 1603800 : fldptr_snowc(g) = cam_out(c)%precsc(i)*1000._r8 * mod2med_areacor(g)
1018 1603800 : fldptr_snowl(g) = cam_out(c)%precsl(i)*1000._r8 * mod2med_areacor(g)
1019 1603800 : fldptr_rainc(g) = (cam_out(c)%precc(i) - cam_out(c)%precsc(i))*1000._r8 * mod2med_areacor(g)
1020 1603800 : fldptr_rainl(g) = (cam_out(c)%precl(i) - cam_out(c)%precsl(i))*1000._r8 * mod2med_areacor(g)
1021 1603800 : fldptr_soll(g) = cam_out(c)%soll(i) * mod2med_areacor(g)
1022 1603800 : fldptr_sols(g) = cam_out(c)%sols(i) * mod2med_areacor(g)
1023 1603800 : fldptr_solld(g) = cam_out(c)%solld(i) * mod2med_areacor(g)
1024 1603800 : fldptr_solsd(g) = cam_out(c)%solsd(i) * mod2med_areacor(g)
1025 1705968 : g = g + 1
1026 : end do
1027 : end do
1028 :
1029 : ! aerosol deposition fluxes
1030 25344 : call state_getfldptr(exportState, 'Faxa_bcph', fldptr2d=fldptr_bcph, rc=rc)
1031 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1032 25344 : call state_getfldptr(exportState, 'Faxa_ocph', fldptr2d=fldptr_ocph, rc=rc)
1033 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1034 25344 : call state_getfldptr(exportState, 'Faxa_dstdry', fldptr2d=fldptr_dstdry, rc=rc)
1035 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1036 25344 : call state_getfldptr(exportState, 'Faxa_dstwet', fldptr2d=fldptr_dstwet, rc=rc)
1037 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1038 : ! (1) => bcphidry, (2) => bcphodry, (3) => bcphiwet
1039 : ! (1) => ocphidry, (2) => ocphodry, (3) => ocphiwet
1040 25344 : g = 1
1041 127512 : do c = begchunk,endchunk
1042 1731312 : do i = 1,get_ncols_p(c)
1043 1603800 : fldptr_bcph(1,g) = cam_out(c)%bcphidry(i) * mod2med_areacor(g)
1044 1603800 : fldptr_bcph(2,g) = cam_out(c)%bcphodry(i) * mod2med_areacor(g)
1045 1603800 : fldptr_bcph(3,g) = cam_out(c)%bcphiwet(i) * mod2med_areacor(g)
1046 1603800 : fldptr_ocph(1,g) = cam_out(c)%ocphidry(i) * mod2med_areacor(g)
1047 1603800 : fldptr_ocph(2,g) = cam_out(c)%ocphodry(i) * mod2med_areacor(g)
1048 1603800 : fldptr_ocph(3,g) = cam_out(c)%ocphiwet(i) * mod2med_areacor(g)
1049 1603800 : fldptr_dstdry(1,g) = cam_out(c)%dstdry1(i) * mod2med_areacor(g)
1050 1603800 : fldptr_dstdry(2,g) = cam_out(c)%dstdry2(i) * mod2med_areacor(g)
1051 1603800 : fldptr_dstdry(3,g) = cam_out(c)%dstdry3(i) * mod2med_areacor(g)
1052 1603800 : fldptr_dstdry(4,g) = cam_out(c)%dstdry4(i) * mod2med_areacor(g)
1053 1603800 : fldptr_dstwet(1,g) = cam_out(c)%dstwet1(i) * mod2med_areacor(g)
1054 1603800 : fldptr_dstwet(2,g) = cam_out(c)%dstwet2(i) * mod2med_areacor(g)
1055 1603800 : fldptr_dstwet(3,g) = cam_out(c)%dstwet3(i) * mod2med_areacor(g)
1056 1603800 : fldptr_dstwet(4,g) = cam_out(c)%dstwet4(i) * mod2med_areacor(g)
1057 1705968 : g = g + 1
1058 : end do
1059 : end do
1060 :
1061 25344 : call state_getfldptr(exportState, 'Sa_o3', fldptr=fldptr_ozone, exists=exists, rc=rc)
1062 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1063 25344 : if (exists) then
1064 25344 : g = 1
1065 127512 : do c = begchunk,endchunk
1066 1731312 : do i = 1,get_ncols_p(c)
1067 1603800 : fldptr_ozone(g) = cam_out(c)%ozone(i) ! atm ozone
1068 1705968 : g = g + 1
1069 : end do
1070 : end do
1071 : end if
1072 :
1073 25344 : call state_getfldptr(exportState, 'Sa_lightning', fldptr=fldptr_lght, exists=exists, rc=rc)
1074 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1075 25344 : if (exists) then
1076 25344 : g = 1
1077 127512 : do c = begchunk,endchunk
1078 1731312 : do i = 1,get_ncols_p(c)
1079 1603800 : fldptr_lght(g) = cam_out(c)%lightning_flash_freq(i) ! cloud-to-ground lightning flash frequency (/min)
1080 1705968 : g = g + 1
1081 : end do
1082 : end do
1083 : end if
1084 :
1085 25344 : call state_getfldptr(exportState, 'Sa_co2prog', fldptr=fldptr_co2prog, exists=exists, rc=rc)
1086 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1087 25344 : if (exists) then
1088 25344 : g = 1
1089 127512 : do c = begchunk,endchunk
1090 1731312 : do i = 1,get_ncols_p(c)
1091 1603800 : fldptr_co2prog(g) = cam_out(c)%co2prog(i) ! atm prognostic co2
1092 1705968 : g = g + 1
1093 : end do
1094 : end do
1095 : end if
1096 :
1097 25344 : call state_getfldptr(exportState, 'Sa_co2diag', fldptr=fldptr_co2diag, exists=exists, rc=rc)
1098 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1099 25344 : if (exists) then
1100 25344 : g = 1
1101 127512 : do c = begchunk,endchunk
1102 1731312 : do i = 1,get_ncols_p(c)
1103 1603800 : fldptr_co2diag(g) = cam_out(c)%co2diag(i) ! atm diagnostic co2
1104 1705968 : g = g + 1
1105 : end do
1106 : end do
1107 : end if
1108 :
1109 25344 : call state_getfldptr(exportState, 'Faxa_ndep', fldptr2d=fldptr_ndep, rc=rc)
1110 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1111 :
1112 4836744 : fldptr_ndep(:,:) = 0._r8
1113 :
1114 25344 : if (.not. (simple_phys .or. aqua_planet)) then
1115 :
1116 : ! The ndep_stream_nl namelist group is read in stream_ndep_init. This sets whether
1117 : ! or not the stream will be used.
1118 25344 : if (.not. stream_ndep_is_initialized) then
1119 2304 : call stream_ndep_init(model_mesh, model_clock, rc)
1120 2304 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1121 2304 : stream_ndep_is_initialized = .true.
1122 : end if
1123 :
1124 25344 : if (ndep_stream_active.or.chem_has_ndep_flx) then
1125 :
1126 : ! Nitrogen dep fluxes are obtained from the ndep input stream if input data is available
1127 : ! otherwise computed by chemistry
1128 25344 : if (ndep_stream_active) then
1129 :
1130 : ! get ndep fluxes from the stream
1131 25344 : call stream_ndep_interp(cam_out, rc)
1132 25344 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1133 :
1134 : end if
1135 :
1136 25344 : g = 1
1137 127512 : do c = begchunk,endchunk
1138 1731312 : do i = 1,get_ncols_p(c)
1139 1603800 : fldptr_ndep(1,g) = cam_out(c)%nhx_nitrogen_flx(i) * mod2med_areacor(g)
1140 1603800 : fldptr_ndep(2,g) = cam_out(c)%noy_nitrogen_flx(i) * mod2med_areacor(g)
1141 1705968 : g = g + 1
1142 : end do
1143 : end do
1144 :
1145 : end if
1146 :
1147 : end if
1148 :
1149 50688 : end subroutine export_fields
1150 :
1151 : !===============================================================================
1152 :
1153 140544 : subroutine fldlist_add(num, fldlist, stdname, ungridded_lbound, ungridded_ubound)
1154 :
1155 : ! input/otuput variables
1156 : integer , intent(inout) :: num
1157 : type(fldlist_type) , intent(inout) :: fldlist(:)
1158 : character(len=*) , intent(in) :: stdname
1159 : integer, optional , intent(in) :: ungridded_lbound
1160 : integer, optional , intent(in) :: ungridded_ubound
1161 :
1162 : ! local variables
1163 : character(len=*), parameter :: subname='(atm_import_export:fldlist_add)'
1164 : !-------------------------------------------------------------------------------
1165 :
1166 : ! Set up a list of field information
1167 :
1168 140544 : num = num + 1
1169 140544 : if (num > fldsMax) then
1170 : call ESMF_LogWrite(trim(subname)//": ERROR num > fldsMax "//trim(stdname), &
1171 0 : ESMF_LOGMSG_ERROR, line=__LINE__, file=__FILE__)
1172 0 : return
1173 : endif
1174 140544 : fldlist(num)%stdname = trim(stdname)
1175 :
1176 140544 : if (present(ungridded_lbound) .and. present(ungridded_ubound)) then
1177 18432 : fldlist(num)%ungridded_lbound = ungridded_lbound
1178 18432 : fldlist(num)%ungridded_ubound = ungridded_ubound
1179 : end if
1180 :
1181 25344 : end subroutine fldlist_add
1182 :
1183 : !===============================================================================
1184 :
1185 4608 : subroutine fldlist_realize(state, fldList, numflds, flds_scalar_name, flds_scalar_num, mesh, tag, rc)
1186 :
1187 : use NUOPC , only : NUOPC_IsConnected, NUOPC_Realize
1188 : use ESMF , only : ESMF_MeshLoc_Element, ESMF_FieldCreate, ESMF_TYPEKIND_R8
1189 : use ESMF , only : ESMF_MAXSTR, ESMF_Field, ESMF_State, ESMF_Mesh, ESMF_StateRemove
1190 : use ESMF , only : ESMF_LogFoundError, ESMF_LOGMSG_INFO, ESMF_SUCCESS
1191 : use ESMF , only : ESMF_LogWrite, ESMF_LOGMSG_ERROR, ESMF_LOGERR_PASSTHRU
1192 :
1193 : ! input/output variables
1194 : type(ESMF_State) , intent(inout) :: state
1195 : type(fldlist_type) , intent(in) :: fldList(:)
1196 : integer , intent(in) :: numflds
1197 : character(len=*) , intent(in) :: flds_scalar_name
1198 : integer , intent(in) :: flds_scalar_num
1199 : character(len=*) , intent(in) :: tag
1200 : type(ESMF_Mesh) , intent(in) :: mesh
1201 : integer , intent(inout) :: rc
1202 :
1203 : ! local variables
1204 : integer :: n
1205 : type(ESMF_Field) :: field
1206 : character(len=80) :: stdname
1207 : character(CL) :: msg
1208 : character(len=*),parameter :: subname='(atm_import_export:fldlist_realize)'
1209 : ! ----------------------------------------------
1210 :
1211 4608 : rc = ESMF_SUCCESS
1212 :
1213 149760 : do n = 1, numflds
1214 140544 : stdname = fldList(n)%stdname
1215 145152 : if (NUOPC_IsConnected(state, fieldName=stdname)) then
1216 140544 : if (stdname == trim(flds_scalar_name)) then
1217 4608 : if (masterproc) then
1218 6 : write(iulog,'(a)') trim(subname)//trim(tag)//" field = "//trim(stdname)//" is connected on root pe"
1219 : end if
1220 : ! Create the scalar field
1221 4608 : call SetScalarField(field, flds_scalar_name, flds_scalar_num, rc=rc)
1222 4608 : if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, line=__LINE__, file=u_FILE_u)) return
1223 : else
1224 : ! Create the field
1225 135936 : if (fldlist(n)%ungridded_lbound > 0 .and. fldlist(n)%ungridded_ubound > 0) then
1226 : field = ESMF_FieldCreate(mesh, ESMF_TYPEKIND_R8, name=stdname, meshloc=ESMF_MESHLOC_ELEMENT, &
1227 : ungriddedLbound=(/fldlist(n)%ungridded_lbound/), &
1228 : ungriddedUbound=(/fldlist(n)%ungridded_ubound/), &
1229 55296 : gridToFieldMap=(/2/), rc=rc)
1230 18432 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1231 18432 : if (masterproc) then
1232 : write(iulog,'(a,i8,a,i8)') trim(subname)// trim(tag)//" Field = "//trim(stdname)// &
1233 24 : " is connected using mesh with lbound ", fldlist(n)%ungridded_lbound,&
1234 48 : " and with ubound ",fldlist(n)%ungridded_ubound
1235 : end if
1236 : else
1237 117504 : field = ESMF_FieldCreate(mesh, ESMF_TYPEKIND_R8, name=stdname, meshloc=ESMF_MESHLOC_ELEMENT, rc=rc)
1238 117504 : if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, line=__LINE__, file=u_FILE_u)) return
1239 117504 : if (masterproc) then
1240 153 : write(iulog,'(a)') trim(subname)// trim(tag)//" Field = "//trim(stdname)// " is connected using mesh "
1241 : end if
1242 : end if
1243 : endif
1244 :
1245 : ! NOW call NUOPC_Realize
1246 140544 : call NUOPC_Realize(state, field=field, rc=rc)
1247 140544 : if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, line=__LINE__, file=u_FILE_u)) return
1248 : else
1249 0 : if (stdname /= trim(flds_scalar_name)) then
1250 0 : if (masterproc) then
1251 0 : write(iulog,'(a)')trim(subname)//trim(tag)//" Field = "//trim(stdname)//" is not connected"
1252 : end if
1253 0 : call ESMF_StateRemove(state, (/stdname/), rc=rc)
1254 0 : if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, line=__LINE__, file=u_FILE_u)) return
1255 : end if
1256 : end if
1257 : end do
1258 :
1259 : contains !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1260 :
1261 4608 : subroutine SetScalarField(field, flds_scalar_name, flds_scalar_num, rc)
1262 : ! ----------------------------------------------
1263 : ! create a field with scalar data on the root pe
1264 : ! ----------------------------------------------
1265 :
1266 4608 : use ESMF, only : ESMF_Field, ESMF_DistGrid, ESMF_Grid
1267 : use ESMF, only : ESMF_DistGridCreate, ESMF_GridCreate, ESMF_LogFoundError, ESMF_LOGERR_PASSTHRU
1268 : use ESMF, only : ESMF_FieldCreate, ESMF_GridCreate, ESMF_TYPEKIND_R8
1269 :
1270 : ! input/output variables
1271 : type(ESMF_Field) , intent(inout) :: field
1272 : character(len=*) , intent(in) :: flds_scalar_name
1273 : integer , intent(in) :: flds_scalar_num
1274 : integer , intent(inout) :: rc
1275 :
1276 : ! local variables
1277 : type(ESMF_Distgrid) :: distgrid
1278 : type(ESMF_Grid) :: grid
1279 : character(len=*), parameter :: subname='(atm_import_export:SetScalarField)'
1280 : ! ----------------------------------------------
1281 :
1282 4608 : rc = ESMF_SUCCESS
1283 :
1284 : ! create a DistGrid with a single index space element, which gets mapped onto DE 0.
1285 4608 : distgrid = ESMF_DistGridCreate(minIndex=(/1/), maxIndex=(/1/), rc=rc)
1286 4608 : if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, line=__LINE__, file=u_FILE_u)) return
1287 :
1288 4608 : grid = ESMF_GridCreate(distgrid, rc=rc)
1289 4608 : if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, line=__LINE__, file=u_FILE_u)) return
1290 :
1291 : field = ESMF_FieldCreate(name=trim(flds_scalar_name), grid=grid, typekind=ESMF_TYPEKIND_R8, &
1292 9216 : ungriddedLBound=(/1/), ungriddedUBound=(/flds_scalar_num/), gridToFieldMap=(/2/), rc=rc)
1293 4608 : if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, line=__LINE__, file=u_FILE_u)) return
1294 :
1295 9216 : end subroutine SetScalarField
1296 :
1297 : end subroutine fldlist_realize
1298 :
1299 : !===============================================================================
1300 1561344 : subroutine state_getfldptr(State, fldname, fldptr, fldptr2d, exists, rc)
1301 :
1302 : ! ----------------------------------------------
1303 : ! Get pointer to a state field
1304 : ! ----------------------------------------------
1305 :
1306 : use ESMF , only : ESMF_State, ESMF_Field, ESMF_Mesh, ESMF_FieldStatus_Flag
1307 : use ESMF , only : ESMF_StateGet, ESMF_FieldGet, ESMF_MeshGet
1308 : use ESMF , only : ESMF_FIELDSTATUS_COMPLETE, ESMF_FAILURE
1309 :
1310 : ! input/output variables
1311 : type(ESMF_State) , intent(in) :: State
1312 : character(len=*) , intent(in) :: fldname
1313 : real(R8), optional, pointer :: fldptr(:)
1314 : real(R8), optional, pointer :: fldptr2d(:,:)
1315 : logical , optional, intent(out) :: exists
1316 : integer , intent(out) :: rc
1317 :
1318 : ! local variables
1319 : type(ESMF_FieldStatus_Flag) :: status
1320 : type(ESMF_StateItem_Flag) :: itemFlag
1321 : type(ESMF_Field) :: lfield
1322 : type(ESMF_Mesh) :: lmesh
1323 : integer :: nnodes, nelements
1324 : logical :: lexists
1325 : character(len=*), parameter :: subname='(atm_import_export:state_getfldptr)'
1326 : ! ----------------------------------------------
1327 :
1328 1561344 : rc = ESMF_SUCCESS
1329 :
1330 1561344 : lexists = .true.
1331 :
1332 : ! Determine if field with name fldname exists in state
1333 1561344 : if (present(exists)) then
1334 470016 : call ESMF_StateGet(state, trim(fldname), itemFlag, rc=rc)
1335 470016 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1336 470016 : if (itemflag == ESMF_STATEITEM_NOTFOUND) then
1337 : lexists = .false.
1338 : end if
1339 470016 : exists = lexists
1340 : end if
1341 :
1342 : if (lexists) then
1343 1423104 : call ESMF_StateGet(State, itemName=trim(fldname), field=lfield, rc=rc)
1344 1423104 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1345 1423104 : if (present(fldptr)) then
1346 1227264 : call ESMF_FieldGet(lfield, farrayPtr=fldptr, rc=rc)
1347 1227264 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1348 195840 : else if (present(fldptr2d)) then
1349 195840 : call ESMF_FieldGet(lfield, farrayPtr=fldptr2d, rc=rc)
1350 195840 : if (ChkErr(rc,__LINE__,u_FILE_u)) return
1351 : end if
1352 : end if
1353 :
1354 3122688 : end subroutine state_getfldptr
1355 :
1356 1561344 : end module atm_import_export
|