1 |
C $Header: /u/u3/gcmpack/MITgcm/model/src/initialise_varia.F,v 1.31.2.2 2003/10/02 18:10:45 edhill Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "PACKAGES_CONFIG.h" |
5 |
#include "CPP_OPTIONS.h" |
6 |
|
7 |
CBOP |
8 |
C !ROUTINE: INITIALISE_VARIA |
9 |
C !INTERFACE: |
10 |
SUBROUTINE INITIALISE_VARIA(myThid) |
11 |
C !DESCRIPTION: \bv |
12 |
C *==========================================================* |
13 |
C | SUBROUTINE INITIALISE_VARIA |
14 |
C | o Set the initial conditions for dynamics variables |
15 |
C | and time dependent arrays |
16 |
C *==========================================================* |
17 |
C | This routine reads/writes data from an input file and |
18 |
C | from various binary files. |
19 |
C | Each thread invokes an instance of this routine as does |
20 |
C | each process in a multi-process parallel environment like |
21 |
C | MPI. |
22 |
C *==========================================================* |
23 |
C \ev |
24 |
|
25 |
C !CALLING SEQUENCE: |
26 |
C INITIALISE_VARIA |
27 |
C | |
28 |
C |-- INI_LINEAR_PHISURF |
29 |
C | |
30 |
C |-- INI_CORI |
31 |
C | |
32 |
C |-- INI_CG2D |
33 |
C | |
34 |
C |-- INI_CG3D |
35 |
C | |
36 |
C |-- INI_MIXING |
37 |
C | |
38 |
C |-- INI_DYNVARS |
39 |
C | |
40 |
C |-- INI_FIELDS |
41 |
C | |
42 |
C |-- INI_TR1 |
43 |
C | |
44 |
C |-- INI_AUTODIFF |
45 |
C | |
46 |
C |-- PACKAGES_INIT_VARIABLES |
47 |
C | |
48 |
C |-- THE_CORRECTION_STEP (if restart from old pickup files) |
49 |
C | |
50 |
C |-- CALL CONVECTIVE_ADJUSTMENT_INI |
51 |
C | |
52 |
C |-- CALC_SURF_DR |
53 |
C | |
54 |
C |-- UPDATE_SURF_DR |
55 |
C | |
56 |
C |-- UPDATE_CG2D |
57 |
C | |
58 |
C |-- INTEGR_CONTINUITY |
59 |
C | |
60 |
C |-- TIMEAVE_STATVARS |
61 |
C | |
62 |
C |-- STATE_SUMMARY |
63 |
|
64 |
C !USES: |
65 |
IMPLICIT NONE |
66 |
C == Global variables == |
67 |
#include "SIZE.h" |
68 |
#include "EEPARAMS.h" |
69 |
#include "PARAMS.h" |
70 |
#include "DYNVARS.h" |
71 |
|
72 |
C !INPUT/OUTPUT PARAMETERS: |
73 |
C == Routine arguments == |
74 |
INTEGER myThid |
75 |
|
76 |
C !LOCAL VARIABLES: |
77 |
C == Local variables == |
78 |
INTEGER bi,bj,K,iMin,iMax,jMin,jMax |
79 |
CEOP |
80 |
|
81 |
#ifdef ALLOW_TAMC_CHECKPOINTING |
82 |
|
83 |
nIter0 = INT( startTime/deltaTClock ) |
84 |
|
85 |
C-- Set Bo_surf => define the Linear Relation: Phi_surf(eta) |
86 |
CALL INI_LINEAR_PHISURF( myThid ) |
87 |
|
88 |
C-- Set coriolis operators |
89 |
CALL INI_CORI( myThid ) |
90 |
|
91 |
C-- Set laplace operators for use in 2D conjugate gradient solver. |
92 |
CALL INI_CG2D( myThid ) |
93 |
|
94 |
#ifdef ALLOW_NONHYDROSTATIC |
95 |
C-- Set laplace operators for use in 3D conjugate gradient solver. |
96 |
ceh3 should add an IF ( useNONHYDROSTATIC ) THEN |
97 |
CALL INI_CG3D( myThid ) |
98 |
#endif |
99 |
|
100 |
#endif /* ALLOW_TAMC_CHECKPOINTING */ |
101 |
_BARRIER |
102 |
|
103 |
C-- Initialise 3-dim. diffusivities |
104 |
CALL INI_MIXING( myThid ) |
105 |
_BARRIER |
106 |
|
107 |
C-- Initialize DYNVARS arrays (state fields + G terms: Gu,Gv,...) to zero [always] |
108 |
CALL INI_DYNVARS( myThid ) |
109 |
|
110 |
C-- Initialise model fields. |
111 |
C Starting values of U, V, W, temp., salt. and tendency terms |
112 |
C are set here. Fields are either set to default or read from |
113 |
C stored files. |
114 |
CALL INI_FIELDS( myThid ) |
115 |
_BARRIER |
116 |
|
117 |
#ifdef ALLOW_PASSIVE_TRACER |
118 |
C-- Initialise passive tracer(s) |
119 |
CALL INI_TR1( myThid ) |
120 |
_BARRIER |
121 |
#endif |
122 |
|
123 |
#ifdef ALLOW_AUTODIFF_TAMC |
124 |
C-- Initialise active fields to help TAMC |
125 |
CALL INI_AUTODIFF( myThid ) |
126 |
_BARRIER |
127 |
#endif |
128 |
|
129 |
C-- Initialize variable data for packages |
130 |
CALL PACKAGES_INIT_VARIABLES( myThid ) |
131 |
|
132 |
#ifndef ALLOW_AUTODIFF_TAMC |
133 |
IF ( usePickupBeforeC35 ) THEN |
134 |
C-- IMPORTANT : Need to activate the following call to restart from |
135 |
C a pickup file written by MITgcmUV_checkpoint34 or earlier. |
136 |
IF ( startTime .NE. 0. ) THEN |
137 |
CALL THE_CORRECTION_STEP(startTime, nIter0, myThid) |
138 |
ENDIF |
139 |
ENDIF |
140 |
#endif |
141 |
|
142 |
C-- Initial conditions are convectively adjusted (for historical reasons) |
143 |
IF ( startTime .EQ. 0. .AND. cAdjFreq .NE. 0. ) THEN |
144 |
CADJ loop = parallel |
145 |
DO bj = myByLo(myThid), myByHi(myThid) |
146 |
CADJ loop = parallel |
147 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
148 |
iMin=1-Olx |
149 |
iMax=sNx+Olx |
150 |
jMin=1-Oly |
151 |
jMax=sNy+Oly |
152 |
CALL CONVECTIVE_ADJUSTMENT_INI( |
153 |
I bi, bj, iMin, iMax, jMin, jMax, |
154 |
I startTime, nIter0, myThid ) |
155 |
ENDDO |
156 |
ENDDO |
157 |
_BARRIER |
158 |
END IF |
159 |
|
160 |
#ifdef NONLIN_FRSURF |
161 |
C-- Compute the surface level thickness <-- function of etaH(n) |
162 |
C and modify hFac(C,W,S) accordingly : |
163 |
IF ( select_rStar.NE.0 ) |
164 |
& CALL CALC_R_STAR(etaH, startTime, -1 , myThid ) |
165 |
IF (nonlinFreeSurf.GT.0) THEN |
166 |
IF ( select_rStar.GT.0 ) THEN |
167 |
CALL UPDATE_R_STAR( startTime, nIter0, myThid ) |
168 |
ELSE |
169 |
CALL CALC_SURF_DR(etaH, startTime, -1 , myThid ) |
170 |
CALL UPDATE_SURF_DR( startTime, nIter0, myThid ) |
171 |
ENDIF |
172 |
ENDIF |
173 |
C- update also CG2D matrix (and preconditioner) |
174 |
IF ( nonlinFreeSurf.GT.2) THEN |
175 |
CALL UPDATE_CG2D( startTime, nIter0, myThid ) |
176 |
ENDIF |
177 |
#endif /* NONLIN_FRSURF */ |
178 |
|
179 |
DO bj=myByLo(myThid),myByHi(myThid) |
180 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
181 |
|
182 |
C-- Integrate continuity vertically for vertical velocity |
183 |
CALL INTEGR_CONTINUITY( bi, bj, uVel, vVel, |
184 |
I startTime, nIter0, myThid ) |
185 |
|
186 |
#ifdef ALLOW_TIMEAVE |
187 |
C-- initialise time-average arrays with initial state values |
188 |
IF (taveFreq.GT.0.) THEN |
189 |
CALL TIMEAVE_STATVARS(startTime, nIter0, bi, bj, myThid) |
190 |
ENDIF |
191 |
cswdptr -- add --- |
192 |
#ifdef ALLOW_PTRACERS |
193 |
che3 needs an IF ( usePTRACERS ) THEN |
194 |
CALL PTRACERS_STATVARS(startTime, nIter0, bi, bj, myThid) |
195 |
#endif |
196 |
cswdptr -- end add --- |
197 |
#endif /* ALLOW_TIMEAVE */ |
198 |
|
199 |
C-- end bi,bj loop. |
200 |
ENDDO |
201 |
ENDDO |
202 |
|
203 |
C-- Fill in overlap regions for wVel : |
204 |
_EXCH_XYZ_R8(wVel,myThid) |
205 |
|
206 |
C-- Finally summarise the model state |
207 |
CALL STATE_SUMMARY( myThid ) |
208 |
|
209 |
RETURN |
210 |
END |