/[MITgcm]/MITgcm/model/src/dynamics.F
ViewVC logotype

Contents of /MITgcm/model/src/dynamics.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.120 - (show annotations) (download)
Mon Jul 11 19:30:42 2005 UTC (18 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint57n_post, checkpoint57l_post, checkpoint57o_post
Changes since 1.119: +8 -32 lines
call diagnostics_scale_fill (instead of diagnostics_fill) and avoid temp array

1 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.119 2005/06/27 12:27:19 jmc Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 CBOP
8 C !ROUTINE: DYNAMICS
9 C !INTERFACE:
10 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
11 C !DESCRIPTION: \bv
12 C *==========================================================*
13 C | SUBROUTINE DYNAMICS
14 C | o Controlling routine for the explicit part of the model
15 C | dynamics.
16 C *==========================================================*
17 C | This routine evaluates the "dynamics" terms for each
18 C | block of ocean in turn. Because the blocks of ocean have
19 C | overlap regions they are independent of one another.
20 C | If terms involving lateral integrals are needed in this
21 C | routine care will be needed. Similarly finite-difference
22 C | operations with stencils wider than the overlap region
23 C | require special consideration.
24 C | The algorithm...
25 C |
26 C | "Correction Step"
27 C | =================
28 C | Here we update the horizontal velocities with the surface
29 C | pressure such that the resulting flow is either consistent
30 C | with the free-surface evolution or the rigid-lid:
31 C | U[n] = U* + dt x d/dx P
32 C | V[n] = V* + dt x d/dy P
33 C |
34 C | "Calculation of Gs"
35 C | ===================
36 C | This is where all the accelerations and tendencies (ie.
37 C | physics, parameterizations etc...) are calculated
38 C | rho = rho ( theta[n], salt[n] )
39 C | b = b(rho, theta)
40 C | K31 = K31 ( rho )
41 C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
42 C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
43 C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
44 C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
45 C |
46 C | "Time-stepping" or "Prediction"
47 C | ================================
48 C | The models variables are stepped forward with the appropriate
49 C | time-stepping scheme (currently we use Adams-Bashforth II)
50 C | - For momentum, the result is always *only* a "prediction"
51 C | in that the flow may be divergent and will be "corrected"
52 C | later with a surface pressure gradient.
53 C | - Normally for tracers the result is the new field at time
54 C | level [n+1} *BUT* in the case of implicit diffusion the result
55 C | is also *only* a prediction.
56 C | - We denote "predictors" with an asterisk (*).
57 C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
58 C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
59 C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
60 C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
61 C | With implicit diffusion:
62 C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
63 C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
64 C | (1 + dt * K * d_zz) theta[n] = theta*
65 C | (1 + dt * K * d_zz) salt[n] = salt*
66 C |
67 C *==========================================================*
68 C \ev
69 C !USES:
70 IMPLICIT NONE
71 C == Global variables ===
72 #include "SIZE.h"
73 #include "EEPARAMS.h"
74 #include "PARAMS.h"
75 #include "DYNVARS.h"
76 #ifdef ALLOW_CD_CODE
77 #include "CD_CODE_VARS.h"
78 #endif
79 #include "GRID.h"
80 #ifdef ALLOW_AUTODIFF_TAMC
81 # include "tamc.h"
82 # include "tamc_keys.h"
83 # include "FFIELDS.h"
84 # include "EOS.h"
85 # ifdef ALLOW_KPP
86 # include "KPP.h"
87 # endif
88 #endif /* ALLOW_AUTODIFF_TAMC */
89
90 C !CALLING SEQUENCE:
91 C DYNAMICS()
92 C |
93 C |-- CALC_GRAD_PHI_SURF
94 C |
95 C |-- CALC_VISCOSITY
96 C |
97 C |-- CALC_PHI_HYD
98 C |
99 C |-- MOM_FLUXFORM
100 C |
101 C |-- MOM_VECINV
102 C |
103 C |-- TIMESTEP
104 C |
105 C |-- OBCS_APPLY_UV
106 C |
107 C |-- IMPLDIFF
108 C |
109 C |-- OBCS_APPLY_UV
110 C |
111 C |-- CALL DEBUG_STATS_RL
112
113 C !INPUT/OUTPUT PARAMETERS:
114 C == Routine arguments ==
115 C myTime - Current time in simulation
116 C myIter - Current iteration number in simulation
117 C myThid - Thread number for this instance of the routine.
118 _RL myTime
119 INTEGER myIter
120 INTEGER myThid
121
122 C !LOCAL VARIABLES:
123 C == Local variables
124 C fVer[UV] o fVer: Vertical flux term - note fVer
125 C is "pipelined" in the vertical
126 C so we need an fVer for each
127 C variable.
128 C phiHydC :: hydrostatic potential anomaly at cell center
129 C In z coords phiHyd is the hydrostatic potential
130 C (=pressure/rho0) anomaly
131 C In p coords phiHyd is the geopotential height anomaly.
132 C phiHydF :: hydrostatic potential anomaly at middle between 2 centers
133 C dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
134 C phiSurfX, :: gradient of Surface potential (Pressure/rho, ocean)
135 C phiSurfY or geopotential (atmos) in X and Y direction
136 C guDissip :: dissipation tendency (all explicit terms), u component
137 C gvDissip :: dissipation tendency (all explicit terms), v component
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 phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
147 _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
148 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
149 _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
150 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
151 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
152 _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
153 _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
154 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
155 _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
156
157 INTEGER iMin, iMax
158 INTEGER jMin, jMax
159 INTEGER bi, bj
160 INTEGER i, j
161 INTEGER k, km1, kp1, kup, kDown
162
163 #ifdef ALLOW_DIAGNOSTICS
164 _RL tmpFac
165 #endif /* ALLOW_DIAGNOSTICS */
166
167
168 C--- The algorithm...
169 C
170 C "Correction Step"
171 C =================
172 C Here we update the horizontal velocities with the surface
173 C pressure such that the resulting flow is either consistent
174 C with the free-surface evolution or the rigid-lid:
175 C U[n] = U* + dt x d/dx P
176 C V[n] = V* + dt x d/dy P
177 C
178 C "Calculation of Gs"
179 C ===================
180 C This is where all the accelerations and tendencies (ie.
181 C physics, parameterizations etc...) are calculated
182 C rho = rho ( theta[n], salt[n] )
183 C b = b(rho, theta)
184 C K31 = K31 ( rho )
185 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
186 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
187 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
188 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
189 C
190 C "Time-stepping" or "Prediction"
191 C ================================
192 C The models variables are stepped forward with the appropriate
193 C time-stepping scheme (currently we use Adams-Bashforth II)
194 C - For momentum, the result is always *only* a "prediction"
195 C in that the flow may be divergent and will be "corrected"
196 C later with a surface pressure gradient.
197 C - Normally for tracers the result is the new field at time
198 C level [n+1} *BUT* in the case of implicit diffusion the result
199 C is also *only* a prediction.
200 C - We denote "predictors" with an asterisk (*).
201 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
202 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
203 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
204 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
205 C With implicit diffusion:
206 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
207 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
208 C (1 + dt * K * d_zz) theta[n] = theta*
209 C (1 + dt * K * d_zz) salt[n] = salt*
210 C---
211 CEOP
212
213 C-- Call to routine for calculation of
214 C Eliassen-Palm-flux-forced U-tendency,
215 C if desired:
216 #ifdef INCLUDE_EP_FORCING_CODE
217 CALL CALC_EP_FORCING(myThid)
218 #endif
219
220 #ifdef ALLOW_AUTODIFF_TAMC
221 C-- HPF directive to help TAMC
222 CHPF$ INDEPENDENT
223 #endif /* ALLOW_AUTODIFF_TAMC */
224
225 DO bj=myByLo(myThid),myByHi(myThid)
226
227 #ifdef ALLOW_AUTODIFF_TAMC
228 C-- HPF directive to help TAMC
229 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
230 CHPF$& ,phiHydF
231 CHPF$& ,KappaRU,KappaRV
232 CHPF$& )
233 #endif /* ALLOW_AUTODIFF_TAMC */
234
235 DO bi=myBxLo(myThid),myBxHi(myThid)
236
237 #ifdef ALLOW_AUTODIFF_TAMC
238 act1 = bi - myBxLo(myThid)
239 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
240 act2 = bj - myByLo(myThid)
241 max2 = myByHi(myThid) - myByLo(myThid) + 1
242 act3 = myThid - 1
243 max3 = nTx*nTy
244 act4 = ikey_dynamics - 1
245 idynkey = (act1 + 1) + act2*max1
246 & + act3*max1*max2
247 & + act4*max1*max2*max3
248 #endif /* ALLOW_AUTODIFF_TAMC */
249
250 C-- Set up work arrays with valid (i.e. not NaN) values
251 C These inital values do not alter the numerical results. They
252 C just ensure that all memory references are to valid floating
253 C point numbers. This prevents spurious hardware signals due to
254 C uninitialised but inert locations.
255
256 DO k=1,Nr
257 DO j=1-OLy,sNy+OLy
258 DO i=1-OLx,sNx+OLx
259 KappaRU(i,j,k) = 0. _d 0
260 KappaRV(i,j,k) = 0. _d 0
261 #ifdef ALLOW_AUTODIFF_TAMC
262 cph(
263 c-- need some re-initialisation here to break dependencies
264 cph)
265 gu(i,j,k,bi,bj) = 0. _d 0
266 gv(i,j,k,bi,bj) = 0. _d 0
267 #endif
268 ENDDO
269 ENDDO
270 ENDDO
271 DO j=1-OLy,sNy+OLy
272 DO i=1-OLx,sNx+OLx
273 fVerU (i,j,1) = 0. _d 0
274 fVerU (i,j,2) = 0. _d 0
275 fVerV (i,j,1) = 0. _d 0
276 fVerV (i,j,2) = 0. _d 0
277 phiHydF (i,j) = 0. _d 0
278 phiHydC (i,j) = 0. _d 0
279 dPhiHydX(i,j) = 0. _d 0
280 dPhiHydY(i,j) = 0. _d 0
281 phiSurfX(i,j) = 0. _d 0
282 phiSurfY(i,j) = 0. _d 0
283 guDissip(i,j) = 0. _d 0
284 gvDissip(i,j) = 0. _d 0
285 ENDDO
286 ENDDO
287
288 C-- Start computation of dynamics
289 iMin = 0
290 iMax = sNx+1
291 jMin = 0
292 jMax = sNy+1
293
294 #ifdef ALLOW_AUTODIFF_TAMC
295 CADJ STORE wvel (:,:,:,bi,bj) =
296 CADJ & comlev1_bibj, key = idynkey, byte = isbyte
297 #endif /* ALLOW_AUTODIFF_TAMC */
298
299 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
300 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
301 IF (implicSurfPress.NE.1.) THEN
302 CALL CALC_GRAD_PHI_SURF(
303 I bi,bj,iMin,iMax,jMin,jMax,
304 I etaN,
305 O phiSurfX,phiSurfY,
306 I myThid )
307 ENDIF
308
309 #ifdef ALLOW_AUTODIFF_TAMC
310 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
311 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
312 #ifdef ALLOW_KPP
313 CADJ STORE KPPviscAz (:,:,:,bi,bj)
314 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
315 #endif /* ALLOW_KPP */
316 #endif /* ALLOW_AUTODIFF_TAMC */
317
318 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
319 C-- Calculate the total vertical diffusivity
320 DO k=1,Nr
321 CALL CALC_VISCOSITY(
322 I bi,bj,iMin,iMax,jMin,jMax,k,
323 O KappaRU,KappaRV,
324 I myThid)
325 ENDDO
326 #endif
327
328 #ifdef ALLOW_AUTODIFF_TAMC
329 CADJ STORE KappaRU(:,:,:)
330 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
331 CADJ STORE KappaRV(:,:,:)
332 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
333 #endif /* ALLOW_AUTODIFF_TAMC */
334
335 C-- Start of dynamics loop
336 DO k=1,Nr
337
338 C-- km1 Points to level above k (=k-1)
339 C-- kup Cycles through 1,2 to point to layer above
340 C-- kDown Cycles through 2,1 to point to current layer
341
342 km1 = MAX(1,k-1)
343 kp1 = MIN(k+1,Nr)
344 kup = 1+MOD(k+1,2)
345 kDown= 1+MOD(k,2)
346
347 #ifdef ALLOW_AUTODIFF_TAMC
348 kkey = (idynkey-1)*Nr + k
349 c
350 CADJ STORE totphihyd (:,:,k,bi,bj)
351 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
352 CADJ STORE theta (:,:,k,bi,bj)
353 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
354 CADJ STORE salt (:,:,k,bi,bj)
355 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
356 #endif /* ALLOW_AUTODIFF_TAMC */
357
358 C-- Integrate hydrostatic balance for phiHyd with BC of
359 C phiHyd(z=0)=0
360 CALL CALC_PHI_HYD(
361 I bi,bj,iMin,iMax,jMin,jMax,k,
362 I theta, salt,
363 U phiHydF,
364 O phiHydC, dPhiHydX, dPhiHydY,
365 I myTime, myIter, myThid )
366
367 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
368 C and step forward storing the result in gU, gV, etc...
369 IF ( momStepping ) THEN
370 #ifdef ALLOW_MOM_FLUXFORM
371 IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
372 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
373 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
374 U fVerU, fVerV,
375 I myTime, myIter, myThid)
376 #endif
377 #ifdef ALLOW_MOM_VECINV
378 IF (vectorInvariantMomentum) CALL MOM_VECINV(
379 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
380 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
381 U fVerU, fVerV,
382 O guDissip, gvDissip,
383 I myTime, myIter, myThid)
384 #endif
385 CALL TIMESTEP(
386 I bi,bj,iMin,iMax,jMin,jMax,k,
387 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
388 I guDissip, gvDissip,
389 I myTime, myIter, myThid)
390
391 #ifdef ALLOW_OBCS
392 C-- Apply open boundary conditions
393 IF (useOBCS) THEN
394 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
395 ENDIF
396 #endif /* ALLOW_OBCS */
397
398 ENDIF
399
400
401 C-- end of dynamics k loop (1:Nr)
402 ENDDO
403
404 C-- Implicit Vertical advection & viscosity
405 #ifdef INCLUDE_IMPLVERTADV_CODE
406 IF ( momImplVertAdv ) THEN
407 CALL MOM_U_IMPLICIT_R( kappaRU,
408 I bi, bj, myTime, myIter, myThid )
409 CALL MOM_V_IMPLICIT_R( kappaRV,
410 I bi, bj, myTime, myIter, myThid )
411 ELSEIF ( implicitViscosity ) THEN
412 #else /* INCLUDE_IMPLVERTADV_CODE */
413 IF ( implicitViscosity ) THEN
414 #endif /* INCLUDE_IMPLVERTADV_CODE */
415 #ifdef ALLOW_AUTODIFF_TAMC
416 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
417 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
418 #endif /* ALLOW_AUTODIFF_TAMC */
419 CALL IMPLDIFF(
420 I bi, bj, iMin, iMax, jMin, jMax,
421 I 0, KappaRU,recip_HFacW,
422 U gU,
423 I myThid )
424 #ifdef ALLOW_AUTODIFF_TAMC
425 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
426 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
427 #endif /* ALLOW_AUTODIFF_TAMC */
428 CALL IMPLDIFF(
429 I bi, bj, iMin, iMax, jMin, jMax,
430 I 0, KappaRV,recip_HFacS,
431 U gV,
432 I myThid )
433 ENDIF
434
435 #ifdef ALLOW_OBCS
436 C-- Apply open boundary conditions
437 IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
438 DO K=1,Nr
439 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
440 ENDDO
441 ENDIF
442 #endif /* ALLOW_OBCS */
443
444 #ifdef ALLOW_CD_CODE
445 IF (implicitViscosity.AND.useCDscheme) THEN
446 #ifdef ALLOW_AUTODIFF_TAMC
447 CADJ STORE vVelD(:,:,:,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 0, KappaRU,recip_HFacW,
452 U vVelD,
453 I myThid )
454 #ifdef ALLOW_AUTODIFF_TAMC
455 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
456 #endif /* ALLOW_AUTODIFF_TAMC */
457 CALL IMPLDIFF(
458 I bi, bj, iMin, iMax, jMin, jMax,
459 I 0, KappaRV,recip_HFacS,
460 U uVelD,
461 I myThid )
462 ENDIF
463 #endif /* ALLOW_CD_CODE */
464 C-- End implicit Vertical advection & viscosity
465
466 ENDDO
467 ENDDO
468
469 #ifdef ALLOW_OBCS
470 IF (useOBCS) THEN
471 CALL OBCS_PRESCRIBE_EXCHANGES(myThid)
472 ENDIF
473 #endif
474
475 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
476
477 Cml(
478 C In order to compare the variance of phiHydLow of a p/z-coordinate
479 C run with etaH of a z/p-coordinate run the drift of phiHydLow
480 C has to be removed by something like the following subroutine:
481 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
482 C & 'phiHydLow', myThid )
483 Cml)
484
485 #ifdef ALLOW_DIAGNOSTICS
486 IF ( usediagnostics ) THEN
487
488 CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
489 CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0, 1,0,1,1,myThid)
490
491 tmpFac = 1. _d 0
492 CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
493 & 'PHIHYDSQ',0,Nr,0,1,1,myThid)
494
495 CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
496 & 'PHIBOTSQ',0, 1,0,1,1,myThid)
497
498 ENDIF
499 #endif /* ALLOW_DIAGNOSTICS */
500
501 #ifdef ALLOW_DEBUG
502 If ( debugLevel .GE. debLevB ) THEN
503 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
504 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
505 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
506 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
507 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
508 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
509 CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
510 CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
511 CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
512 CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
513 #ifndef ALLOW_ADAMSBASHFORTH_3
514 CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
515 CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
516 CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
517 CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
518 #endif
519 ENDIF
520 #endif
521
522 RETURN
523 END

  ViewVC Help
Powered by ViewVC 1.1.22