! $Id$ ! ! Earth System Modeling Framework ! Copyright (c) 2002-2023, University Corporation for Atmospheric Research, ! Massachusetts Institute of Technology, Geophysical Fluid Dynamics ! Laboratory, University of Michigan, National Centers for Environmental ! Prediction, Los Alamos National Laboratory, Argonne National Laboratory, ! NASA Goddard Space Flight Center. ! Licensed under the University of Illinois-NCSA License. ! !============================================================================== #define FILENAME "src/addon/NUOPC/src/NUOPC_ModelBase.F90" !============================================================================== module NUOPC_ModelBase !----------------------------------------------------------------------------- ! Generic Model Component Base implementation !----------------------------------------------------------------------------- use ESMF use NUOPC, SetVM => NUOPC_SetVM implicit none private public & SetVM, & SetServices, & routine_Run, & routine_Nop public & label_Advertise, & label_ModifyAdvertised, & label_RealizeProvided, & label_AcceptTransfer, & label_RealizeAccepted, & label_SetClock, & label_DataInitialize, & label_Advance, & label_AdvanceClock, & label_CheckImport, & label_SetRunClock, & label_TimestampExport, & label_Finalize, & type_InternalStateStruct, & type_InternalState, & label_InternalState character(*), parameter :: & label_InternalState = "ModelBase_InternalState" character(*), parameter :: & label_Advance = "ModelBase_Advance" character(*), parameter :: & label_AdvanceClock = "ModelBase_AdvanceClock" character(*), parameter :: & label_CheckImport = "ModelBase_CheckImport" character(*), parameter :: & label_SetRunClock = "ModelBase_SetRunClock" character(*), parameter :: & label_TimestampExport = "ModelBase_TimestampExport" character(*), parameter :: & label_Finalize = "ModelBase_Finalize" character(*), parameter :: & label_Advertise = "ModelBase_Advertise" character(*), parameter :: & label_ModifyAdvertised = "ModelBase_ModifyAdvertised" character(*), parameter :: & label_RealizeProvided = "ModelBase_RealizeProvided" character(*), parameter :: & label_AcceptTransfer = "ModelBase_AcceptTransfer" character(*), parameter :: & label_RealizeAccepted = "ModelBase_RealizeAccepted" character(*), parameter :: & label_SetClock = "ModelBase_SetClock" character(*), parameter :: & label_DataInitialize = "ModelBase_DataInitialize" type type_InternalStateStruct type(ESMF_Clock) :: driverClock type(ESMF_Time) :: preAdvanceCurrTime type(ESMF_Field), pointer :: cachedExportFieldList(:) end type type type_InternalState type(type_InternalStateStruct), pointer :: wrap end type ! Generic methods public NUOPC_ModelBaseGet !----------------------------------------------------------------------------- contains !----------------------------------------------------------------------------- subroutine SetServices(gcomp, rc) type(ESMF_GridComp) :: gcomp integer, intent(out) :: rc ! local variables character(ESMF_MAXSTR) :: name integer :: stat type(type_InternalState) :: is rc = ESMF_SUCCESS ! query the component for info call NUOPC_CompGet(gcomp, name=name, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! add standard NUOPC GridComp Attribute Package to the Model call NUOPC_CompAttributeInit(gcomp, kind="Model", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out #ifndef NO_COMP_SPECIFIC_COMPLIANCE_CHECK ! set the ESMF compliance checker register Attribute call ESMF_AttributeSet(gcomp, name="ESMF_RUNTIME_COMPLIANCEICREGISTER", & value="NUOPC_Model_ComplianceICR", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out #endif ! Initialize phases ! Phase 0 requires use of ESMF method. call ESMF_GridCompSetEntryPoint(gcomp, ESMF_METHOD_INITIALIZE, & userRoutine=InitializeP0, phase=0, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_INITIALIZE, & phaseLabelList=(/"IPDv00p1", "IPDv01p1", "IPDv02p1", "IPDv03p1", & "IPDv04p1", "IPDv05p1"/), userRoutine=routine_Nop, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_INITIALIZE, & phaseLabelList=(/"IPDv00p2", "IPDv01p3", "IPDv02p3", "IPDv03p3", & "IPDv04p3", "IPDv05p4"/), userRoutine=routine_Nop, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! data initialization call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_INITIALIZE, & phaseLabelList=(/"IPDv00p3", "IPDv01p4", "IPDv02p4", "IPDv03p6", & "IPDv04p6", "IPDv05p7"/), userRoutine=InitializeIPDvXp07, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_INITIALIZE, & phaseLabelList=(/"IPDv00p4", "IPDv01p5"/), & userRoutine=InitializeIPDv01p5, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_INITIALIZE, & phaseLabelList=(/"IPDv02p5", "IPDv03p7", "IPDv04p7", "IPDv05p8"/), & userRoutine=InitializeIPDvXp08, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! Run phases call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_RUN, & phaseLabelList=(/"RunPhase1"/), userRoutine=routine_Run, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! Specialize default Run -> setting the run Clock call NUOPC_CompSpecialize(gcomp, specLabel=label_SetRunClock, & specRoutine=SetRunClock, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! Specialize default Run -> checking import Fields call NUOPC_CompSpecialize(gcomp, specLabel=label_CheckImport, & specRoutine=CheckImport, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! Finalize phases call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_FINALIZE, & phaseLabelList=(/"FinalizePhase1"/), userRoutine=Finalize, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! Set up the internal state allocate(is%wrap, stat=stat) if (ESMF_LogFoundAllocError(statusToCheck=stat, & msg="Allocation of internal state memory failed.", & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out #ifdef ESMF_NO_F2018ASSUMEDTYPE call ESMF_UserCompSetInternalState(gcomp, label_InternalState, is, rc) #else call ESMF_UserCompSetInternalState(gcomp, label_InternalState, is, rc=rc) #endif if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out ! Set IPDvX attribute call NUOPC_CompAttributeSet(gcomp, name="IPDvX", value="true", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out end subroutine !----------------------------------------------------------------------------- subroutine routine_Nop(gcomp, importState, exportState, clock, rc) type(ESMF_GridComp) :: gcomp type(ESMF_State) :: importState, exportState type(ESMF_Clock) :: clock integer, intent(out) :: rc rc = ESMF_SUCCESS end subroutine !----------------------------------------------------------------------------- subroutine handle_verbosity_in(gcomp, verbosity, rc) type(ESMF_GridComp), intent(in) :: gcomp integer, intent(in) :: verbosity integer, intent(out) :: rc ! local variables character(ESMF_MAXSTR) :: name character(ESMF_MAXSTR) :: msgString, pLabel integer :: phase type(ESMF_Clock) :: internalClock logical :: clockIsPresent ! handle verbosity if (btest(verbosity,8)) then call ESMF_GridCompGet(gcomp, name=name, currentPhase=phase, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call NUOPC_CompSearchRevPhaseMap(gcomp, ESMF_METHOD_INITIALIZE, & phaseIndex=phase, phaseLabel=pLabel, rc=rc) if (len_trim(pLabel)==0) pLabel="<none>" if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_GridCompGet(gcomp, clockIsPresent=clockIsPresent, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (clockIsPresent) then call ESMF_GridCompGet(gcomp, clock=internalClock, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call ESMF_ClockPrint(internalClock, options="currTime", & preString=">>>"//trim(name)//& ": entered Initialize (phase="//trim(adjustl(pLabel))// & ") with current time: ", unit=msgString, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out else write(msgString,"(A)") ">>>"//trim(name)//& ": entered Initialize (phase="//trim(adjustl(pLabel))// & ") without valid internal Clock." endif call ESMF_LogWrite(msgString, ESMF_LOGMSG_INFO, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif end subroutine !----------------------------------------------------------------------------- subroutine handle_verbosity_ex(gcomp, verbosity, rc) type(ESMF_GridComp), intent(in) :: gcomp integer, intent(in) :: verbosity integer, intent(out) :: rc ! local variables character(ESMF_MAXSTR) :: name character(ESMF_MAXSTR) :: msgString, pLabel integer :: phase type(ESMF_Clock) :: internalClock logical :: clockIsPresent ! handle verbosity if (btest(verbosity,8)) then call ESMF_GridCompGet(gcomp, name=name, currentPhase=phase, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call NUOPC_CompSearchRevPhaseMap(gcomp, ESMF_METHOD_INITIALIZE, & phaseIndex=phase, phaseLabel=pLabel, rc=rc) if (len_trim(pLabel)==0) pLabel="<none>" if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_GridCompGet(gcomp, clockIsPresent=clockIsPresent, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (clockIsPresent) then call ESMF_GridCompGet(gcomp, clock=internalClock, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call ESMF_ClockPrint(internalClock, options="currTime", & preString="<<<"//trim(name)//& ": leaving Initialize (phase="//trim(adjustl(pLabel))// & ") with current time: ", unit=msgString, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out else write(msgString,"(A)") "<<<"//trim(name)//& ": leaving Initialize (phase="//trim(adjustl(pLabel))// & ") without valid internal Clock." endif call ESMF_LogWrite(msgString, ESMF_LOGMSG_INFO, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif end subroutine !----------------------------------------------------------------------------- subroutine InitializeP0(gcomp, importState, exportState, clock, rc) type(ESMF_GridComp) :: gcomp type(ESMF_State) :: importState, exportState type(ESMF_Clock) :: clock integer, intent(out) :: rc ! local variables character(*), parameter :: rName="InitializeP0" character(ESMF_MAXSTR) :: name integer :: phase type(type_InternalState) :: is integer :: verbosity, diagnostic, profiling type(ESMF_Time) :: currTime character(len=40) :: currTimeString character(ESMF_MAXSTR) :: ipdvxAttr logical :: isPresentAdvertise logical :: isPresentModifyAdvertised logical :: isPresentRealizeProvided logical :: isPresentAcceptTransfer logical :: isPresentRealizeAccepted logical :: isPresentDataInitialize rc = ESMF_SUCCESS ! query the component for info call NUOPC_CompGet(gcomp, name=name, verbosity=verbosity, & diagnostic=diagnostic, profiling=profiling, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,9)) then call ESMF_TraceRegionEnter("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_VMBarrier(rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TraceRegionExit("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(profiling,0)) then call ESMF_TraceRegionEnter(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! intro call NUOPC_LogIntro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle verbosity call handle_verbosity_in(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle diagnostic if (diagnostic>0) then call ESMF_ClockGet(clock, currTime=currTime, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TimePrint(currTime, unit=currTimeString, options="underscore", & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,0)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,1)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! determine whether the component is compatible with IPDvX call NUOPC_CompAttributeGet(gcomp, name="IPDvX", value=ipdvxAttr, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out #if 0 call ESMF_LogWrite("ipdvxAttr: "//ipdvxAttr, ESMF_LOGMSG_DEBUG, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out #endif if (trim(ipdvxAttr)=="true") then ! switch to IPDvX call NUOPC_CompFilterPhaseMap(gcomp, ESMF_METHOD_INITIALIZE, & acceptStringList=(/"IPDvXp"/), rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! Determine the active specializations call ESMF_MethodGet(gcomp, label=label_Advertise, & isPresent=isPresentAdvertise, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_MethodGet(gcomp, label=label_ModifyAdvertised, & isPresent=isPresentModifyAdvertised, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_MethodGet(gcomp, label=label_RealizeProvided, & isPresent=isPresentRealizeProvided, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_MethodGet(gcomp, label=label_AcceptTransfer, & isPresent=isPresentAcceptTransfer, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_MethodGet(gcomp, label=label_RealizeAccepted, & isPresent=isPresentRealizeAccepted, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_MethodGet(gcomp, label=label_DataInitialize, & isPresent=isPresentDataInitialize, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! Advertise: conditionally activate the generic code if (isPresentAdvertise) then call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_INITIALIZE, & phaseLabelList=(/"IPDvXp01"/), userRoutine=InitializeIPDvXp01, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out endif ! Advertise: conditionally activate the generic code if (isPresentModifyAdvertised) then call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_INITIALIZE, & phaseLabelList=(/"IPDvXp02"/), userRoutine=InitializeIPDvXp02, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out endif ! RealizeProvided: conditionally activate the generic code if (isPresentRealizeProvided) then call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_INITIALIZE, & phaseLabelList=(/"IPDvXp04"/), userRoutine=InitializeIPDvXp04, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out endif ! AcceptTransfer: conditionally activate the generic code if (isPresentAcceptTransfer) then call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_INITIALIZE, & phaseLabelList=(/"IPDvXp05"/), userRoutine=InitializeIPDvXp05, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out endif ! RealizeAccepted: conditionally activate the generic code if (isPresentRealizeAccepted) then call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_INITIALIZE, & phaseLabelList=(/"IPDvXp06"/), userRoutine=InitializeIPDvXp06, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out endif ! ConnectedCheck: always activate the generic code call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_INITIALIZE, & phaseLabelList=(/"IPDvXp07"/), userRoutine=InitializeIPDvXp07, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! DataInitialize: choose different generic code depending on spec details if (isPresentDataInitialize) then call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_INITIALIZE, & phaseLabelList=(/"IPDv05p8"/), userRoutine=InitializeIPDvXp08, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out else call NUOPC_CompSetEntryPoint(gcomp, ESMF_METHOD_INITIALIZE, & phaseLabelList=(/"IPDvXp08"/), userRoutine=InitializeIPDv01p5, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out endif else ! switch to IPDv00 (default) call NUOPC_CompFilterPhaseMap(gcomp, ESMF_METHOD_INITIALIZE, & acceptStringList=(/"IPDv00p"/), rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out endif ! query Component for the internal State nullify(is%wrap) #ifdef ESMF_NO_F2018ASSUMEDTYPE call ESMF_UserCompGetInternalState(gcomp, label_InternalState, is, rc) #else call ESMF_UserCompGetInternalState(gcomp, label_InternalState, is, rc=rc) #endif if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out ! store the incoming clock as driverClock in internal state is%wrap%driverClock = clock ! handle diagnostic if (btest(diagnostic,2)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,3)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! handle verbosity call handle_verbosity_ex(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! extro call NUOPC_LogExtro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,0)) then call ESMF_TraceRegionExit(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif end subroutine !----------------------------------------------------------------------------- subroutine InitializeIPDvXp01(gcomp, importState, exportState, clock, rc) type(ESMF_GridComp) :: gcomp type(ESMF_State) :: importState, exportState type(ESMF_Clock) :: clock integer, intent(out) :: rc ! local variables character(*), parameter :: rName="InitializeIPDvXp01" integer :: localrc character(ESMF_MAXSTR) :: name integer :: phase logical :: existflag integer :: verbosity, diagnostic, profiling type(ESMF_Time) :: currTime character(len=40) :: currTimeString rc = ESMF_SUCCESS ! query the component for info call NUOPC_CompGet(gcomp, name=name, verbosity=verbosity, & diagnostic=diagnostic, profiling=profiling, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,9)) then call ESMF_TraceRegionEnter("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_VMBarrier(rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TraceRegionExit("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(profiling,0)) then call ESMF_TraceRegionEnter(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! intro call NUOPC_LogIntro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle verbosity call handle_verbosity_in(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle diagnostic if (diagnostic>0) then call ESMF_ClockGet(clock, currTime=currTime, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TimePrint(currTime, unit=currTimeString, options="underscore", & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,0)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,1)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! SPECIALIZE by calling into attached method if (btest(profiling,1)) then call ESMF_TraceRegionEnter("label_Advertise") endif call ESMF_MethodExecute(gcomp, label=label_Advertise, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) return ! bail out if (btest(profiling,1)) then call ESMF_TraceRegionExit("label_Advertise") endif ! handle diagnostic if (btest(diagnostic,2)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,3)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! handle verbosity call handle_verbosity_ex(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! extro call NUOPC_LogExtro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,0)) then call ESMF_TraceRegionExit(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif end subroutine !----------------------------------------------------------------------------- subroutine InitializeIPDvXp02(gcomp, importState, exportState, clock, rc) type(ESMF_GridComp) :: gcomp type(ESMF_State) :: importState, exportState type(ESMF_Clock) :: clock integer, intent(out) :: rc ! local variables character(*), parameter :: rName="InitializeIPDvXp02" integer :: localrc character(ESMF_MAXSTR) :: name integer :: phase logical :: existflag integer :: verbosity, diagnostic, profiling type(ESMF_Time) :: currTime character(len=40) :: currTimeString rc = ESMF_SUCCESS ! query the component for info call NUOPC_CompGet(gcomp, name=name, verbosity=verbosity, & diagnostic=diagnostic, profiling=profiling, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,9)) then call ESMF_TraceRegionEnter("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_VMBarrier(rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TraceRegionExit("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(profiling,0)) then call ESMF_TraceRegionEnter(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! intro call NUOPC_LogIntro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle verbosity call handle_verbosity_in(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle diagnostic if (diagnostic>0) then call ESMF_ClockGet(clock, currTime=currTime, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TimePrint(currTime, unit=currTimeString, options="underscore", & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,0)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,1)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! SPECIALIZE by calling into attached method if (btest(profiling,1)) then call ESMF_TraceRegionEnter("label_ModifyAdvertised") endif call ESMF_MethodExecute(gcomp, label=label_ModifyAdvertised, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) return ! bail out if (btest(profiling,1)) then call ESMF_TraceRegionExit("label_ModifyAdvertised") endif ! handle diagnostic if (btest(diagnostic,2)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,3)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! handle verbosity call handle_verbosity_ex(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! extro call NUOPC_LogExtro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,0)) then call ESMF_TraceRegionExit(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif end subroutine !----------------------------------------------------------------------------- subroutine InitializeIPDvXp04(gcomp, importState, exportState, clock, rc) type(ESMF_GridComp) :: gcomp type(ESMF_State) :: importState, exportState type(ESMF_Clock) :: clock integer, intent(out) :: rc ! local variables character(*), parameter :: rName="InitializeIPDvXp04" integer :: localrc character(ESMF_MAXSTR) :: name integer :: phase logical :: existflag integer :: verbosity, diagnostic, profiling type(ESMF_Time) :: currTime character(len=40) :: currTimeString rc = ESMF_SUCCESS ! query the component for info call NUOPC_CompGet(gcomp, name=name, verbosity=verbosity, & diagnostic=diagnostic, profiling=profiling, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,9)) then call ESMF_TraceRegionEnter("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_VMBarrier(rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TraceRegionExit("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(profiling,0)) then call ESMF_TraceRegionEnter(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! intro call NUOPC_LogIntro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle verbosity call handle_verbosity_in(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle diagnostic if (diagnostic>0) then call ESMF_ClockGet(clock, currTime=currTime, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TimePrint(currTime, unit=currTimeString, options="underscore", & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,0)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,1)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! SPECIALIZE by calling into attached method if (btest(profiling,1)) then call ESMF_TraceRegionEnter("label_RealizeProvided") endif call ESMF_MethodExecute(gcomp, label=label_RealizeProvided, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) return ! bail out if (btest(profiling,1)) then call ESMF_TraceRegionExit("label_RealizeProvided") endif ! handle diagnostic if (btest(diagnostic,2)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,3)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! handle verbosity call handle_verbosity_ex(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! extro call NUOPC_LogExtro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,0)) then call ESMF_TraceRegionExit(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif end subroutine !----------------------------------------------------------------------------- subroutine InitializeIPDvXp05(gcomp, importState, exportState, clock, rc) type(ESMF_GridComp) :: gcomp type(ESMF_State) :: importState, exportState type(ESMF_Clock) :: clock integer, intent(out) :: rc ! local variables character(*), parameter :: rName="InitializeIPDvXp05" integer :: localrc character(ESMF_MAXSTR) :: name integer :: phase logical :: existflag integer :: verbosity, diagnostic, profiling type(ESMF_Time) :: currTime character(len=40) :: currTimeString rc = ESMF_SUCCESS ! query the component for info call NUOPC_CompGet(gcomp, name=name, verbosity=verbosity, & diagnostic=diagnostic, profiling=profiling, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,9)) then call ESMF_TraceRegionEnter("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_VMBarrier(rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TraceRegionExit("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(profiling,0)) then call ESMF_TraceRegionEnter(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! intro call NUOPC_LogIntro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle verbosity call handle_verbosity_in(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle diagnostic if (diagnostic>0) then call ESMF_ClockGet(clock, currTime=currTime, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TimePrint(currTime, unit=currTimeString, options="underscore", & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,0)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,1)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! SPECIALIZE by calling into attached method if (btest(profiling,1)) then call ESMF_TraceRegionEnter("label_AcceptTransfer") endif call ESMF_MethodExecute(gcomp, label=label_AcceptTransfer, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) return ! bail out if (btest(profiling,1)) then call ESMF_TraceRegionExit("label_AcceptTransfer") endif ! handle diagnostic if (btest(diagnostic,2)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,3)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! handle verbosity call handle_verbosity_ex(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! extro call NUOPC_LogExtro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,0)) then call ESMF_TraceRegionExit(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif end subroutine !----------------------------------------------------------------------------- subroutine InitializeIPDvXp06(gcomp, importState, exportState, clock, rc) type(ESMF_GridComp) :: gcomp type(ESMF_State) :: importState, exportState type(ESMF_Clock) :: clock integer, intent(out) :: rc ! local variables character(*), parameter :: rName="InitializeIPDvXp06" integer :: localrc character(ESMF_MAXSTR) :: name integer :: phase logical :: existflag integer :: verbosity, diagnostic, profiling type(ESMF_Time) :: currTime character(len=40) :: currTimeString rc = ESMF_SUCCESS ! query the component for info call NUOPC_CompGet(gcomp, name=name, verbosity=verbosity, & diagnostic=diagnostic, profiling=profiling, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,9)) then call ESMF_TraceRegionEnter("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_VMBarrier(rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TraceRegionExit("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(profiling,0)) then call ESMF_TraceRegionEnter(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! intro call NUOPC_LogIntro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle verbosity call handle_verbosity_in(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle diagnostic if (diagnostic>0) then call ESMF_ClockGet(clock, currTime=currTime, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TimePrint(currTime, unit=currTimeString, options="underscore", & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,0)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,1)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! SPECIALIZE by calling into attached method if (btest(profiling,1)) then call ESMF_TraceRegionEnter("label_RealizeAccepted") endif call ESMF_MethodExecute(gcomp, label=label_RealizeAccepted, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) return ! bail out if (btest(profiling,1)) then call ESMF_TraceRegionExit("label_RealizeAccepted") endif ! handle diagnostic if (btest(diagnostic,2)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,3)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! handle verbosity call handle_verbosity_ex(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! extro call NUOPC_LogExtro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,0)) then call ESMF_TraceRegionExit(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif end subroutine !----------------------------------------------------------------------------- subroutine InitializeIPDvXp07(gcomp, importState, exportState, clock, rc) type(ESMF_GridComp) :: gcomp type(ESMF_State) :: importState, exportState type(ESMF_Clock) :: clock integer, intent(out) :: rc ! local variables character(*), parameter :: rName="InitializeIPDvXp07" integer :: localrc character(ESMF_MAXSTR) :: name integer :: phase logical :: clockIsPresent logical :: clockIsCreated logical :: existflag logical :: allConnected integer :: i character(ESMF_MAXSTR), pointer :: impStdNameList(:) character(ESMF_MAXSTR), pointer :: impItemNameList(:) character(ESMF_MAXSTR), pointer :: impConnectedList(:) integer :: verbosity, diagnostic, profiling type(ESMF_Time) :: currTime character(len=40) :: currTimeString rc = ESMF_SUCCESS nullify(impStdNameList) nullify(impItemNameList) nullify(impConnectedList) ! query the component for info call NUOPC_CompGet(gcomp, name=name, verbosity=verbosity, & diagnostic=diagnostic, profiling=profiling, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,9)) then call ESMF_TraceRegionEnter("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_VMBarrier(rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TraceRegionExit("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(profiling,0)) then call ESMF_TraceRegionEnter(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! intro call NUOPC_LogIntro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle verbosity call handle_verbosity_in(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle diagnostic if (diagnostic>0) then call ESMF_ClockGet(clock, currTime=currTime, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TimePrint(currTime, unit=currTimeString, options="underscore", & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,0)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,1)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! test whether internal Clock has already been set in the Component call ESMF_GridCompGet(gcomp, clockIsPresent=clockIsPresent, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out clockIsCreated = ESMF_ClockIsCreated(clock, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out if (.not.clockIsPresent .and. clockIsCreated) then ! set the internal Clock as a copy of the incoming Clock by a default call NUOPC_CompSetClock(gcomp, clock, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out endif ! SPECIALIZE by calling into optional attached method to set internal clock if (btest(profiling,1)) then call ESMF_TraceRegionEnter("label_SetClock") endif call ESMF_MethodExecute(gcomp, label=label_SetClock, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) return ! bail out if (btest(profiling,1)) then call ESMF_TraceRegionExit("label_SetClock") endif ! query if all import Fields are connected call NUOPC_GetStateMemberLists(importState, StandardNameList=impStdNameList, & ConnectedList=impConnectedList, itemNameList=impItemNameList, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out allConnected = .true. ! initialize if (associated(impConnectedList)) then do i=1, size(impConnectedList) if (impConnectedList(i) /= "true") then allConnected = .false. call ESMF_LogWrite(trim(name)//": Import Field not connected: "// & trim(impStdNameList(i)), ESMF_LOGMSG_WARNING, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out endif enddo endif ! compatibility check if (.not.allConnected) then !TODO: introduce and use INCOMPATIBILITY return codes!!!! call ESMF_LogSetError(ESMF_RC_ARG_BAD, & msg="NUOPC INCOMPATIBILITY DETECTED: Import Fields not all connected", & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc) return ! bail out endif !TODO: remove Fields that aren't connected from import and export State if (associated(impStdNameList)) deallocate(impStdNameList) if (associated(impItemNameList)) deallocate(impItemNameList) if (associated(impConnectedList)) deallocate(impConnectedList) ! handle diagnostic if (btest(diagnostic,2)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,3)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! handle verbosity call handle_verbosity_ex(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! extro call NUOPC_LogExtro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,0)) then call ESMF_TraceRegionExit(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif end subroutine !----------------------------------------------------------------------------- subroutine InitializeIPDvXp08(gcomp, importState, exportState, clock, rc) type(ESMF_GridComp) :: gcomp type(ESMF_State) :: importState, exportState type(ESMF_Clock) :: clock integer, intent(out) :: rc ! local variables character(*), parameter :: rName="InitializeIPDvXp08" integer :: localrc integer :: phase type(ESMF_Clock) :: internalClock logical :: existflag character(ESMF_MAXSTR) :: oldDataComplete, newDataComplete integer :: oldUpdatedCount, newUpdatedCount logical :: allUpdated character(ESMF_MAXSTR) :: name, valueString1, valueString2 integer :: verbosity, diagnostic, profiling type(ESMF_Time) :: currTime character(len=40) :: currTimeString type(type_InternalState) :: is rc = ESMF_SUCCESS ! query the component for info call NUOPC_CompGet(gcomp, name=name, verbosity=verbosity, & diagnostic=diagnostic, profiling=profiling, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,9)) then call ESMF_TraceRegionEnter("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_VMBarrier(rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TraceRegionExit("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(profiling,0)) then call ESMF_TraceRegionEnter(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! intro call NUOPC_LogIntro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle verbosity call handle_verbosity_in(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle diagnostic if (diagnostic>0) then call ESMF_ClockGet(clock, currTime=currTime, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TimePrint(currTime, unit=currTimeString, options="underscore", & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,0)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,1)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! check how many Fields in the exportState have "Updated" set ! to "true" BEFORE calling the DataInitialize allUpdated = NUOPC_IsUpdated(exportState, count=oldUpdatedCount, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! get the value of the "InitializeDataComplete" attribute call NUOPC_CompAttributeGet(gcomp, name="InitializeDataComplete", & value=oldDataComplete, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! Initialize component data structures, including its export Fields, ! only connected Fields reside in exportState at this time. ! Expect the component to set "InitializeDataComplete" attribute when done. ! SPECIALIZE by calling into attached method to fill initial data if (btest(profiling,1)) then call ESMF_TraceRegionEnter("label_DataInitialize") endif call ESMF_MethodExecute(gcomp, label=label_DataInitialize, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, & rcToReturn=rc)) & return ! bail out if (btest(profiling,1)) then call ESMF_TraceRegionExit("label_DataInitialize") endif ! re-set the "InitializeDataProgress" attribute to "false" call NUOPC_CompAttributeSet(gcomp, & name="InitializeDataProgress", value="false", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=FILENAME)) return ! bail out ! check how many Fields in the exportState have "Updated" set ! to "true" AFTER calling the DataInitialize allUpdated = NUOPC_IsUpdated(exportState, count=newUpdatedCount, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! see if number of updated export fields went up if (newUpdatedCount > oldUpdatedCount) then ! there are more Fields now that have "Updated" set "true" ! -> set "InitializeDataProgress" attribute "true" call NUOPC_CompAttributeSet(gcomp, & name="InitializeDataProgress", value="true", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=FILENAME)) return ! bail out endif ! get the value of the "InitializeDataComplete" attribute call NUOPC_CompAttributeGet(gcomp, name="InitializeDataComplete", & value=newDataComplete, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! see if the "InitializeDataComplete" attribute has changed if (trim(newDataComplete) /= trim(oldDataComplete)) then ! there was a change in the "InitializeDataComplete" attribute setting ! -> set "InitializeDataProgress" attribute "true" call NUOPC_CompAttributeSet(gcomp, & name="InitializeDataProgress", value="true", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=FILENAME)) return ! bail out endif ! query Component for the internal State nullify(is%wrap) #ifdef ESMF_NO_F2018ASSUMEDTYPE call ESMF_UserCompGetInternalState(gcomp, label_InternalState, is, rc) #else call ESMF_UserCompGetInternalState(gcomp, label_InternalState, is, rc=rc) #endif if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! correct setting of timestamps nullify(is%wrap%cachedExportFieldList) call NUOPC_GetStateMemberLists(exportState, & fieldList=is%wrap%cachedExportFieldList, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call ESMF_GridCompGet(gcomp, clock=internalClock, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (allUpdated) then ! update timestamp on all the export Fields call NUOPC_SetTimestamp(is%wrap%cachedExportFieldList, internalClock, & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out else ! update timestamp on only those export Fields that have the ! "Updated" attribute set to "true" call NUOPC_SetTimestamp(is%wrap%cachedExportFieldList, internalClock, & selective=.true., rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out endif ! handle verbosity if (btest(verbosity,11)) then call NUOPC_CompAttributeGet(gcomp, name="InitializeDataProgress", & value=valueString1, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call NUOPC_CompAttributeGet(gcomp, name="InitializeDataComplete", & value=valueString2, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call ESMF_LogWrite(trim(name)//": InitializeDataProgress='"// & trim(valueString1)//"', InitializeDataComplete='"// & trim(valueString2)//"'", ESMF_LOGMSG_INFO, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! handle diagnostic if (btest(diagnostic,2)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,3)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! handle verbosity call handle_verbosity_ex(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! extro call NUOPC_LogExtro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,0)) then call ESMF_TraceRegionExit(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif end subroutine !----------------------------------------------------------------------------- subroutine InitializeIPDv01p5(gcomp, importState, exportState, clock, rc) type(ESMF_GridComp) :: gcomp type(ESMF_State) :: importState, exportState type(ESMF_Clock) :: clock integer, intent(out) :: rc ! local variables character(*), parameter :: rName="InitializeIPDv01p5" integer :: localrc logical :: existflag character(ESMF_MAXSTR) :: name integer :: phase type(ESMF_Clock) :: internalClock logical :: clockIsCreated integer :: verbosity, diagnostic, profiling type(ESMF_Time) :: currTime character(len=40) :: currTimeString type(type_InternalState) :: is rc = ESMF_SUCCESS ! query the component for info call NUOPC_CompGet(gcomp, name=name, verbosity=verbosity, & diagnostic=diagnostic, profiling=profiling, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,9)) then call ESMF_TraceRegionEnter("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_VMBarrier(rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TraceRegionExit("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(profiling,0)) then call ESMF_TraceRegionEnter(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! intro call NUOPC_LogIntro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle verbosity call handle_verbosity_in(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle diagnostic if (diagnostic>0) then call ESMF_ClockGet(clock, currTime=currTime, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TimePrint(currTime, unit=currTimeString, options="underscore", & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,0)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,1)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! if the incoming clock is valid, then use to set currTime on internalClock clockIsCreated = ESMF_ClockIsCreated(clock, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out if (clockIsCreated) then ! reset the currTime of the internalClock call ESMF_ClockGet(clock, currTime=currTime, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_GridCompGet(gcomp, clock=internalClock, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call ESMF_ClockSet(internalClock, currTime=currTime, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out endif ! fill all export Fields with valid initial data for current time ! note that only connected Fields reside in exportState at this time ! SPECIALIZE by calling into attached method to fill initial data if (btest(profiling,1)) then call ESMF_TraceRegionEnter("label_DataInitialize") endif call ESMF_MethodExecute(gcomp, label=label_DataInitialize, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, & rcToReturn=rc)) & return ! bail out if (btest(profiling,1)) then call ESMF_TraceRegionExit("label_DataInitialize") endif ! query Component for the internal State nullify(is%wrap) #ifdef ESMF_NO_F2018ASSUMEDTYPE call ESMF_UserCompGetInternalState(gcomp, label_InternalState, is, rc) #else call ESMF_UserCompGetInternalState(gcomp, label_InternalState, is, rc=rc) #endif if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! update timestamp on all the export Fields nullify(is%wrap%cachedExportFieldList) call NUOPC_GetStateMemberLists(exportState, & fieldList=is%wrap%cachedExportFieldList, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call ESMF_GridCompGet(gcomp, clock=internalClock, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call NUOPC_SetTimestamp(is%wrap%cachedExportFieldList, internalClock, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! indicate data complete because no actual checking for this IPDv/situation call NUOPC_CompAttributeSet(gcomp, name="InitializeDataComplete", & value="true", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! handle diagnostic if (btest(diagnostic,2)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,3)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! handle verbosity call handle_verbosity_ex(gcomp, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! extro call NUOPC_LogExtro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,0)) then call ESMF_TraceRegionExit(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif end subroutine !----------------------------------------------------------------------------- subroutine routine_Run(gcomp, importState, exportState, clock, rc) type(ESMF_GridComp) :: gcomp type(ESMF_State) :: importState, exportState type(ESMF_Clock) :: clock integer, intent(out) :: rc ! local variables character(*), parameter :: rName="Run" integer :: localrc type(type_InternalState) :: is type(ESMF_Clock) :: internalClock logical :: allCurrent logical :: existflag character(ESMF_MAXSTR) :: msgString, pLabel integer :: phase integer :: verbosity, diagnostic, profiling character(ESMF_MAXSTR) :: name type(ESMF_Time) :: currTime character(len=40) :: currTimeString logical :: exportIsCreated rc = ESMF_SUCCESS ! query the component for info call NUOPC_CompGet(gcomp, name=name, verbosity=verbosity, & diagnostic=diagnostic, profiling=profiling, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,10)) then call ESMF_TraceRegionEnter("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_VMBarrier(rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TraceRegionExit("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(profiling,3)) then call ESMF_TraceRegionEnter(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! intro call NUOPC_LogIntro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! get the currentPhase call ESMF_GridCompGet(gcomp, currentPhase=phase, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! handle diagnostic if (diagnostic>0) then call ESMF_ClockGet(clock, currTime=currTime, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TimePrint(currTime, unit=currTimeString, options="underscore", & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,4)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,5)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! query Component for the internal State nullify(is%wrap) #ifdef ESMF_NO_F2018ASSUMEDTYPE call ESMF_UserCompGetInternalState(gcomp, label_InternalState, is, rc) #else call ESMF_UserCompGetInternalState(gcomp, label_InternalState, is, rc=rc) #endif if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out ! store the incoming clock as driverClock in internal state is%wrap%driverClock = clock ! SPECIALIZE required: label_SetRunClock if (btest(profiling,4)) then call ESMF_TraceRegionEnter("label_SetRunClock") endif ! -> first check for the label with phase index call ESMF_MethodExecute(gcomp, label=label_SetRunClock, index=phase, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, & rcToReturn=rc)) & return ! bail out if (.not.existflag) then ! -> next check for the label without phase index call ESMF_MethodExecute(gcomp, label=label_SetRunClock, & userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, & rcToReturn=rc)) & return ! bail out endif if (btest(profiling,4)) then call ESMF_TraceRegionExit("label_SetRunClock") endif ! get the internal clock for the time stepping loop call ESMF_GridCompGet(gcomp, clock=internalClock, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! store the current time of the internalClock before advancing it call ESMF_ClockGet(internalClock, currTime=is%wrap%preAdvanceCurrTime, & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! handle verbosity if (btest(verbosity,9)) then call NUOPC_CompSearchRevPhaseMap(gcomp, ESMF_METHOD_RUN, & phaseIndex=phase, phaseLabel=pLabel, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_ClockPrint(internalClock, options="currTime", & preString=">>>"//trim(name)//& ": entered Run (phase="//trim(adjustl(pLabel))// & ") with current time: ", unit=msgString, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_LogWrite(msgString, ESMF_LOGMSG_INFO, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! SPECIALIZE optionally: label_CheckImport if (btest(profiling,4)) then call ESMF_TraceRegionEnter("label_CheckImport") endif ! -> first check for the label with phase index call ESMF_MethodExecute(gcomp, label=label_CheckImport, index=phase, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, & rcToReturn=rc)) & return ! bail out if (.not.existflag) then ! -> next check for the label without phase index call ESMF_MethodExecute(gcomp, label=label_CheckImport, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, & rcToReturn=rc)) & return ! bail out endif if (btest(profiling,4)) then call ESMF_TraceRegionExit("label_CheckImport") endif ! model time stepping loop do while (.not. ESMF_ClockIsStopTime(internalClock, rc=rc)) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle verbosity if (btest(verbosity,12)) then call ESMF_ClockPrint(internalClock, options="currTime", & preString=trim(name)//": time step-loop starting, current time: ", & unit=msgString, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_LogWrite(msgString, ESMF_LOGMSG_INFO, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! advance the model t->t+dt ! SPECIALIZE required: label_Advance if (btest(profiling,4)) then call ESMF_TraceRegionEnter("label_Advance") endif ! -> first check for the label with phase index call ESMF_MethodExecute(gcomp, label=label_Advance, index=phase, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, & rcToReturn=rc)) & return ! bail out if (.not.existflag) then ! -> next check for the label without phase index call ESMF_MethodExecute(gcomp, label=label_Advance, userRc=localrc, & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, & rcToReturn=rc)) & return ! bail out endif if (btest(profiling,4)) then call ESMF_TraceRegionExit("label_Advance") endif if (btest(profiling,4)) then call ESMF_TraceRegionEnter("label_AdvanceClock") endif ! advance the internalClock to the new current time (optionally specialz) call ESMF_MethodExecute(gcomp, label=label_AdvanceClock, index=phase, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, & rcToReturn=rc)) & return ! bail out if (.not.existflag) then ! -> next check for the label without phase index call ESMF_MethodExecute(gcomp, label=label_AdvanceClock, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, & rcToReturn=rc)) & return ! bail out if (.not.existflag) then ! at last use the DEFAULT implementation to advance the Clock call ESMF_ClockAdvance(internalClock, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out endif endif if (btest(profiling,4)) then call ESMF_TraceRegionExit("label_AdvanceClock") endif ! handle verbosity if (btest(verbosity,12)) then call ESMF_ClockPrint(internalClock, options="currTime", & preString=trim(name)//": time step-loop ending, current time: ", & unit=msgString, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_LogWrite(msgString, ESMF_LOGMSG_INFO, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif enddo ! end of time stepping loop if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out exportIsCreated = ESMF_StateIsCreated(exportState, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (exportIsCreated) then ! SPECIALIZE optionally: label_TimestampExport if (btest(profiling,4)) then call ESMF_TraceRegionEnter("label_TimestampExport") endif ! -> first check for the label with phase index call ESMF_MethodExecute(gcomp, label=label_TimestampExport, index=phase, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, & rcToReturn=rc)) & return ! bail out if (.not.existflag) then ! -> next check for the label without phase index call ESMF_MethodExecute(gcomp, label=label_TimestampExport, & existflag=existflag, userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, & rcToReturn=rc)) & return ! bail out endif if (btest(profiling,4)) then call ESMF_TraceRegionExit("label_TimestampExport") endif endif ! handle diagnostic if (btest(diagnostic,6)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,7)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! handle verbosity if (btest(verbosity,9)) then call ESMF_ClockPrint(internalClock, options="currTime", & preString="<<<"//trim(name)//& ": leaving Run (phase="//trim(adjustl(pLabel))// & ") with current time: ", unit=msgString, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_LogWrite(msgString, ESMF_LOGMSG_INFO, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! extro call NUOPC_LogExtro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,3)) then call ESMF_TraceRegionExit(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif end subroutine !----------------------------------------------------------------------------- subroutine SetRunClock(gcomp, rc) type(ESMF_GridComp) :: gcomp integer, intent(out) :: rc ! Set the model clock according to the incoming driver clock. ! Implement the default explicit timekeeping rules. ! local variables type(type_InternalState) :: is character(ESMF_MAXSTR):: name rc = ESMF_SUCCESS ! query the component for info call NUOPC_CompGet(gcomp, name=name, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! query component for its internal state nullify(is%wrap) #ifdef ESMF_NO_F2018ASSUMEDTYPE call ESMF_UserCompGetInternalState(gcomp, label_InternalState, is, rc) #else call ESMF_UserCompGetInternalState(gcomp, label_InternalState, is, rc=rc) #endif if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! check and set the model clock against the driver clock call NUOPC_CompCheckSetClock(gcomp, is%wrap%driverClock, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out end subroutine !----------------------------------------------------------------------------- subroutine CheckImport(gcomp, rc) type(ESMF_GridComp) :: gcomp integer, intent(out) :: rc ! This is the routine that enforces the explicit time dependence on the ! import fields. This simply means that the timestamps on the Fields in the ! importState are checked against the currentTime on the Component's ! internalClock. Consequenty, this model starts out with forcing fields ! at the current time as it does its forward step from currentTime to ! currentTime + timeStep. ! local variables type(ESMF_Clock) :: clock type(ESMF_Time) :: time type(ESMF_State) :: importState logical :: allCurrent type(ESMF_Field), allocatable :: fieldList(:) integer :: i character(ESMF_MAXSTR) :: fieldName character(ESMF_MAXSTR) :: name character(ESMF_MAXSTR) :: valueString rc = ESMF_SUCCESS ! query the component for info call NUOPC_CompGet(gcomp, name=name, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! query the Component for its clock and importState call ESMF_GridCompGet(gcomp, clock=clock, importState=importState, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! get the current time out of the clock call ESMF_ClockGet(clock, currTime=time, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out ! check that Fields in the importState show correct timestamp allCurrent = NUOPC_IsAtTime(importState, time, fieldList=fieldList, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out if (.not.allCurrent) then !TODO: introduce and use INCOMPATIBILITY return codes!!!! do i=1, size(fieldList) call ESMF_FieldGet(fieldList(i), name=fieldName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call NUOPC_GetAttribute(fieldList(i), name="StandardName", & value=valueString, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call ESMF_LogWrite(trim(name)//": Field '"//trim(fieldName)//& "' in the importState is not at the expected time. StandardName: "& //trim(valueString), ESMF_LOGMSG_WARNING) enddo deallocate(fieldList) call ESMF_LogSetError(ESMF_RC_ARG_BAD, & msg="NUOPC INCOMPATIBILITY DETECTED: Import Fields not at current time", & line=__LINE__, file=trim(name)//":"//FILENAME, & rcToReturn=rc) return ! bail out endif end subroutine !----------------------------------------------------------------------------- recursive subroutine Finalize(gcomp, importState, exportState, clock, rc) type(ESMF_GridComp) :: gcomp type(ESMF_State) :: importState, exportState type(ESMF_Clock) :: clock integer, intent(out) :: rc ! local variables character(*), parameter :: rName="Finalize" integer :: localrc, stat logical :: existflag character(ESMF_MAXSTR) :: name character(ESMF_MAXSTR) :: msgString, pLabel integer :: phase type(ESMF_Clock) :: internalClock integer :: verbosity, diagnostic, profiling type(type_InternalState) :: is type(ESMF_Time) :: currTime character(len=40) :: currTimeString rc = ESMF_SUCCESS ! query the component for info call NUOPC_CompGet(gcomp, name=name, verbosity=verbosity, & diagnostic=diagnostic, profiling=profiling, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,11)) then call ESMF_TraceRegionEnter("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_VMBarrier(rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TraceRegionExit("Leading Barrier", rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(profiling,6)) then call ESMF_TraceRegionEnter(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! intro call NUOPC_LogIntro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle verbosity if (btest(verbosity,10)) then call ESMF_GridCompGet(gcomp, currentPhase=phase, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call NUOPC_CompSearchRevPhaseMap(gcomp, ESMF_METHOD_FINALIZE, & phaseIndex=phase, phaseLabel=pLabel, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_GridCompGet(gcomp, clock=internalClock, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) & return ! bail out call ESMF_ClockPrint(internalClock, options="currTime", & preString=">>>"//trim(name)//& ": entered Finalize (phase="//trim(adjustl(pLabel))// & ") with current time: ", unit=msgString, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_LogWrite(msgString, ESMF_LOGMSG_INFO, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! handle diagnostic if (diagnostic>0) then call ESMF_ClockGet(clock, currTime=currTime, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out call ESMF_TimePrint(currTime, unit=currTimeString, options="underscore", & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,8)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,9)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_enter_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! SPECIALIZE by calling into optional attached method if (btest(profiling,7)) then call ESMF_TraceRegionEnter("label_Finalize") endif call ESMF_MethodExecute(gcomp, label=label_Finalize, existflag=existflag, & userRc=localrc, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out if (btest(profiling,7)) then call ESMF_TraceRegionExit("label_Finalize") endif ! query Component for the internal State nullify(is%wrap) #ifdef ESMF_NO_F2018ASSUMEDTYPE call ESMF_UserCompGetInternalState(gcomp, label_InternalState, is, rc) #else call ESMF_UserCompGetInternalState(gcomp, label_InternalState, is, rc=rc) #endif if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out ! destroy relevant internal State members if (associated(is%wrap%cachedExportFieldList)) & deallocate(is%wrap%cachedExportFieldList) ! deallocate internal state memory deallocate(is%wrap, stat=stat) if (ESMF_LogFoundDeallocError(statusToCheck=stat, & msg="Deallocation of internal state memory failed.", & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out ! handle diagnostic if (btest(diagnostic,10)) then call NUOPC_Write(importState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_import_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif if (btest(diagnostic,11)) then call NUOPC_Write(exportState, fileNamePrefix="diagnostic_"//& trim(name)//"_"//trim(rName)//"_exit_export_"//trim(currTimeString)//& "_", timeslice=1, status=ESMF_FILESTATUS_REPLACE, relaxedFlag=.true., & rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! handle verbosity if (btest(verbosity,10)) then call ESMF_ClockPrint(internalClock, options="currTime", & preString="<<<"//trim(name)//& ": leaving Finalize (phase="//trim(adjustl(pLabel))// & ") with current time: ", unit=msgString, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out call ESMF_LogWrite(msgString, ESMF_LOGMSG_INFO, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif ! extro call NUOPC_LogExtro(name, rName, verbosity, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME)) return ! bail out ! handle profiling if (btest(profiling,6)) then call ESMF_TraceRegionExit(rName, rc=rc) if (ESMF_LogFoundError(rcToCheck=rc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out endif end subroutine !----------------------------------------------------------------------------- !----------------------------------------------------------------------------- !----------------------------------------------------------------------------- !----------------------------------------------------------------------------- !BOP ! !IROUTINE: NUOPC_ModelBaseGet - Get info from a ModelBase ! ! !INTERFACE: subroutine NUOPC_ModelBaseGet(gcomp, driverClock, clock, importState, & exportState, rc) ! !ARGUMENTS: type(ESMF_GridComp) :: gcomp type(ESMF_Clock), intent(out), optional :: driverClock type(ESMF_Clock), intent(out), optional :: clock type(ESMF_State), intent(out), optional :: importState type(ESMF_State), intent(out), optional :: exportState integer, intent(out), optional :: rc ! ! !DESCRIPTION: ! Access ModelBase information. !EOP !----------------------------------------------------------------------------- ! local variables integer :: localrc character(ESMF_MAXSTR) :: name type(type_InternalState) :: is if (present(rc)) rc = ESMF_SUCCESS ! query the component for info call NUOPC_CompGet(gcomp, name=name, rc=localrc) if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) return ! bail out ! driverClock if (present(driverClock)) then ! query Component for the internal State nullify(is%wrap) #ifdef ESMF_NO_F2018ASSUMEDTYPE call ESMF_UserCompGetInternalState(gcomp, label_InternalState, is, localrc) #else call ESMF_UserCompGetInternalState(gcomp, label_InternalState, is, rc=localrc) #endif if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) & return ! bail out driverClock = is%wrap%driverClock endif ! remaining arguments call ESMF_GridCompGet(gcomp, clock=clock, importState=importState, & exportState=exportState, rc=localrc) if (ESMF_LogFoundError(rcToCheck=localrc, msg=ESMF_LOGERR_PASSTHRU, & line=__LINE__, file=trim(name)//":"//FILENAME, rcToReturn=rc)) return ! bail out end subroutine !----------------------------------------------------------------------------- end module