1 |
C $Header$ |
C $Header$ |
2 |
|
C $Name$ |
3 |
|
|
4 |
#include "CPP_OPTIONS.h" |
#include "CPP_OPTIONS.h" |
5 |
|
|
6 |
|
CBOP |
7 |
|
C !ROUTINE: DYNAMICS |
8 |
|
C !INTERFACE: |
9 |
SUBROUTINE DYNAMICS(myTime, myIter, myThid) |
SUBROUTINE DYNAMICS(myTime, myIter, myThid) |
10 |
C /==========================================================\ |
C !DESCRIPTION: \bv |
11 |
C | SUBROUTINE DYNAMICS | |
C *==========================================================* |
12 |
C | o Controlling routine for the explicit part of the model | |
C | SUBROUTINE DYNAMICS |
13 |
C | dynamics. | |
C | o Controlling routine for the explicit part of the model |
14 |
C |==========================================================| |
C | dynamics. |
15 |
C | This routine evaluates the "dynamics" terms for each | |
C *==========================================================* |
16 |
C | block of ocean in turn. Because the blocks of ocean have | |
C | This routine evaluates the "dynamics" terms for each |
17 |
C | overlap regions they are independent of one another. | |
C | block of ocean in turn. Because the blocks of ocean have |
18 |
C | If terms involving lateral integrals are needed in this | |
C | overlap regions they are independent of one another. |
19 |
C | routine care will be needed. Similarly finite-difference | |
C | If terms involving lateral integrals are needed in this |
20 |
C | operations with stencils wider than the overlap region | |
C | routine care will be needed. Similarly finite-difference |
21 |
C | require special consideration. | |
C | operations with stencils wider than the overlap region |
22 |
C | Notes | |
C | require special consideration. |
23 |
C | ===== | |
C | The algorithm... |
24 |
C | C*P* comments indicating place holders for which code is | |
C | |
25 |
C | presently being developed. | |
C | "Correction Step" |
26 |
C \==========================================================/ |
C | ================= |
27 |
|
C | Here we update the horizontal velocities with the surface |
28 |
|
C | pressure such that the resulting flow is either consistent |
29 |
|
C | with the free-surface evolution or the rigid-lid: |
30 |
|
C | U[n] = U* + dt x d/dx P |
31 |
|
C | V[n] = V* + dt x d/dy P |
32 |
|
C | |
33 |
|
C | "Calculation of Gs" |
34 |
|
C | =================== |
35 |
|
C | This is where all the accelerations and tendencies (ie. |
36 |
|
C | physics, parameterizations etc...) are calculated |
37 |
|
C | rho = rho ( theta[n], salt[n] ) |
38 |
|
C | b = b(rho, theta) |
39 |
|
C | K31 = K31 ( rho ) |
40 |
|
C | Gu[n] = Gu( u[n], v[n], wVel, b, ... ) |
41 |
|
C | Gv[n] = Gv( u[n], v[n], wVel, b, ... ) |
42 |
|
C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... ) |
43 |
|
C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... ) |
44 |
|
C | |
45 |
|
C | "Time-stepping" or "Prediction" |
46 |
|
C | ================================ |
47 |
|
C | The models variables are stepped forward with the appropriate |
48 |
|
C | time-stepping scheme (currently we use Adams-Bashforth II) |
49 |
|
C | - For momentum, the result is always *only* a "prediction" |
50 |
|
C | in that the flow may be divergent and will be "corrected" |
51 |
|
C | later with a surface pressure gradient. |
52 |
|
C | - Normally for tracers the result is the new field at time |
53 |
|
C | level [n+1} *BUT* in the case of implicit diffusion the result |
54 |
|
C | is also *only* a prediction. |
55 |
|
C | - We denote "predictors" with an asterisk (*). |
56 |
|
C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] ) |
57 |
|
C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] ) |
58 |
|
C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] ) |
59 |
|
C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] ) |
60 |
|
C | With implicit diffusion: |
61 |
|
C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] ) |
62 |
|
C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] ) |
63 |
|
C | (1 + dt * K * d_zz) theta[n] = theta* |
64 |
|
C | (1 + dt * K * d_zz) salt[n] = salt* |
65 |
|
C | |
66 |
|
C *==========================================================* |
67 |
|
C \ev |
68 |
|
C !USES: |
69 |
|
IMPLICIT NONE |
70 |
C == Global variables === |
C == Global variables === |
71 |
#include "SIZE.h" |
#include "SIZE.h" |
72 |
#include "EEPARAMS.h" |
#include "EEPARAMS.h" |
|
#include "CG2D.h" |
|
73 |
#include "PARAMS.h" |
#include "PARAMS.h" |
74 |
#include "DYNVARS.h" |
#include "DYNVARS.h" |
75 |
|
#include "GRID.h" |
76 |
|
#ifdef ALLOW_PASSIVE_TRACER |
77 |
|
#include "TR1.h" |
78 |
|
#endif |
79 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
80 |
|
# include "tamc.h" |
81 |
|
# include "tamc_keys.h" |
82 |
|
# include "FFIELDS.h" |
83 |
|
# include "EOS.h" |
84 |
|
# ifdef ALLOW_KPP |
85 |
|
# include "KPP.h" |
86 |
|
# endif |
87 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
88 |
|
#ifdef ALLOW_TIMEAVE |
89 |
|
#include "TIMEAVE_STATV.h" |
90 |
|
#endif |
91 |
|
|
92 |
|
C !CALLING SEQUENCE: |
93 |
|
C DYNAMICS() |
94 |
|
C | |
95 |
|
C |-- CALC_GRAD_PHI_SURF |
96 |
|
C | |
97 |
|
C |-- CALC_VISCOSITY |
98 |
|
C | |
99 |
|
C |-- CALC_PHI_HYD |
100 |
|
C | |
101 |
|
C |-- STORE_PRESSURE |
102 |
|
C | |
103 |
|
C |-- MOM_FLUXFORM |
104 |
|
C | |
105 |
|
C |-- MOM_VECINV |
106 |
|
C | |
107 |
|
C |-- TIMESTEP |
108 |
|
C | |
109 |
|
C |-- OBCS_APPLY_UV |
110 |
|
C | |
111 |
|
C |-- IMPLDIFF |
112 |
|
C | |
113 |
|
C |-- OBCS_APPLY_UV |
114 |
|
C | |
115 |
|
C |-- CALL TIMEAVE_CUMUL_1T |
116 |
|
C |-- CALL DEBUG_STATS_RL |
117 |
|
|
118 |
|
C !INPUT/OUTPUT PARAMETERS: |
119 |
C == Routine arguments == |
C == Routine arguments == |
120 |
C myTime - Current time in simulation |
C myTime - Current time in simulation |
121 |
C myIter - Current iteration number in simulation |
C myIter - Current iteration number in simulation |
122 |
C myThid - Thread number for this instance of the routine. |
C myThid - Thread number for this instance of the routine. |
|
INTEGER myThid |
|
123 |
_RL myTime |
_RL myTime |
124 |
INTEGER myIter |
INTEGER myIter |
125 |
|
INTEGER myThid |
126 |
|
|
127 |
|
C !LOCAL VARIABLES: |
128 |
C == Local variables |
C == Local variables |
129 |
C xA, yA - Per block temporaries holding face areas |
C fVer[STUV] o fVer: Vertical flux term - note fVer |
|
C uTrans, vTrans, wTrans - Per block temporaries holding flow transport |
|
|
C wVel o uTrans: Zonal transport |
|
|
C o vTrans: Meridional transport |
|
|
C o wTrans: Vertical transport |
|
|
C o wVel: Vertical velocity at upper and lower |
|
|
C cell faces. |
|
|
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 |
|
130 |
C is "pipelined" in the vertical |
C is "pipelined" in the vertical |
131 |
C so we need an fVer for each |
C so we need an fVer for each |
132 |
C variable. |
C variable. |
133 |
C iMin, iMax - Ranges and sub-block indices on which calculations |
C rhoK, rhoKM1 - Density at current level, and level above |
134 |
C jMin, jMax are applied. |
C phiHyd - Hydrostatic part of the potential. |
135 |
|
C In z coords phiHyd is the hydrostatic |
136 |
|
C Potential (=pressure/rho0) anomaly |
137 |
|
C In p coords phiHyd is the geopotential |
138 |
|
C surface height anomaly. |
139 |
|
C dPhiHydX,Y :: Gradient (X & Y directions) of Hydrostatic Potential |
140 |
|
C phiSurfX, - gradient of Surface potential (Pressure/rho, ocean) |
141 |
|
C phiSurfY or geopotential (atmos) in X and Y direction |
142 |
|
C iMin, iMax - Ranges and sub-block indices on which calculations |
143 |
|
C jMin, jMax are applied. |
144 |
C bi, bj |
C bi, bj |
145 |
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 |
146 |
C are switched with layer to be the appropriate index |
C kDown, km1 are switched with layer to be the appropriate |
147 |
C into fVerTerm |
C index into fVerTerm. |
148 |
_RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
149 |
_RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
150 |
_RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
151 |
_RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
152 |
_RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
153 |
_RL wVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
_RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
154 |
_RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
155 |
_RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
156 |
_RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
157 |
_RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
158 |
_RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
|
_RL mTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL pTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
|
|
_RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
|
|
_RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
|
|
_RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
|
|
_RL pH (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz) |
|
|
_RL rhokm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL rhokp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL rhotmp(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL pSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_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) |
|
159 |
|
|
160 |
INTEGER iMin, iMax |
INTEGER iMin, iMax |
161 |
INTEGER jMin, jMax |
INTEGER jMin, jMax |
162 |
INTEGER bi, bj |
INTEGER bi, bj |
163 |
INTEGER i, j |
INTEGER i, j |
164 |
INTEGER k, kM1, kUp, kDown |
INTEGER k, km1, kp1, kup, kDown |
|
LOGICAL BOTTOM_LAYER |
|
165 |
|
|
166 |
|
LOGICAL DIFFERENT_MULTIPLE |
167 |
|
EXTERNAL DIFFERENT_MULTIPLE |
168 |
|
|
169 |
C--- The algorithm... |
C--- The algorithm... |
170 |
C |
C |
171 |
C "Correction Step" |
C "Correction Step" |
180 |
C =================== |
C =================== |
181 |
C This is where all the accelerations and tendencies (ie. |
C This is where all the accelerations and tendencies (ie. |
182 |
C physics, parameterizations etc...) are calculated |
C physics, parameterizations etc...) are calculated |
|
C w = sum_z ( div. u[n] ) |
|
183 |
C rho = rho ( theta[n], salt[n] ) |
C rho = rho ( theta[n], salt[n] ) |
184 |
|
C b = b(rho, theta) |
185 |
C K31 = K31 ( rho ) |
C K31 = K31 ( rho ) |
186 |
C Gu[n] = Gu( u[n], v[n], w, rho, Ph, ... ) |
C Gu[n] = Gu( u[n], v[n], wVel, b, ... ) |
187 |
C Gv[n] = Gv( u[n], v[n], w, rho, Ph, ... ) |
C Gv[n] = Gv( u[n], v[n], wVel, b, ... ) |
188 |
C Gt[n] = Gt( theta[n], u[n], v[n], w, K31, ... ) |
C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... ) |
189 |
C Gs[n] = Gs( salt[n], u[n], v[n], w, K31, ... ) |
C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... ) |
190 |
C |
C |
191 |
C "Time-stepping" or "Prediction" |
C "Time-stepping" or "Prediction" |
192 |
C ================================ |
C ================================ |
209 |
C (1 + dt * K * d_zz) theta[n] = theta* |
C (1 + dt * K * d_zz) theta[n] = theta* |
210 |
C (1 + dt * K * d_zz) salt[n] = salt* |
C (1 + dt * K * d_zz) salt[n] = salt* |
211 |
C--- |
C--- |
212 |
|
CEOP |
213 |
|
|
214 |
C-- Set up work arrays with valid (i.e. not NaN) values |
C-- Set up work arrays with valid (i.e. not NaN) values |
215 |
C These inital values do not alter the numerical results. They |
C These inital values do not alter the numerical results. They |
218 |
C uninitialised but inert locations. |
C uninitialised but inert locations. |
219 |
DO j=1-OLy,sNy+OLy |
DO j=1-OLy,sNy+OLy |
220 |
DO i=1-OLx,sNx+OLx |
DO i=1-OLx,sNx+OLx |
221 |
xA(i,j) = 0. _d 0 |
rhoKM1 (i,j) = 0. _d 0 |
222 |
yA(i,j) = 0. _d 0 |
rhok (i,j) = 0. _d 0 |
223 |
uTrans(i,j) = 0. _d 0 |
phiSurfX(i,j) = 0. _d 0 |
224 |
vTrans(i,j) = 0. _d 0 |
phiSurfY(i,j) = 0. _d 0 |
|
aTerm(i,j) = 0. _d 0 |
|
|
xTerm(i,j) = 0. _d 0 |
|
|
cTerm(i,j) = 0. _d 0 |
|
|
mTerm(i,j) = 0. _d 0 |
|
|
pTerm(i,j) = 0. _d 0 |
|
|
fZon(i,j) = 0. _d 0 |
|
|
fMer(i,j) = 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 |
|
|
ENDDO |
|
|
rhokm1(i,j) = 0. _d 0 |
|
|
rhok (i,j) = 0. _d 0 |
|
|
rhokp1(i,j) = 0. _d 0 |
|
|
rhotmp(i,j) = 0. _d 0 |
|
|
maskC (i,j) = 0. _d 0 |
|
225 |
ENDDO |
ENDDO |
226 |
ENDDO |
ENDDO |
227 |
|
|
228 |
|
C-- Call to routine for calculation of |
229 |
|
C Eliassen-Palm-flux-forced U-tendency, |
230 |
|
C if desired: |
231 |
|
#ifdef INCLUDE_EP_FORCING_CODE |
232 |
|
CALL CALC_EP_FORCING(myThid) |
233 |
|
#endif |
234 |
|
|
235 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
236 |
|
C-- HPF directive to help TAMC |
237 |
|
CHPF$ INDEPENDENT |
238 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
239 |
|
|
240 |
DO bj=myByLo(myThid),myByHi(myThid) |
DO bj=myByLo(myThid),myByHi(myThid) |
241 |
|
|
242 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
243 |
|
C-- HPF directive to help TAMC |
244 |
|
CHPF$ INDEPENDENT, NEW (fVerU,fVerV |
245 |
|
CHPF$& ,phiHyd |
246 |
|
CHPF$& ,KappaRU,KappaRV |
247 |
|
CHPF$& ) |
248 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
249 |
|
|
250 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
DO bi=myBxLo(myThid),myBxHi(myThid) |
251 |
|
|
252 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
253 |
|
act1 = bi - myBxLo(myThid) |
254 |
|
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
255 |
|
act2 = bj - myByLo(myThid) |
256 |
|
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
257 |
|
act3 = myThid - 1 |
258 |
|
max3 = nTx*nTy |
259 |
|
act4 = ikey_dynamics - 1 |
260 |
|
idynkey = (act1 + 1) + act2*max1 |
261 |
|
& + act3*max1*max2 |
262 |
|
& + act4*max1*max2*max3 |
263 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
264 |
|
|
265 |
C-- Set up work arrays that need valid initial values |
C-- Set up work arrays that need valid initial values |
266 |
DO j=1-OLy,sNy+OLy |
DO j=1-OLy,sNy+OLy |
267 |
DO i=1-OLx,sNx+OLx |
DO i=1-OLx,sNx+OLx |
268 |
wTrans(i,j) = 0. _d 0 |
DO k=1,Nr |
269 |
wVel (i,j,1) = 0. _d 0 |
phiHyd(i,j,k) = 0. _d 0 |
270 |
wVel (i,j,2) = 0. _d 0 |
KappaRU(i,j,k) = 0. _d 0 |
271 |
fVerT(i,j,1) = 0. _d 0 |
KappaRV(i,j,k) = 0. _d 0 |
272 |
fVerT(i,j,2) = 0. _d 0 |
ENDDO |
273 |
fVerS(i,j,1) = 0. _d 0 |
fVerU (i,j,1) = 0. _d 0 |
274 |
fVerS(i,j,2) = 0. _d 0 |
fVerU (i,j,2) = 0. _d 0 |
275 |
fVerU(i,j,1) = 0. _d 0 |
fVerV (i,j,1) = 0. _d 0 |
276 |
fVerU(i,j,2) = 0. _d 0 |
fVerV (i,j,2) = 0. _d 0 |
277 |
fVerV(i,j,1) = 0. _d 0 |
dPhiHydX(i,j) = 0. _d 0 |
278 |
fVerV(i,j,2) = 0. _d 0 |
dPhiHydY(i,j) = 0. _d 0 |
|
pH(i,j,1) = 0. _d 0 |
|
|
K13(i,j,1) = 0. _d 0 |
|
|
K23(i,j,1) = 0. _d 0 |
|
|
K33(i,j,1) = 0. _d 0 |
|
|
KapGM(i,j) = GMkbackground |
|
279 |
ENDDO |
ENDDO |
280 |
ENDDO |
ENDDO |
281 |
|
|
282 |
iMin = 1-OLx+1 |
C-- Start computation of dynamics |
283 |
iMax = sNx+OLx |
iMin = 0 |
284 |
jMin = 1-OLy+1 |
iMax = sNx+1 |
285 |
jMax = sNy+OLy |
jMin = 0 |
286 |
|
jMax = sNy+1 |
287 |
K = 1 |
|
288 |
BOTTOM_LAYER = K .EQ. Nz |
#ifdef ALLOW_AUTODIFF_TAMC |
289 |
|
CADJ STORE wvel (:,:,:,bi,bj) = |
290 |
C-- Calculate gradient of surface pressure |
CADJ & comlev1_bibj, key = idynkey, byte = isbyte |
291 |
CALL GRAD_PSURF( |
#endif /* ALLOW_AUTODIFF_TAMC */ |
292 |
I bi,bj,iMin,iMax,jMin,jMax, |
|
293 |
O pSurfX,pSurfY, |
C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP) |
294 |
I myThid) |
C (note: this loop will be replaced by CALL CALC_GRAD_ETA) |
295 |
|
IF (implicSurfPress.NE.1.) THEN |
296 |
C-- Update fields in top level according to tendency terms |
CALL CALC_GRAD_PHI_SURF( |
297 |
CALL CORRECTION_STEP( |
I bi,bj,iMin,iMax,jMin,jMax, |
298 |
I bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myTime,myThid) |
I etaN, |
299 |
|
O phiSurfX,phiSurfY, |
300 |
IF ( .NOT. BOTTOM_LAYER ) THEN |
I myThid ) |
|
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 1st level (below W(1)) reference to level 1 |
|
|
CALL FIND_RHO( |
|
|
I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType, |
|
|
O rhoKm1, |
|
|
I myThid ) |
|
|
|
|
|
IF ( .NOT. BOTTOM_LAYER ) THEN |
|
|
C-- Check static stability with layer below |
|
|
C and mix as needed. |
|
|
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,rhoKm1,rhoKp1, |
|
|
I myTime,myIter,myThid) |
|
|
C-- Recompute density after mixing |
|
|
CALL FIND_RHO( |
|
|
I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType, |
|
|
O rhoKm1, |
|
|
I myThid ) |
|
301 |
ENDIF |
ENDIF |
302 |
|
|
303 |
C-- Integrate hydrostatic balance for pH with BC of pH(z=0)=0 |
#ifdef ALLOW_AUTODIFF_TAMC |
304 |
CALL CALC_PH( |
CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte |
305 |
I bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoKm1, |
CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte |
306 |
U pH, |
#ifdef ALLOW_KPP |
307 |
I myThid ) |
CADJ STORE KPPviscAz (:,:,:,bi,bj) |
308 |
|
CADJ & = comlev1_bibj, key=idynkey, byte=isbyte |
309 |
DO K=2,Nz |
#endif /* ALLOW_KPP */ |
310 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
|
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-- Update fields in layer below according to tendency terms |
|
|
C CALL CORRECTION_STEP( |
|
|
C I bi,bj,iMin,iMax,jMin,jMax,K,pSurfX,pSurfY,myTime,myThid) |
|
|
|
|
|
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 |
|
|
C-- Integrate hydrostatic balance for pH with BC of pH(z=0)=0 |
|
|
CALL CALC_PH( |
|
|
I bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1,rhoK, |
|
|
U pH, |
|
|
I myThid ) |
|
|
C-- Calculate iso-neutral slopes for the GM/Redi parameterisation |
|
|
CALL FIND_RHO( |
|
|
I bi, bj, iMin, iMax, jMin, jMax, K-1, K, eosType, |
|
|
O rhoTmp, |
|
|
I myThid ) |
|
|
CALL CALC_ISOSLOPES( |
|
|
I bi, bj, iMin, iMax, jMin, jMax, K, |
|
|
I rhoKm1, rhoK, rhotmp, |
|
|
O K13, K23, K33, KapGM, |
|
|
I myThid ) |
|
|
DO J=jMin,jMax |
|
|
DO I=iMin,iMax |
|
|
rhoKm1(I,J)=rhoK(I,J) |
|
|
ENDDO |
|
|
ENDDO |
|
|
|
|
|
ENDDO ! K |
|
|
|
|
|
DO K = Nz, 1, -1 |
|
|
kM1 =max(1,k-1) ! Points to level above k (=k-1) |
|
|
kUp =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above |
|
|
kDown=1+MOD(k,2) ! Cycles through 2,1 to point to current layer |
|
|
iMin = 1-OLx+2 |
|
|
iMax = sNx+OLx-1 |
|
|
jMin = 1-OLy+2 |
|
|
jMax = sNy+OLy-1 |
|
|
|
|
|
C-- Get temporary terms used by tendency routines |
|
|
CALL CALC_COMMON_FACTORS ( |
|
|
I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown, |
|
|
O xA,yA,uTrans,vTrans,wTrans,wVel,maskC,maskUp, |
|
|
I myThid) |
|
311 |
|
|
312 |
|
#ifdef INCLUDE_CALC_DIFFUSIVITY_CALL |
313 |
C-- Calculate the total vertical diffusivity |
C-- Calculate the total vertical diffusivity |
314 |
CALL CALC_DIFFUSIVITY( |
DO k=1,Nr |
315 |
I bi,bj,iMin,iMax,jMin,jMax,K, |
CALL CALC_VISCOSITY( |
316 |
I maskC,maskUp,KapGM,K33, |
I bi,bj,iMin,iMax,jMin,jMax,k, |
317 |
O KappaZT,KappaZS, |
O KappaRU,KappaRV, |
318 |
I myThid) |
I myThid) |
319 |
|
ENDDO |
320 |
|
#endif |
321 |
|
|
322 |
C-- Calculate accelerations in the momentum equations |
C-- Start of dynamics loop |
323 |
IF ( momStepping ) THEN |
DO k=1,Nr |
|
CALL CALC_MOM_RHS( |
|
|
I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown, |
|
|
I xA,yA,uTrans,vTrans,wTrans,wVel,maskC, |
|
|
I pH, |
|
|
U aTerm,xTerm,cTerm,mTerm,pTerm, |
|
|
U fZon, fMer, fVerU, fVerV, |
|
|
I myThid) |
|
|
ENDIF |
|
324 |
|
|
325 |
C-- Calculate active tracer tendencies |
C-- km1 Points to level above k (=k-1) |
326 |
IF ( tempStepping ) THEN |
C-- kup Cycles through 1,2 to point to layer above |
327 |
CALL CALC_GT( |
C-- kDown Cycles through 2,1 to point to current layer |
328 |
I bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown, |
|
329 |
I xA,yA,uTrans,vTrans,wTrans,maskUp,maskC, |
km1 = MAX(1,k-1) |
330 |
I K13,K23,KappaZT,KapGM, |
kp1 = MIN(k+1,Nr) |
331 |
U aTerm,xTerm,fZon,fMer,fVerT, |
kup = 1+MOD(k+1,2) |
332 |
I myThid) |
kDown= 1+MOD(k,2) |
333 |
ENDIF |
|
334 |
IF ( saltStepping ) THEN |
#ifdef ALLOW_AUTODIFF_TAMC |
335 |
CALL CALC_GS( |
kkey = (idynkey-1)*Nr + k |
336 |
I bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown, |
CADJ STORE pressure(:,:,k,bi,bj) = comlev1_bibj_k , |
337 |
I xA,yA,uTrans,vTrans,wTrans,maskUp,maskC, |
CADJ & key=kkey , byte=isbyte |
338 |
I K13,K23,KappaZS,KapGM, |
#endif /* ALLOW_AUTODIFF_TAMC */ |
339 |
U aTerm,xTerm,fZon,fMer,fVerS, |
|
340 |
I myThid) |
C-- Integrate hydrostatic balance for phiHyd with BC of |
341 |
|
C phiHyd(z=0)=0 |
342 |
|
C distinguishe between Stagger and Non Stagger time stepping |
343 |
|
IF (staggerTimeStep) THEN |
344 |
|
CALL CALC_PHI_HYD( |
345 |
|
I bi,bj,iMin,iMax,jMin,jMax,k, |
346 |
|
I gT, gS, |
347 |
|
U phiHyd, |
348 |
|
O dPhiHydX, dPhiHydY, |
349 |
|
I myTime, myIter, myThid ) |
350 |
|
ELSE |
351 |
|
CALL CALC_PHI_HYD( |
352 |
|
I bi,bj,iMin,iMax,jMin,jMax,k, |
353 |
|
I theta, salt, |
354 |
|
U phiHyd, |
355 |
|
O dPhiHydX, dPhiHydY, |
356 |
|
I myTime, myIter, myThid ) |
357 |
ENDIF |
ENDIF |
358 |
|
|
359 |
C-- Prediction step (step forward all model variables) |
C calculate pressure from phiHyd and store it on common block |
360 |
CALL TIMESTEP( |
C variable pressure |
361 |
I bi,bj,iMin,iMax,jMin,jMax,K, |
CALL STORE_PRESSURE( bi, bj, k, phiHyd, myThid ) |
362 |
I myThid) |
|
363 |
|
C-- Calculate accelerations in the momentum equations (gU, gV, ...) |
364 |
C-- Diagnose barotropic divergence of predicted fields |
C and step forward storing the result in gUnm1, gVnm1, etc... |
365 |
CALL DIV_G( |
IF ( momStepping ) THEN |
366 |
I bi,bj,iMin,iMax,jMin,jMax,K, |
#ifndef DISABLE_MOM_FLUXFORM |
367 |
I xA,yA, |
IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM( |
368 |
I myThid) |
I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown, |
369 |
|
I phiHyd,dPhiHydX,dPhiHydY,KappaRU,KappaRV, |
370 |
C-- Cumulative diagnostic calculations (ie. time-averaging) |
U fVerU, fVerV, |
371 |
#ifdef ALLOW_DIAGNOSTICS |
I myTime, myIter, myThid) |
372 |
IF (taveFreq.GT.0.) THEN |
#endif |
373 |
CALL DO_TIME_AVERAGES( |
#ifndef DISABLE_MOM_VECINV |
374 |
I myTime, myIter, bi, bj, K, kUp, kDown, |
IF (vectorInvariantMomentum) CALL MOM_VECINV( |
375 |
I K13, K23, wVel, KapGM, |
I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown, |
376 |
I myThid ) |
I dPhiHydX,dPhiHydY,KappaRU,KappaRV, |
377 |
ENDIF |
U fVerU, fVerV, |
378 |
|
I myTime, myIter, myThid) |
379 |
#endif |
#endif |
380 |
|
CALL TIMESTEP( |
381 |
|
I bi,bj,iMin,iMax,jMin,jMax,k, |
382 |
|
I phiHyd, dPhiHydX,dPhiHydY, phiSurfX, phiSurfY, |
383 |
|
I myIter, myThid) |
384 |
|
|
385 |
|
#ifdef ALLOW_OBCS |
386 |
|
C-- Apply open boundary conditions |
387 |
|
IF (useOBCS) THEN |
388 |
|
CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid ) |
389 |
|
END IF |
390 |
|
#endif /* ALLOW_OBCS */ |
391 |
|
|
392 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
393 |
|
#ifdef INCLUDE_CD_CODE |
394 |
|
ELSE |
395 |
|
DO j=1-OLy,sNy+OLy |
396 |
|
DO i=1-OLx,sNx+OLx |
397 |
|
guCD(i,j,k,bi,bj) = 0.0 |
398 |
|
gvCD(i,j,k,bi,bj) = 0.0 |
399 |
|
END DO |
400 |
|
END DO |
401 |
|
#endif /* INCLUDE_CD_CODE */ |
402 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
403 |
|
ENDIF |
404 |
|
|
|
ENDDO ! K |
|
405 |
|
|
406 |
C-- Implicit diffusion |
C-- end of dynamics k loop (1:Nr) |
407 |
IF (implicitDiffusion) THEN |
ENDDO |
408 |
CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax, |
|
409 |
I KappaZT,KappaZS, |
C-- Implicit viscosity |
410 |
I myThid ) |
IF (implicitViscosity.AND.momStepping) THEN |
411 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
412 |
|
CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte |
413 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
414 |
|
CALL IMPLDIFF( |
415 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
416 |
|
I deltaTmom, KappaRU,recip_HFacW, |
417 |
|
U gUNm1, |
418 |
|
I myThid ) |
419 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
420 |
|
CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte |
421 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
422 |
|
CALL IMPLDIFF( |
423 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
424 |
|
I deltaTmom, KappaRV,recip_HFacS, |
425 |
|
U gVNm1, |
426 |
|
I myThid ) |
427 |
|
|
428 |
|
#ifdef ALLOW_OBCS |
429 |
|
C-- Apply open boundary conditions |
430 |
|
IF (useOBCS) THEN |
431 |
|
DO K=1,Nr |
432 |
|
CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid ) |
433 |
|
ENDDO |
434 |
|
END IF |
435 |
|
#endif /* ALLOW_OBCS */ |
436 |
|
|
437 |
|
#ifdef INCLUDE_CD_CODE |
438 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
439 |
|
CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte |
440 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
441 |
|
CALL IMPLDIFF( |
442 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
443 |
|
I deltaTmom, KappaRU,recip_HFacW, |
444 |
|
U vVelD, |
445 |
|
I myThid ) |
446 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
447 |
|
CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte |
448 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
449 |
|
CALL IMPLDIFF( |
450 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
451 |
|
I deltaTmom, KappaRV,recip_HFacS, |
452 |
|
U uVelD, |
453 |
|
I myThid ) |
454 |
|
#endif /* INCLUDE_CD_CODE */ |
455 |
|
C-- End If implicitViscosity.AND.momStepping |
456 |
ENDIF |
ENDIF |
457 |
|
|
458 |
|
C- jmc: add for diagnostic of phiHyd |
459 |
|
IF ( DIFFERENT_MULTIPLE(diagFreq,myTime+deltaTClock,myTime) |
460 |
|
& .AND. buoyancyRelation .NE. 'OCEANIC' ) THEN |
461 |
|
CALL WRITE_LOCAL_RL('Ph','I10',Nr,phiHyd, |
462 |
|
& bi,bj,1,myIter+1,myThid) |
463 |
|
ENDIF |
464 |
|
|
465 |
|
#ifdef ALLOW_TIMEAVE |
466 |
|
IF (taveFreq.GT.0.) THEN |
467 |
|
CALL TIMEAVE_CUMUL_1T(phiHydtave, phiHyd, Nr, |
468 |
|
I deltaTclock, bi, bj, myThid) |
469 |
|
ENDIF |
470 |
|
#endif /* ALLOW_TIMEAVE */ |
471 |
|
|
472 |
ENDDO |
ENDDO |
473 |
ENDDO |
ENDDO |
474 |
|
|
475 |
C write(0,*) 'dynamics: pS ',minval(cg2d_x(1:sNx,1:sNy,:,:)), |
Cml( |
476 |
C & maxval(cg2d_x(1:sNx,1:sNy,:,:)) |
C In order to compare the variance of phiHydLow of a p/z-coordinate |
477 |
C write(0,*) 'dynamics: U ',minval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.), |
C run with etaH of a z/p-coordinate run the drift of phiHydLow |
478 |
C & maxval(uVel(1:sNx,1:sNy,1,:,:),mask=uVel(1:sNx,1:sNy,1,:,:).NE.0.) |
C has to be removed by something like the following subroutine: |
479 |
C write(0,*) 'dynamics: V ',minval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.), |
C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF, |
480 |
C & maxval(vVel(1:sNx,1:sNy,1,:,:),mask=vVel(1:sNx,1:sNy,1,:,:).NE.0.) |
C & 'phiHydLow', myThid ) |
481 |
C write(0,*) 'dynamics: wVel(1) ', |
Cml) |
482 |
C & minval(wVel(1:sNx,1:sNy,1),mask=wVel(1:sNx,1:sNy,1).NE.0.), |
|
483 |
C & maxval(wVel(1:sNx,1:sNy,1),mask=wVel(1:sNx,1:sNy,1).NE.0.) |
#ifndef DISABLE_DEBUGMODE |
484 |
C write(0,*) 'dynamics: wVel(2) ', |
If (debugMode) THEN |
485 |
C & minval(wVel(1:sNx,1:sNy,2),mask=wVel(1:sNx,1:sNy,2).NE.0.), |
CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid) |
486 |
C & maxval(wVel(1:sNx,1:sNy,2),mask=wVel(1:sNx,1:sNy,2).NE.0.) |
CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid) |
487 |
cblk write(0,*) 'dynamics: K13',minval(K13(1:sNx,1:sNy,:)), |
CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid) |
488 |
cblk & maxval(K13(1:sNx,1:sNy,:)) |
CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid) |
489 |
cblk write(0,*) 'dynamics: K23',minval(K23(1:sNx,1:sNy,:)), |
CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid) |
490 |
cblk & maxval(K23(1:sNx,1:sNy,:)) |
CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid) |
491 |
cblk write(0,*) 'dynamics: K33',minval(K33(1:sNx,1:sNy,:)), |
CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid) |
492 |
cblk & maxval(K33(1:sNx,1:sNy,:)) |
CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid) |
493 |
C write(0,*) 'dynamics: gT ',minval(gT(1:sNx,1:sNy,:,:,:)), |
CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid) |
494 |
C & maxval(gT(1:sNx,1:sNy,:,:,:)) |
CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid) |
495 |
C write(0,*) 'dynamics: T ',minval(Theta(1:sNx,1:sNy,:,:,:)), |
CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid) |
496 |
C & maxval(Theta(1:sNx,1:sNy,:,:,:)) |
CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid) |
497 |
C write(0,*) 'dynamics: gS ',minval(gS(1:sNx,1:sNy,:,:,:)), |
CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid) |
498 |
C & maxval(gS(1:sNx,1:sNy,:,:,:)) |
CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid) |
499 |
C write(0,*) 'dynamics: S ',minval(salt(1:sNx,1:sNy,:,:,:)), |
ENDIF |
500 |
C & maxval(salt(1:sNx,1:sNy,:,:,:)) |
#endif |
|
C write(0,*) 'dynamics: pH ',minval(pH/(Gravity*Rhonil),mask=ph.NE.0.), |
|
|
C & maxval(pH/(Gravity*Rhonil)) |
|
501 |
|
|
502 |
RETURN |
RETURN |
503 |
END |
END |