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" |
29 |
#include "CG2D.h" |
#include "CG2D.h" |
30 |
#include "PARAMS.h" |
#include "PARAMS.h" |
31 |
#include "DYNVARS.h" |
#include "DYNVARS.h" |
32 |
|
#include "GRID.h" |
33 |
|
|
34 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
35 |
|
# include "tamc.h" |
36 |
|
# include "tamc_keys.h" |
37 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
38 |
|
|
39 |
|
#ifdef ALLOW_KPP |
40 |
|
# include "KPP.h" |
41 |
|
#endif |
42 |
|
|
43 |
C == Routine arguments == |
C == Routine arguments == |
44 |
C myTime - Current time in simulation |
C myTime - Current time in simulation |
45 |
C myIter - Current iteration number in simulation |
C myIter - Current iteration number in simulation |
46 |
C myThid - Thread number for this instance of the routine. |
C myThid - Thread number for this instance of the routine. |
|
INTEGER myThid |
|
47 |
_RL myTime |
_RL myTime |
48 |
INTEGER myIter |
INTEGER myIter |
49 |
|
INTEGER myThid |
50 |
|
|
51 |
C == Local variables |
C == Local variables |
52 |
C xA, yA - Per block temporaries holding face areas |
C xA, yA - Per block temporaries holding face areas |
53 |
C uTrans, vTrans, rTrans - Per block temporaries holding flow transport |
C uTrans, vTrans, rTrans - Per block temporaries holding flow |
54 |
|
C transport |
55 |
C rVel o uTrans: Zonal transport |
C rVel o uTrans: Zonal transport |
56 |
C o vTrans: Meridional transport |
C o vTrans: Meridional transport |
57 |
C o rTrans: Vertical transport |
C o rTrans: Vertical transport |
58 |
C o rVel: Vertical velocity at upper and lower |
C o rVel: Vertical velocity at upper and |
59 |
C cell faces. |
C lower cell faces. |
60 |
C maskC,maskUp o maskC: land/water mask for tracer cells |
C maskC,maskUp o maskC: land/water mask for tracer cells |
61 |
C o maskUp: land/water mask for W points |
C o maskUp: land/water mask for W points |
62 |
C aTerm, xTerm, cTerm - Work arrays for holding separate terms in |
C fVer[STUV] o fVer: Vertical flux term - note fVer |
|
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 |
|
63 |
C is "pipelined" in the vertical |
C is "pipelined" in the vertical |
64 |
C so we need an fVer for each |
C so we need an fVer for each |
65 |
C variable. |
C variable. |
66 |
C rhoK, rhoKM1 - Density at current level, level above and level below. |
C rhoK, rhoKM1 - Density at current level, and level above |
|
C rhoKP1 |
|
|
C buoyK, buoyKM1 - Buoyancy at current level and level above. |
|
67 |
C phiHyd - Hydrostatic part of the potential phiHydi. |
C phiHyd - Hydrostatic part of the potential phiHydi. |
68 |
C In z coords phiHydiHyd is the hydrostatic pressure anomaly |
C In z coords phiHydiHyd is the hydrostatic |
69 |
C In p coords phiHydiHyd is the geopotential surface height |
C pressure anomaly |
70 |
|
C In p coords phiHydiHyd is the geopotential |
71 |
|
C surface height |
72 |
C anomaly. |
C anomaly. |
73 |
C etaSurfX, - Holds surface elevation gradient in X and Y. |
C etaSurfX, - Holds surface elevation gradient in X and Y. |
74 |
C etaSurfY |
C etaSurfY |
|
C K13, K23, K33 - Non-zero elements of small-angle approximation |
|
|
C diffusion tensor. |
|
|
C KapGM - Spatially varying Visbeck et. al mixing coeff. |
|
75 |
C KappaRT, - Total diffusion in vertical for T and S. |
C KappaRT, - Total diffusion in vertical for T and S. |
76 |
C KappaRS ( background + spatially varying, isopycnal term). |
C KappaRS (background + spatially varying, isopycnal term). |
77 |
C iMin, iMax - Ranges and sub-block indices on which calculations |
C iMin, iMax - Ranges and sub-block indices on which calculations |
78 |
C jMin, jMax are applied. |
C jMin, jMax are applied. |
79 |
C bi, bj |
C bi, bj |
80 |
C k, kUp, - Index for layer above and below. kUp and kDown |
C k, kup, - Index for layer above and below. kup and kDown |
81 |
C kDown, kM1 are switched with layer to be the appropriate index |
C kDown, km1 are switched with layer to be the appropriate |
82 |
C into fVerTerm |
C index into fVerTerm. |
83 |
_RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
84 |
_RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
85 |
_RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
88 |
_RL rVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
_RL rVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
89 |
_RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
90 |
_RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
_RL aTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL xTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL cTerm (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_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) |
|
91 |
_RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
_RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
92 |
_RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
_RL fVerS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
93 |
_RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
_RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
94 |
_RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
_RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
95 |
_RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
_RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
96 |
_RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
_RL rhokp1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
97 |
_RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
_RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
_RL buoyKM1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL buoyK (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL rhotmp (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL etaSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL etaSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
|
_RL K13 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
|
|
_RL K23 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
|
|
_RL K33 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
|
|
_RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
|
98 |
_RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
_RL KappaRT (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
99 |
_RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
_RL KappaRS (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
100 |
|
_RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
101 |
|
_RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
102 |
|
_RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
103 |
|
_RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
104 |
|
_RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
105 |
|
|
106 |
|
C This is currently also used by IVDC and Diagnostics |
107 |
|
C #ifdef INCLUDE_CONVECT_CALL |
108 |
|
_RL ConvectCount (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
109 |
|
C #endif |
110 |
|
|
111 |
INTEGER iMin, iMax |
INTEGER iMin, iMax |
112 |
INTEGER jMin, jMax |
INTEGER jMin, jMax |
113 |
INTEGER bi, bj |
INTEGER bi, bj |
114 |
INTEGER i, j |
INTEGER i, j |
115 |
INTEGER k, kM1, kUp, kDown |
INTEGER k, km1, kup, kDown |
116 |
LOGICAL BOTTOM_LAYER |
|
117 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
118 |
|
INTEGER isbyte |
119 |
|
PARAMETER( isbyte = 4 ) |
120 |
|
|
121 |
|
INTEGER act1, act2, act3, act4 |
122 |
|
INTEGER max1, max2, max3 |
123 |
|
INTEGER iikey, kkey |
124 |
|
INTEGER maximpl |
125 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
126 |
|
|
127 |
C--- The algorithm... |
C--- The algorithm... |
128 |
C |
C |
137 |
C "Calculation of Gs" |
C "Calculation of Gs" |
138 |
C =================== |
C =================== |
139 |
C This is where all the accelerations and tendencies (ie. |
C This is where all the accelerations and tendencies (ie. |
140 |
C phiHydysics, parameterizations etc...) are calculated |
C physics, parameterizations etc...) are calculated |
141 |
C rVel = sum_r ( div. u[n] ) |
C rVel = sum_r ( div. u[n] ) |
142 |
C rho = rho ( theta[n], salt[n] ) |
C rho = rho ( theta[n], salt[n] ) |
143 |
C b = b(rho, theta) |
C b = b(rho, theta) |
169 |
C (1 + dt * K * d_zz) salt[n] = salt* |
C (1 + dt * K * d_zz) salt[n] = salt* |
170 |
C--- |
C--- |
171 |
|
|
172 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
173 |
|
C-- dummy statement to end declaration part |
174 |
|
ikey = 1 |
175 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
176 |
|
|
177 |
C-- Set up work arrays with valid (i.e. not NaN) values |
C-- Set up work arrays with valid (i.e. not NaN) values |
178 |
C These inital values do not alter the numerical results. They |
C These inital values do not alter the numerical results. They |
179 |
C just ensure that all memory references are to valid floating |
C just ensure that all memory references are to valid floating |
185 |
yA(i,j) = 0. _d 0 |
yA(i,j) = 0. _d 0 |
186 |
uTrans(i,j) = 0. _d 0 |
uTrans(i,j) = 0. _d 0 |
187 |
vTrans(i,j) = 0. _d 0 |
vTrans(i,j) = 0. _d 0 |
188 |
aTerm(i,j) = 0. _d 0 |
DO k=1,Nr |
189 |
xTerm(i,j) = 0. _d 0 |
phiHyd(i,j,k) = 0. _d 0 |
190 |
cTerm(i,j) = 0. _d 0 |
KappaRU(i,j,k) = 0. _d 0 |
191 |
mTerm(i,j) = 0. _d 0 |
KappaRV(i,j,k) = 0. _d 0 |
192 |
pTerm(i,j) = 0. _d 0 |
sigmaX(i,j,k) = 0. _d 0 |
193 |
fZon(i,j) = 0. _d 0 |
sigmaY(i,j,k) = 0. _d 0 |
194 |
fMer(i,j) = 0. _d 0 |
sigmaR(i,j,k) = 0. _d 0 |
|
DO K=1,Nr |
|
|
phiHyd (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 |
|
|
KappaRT(i,j,k) = 0. _d 0 |
|
|
KappaRS(i,j,k) = 0. _d 0 |
|
195 |
ENDDO |
ENDDO |
196 |
rhoKM1 (i,j) = 0. _d 0 |
rhoKM1 (i,j) = 0. _d 0 |
197 |
rhok (i,j) = 0. _d 0 |
rhok (i,j) = 0. _d 0 |
|
rhoKP1 (i,j) = 0. _d 0 |
|
|
rhoTMP (i,j) = 0. _d 0 |
|
|
buoyKM1(i,j) = 0. _d 0 |
|
|
buoyK (i,j) = 0. _d 0 |
|
198 |
maskC (i,j) = 0. _d 0 |
maskC (i,j) = 0. _d 0 |
199 |
ENDDO |
ENDDO |
200 |
ENDDO |
ENDDO |
201 |
|
|
202 |
|
|
203 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
204 |
|
C-- HPF directive to help TAMC |
205 |
|
CHPF$ INDEPENDENT |
206 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
207 |
|
|
208 |
DO bj=myByLo(myThid),myByHi(myThid) |
DO bj=myByLo(myThid),myByHi(myThid) |
209 |
|
|
210 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
211 |
|
C-- HPF directive to help TAMC |
212 |
|
CHPF$ INDEPENDENT, NEW (rTrans,rVel,fVerT,fVerS,fVerU,fVerV |
213 |
|
CHPF$& ,phiHyd,utrans,vtrans,maskc,xA,yA |
214 |
|
CHPF$& ,KappaRT,KappaRS,KappaRU,KappaRV |
215 |
|
CHPF$& ) |
216 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
217 |
|
|
218 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
DO bi=myBxLo(myThid),myBxHi(myThid) |
219 |
|
|
220 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
221 |
|
act1 = bi - myBxLo(myThid) |
222 |
|
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
223 |
|
|
224 |
|
act2 = bj - myByLo(myThid) |
225 |
|
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
226 |
|
|
227 |
|
act3 = myThid - 1 |
228 |
|
max3 = nTx*nTy |
229 |
|
|
230 |
|
act4 = ikey_dynamics - 1 |
231 |
|
|
232 |
|
ikey = (act1 + 1) + act2*max1 |
233 |
|
& + act3*max1*max2 |
234 |
|
& + act4*max1*max2*max3 |
235 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
236 |
|
|
237 |
C-- Set up work arrays that need valid initial values |
C-- Set up work arrays that need valid initial values |
238 |
DO j=1-OLy,sNy+OLy |
DO j=1-OLy,sNy+OLy |
239 |
DO i=1-OLx,sNx+OLx |
DO i=1-OLx,sNx+OLx |
248 |
fVerU (i,j,2) = 0. _d 0 |
fVerU (i,j,2) = 0. _d 0 |
249 |
fVerV (i,j,1) = 0. _d 0 |
fVerV (i,j,1) = 0. _d 0 |
250 |
fVerV (i,j,2) = 0. _d 0 |
fVerV (i,j,2) = 0. _d 0 |
251 |
phiHyd(i,j,1) = 0. _d 0 |
ENDDO |
252 |
K13 (i,j,1) = 0. _d 0 |
ENDDO |
253 |
K23 (i,j,1) = 0. _d 0 |
|
254 |
K33 (i,j,1) = 0. _d 0 |
DO k=1,Nr |
255 |
KapGM (i,j) = GMkbackground |
DO j=1-OLy,sNy+OLy |
256 |
|
DO i=1-OLx,sNx+OLx |
257 |
|
#ifdef INCLUDE_CONVECT_CALL |
258 |
|
ConvectCount(i,j,k) = 0. |
259 |
|
#endif |
260 |
|
KappaRT(i,j,k) = 0. _d 0 |
261 |
|
KappaRS(i,j,k) = 0. _d 0 |
262 |
|
ENDDO |
263 |
ENDDO |
ENDDO |
264 |
ENDDO |
ENDDO |
265 |
|
|
268 |
jMin = 1-OLy+1 |
jMin = 1-OLy+1 |
269 |
jMax = sNy+OLy |
jMax = sNy+OLy |
270 |
|
|
|
K = 1 |
|
|
BOTTOM_LAYER = K .EQ. Nr |
|
271 |
|
|
272 |
C-- Calculate gradient of surface pressure |
C-- Start of diagnostic loop |
273 |
CALL CALC_GRAD_ETA_SURF( |
DO k=Nr,1,-1 |
274 |
I bi,bj,iMin,iMax,jMin,jMax, |
|
275 |
O etaSurfX,etaSurfY, |
#ifdef ALLOW_AUTODIFF_TAMC |
276 |
I myThid) |
C? Patrick, is this formula correct now that we change the loop range? |
277 |
C-- Update fields in top level according to tendency terms |
C? Do we still need this? |
278 |
CALL CORRECTION_STEP( |
kkey = (ikey-1)*(Nr-2+1) + (k-2) + 1 |
279 |
I bi,bj,iMin,iMax,jMin,jMax,K, |
#endif /* ALLOW_AUTODIFF_TAMC */ |
280 |
I etaSurfX,etaSurfY,myTime,myThid) |
|
281 |
IF ( .NOT. BOTTOM_LAYER ) THEN |
C-- Integrate continuity vertically for vertical velocity |
282 |
C-- Update fields in layer below according to tendency terms |
CALL INTEGRATE_FOR_W( |
283 |
CALL CORRECTION_STEP( |
I bi, bj, k, uVel, vVel, |
284 |
I bi,bj,iMin,iMax,jMin,jMax,K+1, |
O wVel, |
285 |
I etaSurfX,etaSurfY,myTime,myThid) |
I myThid ) |
286 |
ENDIF |
|
287 |
C-- Density of 1st level (below W(1)) reference to level 1 |
#ifdef ALLOW_OBCS |
288 |
CALL FIND_RHO( |
#ifdef ALLOW_NONHYDROSTATIC |
289 |
I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType, |
C-- Apply OBC to W if in N-H mode |
290 |
O rhoKm1, |
IF (useOBCS.AND.nonHydrostatic) THEN |
291 |
I myThid ) |
CALL OBCS_APPLY_W( bi, bj, k, wVel, myThid ) |
292 |
|
ENDIF |
293 |
IF ( .NOT. BOTTOM_LAYER ) THEN |
#endif /* ALLOW_NONHYDROSTATIC */ |
294 |
C-- Check static stability with layer below |
#endif /* ALLOW_OBCS */ |
295 |
C-- and mix as needed. |
|
296 |
CALL FIND_RHO( |
C-- Calculate gradients of potential density for isoneutral |
297 |
I bi, bj, iMin, iMax, jMin, jMax, K+1, K, eosType, |
C slope terms (e.g. GM/Redi tensor or IVDC diffusivity) |
298 |
O rhoKp1, |
c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN |
299 |
I myThid ) |
IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) ) THEN |
300 |
CALL CONVECT( |
CALL FIND_RHO( |
301 |
I bi,bj,iMin,iMax,jMin,jMax,K+1,rhoKm1,rhoKp1, |
I bi, bj, iMin, iMax, jMin, jMax, k, k, eosType, |
302 |
I myTime,myIter,myThid) |
I theta, salt, |
303 |
C-- Recompute density after mixing |
O rhoK, |
|
CALL FIND_RHO( |
|
|
I bi, bj, iMin, iMax, jMin, jMax, K, K, eosType, |
|
|
O rhoKm1, |
|
|
I myThid ) |
|
|
ENDIF |
|
|
C-- Calculate buoyancy |
|
|
CALL CALC_BUOYANCY( |
|
|
I bi,bj,iMin,iMax,jMin,jMax,K,rhoKm1, |
|
|
O buoyKm1, |
|
|
I myThid ) |
|
|
C-- Integrate hydrostatic balance for phiHyd with BC of phiHyd(z=0)=0 |
|
|
CALL CALC_PHI_HYD( |
|
|
I bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyKm1, |
|
|
U phiHyd, |
|
|
I myThid ) |
|
|
|
|
|
DO K=2,Nr |
|
|
BOTTOM_LAYER = K .EQ. Nr |
|
|
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, |
|
|
I etaSurfX,etaSurfY,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 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-- Calculate buoyancy |
|
|
CALL CALC_BUOYANCY( |
|
|
I bi,bj,iMin,iMax,jMin,jMax,K,rhoK, |
|
|
O buoyK, |
|
|
I myThid ) |
|
|
C-- Integrate hydrostatic balance for phiHyd with BC of phiHyd(z=0)=0 |
|
|
CALL CALC_PHI_HYD( |
|
|
I bi,bj,iMin,iMax,jMin,jMax,K,buoyKm1,buoyK, |
|
|
U phiHyd, |
|
|
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, |
|
304 |
I myThid ) |
I myThid ) |
305 |
CALL CALC_ISOSLOPES( |
IF (k.GT.1) CALL FIND_RHO( |
306 |
I bi, bj, iMin, iMax, jMin, jMax, K, |
I bi, bj, iMin, iMax, jMin, jMax, k-1, k, eosType, |
307 |
I rhoKm1, rhoK, rhotmp, |
I theta, salt, |
308 |
O K13, K23, K33, KapGM, |
O rhoKm1, |
309 |
I myThid ) |
I myThid ) |
310 |
DO J=jMin,jMax |
CALL GRAD_SIGMA( |
311 |
DO I=iMin,iMax |
I bi, bj, iMin, iMax, jMin, jMax, k, |
312 |
rhoKm1 (I,J) = rhoK(I,J) |
I rhoK, rhoKm1, rhoK, |
313 |
buoyKm1(I,J) = buoyK(I,J) |
O sigmaX, sigmaY, sigmaR, |
314 |
|
I myThid ) |
315 |
|
ENDIF |
316 |
|
|
317 |
|
C-- Implicit Vertical Diffusion for Convection |
318 |
|
c ==> should use sigmaR !!! |
319 |
|
IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN |
320 |
|
CALL CALC_IVDC( |
321 |
|
I bi, bj, iMin, iMax, jMin, jMax, k, |
322 |
|
I rhoKm1, rhoK, |
323 |
|
U ConvectCount, KappaRT, KappaRS, |
324 |
|
I myTime, myIter, myThid) |
325 |
|
END IF |
326 |
|
|
327 |
|
C-- end of diagnostic k loop (Nr:1) |
328 |
|
ENDDO |
329 |
|
|
330 |
|
#ifdef ALLOW_OBCS |
331 |
|
C-- Calculate future values on open boundaries |
332 |
|
IF (useOBCS) THEN |
333 |
|
CALL OBCS_CALC( bi, bj, myTime+deltaT, |
334 |
|
I uVel, vVel, wVel, theta, salt, |
335 |
|
I myThid ) |
336 |
|
ENDIF |
337 |
|
#endif /* ALLOW_OBCS */ |
338 |
|
|
339 |
|
C-- Determines forcing terms based on external fields |
340 |
|
C relaxation terms, etc. |
341 |
|
CALL EXTERNAL_FORCING_SURF( |
342 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
343 |
|
I myThid ) |
344 |
|
|
345 |
|
#ifdef ALLOW_GMREDI |
346 |
|
C-- Calculate iso-neutral slopes for the GM/Redi parameterisation |
347 |
|
IF (useGMRedi) THEN |
348 |
|
DO k=1,Nr |
349 |
|
CALL GMREDI_CALC_TENSOR( |
350 |
|
I bi, bj, iMin, iMax, jMin, jMax, k, |
351 |
|
I sigmaX, sigmaY, sigmaR, |
352 |
|
I myThid ) |
353 |
ENDDO |
ENDDO |
354 |
ENDDO |
#ifdef ALLOW_AUTODIFF_TAMC |
355 |
ENDDO ! K |
ELSE |
356 |
|
DO k=1, Nr |
357 |
|
CALL GMREDI_CALC_TENSOR_DUMMY( |
358 |
|
I bi, bj, iMin, iMax, jMin, jMax, k, |
359 |
|
I sigmaX, sigmaY, sigmaR, |
360 |
|
I myThid ) |
361 |
|
ENDDO |
362 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
363 |
|
ENDIF |
364 |
|
#endif /* ALLOW_GMREDI */ |
365 |
|
|
366 |
|
#ifdef ALLOW_KPP |
367 |
|
C-- Compute KPP mixing coefficients |
368 |
|
IF (useKPP) THEN |
369 |
|
CALL KPP_CALC( |
370 |
|
I bi, bj, myTime, myThid ) |
371 |
|
ENDIF |
372 |
|
#endif /* ALLOW_KPP */ |
373 |
|
|
374 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
375 |
|
CADJ STORE KappaRT(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte |
376 |
|
CADJ STORE KappaRS(:,:,:) = comlev1_bibj, key = ikey, byte = isbyte |
377 |
|
CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte |
378 |
|
CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte |
379 |
|
CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte |
380 |
|
CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte |
381 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
382 |
|
|
383 |
|
#ifdef ALLOW_AIM |
384 |
|
C AIM - atmospheric intermediate model, physics package code. |
385 |
|
C note(jmc) : phiHyd=0 at this point but is not really used in Molteni Physics |
386 |
|
IF ( useAIM ) THEN |
387 |
|
CALL TIMER_START('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid) |
388 |
|
CALL AIM_DO_ATMOS_PHYSICS( phiHyd, myTime, myThid ) |
389 |
|
CALL TIMER_STOP ('AIM_DO_ATMOS_PHYS [DYNAMICS]', myThid) |
390 |
|
ENDIF |
391 |
|
#endif /* ALLOW_AIM */ |
392 |
|
|
393 |
|
|
394 |
DO K = Nr, 1, -1 |
C-- Start of thermodynamics loop |
395 |
|
DO k=Nr,1,-1 |
396 |
|
|
397 |
kM1 =max(1,k-1) ! Points to level above k (=k-1) |
C-- km1 Points to level above k (=k-1) |
398 |
kUp =1+MOD(k+1,2) ! Cycles through 1,2 to point to layer above |
C-- kup Cycles through 1,2 to point to layer above |
399 |
kDown=1+MOD(k,2) ! Cycles through 2,1 to point to current layer |
C-- kDown Cycles through 2,1 to point to current layer |
400 |
iMin = 1-OLx+2 |
|
401 |
iMax = sNx+OLx-1 |
km1 = MAX(1,k-1) |
402 |
jMin = 1-OLy+2 |
kup = 1+MOD(k+1,2) |
403 |
jMax = sNy+OLy-1 |
kDown= 1+MOD(k,2) |
404 |
|
|
405 |
|
iMin = 1-OLx+2 |
406 |
|
iMax = sNx+OLx-1 |
407 |
|
jMin = 1-OLy+2 |
408 |
|
jMax = sNy+OLy-1 |
409 |
|
|
410 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
411 |
|
CPatrick Is this formula correct? |
412 |
|
kkey = (ikey-1)*(Nr-1+1) + (k-1) + 1 |
413 |
|
CADJ STORE rvel (:,:,kDown) = comlev1_bibj_k, key = kkey, byte = isbyte |
414 |
|
CADJ STORE rTrans(:,:) = comlev1_bibj_k, key = kkey, byte = isbyte |
415 |
|
CADJ STORE KappaRT(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte |
416 |
|
CADJ STORE KappaRS(:,:,:) = comlev1_bibj_k, key = kkey, byte = isbyte |
417 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
418 |
|
|
419 |
C-- Get temporary terms used by tendency routines |
C-- Get temporary terms used by tendency routines |
420 |
CALL CALC_COMMON_FACTORS ( |
CALL CALC_COMMON_FACTORS ( |
421 |
I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown, |
I bi,bj,iMin,iMax,jMin,jMax,k,km1,kup,kDown, |
422 |
O xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp, |
O xA,yA,uTrans,vTrans,rTrans,rVel,maskC,maskUp, |
423 |
I myThid) |
I myThid) |
424 |
|
|
425 |
|
#ifdef INCLUDE_CALC_DIFFUSIVITY_CALL |
426 |
C-- Calculate the total vertical diffusivity |
C-- Calculate the total vertical diffusivity |
427 |
CALL CALC_DIFFUSIVITY( |
CALL CALC_DIFFUSIVITY( |
428 |
I bi,bj,iMin,iMax,jMin,jMax,K, |
I bi,bj,iMin,iMax,jMin,jMax,k, |
429 |
I maskC,maskUp,KapGM,K33, |
I maskC,maskup, |
430 |
O KappaRT,KappaRS, |
O KappaRT,KappaRS,KappaRU,KappaRV, |
431 |
I myThid) |
I myThid) |
432 |
C-- Calculate accelerations in the momentum equations |
#endif |
433 |
IF ( momStepping ) THEN |
|
434 |
CALL CALC_MOM_RHS( |
C-- Calculate active tracer tendencies (gT,gS,...) |
435 |
I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown, |
C and step forward storing result in gTnm1, gSnm1, etc. |
|
I xA,yA,uTrans,vTrans,rTrans,rVel,maskC, |
|
|
I phiHyd, |
|
|
U aTerm,xTerm,cTerm,mTerm,pTerm, |
|
|
U fZon, fMer, fVerU, fVerV, |
|
|
I myThid) |
|
|
ENDIF |
|
|
C-- Calculate active tracer tendencies |
|
436 |
IF ( tempStepping ) THEN |
IF ( tempStepping ) THEN |
437 |
CALL CALC_GT( |
CALL CALC_GT( |
438 |
I bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown, |
I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown, |
439 |
I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC, |
I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC, |
440 |
I K13,K23,KappaRT,KapGM, |
I KappaRT, |
441 |
U aTerm,xTerm,fZon,fMer,fVerT, |
U fVerT, |
442 |
I myThid) |
I myTime, myThid) |
443 |
|
CALL TIMESTEP_TRACER( |
444 |
|
I bi,bj,iMin,iMax,jMin,jMax,k, |
445 |
|
I theta, gT, |
446 |
|
U gTnm1, |
447 |
|
I myIter, myThid) |
448 |
ENDIF |
ENDIF |
449 |
IF ( saltStepping ) THEN |
IF ( saltStepping ) THEN |
450 |
CALL CALC_GS( |
CALL CALC_GS( |
451 |
I bi,bj,iMin,iMax,jMin,jMax, k,kM1,kUp,kDown, |
I bi,bj,iMin,iMax,jMin,jMax, k,km1,kup,kDown, |
452 |
I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC, |
I xA,yA,uTrans,vTrans,rTrans,maskUp,maskC, |
453 |
I K13,K23,KappaRS,KapGM, |
I KappaRS, |
454 |
U aTerm,xTerm,fZon,fMer,fVerS, |
U fVerS, |
455 |
I myThid) |
I myTime, myThid) |
456 |
|
CALL TIMESTEP_TRACER( |
457 |
|
I bi,bj,iMin,iMax,jMin,jMax,k, |
458 |
|
I salt, gS, |
459 |
|
U gSnm1, |
460 |
|
I myIter, myThid) |
461 |
ENDIF |
ENDIF |
|
C-- Prediction step (step forward all model variables) |
|
|
CALL TIMESTEP( |
|
|
I bi,bj,iMin,iMax,jMin,jMax,K, |
|
|
I myThid) |
|
|
C-- Diagnose barotropic divergence of predicted fields |
|
|
CALL CALC_DIV_GHAT( |
|
|
I bi,bj,iMin,iMax,jMin,jMax,K, |
|
|
I xA,yA, |
|
|
I myThid) |
|
|
|
|
|
C-- Cumulative diagnostic calculations (ie. time-averaging) |
|
|
#ifdef ALLOW_DIAGNOSTICS |
|
|
IF (taveFreq.GT.0.) THEN |
|
|
CALL DO_TIME_AVERAGES( |
|
|
I myTime, myIter, bi, bj, K, kUp, kDown, |
|
|
I K13, K23, rVel, KapGM, |
|
|
I myThid ) |
|
|
ENDIF |
|
|
#endif |
|
462 |
|
|
463 |
ENDDO ! K |
#ifdef ALLOW_OBCS |
464 |
|
C-- Apply open boundary conditions |
465 |
|
IF (useOBCS) THEN |
466 |
|
CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid ) |
467 |
|
END IF |
468 |
|
#endif /* ALLOW_OBCS */ |
469 |
|
|
470 |
|
C-- Freeze water |
471 |
|
IF (allowFreezing) THEN |
472 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
473 |
|
CADJ STORE gTNm1(:,:,k,bi,bj) = comlev1_bibj_k |
474 |
|
CADJ & , key = kkey, byte = isbyte |
475 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
476 |
|
CALL FREEZE( bi, bj, iMin, iMax, jMin, jMax, k, myThid ) |
477 |
|
END IF |
478 |
|
|
479 |
|
C-- end of thermodynamic k loop (Nr:1) |
480 |
|
ENDDO |
481 |
|
|
482 |
|
|
483 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
484 |
|
CPatrick? What about this one? |
485 |
|
maximpl = 6 |
486 |
|
iikey = (ikey-1)*maximpl |
487 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
488 |
|
|
489 |
C-- Implicit diffusion |
C-- Implicit diffusion |
490 |
IF (implicitDiffusion) THEN |
IF (implicitDiffusion) THEN |
491 |
CALL IMPLDIFF( bi, bj, iMin, iMax, jMin, jMax, |
|
492 |
I KappaRT,KappaRS, |
IF (tempStepping) THEN |
493 |
I myThid ) |
#ifdef ALLOW_AUTODIFF_TAMC |
494 |
|
idkey = iikey + 1 |
495 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
496 |
|
CALL IMPLDIFF( |
497 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
498 |
|
I deltaTtracer, KappaRT, recip_HFacC, |
499 |
|
U gTNm1, |
500 |
|
I myThid ) |
501 |
|
ENDIF |
502 |
|
|
503 |
|
IF (saltStepping) THEN |
504 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
505 |
|
idkey = iikey + 2 |
506 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
507 |
|
CALL IMPLDIFF( |
508 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
509 |
|
I deltaTtracer, KappaRS, recip_HFacC, |
510 |
|
U gSNm1, |
511 |
|
I myThid ) |
512 |
|
ENDIF |
513 |
|
|
514 |
|
#ifdef ALLOW_OBCS |
515 |
|
C-- Apply open boundary conditions |
516 |
|
IF (useOBCS) THEN |
517 |
|
DO K=1,Nr |
518 |
|
CALL OBCS_APPLY_TS( bi, bj, k, gTnm1, gSnm1, myThid ) |
519 |
|
ENDDO |
520 |
|
END IF |
521 |
|
#endif /* ALLOW_OBCS */ |
522 |
|
|
523 |
|
C-- End If implicitDiffusion |
524 |
|
ENDIF |
525 |
|
|
526 |
|
|
527 |
|
|
528 |
|
C-- Start of dynamics loop |
529 |
|
DO k=1,Nr |
530 |
|
|
531 |
|
C-- km1 Points to level above k (=k-1) |
532 |
|
C-- kup Cycles through 1,2 to point to layer above |
533 |
|
C-- kDown Cycles through 2,1 to point to current layer |
534 |
|
|
535 |
|
km1 = MAX(1,k-1) |
536 |
|
kup = 1+MOD(k+1,2) |
537 |
|
kDown= 1+MOD(k,2) |
538 |
|
|
539 |
|
iMin = 1-OLx+2 |
540 |
|
iMax = sNx+OLx-1 |
541 |
|
jMin = 1-OLy+2 |
542 |
|
jMax = sNy+OLy-1 |
543 |
|
|
544 |
|
C-- Integrate hydrostatic balance for phiHyd with BC of |
545 |
|
C phiHyd(z=0)=0 |
546 |
|
C distinguishe between Stagger and Non Stagger time stepping |
547 |
|
IF (staggerTimeStep) THEN |
548 |
|
CALL CALC_PHI_HYD( |
549 |
|
I bi,bj,iMin,iMax,jMin,jMax,k, |
550 |
|
I gTnm1, gSnm1, |
551 |
|
U phiHyd, |
552 |
|
I myThid ) |
553 |
|
ELSE |
554 |
|
CALL CALC_PHI_HYD( |
555 |
|
I bi,bj,iMin,iMax,jMin,jMax,k, |
556 |
|
I theta, salt, |
557 |
|
U phiHyd, |
558 |
|
I myThid ) |
559 |
|
ENDIF |
560 |
|
|
561 |
|
C-- Calculate accelerations in the momentum equations (gU, gV, ...) |
562 |
|
C and step forward storing the result in gUnm1, gVnm1, etc... |
563 |
|
IF ( momStepping ) THEN |
564 |
|
CALL CALC_MOM_RHS( |
565 |
|
I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown, |
566 |
|
I phiHyd,KappaRU,KappaRV, |
567 |
|
U fVerU, fVerV, |
568 |
|
I myTime, myThid) |
569 |
|
CALL TIMESTEP( |
570 |
|
I bi,bj,iMin,iMax,jMin,jMax,k,phiHyd, |
571 |
|
I myIter, myThid) |
572 |
|
|
573 |
|
#ifdef ALLOW_OBCS |
574 |
|
C-- Apply open boundary conditions |
575 |
|
IF (useOBCS) THEN |
576 |
|
CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid ) |
577 |
|
END IF |
578 |
|
#endif /* ALLOW_OBCS */ |
579 |
|
|
580 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
581 |
|
#ifdef INCLUDE_CD_CODE |
582 |
|
ELSE |
583 |
|
DO j=1-OLy,sNy+OLy |
584 |
|
DO i=1-OLx,sNx+OLx |
585 |
|
guCD(i,j,k,bi,bj) = 0.0 |
586 |
|
gvCD(i,j,k,bi,bj) = 0.0 |
587 |
|
END DO |
588 |
|
END DO |
589 |
|
#endif /* INCLUDE_CD_CODE */ |
590 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
591 |
|
ENDIF |
592 |
|
|
593 |
|
|
594 |
|
C-- end of dynamics k loop (1:Nr) |
595 |
|
ENDDO |
596 |
|
|
597 |
|
|
598 |
|
|
599 |
|
C-- Implicit viscosity |
600 |
|
IF (implicitViscosity.AND.momStepping) THEN |
601 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
602 |
|
idkey = iikey + 3 |
603 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
604 |
|
CALL IMPLDIFF( |
605 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
606 |
|
I deltaTmom, KappaRU,recip_HFacW, |
607 |
|
U gUNm1, |
608 |
|
I myThid ) |
609 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
610 |
|
idkey = iikey + 4 |
611 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
612 |
|
CALL IMPLDIFF( |
613 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
614 |
|
I deltaTmom, KappaRV,recip_HFacS, |
615 |
|
U gVNm1, |
616 |
|
I myThid ) |
617 |
|
|
618 |
|
#ifdef ALLOW_OBCS |
619 |
|
C-- Apply open boundary conditions |
620 |
|
IF (useOBCS) THEN |
621 |
|
DO K=1,Nr |
622 |
|
CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid ) |
623 |
|
ENDDO |
624 |
|
END IF |
625 |
|
#endif /* ALLOW_OBCS */ |
626 |
|
|
627 |
|
#ifdef INCLUDE_CD_CODE |
628 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
629 |
|
idkey = iikey + 5 |
630 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
631 |
|
CALL IMPLDIFF( |
632 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
633 |
|
I deltaTmom, KappaRU,recip_HFacW, |
634 |
|
U vVelD, |
635 |
|
I myThid ) |
636 |
|
#ifdef ALLOW_AUTODIFF_TAMC |
637 |
|
idkey = iikey + 6 |
638 |
|
#endif /* ALLOW_AUTODIFF_TAMC */ |
639 |
|
CALL IMPLDIFF( |
640 |
|
I bi, bj, iMin, iMax, jMin, jMax, |
641 |
|
I deltaTmom, KappaRV,recip_HFacS, |
642 |
|
U uVelD, |
643 |
|
I myThid ) |
644 |
|
#endif /* INCLUDE_CD_CODE */ |
645 |
|
C-- End If implicitViscosity.AND.momStepping |
646 |
ENDIF |
ENDIF |
647 |
|
|
648 |
ENDDO |
ENDDO |
649 |
ENDDO |
ENDDO |
650 |
|
|
|
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: rVel(1) ', |
|
|
C & minval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.), |
|
|
C & maxval(rVel(1:sNx,1:sNy,1),mask=rVel(1:sNx,1:sNy,1).NE.0.) |
|
|
C write(0,*) 'dynamics: rVel(2) ', |
|
|
C & minval(rVel(1:sNx,1:sNy,2),mask=rVel(1:sNx,1:sNy,2).NE.0.), |
|
|
C & maxval(rVel(1:sNx,1:sNy,2),mask=rVel(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: phiHyd ',minval(phiHyd/(Gravity*Rhonil),mask=phiHyd.NE.0.), |
|
|
C & maxval(phiHyd/(Gravity*Rhonil)) |
|
|
|
|
651 |
RETURN |
RETURN |
652 |
END |
END |