/[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.117 - (hide annotations) (download)
Mon May 23 20:49:37 2005 UTC (19 years ago) by molod
Branch: MAIN
CVS Tags: checkpoint57i_post
Changes since 1.116: +3 -3 lines
Bug fix with new diag

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

  ViewVC Help
Powered by ViewVC 1.1.22