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

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

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


Revision 1.83.4.2 - (hide annotations) (download)
Sun Mar 24 17:25:19 2002 UTC (22 years, 3 months ago) by heimbach
Branch: ecco-branch
CVS Tags: ecco_c44_e19, ecco_c44_e22, ecco_c44_e20, ecco_c44_e21
Changes since 1.83.4.1: +0 -5 lines
Merged changes to enable stable adjoint of GM from release1_final.

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

  ViewVC Help
Powered by ViewVC 1.1.22