1 |
C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.83.2.5 2002/07/11 14:24:26 heimbach Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CPP_OPTIONS.h" |
5 |
|
6 |
CBOP |
7 |
C !ROUTINE: DYNAMICS |
8 |
C !INTERFACE: |
9 |
SUBROUTINE DYNAMICS(myTime, myIter, myThid) |
10 |
C !DESCRIPTION: \bv |
11 |
C *==========================================================* |
12 |
C | SUBROUTINE DYNAMICS |
13 |
C | o Controlling routine for the explicit part of the model |
14 |
C | dynamics. |
15 |
C *==========================================================* |
16 |
C | This routine evaluates the "dynamics" terms for each |
17 |
C | block of ocean in turn. Because the blocks of ocean have |
18 |
C | overlap regions they are independent of one another. |
19 |
C | If terms involving lateral integrals are needed in this |
20 |
C | routine care will be needed. Similarly finite-difference |
21 |
C | operations with stencils wider than the overlap region |
22 |
C | require special consideration. |
23 |
C | The algorithm... |
24 |
C | |
25 |
C | "Correction Step" |
26 |
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 === |
71 |
#include "SIZE.h" |
72 |
#include "EEPARAMS.h" |
73 |
#include "PARAMS.h" |
74 |
#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 |
# ifdef ALLOW_KPP |
84 |
# include "KPP.h" |
85 |
# endif |
86 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
87 |
#ifdef ALLOW_TIMEAVE |
88 |
#include "TIMEAVE_STATV.h" |
89 |
#endif |
90 |
|
91 |
C !CALLING SEQUENCE: |
92 |
C DYNAMICS() |
93 |
C | |
94 |
C |-- CALC_GRAD_PHI_SURF |
95 |
C | |
96 |
C |-- CALC_VISCOSITY |
97 |
C | |
98 |
C |-- CALC_PHI_HYD |
99 |
C | |
100 |
C |-- MOM_FLUXFORM |
101 |
C | |
102 |
C |-- MOM_VECINV |
103 |
C | |
104 |
C |-- TIMESTEP |
105 |
C | |
106 |
C |-- OBCS_APPLY_UV |
107 |
C | |
108 |
C |-- IMPLDIFF |
109 |
C | |
110 |
C |-- OBCS_APPLY_UV |
111 |
C | |
112 |
C |-- CALL TIMEAVE_CUMUL_1T |
113 |
C |-- CALL DEBUG_STATS_RL |
114 |
|
115 |
C !INPUT/OUTPUT PARAMETERS: |
116 |
C == Routine arguments == |
117 |
C myTime - Current time in simulation |
118 |
C myIter - Current iteration number in simulation |
119 |
C myThid - Thread number for this instance of the routine. |
120 |
_RL myTime |
121 |
INTEGER myIter |
122 |
INTEGER myThid |
123 |
|
124 |
C !LOCAL VARIABLES: |
125 |
C == Local variables |
126 |
C fVer[STUV] o fVer: Vertical flux term - note fVer |
127 |
C is "pipelined" in the vertical |
128 |
C so we need an fVer for each |
129 |
C variable. |
130 |
C rhoK, rhoKM1 - Density at current level, and level above |
131 |
C phiHyd - Hydrostatic part of the potential phiHydi. |
132 |
C In z coords phiHydiHyd is the hydrostatic |
133 |
C Potential (=pressure/rho0) anomaly |
134 |
C In p coords phiHydiHyd is the geopotential |
135 |
C surface height anomaly. |
136 |
C phiSurfX, - gradient of Surface potentiel (Pressure/rho, ocean) |
137 |
C phiSurfY or geopotentiel (atmos) in X and Y direction |
138 |
C iMin, iMax - Ranges and sub-block indices on which calculations |
139 |
C jMin, jMax are applied. |
140 |
C bi, bj |
141 |
C k, kup, - Index for layer above and below. kup and kDown |
142 |
C kDown, km1 are switched with layer to be the appropriate |
143 |
C index into fVerTerm. |
144 |
_RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
145 |
_RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2) |
146 |
_RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
147 |
_RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
148 |
_RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
149 |
_RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
150 |
_RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
151 |
_RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
152 |
_RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr) |
153 |
|
154 |
INTEGER iMin, iMax |
155 |
INTEGER jMin, jMax |
156 |
INTEGER bi, bj |
157 |
INTEGER i, j |
158 |
INTEGER k, km1, kp1, kup, kDown |
159 |
|
160 |
Cjmc : add for phiHyd output <- but not working if multi tile per CPU |
161 |
c CHARACTER*(MAX_LEN_MBUF) suff |
162 |
c LOGICAL DIFFERENT_MULTIPLE |
163 |
c EXTERNAL DIFFERENT_MULTIPLE |
164 |
Cjmc(end) |
165 |
|
166 |
C--- The algorithm... |
167 |
C |
168 |
C "Correction Step" |
169 |
C ================= |
170 |
C Here we update the horizontal velocities with the surface |
171 |
C pressure such that the resulting flow is either consistent |
172 |
C with the free-surface evolution or the rigid-lid: |
173 |
C U[n] = U* + dt x d/dx P |
174 |
C V[n] = V* + dt x d/dy P |
175 |
C |
176 |
C "Calculation of Gs" |
177 |
C =================== |
178 |
C This is where all the accelerations and tendencies (ie. |
179 |
C physics, parameterizations etc...) are calculated |
180 |
C rho = rho ( theta[n], salt[n] ) |
181 |
C b = b(rho, theta) |
182 |
C K31 = K31 ( rho ) |
183 |
C Gu[n] = Gu( u[n], v[n], wVel, b, ... ) |
184 |
C Gv[n] = Gv( u[n], v[n], wVel, b, ... ) |
185 |
C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... ) |
186 |
C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... ) |
187 |
C |
188 |
C "Time-stepping" or "Prediction" |
189 |
C ================================ |
190 |
C The models variables are stepped forward with the appropriate |
191 |
C time-stepping scheme (currently we use Adams-Bashforth II) |
192 |
C - For momentum, the result is always *only* a "prediction" |
193 |
C in that the flow may be divergent and will be "corrected" |
194 |
C later with a surface pressure gradient. |
195 |
C - Normally for tracers the result is the new field at time |
196 |
C level [n+1} *BUT* in the case of implicit diffusion the result |
197 |
C is also *only* a prediction. |
198 |
C - We denote "predictors" with an asterisk (*). |
199 |
C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] ) |
200 |
C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] ) |
201 |
C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] ) |
202 |
C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] ) |
203 |
C With implicit diffusion: |
204 |
C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] ) |
205 |
C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] ) |
206 |
C (1 + dt * K * d_zz) theta[n] = theta* |
207 |
C (1 + dt * K * d_zz) salt[n] = salt* |
208 |
C--- |
209 |
CEOP |
210 |
|
211 |
C-- Set up work arrays with valid (i.e. not NaN) values |
212 |
C These inital values do not alter the numerical results. They |
213 |
C just ensure that all memory references are to valid floating |
214 |
C point numbers. This prevents spurious hardware signals due to |
215 |
C uninitialised but inert locations. |
216 |
DO j=1-OLy,sNy+OLy |
217 |
DO i=1-OLx,sNx+OLx |
218 |
rhoKM1 (i,j) = 0. _d 0 |
219 |
rhok (i,j) = 0. _d 0 |
220 |
phiSurfX(i,j) = 0. _d 0 |
221 |
phiSurfY(i,j) = 0. _d 0 |
222 |
ENDDO |
223 |
ENDDO |
224 |
|
225 |
C-- Call to routine for calculation of |
226 |
C Eliassen-Palm-flux-forced U-tendency, |
227 |
C if desired: |
228 |
#ifdef INCLUDE_EP_FORCING_CODE |
229 |
CALL CALC_EP_FORCING(myThid) |
230 |
#endif |
231 |
|
232 |
#ifdef ALLOW_AUTODIFF_TAMC |
233 |
C-- HPF directive to help TAMC |
234 |
CHPF$ INDEPENDENT |
235 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
236 |
|
237 |
DO bj=myByLo(myThid),myByHi(myThid) |
238 |
|
239 |
#ifdef ALLOW_AUTODIFF_TAMC |
240 |
C-- HPF directive to help TAMC |
241 |
CHPF$ INDEPENDENT, NEW (fVerU,fVerV |
242 |
CHPF$& ,phiHyd |
243 |
CHPF$& ,KappaRU,KappaRV |
244 |
CHPF$& ) |
245 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
246 |
|
247 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
248 |
|
249 |
#ifdef ALLOW_AUTODIFF_TAMC |
250 |
act1 = bi - myBxLo(myThid) |
251 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
252 |
act2 = bj - myByLo(myThid) |
253 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
254 |
act3 = myThid - 1 |
255 |
max3 = nTx*nTy |
256 |
act4 = ikey_dynamics - 1 |
257 |
ikey = (act1 + 1) + act2*max1 |
258 |
& + act3*max1*max2 |
259 |
& + act4*max1*max2*max3 |
260 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
261 |
|
262 |
C-- Set up work arrays that need valid initial values |
263 |
DO j=1-OLy,sNy+OLy |
264 |
DO i=1-OLx,sNx+OLx |
265 |
DO k=1,Nr |
266 |
phiHyd(i,j,k) = 0. _d 0 |
267 |
KappaRU(i,j,k) = 0. _d 0 |
268 |
KappaRV(i,j,k) = 0. _d 0 |
269 |
ENDDO |
270 |
fVerU (i,j,1) = 0. _d 0 |
271 |
fVerU (i,j,2) = 0. _d 0 |
272 |
fVerV (i,j,1) = 0. _d 0 |
273 |
fVerV (i,j,2) = 0. _d 0 |
274 |
ENDDO |
275 |
ENDDO |
276 |
|
277 |
C-- Start computation of dynamics |
278 |
iMin = 1-OLx+2 |
279 |
iMax = sNx+OLx-1 |
280 |
jMin = 1-OLy+2 |
281 |
jMax = sNy+OLy-1 |
282 |
|
283 |
#ifdef ALLOW_AUTODIFF_TAMC |
284 |
CADJ STORE wvel (:,:,:,bi,bj) = comlev1_bibj, key = ikey, byte = isbyte |
285 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
286 |
|
287 |
C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP) |
288 |
C (note: this loop will be replaced by CALL CALC_GRAD_ETA) |
289 |
IF (implicSurfPress.NE.1.) THEN |
290 |
CALL CALC_GRAD_PHI_SURF( |
291 |
I bi,bj,iMin,iMax,jMin,jMax, |
292 |
I etaN, |
293 |
O phiSurfX,phiSurfY, |
294 |
I myThid ) |
295 |
ENDIF |
296 |
|
297 |
#ifdef ALLOW_AUTODIFF_TAMC |
298 |
CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
299 |
CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte |
300 |
#ifdef ALLOW_KPP |
301 |
CADJ STORE KPPviscAz (:,:,:,bi,bj) |
302 |
CADJ & = comlev1_bibj, key=ikey, byte=isbyte |
303 |
#endif /* ALLOW_KPP */ |
304 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
305 |
|
306 |
#ifdef INCLUDE_CALC_DIFFUSIVITY_CALL |
307 |
C-- Calculate the total vertical diffusivity |
308 |
DO k=1,Nr |
309 |
CALL CALC_VISCOSITY( |
310 |
I bi,bj,iMin,iMax,jMin,jMax,k, |
311 |
O KappaRU,KappaRV, |
312 |
I myThid) |
313 |
ENDDO |
314 |
#endif |
315 |
|
316 |
C-- Start of dynamics loop |
317 |
DO k=1,Nr |
318 |
|
319 |
C-- km1 Points to level above k (=k-1) |
320 |
C-- kup Cycles through 1,2 to point to layer above |
321 |
C-- kDown Cycles through 2,1 to point to current layer |
322 |
|
323 |
km1 = MAX(1,k-1) |
324 |
kp1 = MIN(k+1,Nr) |
325 |
kup = 1+MOD(k+1,2) |
326 |
kDown= 1+MOD(k,2) |
327 |
|
328 |
#ifdef ALLOW_AUTODIFF_TAMC |
329 |
kkey = (ikey-1)*Nr + k |
330 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
331 |
|
332 |
C-- Integrate hydrostatic balance for phiHyd with BC of |
333 |
C phiHyd(z=0)=0 |
334 |
C distinguishe between Stagger and Non Stagger time stepping |
335 |
IF (staggerTimeStep) THEN |
336 |
CALL CALC_PHI_HYD( |
337 |
I bi,bj,iMin,iMax,jMin,jMax,k, |
338 |
I gT, gS, |
339 |
U phiHyd, |
340 |
I myThid ) |
341 |
ELSE |
342 |
CALL CALC_PHI_HYD( |
343 |
I bi,bj,iMin,iMax,jMin,jMax,k, |
344 |
I theta, salt, |
345 |
U phiHyd, |
346 |
I myThid ) |
347 |
ENDIF |
348 |
|
349 |
C-- Calculate accelerations in the momentum equations (gU, gV, ...) |
350 |
C and step forward storing the result in gUnm1, gVnm1, etc... |
351 |
IF ( momStepping ) THEN |
352 |
#ifndef DISABLE_MOM_FLUXFORM |
353 |
IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM( |
354 |
I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown, |
355 |
I phiHyd,KappaRU,KappaRV, |
356 |
U fVerU, fVerV, |
357 |
I myTime, myIter, myThid) |
358 |
#endif |
359 |
#ifndef DISABLE_MOM_VECINV |
360 |
IF (vectorInvariantMomentum) CALL MOM_VECINV( |
361 |
I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown, |
362 |
I phiHyd,KappaRU,KappaRV, |
363 |
U fVerU, fVerV, |
364 |
I myTime, myIter, myThid) |
365 |
#endif |
366 |
CALL TIMESTEP( |
367 |
I bi,bj,iMin,iMax,jMin,jMax,k, |
368 |
I phiHyd, phiSurfX, phiSurfY, |
369 |
I myIter, myThid) |
370 |
|
371 |
#ifdef ALLOW_OBCS |
372 |
C-- Apply open boundary conditions |
373 |
IF (useOBCS) THEN |
374 |
CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid ) |
375 |
END IF |
376 |
#endif /* ALLOW_OBCS */ |
377 |
|
378 |
#ifdef ALLOW_AUTODIFF_TAMC |
379 |
#ifdef INCLUDE_CD_CODE |
380 |
ELSE |
381 |
DO j=1-OLy,sNy+OLy |
382 |
DO i=1-OLx,sNx+OLx |
383 |
guCD(i,j,k,bi,bj) = 0.0 |
384 |
gvCD(i,j,k,bi,bj) = 0.0 |
385 |
END DO |
386 |
END DO |
387 |
#endif /* INCLUDE_CD_CODE */ |
388 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
389 |
ENDIF |
390 |
|
391 |
|
392 |
C-- end of dynamics k loop (1:Nr) |
393 |
ENDDO |
394 |
|
395 |
C-- Implicit viscosity |
396 |
IF (implicitViscosity.AND.momStepping) THEN |
397 |
#ifdef ALLOW_AUTODIFF_TAMC |
398 |
CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte |
399 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
400 |
CALL IMPLDIFF( |
401 |
I bi, bj, iMin, iMax, jMin, jMax, |
402 |
I deltaTmom, KappaRU,recip_HFacW, |
403 |
U gUNm1, |
404 |
I myThid ) |
405 |
#ifdef ALLOW_AUTODIFF_TAMC |
406 |
CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte |
407 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
408 |
CALL IMPLDIFF( |
409 |
I bi, bj, iMin, iMax, jMin, jMax, |
410 |
I deltaTmom, KappaRV,recip_HFacS, |
411 |
U gVNm1, |
412 |
I myThid ) |
413 |
|
414 |
#ifdef ALLOW_OBCS |
415 |
C-- Apply open boundary conditions |
416 |
IF (useOBCS) THEN |
417 |
DO K=1,Nr |
418 |
CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid ) |
419 |
ENDDO |
420 |
END IF |
421 |
#endif /* ALLOW_OBCS */ |
422 |
|
423 |
#ifdef INCLUDE_CD_CODE |
424 |
#ifdef ALLOW_AUTODIFF_TAMC |
425 |
CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte |
426 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
427 |
CALL IMPLDIFF( |
428 |
I bi, bj, iMin, iMax, jMin, jMax, |
429 |
I deltaTmom, KappaRU,recip_HFacW, |
430 |
U vVelD, |
431 |
I myThid ) |
432 |
#ifdef ALLOW_AUTODIFF_TAMC |
433 |
CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=ikey, byte=isbyte |
434 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
435 |
CALL IMPLDIFF( |
436 |
I bi, bj, iMin, iMax, jMin, jMax, |
437 |
I deltaTmom, KappaRV,recip_HFacS, |
438 |
U uVelD, |
439 |
I myThid ) |
440 |
#endif /* INCLUDE_CD_CODE */ |
441 |
C-- End If implicitViscosity.AND.momStepping |
442 |
ENDIF |
443 |
|
444 |
Cjmc : add for phiHyd output <- but not working if multi tile per CPU |
445 |
c IF ( DIFFERENT_MULTIPLE(dumpFreq,myTime+deltaTClock,myTime) |
446 |
c & .AND. buoyancyRelation .eq. 'ATMOSPHERIC' ) THEN |
447 |
c WRITE(suff,'(I10.10)') myIter+1 |
448 |
c CALL WRITE_FLD_XYZ_RL('PH.',suff,phiHyd,myIter+1,myThid) |
449 |
c ENDIF |
450 |
Cjmc(end) |
451 |
|
452 |
#ifdef ALLOW_TIMEAVE |
453 |
IF (taveFreq.GT.0.) THEN |
454 |
CALL TIMEAVE_CUMUL_1T(phiHydtave, phiHyd, Nr, |
455 |
I deltaTclock, bi, bj, myThid) |
456 |
ENDIF |
457 |
#endif /* ALLOW_TIMEAVE */ |
458 |
|
459 |
ENDDO |
460 |
ENDDO |
461 |
|
462 |
#ifndef DISABLE_DEBUGMODE |
463 |
If (debugMode) THEN |
464 |
CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid) |
465 |
CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid) |
466 |
CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid) |
467 |
CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid) |
468 |
CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid) |
469 |
CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid) |
470 |
CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid) |
471 |
CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid) |
472 |
CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid) |
473 |
CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid) |
474 |
CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid) |
475 |
CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid) |
476 |
CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid) |
477 |
CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid) |
478 |
ENDIF |
479 |
#endif |
480 |
|
481 |
RETURN |
482 |
END |