/[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.105 - (hide annotations) (download)
Tue Nov 4 19:51:53 2003 UTC (20 years, 6 months ago) by edhill
Branch: MAIN
CVS Tags: checkpoint52d_pre, checkpoint52, checkpoint51t_post, checkpoint51s_post, checkpoint52e_pre, checkpoint52b_pre, checkpoint52b_post, checkpoint52c_post, checkpoint52d_post, checkpoint52a_pre, branch-netcdf, checkpoint52a_post, ecco_c52_e35, checkpoint51u_post
Branch point for: netcdf-sm0
Changes since 1.104: +3 -3 lines
 o cleanup: convert DISABLE_MOM_FLUXFORM & DISABLE_MOM_VECINV to the
   newer ALLOW_${pkg} form
   - the only remaining package-based "special case" within genmake2
     is the one for AIM vs. AIM_V23

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

  ViewVC Help
Powered by ViewVC 1.1.22