/[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.89 - (hide annotations) (download)
Wed Aug 7 16:55:52 2002 UTC (21 years, 9 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint46g_pre, checkpoint46f_post, checkpoint46b_post, checkpoint46d_pre, checkpoint46e_pre, checkpoint46c_pre, checkpoint46c_post, checkpoint46e_post, checkpoint46d_post
Changes since 1.88: +8 -1 lines
o Added new equation of state -> JMD95Z and JMD95P
  - EOS of Jackett and McDougall, 1995, JPO
  - moved all EOS parameters into EOS.h
  - new routines ini_eos.F, store_pressure.F
o Added UNESCO EOS, but not recommended because it requires
  in-situ temperature (see JMD95)
o Modified formatting for knudsen2.f in utils/knudsen2 and added
  unesco.f to be used with POLY3

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

  ViewVC Help
Powered by ViewVC 1.1.22