/[MITgcm]/MITgcm/model/src/initialise_varia.F
ViewVC logotype

Annotation of /MITgcm/model/src/initialise_varia.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.17 - (hide annotations) (download)
Mon Aug 27 18:50:41 2001 UTC (22 years, 8 months ago) by jmc
Branch: MAIN
Changes since 1.16: +15 -2 lines
modified to incorporate NonLin-FreeSurf

1 jmc 1.17 C $Header: /u/gcmpack/models/MITgcmUV/model/src/initialise_varia.F,v 1.16 2001/08/13 18:07:35 heimbach Exp $
2     C $Name: $
3 adcroft 1.1
4     #include "CPP_OPTIONS.h"
5    
6     CStartOfInterface
7     SUBROUTINE INITIALISE_VARIA(myThid)
8     C /==========================================================\
9     C | SUBROUTINE INITIALISE_VARIA |
10     C | o Set the initial conditions for dynamics variables |
11     C | and time dependent arrays |
12     C |==========================================================|
13     C | This routine reads/writes data from an input file and |
14     C | from various binary files. |
15     C | Each thread invokes an instance of this routine as does |
16     C | each process in a multi-process parallel environment like|
17     C | MPI. |
18     C \==========================================================/
19     IMPLICIT NONE
20    
21     C == Global variables ==
22     #include "SIZE.h"
23     #include "EEPARAMS.h"
24     #include "PARAMS.h"
25 adcroft 1.5 #include "DYNVARS.h"
26 adcroft 1.1
27     C == Routine arguments ==
28     INTEGER myThid
29     CEndOfInterface
30    
31     C == Local variables ==
32 adcroft 1.5 INTEGER bi,bj,K,iMin,iMax,jMin,jMax
33 heimbach 1.14
34     #ifdef ALLOW_TAMC_CHECKPOINTING
35 heimbach 1.15
36     nIter0 = INT( startTime/deltaTClock )
37    
38     C-- Set Bo_surf => define the Linear Relation: Phi_surf(eta)
39     CALL INI_LINEAR_PHISURF( myThid )
40    
41     C-- Set coriolis operators
42     CALL INI_CORI( myThid )
43    
44     C-- Set laplace operators for use in 2D conjugate gradient solver.
45     CALL INI_CG2D( myThid )
46    
47     #ifdef ALLOW_NONHYDROSTATIC
48     C-- Set laplace operators for use in 3D conjugate gradient solver.
49     CALL INI_CG3D( myThid )
50 heimbach 1.14 #endif
51 adcroft 1.1
52 heimbach 1.15 #endif /* ALLOW_TAMC_CHECKPOINTING */
53 heimbach 1.16 _BARRIER
54 heimbach 1.15
55 heimbach 1.16 C-- Initialise 3-dim. diffusivities
56     CALL INI_MIXING( myThid )
57 adcroft 1.5 _BARRIER
58 heimbach 1.16
59 adcroft 1.12 C-- Initialize DYNVARS arrays (state fields + G terms: Gu,Gv,...) to zero [always]
60     CALL INI_DYNVARS( myThid )
61    
62 adcroft 1.1 C-- Initialise model fields.
63     C Starting values of U, V, W, temp., salt. and tendency terms
64     C are set here. Fields are either set to default or read from
65     C stored files.
66     CALL INI_FIELDS( myThid )
67     _BARRIER
68 jmc 1.8
69 heimbach 1.15 #ifdef ALLOW_PASSIVE_TRACER
70 heimbach 1.13 C-- Initialise passive tracer(s)
71     CALL INI_TR1( myThid )
72     _BARRIER
73 heimbach 1.15 #endif
74 heimbach 1.13
75 heimbach 1.11 IF ( usePickupBeforeC35 ) THEN
76 jmc 1.8 C-- IMPORTANT : Need to activate the following call to restart from
77     C a pickup file written by MITgcmUV_checkpoint34 or earlier.
78 heimbach 1.11 IF ( startTime .NE. 0. ) THEN
79     CALL THE_CORRECTION_STEP(startTime, nIter0, myThid)
80     ENDIF
81     ENDIF
82 adcroft 1.1
83 adcroft 1.5 C-- Initial conditions are convectively adjusted (for historical reasons)
84 heimbach 1.3 IF ( startTime .EQ. 0. ) THEN
85 heimbach 1.13 CADJ loop = parallel
86 adcroft 1.1 DO bj = myByLo(myThid), myByHi(myThid)
87 heimbach 1.13 CADJ loop = parallel
88 adcroft 1.1 DO bi = myBxLo(myThid), myBxHi(myThid)
89 adcroft 1.5 iMin=1-Olx
90     iMax=sNx+Olx
91     jMin=1-Oly
92     jMax=sNy+Oly
93 heimbach 1.11 CALL CONVECTIVE_ADJUSTMENT_INI(
94 adcroft 1.5 I bi, bj, iMin, iMax, jMin, jMax,
95 jmc 1.7 I startTime, nIter0, myThid )
96 adcroft 1.1 ENDDO
97     ENDDO
98 adcroft 1.5 _BARRIER
99     END IF
100    
101     C-- Initialize variable data for packages
102     CALL PACKAGES_INIT_VARIABLES( myThid )
103 jmc 1.17
104     #ifdef NONLIN_FRSURF
105     C-- Compute the surface level thickness, function of eta_n-1
106     C and modify hFac(C,W,S) accordingly :
107     IF (nonlinFreeSurf.GT.0) THEN
108     CALL CALC_SURF_DR(etaNm1, startTime, nIter0, myThid )
109     CALL UPDATE_SURF_DR( startTime, nIter0, myThid )
110     ENDIF
111     C- update also CG2D matrix (and preconditioner)
112     IF ( nonlinFreeSurf.GT.2) THEN
113     CALL UPDATE_CG2D( startTime, nIter0, myThid )
114     ENDIF
115     #endif
116 adcroft 1.1
117     C-- Finally summarise the model state
118     CALL STATE_SUMMARY( myThid )
119    
120 jmc 1.9 #ifdef ALLOW_TIMEAVE
121 jmc 1.7 C-- initialise time-average arrays with initial state values
122     IF (taveFreq.GT.0.) THEN
123     DO bj=myByLo(myThid),myByHi(myThid)
124     DO bi=myBxLo(myThid),myBxHi(myThid)
125 jmc 1.9 CALL TIMEAVE_STATVARS(startTime, nIter0, bi, bj, myThid)
126 jmc 1.7 ENDDO
127     ENDDO
128     ENDIF
129 jmc 1.9 #endif /* ALLOW_TIMEAVE */
130 jmc 1.7
131 adcroft 1.12 #ifdef EXACT_CONSERV
132     IF (exactConserv) THEN
133     C-- Compute the initial Barotropic Flow Divergence :
134     DO bj=myByLo(myThid),myByHi(myThid)
135     DO bi=myBxLo(myThid),myBxHi(myThid)
136     CALL CALC_EXACT_ETA( bi,bj, uVel,vVel,
137     I startTime, nIter0, myThid )
138     ENDDO
139     ENDDO
140     ENDIF
141     #endif /* EXACT_CONSERV */
142    
143 jmc 1.7 RETURN
144 adcroft 1.1 END

  ViewVC Help
Powered by ViewVC 1.1.22