/[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.112 - (hide annotations) (download)
Mon Jan 24 17:00:17 2005 UTC (19 years, 4 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57c_post, checkpoint57c_pre
Changes since 1.111: +1 -4 lines
remove re-initialisation of totPhiHyd
o now obsolete for z-coord.
o and wrong for p-coord. with pickup

1 heimbach 1.112 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.111 2004/12/16 23:20:06 jmc 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.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 jmc 1.110 C guDissip :: dissipation tendency (all explicit terms), u component
138     C gvDissip :: dissipation tendency (all explicit terms), v component
139 cnh 1.30 C iMin, iMax - Ranges and sub-block indices on which calculations
140     C jMin, jMax are applied.
141 cnh 1.1 C bi, bj
142 heimbach 1.53 C k, kup, - Index for layer above and below. kup and kDown
143     C kDown, km1 are switched with layer to be the appropriate
144 cnh 1.38 C index into fVerTerm.
145 cnh 1.30 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
146     _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
147 jmc 1.94 _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
148     _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
149 jmc 1.92 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
150     _RL dPhiHydY(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 jmc 1.110 _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
154     _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
155 adcroft 1.42 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
156     _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
157 adcroft 1.12
158 cnh 1.1 INTEGER iMin, iMax
159     INTEGER jMin, jMax
160     INTEGER bi, bj
161     INTEGER i, j
162 heimbach 1.77 INTEGER k, km1, kp1, kup, kDown
163 cnh 1.1
164 jmc 1.92 LOGICAL DIFFERENT_MULTIPLE
165     EXTERNAL DIFFERENT_MULTIPLE
166 jmc 1.62
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.88 C-- Call to routine for calculation of
213     C Eliassen-Palm-flux-forced U-tendency,
214     C if desired:
215     #ifdef INCLUDE_EP_FORCING_CODE
216     CALL CALC_EP_FORCING(myThid)
217     #endif
218    
219 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
220     C-- HPF directive to help TAMC
221     CHPF$ INDEPENDENT
222     #endif /* ALLOW_AUTODIFF_TAMC */
223    
224 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
225 heimbach 1.76
226     #ifdef ALLOW_AUTODIFF_TAMC
227     C-- HPF directive to help TAMC
228     CHPF$ INDEPENDENT, NEW (fVerU,fVerV
229 jmc 1.94 CHPF$& ,phiHydF
230 heimbach 1.76 CHPF$& ,KappaRU,KappaRV
231     CHPF$& )
232     #endif /* ALLOW_AUTODIFF_TAMC */
233    
234 cnh 1.1 DO bi=myBxLo(myThid),myBxHi(myThid)
235 heimbach 1.76
236     #ifdef ALLOW_AUTODIFF_TAMC
237     act1 = bi - myBxLo(myThid)
238     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
239     act2 = bj - myByLo(myThid)
240     max2 = myByHi(myThid) - myByLo(myThid) + 1
241     act3 = myThid - 1
242     max3 = nTx*nTy
243     act4 = ikey_dynamics - 1
244 heimbach 1.91 idynkey = (act1 + 1) + act2*max1
245 heimbach 1.76 & + act3*max1*max2
246     & + act4*max1*max2*max3
247     #endif /* ALLOW_AUTODIFF_TAMC */
248    
249 heimbach 1.97 C-- Set up work arrays with valid (i.e. not NaN) values
250     C These inital values do not alter the numerical results. They
251     C just ensure that all memory references are to valid floating
252     C point numbers. This prevents spurious hardware signals due to
253     C uninitialised but inert locations.
254    
255 jmc 1.94 DO k=1,Nr
256     DO j=1-OLy,sNy+OLy
257     DO i=1-OLx,sNx+OLx
258 heimbach 1.87 KappaRU(i,j,k) = 0. _d 0
259     KappaRV(i,j,k) = 0. _d 0
260 heimbach 1.97 #ifdef ALLOW_AUTODIFF_TAMC
261     cph(
262     c-- need some re-initialisation here to break dependencies
263     cph)
264     gu(i,j,k,bi,bj) = 0. _d 0
265     gv(i,j,k,bi,bj) = 0. _d 0
266     #endif
267 heimbach 1.87 ENDDO
268 jmc 1.94 ENDDO
269     ENDDO
270     DO j=1-OLy,sNy+OLy
271     DO i=1-OLx,sNx+OLx
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 jmc 1.94 phiHydF (i,j) = 0. _d 0
277     phiHydC (i,j) = 0. _d 0
278 jmc 1.92 dPhiHydX(i,j) = 0. _d 0
279     dPhiHydY(i,j) = 0. _d 0
280 heimbach 1.97 phiSurfX(i,j) = 0. _d 0
281     phiSurfY(i,j) = 0. _d 0
282 jmc 1.110 guDissip(i,j) = 0. _d 0
283     gvDissip(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 theta (:,:,k,bi,bj)
352     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
353     CADJ STORE salt (:,:,k,bi,bj)
354 heimbach 1.95 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
355 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
356    
357 adcroft 1.58 C-- Integrate hydrostatic balance for phiHyd with BC of
358     C phiHyd(z=0)=0
359 jmc 1.107 CALL CALC_PHI_HYD(
360 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,
361     I theta, salt,
362 jmc 1.94 U phiHydF,
363     O phiHydC, dPhiHydX, dPhiHydY,
364 jmc 1.92 I myTime, myIter, myThid )
365 mlosch 1.89
366 adcroft 1.58 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
367 jmc 1.96 C and step forward storing the result in gU, gV, etc...
368 adcroft 1.58 IF ( momStepping ) THEN
369 edhill 1.105 #ifdef ALLOW_MOM_FLUXFORM
370 adcroft 1.79 IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
371 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
372 jmc 1.94 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
373 adcroft 1.58 U fVerU, fVerV,
374 adcroft 1.80 I myTime, myIter, myThid)
375 adcroft 1.79 #endif
376 edhill 1.105 #ifdef ALLOW_MOM_VECINV
377 adcroft 1.79 IF (vectorInvariantMomentum) CALL MOM_VECINV(
378     I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
379 jmc 1.92 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
380 adcroft 1.79 U fVerU, fVerV,
381 jmc 1.110 O guDissip, gvDissip,
382 adcroft 1.80 I myTime, myIter, myThid)
383 adcroft 1.79 #endif
384 adcroft 1.58 CALL TIMESTEP(
385 jmc 1.63 I bi,bj,iMin,iMax,jMin,jMax,k,
386 jmc 1.94 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
387 jmc 1.110 I guDissip, gvDissip,
388 jmc 1.96 I myTime, myIter, myThid)
389 adcroft 1.58
390     #ifdef ALLOW_OBCS
391     C-- Apply open boundary conditions
392 jmc 1.96 IF (useOBCS) THEN
393     CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
394     ENDIF
395 adcroft 1.58 #endif /* ALLOW_OBCS */
396    
397     ENDIF
398    
399    
400     C-- end of dynamics k loop (1:Nr)
401     ENDDO
402    
403 jmc 1.106 C-- Implicit Vertical advection & viscosity
404     #ifdef INCLUDE_IMPLVERTADV_CODE
405     IF ( momImplVertAdv ) THEN
406     CALL MOM_U_IMPLICIT_R( kappaRU,
407     I bi, bj, myTime, myIter, myThid )
408     CALL MOM_V_IMPLICIT_R( kappaRV,
409     I bi, bj, myTime, myIter, myThid )
410     ELSEIF ( implicitViscosity ) THEN
411     #else /* INCLUDE_IMPLVERTADV_CODE */
412     IF ( implicitViscosity ) THEN
413     #endif /* INCLUDE_IMPLVERTADV_CODE */
414 adcroft 1.58 #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 jmc 1.111 I 0, 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 jmc 1.111 I 0, KappaRV,recip_HFacS,
430 jmc 1.96 U gV,
431 adcroft 1.42 I myThid )
432 jmc 1.106 ENDIF
433 heimbach 1.49
434 adcroft 1.58 #ifdef ALLOW_OBCS
435     C-- Apply open boundary conditions
436 jmc 1.106 IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
437 adcroft 1.58 DO K=1,Nr
438 jmc 1.96 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
439 adcroft 1.58 ENDDO
440 jmc 1.106 ENDIF
441 adcroft 1.58 #endif /* ALLOW_OBCS */
442 heimbach 1.49
443 edhill 1.102 #ifdef ALLOW_CD_CODE
444 jmc 1.106 IF (implicitViscosity.AND.useCDscheme) THEN
445 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
446 heimbach 1.91 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
447 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
448 adcroft 1.42 CALL IMPLDIFF(
449     I bi, bj, iMin, iMax, jMin, jMax,
450 jmc 1.111 I 0, KappaRU,recip_HFacW,
451 adcroft 1.42 U vVelD,
452     I myThid )
453 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
454 heimbach 1.91 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
455 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
456 adcroft 1.42 CALL IMPLDIFF(
457     I bi, bj, iMin, iMax, jMin, jMax,
458 jmc 1.111 I 0, KappaRV,recip_HFacS,
459 adcroft 1.42 U uVelD,
460     I myThid )
461 jmc 1.106 ENDIF
462 edhill 1.102 #endif /* ALLOW_CD_CODE */
463 jmc 1.106 C-- End implicit Vertical advection & viscosity
464 cnh 1.1
465     ENDDO
466     ENDDO
467 mlosch 1.90
468 heimbach 1.109 #ifdef ALLOW_OBCS
469     IF (useOBCS) THEN
470     CALL OBCS_PRESCRIBE_EXCHANGES(myThid)
471     ENDIF
472     #endif
473    
474 mlosch 1.90 Cml(
475     C In order to compare the variance of phiHydLow of a p/z-coordinate
476     C run with etaH of a z/p-coordinate run the drift of phiHydLow
477     C has to be removed by something like the following subroutine:
478     C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
479     C & 'phiHydLow', myThid )
480     Cml)
481 adcroft 1.69
482 edhill 1.104 #ifdef ALLOW_DEBUG
483 heimbach 1.98 If ( debugLevel .GE. debLevB ) THEN
484 adcroft 1.69 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
485 adcroft 1.73 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
486 adcroft 1.69 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
487     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
488     CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
489     CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
490     CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
491     CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
492     CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
493     CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
494     CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
495     CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
496     CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
497     CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
498 adcroft 1.70 ENDIF
499 adcroft 1.69 #endif
500 cnh 1.1
501     RETURN
502     END

  ViewVC Help
Powered by ViewVC 1.1.22