C $Header: /home/ubuntu/mnt/e9_copy/MITgcm/model/src/initialise_varia.F,v 1.17 2001/08/27 18:50:41 jmc Exp $ C $Name: $ #include "CPP_OPTIONS.h" CStartOfInterface SUBROUTINE INITIALISE_VARIA(myThid) C /==========================================================\ C | SUBROUTINE INITIALISE_VARIA | C | o Set the initial conditions for dynamics variables | C | and time dependent arrays | C |==========================================================| C | This routine reads/writes data from an input file and | C | from various binary files. | C | Each thread invokes an instance of this routine as does | C | each process in a multi-process parallel environment like| C | MPI. | C \==========================================================/ IMPLICIT NONE C == Global variables == #include "SIZE.h" #include "EEPARAMS.h" #include "PARAMS.h" #include "DYNVARS.h" C == Routine arguments == INTEGER myThid CEndOfInterface C == Local variables == INTEGER bi,bj,K,iMin,iMax,jMin,jMax #ifdef ALLOW_TAMC_CHECKPOINTING nIter0 = INT( startTime/deltaTClock ) C-- Set Bo_surf => define the Linear Relation: Phi_surf(eta) CALL INI_LINEAR_PHISURF( myThid ) C-- Set coriolis operators CALL INI_CORI( myThid ) C-- Set laplace operators for use in 2D conjugate gradient solver. CALL INI_CG2D( myThid ) #ifdef ALLOW_NONHYDROSTATIC C-- Set laplace operators for use in 3D conjugate gradient solver. CALL INI_CG3D( myThid ) #endif #endif /* ALLOW_TAMC_CHECKPOINTING */ _BARRIER C-- Initialise 3-dim. diffusivities CALL INI_MIXING( myThid ) _BARRIER C-- Initialize DYNVARS arrays (state fields + G terms: Gu,Gv,...) to zero [always] CALL INI_DYNVARS( myThid ) C-- Initialise model fields. C Starting values of U, V, W, temp., salt. and tendency terms C are set here. Fields are either set to default or read from C stored files. CALL INI_FIELDS( myThid ) _BARRIER #ifdef ALLOW_PASSIVE_TRACER C-- Initialise passive tracer(s) CALL INI_TR1( myThid ) _BARRIER #endif IF ( usePickupBeforeC35 ) THEN C-- IMPORTANT : Need to activate the following call to restart from C a pickup file written by MITgcmUV_checkpoint34 or earlier. IF ( startTime .NE. 0. ) THEN CALL THE_CORRECTION_STEP(startTime, nIter0, myThid) ENDIF ENDIF C-- Initial conditions are convectively adjusted (for historical reasons) IF ( startTime .EQ. 0. ) THEN CADJ loop = parallel DO bj = myByLo(myThid), myByHi(myThid) CADJ loop = parallel DO bi = myBxLo(myThid), myBxHi(myThid) iMin=1-Olx iMax=sNx+Olx jMin=1-Oly jMax=sNy+Oly CALL CONVECTIVE_ADJUSTMENT_INI( I bi, bj, iMin, iMax, jMin, jMax, I startTime, nIter0, myThid ) ENDDO ENDDO _BARRIER END IF C-- Initialize variable data for packages CALL PACKAGES_INIT_VARIABLES( myThid ) #ifdef NONLIN_FRSURF C-- Compute the surface level thickness, function of eta_n-1 C and modify hFac(C,W,S) accordingly : IF (nonlinFreeSurf.GT.0) THEN CALL CALC_SURF_DR(etaNm1, startTime, nIter0, myThid ) CALL UPDATE_SURF_DR( startTime, nIter0, myThid ) ENDIF C- update also CG2D matrix (and preconditioner) IF ( nonlinFreeSurf.GT.2) THEN CALL UPDATE_CG2D( startTime, nIter0, myThid ) ENDIF #endif C-- Finally summarise the model state CALL STATE_SUMMARY( myThid ) #ifdef ALLOW_TIMEAVE C-- initialise time-average arrays with initial state values IF (taveFreq.GT.0.) THEN DO bj=myByLo(myThid),myByHi(myThid) DO bi=myBxLo(myThid),myBxHi(myThid) CALL TIMEAVE_STATVARS(startTime, nIter0, bi, bj, myThid) ENDDO ENDDO ENDIF #endif /* ALLOW_TIMEAVE */ #ifdef EXACT_CONSERV IF (exactConserv) THEN C-- Compute the initial Barotropic Flow Divergence : DO bj=myByLo(myThid),myByHi(myThid) DO bi=myBxLo(myThid),myBxHi(myThid) CALL CALC_EXACT_ETA( bi,bj, uVel,vVel, I startTime, nIter0, myThid ) ENDDO ENDDO ENDIF #endif /* EXACT_CONSERV */ RETURN END