clm5.0/src_clm40/main/lnd_comp_mct.F90
2025-01-12 20:48:10 +08:00

650 lines
26 KiB
Fortran

module lnd_comp_mct
!---------------------------------------------------------------------------
! !DESCRIPTION:
! Interface of the active land model component of CESM the CLM (Community Land Model)
! with the main CESM driver. This is a thin interface taking CESM driver information
! in MCT (Model Coupling Toolkit) format and converting it to use by CLM.
!
! !uses:
use shr_kind_mod , only : r8 => shr_kind_r8
use shr_sys_mod , only : shr_sys_flush
use mct_mod , only : mct_avect, mct_gsmap, mct_gGrid
use decompmod , only : bounds_type, ldecomp
use lnd_import_export
!
! !public member functions:
implicit none
save
private ! by default make data private
!
! !public member functions:
public :: lnd_init_mct ! clm initialization
public :: lnd_run_mct ! clm run phase
public :: lnd_final_mct ! clm finalization/cleanup
!
! !private member functions:
private :: lnd_setgsmap_mct ! set the land model mct gs map
private :: lnd_domain_mct ! set the land model domain information
!---------------------------------------------------------------------------
contains
!====================================================================================
subroutine lnd_init_mct( EClock, cdata_l, x2l_l, l2x_l, NLFilename )
!
! !DESCRIPTION:
! Initialize land surface model and obtain relevant atmospheric model arrays
! back from (i.e. albedos, surface temperature and snow cover over land).
!
! !USES:
use shr_kind_mod , only : shr_kind_cl
use abortutils , only : endrun
use clm_time_manager , only : get_nstep, get_step_size, set_timemgr_init, &
set_nextsw_cday
use clm_atmlnd , only : clm_l2a
use clm_glclnd , only : clm_s2x
use clm_initializeMod, only : initialize1, initialize2
use clm_varctl , only : finidat,single_column, clm_varctl_set, iulog, noland, &
inst_index, inst_suffix, inst_name
use clm_varorb , only : eccen, obliqr, lambm0, mvelpp
use controlMod , only : control_setNL
use decompMod , only : get_proc_bounds
use domainMod , only : ldomain
use shr_file_mod , only : shr_file_setLogUnit, shr_file_setLogLevel, &
shr_file_getLogUnit, shr_file_getLogLevel, &
shr_file_getUnit, shr_file_setIO
use seq_cdata_mod , only : seq_cdata, seq_cdata_setptrs
use seq_timemgr_mod , only : seq_timemgr_EClockGetData
use seq_infodata_mod , only : seq_infodata_type, seq_infodata_GetData, seq_infodata_PutData, &
seq_infodata_start_type_start, seq_infodata_start_type_cont, &
seq_infodata_start_type_brnch
use seq_comm_mct , only : seq_comm_suffix, seq_comm_inst, seq_comm_name
use seq_flds_mod , only : seq_flds_x2l_fields, seq_flds_l2x_fields
use spmdMod , only : masterproc, spmd_init
use clm_varctl , only : nsrStartup, nsrContinue, nsrBranch
use clm_cpl_indices , only : clm_cpl_indices_set
use mct_mod , only : mct_aVect_init, mct_aVect_zero, mct_gsMap_lsize
use ESMF
!
! !ARGUMENTS:
type(ESMF_Clock), intent(inout) :: EClock ! Input synchronization clock
type(seq_cdata), intent(inout) :: cdata_l ! Input land-model driver data
type(mct_aVect), intent(inout) :: x2l_l, l2x_l ! land model import and export states
character(len=*), optional, intent(in) :: NLFilename ! Namelist filename to read
!
! !LOCAL VARIABLES:
integer :: LNDID ! Land identifyer
integer :: mpicom_lnd ! MPI communicator
type(mct_gsMap), pointer :: GSMap_lnd ! Land model MCT GS map
type(mct_gGrid), pointer :: dom_l ! Land model domain
type(seq_infodata_type), pointer :: infodata ! CESM driver level info data
integer :: lsize ! size of attribute vector
integer :: g,i,j ! indices
integer :: dtime_sync ! coupling time-step from the input synchronization clock
integer :: dtime_clm ! clm time-step
logical :: exists ! true if file exists
logical :: atm_aero ! Flag if aerosol data sent from atm model
real(r8) :: scmlat ! single-column latitude
real(r8) :: scmlon ! single-column longitude
real(r8) :: nextsw_cday ! calday from clock of next radiation computation
character(len=SHR_KIND_CL) :: caseid ! case identifier name
character(len=SHR_KIND_CL) :: ctitle ! case description title
character(len=SHR_KIND_CL) :: starttype ! start-type (startup, continue, branch, hybrid)
character(len=SHR_KIND_CL) :: calendar ! calendar type name
character(len=SHR_KIND_CL) :: hostname ! hostname of machine running on
character(len=SHR_KIND_CL) :: version ! Model version
character(len=SHR_KIND_CL) :: username ! user running the model
integer :: nsrest ! clm restart type
integer :: ref_ymd ! reference date (YYYYMMDD)
integer :: ref_tod ! reference time of day (sec)
integer :: start_ymd ! start date (YYYYMMDD)
integer :: start_tod ! start time of day (sec)
integer :: stop_ymd ! stop date (YYYYMMDD)
integer :: stop_tod ! stop time of day (sec)
logical :: brnch_retain_casename ! flag if should retain the case name on a branch start type
integer :: lbnum ! input to memory diagnostic
integer :: shrlogunit,shrloglev ! old values for log unit and log level
type(bounds_type) :: bounds ! bounds
character(len=32), parameter :: sub = 'lnd_init_mct'
character(len=*), parameter :: format = "('("//trim(sub)//") :',A)"
!-----------------------------------------------------------------------
! Set cdata data
call seq_cdata_setptrs(cdata_l, ID=LNDID, mpicom=mpicom_lnd, &
gsMap=GSMap_lnd, dom=dom_l, infodata=infodata)
! Determine attriute vector indices
call clm_cpl_indices_set()
! Initialize clm MPI communicator
call spmd_init( mpicom_lnd, LNDID )
#if (defined _MEMTRACE)
if(masterproc) then
lbnum=1
call memmon_dump_fort('memmon.out','lnd_init_mct:start::',lbnum)
endif
#endif
inst_name = seq_comm_name(LNDID)
inst_index = seq_comm_inst(LNDID)
inst_suffix = seq_comm_suffix(LNDID)
! Initialize io log unit
call shr_file_getLogUnit (shrlogunit)
if (masterproc) then
inquire(file='lnd_modelio.nml'//trim(inst_suffix),exist=exists)
if (exists) then
iulog = shr_file_getUnit()
call shr_file_setIO('lnd_modelio.nml'//trim(inst_suffix),iulog)
end if
write(iulog,format) "CLM land model initialization"
else
iulog = shrlogunit
end if
call shr_file_getLogLevel(shrloglev)
call shr_file_setLogUnit (iulog)
! Use infodata to set orbital values
call seq_infodata_GetData( infodata, orb_eccen=eccen, orb_mvelpp=mvelpp, &
orb_lambm0=lambm0, orb_obliqr=obliqr )
! Consistency check on namelist filename
call control_setNL("lnd_in"//trim(inst_suffix))
! Initialize clm
! initialize1 reads namelist, grid and surface data (need this to initialize gsmap)
! initialize2 performs rest of initialization
call seq_timemgr_EClockGetData(EClock, &
start_ymd=start_ymd, &
start_tod=start_tod, ref_ymd=ref_ymd, &
ref_tod=ref_tod, stop_ymd=stop_ymd, &
stop_tod=stop_tod, &
calendar=calendar )
call seq_infodata_GetData(infodata, case_name=caseid, &
case_desc=ctitle, single_column=single_column, &
scmlat=scmlat, scmlon=scmlon, &
brnch_retain_casename=brnch_retain_casename, &
start_type=starttype, model_version=version, &
hostname=hostname, username=username )
call set_timemgr_init( calendar_in=calendar, start_ymd_in=start_ymd, start_tod_in=start_tod, &
ref_ymd_in=ref_ymd, ref_tod_in=ref_tod, stop_ymd_in=stop_ymd, &
stop_tod_in=stop_tod)
if ( trim(starttype) == trim(seq_infodata_start_type_start)) then
nsrest = nsrStartup
else if (trim(starttype) == trim(seq_infodata_start_type_cont) ) then
nsrest = nsrContinue
else if (trim(starttype) == trim(seq_infodata_start_type_brnch)) then
nsrest = nsrBranch
else
call endrun( sub//' ERROR: unknown starttype' )
end if
call clm_varctl_set(caseid_in=caseid, ctitle_in=ctitle, &
brnch_retain_casename_in=brnch_retain_casename, &
single_column_in=single_column, scmlat_in=scmlat, &
scmlon_in=scmlon, nsrest_in=nsrest, version_in=version, &
hostname_in=hostname, username_in=username)
! Read namelist, grid and surface data
call initialize1( )
! If no land then exit out of initialization
if ( noland ) then
call seq_infodata_PutData( infodata, lnd_present =.false.)
call seq_infodata_PutData( infodata, lnd_prognostic=.false.)
return
end if
! Determine if aerosol and dust deposition come from atmosphere component
call seq_infodata_GetData(infodata, atm_aero=atm_aero )
if ( .not. atm_aero )then
call endrun( sub//' ERROR: atmosphere model MUST send aerosols to CLM' )
end if
! Initialize clm gsMap, clm domain and clm attribute vectors
call get_proc_bounds( bounds )
call lnd_SetgsMap_mct( bounds, mpicom_lnd, LNDID, gsMap_lnd )
lsize = mct_gsMap_lsize(gsMap_lnd, mpicom_lnd)
call lnd_domain_mct( bounds, lsize, gsMap_lnd, dom_l )
call mct_aVect_init(x2l_l, rList=seq_flds_x2l_fields, lsize=lsize)
call mct_aVect_zero(x2l_l)
call mct_aVect_init(l2x_l, rList=seq_flds_l2x_fields, lsize=lsize)
call mct_aVect_zero(l2x_l)
! Finish initializing clm
call initialize2()
! Check that clm internal dtime aligns with clm coupling interval
call seq_timemgr_EClockGetData(EClock, dtime=dtime_sync )
dtime_clm = get_step_size()
if (masterproc) then
write(iulog,*)'dtime_sync= ',dtime_sync,&
' dtime_clm= ',dtime_clm,' mod = ',mod(dtime_sync,dtime_clm)
end if
if (mod(dtime_sync,dtime_clm) /= 0) then
write(iulog,*)'clm dtime ',dtime_clm,' and Eclock dtime ',&
dtime_sync,' never align'
call endrun( sub//' ERROR: time out of sync' )
end if
! Create land export state
call lnd_export(bounds, clm_l2a, clm_s2x, l2x_l%rattr)
! Fill in infodata settings
call seq_infodata_PutData(infodata, lnd_prognostic=.true.)
call seq_infodata_PutData(infodata, lnd_nx=ldomain%ni, lnd_ny=ldomain%nj)
! Get infodata info
call seq_infodata_GetData(infodata, nextsw_cday=nextsw_cday )
call set_nextsw_cday(nextsw_cday)
! Reset shr logging to original values
call shr_file_setLogUnit (shrlogunit)
call shr_file_setLogLevel(shrloglev)
#if (defined _MEMTRACE)
if(masterproc) then
write(iulog,*) TRIM(Sub) // ':end::'
lbnum=1
call memmon_dump_fort('memmon.out','lnd_int_mct:end::',lbnum)
call memmon_reset_addr()
endif
#endif
end subroutine lnd_init_mct
!====================================================================================
subroutine lnd_run_mct(EClock, cdata_l, x2l_l, l2x_l)
!
! !DESCRIPTION:
! Run clm model
!
! !USES:
use shr_kind_mod ,only : r8 => shr_kind_r8
use clmtype
use clm_atmlnd ,only : clm_l2a, clm_a2l, a2l_not_downscaled_gcell
use clm_glclnd ,only : clm_s2x, clm_x2s
use clm_driver ,only : clm_drv
use clm_time_manager,only : get_curr_date, get_nstep, get_curr_calday, get_step_size, &
advance_timestep, set_nextsw_cday,update_rad_dtime
use decompMod ,only : get_proc_bounds
use abortutils ,only : endrun
use clm_varctl ,only : iulog
use clm_varorb ,only : eccen, obliqr, lambm0, mvelpp
use shr_file_mod ,only : shr_file_setLogUnit, shr_file_setLogLevel, &
shr_file_getLogUnit, shr_file_getLogLevel
use seq_cdata_mod ,only : seq_cdata, seq_cdata_setptrs
use seq_timemgr_mod ,only : seq_timemgr_EClockGetData, seq_timemgr_StopAlarmIsOn, &
seq_timemgr_RestartAlarmIsOn, seq_timemgr_EClockDateInSync
use seq_infodata_mod,only : seq_infodata_type, seq_infodata_GetData
use spmdMod ,only : masterproc, mpicom
use perf_mod ,only : t_startf, t_stopf, t_barrierf
use shr_orb_mod ,only : shr_orb_decl
use ESMF
!
! !ARGUMENTS:
type(ESMF_Clock) , intent(inout) :: EClock ! Input synchronization clock from driver
type(seq_cdata) , intent(inout) :: cdata_l ! Input driver data for land model
type(mct_aVect) , intent(inout) :: x2l_l ! Import state to land model
type(mct_aVect) , intent(inout) :: l2x_l ! Export state from land model
!
! !LOCAL VARIABLES:
integer :: ymd_sync ! Sync date (YYYYMMDD)
integer :: yr_sync ! Sync current year
integer :: mon_sync ! Sync current month
integer :: day_sync ! Sync current day
integer :: tod_sync ! Sync current time of day (sec)
integer :: ymd ! CLM current date (YYYYMMDD)
integer :: yr ! CLM current year
integer :: mon ! CLM current month
integer :: day ! CLM current day
integer :: tod ! CLM current time of day (sec)
integer :: dtime ! time step increment (sec)
integer :: nstep ! time step index
logical :: rstwr_sync ! .true. ==> write restart file before returning
logical :: rstwr ! .true. ==> write restart file before returning
logical :: nlend_sync ! Flag signaling last time-step
logical :: nlend ! .true. ==> last time-step
logical :: dosend ! true => send data back to driver
logical :: doalb ! .true. ==> do albedo calculation on this time step
real(r8) :: nextsw_cday ! calday from clock of next radiation computation
real(r8) :: caldayp1 ! clm calday plus dtime offset
integer :: shrlogunit,shrloglev ! old values for share log unit and log level
integer :: lbnum ! input to memory diagnostic
integer :: g,i,lsize ! counters
real(r8) :: calday ! calendar day for nstep
real(r8) :: declin ! solar declination angle in radians for nstep
real(r8) :: declinp1 ! solar declination angle in radians for nstep+1
real(r8) :: eccf ! earth orbit eccentricity factor
real(r8) :: recip ! reciprical
logical,save :: first_call = .true. ! first call work
type(seq_infodata_type),pointer :: infodata ! CESM information from the driver
type(mct_gGrid), pointer :: dom_l ! Land model domain data
type(bounds_type) :: bounds ! bounds
character(len=32) :: rdate ! date char string for restart file names
character(len=32), parameter :: sub = "lnd_run_mct"
!---------------------------------------------------------------------------
! Determine processor bounds
call get_proc_bounds(bounds)
#if (defined _MEMTRACE)
if(masterproc) then
lbnum=1
call memmon_dump_fort('memmon.out','lnd_run_mct:start::',lbnum)
endif
#endif
! Reset shr logging to my log file
call shr_file_getLogUnit (shrlogunit)
call shr_file_getLogLevel(shrloglev)
call shr_file_setLogUnit (iulog)
! Determine time of next atmospheric shortwave calculation
call seq_cdata_setptrs(cdata_l, infodata=infodata, dom=dom_l)
call seq_timemgr_EClockGetData(EClock, &
curr_ymd=ymd, curr_tod=tod_sync, &
curr_yr=yr_sync, curr_mon=mon_sync, curr_day=day_sync)
call seq_infodata_GetData(infodata, nextsw_cday=nextsw_cday )
call set_nextsw_cday( nextsw_cday )
dtime = get_step_size()
write(rdate,'(i4.4,"-",i2.2,"-",i2.2,"-",i5.5)') yr_sync,mon_sync,day_sync,tod_sync
nlend_sync = seq_timemgr_StopAlarmIsOn( EClock )
rstwr_sync = seq_timemgr_RestartAlarmIsOn( EClock )
! Map MCT to land data type
! Perform downscaling if appropriate
! Map to clm (only when state and/or fluxes need to be updated)
call t_startf ('lc_lnd_import')
call lnd_import( bounds, x2l_l%rattr, clm_a2l, a2l_not_downscaled_gcell, clm_x2s)
call t_stopf ('lc_lnd_import')
! Use infodata to set orbital values if updated mid-run
call seq_infodata_GetData( infodata, orb_eccen=eccen, orb_mvelpp=mvelpp, &
orb_lambm0=lambm0, orb_obliqr=obliqr )
! Loop over time steps in coupling interval
dosend = .false.
do while(.not. dosend)
! Determine if dosend
! When time is not updated at the beginning of the loop - then return only if
! are in sync with clock before time is updated
call get_curr_date( yr, mon, day, tod )
ymd = yr*10000 + mon*100 + day
tod = tod
dosend = (seq_timemgr_EClockDateInSync( EClock, ymd, tod))
! Determine doalb based on nextsw_cday sent from atm model
nstep = get_nstep()
caldayp1 = get_curr_calday(offset=dtime)
if (nstep == 0) then
doalb = .false.
else if (nstep == 1) then
doalb = (abs(nextsw_cday- caldayp1) < 1.e-10_r8)
else
doalb = (nextsw_cday >= -0.5_r8)
end if
call update_rad_dtime(doalb)
! Determine if time to write cam restart and stop
rstwr = .false.
if (rstwr_sync .and. dosend) rstwr = .true.
nlend = .false.
if (nlend_sync .and. dosend) nlend = .true.
! Run clm
call t_barrierf('sync_clm_run1', mpicom)
call t_startf ('clm_run')
call t_startf ('shr_orb_decl')
calday = get_curr_calday()
call shr_orb_decl( calday , eccen, mvelpp, lambm0, obliqr, declin , eccf )
call shr_orb_decl( nextsw_cday, eccen, mvelpp, lambm0, obliqr, declinp1, eccf )
call t_stopf ('shr_orb_decl')
call clm_drv(doalb, nextsw_cday, declinp1, declin, rstwr, nlend, rdate)
call t_stopf ('clm_run')
! Create l2x_l export state - add river runoff input to l2x_l if appropriate
call t_startf ('lc_lnd_export')
call lnd_export(bounds, clm_l2a, clm_s2x, l2x_l%rattr)
call t_stopf ('lc_lnd_export')
! Advance clm time step
call t_startf ('lc_clm2_adv_timestep')
call advance_timestep()
call t_stopf ('lc_clm2_adv_timestep')
end do
! Check that internal clock is in sync with master clock
call get_curr_date( yr, mon, day, tod, offset=-dtime )
ymd = yr*10000 + mon*100 + day
tod = tod
if ( .not. seq_timemgr_EClockDateInSync( EClock, ymd, tod ) )then
call seq_timemgr_EclockGetData( EClock, curr_ymd=ymd_sync, curr_tod=tod_sync )
write(iulog,*)' clm ymd=',ymd ,' clm tod= ',tod
write(iulog,*)'sync ymd=',ymd_sync,' sync tod= ',tod_sync
call endrun( sub//":: CLM clock not in sync with Master Sync clock" )
end if
! Reset shr logging to my original values
call shr_file_setLogUnit (shrlogunit)
call shr_file_setLogLevel(shrloglev)
#if (defined _MEMTRACE)
if(masterproc) then
lbnum=1
call memmon_dump_fort('memmon.out','lnd_run_mct:end::',lbnum)
call memmon_reset_addr()
endif
#endif
first_call = .false.
end subroutine lnd_run_mct
!====================================================================================
subroutine lnd_final_mct( EClock, cdata_l, x2l_l, l2x_l)
!
! !DESCRIPTION:
! Finalize land surface model
use seq_cdata_mod ,only : seq_cdata, seq_cdata_setptrs
use seq_timemgr_mod ,only : seq_timemgr_EClockGetData, seq_timemgr_StopAlarmIsOn, &
seq_timemgr_RestartAlarmIsOn, seq_timemgr_EClockDateInSync
use mct_mod
use esmf
!
! !ARGUMENTS:
type(ESMF_Clock) , intent(inout) :: EClock ! Input synchronization clock from driver
type(seq_cdata) , intent(inout) :: cdata_l ! Input driver data for land model
type(mct_aVect) , intent(inout) :: x2l_l ! Import state to land model
type(mct_aVect) , intent(inout) :: l2x_l ! Export state from land model
!---------------------------------------------------------------------------
! fill this in
end subroutine lnd_final_mct
!====================================================================================
subroutine lnd_setgsmap_mct( bounds, mpicom_lnd, LNDID, gsMap_lnd )
!
! !DESCRIPTION:
! Set the MCT GS map for the land model
!
! !USES:
use shr_kind_mod , only : r8 => shr_kind_r8
use domainMod , only : ldomain
use mct_mod , only : mct_gsMap, mct_gsMap_init
implicit none
!
! !ARGUMENTS:
type(bounds_type) , intent(in) :: bounds ! bounds
integer , intent(in) :: mpicom_lnd ! MPI communicator for the clm land model
integer , intent(in) :: LNDID ! Land model identifyer number
type(mct_gsMap) , intent(out) :: gsMap_lnd ! Resulting MCT GS map for the land model
!
! !LOCAL VARIABLES:
integer,allocatable :: gindex(:) ! Number the local grid points
integer :: i, j, n, gi ! Indices
integer :: lsize,gsize ! GS Map size
integer :: ier ! Error code
!---------------------------------------------------------------------------
! Build the land grid numbering for MCT
! NOTE: Numbering scheme is: West to East and South to North
! starting at south pole. Should be the same as what's used in SCRIP
allocate(gindex(bounds%begg:bounds%endg),stat=ier)
! number the local grid
do n = bounds%begg, bounds%endg
gindex(n) = ldecomp%gdc2glo(n)
end do
lsize = bounds%endg - bounds%begg + 1
gsize = ldomain%ni * ldomain%nj
call mct_gsMap_init( gsMap_lnd, gindex, mpicom_lnd, LNDID, lsize, gsize )
deallocate(gindex)
end subroutine lnd_SetgsMap_mct
!====================================================================================
subroutine lnd_domain_mct( bounds, lsize, gsMap_l, dom_l )
!
! !DESCRIPTION:
! Send the land model domain information to the coupler
!
! !USES:
use clm_varcon , only: re
use domainMod , only: ldomain
use spmdMod , only: iam
use mct_mod , only: mct_gGrid_importIAttr
use mct_mod , only: mct_gGrid_importRAttr, mct_gGrid_init, mct_gsMap_orderedPoints
use seq_flds_mod, only: seq_flds_dom_coord, seq_flds_dom_other
!
! !ARGUMENTS:
type(bounds_type), intent(in) :: bounds ! bounds
integer , intent(in) :: lsize ! land model domain data size
type(mct_gsMap), intent(inout) :: gsMap_l ! Output land model MCT GS map
type(mct_ggrid), intent(out) :: dom_l ! Output domain information for land model
!
! Local Variables
integer :: g,i,j ! index
real(r8), pointer :: data(:) ! temporary
integer , pointer :: idata(:) ! temporary
!---------------------------------------------------------------------------
!
! Initialize mct domain type
! lat/lon in degrees, area in radians^2, mask is 1 (land), 0 (non-land)
! Note that in addition land carries around landfrac for the purposes of domain checking
!
call mct_gGrid_init( GGrid=dom_l, CoordChars=trim(seq_flds_dom_coord), &
OtherChars=trim(seq_flds_dom_other), lsize=lsize )
!
! Allocate memory
!
allocate(data(lsize))
!
! Determine global gridpoint number attribute, GlobGridNum, which is set automatically by MCT
!
call mct_gsMap_orderedPoints(gsMap_l, iam, idata)
call mct_gGrid_importIAttr(dom_l,'GlobGridNum',idata,lsize)
!
! Determine domain (numbering scheme is: West to East and South to North to South pole)
! Initialize attribute vector with special value
!
data(:) = -9999.0_R8
call mct_gGrid_importRAttr(dom_l,"lat" ,data,lsize)
call mct_gGrid_importRAttr(dom_l,"lon" ,data,lsize)
call mct_gGrid_importRAttr(dom_l,"area" ,data,lsize)
call mct_gGrid_importRAttr(dom_l,"aream",data,lsize)
data(:) = 0.0_R8
call mct_gGrid_importRAttr(dom_l,"mask" ,data,lsize)
!
! Fill in correct values for domain components
! Note aream will be filled in in the atm-lnd mapper
!
do g = bounds%begg,bounds%endg
i = 1 + (g - bounds%begg)
data(i) = ldomain%lonc(g)
end do
call mct_gGrid_importRattr(dom_l,"lon",data,lsize)
do g = bounds%begg,bounds%endg
i = 1 + (g - bounds%begg)
data(i) = ldomain%latc(g)
end do
call mct_gGrid_importRattr(dom_l,"lat",data,lsize)
do g = bounds%begg,bounds%endg
i = 1 + (g - bounds%begg)
data(i) = ldomain%area(g)/(re*re)
end do
call mct_gGrid_importRattr(dom_l,"area",data,lsize)
do g = bounds%begg,bounds%endg
i = 1 + (g - bounds%begg)
data(i) = real(ldomain%mask(g), r8)
end do
call mct_gGrid_importRattr(dom_l,"mask",data,lsize)
do g = bounds%begg,bounds%endg
i = 1 + (g - bounds%begg)
data(i) = real(ldomain%frac(g), r8)
end do
call mct_gGrid_importRattr(dom_l,"frac",data,lsize)
deallocate(data)
deallocate(idata)
end subroutine lnd_domain_mct
end module lnd_comp_mct