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