/[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.123 - (hide annotations) (download)
Tue Aug 16 22:52:06 2005 UTC (18 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57r_post, checkpoint57q_post
Changes since 1.122: +13 -3 lines
add debug messages.

1 jmc 1.123 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.122 2005/07/30 23:39:48 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 jmc 1.122 C | W[n] = W* + dt x d/dz P (NH mode)
34 cnh 1.82 C |
35     C | "Calculation of Gs"
36     C | ===================
37     C | This is where all the accelerations and tendencies (ie.
38     C | physics, parameterizations etc...) are calculated
39     C | rho = rho ( theta[n], salt[n] )
40     C | b = b(rho, theta)
41     C | K31 = K31 ( rho )
42     C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
43     C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
44     C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
45     C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
46     C |
47     C | "Time-stepping" or "Prediction"
48     C | ================================
49     C | The models variables are stepped forward with the appropriate
50     C | time-stepping scheme (currently we use Adams-Bashforth II)
51     C | - For momentum, the result is always *only* a "prediction"
52     C | in that the flow may be divergent and will be "corrected"
53     C | later with a surface pressure gradient.
54     C | - Normally for tracers the result is the new field at time
55     C | level [n+1} *BUT* in the case of implicit diffusion the result
56     C | is also *only* a prediction.
57     C | - We denote "predictors" with an asterisk (*).
58     C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
59     C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
60     C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
61     C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
62     C | With implicit diffusion:
63     C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
64     C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
65     C | (1 + dt * K * d_zz) theta[n] = theta*
66     C | (1 + dt * K * d_zz) salt[n] = salt*
67     C |
68     C *==========================================================*
69     C \ev
70     C !USES:
71 adcroft 1.40 IMPLICIT NONE
72 cnh 1.1 C == Global variables ===
73     #include "SIZE.h"
74     #include "EEPARAMS.h"
75 adcroft 1.6 #include "PARAMS.h"
76 adcroft 1.3 #include "DYNVARS.h"
77 edhill 1.103 #ifdef ALLOW_CD_CODE
78     #include "CD_CODE_VARS.h"
79     #endif
80 adcroft 1.42 #include "GRID.h"
81 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
82 heimbach 1.53 # include "tamc.h"
83     # include "tamc_keys.h"
84 heimbach 1.67 # include "FFIELDS.h"
85 heimbach 1.91 # include "EOS.h"
86 heimbach 1.67 # ifdef ALLOW_KPP
87     # include "KPP.h"
88     # endif
89 heimbach 1.53 #endif /* ALLOW_AUTODIFF_TAMC */
90 jmc 1.62
91 cnh 1.82 C !CALLING SEQUENCE:
92     C DYNAMICS()
93     C |
94 jmc 1.122 C |-- CALC_EP_FORCING
95     C |
96 cnh 1.82 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 jmc 1.122 C |-- MOM_U_IMPLICIT_R
111     C |-- MOM_V_IMPLICIT_R
112     C |
113 cnh 1.82 C |-- IMPLDIFF
114     C |
115     C |-- OBCS_APPLY_UV
116     C |
117 jmc 1.122 C |-- CALC_GW
118     C |
119     C |-- DIAGNOSTICS_FILL
120     C |-- DEBUG_STATS_RL
121 cnh 1.82
122     C !INPUT/OUTPUT PARAMETERS:
123 cnh 1.1 C == Routine arguments ==
124 cnh 1.8 C myTime - Current time in simulation
125     C myIter - Current iteration number in simulation
126 cnh 1.1 C myThid - Thread number for this instance of the routine.
127 cnh 1.8 _RL myTime
128     INTEGER myIter
129 adcroft 1.47 INTEGER myThid
130 cnh 1.1
131 cnh 1.82 C !LOCAL VARIABLES:
132 cnh 1.1 C == Local variables
133 jmc 1.113 C fVer[UV] o fVer: Vertical flux term - note fVer
134     C is "pipelined" in the vertical
135     C so we need an fVer for each
136     C variable.
137 jmc 1.94 C phiHydC :: hydrostatic potential anomaly at cell center
138     C In z coords phiHyd is the hydrostatic potential
139     C (=pressure/rho0) anomaly
140     C In p coords phiHyd is the geopotential height anomaly.
141     C phiHydF :: hydrostatic potential anomaly at middle between 2 centers
142     C dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
143     C phiSurfX, :: gradient of Surface potential (Pressure/rho, ocean)
144 jmc 1.92 C phiSurfY or geopotential (atmos) in X and Y direction
145 jmc 1.110 C guDissip :: dissipation tendency (all explicit terms), u component
146     C gvDissip :: dissipation tendency (all explicit terms), v component
147 cnh 1.30 C iMin, iMax - Ranges and sub-block indices on which calculations
148     C jMin, jMax are applied.
149 cnh 1.1 C bi, bj
150 heimbach 1.53 C k, kup, - Index for layer above and below. kup and kDown
151     C kDown, km1 are switched with layer to be the appropriate
152 cnh 1.38 C index into fVerTerm.
153 cnh 1.30 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
154     _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
155 jmc 1.94 _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
156     _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
157 jmc 1.92 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
158     _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
159 jmc 1.63 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
160     _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
161 jmc 1.110 _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
162     _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
163 adcroft 1.42 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
164     _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
165 adcroft 1.12
166 cnh 1.1 INTEGER iMin, iMax
167     INTEGER jMin, jMax
168     INTEGER bi, bj
169     INTEGER i, j
170 heimbach 1.77 INTEGER k, km1, kp1, kup, kDown
171 cnh 1.1
172 jmc 1.113 #ifdef ALLOW_DIAGNOSTICS
173 jmc 1.120 _RL tmpFac
174 jmc 1.113 #endif /* ALLOW_DIAGNOSTICS */
175    
176 jmc 1.62
177 adcroft 1.11 C--- The algorithm...
178     C
179     C "Correction Step"
180     C =================
181     C Here we update the horizontal velocities with the surface
182     C pressure such that the resulting flow is either consistent
183     C with the free-surface evolution or the rigid-lid:
184     C U[n] = U* + dt x d/dx P
185     C V[n] = V* + dt x d/dy P
186     C
187     C "Calculation of Gs"
188     C ===================
189     C This is where all the accelerations and tendencies (ie.
190 heimbach 1.53 C physics, parameterizations etc...) are calculated
191 adcroft 1.11 C rho = rho ( theta[n], salt[n] )
192 cnh 1.27 C b = b(rho, theta)
193 adcroft 1.11 C K31 = K31 ( rho )
194 jmc 1.61 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
195     C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
196     C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
197     C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
198 adcroft 1.11 C
199 adcroft 1.12 C "Time-stepping" or "Prediction"
200 adcroft 1.11 C ================================
201     C The models variables are stepped forward with the appropriate
202     C time-stepping scheme (currently we use Adams-Bashforth II)
203     C - For momentum, the result is always *only* a "prediction"
204     C in that the flow may be divergent and will be "corrected"
205     C later with a surface pressure gradient.
206     C - Normally for tracers the result is the new field at time
207     C level [n+1} *BUT* in the case of implicit diffusion the result
208     C is also *only* a prediction.
209     C - We denote "predictors" with an asterisk (*).
210     C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
211     C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
212     C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
213     C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
214 adcroft 1.12 C With implicit diffusion:
215 adcroft 1.11 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
216     C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
217 adcroft 1.12 C (1 + dt * K * d_zz) theta[n] = theta*
218     C (1 + dt * K * d_zz) salt[n] = salt*
219 adcroft 1.11 C---
220 cnh 1.82 CEOP
221 adcroft 1.11
222 jmc 1.123 #ifdef ALLOW_DEBUG
223     IF ( debugLevel .GE. debLevB )
224     & CALL DEBUG_ENTER( 'DYNAMICS', myThid )
225     #endif
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 jmc 1.94 CHPF$& ,phiHydF
245 heimbach 1.76 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 heimbach 1.91 idynkey = (act1 + 1) + act2*max1
260 heimbach 1.76 & + act3*max1*max2
261     & + act4*max1*max2*max3
262     #endif /* ALLOW_AUTODIFF_TAMC */
263    
264 heimbach 1.97 C-- Set up work arrays with valid (i.e. not NaN) values
265     C These inital values do not alter the numerical results. They
266     C just ensure that all memory references are to valid floating
267     C point numbers. This prevents spurious hardware signals due to
268     C uninitialised but inert locations.
269    
270 jmc 1.94 DO k=1,Nr
271     DO j=1-OLy,sNy+OLy
272     DO i=1-OLx,sNx+OLx
273 heimbach 1.87 KappaRU(i,j,k) = 0. _d 0
274     KappaRV(i,j,k) = 0. _d 0
275 heimbach 1.97 #ifdef ALLOW_AUTODIFF_TAMC
276     cph(
277     c-- need some re-initialisation here to break dependencies
278     cph)
279 jmc 1.122 gU(i,j,k,bi,bj) = 0. _d 0
280     gV(i,j,k,bi,bj) = 0. _d 0
281 heimbach 1.97 #endif
282 heimbach 1.87 ENDDO
283 jmc 1.94 ENDDO
284     ENDDO
285     DO j=1-OLy,sNy+OLy
286     DO i=1-OLx,sNx+OLx
287 heimbach 1.76 fVerU (i,j,1) = 0. _d 0
288     fVerU (i,j,2) = 0. _d 0
289     fVerV (i,j,1) = 0. _d 0
290     fVerV (i,j,2) = 0. _d 0
291 jmc 1.94 phiHydF (i,j) = 0. _d 0
292     phiHydC (i,j) = 0. _d 0
293 jmc 1.92 dPhiHydX(i,j) = 0. _d 0
294     dPhiHydY(i,j) = 0. _d 0
295 heimbach 1.97 phiSurfX(i,j) = 0. _d 0
296     phiSurfY(i,j) = 0. _d 0
297 jmc 1.110 guDissip(i,j) = 0. _d 0
298     gvDissip(i,j) = 0. _d 0
299 heimbach 1.76 ENDDO
300     ENDDO
301 heimbach 1.49
302 jmc 1.63 C-- Start computation of dynamics
303 jmc 1.93 iMin = 0
304     iMax = sNx+1
305     jMin = 0
306     jMax = sNy+1
307 jmc 1.63
308 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
309 heimbach 1.91 CADJ STORE wvel (:,:,:,bi,bj) =
310     CADJ & comlev1_bibj, key = idynkey, byte = isbyte
311 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
312    
313 jmc 1.65 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
314 jmc 1.63 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
315     IF (implicSurfPress.NE.1.) THEN
316 jmc 1.65 CALL CALC_GRAD_PHI_SURF(
317     I bi,bj,iMin,iMax,jMin,jMax,
318     I etaN,
319     O phiSurfX,phiSurfY,
320     I myThid )
321 jmc 1.63 ENDIF
322 heimbach 1.83
323     #ifdef ALLOW_AUTODIFF_TAMC
324 heimbach 1.91 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
325     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
326 heimbach 1.83 #ifdef ALLOW_KPP
327     CADJ STORE KPPviscAz (:,:,:,bi,bj)
328 heimbach 1.91 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
329 heimbach 1.83 #endif /* ALLOW_KPP */
330     #endif /* ALLOW_AUTODIFF_TAMC */
331 adcroft 1.58
332 heimbach 1.77 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
333     C-- Calculate the total vertical diffusivity
334     DO k=1,Nr
335     CALL CALC_VISCOSITY(
336     I bi,bj,iMin,iMax,jMin,jMax,k,
337     O KappaRU,KappaRV,
338     I myThid)
339     ENDDO
340     #endif
341    
342 heimbach 1.101 #ifdef ALLOW_AUTODIFF_TAMC
343     CADJ STORE KappaRU(:,:,:)
344     CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
345     CADJ STORE KappaRV(:,:,:)
346     CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
347     #endif /* ALLOW_AUTODIFF_TAMC */
348    
349 adcroft 1.58 C-- Start of dynamics loop
350     DO k=1,Nr
351    
352     C-- km1 Points to level above k (=k-1)
353     C-- kup Cycles through 1,2 to point to layer above
354     C-- kDown Cycles through 2,1 to point to current layer
355    
356     km1 = MAX(1,k-1)
357 heimbach 1.77 kp1 = MIN(k+1,Nr)
358 adcroft 1.58 kup = 1+MOD(k+1,2)
359     kDown= 1+MOD(k,2)
360    
361 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
362 heimbach 1.91 kkey = (idynkey-1)*Nr + k
363 heimbach 1.99 c
364 heimbach 1.95 CADJ STORE totphihyd (:,:,k,bi,bj)
365 heimbach 1.99 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
366     CADJ STORE theta (:,:,k,bi,bj)
367     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
368     CADJ STORE salt (:,:,k,bi,bj)
369 heimbach 1.95 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
370 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
371    
372 adcroft 1.58 C-- Integrate hydrostatic balance for phiHyd with BC of
373     C phiHyd(z=0)=0
374 jmc 1.107 CALL CALC_PHI_HYD(
375 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,
376     I theta, salt,
377 jmc 1.94 U phiHydF,
378     O phiHydC, dPhiHydX, dPhiHydY,
379 jmc 1.92 I myTime, myIter, myThid )
380 mlosch 1.89
381 adcroft 1.58 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
382 jmc 1.96 C and step forward storing the result in gU, gV, etc...
383 adcroft 1.58 IF ( momStepping ) THEN
384 edhill 1.105 #ifdef ALLOW_MOM_FLUXFORM
385 adcroft 1.79 IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
386 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
387 jmc 1.121 I KappaRU, KappaRV,
388 adcroft 1.58 U fVerU, fVerV,
389 jmc 1.121 O guDissip, gvDissip,
390 adcroft 1.80 I myTime, myIter, myThid)
391 adcroft 1.79 #endif
392 edhill 1.105 #ifdef ALLOW_MOM_VECINV
393 adcroft 1.79 IF (vectorInvariantMomentum) CALL MOM_VECINV(
394     I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
395 jmc 1.121 I KappaRU, KappaRV,
396 adcroft 1.79 U fVerU, fVerV,
397 jmc 1.110 O guDissip, gvDissip,
398 adcroft 1.80 I myTime, myIter, myThid)
399 adcroft 1.79 #endif
400 adcroft 1.58 CALL TIMESTEP(
401 jmc 1.63 I bi,bj,iMin,iMax,jMin,jMax,k,
402 jmc 1.94 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
403 jmc 1.110 I guDissip, gvDissip,
404 jmc 1.96 I myTime, myIter, myThid)
405 adcroft 1.58
406     #ifdef ALLOW_OBCS
407     C-- Apply open boundary conditions
408 jmc 1.96 IF (useOBCS) THEN
409     CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
410     ENDIF
411 adcroft 1.58 #endif /* ALLOW_OBCS */
412    
413     ENDIF
414    
415    
416     C-- end of dynamics k loop (1:Nr)
417     ENDDO
418    
419 jmc 1.106 C-- Implicit Vertical advection & viscosity
420     #ifdef INCLUDE_IMPLVERTADV_CODE
421     IF ( momImplVertAdv ) THEN
422     CALL MOM_U_IMPLICIT_R( kappaRU,
423     I bi, bj, myTime, myIter, myThid )
424     CALL MOM_V_IMPLICIT_R( kappaRV,
425     I bi, bj, myTime, myIter, myThid )
426     ELSEIF ( implicitViscosity ) THEN
427     #else /* INCLUDE_IMPLVERTADV_CODE */
428     IF ( implicitViscosity ) THEN
429     #endif /* INCLUDE_IMPLVERTADV_CODE */
430 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
431 heimbach 1.101 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
432 jmc 1.96 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, 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 jmc 1.111 I 0, KappaRU,recip_HFacW,
437 jmc 1.96 U gU,
438 adcroft 1.42 I myThid )
439 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
440 heimbach 1.101 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
441 heimbach 1.97 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
442 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
443 adcroft 1.42 CALL IMPLDIFF(
444     I bi, bj, iMin, iMax, jMin, jMax,
445 jmc 1.111 I 0, KappaRV,recip_HFacS,
446 jmc 1.96 U gV,
447 adcroft 1.42 I myThid )
448 jmc 1.106 ENDIF
449 heimbach 1.49
450 adcroft 1.58 #ifdef ALLOW_OBCS
451     C-- Apply open boundary conditions
452 jmc 1.106 IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
453 adcroft 1.58 DO K=1,Nr
454 jmc 1.96 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
455 adcroft 1.58 ENDDO
456 jmc 1.106 ENDIF
457 adcroft 1.58 #endif /* ALLOW_OBCS */
458 heimbach 1.49
459 edhill 1.102 #ifdef ALLOW_CD_CODE
460 jmc 1.106 IF (implicitViscosity.AND.useCDscheme) THEN
461 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
462 heimbach 1.91 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
463 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
464 adcroft 1.42 CALL IMPLDIFF(
465     I bi, bj, iMin, iMax, jMin, jMax,
466 jmc 1.111 I 0, KappaRU,recip_HFacW,
467 adcroft 1.42 U vVelD,
468     I myThid )
469 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
470 heimbach 1.91 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
471 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
472 adcroft 1.42 CALL IMPLDIFF(
473     I bi, bj, iMin, iMax, jMin, jMax,
474 jmc 1.111 I 0, KappaRV,recip_HFacS,
475 adcroft 1.42 U uVelD,
476     I myThid )
477 jmc 1.106 ENDIF
478 edhill 1.102 #endif /* ALLOW_CD_CODE */
479 jmc 1.106 C-- End implicit Vertical advection & viscosity
480 cnh 1.1
481     ENDDO
482     ENDDO
483 mlosch 1.90
484 heimbach 1.109 #ifdef ALLOW_OBCS
485     IF (useOBCS) THEN
486     CALL OBCS_PRESCRIBE_EXCHANGES(myThid)
487     ENDIF
488     #endif
489    
490 jmc 1.113 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
491    
492 jmc 1.122 #ifdef ALLOW_NONHYDROSTATIC
493     C-- Step forward W field in N-H algorithm
494     IF ( momStepping .AND. nonHydrostatic ) THEN
495     #ifdef ALLOW_DEBUG
496 jmc 1.123 IF ( debugLevel .GE. debLevB )
497     & CALL DEBUG_CALL('CALC_GW', myThid )
498 jmc 1.122 #endif
499     CALL TIMER_START('CALC_GW [DYNAMICS]',myThid)
500     CALL CALC_GW( myTime, myIter, myThid )
501     CALL TIMER_STOP ('CALC_GW [DYNAMICS]',myThid)
502     ENDIF
503     #endif
504    
505     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
506    
507 mlosch 1.90 Cml(
508     C In order to compare the variance of phiHydLow of a p/z-coordinate
509     C run with etaH of a z/p-coordinate run the drift of phiHydLow
510     C has to be removed by something like the following subroutine:
511     C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
512     C & 'phiHydLow', myThid )
513     Cml)
514 adcroft 1.69
515 jmc 1.113 #ifdef ALLOW_DIAGNOSTICS
516     IF ( usediagnostics ) THEN
517    
518     CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
519 jmc 1.120 CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0, 1,0,1,1,myThid)
520 molod 1.116
521 jmc 1.120 tmpFac = 1. _d 0
522     CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
523     & 'PHIHYDSQ',0,Nr,0,1,1,myThid)
524 molod 1.116
525 jmc 1.120 CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
526     & 'PHIBOTSQ',0, 1,0,1,1,myThid)
527 jmc 1.113
528     ENDIF
529     #endif /* ALLOW_DIAGNOSTICS */
530    
531 edhill 1.104 #ifdef ALLOW_DEBUG
532 heimbach 1.98 If ( debugLevel .GE. debLevB ) THEN
533 adcroft 1.69 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
534 adcroft 1.73 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
535 adcroft 1.69 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
536     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
537     CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
538     CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
539 jmc 1.115 CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
540     CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
541     CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
542     CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
543     #ifndef ALLOW_ADAMSBASHFORTH_3
544     CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
545     CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
546     CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
547     CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
548     #endif
549 adcroft 1.70 ENDIF
550 adcroft 1.69 #endif
551 cnh 1.1
552 jmc 1.123 #ifdef ALLOW_DEBUG
553     IF ( debugLevel .GE. debLevB )
554     & CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
555     #endif
556    
557 cnh 1.1 RETURN
558     END

  ViewVC Help
Powered by ViewVC 1.1.22