/[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.110 - (show annotations) (download)
Wed Nov 10 03:02:00 2004 UTC (19 years, 6 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint56b_post, checkpoint57, checkpoint56, checkpoint57a_post, checkpoint55j_post, checkpoint56a_post, checkpoint56c_post, checkpoint57a_pre
Changes since 1.109: +9 -1 lines
isolate dissipation tendency (allow to keep it out off AB)
 note: only implemented in vector-invariant form.

1 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.109 2004/09/23 17:48:24 heimbach 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 TIMEAVE_CUMUL_1T
112 C |-- CALL DEBUG_STATS_RL
113
114 C !INPUT/OUTPUT PARAMETERS:
115 C == Routine arguments ==
116 C myTime - Current time in simulation
117 C myIter - Current iteration number in simulation
118 C myThid - Thread number for this instance of the routine.
119 _RL myTime
120 INTEGER myIter
121 INTEGER myThid
122
123 C !LOCAL VARIABLES:
124 C == Local variables
125 C fVer[STUV] o fVer: Vertical flux term - note fVer
126 C is "pipelined" in the vertical
127 C so we need an fVer for each
128 C variable.
129 C phiHydC :: hydrostatic potential anomaly at cell center
130 C In z coords phiHyd is the hydrostatic potential
131 C (=pressure/rho0) anomaly
132 C In p coords phiHyd is the geopotential height anomaly.
133 C phiHydF :: hydrostatic potential anomaly at middle between 2 centers
134 C dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
135 C phiSurfX, :: gradient of Surface potential (Pressure/rho, ocean)
136 C phiSurfY or geopotential (atmos) in X and Y direction
137 C guDissip :: dissipation tendency (all explicit terms), u component
138 C gvDissip :: dissipation tendency (all explicit terms), v component
139 C iMin, iMax - Ranges and sub-block indices on which calculations
140 C jMin, jMax are applied.
141 C bi, bj
142 C k, kup, - Index for layer above and below. kup and kDown
143 C kDown, km1 are switched with layer to be the appropriate
144 C index into fVerTerm.
145 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
146 _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
147 _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
148 _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
149 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
150 _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
151 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
152 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
153 _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
154 _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
155 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
156 _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
157
158 INTEGER iMin, iMax
159 INTEGER jMin, jMax
160 INTEGER bi, bj
161 INTEGER i, j
162 INTEGER k, km1, kp1, kup, kDown
163
164 LOGICAL DIFFERENT_MULTIPLE
165 EXTERNAL DIFFERENT_MULTIPLE
166
167 C--- The algorithm...
168 C
169 C "Correction Step"
170 C =================
171 C Here we update the horizontal velocities with the surface
172 C pressure such that the resulting flow is either consistent
173 C with the free-surface evolution or the rigid-lid:
174 C U[n] = U* + dt x d/dx P
175 C V[n] = V* + dt x d/dy P
176 C
177 C "Calculation of Gs"
178 C ===================
179 C This is where all the accelerations and tendencies (ie.
180 C physics, parameterizations etc...) are calculated
181 C rho = rho ( theta[n], salt[n] )
182 C b = b(rho, theta)
183 C K31 = K31 ( rho )
184 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
185 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
186 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
187 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
188 C
189 C "Time-stepping" or "Prediction"
190 C ================================
191 C The models variables are stepped forward with the appropriate
192 C time-stepping scheme (currently we use Adams-Bashforth II)
193 C - For momentum, the result is always *only* a "prediction"
194 C in that the flow may be divergent and will be "corrected"
195 C later with a surface pressure gradient.
196 C - Normally for tracers the result is the new field at time
197 C level [n+1} *BUT* in the case of implicit diffusion the result
198 C is also *only* a prediction.
199 C - We denote "predictors" with an asterisk (*).
200 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
201 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
202 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
203 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
204 C With implicit diffusion:
205 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
206 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
207 C (1 + dt * K * d_zz) theta[n] = theta*
208 C (1 + dt * K * d_zz) salt[n] = salt*
209 C---
210 CEOP
211
212 C-- Call to routine for calculation of
213 C Eliassen-Palm-flux-forced U-tendency,
214 C if desired:
215 #ifdef INCLUDE_EP_FORCING_CODE
216 CALL CALC_EP_FORCING(myThid)
217 #endif
218
219 #ifdef ALLOW_AUTODIFF_TAMC
220 C-- HPF directive to help TAMC
221 CHPF$ INDEPENDENT
222 #endif /* ALLOW_AUTODIFF_TAMC */
223
224 DO bj=myByLo(myThid),myByHi(myThid)
225
226 #ifdef ALLOW_AUTODIFF_TAMC
227 C-- HPF directive to help TAMC
228 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
229 CHPF$& ,phiHydF
230 CHPF$& ,KappaRU,KappaRV
231 CHPF$& )
232 #endif /* ALLOW_AUTODIFF_TAMC */
233
234 DO bi=myBxLo(myThid),myBxHi(myThid)
235
236 #ifdef ALLOW_AUTODIFF_TAMC
237 act1 = bi - myBxLo(myThid)
238 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
239 act2 = bj - myByLo(myThid)
240 max2 = myByHi(myThid) - myByLo(myThid) + 1
241 act3 = myThid - 1
242 max3 = nTx*nTy
243 act4 = ikey_dynamics - 1
244 idynkey = (act1 + 1) + act2*max1
245 & + act3*max1*max2
246 & + act4*max1*max2*max3
247 #endif /* ALLOW_AUTODIFF_TAMC */
248
249 C-- Set up work arrays with valid (i.e. not NaN) values
250 C These inital values do not alter the numerical results. They
251 C just ensure that all memory references are to valid floating
252 C point numbers. This prevents spurious hardware signals due to
253 C uninitialised but inert locations.
254
255 DO k=1,Nr
256 DO j=1-OLy,sNy+OLy
257 DO i=1-OLx,sNx+OLx
258 KappaRU(i,j,k) = 0. _d 0
259 KappaRV(i,j,k) = 0. _d 0
260 #ifdef ALLOW_AUTODIFF_TAMC
261 cph(
262 c-- need some re-initialisation here to break dependencies
263 c-- totphihyd is assumed zero from ini_pressure, i.e.
264 c-- avoiding iterate pressure p = integral of (g*rho(p)*dz)
265 cph)
266 totPhiHyd(i,j,k,bi,bj) = 0. _d 0
267 gu(i,j,k,bi,bj) = 0. _d 0
268 gv(i,j,k,bi,bj) = 0. _d 0
269 #endif
270 ENDDO
271 ENDDO
272 ENDDO
273 DO j=1-OLy,sNy+OLy
274 DO i=1-OLx,sNx+OLx
275 fVerU (i,j,1) = 0. _d 0
276 fVerU (i,j,2) = 0. _d 0
277 fVerV (i,j,1) = 0. _d 0
278 fVerV (i,j,2) = 0. _d 0
279 phiHydF (i,j) = 0. _d 0
280 phiHydC (i,j) = 0. _d 0
281 dPhiHydX(i,j) = 0. _d 0
282 dPhiHydY(i,j) = 0. _d 0
283 phiSurfX(i,j) = 0. _d 0
284 phiSurfY(i,j) = 0. _d 0
285 guDissip(i,j) = 0. _d 0
286 gvDissip(i,j) = 0. _d 0
287 ENDDO
288 ENDDO
289
290 C-- Start computation of dynamics
291 iMin = 0
292 iMax = sNx+1
293 jMin = 0
294 jMax = sNy+1
295
296 #ifdef ALLOW_AUTODIFF_TAMC
297 CADJ STORE wvel (:,:,:,bi,bj) =
298 CADJ & comlev1_bibj, key = idynkey, byte = isbyte
299 #endif /* ALLOW_AUTODIFF_TAMC */
300
301 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
302 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
303 IF (implicSurfPress.NE.1.) THEN
304 CALL CALC_GRAD_PHI_SURF(
305 I bi,bj,iMin,iMax,jMin,jMax,
306 I etaN,
307 O phiSurfX,phiSurfY,
308 I myThid )
309 ENDIF
310
311 #ifdef ALLOW_AUTODIFF_TAMC
312 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
313 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
314 #ifdef ALLOW_KPP
315 CADJ STORE KPPviscAz (:,:,:,bi,bj)
316 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
317 #endif /* ALLOW_KPP */
318 #endif /* ALLOW_AUTODIFF_TAMC */
319
320 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
321 C-- Calculate the total vertical diffusivity
322 DO k=1,Nr
323 CALL CALC_VISCOSITY(
324 I bi,bj,iMin,iMax,jMin,jMax,k,
325 O KappaRU,KappaRV,
326 I myThid)
327 ENDDO
328 #endif
329
330 #ifdef ALLOW_AUTODIFF_TAMC
331 CADJ STORE KappaRU(:,:,:)
332 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
333 CADJ STORE KappaRV(:,:,:)
334 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
335 #endif /* ALLOW_AUTODIFF_TAMC */
336
337 C-- Start of dynamics loop
338 DO k=1,Nr
339
340 C-- km1 Points to level above k (=k-1)
341 C-- kup Cycles through 1,2 to point to layer above
342 C-- kDown Cycles through 2,1 to point to current layer
343
344 km1 = MAX(1,k-1)
345 kp1 = MIN(k+1,Nr)
346 kup = 1+MOD(k+1,2)
347 kDown= 1+MOD(k,2)
348
349 #ifdef ALLOW_AUTODIFF_TAMC
350 kkey = (idynkey-1)*Nr + k
351 c
352 CADJ STORE totphihyd (:,:,k,bi,bj)
353 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
354 CADJ STORE theta (:,:,k,bi,bj)
355 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
356 CADJ STORE salt (:,:,k,bi,bj)
357 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
358 #endif /* ALLOW_AUTODIFF_TAMC */
359
360 C-- Integrate hydrostatic balance for phiHyd with BC of
361 C phiHyd(z=0)=0
362 CALL CALC_PHI_HYD(
363 I bi,bj,iMin,iMax,jMin,jMax,k,
364 I theta, salt,
365 U phiHydF,
366 O phiHydC, dPhiHydX, dPhiHydY,
367 I myTime, myIter, myThid )
368
369 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
370 C and step forward storing the result in gU, gV, etc...
371 IF ( momStepping ) THEN
372 #ifdef ALLOW_MOM_FLUXFORM
373 IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
374 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
375 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
376 U fVerU, fVerV,
377 I myTime, myIter, myThid)
378 #endif
379 #ifdef ALLOW_MOM_VECINV
380 IF (vectorInvariantMomentum) CALL MOM_VECINV(
381 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
382 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
383 U fVerU, fVerV,
384 O guDissip, gvDissip,
385 I myTime, myIter, myThid)
386 #endif
387 CALL TIMESTEP(
388 I bi,bj,iMin,iMax,jMin,jMax,k,
389 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
390 I guDissip, gvDissip,
391 I myTime, myIter, myThid)
392
393 #ifdef ALLOW_OBCS
394 C-- Apply open boundary conditions
395 IF (useOBCS) THEN
396 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
397 ENDIF
398 #endif /* ALLOW_OBCS */
399
400 ENDIF
401
402
403 C-- end of dynamics k loop (1:Nr)
404 ENDDO
405
406 C-- Implicit Vertical advection & viscosity
407 #ifdef INCLUDE_IMPLVERTADV_CODE
408 IF ( momImplVertAdv ) THEN
409 CALL MOM_U_IMPLICIT_R( kappaRU,
410 I bi, bj, myTime, myIter, myThid )
411 CALL MOM_V_IMPLICIT_R( kappaRV,
412 I bi, bj, myTime, myIter, myThid )
413 ELSEIF ( implicitViscosity ) THEN
414 #else /* INCLUDE_IMPLVERTADV_CODE */
415 IF ( implicitViscosity ) THEN
416 #endif /* INCLUDE_IMPLVERTADV_CODE */
417 #ifdef ALLOW_AUTODIFF_TAMC
418 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
419 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
420 #endif /* ALLOW_AUTODIFF_TAMC */
421 CALL IMPLDIFF(
422 I bi, bj, iMin, iMax, jMin, jMax,
423 I deltaTmom, KappaRU,recip_HFacW,
424 U gU,
425 I myThid )
426 #ifdef ALLOW_AUTODIFF_TAMC
427 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
428 CADJ STORE gV(:,:,:,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, KappaRV,recip_HFacS,
433 U gV,
434 I myThid )
435 ENDIF
436
437 #ifdef ALLOW_OBCS
438 C-- Apply open boundary conditions
439 IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
440 DO K=1,Nr
441 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
442 ENDDO
443 ENDIF
444 #endif /* ALLOW_OBCS */
445
446 #ifdef ALLOW_CD_CODE
447 IF (implicitViscosity.AND.useCDscheme) THEN
448 #ifdef ALLOW_AUTODIFF_TAMC
449 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
450 #endif /* ALLOW_AUTODIFF_TAMC */
451 CALL IMPLDIFF(
452 I bi, bj, iMin, iMax, jMin, jMax,
453 I deltaTmom, KappaRU,recip_HFacW,
454 U vVelD,
455 I myThid )
456 #ifdef ALLOW_AUTODIFF_TAMC
457 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
458 #endif /* ALLOW_AUTODIFF_TAMC */
459 CALL IMPLDIFF(
460 I bi, bj, iMin, iMax, jMin, jMax,
461 I deltaTmom, KappaRV,recip_HFacS,
462 U uVelD,
463 I myThid )
464 ENDIF
465 #endif /* ALLOW_CD_CODE */
466 C-- End implicit Vertical advection & viscosity
467
468 ENDDO
469 ENDDO
470
471 #ifdef ALLOW_OBCS
472 IF (useOBCS) THEN
473 CALL OBCS_PRESCRIBE_EXCHANGES(myThid)
474 ENDIF
475 #endif
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_DEBUG
486 If ( debugLevel .GE. debLevB ) THEN
487 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
488 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
489 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
490 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
491 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
492 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
493 CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
494 CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
495 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
496 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
497 CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
498 CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
499 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
500 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
501 ENDIF
502 #endif
503
504 RETURN
505 END

  ViewVC Help
Powered by ViewVC 1.1.22