/[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.113 - (hide annotations) (download)
Fri Jan 28 01:00:13 2005 UTC (19 years, 4 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57d_post, checkpoint57e_post, eckpoint57e_pre, checkpoint57f_pre
Changes since 1.112: +36 -6 lines
move state variable diagnostics to the beginning of the time step.

1 jmc 1.113 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.112 2005/01/24 17:00:17 heimbach 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 DEBUG_STATS_RL
112    
113     C !INPUT/OUTPUT PARAMETERS:
114 cnh 1.1 C == Routine arguments ==
115 cnh 1.8 C myTime - Current time in simulation
116     C myIter - Current iteration number in simulation
117 cnh 1.1 C myThid - Thread number for this instance of the routine.
118 cnh 1.8 _RL myTime
119     INTEGER myIter
120 adcroft 1.47 INTEGER myThid
121 cnh 1.1
122 cnh 1.82 C !LOCAL VARIABLES:
123 cnh 1.1 C == Local variables
124 jmc 1.113 C fVer[UV] o fVer: Vertical flux term - note fVer
125     C is "pipelined" in the vertical
126     C so we need an fVer for each
127     C variable.
128 jmc 1.94 C phiHydC :: hydrostatic potential anomaly at cell center
129     C In z coords phiHyd is the hydrostatic potential
130     C (=pressure/rho0) anomaly
131     C In p coords phiHyd is the geopotential height anomaly.
132     C phiHydF :: hydrostatic potential anomaly at middle between 2 centers
133     C dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
134     C phiSurfX, :: gradient of Surface potential (Pressure/rho, ocean)
135 jmc 1.92 C phiSurfY or geopotential (atmos) in X and Y direction
136 jmc 1.110 C guDissip :: dissipation tendency (all explicit terms), u component
137     C gvDissip :: dissipation tendency (all explicit terms), v component
138 cnh 1.30 C iMin, iMax - Ranges and sub-block indices on which calculations
139     C jMin, jMax are applied.
140 cnh 1.1 C bi, bj
141 heimbach 1.53 C k, kup, - Index for layer above and below. kup and kDown
142     C kDown, km1 are switched with layer to be the appropriate
143 cnh 1.38 C index into fVerTerm.
144 cnh 1.30 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
145     _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
146 jmc 1.94 _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
147     _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
148 jmc 1.92 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
149     _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
150 jmc 1.63 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
151     _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
152 jmc 1.110 _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
153     _RL gvDissip(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.113
166     #ifdef ALLOW_DIAGNOSTICS
167     _RL tmpFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
168     LOGICAL DIAGNOSTICS_IS_ON
169     EXTERNAL DIAGNOSTICS_IS_ON
170     #endif /* ALLOW_DIAGNOSTICS */
171    
172 jmc 1.62
173 adcroft 1.11 C--- The algorithm...
174     C
175     C "Correction Step"
176     C =================
177     C Here we update the horizontal velocities with the surface
178     C pressure such that the resulting flow is either consistent
179     C with the free-surface evolution or the rigid-lid:
180     C U[n] = U* + dt x d/dx P
181     C V[n] = V* + dt x d/dy P
182     C
183     C "Calculation of Gs"
184     C ===================
185     C This is where all the accelerations and tendencies (ie.
186 heimbach 1.53 C physics, parameterizations etc...) are calculated
187 adcroft 1.11 C rho = rho ( theta[n], salt[n] )
188 cnh 1.27 C b = b(rho, theta)
189 adcroft 1.11 C K31 = K31 ( rho )
190 jmc 1.61 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
191     C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
192     C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
193     C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
194 adcroft 1.11 C
195 adcroft 1.12 C "Time-stepping" or "Prediction"
196 adcroft 1.11 C ================================
197     C The models variables are stepped forward with the appropriate
198     C time-stepping scheme (currently we use Adams-Bashforth II)
199     C - For momentum, the result is always *only* a "prediction"
200     C in that the flow may be divergent and will be "corrected"
201     C later with a surface pressure gradient.
202     C - Normally for tracers the result is the new field at time
203     C level [n+1} *BUT* in the case of implicit diffusion the result
204     C is also *only* a prediction.
205     C - We denote "predictors" with an asterisk (*).
206     C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
207     C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
208     C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
209     C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
210 adcroft 1.12 C With implicit diffusion:
211 adcroft 1.11 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
212     C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
213 adcroft 1.12 C (1 + dt * K * d_zz) theta[n] = theta*
214     C (1 + dt * K * d_zz) salt[n] = salt*
215 adcroft 1.11 C---
216 cnh 1.82 CEOP
217 adcroft 1.11
218 heimbach 1.88 C-- Call to routine for calculation of
219     C Eliassen-Palm-flux-forced U-tendency,
220     C if desired:
221     #ifdef INCLUDE_EP_FORCING_CODE
222     CALL CALC_EP_FORCING(myThid)
223     #endif
224    
225 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
226     C-- HPF directive to help TAMC
227     CHPF$ INDEPENDENT
228     #endif /* ALLOW_AUTODIFF_TAMC */
229    
230 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
231 heimbach 1.76
232     #ifdef ALLOW_AUTODIFF_TAMC
233     C-- HPF directive to help TAMC
234     CHPF$ INDEPENDENT, NEW (fVerU,fVerV
235 jmc 1.94 CHPF$& ,phiHydF
236 heimbach 1.76 CHPF$& ,KappaRU,KappaRV
237     CHPF$& )
238     #endif /* ALLOW_AUTODIFF_TAMC */
239    
240 cnh 1.1 DO bi=myBxLo(myThid),myBxHi(myThid)
241 heimbach 1.76
242     #ifdef ALLOW_AUTODIFF_TAMC
243     act1 = bi - myBxLo(myThid)
244     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
245     act2 = bj - myByLo(myThid)
246     max2 = myByHi(myThid) - myByLo(myThid) + 1
247     act3 = myThid - 1
248     max3 = nTx*nTy
249     act4 = ikey_dynamics - 1
250 heimbach 1.91 idynkey = (act1 + 1) + act2*max1
251 heimbach 1.76 & + act3*max1*max2
252     & + act4*max1*max2*max3
253     #endif /* ALLOW_AUTODIFF_TAMC */
254    
255 heimbach 1.97 C-- Set up work arrays with valid (i.e. not NaN) values
256     C These inital values do not alter the numerical results. They
257     C just ensure that all memory references are to valid floating
258     C point numbers. This prevents spurious hardware signals due to
259     C uninitialised but inert locations.
260    
261 jmc 1.94 DO k=1,Nr
262     DO j=1-OLy,sNy+OLy
263     DO i=1-OLx,sNx+OLx
264 heimbach 1.87 KappaRU(i,j,k) = 0. _d 0
265     KappaRV(i,j,k) = 0. _d 0
266 heimbach 1.97 #ifdef ALLOW_AUTODIFF_TAMC
267     cph(
268     c-- need some re-initialisation here to break dependencies
269     cph)
270     gu(i,j,k,bi,bj) = 0. _d 0
271     gv(i,j,k,bi,bj) = 0. _d 0
272     #endif
273 heimbach 1.87 ENDDO
274 jmc 1.94 ENDDO
275     ENDDO
276     DO j=1-OLy,sNy+OLy
277     DO i=1-OLx,sNx+OLx
278 heimbach 1.76 fVerU (i,j,1) = 0. _d 0
279     fVerU (i,j,2) = 0. _d 0
280     fVerV (i,j,1) = 0. _d 0
281     fVerV (i,j,2) = 0. _d 0
282 jmc 1.94 phiHydF (i,j) = 0. _d 0
283     phiHydC (i,j) = 0. _d 0
284 jmc 1.92 dPhiHydX(i,j) = 0. _d 0
285     dPhiHydY(i,j) = 0. _d 0
286 heimbach 1.97 phiSurfX(i,j) = 0. _d 0
287     phiSurfY(i,j) = 0. _d 0
288 jmc 1.110 guDissip(i,j) = 0. _d 0
289     gvDissip(i,j) = 0. _d 0
290 heimbach 1.76 ENDDO
291     ENDDO
292 heimbach 1.49
293 jmc 1.63 C-- Start computation of dynamics
294 jmc 1.93 iMin = 0
295     iMax = sNx+1
296     jMin = 0
297     jMax = sNy+1
298 jmc 1.63
299 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
300 heimbach 1.91 CADJ STORE wvel (:,:,:,bi,bj) =
301     CADJ & comlev1_bibj, key = idynkey, byte = isbyte
302 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
303    
304 jmc 1.65 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
305 jmc 1.63 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
306     IF (implicSurfPress.NE.1.) THEN
307 jmc 1.65 CALL CALC_GRAD_PHI_SURF(
308     I bi,bj,iMin,iMax,jMin,jMax,
309     I etaN,
310     O phiSurfX,phiSurfY,
311     I myThid )
312 jmc 1.63 ENDIF
313 heimbach 1.83
314     #ifdef ALLOW_AUTODIFF_TAMC
315 heimbach 1.91 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
316     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
317 heimbach 1.83 #ifdef ALLOW_KPP
318     CADJ STORE KPPviscAz (:,:,:,bi,bj)
319 heimbach 1.91 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
320 heimbach 1.83 #endif /* ALLOW_KPP */
321     #endif /* ALLOW_AUTODIFF_TAMC */
322 adcroft 1.58
323 heimbach 1.77 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
324     C-- Calculate the total vertical diffusivity
325     DO k=1,Nr
326     CALL CALC_VISCOSITY(
327     I bi,bj,iMin,iMax,jMin,jMax,k,
328     O KappaRU,KappaRV,
329     I myThid)
330     ENDDO
331     #endif
332    
333 heimbach 1.101 #ifdef ALLOW_AUTODIFF_TAMC
334     CADJ STORE KappaRU(:,:,:)
335     CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
336     CADJ STORE KappaRV(:,:,:)
337     CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
338     #endif /* ALLOW_AUTODIFF_TAMC */
339    
340 adcroft 1.58 C-- Start of dynamics loop
341     DO k=1,Nr
342    
343     C-- km1 Points to level above k (=k-1)
344     C-- kup Cycles through 1,2 to point to layer above
345     C-- kDown Cycles through 2,1 to point to current layer
346    
347     km1 = MAX(1,k-1)
348 heimbach 1.77 kp1 = MIN(k+1,Nr)
349 adcroft 1.58 kup = 1+MOD(k+1,2)
350     kDown= 1+MOD(k,2)
351    
352 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
353 heimbach 1.91 kkey = (idynkey-1)*Nr + k
354 heimbach 1.99 c
355 heimbach 1.95 CADJ STORE totphihyd (:,:,k,bi,bj)
356 heimbach 1.99 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
357     CADJ STORE theta (:,:,k,bi,bj)
358     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
359     CADJ STORE salt (:,:,k,bi,bj)
360 heimbach 1.95 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
361 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
362    
363 adcroft 1.58 C-- Integrate hydrostatic balance for phiHyd with BC of
364     C phiHyd(z=0)=0
365 jmc 1.107 CALL CALC_PHI_HYD(
366 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,
367     I theta, salt,
368 jmc 1.94 U phiHydF,
369     O phiHydC, dPhiHydX, dPhiHydY,
370 jmc 1.92 I myTime, myIter, myThid )
371 mlosch 1.89
372 adcroft 1.58 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
373 jmc 1.96 C and step forward storing the result in gU, gV, etc...
374 adcroft 1.58 IF ( momStepping ) THEN
375 edhill 1.105 #ifdef ALLOW_MOM_FLUXFORM
376 adcroft 1.79 IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
377 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
378 jmc 1.94 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
379 adcroft 1.58 U fVerU, fVerV,
380 adcroft 1.80 I myTime, myIter, myThid)
381 adcroft 1.79 #endif
382 edhill 1.105 #ifdef ALLOW_MOM_VECINV
383 adcroft 1.79 IF (vectorInvariantMomentum) CALL MOM_VECINV(
384     I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
385 jmc 1.92 I dPhiHydX,dPhiHydY,KappaRU,KappaRV,
386 adcroft 1.79 U fVerU, fVerV,
387 jmc 1.110 O guDissip, gvDissip,
388 adcroft 1.80 I myTime, myIter, myThid)
389 adcroft 1.79 #endif
390 adcroft 1.58 CALL TIMESTEP(
391 jmc 1.63 I bi,bj,iMin,iMax,jMin,jMax,k,
392 jmc 1.94 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
393 jmc 1.110 I guDissip, gvDissip,
394 jmc 1.96 I myTime, myIter, myThid)
395 adcroft 1.58
396     #ifdef ALLOW_OBCS
397     C-- Apply open boundary conditions
398 jmc 1.96 IF (useOBCS) THEN
399     CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
400     ENDIF
401 adcroft 1.58 #endif /* ALLOW_OBCS */
402    
403     ENDIF
404    
405    
406     C-- end of dynamics k loop (1:Nr)
407     ENDDO
408    
409 jmc 1.106 C-- Implicit Vertical advection & viscosity
410     #ifdef INCLUDE_IMPLVERTADV_CODE
411     IF ( momImplVertAdv ) THEN
412     CALL MOM_U_IMPLICIT_R( kappaRU,
413     I bi, bj, myTime, myIter, myThid )
414     CALL MOM_V_IMPLICIT_R( kappaRV,
415     I bi, bj, myTime, myIter, myThid )
416     ELSEIF ( implicitViscosity ) THEN
417     #else /* INCLUDE_IMPLVERTADV_CODE */
418     IF ( implicitViscosity ) THEN
419     #endif /* INCLUDE_IMPLVERTADV_CODE */
420 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
421 heimbach 1.101 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
422 jmc 1.96 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
423 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
424 adcroft 1.42 CALL IMPLDIFF(
425     I bi, bj, iMin, iMax, jMin, jMax,
426 jmc 1.111 I 0, KappaRU,recip_HFacW,
427 jmc 1.96 U gU,
428 adcroft 1.42 I myThid )
429 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
430 heimbach 1.101 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
431 heimbach 1.97 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
432 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
433 adcroft 1.42 CALL IMPLDIFF(
434     I bi, bj, iMin, iMax, jMin, jMax,
435 jmc 1.111 I 0, KappaRV,recip_HFacS,
436 jmc 1.96 U gV,
437 adcroft 1.42 I myThid )
438 jmc 1.106 ENDIF
439 heimbach 1.49
440 adcroft 1.58 #ifdef ALLOW_OBCS
441     C-- Apply open boundary conditions
442 jmc 1.106 IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
443 adcroft 1.58 DO K=1,Nr
444 jmc 1.96 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
445 adcroft 1.58 ENDDO
446 jmc 1.106 ENDIF
447 adcroft 1.58 #endif /* ALLOW_OBCS */
448 heimbach 1.49
449 edhill 1.102 #ifdef ALLOW_CD_CODE
450 jmc 1.106 IF (implicitViscosity.AND.useCDscheme) THEN
451 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
452 heimbach 1.91 CADJ STORE vVelD(:,:,:,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 jmc 1.111 I 0, KappaRU,recip_HFacW,
457 adcroft 1.42 U vVelD,
458     I myThid )
459 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
460 heimbach 1.91 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
461 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
462 adcroft 1.42 CALL IMPLDIFF(
463     I bi, bj, iMin, iMax, jMin, jMax,
464 jmc 1.111 I 0, KappaRV,recip_HFacS,
465 adcroft 1.42 U uVelD,
466     I myThid )
467 jmc 1.106 ENDIF
468 edhill 1.102 #endif /* ALLOW_CD_CODE */
469 jmc 1.106 C-- End implicit Vertical advection & viscosity
470 cnh 1.1
471     ENDDO
472     ENDDO
473 mlosch 1.90
474 heimbach 1.109 #ifdef ALLOW_OBCS
475     IF (useOBCS) THEN
476     CALL OBCS_PRESCRIBE_EXCHANGES(myThid)
477     ENDIF
478     #endif
479    
480 jmc 1.113 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
481    
482 mlosch 1.90 Cml(
483     C In order to compare the variance of phiHydLow of a p/z-coordinate
484     C run with etaH of a z/p-coordinate run the drift of phiHydLow
485     C has to be removed by something like the following subroutine:
486     C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
487     C & 'phiHydLow', myThid )
488     Cml)
489 adcroft 1.69
490 jmc 1.113 #ifdef ALLOW_DIAGNOSTICS
491     IF ( usediagnostics ) THEN
492    
493     CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
494     CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0,1,0,1,1,myThid)
495    
496     IF ( DIAGNOSTICS_IS_ON('PHIBOTSQ',myThid) ) THEN
497     DO bj = myByLo(myThid), myByHi(myThid)
498     DO bi = myBxLo(myThid), myBxHi(myThid)
499     DO j = 1,sNy
500     DO i = 1,sNx
501     tmpFld(i,j) = phiHydLow(i,j,bi,bj)*phiHydLow(i,j,bi,bj)
502     ENDDO
503     ENDDO
504     CALL DIAGNOSTICS_FILL(tmpFld,'PHIBOTSQ',0,1,2,bi,bj,myThid)
505     ENDDO
506     ENDDO
507     ENDIF
508    
509     ENDIF
510     #endif /* ALLOW_DIAGNOSTICS */
511    
512 edhill 1.104 #ifdef ALLOW_DEBUG
513 heimbach 1.98 If ( debugLevel .GE. debLevB ) THEN
514 adcroft 1.69 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
515 adcroft 1.73 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
516 adcroft 1.69 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
517     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
518     CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
519     CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
520     CALL DEBUG_STATS_RL(Nr,Gu,'Gu (DYNAMICS)',myThid)
521     CALL DEBUG_STATS_RL(Nr,Gv,'Gv (DYNAMICS)',myThid)
522     CALL DEBUG_STATS_RL(Nr,Gt,'Gt (DYNAMICS)',myThid)
523     CALL DEBUG_STATS_RL(Nr,Gs,'Gs (DYNAMICS)',myThid)
524     CALL DEBUG_STATS_RL(Nr,GuNm1,'GuNm1 (DYNAMICS)',myThid)
525     CALL DEBUG_STATS_RL(Nr,GvNm1,'GvNm1 (DYNAMICS)',myThid)
526     CALL DEBUG_STATS_RL(Nr,GtNm1,'GtNm1 (DYNAMICS)',myThid)
527     CALL DEBUG_STATS_RL(Nr,GsNm1,'GsNm1 (DYNAMICS)',myThid)
528 adcroft 1.70 ENDIF
529 adcroft 1.69 #endif
530 cnh 1.1
531     RETURN
532     END

  ViewVC Help
Powered by ViewVC 1.1.22