/[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.106 - (show annotations) (download)
Sat Jan 3 01:01:34 2004 UTC (20 years, 4 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint52l_pre, hrcube4, hrcube5, checkpoint52j_pre, checkpoint52l_post, checkpoint52k_post, checkpoint54, checkpoint53, checkpoint52f_post, checkpoint52i_pre, hrcube_1, hrcube_2, hrcube_3, checkpoint52e_post, checkpoint53d_post, checkpoint52m_post, checkpoint52f_pre, checkpoint54a_pre, checkpoint53c_post, checkpoint53a_post, checkpoint53g_post, checkpoint52i_post, checkpoint52h_pre, checkpoint53f_post, checkpoint52j_post, checkpoint52n_post, checkpoint53b_pre, checkpoint53b_post, checkpoint53d_pre
Changes since 1.105: +18 -7 lines
add calls for implicit vertical direction (advection & diffusion)
    but keep impldiff for implicit diffusion & viscosity only.

1 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.105 2003/11/04 19:51:53 edhill 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_PASSIVE_TRACER
81 #include "TR1.h"
82 #endif
83 #ifdef ALLOW_AUTODIFF_TAMC
84 # include "tamc.h"
85 # include "tamc_keys.h"
86 # include "FFIELDS.h"
87 # include "EOS.h"
88 # ifdef ALLOW_KPP
89 # include "KPP.h"
90 # endif
91 #endif /* ALLOW_AUTODIFF_TAMC */
92
93 C !CALLING SEQUENCE:
94 C DYNAMICS()
95 C |
96 C |-- CALC_GRAD_PHI_SURF
97 C |
98 C |-- CALC_VISCOSITY
99 C |
100 C |-- CALC_PHI_HYD
101 C |
102 C |-- MOM_FLUXFORM
103 C |
104 C |-- MOM_VECINV
105 C |
106 C |-- TIMESTEP
107 C |
108 C |-- OBCS_APPLY_UV
109 C |
110 C |-- IMPLDIFF
111 C |
112 C |-- OBCS_APPLY_UV
113 C |
114 C |-- CALL TIMEAVE_CUMUL_1T
115 C |-- CALL DEBUG_STATS_RL
116
117 C !INPUT/OUTPUT PARAMETERS:
118 C == Routine arguments ==
119 C myTime - Current time in simulation
120 C myIter - Current iteration number in simulation
121 C myThid - Thread number for this instance of the routine.
122 _RL myTime
123 INTEGER myIter
124 INTEGER myThid
125
126 C !LOCAL VARIABLES:
127 C == Local variables
128 C fVer[STUV] o fVer: Vertical flux term - note fVer
129 C is "pipelined" in the vertical
130 C so we need an fVer for each
131 C variable.
132 C phiHydC :: hydrostatic potential anomaly at cell center
133 C In z coords phiHyd is the hydrostatic potential
134 C (=pressure/rho0) anomaly
135 C In p coords phiHyd is the geopotential height anomaly.
136 C phiHydF :: hydrostatic potential anomaly at middle between 2 centers
137 C dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
138 C phiSurfX, :: gradient of Surface potential (Pressure/rho, ocean)
139 C phiSurfY or geopotential (atmos) in X and Y direction
140 C iMin, iMax - Ranges and sub-block indices on which calculations
141 C jMin, jMax are applied.
142 C bi, bj
143 C k, kup, - Index for layer above and below. kup and kDown
144 C kDown, km1 are switched with layer to be the appropriate
145 C index into fVerTerm.
146 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
147 _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
148 _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
149 _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
150 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
151 _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
152 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
153 _RL phiSurfY(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 LOGICAL DIFFERENT_MULTIPLE
164 EXTERNAL DIFFERENT_MULTIPLE
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-- Call to routine for calculation of
212 C Eliassen-Palm-flux-forced U-tendency,
213 C if desired:
214 #ifdef INCLUDE_EP_FORCING_CODE
215 CALL CALC_EP_FORCING(myThid)
216 #endif
217
218 #ifdef ALLOW_AUTODIFF_TAMC
219 C-- HPF directive to help TAMC
220 CHPF$ INDEPENDENT
221 #endif /* ALLOW_AUTODIFF_TAMC */
222
223 DO bj=myByLo(myThid),myByHi(myThid)
224
225 #ifdef ALLOW_AUTODIFF_TAMC
226 C-- HPF directive to help TAMC
227 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
228 CHPF$& ,phiHydF
229 CHPF$& ,KappaRU,KappaRV
230 CHPF$& )
231 #endif /* ALLOW_AUTODIFF_TAMC */
232
233 DO bi=myBxLo(myThid),myBxHi(myThid)
234
235 #ifdef ALLOW_AUTODIFF_TAMC
236 act1 = bi - myBxLo(myThid)
237 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
238 act2 = bj - myByLo(myThid)
239 max2 = myByHi(myThid) - myByLo(myThid) + 1
240 act3 = myThid - 1
241 max3 = nTx*nTy
242 act4 = ikey_dynamics - 1
243 idynkey = (act1 + 1) + act2*max1
244 & + act3*max1*max2
245 & + act4*max1*max2*max3
246 #endif /* ALLOW_AUTODIFF_TAMC */
247
248 C-- Set up work arrays with valid (i.e. not NaN) values
249 C These inital values do not alter the numerical results. They
250 C just ensure that all memory references are to valid floating
251 C point numbers. This prevents spurious hardware signals due to
252 C uninitialised but inert locations.
253
254 DO k=1,Nr
255 DO j=1-OLy,sNy+OLy
256 DO i=1-OLx,sNx+OLx
257 KappaRU(i,j,k) = 0. _d 0
258 KappaRV(i,j,k) = 0. _d 0
259 #ifdef ALLOW_AUTODIFF_TAMC
260 cph(
261 c-- need some re-initialisation here to break dependencies
262 c-- totphihyd is assumed zero from ini_pressure, i.e.
263 c-- avoiding iterate pressure p = integral of (g*rho(p)*dz)
264 cph)
265 totPhiHyd(i,j,k,bi,bj) = 0. _d 0
266 gu(i,j,k,bi,bj) = 0. _d 0
267 gv(i,j,k,bi,bj) = 0. _d 0
268 #endif
269 ENDDO
270 ENDDO
271 ENDDO
272 DO j=1-OLy,sNy+OLy
273 DO i=1-OLx,sNx+OLx
274 fVerU (i,j,1) = 0. _d 0
275 fVerU (i,j,2) = 0. _d 0
276 fVerV (i,j,1) = 0. _d 0
277 fVerV (i,j,2) = 0. _d 0
278 phiHydF (i,j) = 0. _d 0
279 phiHydC (i,j) = 0. _d 0
280 dPhiHydX(i,j) = 0. _d 0
281 dPhiHydY(i,j) = 0. _d 0
282 phiSurfX(i,j) = 0. _d 0
283 phiSurfY(i,j) = 0. _d 0
284 ENDDO
285 ENDDO
286
287 C-- Start computation of dynamics
288 iMin = 0
289 iMax = sNx+1
290 jMin = 0
291 jMax = sNy+1
292
293 #ifdef ALLOW_AUTODIFF_TAMC
294 CADJ STORE wvel (:,:,:,bi,bj) =
295 CADJ & comlev1_bibj, key = idynkey, byte = isbyte
296 #endif /* ALLOW_AUTODIFF_TAMC */
297
298 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
299 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
300 IF (implicSurfPress.NE.1.) THEN
301 CALL CALC_GRAD_PHI_SURF(
302 I bi,bj,iMin,iMax,jMin,jMax,
303 I etaN,
304 O phiSurfX,phiSurfY,
305 I myThid )
306 ENDIF
307
308 #ifdef ALLOW_AUTODIFF_TAMC
309 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
310 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
311 #ifdef ALLOW_KPP
312 CADJ STORE KPPviscAz (:,:,:,bi,bj)
313 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
314 #endif /* ALLOW_KPP */
315 #endif /* ALLOW_AUTODIFF_TAMC */
316
317 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
318 C-- Calculate the total vertical diffusivity
319 DO k=1,Nr
320 CALL CALC_VISCOSITY(
321 I bi,bj,iMin,iMax,jMin,jMax,k,
322 O KappaRU,KappaRV,
323 I myThid)
324 ENDDO
325 #endif
326
327 #ifdef ALLOW_AUTODIFF_TAMC
328 CADJ STORE KappaRU(:,:,:)
329 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
330 CADJ STORE KappaRV(:,:,:)
331 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
332 #endif /* ALLOW_AUTODIFF_TAMC */
333
334 C-- Start of dynamics loop
335 DO k=1,Nr
336
337 C-- km1 Points to level above k (=k-1)
338 C-- kup Cycles through 1,2 to point to layer above
339 C-- kDown Cycles through 2,1 to point to current layer
340
341 km1 = MAX(1,k-1)
342 kp1 = MIN(k+1,Nr)
343 kup = 1+MOD(k+1,2)
344 kDown= 1+MOD(k,2)
345
346 #ifdef ALLOW_AUTODIFF_TAMC
347 kkey = (idynkey-1)*Nr + k
348 c
349 CADJ STORE totphihyd (:,:,k,bi,bj)
350 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
351 CADJ STORE gt (:,:,k,bi,bj)
352 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
353 CADJ STORE gs (:,:,k,bi,bj)
354 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
355 CADJ STORE theta (:,:,k,bi,bj)
356 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
357 CADJ STORE salt (:,:,k,bi,bj)
358 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
359 #endif /* ALLOW_AUTODIFF_TAMC */
360
361 C-- Integrate hydrostatic balance for phiHyd with BC of
362 C phiHyd(z=0)=0
363 C distinguishe between Stagger and Non Stagger time stepping
364 IF (staggerTimeStep) THEN
365 CALL CALC_PHI_HYD(
366 I bi,bj,iMin,iMax,jMin,jMax,k,
367 I gT, gS,
368 U phiHydF,
369 O phiHydC, dPhiHydX, dPhiHydY,
370 I myTime, myIter, myThid )
371 ELSE
372 CALL CALC_PHI_HYD(
373 I bi,bj,iMin,iMax,jMin,jMax,k,
374 I theta, salt,
375 U phiHydF,
376 O phiHydC, dPhiHydX, dPhiHydY,
377 I myTime, myIter, myThid )
378 ENDIF
379
380 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
381 C and step forward storing the result in gU, gV, etc...
382 IF ( momStepping ) THEN
383 #ifdef ALLOW_MOM_FLUXFORM
384 IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
385 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
386 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
387 U fVerU, fVerV,
388 I myTime, myIter, myThid)
389 #endif
390 #ifdef ALLOW_MOM_VECINV
391 IF (vectorInvariantMomentum) CALL MOM_VECINV(
392 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
393 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
394 U fVerU, fVerV,
395 I myTime, myIter, myThid)
396 #endif
397 CALL TIMESTEP(
398 I bi,bj,iMin,iMax,jMin,jMax,k,
399 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
400 I myTime, myIter, myThid)
401
402 #ifdef ALLOW_OBCS
403 C-- Apply open boundary conditions
404 IF (useOBCS) THEN
405 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
406 ENDIF
407 #endif /* ALLOW_OBCS */
408
409 ENDIF
410
411
412 C-- end of dynamics k loop (1:Nr)
413 ENDDO
414
415 C-- Implicit Vertical advection & viscosity
416 #ifdef INCLUDE_IMPLVERTADV_CODE
417 IF ( momImplVertAdv ) THEN
418 CALL MOM_U_IMPLICIT_R( kappaRU,
419 I bi, bj, myTime, myIter, myThid )
420 CALL MOM_V_IMPLICIT_R( kappaRV,
421 I bi, bj, myTime, myIter, myThid )
422 ELSEIF ( implicitViscosity ) THEN
423 #else /* INCLUDE_IMPLVERTADV_CODE */
424 IF ( implicitViscosity ) THEN
425 #endif /* INCLUDE_IMPLVERTADV_CODE */
426 #ifdef ALLOW_AUTODIFF_TAMC
427 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
428 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
429 #endif /* ALLOW_AUTODIFF_TAMC */
430 CALL IMPLDIFF(
431 I bi, bj, iMin, iMax, jMin, jMax,
432 I deltaTmom, KappaRU,recip_HFacW,
433 U gU,
434 I myThid )
435 #ifdef ALLOW_AUTODIFF_TAMC
436 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
437 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
438 #endif /* ALLOW_AUTODIFF_TAMC */
439 CALL IMPLDIFF(
440 I bi, bj, iMin, iMax, jMin, jMax,
441 I deltaTmom, KappaRV,recip_HFacS,
442 U gV,
443 I myThid )
444 ENDIF
445
446 #ifdef ALLOW_OBCS
447 C-- Apply open boundary conditions
448 IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
449 DO K=1,Nr
450 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
451 ENDDO
452 ENDIF
453 #endif /* ALLOW_OBCS */
454
455 #ifdef ALLOW_CD_CODE
456 IF (implicitViscosity.AND.useCDscheme) THEN
457 #ifdef ALLOW_AUTODIFF_TAMC
458 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
459 #endif /* ALLOW_AUTODIFF_TAMC */
460 CALL IMPLDIFF(
461 I bi, bj, iMin, iMax, jMin, jMax,
462 I deltaTmom, KappaRU,recip_HFacW,
463 U vVelD,
464 I myThid )
465 #ifdef ALLOW_AUTODIFF_TAMC
466 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
467 #endif /* ALLOW_AUTODIFF_TAMC */
468 CALL IMPLDIFF(
469 I bi, bj, iMin, iMax, jMin, jMax,
470 I deltaTmom, KappaRV,recip_HFacS,
471 U uVelD,
472 I myThid )
473 ENDIF
474 #endif /* ALLOW_CD_CODE */
475 C-- End implicit Vertical advection & viscosity
476
477 ENDDO
478 ENDDO
479
480 Cml(
481 C In order to compare the variance of phiHydLow of a p/z-coordinate
482 C run with etaH of a z/p-coordinate run the drift of phiHydLow
483 C has to be removed by something like the following subroutine:
484 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
485 C & 'phiHydLow', myThid )
486 Cml)
487
488 #ifdef ALLOW_DEBUG
489 If ( debugLevel .GE. debLevB ) THEN
490 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
491 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
492 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
493 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
494 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
495 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
496 CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
497 CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
498 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
499 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
500 CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
501 CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
502 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
503 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
504 ENDIF
505 #endif
506
507 RETURN
508 END

  ViewVC Help
Powered by ViewVC 1.1.22