/[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.101 - (hide annotations) (download)
Fri Oct 10 22:56:08 2003 UTC (20 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint51k_post, checkpoint51o_pre, checkpoint51l_post, checkpoint51n_post, checkpoint51j_post, checkpoint51n_pre, checkpoint51l_pre, checkpoint51o_post, checkpoint51m_post
Branch point for: tg2-branch, checkpoint51n_branch
Changes since 1.100: +10 -1 lines
adjusted some flow directives

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

  ViewVC Help
Powered by ViewVC 1.1.22