/[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.144 - (hide annotations) (download)
Sat Jan 16 22:55:53 2010 UTC (14 years, 3 months ago) by jmc
Branch: MAIN
Changes since 1.143: +19 -19 lines
replace maskH by maskInC

1 jmc 1.144 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.143 2009/10/26 21:48:43 gforget 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 heimbach 1.131 #ifdef ALLOW_OBCS
7     # include "OBCS_OPTIONS.h"
8     #endif
9    
10 jmc 1.125 #undef DYNAMICS_GUGV_EXCH_CHECK
11 cnh 1.1
12 cnh 1.82 CBOP
13     C !ROUTINE: DYNAMICS
14     C !INTERFACE:
15 cnh 1.8 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
16 cnh 1.82 C !DESCRIPTION: \bv
17     C *==========================================================*
18 jmc 1.144 C | SUBROUTINE DYNAMICS
19     C | o Controlling routine for the explicit part of the model
20     C | dynamics.
21 cnh 1.82 C *==========================================================*
22 jmc 1.144 C | This routine evaluates the "dynamics" terms for each
23     C | block of ocean in turn. Because the blocks of ocean have
24     C | overlap regions they are independent of one another.
25     C | If terms involving lateral integrals are needed in this
26     C | routine care will be needed. Similarly finite-difference
27     C | operations with stencils wider than the overlap region
28     C | require special consideration.
29 cnh 1.82 C | The algorithm...
30     C |
31     C | "Correction Step"
32     C | =================
33     C | Here we update the horizontal velocities with the surface
34     C | pressure such that the resulting flow is either consistent
35     C | with the free-surface evolution or the rigid-lid:
36     C | U[n] = U* + dt x d/dx P
37     C | V[n] = V* + dt x d/dy P
38 jmc 1.122 C | W[n] = W* + dt x d/dz P (NH mode)
39 cnh 1.82 C |
40     C | "Calculation of Gs"
41     C | ===================
42     C | This is where all the accelerations and tendencies (ie.
43     C | physics, parameterizations etc...) are calculated
44     C | rho = rho ( theta[n], salt[n] )
45     C | b = b(rho, theta)
46     C | K31 = K31 ( rho )
47     C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
48     C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
49     C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
50     C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
51     C |
52     C | "Time-stepping" or "Prediction"
53     C | ================================
54     C | The models variables are stepped forward with the appropriate
55     C | time-stepping scheme (currently we use Adams-Bashforth II)
56     C | - For momentum, the result is always *only* a "prediction"
57     C | in that the flow may be divergent and will be "corrected"
58     C | later with a surface pressure gradient.
59     C | - Normally for tracers the result is the new field at time
60     C | level [n+1} *BUT* in the case of implicit diffusion the result
61     C | is also *only* a prediction.
62     C | - We denote "predictors" with an asterisk (*).
63     C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
64     C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
65     C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
66     C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
67     C | With implicit diffusion:
68     C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
69     C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
70     C | (1 + dt * K * d_zz) theta[n] = theta*
71     C | (1 + dt * K * d_zz) salt[n] = salt*
72     C |
73     C *==========================================================*
74     C \ev
75     C !USES:
76 adcroft 1.40 IMPLICIT NONE
77 cnh 1.1 C == Global variables ===
78     #include "SIZE.h"
79     #include "EEPARAMS.h"
80 adcroft 1.6 #include "PARAMS.h"
81 adcroft 1.3 #include "DYNVARS.h"
82 edhill 1.103 #ifdef ALLOW_CD_CODE
83     #include "CD_CODE_VARS.h"
84     #endif
85 adcroft 1.42 #include "GRID.h"
86 heimbach 1.49 #ifdef ALLOW_AUTODIFF_TAMC
87 heimbach 1.53 # include "tamc.h"
88     # include "tamc_keys.h"
89 heimbach 1.67 # include "FFIELDS.h"
90 heimbach 1.91 # include "EOS.h"
91 heimbach 1.67 # ifdef ALLOW_KPP
92     # include "KPP.h"
93     # endif
94 heimbach 1.131 # ifdef ALLOW_PTRACERS
95     # include "PTRACERS_SIZE.h"
96 jmc 1.139 # include "PTRACERS_FIELDS.h"
97 heimbach 1.131 # endif
98     # ifdef ALLOW_OBCS
99     # include "OBCS.h"
100     # ifdef ALLOW_PTRACERS
101     # include "OBCS_PTRACERS.h"
102     # endif
103     # endif
104 heimbach 1.133 # ifdef ALLOW_MOM_FLUXFORM
105     # include "MOM_FLUXFORM.h"
106     # endif
107 heimbach 1.53 #endif /* ALLOW_AUTODIFF_TAMC */
108 jmc 1.62
109 cnh 1.82 C !CALLING SEQUENCE:
110     C DYNAMICS()
111     C |
112 jmc 1.122 C |-- CALC_EP_FORCING
113     C |
114 cnh 1.82 C |-- CALC_GRAD_PHI_SURF
115     C |
116     C |-- CALC_VISCOSITY
117     C |
118 jmc 1.136 C |-- CALC_PHI_HYD
119 cnh 1.82 C |
120 jmc 1.136 C |-- MOM_FLUXFORM
121 cnh 1.82 C |
122 jmc 1.136 C |-- MOM_VECINV
123 cnh 1.82 C |
124 jmc 1.136 C |-- TIMESTEP
125 cnh 1.82 C |
126 jmc 1.136 C |-- OBCS_APPLY_UV
127 cnh 1.82 C |
128 jmc 1.136 C |-- MOM_U_IMPLICIT_R
129     C |-- MOM_V_IMPLICIT_R
130 jmc 1.122 C |
131 jmc 1.136 C |-- IMPLDIFF
132 cnh 1.82 C |
133     C |-- OBCS_APPLY_UV
134     C |
135 jmc 1.122 C |-- CALC_GW
136     C |
137     C |-- DIAGNOSTICS_FILL
138     C |-- DEBUG_STATS_RL
139 cnh 1.82
140     C !INPUT/OUTPUT PARAMETERS:
141 cnh 1.1 C == Routine arguments ==
142 jmc 1.140 C myTime :: Current time in simulation
143     C myIter :: Current iteration number in simulation
144     C myThid :: Thread number for this instance of the routine.
145 cnh 1.8 _RL myTime
146     INTEGER myIter
147 adcroft 1.47 INTEGER myThid
148 cnh 1.1
149 cnh 1.82 C !LOCAL VARIABLES:
150 cnh 1.1 C == Local variables
151 jmc 1.113 C fVer[UV] o fVer: Vertical flux term - note fVer
152     C is "pipelined" in the vertical
153     C so we need an fVer for each
154     C variable.
155 jmc 1.94 C phiHydC :: hydrostatic potential anomaly at cell center
156     C In z coords phiHyd is the hydrostatic potential
157     C (=pressure/rho0) anomaly
158     C In p coords phiHyd is the geopotential height anomaly.
159     C phiHydF :: hydrostatic potential anomaly at middle between 2 centers
160     C dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
161     C phiSurfX, :: gradient of Surface potential (Pressure/rho, ocean)
162 jmc 1.92 C phiSurfY or geopotential (atmos) in X and Y direction
163 jmc 1.110 C guDissip :: dissipation tendency (all explicit terms), u component
164     C gvDissip :: dissipation tendency (all explicit terms), v component
165 jmc 1.140 C KappaRU :: vertical viscosity
166     C KappaRV :: vertical viscosity
167 cnh 1.30 C iMin, iMax - Ranges and sub-block indices on which calculations
168     C jMin, jMax are applied.
169 cnh 1.1 C bi, bj
170 heimbach 1.53 C k, kup, - Index for layer above and below. kup and kDown
171 jmc 1.144 C kDown, km1 are switched with layer to be the appropriate
172 cnh 1.38 C index into fVerTerm.
173 cnh 1.30 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
174     _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
175 jmc 1.94 _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
176     _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
177 jmc 1.92 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
178     _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
179 jmc 1.63 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
180     _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
181 jmc 1.110 _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
182     _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
183 adcroft 1.42 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
184     _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
185 adcroft 1.12
186 cnh 1.1 INTEGER iMin, iMax
187     INTEGER jMin, jMax
188     INTEGER bi, bj
189     INTEGER i, j
190 heimbach 1.77 INTEGER k, km1, kp1, kup, kDown
191 cnh 1.1
192 jmc 1.113 #ifdef ALLOW_DIAGNOSTICS
193 jmc 1.120 _RL tmpFac
194 jmc 1.113 #endif /* ALLOW_DIAGNOSTICS */
195    
196 jmc 1.136
197 adcroft 1.11 C--- The algorithm...
198     C
199     C "Correction Step"
200     C =================
201     C Here we update the horizontal velocities with the surface
202     C pressure such that the resulting flow is either consistent
203     C with the free-surface evolution or the rigid-lid:
204     C U[n] = U* + dt x d/dx P
205     C V[n] = V* + dt x d/dy P
206     C
207     C "Calculation of Gs"
208     C ===================
209     C This is where all the accelerations and tendencies (ie.
210 heimbach 1.53 C physics, parameterizations etc...) are calculated
211 adcroft 1.11 C rho = rho ( theta[n], salt[n] )
212 cnh 1.27 C b = b(rho, theta)
213 adcroft 1.11 C K31 = K31 ( rho )
214 jmc 1.61 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
215     C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
216     C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
217     C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
218 adcroft 1.11 C
219 adcroft 1.12 C "Time-stepping" or "Prediction"
220 adcroft 1.11 C ================================
221     C The models variables are stepped forward with the appropriate
222     C time-stepping scheme (currently we use Adams-Bashforth II)
223     C - For momentum, the result is always *only* a "prediction"
224     C in that the flow may be divergent and will be "corrected"
225     C later with a surface pressure gradient.
226     C - Normally for tracers the result is the new field at time
227     C level [n+1} *BUT* in the case of implicit diffusion the result
228     C is also *only* a prediction.
229     C - We denote "predictors" with an asterisk (*).
230     C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
231     C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
232     C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
233     C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
234 adcroft 1.12 C With implicit diffusion:
235 adcroft 1.11 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
236     C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
237 adcroft 1.12 C (1 + dt * K * d_zz) theta[n] = theta*
238     C (1 + dt * K * d_zz) salt[n] = salt*
239 adcroft 1.11 C---
240 cnh 1.82 CEOP
241 adcroft 1.11
242 jmc 1.123 #ifdef ALLOW_DEBUG
243     IF ( debugLevel .GE. debLevB )
244     & CALL DEBUG_ENTER( 'DYNAMICS', myThid )
245     #endif
246    
247 heimbach 1.88 C-- Call to routine for calculation of
248     C Eliassen-Palm-flux-forced U-tendency,
249     C if desired:
250     #ifdef INCLUDE_EP_FORCING_CODE
251     CALL CALC_EP_FORCING(myThid)
252     #endif
253    
254 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
255     C-- HPF directive to help TAMC
256     CHPF$ INDEPENDENT
257     #endif /* ALLOW_AUTODIFF_TAMC */
258    
259 cnh 1.1 DO bj=myByLo(myThid),myByHi(myThid)
260 heimbach 1.76
261     #ifdef ALLOW_AUTODIFF_TAMC
262     C-- HPF directive to help TAMC
263     CHPF$ INDEPENDENT, NEW (fVerU,fVerV
264 jmc 1.94 CHPF$& ,phiHydF
265 heimbach 1.76 CHPF$& ,KappaRU,KappaRV
266     CHPF$& )
267     #endif /* ALLOW_AUTODIFF_TAMC */
268    
269 cnh 1.1 DO bi=myBxLo(myThid),myBxHi(myThid)
270 heimbach 1.76
271     #ifdef ALLOW_AUTODIFF_TAMC
272     act1 = bi - myBxLo(myThid)
273     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
274     act2 = bj - myByLo(myThid)
275     max2 = myByHi(myThid) - myByLo(myThid) + 1
276     act3 = myThid - 1
277     max3 = nTx*nTy
278     act4 = ikey_dynamics - 1
279 heimbach 1.91 idynkey = (act1 + 1) + act2*max1
280 heimbach 1.76 & + act3*max1*max2
281     & + act4*max1*max2*max3
282     #endif /* ALLOW_AUTODIFF_TAMC */
283    
284 heimbach 1.97 C-- Set up work arrays with valid (i.e. not NaN) values
285     C These inital values do not alter the numerical results. They
286     C just ensure that all memory references are to valid floating
287     C point numbers. This prevents spurious hardware signals due to
288     C uninitialised but inert locations.
289    
290 jmc 1.140 #ifdef ALLOW_AUTODIFF_TAMC
291 jmc 1.94 DO k=1,Nr
292     DO j=1-OLy,sNy+OLy
293     DO i=1-OLx,sNx+OLx
294 heimbach 1.87 KappaRU(i,j,k) = 0. _d 0
295     KappaRV(i,j,k) = 0. _d 0
296 heimbach 1.97 cph(
297     c-- need some re-initialisation here to break dependencies
298     cph)
299 jmc 1.122 gU(i,j,k,bi,bj) = 0. _d 0
300     gV(i,j,k,bi,bj) = 0. _d 0
301 heimbach 1.87 ENDDO
302 jmc 1.94 ENDDO
303     ENDDO
304 jmc 1.140 #endif /* ALLOW_AUTODIFF_TAMC */
305 jmc 1.94 DO j=1-OLy,sNy+OLy
306     DO i=1-OLx,sNx+OLx
307 heimbach 1.76 fVerU (i,j,1) = 0. _d 0
308     fVerU (i,j,2) = 0. _d 0
309     fVerV (i,j,1) = 0. _d 0
310     fVerV (i,j,2) = 0. _d 0
311 jmc 1.136 phiHydF (i,j) = 0. _d 0
312     phiHydC (i,j) = 0. _d 0
313 jmc 1.92 dPhiHydX(i,j) = 0. _d 0
314 jmc 1.136 dPhiHydY(i,j) = 0. _d 0
315 heimbach 1.97 phiSurfX(i,j) = 0. _d 0
316     phiSurfY(i,j) = 0. _d 0
317 jmc 1.110 guDissip(i,j) = 0. _d 0
318     gvDissip(i,j) = 0. _d 0
319 heimbach 1.138 #ifdef ALLOW_AUTODIFF_TAMC
320 gforget 1.143 phiHydLow(i,j,bi,bj) = 0. _d 0
321 heimbach 1.138 # ifdef NONLIN_FRSURF
322     # ifndef DISABLE_RSTAR_CODE
323     dWtransC(i,j,bi,bj) = 0. _d 0
324     dWtransU(i,j,bi,bj) = 0. _d 0
325     dWtransV(i,j,bi,bj) = 0. _d 0
326     # endif
327     # endif
328     #endif
329 heimbach 1.76 ENDDO
330     ENDDO
331 heimbach 1.49
332 jmc 1.63 C-- Start computation of dynamics
333 jmc 1.93 iMin = 0
334     iMax = sNx+1
335     jMin = 0
336     jMax = sNy+1
337 jmc 1.63
338 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
339 heimbach 1.91 CADJ STORE wvel (:,:,:,bi,bj) =
340 heimbach 1.141 CADJ & comlev1_bibj, key=idynkey, byte=isbyte
341 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
342    
343 jmc 1.65 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
344 jmc 1.63 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
345     IF (implicSurfPress.NE.1.) THEN
346 jmc 1.65 CALL CALC_GRAD_PHI_SURF(
347     I bi,bj,iMin,iMax,jMin,jMax,
348     I etaN,
349     O phiSurfX,phiSurfY,
350 jmc 1.136 I myThid )
351 jmc 1.63 ENDIF
352 heimbach 1.83
353     #ifdef ALLOW_AUTODIFF_TAMC
354 heimbach 1.91 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
355     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
356 heimbach 1.83 #ifdef ALLOW_KPP
357     CADJ STORE KPPviscAz (:,:,:,bi,bj)
358 heimbach 1.91 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
359 heimbach 1.83 #endif /* ALLOW_KPP */
360     #endif /* ALLOW_AUTODIFF_TAMC */
361 adcroft 1.58
362 heimbach 1.77 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
363 jmc 1.140 C-- Calculate the total vertical viscosity
364     CALL CALC_VISCOSITY(
365     I bi,bj, iMin,iMax,jMin,jMax,
366     O KappaRU, KappaRV,
367     I myThid )
368     #else
369 heimbach 1.77 DO k=1,Nr
370 jmc 1.140 DO j=1-OLy,sNy+OLy
371     DO i=1-OLx,sNx+OLx
372     KappaRU(i,j,k) = 0. _d 0
373     KappaRV(i,j,k) = 0. _d 0
374     ENDDO
375     ENDDO
376     ENDDO
377 heimbach 1.77 #endif
378    
379 heimbach 1.101 #ifdef ALLOW_AUTODIFF_TAMC
380     CADJ STORE KappaRU(:,:,:)
381 heimbach 1.132 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
382 heimbach 1.101 CADJ STORE KappaRV(:,:,:)
383 heimbach 1.132 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
384 heimbach 1.101 #endif /* ALLOW_AUTODIFF_TAMC */
385    
386 adcroft 1.58 C-- Start of dynamics loop
387     DO k=1,Nr
388    
389     C-- km1 Points to level above k (=k-1)
390     C-- kup Cycles through 1,2 to point to layer above
391     C-- kDown Cycles through 2,1 to point to current layer
392    
393     km1 = MAX(1,k-1)
394 heimbach 1.77 kp1 = MIN(k+1,Nr)
395 adcroft 1.58 kup = 1+MOD(k+1,2)
396     kDown= 1+MOD(k,2)
397    
398 jmc 1.144 #ifdef ALLOW_AUTODIFF_TAMC
399 heimbach 1.91 kkey = (idynkey-1)*Nr + k
400 heimbach 1.99 c
401 heimbach 1.95 CADJ STORE totphihyd (:,:,k,bi,bj)
402 heimbach 1.99 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
403 gforget 1.143 CADJ STORE phihydlow (:,:,bi,bj)
404     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
405 heimbach 1.99 CADJ STORE theta (:,:,k,bi,bj)
406     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
407     CADJ STORE salt (:,:,k,bi,bj)
408 heimbach 1.95 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
409 heimbach 1.129 CADJ STORE gt(:,:,k,bi,bj)
410     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
411     CADJ STORE gs(:,:,k,bi,bj)
412     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
413 heimbach 1.126 # ifdef NONLIN_FRSURF
414     cph-test
415     CADJ STORE phiHydC (:,:)
416     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
417     CADJ STORE phiHydF (:,:)
418     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
419     CADJ STORE gudissip (:,:)
420     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
421     CADJ STORE gvdissip (:,:)
422     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
423     CADJ STORE fVerU (:,:,:)
424     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
425     CADJ STORE fVerV (:,:,:)
426     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
427     CADJ STORE gu(:,:,k,bi,bj)
428     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
429     CADJ STORE gv(:,:,k,bi,bj)
430     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
431     CADJ STORE gunm1(:,:,k,bi,bj)
432     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
433     CADJ STORE gvnm1(:,:,k,bi,bj)
434     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
435     # ifdef ALLOW_CD_CODE
436     CADJ STORE unm1(:,:,k,bi,bj)
437     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
438     CADJ STORE vnm1(:,:,k,bi,bj)
439     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
440     CADJ STORE uVelD(:,:,k,bi,bj)
441     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
442     CADJ STORE vVelD(:,:,k,bi,bj)
443     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
444     # endif
445     # endif
446 heimbach 1.134 # ifdef ALLOW_DEPTH_CONTROL
447     CADJ STORE fVerU (:,:,:)
448     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
449     CADJ STORE fVerV (:,:,:)
450     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
451     # endif
452 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
453    
454 jmc 1.144 C-- Integrate hydrostatic balance for phiHyd with BC of
455 adcroft 1.58 C phiHyd(z=0)=0
456 jmc 1.128 IF ( implicitIntGravWave ) THEN
457     CALL CALC_PHI_HYD(
458     I bi,bj,iMin,iMax,jMin,jMax,k,
459     I gT, gS,
460     U phiHydF,
461     O phiHydC, dPhiHydX, dPhiHydY,
462     I myTime, myIter, myThid )
463     ELSE
464     CALL CALC_PHI_HYD(
465 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,
466     I theta, salt,
467 jmc 1.94 U phiHydF,
468     O phiHydC, dPhiHydX, dPhiHydY,
469 jmc 1.92 I myTime, myIter, myThid )
470 jmc 1.128 ENDIF
471 mlosch 1.89
472 adcroft 1.58 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
473 jmc 1.96 C and step forward storing the result in gU, gV, etc...
474 adcroft 1.58 IF ( momStepping ) THEN
475 heimbach 1.138 #ifdef ALLOW_AUTODIFF_TAMC
476     # ifdef NONLIN_FRSURF
477     # ifndef DISABLE_RSTAR_CODE
478     CADJ STORE dWtransC(:,:,bi,bj)
479     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
480     CADJ STORE dWtransU(:,:,bi,bj)
481     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
482     CADJ STORE dWtransV(:,:,bi,bj)
483     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
484     # endif
485     # endif
486     #endif
487 heimbach 1.132 IF (.NOT. vectorInvariantMomentum) THEN
488 edhill 1.105 #ifdef ALLOW_MOM_FLUXFORM
489 heimbach 1.132 C
490     CALL MOM_FLUXFORM(
491 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
492 jmc 1.121 I KappaRU, KappaRV,
493 adcroft 1.58 U fVerU, fVerV,
494 jmc 1.121 O guDissip, gvDissip,
495 adcroft 1.80 I myTime, myIter, myThid)
496 adcroft 1.79 #endif
497 heimbach 1.132 ELSE
498 edhill 1.105 #ifdef ALLOW_MOM_VECINV
499 heimbach 1.126 C
500 jmc 1.144 # ifdef ALLOW_AUTODIFF_TAMC
501 heimbach 1.126 # ifdef NONLIN_FRSURF
502     CADJ STORE fVerU(:,:,:)
503     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
504     CADJ STORE fVerV(:,:,:)
505     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
506     # endif
507     # endif /* ALLOW_AUTODIFF_TAMC */
508     C
509     CALL MOM_VECINV(
510 adcroft 1.79 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
511 jmc 1.121 I KappaRU, KappaRV,
512 adcroft 1.79 U fVerU, fVerV,
513 jmc 1.110 O guDissip, gvDissip,
514 adcroft 1.80 I myTime, myIter, myThid)
515 heimbach 1.132 #endif
516 heimbach 1.126 ENDIF
517 heimbach 1.132 C
518 adcroft 1.58 CALL TIMESTEP(
519 jmc 1.63 I bi,bj,iMin,iMax,jMin,jMax,k,
520 jmc 1.94 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
521 jmc 1.110 I guDissip, gvDissip,
522 jmc 1.96 I myTime, myIter, myThid)
523 adcroft 1.58
524     #ifdef ALLOW_OBCS
525     C-- Apply open boundary conditions
526 jmc 1.96 IF (useOBCS) THEN
527     CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
528     ENDIF
529 adcroft 1.58 #endif /* ALLOW_OBCS */
530    
531     ENDIF
532    
533    
534     C-- end of dynamics k loop (1:Nr)
535     ENDDO
536    
537 jmc 1.106 C-- Implicit Vertical advection & viscosity
538 jmc 1.130 #if (defined (INCLUDE_IMPLVERTADV_CODE) && defined (ALLOW_MOM_COMMON))
539 jmc 1.106 IF ( momImplVertAdv ) THEN
540 jmc 1.136 CALL MOM_U_IMPLICIT_R( kappaRU,
541 jmc 1.106 I bi, bj, myTime, myIter, myThid )
542     CALL MOM_V_IMPLICIT_R( kappaRV,
543     I bi, bj, myTime, myIter, myThid )
544     ELSEIF ( implicitViscosity ) THEN
545     #else /* INCLUDE_IMPLVERTADV_CODE */
546     IF ( implicitViscosity ) THEN
547     #endif /* INCLUDE_IMPLVERTADV_CODE */
548 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
549 heimbach 1.101 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
550 jmc 1.96 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
551 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
552 adcroft 1.42 CALL IMPLDIFF(
553     I bi, bj, iMin, iMax, jMin, jMax,
554 jmc 1.124 I -1, KappaRU,recip_HFacW,
555 jmc 1.96 U gU,
556 adcroft 1.42 I myThid )
557 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
558 heimbach 1.101 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
559 heimbach 1.97 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
560 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
561 adcroft 1.42 CALL IMPLDIFF(
562     I bi, bj, iMin, iMax, jMin, jMax,
563 jmc 1.124 I -2, KappaRV,recip_HFacS,
564 jmc 1.96 U gV,
565 adcroft 1.42 I myThid )
566 jmc 1.106 ENDIF
567 heimbach 1.49
568 adcroft 1.58 #ifdef ALLOW_OBCS
569     C-- Apply open boundary conditions
570 jmc 1.106 IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
571 adcroft 1.58 DO K=1,Nr
572 jmc 1.96 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
573 adcroft 1.58 ENDDO
574 jmc 1.106 ENDIF
575 adcroft 1.58 #endif /* ALLOW_OBCS */
576 heimbach 1.49
577 edhill 1.102 #ifdef ALLOW_CD_CODE
578 jmc 1.106 IF (implicitViscosity.AND.useCDscheme) THEN
579 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
580 heimbach 1.91 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
581 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
582 adcroft 1.42 CALL IMPLDIFF(
583     I bi, bj, iMin, iMax, jMin, jMax,
584 jmc 1.111 I 0, KappaRU,recip_HFacW,
585 adcroft 1.42 U vVelD,
586     I myThid )
587 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
588 heimbach 1.91 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
589 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
590 adcroft 1.42 CALL IMPLDIFF(
591     I bi, bj, iMin, iMax, jMin, jMax,
592 jmc 1.111 I 0, KappaRV,recip_HFacS,
593 adcroft 1.42 U uVelD,
594     I myThid )
595 jmc 1.106 ENDIF
596 edhill 1.102 #endif /* ALLOW_CD_CODE */
597 jmc 1.106 C-- End implicit Vertical advection & viscosity
598 heimbach 1.109
599 jmc 1.113 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
600    
601 jmc 1.122 #ifdef ALLOW_NONHYDROSTATIC
602     C-- Step forward W field in N-H algorithm
603 jmc 1.136 IF ( nonHydrostatic ) THEN
604 jmc 1.122 #ifdef ALLOW_DEBUG
605 jmc 1.123 IF ( debugLevel .GE. debLevB )
606     & CALL DEBUG_CALL('CALC_GW', myThid )
607 jmc 1.122 #endif
608     CALL TIMER_START('CALC_GW [DYNAMICS]',myThid)
609 baylor 1.135 CALL CALC_GW(
610 jmc 1.136 I bi,bj, KappaRU, KappaRV,
611     I myTime, myIter, myThid )
612     ENDIF
613     IF ( nonHydrostatic.OR.implicitIntGravWave )
614     & CALL TIMESTEP_WVEL( bi,bj, myTime, myIter, myThid )
615     IF ( nonHydrostatic )
616 jmc 1.128 & CALL TIMER_STOP ('CALC_GW [DYNAMICS]',myThid)
617 jmc 1.122 #endif
618    
619     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
620    
621 jmc 1.136 C- end of bi,bj loops
622     ENDDO
623     ENDDO
624    
625     #ifdef ALLOW_OBCS
626     IF (useOBCS) THEN
627     CALL OBCS_PRESCRIBE_EXCHANGES(myThid)
628     ENDIF
629     #endif
630    
631 mlosch 1.90 Cml(
632     C In order to compare the variance of phiHydLow of a p/z-coordinate
633     C run with etaH of a z/p-coordinate run the drift of phiHydLow
634     C has to be removed by something like the following subroutine:
635 jmc 1.144 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskInC, maskInC, rA, drF,
636     C & 'phiHydLow', myTime, myThid )
637 mlosch 1.90 Cml)
638 adcroft 1.69
639 jmc 1.113 #ifdef ALLOW_DIAGNOSTICS
640 jmc 1.130 IF ( useDiagnostics ) THEN
641 jmc 1.113
642     CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
643 jmc 1.120 CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0, 1,0,1,1,myThid)
644 molod 1.116
645 jmc 1.120 tmpFac = 1. _d 0
646     CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
647     & 'PHIHYDSQ',0,Nr,0,1,1,myThid)
648 molod 1.116
649 jmc 1.120 CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
650     & 'PHIBOTSQ',0, 1,0,1,1,myThid)
651 jmc 1.113
652     ENDIF
653     #endif /* ALLOW_DIAGNOSTICS */
654 jmc 1.136
655 edhill 1.104 #ifdef ALLOW_DEBUG
656 heimbach 1.98 If ( debugLevel .GE. debLevB ) THEN
657 adcroft 1.69 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
658 adcroft 1.73 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
659 adcroft 1.69 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
660     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
661     CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
662     CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
663 jmc 1.115 CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
664     CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
665     CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
666     CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
667     #ifndef ALLOW_ADAMSBASHFORTH_3
668     CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
669     CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
670     CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
671     CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
672     #endif
673 adcroft 1.70 ENDIF
674 adcroft 1.69 #endif
675 cnh 1.1
676 jmc 1.125 #ifdef DYNAMICS_GUGV_EXCH_CHECK
677 jmc 1.144 C- jmc: For safety checking only: This Exchange here should not change
678     C the solution. If solution changes, it means something is wrong,
679 jmc 1.125 C but it does not mean that it is less wrong with this exchange.
680     IF ( debugLevel .GT. debLevB ) THEN
681     CALL EXCH_UV_XYZ_RL(gU,gV,.TRUE.,myThid)
682     ENDIF
683     #endif
684    
685 jmc 1.123 #ifdef ALLOW_DEBUG
686     IF ( debugLevel .GE. debLevB )
687     & CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
688     #endif
689    
690 cnh 1.1 RETURN
691     END

  ViewVC Help
Powered by ViewVC 1.1.22