11! $Id$
22
3+ #define _DEALLOC(A) \
4+ if (associated (A))then ; \
5+ if (MAPL_ShmInitialized)then ; \
6+ call MAPL_SyncSharedMemory(rc= STATUS); \
7+ call MAPL_DeAllocNodeArray(A,rc= STATUS); \
8+ else ; \
9+ deallocate (A,stat= STATUS); \
10+ endif ; \
11+ _VERIFY(STATUS); \
12+ NULLIFY(A); \
13+ endif
14+
315#include " MAPL_Generic.h"
416
517
@@ -25,7 +37,13 @@ module GEOS_LakeGridCompMod
2537 integer , parameter :: NUM_SUBTILES = 2
2638
2739 type lake_state
28- integer :: CHOOSEMOSFC
40+ private
41+ integer :: CHOOSEMOSFC
42+ logical :: InitDone
43+ logical , pointer :: mask(:) = > null ()
44+ real :: tol_frice
45+ character (len= ESMF_MAXSTR) :: sstfile
46+ character (len= ESMF_MAXSTR) :: DataFrtFile
2947 end type lake_state
3048
3149 type lake_state_wrap
@@ -854,6 +872,7 @@ subroutine SetServices ( GC, RC )
854872
855873 allocate (mystate,stat= status)
856874 VERIFY_(status)
875+ mystate% InitDone = .false.
857876 call MAPL_GetResource (MAPL, SURFRC, label = ' SURFRC:' , default = ' GEOS_SurfaceGridComp.rc' , RC= STATUS) ; VERIFY_(STATUS)
858877 SCF = ESMF_ConfigCreate(rc= status) ; VERIFY_(STATUS)
859878 call ESMF_ConfigLoadFile (SCF,SURFRC,rc= status) ; VERIFY_(STATUS)
@@ -914,8 +933,6 @@ subroutine RUN1 ( GC, IMPORT, EXPORT, CLOCK, RC )
914933
915934 type (MAPL_MetaComp), pointer :: MAPL
916935 type (ESMF_State ) :: INTERNAL
917- type (ESMF_Alarm ) :: ALARM
918- type (ESMF_Config ) :: CF
919936
920937! pointers to export
921938
@@ -1553,6 +1570,13 @@ subroutine LAKECORE(NT,RC)
15531570 real , parameter :: LAKECAP = (MAPL_RHOWTR* MAPL_CAPWTR* LAKEDEPTH )
15541571 real , parameter :: LAKEICECAP = (MAPL_RHOWTR* MAPL_CAPWTR* LAKEICEDEPTH )
15551572
1573+ logical :: datalake
1574+ integer :: dlk
1575+ real , allocatable :: DATA_SST(:), DATA_FR(:)
1576+ character (len= ESMF_MAXSTR) :: maskfile
1577+ real , parameter :: Tfreeze = MAPL_TICE - 1.8
1578+ type (lake_state_wrap) :: wrap
1579+ type (lake_state), pointer :: mystate
15561580
15571581! Begin...
15581582!- ---------
@@ -1713,6 +1737,80 @@ subroutine LAKECORE(NT,RC)
17131737 if (associated (QST )) QST = 0.0
17141738 if (associated (HLWUP )) HLWUP = ALW
17151739 if (associated (LWNDSRF)) LWNDSRF = LWDNSRF - ALW
1740+ ! datalake addition
1741+ ! ==================
1742+
1743+ ! check resource if datalake is enabled
1744+ call MAPL_GetResource ( MAPL, DLK, Label= " DATALAKE:" , DEFAULT= 0 , RC= STATUS)
1745+ VERIFY_(STATUS)
1746+ DATALAKE = (DLK /= 0 )
1747+
1748+ call ESMF_UserCompGetInternalState(gc,' lake_private' ,wrap,status)
1749+ VERIFY_(status)
1750+ mystate = > wrap% ptr
1751+
1752+ ! Initialize a mask where we could apply the SST/FR from Reynolds/Ostia
1753+ if (.not. associated (mystate% mask)) then
1754+ allocate (mystate% mask(NT), stat= status); VERIFY_(STATUS)
1755+ mystate% mask = .false.
1756+ end if
1757+
1758+ if (datalake) then
1759+
1760+ ! next section is done only once. We do it here since the Initalize
1761+ ! method of this component defaults to MAPL_GenericInitialize
1762+
1763+ if (.not. mystate% InitDone) then
1764+ mystate% InitDone = .true.
1765+ call MAPL_GetResource ( MAPL, mystate% sstfile, &
1766+ Label= " LAKE_SST_FILE:" , DEFAULT= " sst.data" , RC= STATUS)
1767+ VERIFY_(STATUS)
1768+ call MAPL_GetResource ( MAPL, mystate% dataFrtFile, &
1769+ Label= " LAKE_FRT_FILE:" , DEFAULT= " fraci.data" , RC= STATUS)
1770+ VERIFY_(STATUS)
1771+ call MAPL_GetResource ( MAPL, mystate% tol_frice, &
1772+ Label= " LAKE_TOL_FRICE:" , DEFAULT= 1.0e-2 , RC= STATUS)
1773+ VERIFY_(STATUS)
1774+
1775+ call MAPL_GetResource ( MAPL, MASKFILE, &
1776+ Label= " DATALAKE_MASK_FILE:" , DEFAULT= " DataLakeMask.data" , RC= STATUS)
1777+ VERIFY_(STATUS)
1778+
1779+ ! ALT: For now we are reading the entire mask from a file
1780+ ! we might decide later to use a different strategy (for example
1781+ ! Santha suggested lat-lon based description)
1782+
1783+ call DataLakeReadMask(MAPL, mystate% mask, &
1784+ maskfile, RC= status)
1785+ VERIFY_(STATUS)
1786+ end if
1787+
1788+
1789+ allocate (DATA_SST(NT), DATA_FR(NT), stat= status); VERIFY_(STATUS)
1790+ call MAPL_ReadForcing(MAPL, ' SST' , mystate% sstfile,&
1791+ CURRENT_TIME, DATA_SST, RC= STATUS)
1792+ VERIFY_(STATUS)
1793+ call MAPL_ReadForcing(MAPL, ' FR' , mystate% dataFrtFile,&
1794+ CURRENT_TIME, DATA_FR, RC= STATUS)
1795+ VERIFY_(STATUS)
1796+
1797+ do I= 1 ,NT
1798+ if (mystate% mask(i)) then
1799+ ! we are operating over the 'observed' lake: Great Lakes and Caspian Sea (set by mask)
1800+ TS(I,WATER) = DATA_SST(I)
1801+ TS(I,ICE) = Tfreeze
1802+ if (data_fr(i) > mystate% tol_frice) then ! have lake ice
1803+ FR(I,WATER) = 1.0 - DATA_FR(I)
1804+ FR(I,ICE) = DATA_FR(I)
1805+ else
1806+ FR(I,WATER) = 1.0
1807+ FR(I,ICE) = 0.0
1808+ end if
1809+ end if
1810+ end do
1811+
1812+ deallocate (DATA_SST, DATA_FR)
1813+ end if
17161814
17171815 do N= 1 ,NUM_SUBTILES
17181816 CFT = (CH(:,N)/ CTATM)
@@ -1749,7 +1847,7 @@ subroutine LAKECORE(NT,RC)
17491847! Update surface temperature and moisture
17501848!- ---------------------------------------
17511849
1752- TS(:,N) = TS(:,N) + DTS
1850+ where ( .not. mystate % mask) TS(:,N) = TS(:,N) + DTS
17531851 DQS = GEOS_QSAT(TS(:,N), PS, RAMP= 0.0 , PASCALS= .TRUE. ) - QS(:,N)
17541852 QS(:,N) = QS(:,N) + DQS
17551853
@@ -1773,6 +1871,10 @@ subroutine LAKECORE(NT,RC)
17731871! Update Ice fraction
17741872!- -------------------
17751873 do I= 1 ,NT
1874+ ! $ if(mystate%mask(i)) cycle
1875+ if (mystate% mask(i)) then
1876+ cycle
1877+ end if
17761878 if (TS(I,ICE)>MAPL_TICE .and. FR(I,ICE)>0.0 ) then
17771879 ! MELT
17781880 FR(I,WATER) = 1.0
@@ -1919,6 +2021,54 @@ end subroutine ALBLAKEICE
19192021
19202022end subroutine RUN2
19212023
2024+ subroutine DataLakeReadMask (mapl , msk , maskfile , rc )
2025+ implicit none
2026+ ! arguments
2027+ type (MAPL_MetaComp), pointer :: MAPL
2028+ logical , intent (INOUT ) :: msk(:)
2029+ character (len=* ), intent (IN ) :: maskfile
2030+ integer , optional , intent (OUT ) :: rc
2031+
2032+ ! errlog vars
2033+ integer :: status
2034+ character (len= ESMF_MAXSTR), parameter :: Iam= ' DataLakeReadMask'
2035+
2036+ ! local vars
2037+ integer :: unit
2038+ integer , pointer :: tilemask(:) = > null ()
2039+ type (ESMF_Grid) :: TILEGRID
2040+ type (MAPL_LocStream) :: LOCSTREAM
2041+ integer :: NT
2042+ real , allocatable :: imask(:)
2043+
2044+ ! this is the first attempt to read the mask. For now we support binary
2045+ call MAPL_Get(MAPL, LocStream= LOCSTREAM, RC= STATUS)
2046+ VERIFY_(STATUS)
2047+ call MAPL_LocStreamGet(LOCSTREAM, TILEGRID= TILEGRID, RC= STATUS)
2048+ VERIFY_(STATUS)
2049+
2050+ call MAPL_TileMaskGet(tilegrid, tilemask, rc= status)
2051+ VERIFY_(STATUS)
2052+
2053+ nt = size (msk)
2054+ allocate (imask(nt), stat= status)
2055+
2056+ unit = GETFILE( maskfile, form= " unformatted" , RC= STATUS )
2057+ VERIFY_(STATUS)
2058+
2059+ call MAPL_VarRead(unit, tilegrid, imask, mask= tilemask, rc= status)
2060+ VERIFY_(STATUS)
2061+
2062+ call FREE_FILE(unit, RC= STATUS)
2063+ VERIFY_(STATUS)
2064+
2065+ msk = (imask /= 0.0 )
2066+ deallocate (imask)
2067+ _DEALLOC(tilemask)
2068+
2069+ RETURN_(ESMF_SUCCESS)
2070+ end subroutine DataLakeReadMask
2071+
19222072! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
19232073
19242074end module GEOS_LakeGridCompMod
0 commit comments