/[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.113 - (show annotations) (download)
Fri Jan 28 01:00:13 2005 UTC (19 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57d_post, checkpoint57e_post, eckpoint57e_pre, checkpoint57f_pre
Changes since 1.112: +36 -6 lines
move state variable diagnostics to the beginning of the time step.

1 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.112 2005/01/24 17:00:17 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 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 LOGICAL DIFFERENT_MULTIPLE
164 EXTERNAL DIFFERENT_MULTIPLE
165
166 #ifdef ALLOW_DIAGNOSTICS
167 _RL tmpFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
168 LOGICAL DIAGNOSTICS_IS_ON
169 EXTERNAL DIAGNOSTICS_IS_ON
170 #endif /* ALLOW_DIAGNOSTICS */
171
172
173 C--- The algorithm...
174 C
175 C "Correction Step"
176 C =================
177 C Here we update the horizontal velocities with the surface
178 C pressure such that the resulting flow is either consistent
179 C with the free-surface evolution or the rigid-lid:
180 C U[n] = U* + dt x d/dx P
181 C V[n] = V* + dt x d/dy P
182 C
183 C "Calculation of Gs"
184 C ===================
185 C This is where all the accelerations and tendencies (ie.
186 C physics, parameterizations etc...) are calculated
187 C rho = rho ( theta[n], salt[n] )
188 C b = b(rho, theta)
189 C K31 = K31 ( rho )
190 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
191 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
192 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
193 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
194 C
195 C "Time-stepping" or "Prediction"
196 C ================================
197 C The models variables are stepped forward with the appropriate
198 C time-stepping scheme (currently we use Adams-Bashforth II)
199 C - For momentum, the result is always *only* a "prediction"
200 C in that the flow may be divergent and will be "corrected"
201 C later with a surface pressure gradient.
202 C - Normally for tracers the result is the new field at time
203 C level [n+1} *BUT* in the case of implicit diffusion the result
204 C is also *only* a prediction.
205 C - We denote "predictors" with an asterisk (*).
206 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
207 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
208 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
209 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
210 C With implicit diffusion:
211 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
212 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
213 C (1 + dt * K * d_zz) theta[n] = theta*
214 C (1 + dt * K * d_zz) salt[n] = salt*
215 C---
216 CEOP
217
218 C-- Call to routine for calculation of
219 C Eliassen-Palm-flux-forced U-tendency,
220 C if desired:
221 #ifdef INCLUDE_EP_FORCING_CODE
222 CALL CALC_EP_FORCING(myThid)
223 #endif
224
225 #ifdef ALLOW_AUTODIFF_TAMC
226 C-- HPF directive to help TAMC
227 CHPF$ INDEPENDENT
228 #endif /* ALLOW_AUTODIFF_TAMC */
229
230 DO bj=myByLo(myThid),myByHi(myThid)
231
232 #ifdef ALLOW_AUTODIFF_TAMC
233 C-- HPF directive to help TAMC
234 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
235 CHPF$& ,phiHydF
236 CHPF$& ,KappaRU,KappaRV
237 CHPF$& )
238 #endif /* ALLOW_AUTODIFF_TAMC */
239
240 DO bi=myBxLo(myThid),myBxHi(myThid)
241
242 #ifdef ALLOW_AUTODIFF_TAMC
243 act1 = bi - myBxLo(myThid)
244 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
245 act2 = bj - myByLo(myThid)
246 max2 = myByHi(myThid) - myByLo(myThid) + 1
247 act3 = myThid - 1
248 max3 = nTx*nTy
249 act4 = ikey_dynamics - 1
250 idynkey = (act1 + 1) + act2*max1
251 & + act3*max1*max2
252 & + act4*max1*max2*max3
253 #endif /* ALLOW_AUTODIFF_TAMC */
254
255 C-- Set up work arrays with valid (i.e. not NaN) values
256 C These inital values do not alter the numerical results. They
257 C just ensure that all memory references are to valid floating
258 C point numbers. This prevents spurious hardware signals due to
259 C uninitialised but inert locations.
260
261 DO k=1,Nr
262 DO j=1-OLy,sNy+OLy
263 DO i=1-OLx,sNx+OLx
264 KappaRU(i,j,k) = 0. _d 0
265 KappaRV(i,j,k) = 0. _d 0
266 #ifdef ALLOW_AUTODIFF_TAMC
267 cph(
268 c-- need some re-initialisation here to break dependencies
269 cph)
270 gu(i,j,k,bi,bj) = 0. _d 0
271 gv(i,j,k,bi,bj) = 0. _d 0
272 #endif
273 ENDDO
274 ENDDO
275 ENDDO
276 DO j=1-OLy,sNy+OLy
277 DO i=1-OLx,sNx+OLx
278 fVerU (i,j,1) = 0. _d 0
279 fVerU (i,j,2) = 0. _d 0
280 fVerV (i,j,1) = 0. _d 0
281 fVerV (i,j,2) = 0. _d 0
282 phiHydF (i,j) = 0. _d 0
283 phiHydC (i,j) = 0. _d 0
284 dPhiHydX(i,j) = 0. _d 0
285 dPhiHydY(i,j) = 0. _d 0
286 phiSurfX(i,j) = 0. _d 0
287 phiSurfY(i,j) = 0. _d 0
288 guDissip(i,j) = 0. _d 0
289 gvDissip(i,j) = 0. _d 0
290 ENDDO
291 ENDDO
292
293 C-- Start computation of dynamics
294 iMin = 0
295 iMax = sNx+1
296 jMin = 0
297 jMax = sNy+1
298
299 #ifdef ALLOW_AUTODIFF_TAMC
300 CADJ STORE wvel (:,:,:,bi,bj) =
301 CADJ & comlev1_bibj, key = idynkey, byte = isbyte
302 #endif /* ALLOW_AUTODIFF_TAMC */
303
304 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
305 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
306 IF (implicSurfPress.NE.1.) THEN
307 CALL CALC_GRAD_PHI_SURF(
308 I bi,bj,iMin,iMax,jMin,jMax,
309 I etaN,
310 O phiSurfX,phiSurfY,
311 I myThid )
312 ENDIF
313
314 #ifdef ALLOW_AUTODIFF_TAMC
315 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
316 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
317 #ifdef ALLOW_KPP
318 CADJ STORE KPPviscAz (:,:,:,bi,bj)
319 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
320 #endif /* ALLOW_KPP */
321 #endif /* ALLOW_AUTODIFF_TAMC */
322
323 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
324 C-- Calculate the total vertical diffusivity
325 DO k=1,Nr
326 CALL CALC_VISCOSITY(
327 I bi,bj,iMin,iMax,jMin,jMax,k,
328 O KappaRU,KappaRV,
329 I myThid)
330 ENDDO
331 #endif
332
333 #ifdef ALLOW_AUTODIFF_TAMC
334 CADJ STORE KappaRU(:,:,:)
335 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
336 CADJ STORE KappaRV(:,:,:)
337 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
338 #endif /* ALLOW_AUTODIFF_TAMC */
339
340 C-- Start of dynamics loop
341 DO k=1,Nr
342
343 C-- km1 Points to level above k (=k-1)
344 C-- kup Cycles through 1,2 to point to layer above
345 C-- kDown Cycles through 2,1 to point to current layer
346
347 km1 = MAX(1,k-1)
348 kp1 = MIN(k+1,Nr)
349 kup = 1+MOD(k+1,2)
350 kDown= 1+MOD(k,2)
351
352 #ifdef ALLOW_AUTODIFF_TAMC
353 kkey = (idynkey-1)*Nr + k
354 c
355 CADJ STORE totphihyd (:,:,k,bi,bj)
356 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
357 CADJ STORE theta (:,:,k,bi,bj)
358 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
359 CADJ STORE salt (:,:,k,bi,bj)
360 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
361 #endif /* ALLOW_AUTODIFF_TAMC */
362
363 C-- Integrate hydrostatic balance for phiHyd with BC of
364 C phiHyd(z=0)=0
365 CALL CALC_PHI_HYD(
366 I bi,bj,iMin,iMax,jMin,jMax,k,
367 I theta, salt,
368 U phiHydF,
369 O phiHydC, dPhiHydX, dPhiHydY,
370 I myTime, myIter, myThid )
371
372 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
373 C and step forward storing the result in gU, gV, etc...
374 IF ( momStepping ) THEN
375 #ifdef ALLOW_MOM_FLUXFORM
376 IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
377 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
378 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
379 U fVerU, fVerV,
380 I myTime, myIter, myThid)
381 #endif
382 #ifdef ALLOW_MOM_VECINV
383 IF (vectorInvariantMomentum) CALL MOM_VECINV(
384 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
385 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
386 U fVerU, fVerV,
387 O guDissip, gvDissip,
388 I myTime, myIter, myThid)
389 #endif
390 CALL TIMESTEP(
391 I bi,bj,iMin,iMax,jMin,jMax,k,
392 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
393 I guDissip, gvDissip,
394 I myTime, myIter, myThid)
395
396 #ifdef ALLOW_OBCS
397 C-- Apply open boundary conditions
398 IF (useOBCS) THEN
399 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
400 ENDIF
401 #endif /* ALLOW_OBCS */
402
403 ENDIF
404
405
406 C-- end of dynamics k loop (1:Nr)
407 ENDDO
408
409 C-- Implicit Vertical advection & viscosity
410 #ifdef INCLUDE_IMPLVERTADV_CODE
411 IF ( momImplVertAdv ) THEN
412 CALL MOM_U_IMPLICIT_R( kappaRU,
413 I bi, bj, myTime, myIter, myThid )
414 CALL MOM_V_IMPLICIT_R( kappaRV,
415 I bi, bj, myTime, myIter, myThid )
416 ELSEIF ( implicitViscosity ) THEN
417 #else /* INCLUDE_IMPLVERTADV_CODE */
418 IF ( implicitViscosity ) THEN
419 #endif /* INCLUDE_IMPLVERTADV_CODE */
420 #ifdef ALLOW_AUTODIFF_TAMC
421 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
422 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
423 #endif /* ALLOW_AUTODIFF_TAMC */
424 CALL IMPLDIFF(
425 I bi, bj, iMin, iMax, jMin, jMax,
426 I 0, KappaRU,recip_HFacW,
427 U gU,
428 I myThid )
429 #ifdef ALLOW_AUTODIFF_TAMC
430 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
431 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
432 #endif /* ALLOW_AUTODIFF_TAMC */
433 CALL IMPLDIFF(
434 I bi, bj, iMin, iMax, jMin, jMax,
435 I 0, KappaRV,recip_HFacS,
436 U gV,
437 I myThid )
438 ENDIF
439
440 #ifdef ALLOW_OBCS
441 C-- Apply open boundary conditions
442 IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
443 DO K=1,Nr
444 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
445 ENDDO
446 ENDIF
447 #endif /* ALLOW_OBCS */
448
449 #ifdef ALLOW_CD_CODE
450 IF (implicitViscosity.AND.useCDscheme) THEN
451 #ifdef ALLOW_AUTODIFF_TAMC
452 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
453 #endif /* ALLOW_AUTODIFF_TAMC */
454 CALL IMPLDIFF(
455 I bi, bj, iMin, iMax, jMin, jMax,
456 I 0, KappaRU,recip_HFacW,
457 U vVelD,
458 I myThid )
459 #ifdef ALLOW_AUTODIFF_TAMC
460 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
461 #endif /* ALLOW_AUTODIFF_TAMC */
462 CALL IMPLDIFF(
463 I bi, bj, iMin, iMax, jMin, jMax,
464 I 0, KappaRV,recip_HFacS,
465 U uVelD,
466 I myThid )
467 ENDIF
468 #endif /* ALLOW_CD_CODE */
469 C-- End implicit Vertical advection & viscosity
470
471 ENDDO
472 ENDDO
473
474 #ifdef ALLOW_OBCS
475 IF (useOBCS) THEN
476 CALL OBCS_PRESCRIBE_EXCHANGES(myThid)
477 ENDIF
478 #endif
479
480 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
481
482 Cml(
483 C In order to compare the variance of phiHydLow of a p/z-coordinate
484 C run with etaH of a z/p-coordinate run the drift of phiHydLow
485 C has to be removed by something like the following subroutine:
486 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
487 C & 'phiHydLow', myThid )
488 Cml)
489
490 #ifdef ALLOW_DIAGNOSTICS
491 IF ( usediagnostics ) THEN
492
493 CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
494 CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0,1,0,1,1,myThid)
495
496 IF ( DIAGNOSTICS_IS_ON('PHIBOTSQ',myThid) ) THEN
497 DO bj = myByLo(myThid), myByHi(myThid)
498 DO bi = myBxLo(myThid), myBxHi(myThid)
499 DO j = 1,sNy
500 DO i = 1,sNx
501 tmpFld(i,j) = phiHydLow(i,j,bi,bj)*phiHydLow(i,j,bi,bj)
502 ENDDO
503 ENDDO
504 CALL DIAGNOSTICS_FILL(tmpFld,'PHIBOTSQ',0,1,2,bi,bj,myThid)
505 ENDDO
506 ENDDO
507 ENDIF
508
509 ENDIF
510 #endif /* ALLOW_DIAGNOSTICS */
511
512 #ifdef ALLOW_DEBUG
513 If ( debugLevel .GE. debLevB ) THEN
514 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
515 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
516 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
517 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
518 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
519 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
520 CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
521 CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
522 CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
523 CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
524 CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
525 CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
526 CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
527 CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
528 ENDIF
529 #endif
530
531 RETURN
532 END

  ViewVC Help
Powered by ViewVC 1.1.22