1 |
C $Header$ |
C $Header$ |
2 |
|
C $Name$ |
3 |
|
|
4 |
#include "CPP_OPTIONS.h" |
#include "CPP_OPTIONS.h" |
5 |
|
|
21 |
C | C*P* comments indicating place holders for which code is | |
C | C*P* comments indicating place holders for which code is | |
22 |
C | presently being developed. | |
C | presently being developed. | |
23 |
C \==========================================================/ |
C \==========================================================/ |
24 |
|
IMPLICIT NONE |
25 |
|
|
26 |
C == Global variables === |
C == Global variables === |
27 |
#include "SIZE.h" |
#include "SIZE.h" |
28 |
#include "EEPARAMS.h" |
#include "EEPARAMS.h" |
|
#include "CG2D.h" |
|
29 |
#include "PARAMS.h" |
#include "PARAMS.h" |
30 |
#include "DYNVARS.h" |
#include "DYNVARS.h" |
31 |
|
#include "GRID.h" |
32 |
|
|
33 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
34 |
|
# include "tamc.h" |
35 |
|
# include "tamc_keys.h" |
36 |
|
# include "FFIELDS.h" |
37 |
|
# ifdef ALLOW_KPP |
38 |
|
# include "KPP.h" |
39 |
|
# endif |
40 |
|
# ifdef ALLOW_GMREDI |
41 |
|
# include "GMREDI.h" |
42 |
|
# endif |
43 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
44 |
|
|
45 |
|
#ifdef ALLOW_TIMEAVE |
46 |
|
#include "TIMEAVE_STATV.h" |
47 |
|
#endif |
48 |
|
|
49 |
C == Routine arguments == |
C == Routine arguments == |
50 |
C myTime - Current time in simulation |
C myTime - Current time in simulation |
51 |
C myIter - Current iteration number in simulation |
C myIter - Current iteration number in simulation |
52 |
C myThid - Thread number for this instance of the routine. |
C myThid - Thread number for this instance of the routine. |
|
INTEGER myThid |
|
53 |
_RL myTime |
_RL myTime |
54 |
INTEGER myIter |
INTEGER myIter |
55 |
|
INTEGER myThid |
56 |
|
|
57 |
C == Local variables |
C == Local variables |
58 |
C xA, yA - Per block temporaries holding face areas |
C xA, yA - Per block temporaries holding face areas |
59 |
C uTrans, vTrans, wTrans - Per block temporaries holding flow transport |
C uTrans, vTrans, rTrans - Per block temporaries holding flow |
60 |
C wVel o uTrans: Zonal transport |
C transport |
61 |
|
C o uTrans: Zonal transport |
62 |
C o vTrans: Meridional transport |
C o vTrans: Meridional transport |
63 |
C o wTrans: Vertical transport |
C o rTrans: Vertical transport |
64 |
C o wVel: Vertical velocity at upper and lower |
C maskUp o maskUp: land/water mask for W points |
65 |
C cell faces. |
C fVer[STUV] o fVer: Vertical flux term - note fVer |
|
C maskC,maskUp o maskC: land/water mask for tracer cells |
|
|
C o maskUp: land/water mask for W points |
|
|
C aTerm, xTerm, cTerm - Work arrays for holding separate terms in |
|
|
C mTerm, pTerm, tendency equations. |
|
|
C fZon, fMer, fVer[STUV] o aTerm: Advection term |
|
|
C o xTerm: Mixing term |
|
|
C o cTerm: Coriolis term |
|
|
C o mTerm: Metric term |
|
|
C o pTerm: Pressure term |
|
|
C o fZon: Zonal flux term |
|
|
C o fMer: Meridional flux term |
|
|
C o fVer: Vertical flux term - note fVer |
|
66 |
C is "pipelined" in the vertical |
C is "pipelined" in the vertical |
67 |
C so we need an fVer for each |
C so we need an fVer for each |
68 |
C variable. |
C variable. |
69 |
C rhoK, rhoKM1 - Density at current level, level above and level below. |
C rhoK, rhoKM1 - Density at current level, and level above |
70 |
C rhoKP1 |
C phiHyd - Hydrostatic part of the potential phiHydi. |
71 |
C buoyK, buoyKM1 - Buoyancy at current level and level above. |
C In z coords phiHydiHyd is the hydrostatic |
72 |
C phiHyd - Hydrostatic part of the potential phi. |
C Potential (=pressure/rho0) anomaly |
73 |
C In z coords phiHyd is the hydrostatic pressure anomaly |
C In p coords phiHydiHyd is the geopotential |
74 |
C In p coords phiHyd is the geopotential surface height anomaly. |
C surface height anomaly. |
75 |
C iMin, iMax - Ranges and sub-block indices on which calculations |
C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean) |
76 |
C jMin, jMax are applied. |
C phiSurfY or geopotentiel (atmos) in X and Y direction |
77 |
|
C KappaRT, - Total diffusion in vertical for T and S. |
78 |
|
C KappaRS (background + spatially varying, isopycnal term). |
79 |
|
C iMin, iMax - Ranges and sub-block indices on which calculations |
80 |
|
C jMin, jMax are applied. |
81 |
C bi, bj |
C bi, bj |
82 |
C k, kUp, kDown, kM1 - Index for layer above and below. kUp and kDown |
C k, kup, - Index for layer above and below. kup and kDown |
83 |
C are switched with layer to be the appropriate index |
C kDown, km1 are switched with layer to be the appropriate |
84 |
C into fVerTerm |
C index into fVerTerm. |
85 |
_RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
C tauAB - Adams-Bashforth timestepping weight: 0=forward ; 1/2=Adams-Bashf. |
86 |
_RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
87 |
_RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
88 |
_RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
89 |
_RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
90 |
_RL wVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
_RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
91 |
_RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
92 |
_RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
93 |
_RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
94 |
_RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
95 |
_RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
96 |
_RL mTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
97 |
_RL pTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
98 |
_RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
99 |
_RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
100 |
_RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
_RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
101 |
_RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
_RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
102 |
_RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
_RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
103 |
_RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
_RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
104 |
_RL phiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz) |
_RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
105 |
_RL rhokm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
106 |
_RL rhokp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
107 |
_RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
108 |
_RL buoyKM1(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL tauAB |
109 |
_RL buoyK (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
110 |
_RL rhotmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
C This is currently used by IVDC and Diagnostics |
111 |
_RL pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
|
_RL pSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL K13 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz) |
|
|
_RL K23 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz) |
|
|
_RL K33 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz) |
|
|
_RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL KappaZT(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz) |
|
|
_RL KappaZS(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz) |
|
112 |
|
|
113 |
INTEGER iMin, iMax |
INTEGER iMin, iMax |
114 |
INTEGER jMin, jMax |
INTEGER jMin, jMax |
115 |
INTEGER bi, bj |
INTEGER bi, bj |
116 |
INTEGER i, j |
INTEGER i, j |
117 |
INTEGER k, kM1, kUp, kDown |
INTEGER k, km1, kup, kDown |
|
LOGICAL BOTTOM_LAYER |
|
118 |
|
|
119 |
|
Cjmc : add for phiHyd output <- but not working if multi tile per CPU |
120 |
|
c CHARACTER*(MAX_LEN_MBUF) suff |
121 |
|
c LOGICAL DIFFERENT_MULTIPLE |
122 |
|
c EXTERNAL DIFFERENT_MULTIPLE |
123 |
|
Cjmc(end) |
124 |
|
|
125 |
C--- The algorithm... |
C--- The algorithm... |
126 |
C |
C |
127 |
C "Correction Step" |
C "Correction Step" |
136 |
C =================== |
C =================== |
137 |
C This is where all the accelerations and tendencies (ie. |
C This is where all the accelerations and tendencies (ie. |
138 |
C physics, parameterizations etc...) are calculated |
C physics, parameterizations etc...) are calculated |
|
C w = sum_z ( div. u[n] ) |
|
139 |
C rho = rho ( theta[n], salt[n] ) |
C rho = rho ( theta[n], salt[n] ) |
140 |
|
C b = b(rho, theta) |
141 |
C K31 = K31 ( rho ) |
C K31 = K31 ( rho ) |
142 |
C Gu[n] = Gu( u[n], v[n], w, rho, Ph, ... ) |
C Gu[n] = Gu( u[n], v[n], wVel, b, ... ) |
143 |
C Gv[n] = Gv( u[n], v[n], w, rho, Ph, ... ) |
C Gv[n] = Gv( u[n], v[n], wVel, b, ... ) |
144 |
C Gt[n] = Gt( theta[n], u[n], v[n], w, K31, ... ) |
C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... ) |
145 |
C Gs[n] = Gs( salt[n], u[n], v[n], w, K31, ... ) |
C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... ) |
146 |
C |
C |
147 |
C "Time-stepping" or "Prediction" |
C "Time-stepping" or "Prediction" |
148 |
C ================================ |
C ================================ |
166 |
C (1 + dt * K * d_zz) salt[n] = salt* |
C (1 + dt * K * d_zz) salt[n] = salt* |
167 |
C--- |
C--- |
168 |
|
|
169 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
170 |
|
C-- dummy statement to end declaration part |
171 |
|
ikey = 1 |
172 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
173 |
|
|
174 |
C-- Set up work arrays with valid (i.e. not NaN) values |
C-- Set up work arrays with valid (i.e. not NaN) values |
175 |
C These inital values do not alter the numerical results. They |
C These inital values do not alter the numerical results. They |
176 |
C just ensure that all memory references are to valid floating |
C just ensure that all memory references are to valid floating |
182 |
yA(i,j) = 0. _d 0 |
yA(i,j) = 0. _d 0 |
183 |
uTrans(i,j) = 0. _d 0 |
uTrans(i,j) = 0. _d 0 |
184 |
vTrans(i,j) = 0. _d 0 |
vTrans(i,j) = 0. _d 0 |
185 |
aTerm(i,j) = 0. _d 0 |
DO k=1,Nr |
186 |
xTerm(i,j) = 0. _d 0 |
phiHyd(i,j,k) = 0. _d 0 |
187 |
cTerm(i,j) = 0. _d 0 |
KappaRU(i,j,k) = 0. _d 0 |
188 |
mTerm(i,j) = 0. _d 0 |
KappaRV(i,j,k) = 0. _d 0 |
189 |
pTerm(i,j) = 0. _d 0 |
sigmaX(i,j,k) = 0. _d 0 |
190 |
fZon(i,j) = 0. _d 0 |
sigmaY(i,j,k) = 0. _d 0 |
191 |
fMer(i,j) = 0. _d 0 |
sigmaR(i,j,k) = 0. _d 0 |
|
DO K=1,nZ |
|
|
pH (i,j,k) = 0. _d 0 |
|
|
K13(i,j,k) = 0. _d 0 |
|
|
K23(i,j,k) = 0. _d 0 |
|
|
K33(i,j,k) = 0. _d 0 |
|
|
KappaZT(i,j,k) = 0. _d 0 |
|
192 |
ENDDO |
ENDDO |
193 |
rhokm1(i,j) = 0. _d 0 |
rhoKM1 (i,j) = 0. _d 0 |
194 |
rhok (i,j) = 0. _d 0 |
rhok (i,j) = 0. _d 0 |
195 |
rhokp1(i,j) = 0. _d 0 |
phiSurfX(i,j) = 0. _d 0 |
196 |
rhotmp(i,j) = 0. _d 0 |
phiSurfY(i,j) = 0. _d 0 |
|
buoyKM1(i,j) = 0. _d 0 |
|
|
buoyK (i,j) = 0. _d 0 |
|
|
maskC (i,j) = 0. _d 0 |
|
197 |
ENDDO |
ENDDO |
198 |
ENDDO |
ENDDO |
199 |
|
|
200 |
|
|
201 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
202 |
|
C-- HPF directive to help TAMC |
203 |
|
CHPF$ INDEPENDENT |
204 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
205 |
|
|
206 |
DO bj=myByLo(myThid),myByHi(myThid) |
DO bj=myByLo(myThid),myByHi(myThid) |
207 |
|
|
208 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
209 |
|
C-- HPF directive to help TAMC |
210 |
|
CHPF$ INDEPENDENT, NEW (rTrans,fVerT,fVerS,fVerU,fVerV |
211 |
|
CHPF$& ,phiHyd,utrans,vtrans,xA,yA |
212 |
|
CHPF$& ,KappaRT,KappaRS,KappaRU,KappaRV |
213 |
|
CHPF$& ) |
214 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
215 |
|
|
216 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
DO bi=myBxLo(myThid),myBxHi(myThid) |
217 |
|
|
218 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
219 |
|
act1 = bi - myBxLo(myThid) |
220 |
|
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
221 |
|
|
222 |
|
act2 = bj - myByLo(myThid) |
223 |
|
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
224 |
|
|
225 |
|
act3 = myThid - 1 |
226 |
|
max3 = nTx*nTy |
227 |
|
|
228 |
|
act4 = ikey_dynamics - 1 |
229 |
|
|
230 |
|
ikey = (act1 + 1) + act2*max1 |
231 |
|
& + act3*max1*max2 |
232 |
|
& + act4*max1*max2*max3 |
233 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
234 |
|
|
235 |
C-- Set up work arrays that need valid initial values |
C-- Set up work arrays that need valid initial values |
236 |
DO j=1-OLy,sNy+OLy |
DO j=1-OLy,sNy+OLy |
237 |
DO i=1-OLx,sNx+OLx |
DO i=1-OLx,sNx+OLx |
238 |
wTrans(i,j) = 0. _d 0 |
rTrans(i,j) = 0. _d 0 |
239 |
wVel (i,j,1) = 0. _d 0 |
fVerT (i,j,1) = 0. _d 0 |
240 |
wVel (i,j,2) = 0. _d 0 |
fVerT (i,j,2) = 0. _d 0 |
241 |
fVerT(i,j,1) = 0. _d 0 |
fVerS (i,j,1) = 0. _d 0 |
242 |
fVerT(i,j,2) = 0. _d 0 |
fVerS (i,j,2) = 0. _d 0 |
243 |
fVerS(i,j,1) = 0. _d 0 |
fVerU (i,j,1) = 0. _d 0 |
244 |
fVerS(i,j,2) = 0. _d 0 |
fVerU (i,j,2) = 0. _d 0 |
245 |
fVerU(i,j,1) = 0. _d 0 |
fVerV (i,j,1) = 0. _d 0 |
246 |
fVerU(i,j,2) = 0. _d 0 |
fVerV (i,j,2) = 0. _d 0 |
247 |
fVerV(i,j,1) = 0. _d 0 |
ENDDO |
248 |
fVerV(i,j,2) = 0. _d 0 |
ENDDO |
249 |
pH(i,j,1) = 0. _d 0 |
|
250 |
K13(i,j,1) = 0. _d 0 |
DO k=1,Nr |
251 |
K23(i,j,1) = 0. _d 0 |
DO j=1-OLy,sNy+OLy |
252 |
K33(i,j,1) = 0. _d 0 |
DO i=1-OLx,sNx+OLx |
253 |
KapGM(i,j) = GMkbackground |
C This is currently also used by IVDC and Diagnostics |
254 |
|
ConvectCount(i,j,k) = 0. |
255 |
|
KappaRT(i,j,k) = 0. _d 0 |
256 |
|
KappaRS(i,j,k) = 0. _d 0 |
257 |
|
ENDDO |
258 |
ENDDO |
ENDDO |
259 |
ENDDO |
ENDDO |
260 |
|
|
263 |
jMin = 1-OLy+1 |
jMin = 1-OLy+1 |
264 |
jMax = sNy+OLy |
jMax = sNy+OLy |
265 |
|
|
|
K = 1 |
|
|
BOTTOM_LAYER = K .EQ. Nz |
|
266 |
|
|
267 |
C-- Calculate gradient of surface pressure |
#ifdef ALLOW_AUTODIFF_TAMC |
268 |
CALL GRAD_PSURF( |
CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte |
269 |
I bi,bj,iMin,iMax,jMin,jMax, |
CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte |
270 |
O pSurfX,pSurfY, |
CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte |
271 |
I myThid) |
CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte |
272 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
273 |
C-- Update fields in top level according to tendency terms |
|
274 |
CALL CORRECTION_STEP( |
C-- Start of diagnostic loop |
275 |
I bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myTime,myThid) |
DO k=Nr,1,-1 |
276 |
|
|
277 |
IF ( .NOT. BOTTOM_LAYER ) THEN |
#ifdef ALLOW_AUTODIFF_TAMC |
278 |
C-- Update fields in layer below according to tendency terms |
C? Patrick, is this formula correct now that we change the loop range? |
279 |
CALL CORRECTION_STEP( |
C? Do we still need this? |
280 |
I bi,bj,iMin,iMax,jMin,jMax,K+1,pSurfX,pSurfY,myTime,myThid) |
cph kkey formula corrected. |
281 |
ENDIF |
cph Needed for rhok, rhokm1, in the case useGMREDI. |
282 |
C-- Density of 1st level (below W(1)) reference to level 1 |
kkey = (ikey-1)*Nr + k |
283 |
CALL FIND_RHO( |
CADJ STORE rhokm1(:,:) = comlev1_bibj_k , key=kkey, byte=isbyte |
284 |
I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType, |
CADJ STORE rhok (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte |
285 |
O rhoKm1, |
#endif /* ALLOW_AUTODIFF_TAMC */ |
286 |
I myThid ) |
|
287 |
|
C-- Integrate continuity vertically for vertical velocity |
288 |
IF ( .NOT. BOTTOM_LAYER ) THEN |
CALL INTEGRATE_FOR_W( |
289 |
|
I bi, bj, k, uVel, vVel, |
290 |
C-- Check static stability with layer below |
O wVel, |
291 |
C and mix as needed. |
I myThid ) |
292 |
CALL FIND_RHO( |
|
293 |
I bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType, |
#ifdef ALLOW_OBCS |
294 |
O rhoKp1, |
#ifdef ALLOW_NONHYDROSTATIC |
295 |
I myThid ) |
C-- Apply OBC to W if in N-H mode |
296 |
CALL CONVECT( |
IF (useOBCS.AND.nonHydrostatic) THEN |
297 |
I bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1, |
CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid ) |
298 |
I myTime,myIter,myThid) |
ENDIF |
299 |
C-- Recompute density after mixing |
#endif /* ALLOW_NONHYDROSTATIC */ |
300 |
CALL FIND_RHO( |
#endif /* ALLOW_OBCS */ |
301 |
I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType, |
|
302 |
O rhoKm1, |
C-- Calculate gradients of potential density for isoneutral |
303 |
I myThid ) |
C slope terms (e.g. GM/Redi tensor or IVDC diffusivity) |
304 |
ENDIF |
c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN |
305 |
|
IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN |
306 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
307 |
|
CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
308 |
|
CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
309 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
310 |
|
CALL FIND_RHO( |
311 |
|
I bi, bj, iMin, iMax, jMin, jMax, k, k, eosType, |
312 |
|
I theta, salt, |
313 |
|
O rhoK, |
314 |
|
I myThid ) |
315 |
|
IF (k.GT.1) THEN |
316 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
317 |
|
CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
318 |
|
CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
319 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
320 |
|
CALL FIND_RHO( |
321 |
|
I bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType, |
322 |
|
I theta, salt, |
323 |
|
O rhoKm1, |
324 |
|
I myThid ) |
325 |
|
ENDIF |
326 |
|
CALL GRAD_SIGMA( |
327 |
|
I bi, bj, iMin, iMax, jMin, jMax, k, |
328 |
|
I rhoK, rhoKm1, rhoK, |
329 |
|
O sigmaX, sigmaY, sigmaR, |
330 |
|
I myThid ) |
331 |
|
ENDIF |
332 |
|
|
333 |
C-- Calculate buoyancy |
C-- Implicit Vertical Diffusion for Convection |
334 |
CALL CALC_BUOY( |
c ==> should use sigmaR !!! |
335 |
I bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1, |
IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN |
336 |
O buoyKm1, |
CALL CALC_IVDC( |
337 |
I myThid ) |
I bi, bj, iMin, iMax, jMin, jMax, k, |
338 |
|
I rhoKm1, rhoK, |
339 |
C-- Integrate hydrostatic balance for pH with BC of pH(z=0)=0 |
U ConvectCount, KappaRT, KappaRS, |
340 |
CALL CALC_PHI_HYD( |
I myTime, myIter, myThid) |
341 |
I bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyKm1, |
ENDIF |
|
U phiHyd, |
|
|
I myThid ) |
|
|
|
|
|
DO K=2,Nz |
|
|
|
|
|
BOTTOM_LAYER = K .EQ. Nz |
|
|
IF ( .NOT. BOTTOM_LAYER ) THEN |
|
|
C-- Update fields in layer below according to tendency terms |
|
|
CALL CORRECTION_STEP( |
|
|
I bi,bj,iMin,iMax,jMin,jMax,K+1,pSurfX,pSurfY,myTime,myThid) |
|
|
ENDIF |
|
|
C-- Density of K level (below W(K)) reference to K level |
|
|
CALL FIND_RHO( |
|
|
I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType, |
|
|
O rhoK, |
|
|
I myThid ) |
|
|
IF ( .NOT. BOTTOM_LAYER ) THEN |
|
|
C-- Check static stability with layer below |
|
|
C and mix as needed. |
|
|
C-- Density of K+1 level (below W(K+1)) reference to K level |
|
|
CALL FIND_RHO( |
|
|
I bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType, |
|
|
O rhoKp1, |
|
|
I myThid ) |
|
|
CALL CONVECT( |
|
|
I bi,bj,iMin,iMax,jMin,jMax,K+1,rhoK,rhoKp1, |
|
|
I myTime,myIter,myThid) |
|
|
C-- Recompute density after mixing |
|
|
CALL FIND_RHO( |
|
|
I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType, |
|
|
O rhoK, |
|
|
I myThid ) |
|
|
ENDIF |
|
342 |
|
|
343 |
C-- Calculate buoyancy |
C-- end of diagnostic k loop (Nr:1) |
344 |
CALL CALC_BUOY( |
ENDDO |
345 |
I bi,bj,iMin,iMax,jMin,jMax,K,rhoK, |
|
346 |
O buoyK, |
#ifdef ALLOW_AUTODIFF_TAMC |
347 |
I myThid ) |
cph avoids recomputation of integrate_for_w |
348 |
|
CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte |
349 |
C-- Integrate hydrostatic balance for pH with BC of pH(z=0)=0 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
350 |
CALL CALC_PHI_HYD( |
|
351 |
I bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyK, |
#ifdef ALLOW_OBCS |
352 |
U phiHyd, |
C-- Calculate future values on open boundaries |
353 |
I myThid ) |
IF (useOBCS) THEN |
354 |
C-- Calculate iso-neutral slopes for the GM/Redi parameterisation |
CALL OBCS_CALC( bi, bj, myTime+deltaT, |
355 |
CALL FIND_RHO( |
I uVel, vVel, wVel, theta, salt, |
356 |
I bi, bj, iMin, iMax, jMin, jMax, K-1, K, eosType, |
I myThid ) |
357 |
O rhoTmp, |
ENDIF |
358 |
I myThid ) |
#endif /* ALLOW_OBCS */ |
359 |
CALL CALC_ISOSLOPES( |
|
360 |
I bi, bj, iMin, iMax, jMin, jMax, K, |
C-- Determines forcing terms based on external fields |
361 |
I rhoKm1, rhoK, rhotmp, |
C relaxation terms, etc. |
362 |
O K13, K23, K33, KapGM, |
CALL EXTERNAL_FORCING_SURF( |
363 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
364 |
|
I myThid ) |
365 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
366 |
|
cph needed for KPP |
367 |
|
CADJ STORE surfacetendencyU(:,:,bi,bj) |
368 |
|
CADJ & = comlev1_bibj, key=ikey, byte=isbyte |
369 |
|
CADJ STORE surfacetendencyV(:,:,bi,bj) |
370 |
|
CADJ & = comlev1_bibj, key=ikey, byte=isbyte |
371 |
|
CADJ STORE surfacetendencyS(:,:,bi,bj) |
372 |
|
CADJ & = comlev1_bibj, key=ikey, byte=isbyte |
373 |
|
CADJ STORE surfacetendencyT(:,:,bi,bj) |
374 |
|
CADJ & = comlev1_bibj, key=ikey, byte=isbyte |
375 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
376 |
|
|
377 |
|
#ifdef ALLOW_GMREDI |
378 |
|
|
379 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
380 |
|
CADJ STORE sigmaX(:,:,:) = comlev1, key=ikey, byte=isbyte |
381 |
|
CADJ STORE sigmaY(:,:,:) = comlev1, key=ikey, byte=isbyte |
382 |
|
CADJ STORE sigmaR(:,:,:) = comlev1, key=ikey, byte=isbyte |
383 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
384 |
|
C-- Calculate iso-neutral slopes for the GM/Redi parameterisation |
385 |
|
IF (useGMRedi) THEN |
386 |
|
DO k=1,Nr |
387 |
|
CALL GMREDI_CALC_TENSOR( |
388 |
|
I bi, bj, iMin, iMax, jMin, jMax, k, |
389 |
|
I sigmaX, sigmaY, sigmaR, |
390 |
I myThid ) |
I myThid ) |
|
DO J=jMin,jMax |
|
|
DO I=iMin,iMax |
|
|
rhoKm1(I,J) =rhoK(I,J) |
|
|
buoyKm1(I,J)=buoyK(I,J) |
|
391 |
ENDDO |
ENDDO |
392 |
ENDDO |
#ifdef ALLOW_AUTODIFF_TAMC |
393 |
|
ELSE |
394 |
|
DO k=1, Nr |
395 |
|
CALL GMREDI_CALC_TENSOR_DUMMY( |
396 |
|
I bi, bj, iMin, iMax, jMin, jMax, k, |
397 |
|
I sigmaX, sigmaY, sigmaR, |
398 |
|
I myThid ) |
399 |
|
ENDDO |
400 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
401 |
|
ENDIF |
402 |
|
|
403 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
404 |
|
CADJ STORE Kwx(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
405 |
|
CADJ STORE Kwy(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
406 |
|
CADJ STORE Kwz(:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
407 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
408 |
|
|
409 |
|
#endif /* ALLOW_GMREDI */ |
410 |
|
|
411 |
|
#ifdef ALLOW_KPP |
412 |
|
C-- Compute KPP mixing coefficients |
413 |
|
IF (useKPP) THEN |
414 |
|
CALL KPP_CALC( |
415 |
|
I bi, bj, myTime, myThid ) |
416 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
417 |
|
ELSE |
418 |
|
CALL KPP_CALC_DUMMY( |
419 |
|
I bi, bj, myTime, myThid ) |
420 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
421 |
|
ENDIF |
422 |
|
|
423 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
424 |
|
CADJ STORE KPPghat (:,:,:,bi,bj) |
425 |
|
CADJ & , KPPviscAz (:,:,:,bi,bj) |
426 |
|
CADJ & , KPPdiffKzT(:,:,:,bi,bj) |
427 |
|
CADJ & , KPPdiffKzS(:,:,:,bi,bj) |
428 |
|
CADJ & , KPPfrac (:,: ,bi,bj) |
429 |
|
CADJ & = comlev1_bibj, key=ikey, byte=isbyte |
430 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
431 |
|
|
432 |
|
#endif /* ALLOW_KPP */ |
433 |
|
|
434 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
435 |
|
CADJ STORE KappaRT(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte |
436 |
|
CADJ STORE KappaRS(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte |
437 |
|
CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte |
438 |
|
CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte |
439 |
|
CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte |
440 |
|
CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte |
441 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
442 |
|
|
443 |
|
#ifdef ALLOW_AIM |
444 |
|
C AIM - atmospheric intermediate model, physics package code. |
445 |
|
C note(jmc) : phiHyd=0 at this point but is not really used in Molteni Physics |
446 |
|
IF ( useAIM ) THEN |
447 |
|
CALL TIMER_START('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid) |
448 |
|
CALL AIM_DO_ATMOS_PHYSICS( phiHyd, myTime, myThid ) |
449 |
|
CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid) |
450 |
|
ENDIF |
451 |
|
#endif /* ALLOW_AIM */ |
452 |
|
|
|
ENDDO ! K |
|
453 |
|
|
454 |
DO K = Nz, 1, -1 |
C-- Start of thermodynamics loop |
455 |
kM1 =max(1,k-1) ! Points to level above k (=k-1) |
DO k=Nr,1,-1 |
456 |
kUp =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above |
#ifdef ALLOW_AUTODIFF_TAMC |
457 |
kDown=1+MOD(k,2) ! Cycles through 2,1 to point to current layer |
C? Patrick Is this formula correct? |
458 |
iMin = 1-OLx+2 |
cph Yes, but I rewrote it. |
459 |
iMax = sNx+OLx-1 |
cph Also, the KappaR? need the index and subscript k! |
460 |
jMin = 1-OLy+2 |
kkey = (ikey-1)*Nr + k |
461 |
jMax = sNy+OLy-1 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
462 |
|
|
463 |
|
C-- km1 Points to level above k (=k-1) |
464 |
|
C-- kup Cycles through 1,2 to point to layer above |
465 |
|
C-- kDown Cycles through 2,1 to point to current layer |
466 |
|
|
467 |
|
km1 = MAX(1,k-1) |
468 |
|
kup = 1+MOD(k+1,2) |
469 |
|
kDown= 1+MOD(k,2) |
470 |
|
|
471 |
|
iMin = 1-OLx+2 |
472 |
|
iMax = sNx+OLx-1 |
473 |
|
jMin = 1-OLy+2 |
474 |
|
jMax = sNy+OLy-1 |
475 |
|
|
476 |
C-- Get temporary terms used by tendency routines |
C-- Get temporary terms used by tendency routines |
477 |
CALL CALC_COMMON_FACTORS ( |
CALL CALC_COMMON_FACTORS ( |
478 |
I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown, |
I bi,bj,iMin,iMax,jMin,jMax,k, |
479 |
O xA,yA,uTrans,vTrans,wTrans,wVel,maskC,maskUp, |
O xA,yA,uTrans,vTrans,rTrans,maskUp, |
480 |
I myThid) |
I myThid) |
481 |
|
|
482 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
483 |
|
CADJ STORE KappaRT(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte |
484 |
|
CADJ STORE KappaRS(:,:,k) = comlev1_bibj_k, key=kkey, byte=isbyte |
485 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
486 |
|
|
487 |
|
#ifdef INCLUDE_CALC_DIFFUSIVITY_CALL |
488 |
C-- Calculate the total vertical diffusivity |
C-- Calculate the total vertical diffusivity |
489 |
CALL CALC_DIFFUSIVITY( |
CALL CALC_DIFFUSIVITY( |
490 |
I bi,bj,iMin,iMax,jMin,jMax,K, |
I bi,bj,iMin,iMax,jMin,jMax,k, |
491 |
I maskC,maskUp,KapGM,K33, |
I maskUp, |
492 |
O KappaZT,KappaZS, |
O KappaRT,KappaRS,KappaRU,KappaRV, |
493 |
I myThid) |
I myThid) |
494 |
|
#endif |
495 |
|
|
496 |
C-- Calculate accelerations in the momentum equations |
C-- Calculate active tracer tendencies (gT,gS,...) |
497 |
IF ( momStepping ) THEN |
C and step forward storing result in gTnm1, gSnm1, etc. |
|
CALL CALC_MOM_RHS( |
|
|
I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown, |
|
|
I xA,yA,uTrans,vTrans,wTrans,wVel,maskC, |
|
|
I phiHyd, |
|
|
U aTerm,xTerm,cTerm,mTerm,pTerm, |
|
|
U fZon, fMer, fVerU, fVerV, |
|
|
I myThid) |
|
|
ENDIF |
|
|
|
|
|
C-- Calculate active tracer tendencies |
|
498 |
IF ( tempStepping ) THEN |
IF ( tempStepping ) THEN |
499 |
CALL CALC_GT( |
CALL CALC_GT( |
500 |
I bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown, |
I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown, |
501 |
I xA,yA,uTrans,vTrans,wTrans,maskUp,maskC, |
I xA,yA,uTrans,vTrans,rTrans,maskUp, |
502 |
I K13,K23,KappaZT,KapGM, |
I KappaRT, |
503 |
U aTerm,xTerm,fZon,fMer,fVerT, |
U fVerT, |
504 |
I myThid) |
I myTime, myThid) |
505 |
|
tauAB = 0.5d0 + abEps |
506 |
|
CALL TIMESTEP_TRACER( |
507 |
|
I bi,bj,iMin,iMax,jMin,jMax,k,tauAB, |
508 |
|
I theta, gT, |
509 |
|
U gTnm1, |
510 |
|
I myIter, myThid) |
511 |
ENDIF |
ENDIF |
512 |
IF ( saltStepping ) THEN |
IF ( saltStepping ) THEN |
513 |
CALL CALC_GS( |
CALL CALC_GS( |
514 |
I bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown, |
I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown, |
515 |
I xA,yA,uTrans,vTrans,wTrans,maskUp,maskC, |
I xA,yA,uTrans,vTrans,rTrans,maskUp, |
516 |
I K13,K23,KappaZS,KapGM, |
I KappaRS, |
517 |
U aTerm,xTerm,fZon,fMer,fVerS, |
U fVerS, |
518 |
I myThid) |
I myTime, myThid) |
519 |
|
tauAB = 0.5d0 + abEps |
520 |
|
CALL TIMESTEP_TRACER( |
521 |
|
I bi,bj,iMin,iMax,jMin,jMax,k,tauAB, |
522 |
|
I salt, gS, |
523 |
|
U gSnm1, |
524 |
|
I myIter, myThid) |
525 |
ENDIF |
ENDIF |
526 |
|
|
527 |
C-- Prediction step (step forward all model variables) |
#ifdef ALLOW_OBCS |
528 |
CALL TIMESTEP( |
C-- Apply open boundary conditions |
529 |
I bi,bj,iMin,iMax,jMin,jMax,K, |
IF (useOBCS) THEN |
530 |
I myThid) |
CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid ) |
531 |
|
END IF |
532 |
C-- Diagnose barotropic divergence of predicted fields |
#endif /* ALLOW_OBCS */ |
533 |
CALL DIV_G( |
|
534 |
I bi,bj,iMin,iMax,jMin,jMax,K, |
C-- Freeze water |
535 |
I xA,yA, |
IF (allowFreezing) THEN |
536 |
I myThid) |
#ifdef ALLOW_AUTODIFF_TAMC |
537 |
|
CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k |
538 |
C-- Cumulative diagnostic calculations (ie. time-averaging) |
CADJ & , key = kkey, byte = isbyte |
539 |
#ifdef ALLOW_DIAGNOSTICS |
#endif /* ALLOW_AUTODIFF_TAMC */ |
540 |
IF (taveFreq.GT.0.) THEN |
CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid ) |
541 |
CALL DO_TIME_AVERAGES( |
END IF |
542 |
I myTime, myIter, bi, bj, K, kUp, kDown, |
|
543 |
I K13, K23, wVel, KapGM, |
C-- end of thermodynamic k loop (Nr:1) |
544 |
I myThid ) |
ENDDO |
|
ENDIF |
|
|
#endif |
|
545 |
|
|
546 |
ENDDO ! K |
|
547 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
548 |
|
C? Patrick? What about this one? |
549 |
|
cph Keys iikey and idkey don't seem to be needed |
550 |
|
cph since storing occurs on different tape for each |
551 |
|
cph impldiff call anyways. |
552 |
|
cph Thus, common block comlev1_impl isn't needed either. |
553 |
|
cph Storing below needed in the case useGMREDI. |
554 |
|
iikey = (ikey-1)*maximpl |
555 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
556 |
|
|
557 |
C-- Implicit diffusion |
C-- Implicit diffusion |
558 |
IF (implicitDiffusion) THEN |
IF (implicitDiffusion) THEN |
559 |
CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax, |
|
560 |
I KappaZT,KappaZS, |
IF (tempStepping) THEN |
561 |
I myThid ) |
#ifdef ALLOW_AUTODIFF_TAMC |
562 |
|
idkey = iikey + 1 |
563 |
|
CADJ STORE gTNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte |
564 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
565 |
|
CALL IMPLDIFF( |
566 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
567 |
|
I deltaTtracer, KappaRT, recip_HFacC, |
568 |
|
U gTNm1, |
569 |
|
I myThid ) |
570 |
|
ENDIF |
571 |
|
|
572 |
|
IF (saltStepping) THEN |
573 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
574 |
|
idkey = iikey + 2 |
575 |
|
CADJ STORE gSNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte |
576 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
577 |
|
CALL IMPLDIFF( |
578 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
579 |
|
I deltaTtracer, KappaRS, recip_HFacC, |
580 |
|
U gSNm1, |
581 |
|
I myThid ) |
582 |
|
ENDIF |
583 |
|
|
584 |
|
#ifdef ALLOW_OBCS |
585 |
|
C-- Apply open boundary conditions |
586 |
|
IF (useOBCS) THEN |
587 |
|
DO K=1,Nr |
588 |
|
CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid ) |
589 |
|
ENDDO |
590 |
|
END IF |
591 |
|
#endif /* ALLOW_OBCS */ |
592 |
|
|
593 |
|
C-- End If implicitDiffusion |
594 |
|
ENDIF |
595 |
|
|
596 |
|
C-- Start computation of dynamics |
597 |
|
iMin = 1-OLx+2 |
598 |
|
iMax = sNx+OLx-1 |
599 |
|
jMin = 1-OLy+2 |
600 |
|
jMax = sNy+OLy-1 |
601 |
|
|
602 |
|
C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP) |
603 |
|
C (note: this loop will be replaced by CALL CALC_GRAD_ETA) |
604 |
|
IF (implicSurfPress.NE.1.) THEN |
605 |
|
CALL CALC_GRAD_PHI_SURF( |
606 |
|
I bi,bj,iMin,iMax,jMin,jMax, |
607 |
|
I etaN, |
608 |
|
O phiSurfX,phiSurfY, |
609 |
|
I myThid ) |
610 |
|
ENDIF |
611 |
|
|
612 |
|
C-- Start of dynamics loop |
613 |
|
DO k=1,Nr |
614 |
|
|
615 |
|
C-- km1 Points to level above k (=k-1) |
616 |
|
C-- kup Cycles through 1,2 to point to layer above |
617 |
|
C-- kDown Cycles through 2,1 to point to current layer |
618 |
|
|
619 |
|
km1 = MAX(1,k-1) |
620 |
|
kup = 1+MOD(k+1,2) |
621 |
|
kDown= 1+MOD(k,2) |
622 |
|
|
623 |
|
C-- Integrate hydrostatic balance for phiHyd with BC of |
624 |
|
C phiHyd(z=0)=0 |
625 |
|
C distinguishe between Stagger and Non Stagger time stepping |
626 |
|
IF (staggerTimeStep) THEN |
627 |
|
CALL CALC_PHI_HYD( |
628 |
|
I bi,bj,iMin,iMax,jMin,jMax,k, |
629 |
|
I gTnm1, gSnm1, |
630 |
|
U phiHyd, |
631 |
|
I myThid ) |
632 |
|
ELSE |
633 |
|
CALL CALC_PHI_HYD( |
634 |
|
I bi,bj,iMin,iMax,jMin,jMax,k, |
635 |
|
I theta, salt, |
636 |
|
U phiHyd, |
637 |
|
I myThid ) |
638 |
|
ENDIF |
639 |
|
|
640 |
|
C-- Calculate accelerations in the momentum equations (gU, gV, ...) |
641 |
|
C and step forward storing the result in gUnm1, gVnm1, etc... |
642 |
|
IF ( momStepping ) THEN |
643 |
|
CALL CALC_MOM_RHS( |
644 |
|
I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown, |
645 |
|
I phiHyd,KappaRU,KappaRV, |
646 |
|
U fVerU, fVerV, |
647 |
|
I myTime, myThid) |
648 |
|
CALL TIMESTEP( |
649 |
|
I bi,bj,iMin,iMax,jMin,jMax,k, |
650 |
|
I phiHyd, phiSurfX, phiSurfY, |
651 |
|
I myIter, myThid) |
652 |
|
|
653 |
|
#ifdef ALLOW_OBCS |
654 |
|
C-- Apply open boundary conditions |
655 |
|
IF (useOBCS) THEN |
656 |
|
CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid ) |
657 |
|
END IF |
658 |
|
#endif /* ALLOW_OBCS */ |
659 |
|
|
660 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
661 |
|
#ifdef INCLUDE_CD_CODE |
662 |
|
ELSE |
663 |
|
DO j=1-OLy,sNy+OLy |
664 |
|
DO i=1-OLx,sNx+OLx |
665 |
|
guCD(i,j,k,bi,bj) = 0.0 |
666 |
|
gvCD(i,j,k,bi,bj) = 0.0 |
667 |
|
END DO |
668 |
|
END DO |
669 |
|
#endif /* INCLUDE_CD_CODE */ |
670 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
671 |
|
ENDIF |
672 |
|
|
673 |
|
|
674 |
|
C-- end of dynamics k loop (1:Nr) |
675 |
|
ENDDO |
676 |
|
|
677 |
|
|
678 |
|
|
679 |
|
C-- Implicit viscosity |
680 |
|
IF (implicitViscosity.AND.momStepping) THEN |
681 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
682 |
|
idkey = iikey + 3 |
683 |
|
CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte |
684 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
685 |
|
CALL IMPLDIFF( |
686 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
687 |
|
I deltaTmom, KappaRU,recip_HFacW, |
688 |
|
U gUNm1, |
689 |
|
I myThid ) |
690 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
691 |
|
idkey = iikey + 4 |
692 |
|
CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte |
693 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
694 |
|
CALL IMPLDIFF( |
695 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
696 |
|
I deltaTmom, KappaRV,recip_HFacS, |
697 |
|
U gVNm1, |
698 |
|
I myThid ) |
699 |
|
|
700 |
|
#ifdef ALLOW_OBCS |
701 |
|
C-- Apply open boundary conditions |
702 |
|
IF (useOBCS) THEN |
703 |
|
DO K=1,Nr |
704 |
|
CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid ) |
705 |
|
ENDDO |
706 |
|
END IF |
707 |
|
#endif /* ALLOW_OBCS */ |
708 |
|
|
709 |
|
#ifdef INCLUDE_CD_CODE |
710 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
711 |
|
idkey = iikey + 5 |
712 |
|
CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte |
713 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
714 |
|
CALL IMPLDIFF( |
715 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
716 |
|
I deltaTmom, KappaRU,recip_HFacW, |
717 |
|
U vVelD, |
718 |
|
I myThid ) |
719 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
720 |
|
idkey = iikey + 6 |
721 |
|
CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte |
722 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
723 |
|
CALL IMPLDIFF( |
724 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
725 |
|
I deltaTmom, KappaRV,recip_HFacS, |
726 |
|
U uVelD, |
727 |
|
I myThid ) |
728 |
|
#endif /* INCLUDE_CD_CODE */ |
729 |
|
C-- End If implicitViscosity.AND.momStepping |
730 |
ENDIF |
ENDIF |
731 |
|
|
732 |
|
Cjmc : add for phiHyd output <- but not working if multi tile per CPU |
733 |
|
c IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime) |
734 |
|
c & .AND. buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN |
735 |
|
c WRITE(suff,'(I10.10)') myIter+1 |
736 |
|
c CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid) |
737 |
|
c ENDIF |
738 |
|
Cjmc(end) |
739 |
|
|
740 |
|
#ifdef ALLOW_TIMEAVE |
741 |
|
IF (taveFreq.GT.0.) THEN |
742 |
|
CALL TIMEAVE_CUMUL_1T(phiHydtave, phiHyd, Nr, |
743 |
|
I deltaTclock, bi, bj, myThid) |
744 |
|
IF (ivdc_kappa.NE.0.) THEN |
745 |
|
CALL TIMEAVE_CUMULATE(ConvectCountTave, ConvectCount, Nr, |
746 |
|
I deltaTclock, bi, bj, myThid) |
747 |
|
ENDIF |
748 |
|
ENDIF |
749 |
|
#endif /* ALLOW_TIMEAVE */ |
750 |
|
|
751 |
ENDDO |
ENDDO |
752 |
ENDDO |
ENDDO |
753 |
|
|
|
C write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)), |
|
|
C & maxval(cg2d_x(1:sNx,1:sNy,:,:)) |
|
|
C write(0,*) 'dynamics: U ',minval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.), |
|
|
C & maxval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.) |
|
|
C write(0,*) 'dynamics: V ',minval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.), |
|
|
C & maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.) |
|
|
C write(0,*) 'dynamics: wVel(1) ', |
|
|
C & minval(wVel(1:sNx,1:sNy,1),mask=wVel(1:sNx,1:sNy,1).NE.0.), |
|
|
C & maxval(wVel(1:sNx,1:sNy,1),mask=wVel(1:sNx,1:sNy,1).NE.0.) |
|
|
C write(0,*) 'dynamics: wVel(2) ', |
|
|
C & minval(wVel(1:sNx,1:sNy,2),mask=wVel(1:sNx,1:sNy,2).NE.0.), |
|
|
C & maxval(wVel(1:sNx,1:sNy,2),mask=wVel(1:sNx,1:sNy,2).NE.0.) |
|
|
cblk write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)), |
|
|
cblk & maxval(K13(1:sNx,1:sNy,:)) |
|
|
cblk write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)), |
|
|
cblk & maxval(K23(1:sNx,1:sNy,:)) |
|
|
cblk write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)), |
|
|
cblk & maxval(K33(1:sNx,1:sNy,:)) |
|
|
C write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)), |
|
|
C & maxval(gT(1:sNx,1:sNy,:,:,:)) |
|
|
C write(0,*) 'dynamics: T ',minval(Theta(1:sNx,1:sNy,:,:,:)), |
|
|
C & maxval(Theta(1:sNx,1:sNy,:,:,:)) |
|
|
C write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)), |
|
|
C & maxval(gS(1:sNx,1:sNy,:,:,:)) |
|
|
C write(0,*) 'dynamics: S ',minval(salt(1:sNx,1:sNy,:,:,:)), |
|
|
C & maxval(salt(1:sNx,1:sNy,:,:,:)) |
|
|
C write(0,*) 'dynamics: pH ',minval(pH/(Gravity*Rhonil),mask=ph.NE.0.), |
|
|
C & maxval(pH/(Gravity*Rhonil)) |
|
|
|
|
754 |
RETURN |
RETURN |
755 |
END |
END |