/[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.93 - (hide annotations) (download)
Tue Feb 11 04:05:32 2003 UTC (21 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint48e_post, checkpoint48f_post
Changes since 1.92: +5 -5 lines
dynamics: change definition of computational domain & adapt mom_fluxform
 accordingly ; when viscA4=0, allows to run the dynamics with Olx=Oly=2.

1 jmc 1.93 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.92 2003/02/08 02:09:20 jmc 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 heimbach 1.91 # include "EOS.h"
84 heimbach 1.67 # ifdef ALLOW_KPP
85     # include "KPP.h"
86     # endif
87 heimbach 1.53 #endif /* ALLOW_AUTODIFF_TAMC */
88 jmc 1.64 #ifdef ALLOW_TIMEAVE
89     #include "TIMEAVE_STATV.h"
90 jmc 1.62 #endif
91    
92 cnh 1.82 C !CALLING SEQUENCE:
93     C DYNAMICS()
94     C |
95     C |-- CALC_GRAD_PHI_SURF
96     C |
97     C |-- CALC_VISCOSITY
98     C |
99     C |-- CALC_PHI_HYD
100     C |
101 mlosch 1.89 C |-- STORE_PRESSURE
102     C |
103 cnh 1.82 C |-- MOM_FLUXFORM
104     C |
105     C |-- MOM_VECINV
106     C |
107     C |-- TIMESTEP
108     C |
109     C |-- OBCS_APPLY_UV
110     C |
111     C |-- IMPLDIFF
112     C |
113     C |-- OBCS_APPLY_UV
114     C |
115     C |-- CALL TIMEAVE_CUMUL_1T
116     C |-- CALL DEBUG_STATS_RL
117    
118     C !INPUT/OUTPUT PARAMETERS:
119 cnh 1.1 C == Routine arguments ==
120 cnh 1.8 C myTime - Current time in simulation
121     C myIter - Current iteration number in simulation
122 cnh 1.1 C myThid - Thread number for this instance of the routine.
123 cnh 1.8 _RL myTime
124     INTEGER myIter
125 adcroft 1.47 INTEGER myThid
126 cnh 1.1
127 cnh 1.82 C !LOCAL VARIABLES:
128 cnh 1.1 C == Local variables
129 adcroft 1.58 C fVer[STUV] o fVer: Vertical flux term - note fVer
130 cnh 1.1 C is "pipelined" in the vertical
131     C so we need an fVer for each
132     C variable.
133 adcroft 1.58 C rhoK, rhoKM1 - Density at current level, and level above
134 jmc 1.92 C phiHyd - Hydrostatic part of the potential.
135     C In z coords phiHyd is the hydrostatic
136 jmc 1.65 C Potential (=pressure/rho0) anomaly
137 jmc 1.92 C In p coords phiHyd is the geopotential
138 jmc 1.65 C surface height anomaly.
139 jmc 1.92 C dPhiHydX,Y :: Gradient (X & Y directions) of Hydrostatic Potential
140     C phiSurfX, - gradient of Surface potential (Pressure/rho, ocean)
141     C phiSurfY or geopotential (atmos) in X and Y direction
142 cnh 1.30 C iMin, iMax - Ranges and sub-block indices on which calculations
143     C jMin, jMax are applied.
144 cnh 1.1 C bi, bj
145 heimbach 1.53 C k, kup, - Index for layer above and below. kup and kDown
146     C kDown, km1 are switched with layer to be the appropriate
147 cnh 1.38 C index into fVerTerm.
148 cnh 1.30 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
149     _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
150 cnh 1.31 _RL phiHyd (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
151 jmc 1.92 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
152     _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
153 cnh 1.30 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
154     _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
155 jmc 1.63 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
156     _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
157 adcroft 1.42 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
158     _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
159 adcroft 1.12
160 cnh 1.1 INTEGER iMin, iMax
161     INTEGER jMin, jMax
162     INTEGER bi, bj
163     INTEGER i, j
164 heimbach 1.77 INTEGER k, km1, kp1, kup, kDown
165 cnh 1.1
166 jmc 1.92 LOGICAL DIFFERENT_MULTIPLE
167     EXTERNAL DIFFERENT_MULTIPLE
168 jmc 1.62
169 adcroft 1.11 C--- The algorithm...
170     C
171     C "Correction Step"
172     C =================
173     C Here we update the horizontal velocities with the surface
174     C pressure such that the resulting flow is either consistent
175     C with the free-surface evolution or the rigid-lid:
176     C U[n] = U* + dt x d/dx P
177     C V[n] = V* + dt x d/dy P
178     C
179     C "Calculation of Gs"
180     C ===================
181     C This is where all the accelerations and tendencies (ie.
182 heimbach 1.53 C physics, parameterizations etc...) are calculated
183 adcroft 1.11 C rho = rho ( theta[n], salt[n] )
184 cnh 1.27 C b = b(rho, theta)
185 adcroft 1.11 C K31 = K31 ( rho )
186 jmc 1.61 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
187     C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
188     C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
189     C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
190 adcroft 1.11 C
191 adcroft 1.12 C "Time-stepping" or "Prediction"
192 adcroft 1.11 C ================================
193     C The models variables are stepped forward with the appropriate
194     C time-stepping scheme (currently we use Adams-Bashforth II)
195     C - For momentum, the result is always *only* a "prediction"
196     C in that the flow may be divergent and will be "corrected"
197     C later with a surface pressure gradient.
198     C - Normally for tracers the result is the new field at time
199     C level [n+1} *BUT* in the case of implicit diffusion the result
200     C is also *only* a prediction.
201     C - We denote "predictors" with an asterisk (*).
202     C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
203     C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
204     C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
205     C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
206 adcroft 1.12 C With implicit diffusion:
207 adcroft 1.11 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
208     C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
209 adcroft 1.12 C (1 + dt * K * d_zz) theta[n] = theta*
210     C (1 + dt * K * d_zz) salt[n] = salt*
211 adcroft 1.11 C---
212 cnh 1.82 CEOP
213 adcroft 1.11
214 heimbach 1.76 C-- Set up work arrays with valid (i.e. not NaN) values
215     C These inital values do not alter the numerical results. They
216     C just ensure that all memory references are to valid floating
217     C point numbers. This prevents spurious hardware signals due to
218     C uninitialised but inert locations.
219     DO j=1-OLy,sNy+OLy
220     DO i=1-OLx,sNx+OLx
221     rhoKM1 (i,j) = 0. _d 0
222     rhok (i,j) = 0. _d 0
223     phiSurfX(i,j) = 0. _d 0
224     phiSurfY(i,j) = 0. _d 0
225     ENDDO
226     ENDDO
227    
228 heimbach 1.88 C-- Call to routine for calculation of
229     C Eliassen-Palm-flux-forced U-tendency,
230     C if desired:
231     #ifdef INCLUDE_EP_FORCING_CODE
232     CALL CALC_EP_FORCING(myThid)
233     #endif
234    
235 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
236     C-- HPF directive to help TAMC
237     CHPF$ INDEPENDENT
238     #endif /* ALLOW_AUTODIFF_TAMC */
239    
240 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
241 heimbach 1.76
242     #ifdef ALLOW_AUTODIFF_TAMC
243     C-- HPF directive to help TAMC
244     CHPF$ INDEPENDENT, NEW (fVerU,fVerV
245     CHPF$& ,phiHyd
246     CHPF$& ,KappaRU,KappaRV
247     CHPF$& )
248     #endif /* ALLOW_AUTODIFF_TAMC */
249    
250 cnh 1.1 DO bi=myBxLo(myThid),myBxHi(myThid)
251 heimbach 1.76
252     #ifdef ALLOW_AUTODIFF_TAMC
253     act1 = bi - myBxLo(myThid)
254     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
255     act2 = bj - myByLo(myThid)
256     max2 = myByHi(myThid) - myByLo(myThid) + 1
257     act3 = myThid - 1
258     max3 = nTx*nTy
259     act4 = ikey_dynamics - 1
260 heimbach 1.91 idynkey = (act1 + 1) + act2*max1
261 heimbach 1.76 & + act3*max1*max2
262     & + act4*max1*max2*max3
263     #endif /* ALLOW_AUTODIFF_TAMC */
264    
265     C-- Set up work arrays that need valid initial values
266     DO j=1-OLy,sNy+OLy
267     DO i=1-OLx,sNx+OLx
268 heimbach 1.87 DO k=1,Nr
269     phiHyd(i,j,k) = 0. _d 0
270     KappaRU(i,j,k) = 0. _d 0
271     KappaRV(i,j,k) = 0. _d 0
272     ENDDO
273 heimbach 1.76 fVerU (i,j,1) = 0. _d 0
274     fVerU (i,j,2) = 0. _d 0
275     fVerV (i,j,1) = 0. _d 0
276     fVerV (i,j,2) = 0. _d 0
277 jmc 1.92 dPhiHydX(i,j) = 0. _d 0
278     dPhiHydY(i,j) = 0. _d 0
279 heimbach 1.76 ENDDO
280     ENDDO
281 heimbach 1.49
282 jmc 1.63 C-- Start computation of dynamics
283 jmc 1.93 iMin = 0
284     iMax = sNx+1
285     jMin = 0
286     jMax = sNy+1
287 jmc 1.63
288 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
289 heimbach 1.91 CADJ STORE wvel (:,:,:,bi,bj) =
290     CADJ & comlev1_bibj, key = idynkey, byte = isbyte
291 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
292    
293 jmc 1.65 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
294 jmc 1.63 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
295     IF (implicSurfPress.NE.1.) THEN
296 jmc 1.65 CALL CALC_GRAD_PHI_SURF(
297     I bi,bj,iMin,iMax,jMin,jMax,
298     I etaN,
299     O phiSurfX,phiSurfY,
300     I myThid )
301 jmc 1.63 ENDIF
302 heimbach 1.83
303     #ifdef ALLOW_AUTODIFF_TAMC
304 heimbach 1.91 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
305     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
306 heimbach 1.83 #ifdef ALLOW_KPP
307     CADJ STORE KPPviscAz (:,:,:,bi,bj)
308 heimbach 1.91 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
309 heimbach 1.83 #endif /* ALLOW_KPP */
310     #endif /* ALLOW_AUTODIFF_TAMC */
311 adcroft 1.58
312 heimbach 1.77 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
313     C-- Calculate the total vertical diffusivity
314     DO k=1,Nr
315     CALL CALC_VISCOSITY(
316     I bi,bj,iMin,iMax,jMin,jMax,k,
317     O KappaRU,KappaRV,
318     I myThid)
319     ENDDO
320     #endif
321    
322 adcroft 1.58 C-- Start of dynamics loop
323     DO k=1,Nr
324    
325     C-- km1 Points to level above k (=k-1)
326     C-- kup Cycles through 1,2 to point to layer above
327     C-- kDown Cycles through 2,1 to point to current layer
328    
329     km1 = MAX(1,k-1)
330 heimbach 1.77 kp1 = MIN(k+1,Nr)
331 adcroft 1.58 kup = 1+MOD(k+1,2)
332     kDown= 1+MOD(k,2)
333    
334 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
335 heimbach 1.91 kkey = (idynkey-1)*Nr + k
336     CADJ STORE pressure(:,:,k,bi,bj) = comlev1_bibj_k ,
337     CADJ & key=kkey , byte=isbyte
338 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
339    
340 adcroft 1.58 C-- Integrate hydrostatic balance for phiHyd with BC of
341     C phiHyd(z=0)=0
342     C distinguishe between Stagger and Non Stagger time stepping
343     IF (staggerTimeStep) THEN
344     CALL CALC_PHI_HYD(
345     I bi,bj,iMin,iMax,jMin,jMax,k,
346 adcroft 1.81 I gT, gS,
347 adcroft 1.58 U phiHyd,
348 jmc 1.92 O dPhiHydX, dPhiHydY,
349     I myTime, myIter, myThid )
350 adcroft 1.58 ELSE
351     CALL CALC_PHI_HYD(
352     I bi,bj,iMin,iMax,jMin,jMax,k,
353     I theta, salt,
354     U phiHyd,
355 jmc 1.92 O dPhiHydX, dPhiHydY,
356     I myTime, myIter, myThid )
357 adcroft 1.58 ENDIF
358 mlosch 1.89
359     C calculate pressure from phiHyd and store it on common block
360     C variable pressure
361     CALL STORE_PRESSURE( bi, bj, k, phiHyd, myThid )
362    
363 adcroft 1.58 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
364     C and step forward storing the result in gUnm1, gVnm1, etc...
365     IF ( momStepping ) THEN
366 adcroft 1.79 #ifndef DISABLE_MOM_FLUXFORM
367     IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
368 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
369 jmc 1.92 I phiHyd,dPhiHydX,dPhiHydY,KappaRU,KappaRV,
370 adcroft 1.58 U fVerU, fVerV,
371 adcroft 1.80 I myTime, myIter, myThid)
372 adcroft 1.79 #endif
373     #ifndef DISABLE_MOM_VECINV
374     IF (vectorInvariantMomentum) CALL MOM_VECINV(
375     I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
376 jmc 1.92 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
377 adcroft 1.79 U fVerU, fVerV,
378 adcroft 1.80 I myTime, myIter, myThid)
379 adcroft 1.79 #endif
380 adcroft 1.58 CALL TIMESTEP(
381 jmc 1.63 I bi,bj,iMin,iMax,jMin,jMax,k,
382 jmc 1.92 I phiHyd, dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
383 adcroft 1.58 I myIter, myThid)
384    
385     #ifdef ALLOW_OBCS
386     C-- Apply open boundary conditions
387     IF (useOBCS) THEN
388     CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
389     END IF
390     #endif /* ALLOW_OBCS */
391    
392     #ifdef ALLOW_AUTODIFF_TAMC
393     #ifdef INCLUDE_CD_CODE
394     ELSE
395     DO j=1-OLy,sNy+OLy
396     DO i=1-OLx,sNx+OLx
397     guCD(i,j,k,bi,bj) = 0.0
398     gvCD(i,j,k,bi,bj) = 0.0
399     END DO
400     END DO
401     #endif /* INCLUDE_CD_CODE */
402     #endif /* ALLOW_AUTODIFF_TAMC */
403     ENDIF
404    
405    
406     C-- end of dynamics k loop (1:Nr)
407     ENDDO
408    
409 adcroft 1.44 C-- Implicit viscosity
410 adcroft 1.58 IF (implicitViscosity.AND.momStepping) THEN
411     #ifdef ALLOW_AUTODIFF_TAMC
412 heimbach 1.91 CADJ STORE gUNm1(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
413 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
414 adcroft 1.42 CALL IMPLDIFF(
415     I bi, bj, iMin, iMax, jMin, jMax,
416     I deltaTmom, KappaRU,recip_HFacW,
417     U gUNm1,
418     I myThid )
419 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
420 heimbach 1.91 CADJ STORE gVNm1(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
421 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
422 adcroft 1.42 CALL IMPLDIFF(
423     I bi, bj, iMin, iMax, jMin, jMax,
424     I deltaTmom, KappaRV,recip_HFacS,
425     U gVNm1,
426     I myThid )
427 heimbach 1.49
428 adcroft 1.58 #ifdef ALLOW_OBCS
429     C-- Apply open boundary conditions
430     IF (useOBCS) THEN
431     DO K=1,Nr
432     CALL OBCS_APPLY_UV( bi, bj, k, gUnm1, gVnm1, myThid )
433     ENDDO
434     END IF
435     #endif /* ALLOW_OBCS */
436 heimbach 1.49
437 adcroft 1.58 #ifdef INCLUDE_CD_CODE
438     #ifdef ALLOW_AUTODIFF_TAMC
439 heimbach 1.91 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
440 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
441 adcroft 1.42 CALL IMPLDIFF(
442     I bi, bj, iMin, iMax, jMin, jMax,
443     I deltaTmom, KappaRU,recip_HFacW,
444     U vVelD,
445     I myThid )
446 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
447 heimbach 1.91 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
448 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
449 adcroft 1.42 CALL IMPLDIFF(
450     I bi, bj, iMin, iMax, jMin, jMax,
451     I deltaTmom, KappaRV,recip_HFacS,
452     U uVelD,
453     I myThid )
454 adcroft 1.58 #endif /* INCLUDE_CD_CODE */
455     C-- End If implicitViscosity.AND.momStepping
456 heimbach 1.53 ENDIF
457 cnh 1.1
458 jmc 1.92 C- jmc: add for diagnostic of phiHyd
459     IF ( DIFFERENT_MULTIPLE(diagFreq,myTime+deltaTClock,myTime)
460     & .AND. buoyancyRelation .NE. 'OCEANIC' ) THEN
461     CALL WRITE_LOCAL_RL('Ph','I10',Nr,phiHyd,
462     & bi,bj,1,myIter+1,myThid)
463     ENDIF
464 jmc 1.62
465 jmc 1.64 #ifdef ALLOW_TIMEAVE
466 jmc 1.62 IF (taveFreq.GT.0.) THEN
467 adcroft 1.68 CALL TIMEAVE_CUMUL_1T(phiHydtave, phiHyd, Nr,
468 jmc 1.64 I deltaTclock, bi, bj, myThid)
469 jmc 1.62 ENDIF
470 jmc 1.64 #endif /* ALLOW_TIMEAVE */
471 jmc 1.62
472 cnh 1.1 ENDDO
473     ENDDO
474 mlosch 1.90
475     Cml(
476     C In order to compare the variance of phiHydLow of a p/z-coordinate
477     C run with etaH of a z/p-coordinate run the drift of phiHydLow
478     C has to be removed by something like the following subroutine:
479     C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
480     C & 'phiHydLow', myThid )
481     Cml)
482 adcroft 1.69
483 adcroft 1.79 #ifndef DISABLE_DEBUGMODE
484 adcroft 1.70 If (debugMode) THEN
485 adcroft 1.69 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
486 adcroft 1.73 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
487 adcroft 1.69 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
488     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
489     CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
490     CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
491     CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
492     CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
493     CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
494     CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
495     CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
496     CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
497     CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
498     CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
499 adcroft 1.70 ENDIF
500 adcroft 1.69 #endif
501 cnh 1.1
502     RETURN
503     END

  ViewVC Help
Powered by ViewVC 1.1.22