LCOV - code coverage report
Current view: top level - cpl/nuopc - atm_import_export.F90 (source / functions) Hit Total Coverage
Test: coverage.info Lines: 563 688 81.8 %
Date: 2025-01-13 21:54:50 Functions: 9 10 90.0 %

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

Generated by: LCOV version 1.14