/[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.141 - (hide annotations) (download)
Fri Feb 13 21:56:48 2009 UTC (15 years, 3 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint61l, checkpoint61j, checkpoint61k
Changes since 1.140: +2 -2 lines
Add TAF option "kind" (or adjust "byte") to enable real*4 common blocks

1 heimbach 1.141 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.140 2008/10/20 23:51:39 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 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     C | SUBROUTINE DYNAMICS
19     C | o Controlling routine for the explicit part of the model
20     C | dynamics.
21     C *==========================================================*
22     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     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     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     # ifdef NONLIN_FRSURF
321     # ifndef DISABLE_RSTAR_CODE
322     dWtransC(i,j,bi,bj) = 0. _d 0
323     dWtransU(i,j,bi,bj) = 0. _d 0
324     dWtransV(i,j,bi,bj) = 0. _d 0
325     # endif
326     # endif
327     #endif
328 heimbach 1.76 ENDDO
329     ENDDO
330 heimbach 1.49
331 jmc 1.63 C-- Start computation of dynamics
332 jmc 1.93 iMin = 0
333     iMax = sNx+1
334     jMin = 0
335     jMax = sNy+1
336 jmc 1.63
337 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
338 heimbach 1.91 CADJ STORE wvel (:,:,:,bi,bj) =
339 heimbach 1.141 CADJ & comlev1_bibj, key=idynkey, byte=isbyte
340 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
341    
342 jmc 1.65 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
343 jmc 1.63 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
344     IF (implicSurfPress.NE.1.) THEN
345 jmc 1.65 CALL CALC_GRAD_PHI_SURF(
346     I bi,bj,iMin,iMax,jMin,jMax,
347     I etaN,
348     O phiSurfX,phiSurfY,
349 jmc 1.136 I myThid )
350 jmc 1.63 ENDIF
351 heimbach 1.83
352     #ifdef ALLOW_AUTODIFF_TAMC
353 heimbach 1.91 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
354     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
355 heimbach 1.83 #ifdef ALLOW_KPP
356     CADJ STORE KPPviscAz (:,:,:,bi,bj)
357 heimbach 1.91 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
358 heimbach 1.83 #endif /* ALLOW_KPP */
359     #endif /* ALLOW_AUTODIFF_TAMC */
360 adcroft 1.58
361 heimbach 1.77 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
362 jmc 1.140 C-- Calculate the total vertical viscosity
363     CALL CALC_VISCOSITY(
364     I bi,bj, iMin,iMax,jMin,jMax,
365     O KappaRU, KappaRV,
366     I myThid )
367     #else
368 heimbach 1.77 DO k=1,Nr
369 jmc 1.140 DO j=1-OLy,sNy+OLy
370     DO i=1-OLx,sNx+OLx
371     KappaRU(i,j,k) = 0. _d 0
372     KappaRV(i,j,k) = 0. _d 0
373     ENDDO
374     ENDDO
375     ENDDO
376 heimbach 1.77 #endif
377    
378 heimbach 1.101 #ifdef ALLOW_AUTODIFF_TAMC
379     CADJ STORE KappaRU(:,:,:)
380 heimbach 1.132 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
381 heimbach 1.101 CADJ STORE KappaRV(:,:,:)
382 heimbach 1.132 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
383 heimbach 1.101 #endif /* ALLOW_AUTODIFF_TAMC */
384    
385 adcroft 1.58 C-- Start of dynamics loop
386     DO k=1,Nr
387    
388     C-- km1 Points to level above k (=k-1)
389     C-- kup Cycles through 1,2 to point to layer above
390     C-- kDown Cycles through 2,1 to point to current layer
391    
392     km1 = MAX(1,k-1)
393 heimbach 1.77 kp1 = MIN(k+1,Nr)
394 adcroft 1.58 kup = 1+MOD(k+1,2)
395     kDown= 1+MOD(k,2)
396    
397 heimbach 1.76 #ifdef ALLOW_AUTODIFF_TAMC
398 heimbach 1.91 kkey = (idynkey-1)*Nr + k
399 heimbach 1.99 c
400 heimbach 1.95 CADJ STORE totphihyd (:,:,k,bi,bj)
401 heimbach 1.99 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
402     CADJ STORE theta (:,:,k,bi,bj)
403     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
404     CADJ STORE salt (:,:,k,bi,bj)
405 heimbach 1.95 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
406 heimbach 1.129 CADJ STORE gt(:,:,k,bi,bj)
407     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
408     CADJ STORE gs(:,:,k,bi,bj)
409     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
410 heimbach 1.126 # ifdef NONLIN_FRSURF
411     cph-test
412     CADJ STORE phiHydC (:,:)
413     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
414     CADJ STORE phiHydF (:,:)
415     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
416     CADJ STORE gudissip (:,:)
417     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
418     CADJ STORE gvdissip (:,:)
419     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
420     CADJ STORE fVerU (:,:,:)
421     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
422     CADJ STORE fVerV (:,:,:)
423     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
424     CADJ STORE gu(:,:,k,bi,bj)
425     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
426     CADJ STORE gv(:,:,k,bi,bj)
427     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
428     CADJ STORE gunm1(:,:,k,bi,bj)
429     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
430     CADJ STORE gvnm1(:,:,k,bi,bj)
431     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
432     # ifdef ALLOW_CD_CODE
433     CADJ STORE unm1(:,:,k,bi,bj)
434     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
435     CADJ STORE vnm1(:,:,k,bi,bj)
436     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
437     CADJ STORE uVelD(:,:,k,bi,bj)
438     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
439     CADJ STORE vVelD(:,:,k,bi,bj)
440     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
441     # endif
442     # endif
443 heimbach 1.134 # ifdef ALLOW_DEPTH_CONTROL
444     CADJ STORE fVerU (:,:,:)
445     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
446     CADJ STORE fVerV (:,:,:)
447     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
448     # endif
449 heimbach 1.76 #endif /* ALLOW_AUTODIFF_TAMC */
450    
451 adcroft 1.58 C-- Integrate hydrostatic balance for phiHyd with BC of
452     C phiHyd(z=0)=0
453 jmc 1.128 IF ( implicitIntGravWave ) THEN
454     CALL CALC_PHI_HYD(
455     I bi,bj,iMin,iMax,jMin,jMax,k,
456     I gT, gS,
457     U phiHydF,
458     O phiHydC, dPhiHydX, dPhiHydY,
459     I myTime, myIter, myThid )
460     ELSE
461     CALL CALC_PHI_HYD(
462 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,
463     I theta, salt,
464 jmc 1.94 U phiHydF,
465     O phiHydC, dPhiHydX, dPhiHydY,
466 jmc 1.92 I myTime, myIter, myThid )
467 jmc 1.128 ENDIF
468 mlosch 1.89
469 adcroft 1.58 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
470 jmc 1.96 C and step forward storing the result in gU, gV, etc...
471 adcroft 1.58 IF ( momStepping ) THEN
472 heimbach 1.138 #ifdef ALLOW_AUTODIFF_TAMC
473     # ifdef NONLIN_FRSURF
474     # ifndef DISABLE_RSTAR_CODE
475     CADJ STORE dWtransC(:,:,bi,bj)
476     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
477     CADJ STORE dWtransU(:,:,bi,bj)
478     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
479     CADJ STORE dWtransV(:,:,bi,bj)
480     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
481     # endif
482     # endif
483     #endif
484 heimbach 1.132 IF (.NOT. vectorInvariantMomentum) THEN
485 edhill 1.105 #ifdef ALLOW_MOM_FLUXFORM
486 heimbach 1.132 C
487     CALL MOM_FLUXFORM(
488 adcroft 1.58 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
489 jmc 1.121 I KappaRU, KappaRV,
490 adcroft 1.58 U fVerU, fVerV,
491 jmc 1.121 O guDissip, gvDissip,
492 adcroft 1.80 I myTime, myIter, myThid)
493 adcroft 1.79 #endif
494 heimbach 1.132 ELSE
495 edhill 1.105 #ifdef ALLOW_MOM_VECINV
496 heimbach 1.126 C
497     # ifdef ALLOW_AUTODIFF_TAMC
498     # ifdef NONLIN_FRSURF
499     CADJ STORE fVerU(:,:,:)
500     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
501     CADJ STORE fVerV(:,:,:)
502     CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
503     # endif
504     # endif /* ALLOW_AUTODIFF_TAMC */
505     C
506     CALL MOM_VECINV(
507 adcroft 1.79 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
508 jmc 1.121 I KappaRU, KappaRV,
509 adcroft 1.79 U fVerU, fVerV,
510 jmc 1.110 O guDissip, gvDissip,
511 adcroft 1.80 I myTime, myIter, myThid)
512 heimbach 1.132 #endif
513 heimbach 1.126 ENDIF
514 heimbach 1.132 C
515 adcroft 1.58 CALL TIMESTEP(
516 jmc 1.63 I bi,bj,iMin,iMax,jMin,jMax,k,
517 jmc 1.94 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
518 jmc 1.110 I guDissip, gvDissip,
519 jmc 1.96 I myTime, myIter, myThid)
520 adcroft 1.58
521     #ifdef ALLOW_OBCS
522     C-- Apply open boundary conditions
523 jmc 1.96 IF (useOBCS) THEN
524     CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
525     ENDIF
526 adcroft 1.58 #endif /* ALLOW_OBCS */
527    
528     ENDIF
529    
530    
531     C-- end of dynamics k loop (1:Nr)
532     ENDDO
533    
534 jmc 1.106 C-- Implicit Vertical advection & viscosity
535 jmc 1.130 #if (defined (INCLUDE_IMPLVERTADV_CODE) && defined (ALLOW_MOM_COMMON))
536 jmc 1.106 IF ( momImplVertAdv ) THEN
537 jmc 1.136 CALL MOM_U_IMPLICIT_R( kappaRU,
538 jmc 1.106 I bi, bj, myTime, myIter, myThid )
539     CALL MOM_V_IMPLICIT_R( kappaRV,
540     I bi, bj, myTime, myIter, myThid )
541     ELSEIF ( implicitViscosity ) THEN
542     #else /* INCLUDE_IMPLVERTADV_CODE */
543     IF ( implicitViscosity ) THEN
544     #endif /* INCLUDE_IMPLVERTADV_CODE */
545 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
546 heimbach 1.101 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
547 jmc 1.96 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
548 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
549 adcroft 1.42 CALL IMPLDIFF(
550     I bi, bj, iMin, iMax, jMin, jMax,
551 jmc 1.124 I -1, KappaRU,recip_HFacW,
552 jmc 1.96 U gU,
553 adcroft 1.42 I myThid )
554 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
555 heimbach 1.101 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
556 heimbach 1.97 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
557 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
558 adcroft 1.42 CALL IMPLDIFF(
559     I bi, bj, iMin, iMax, jMin, jMax,
560 jmc 1.124 I -2, KappaRV,recip_HFacS,
561 jmc 1.96 U gV,
562 adcroft 1.42 I myThid )
563 jmc 1.106 ENDIF
564 heimbach 1.49
565 adcroft 1.58 #ifdef ALLOW_OBCS
566     C-- Apply open boundary conditions
567 jmc 1.106 IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
568 adcroft 1.58 DO K=1,Nr
569 jmc 1.96 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
570 adcroft 1.58 ENDDO
571 jmc 1.106 ENDIF
572 adcroft 1.58 #endif /* ALLOW_OBCS */
573 heimbach 1.49
574 edhill 1.102 #ifdef ALLOW_CD_CODE
575 jmc 1.106 IF (implicitViscosity.AND.useCDscheme) THEN
576 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
577 heimbach 1.91 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
578 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
579 adcroft 1.42 CALL IMPLDIFF(
580     I bi, bj, iMin, iMax, jMin, jMax,
581 jmc 1.111 I 0, KappaRU,recip_HFacW,
582 adcroft 1.42 U vVelD,
583     I myThid )
584 adcroft 1.58 #ifdef ALLOW_AUTODIFF_TAMC
585 heimbach 1.91 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
586 adcroft 1.58 #endif /* ALLOW_AUTODIFF_TAMC */
587 adcroft 1.42 CALL IMPLDIFF(
588     I bi, bj, iMin, iMax, jMin, jMax,
589 jmc 1.111 I 0, KappaRV,recip_HFacS,
590 adcroft 1.42 U uVelD,
591     I myThid )
592 jmc 1.106 ENDIF
593 edhill 1.102 #endif /* ALLOW_CD_CODE */
594 jmc 1.106 C-- End implicit Vertical advection & viscosity
595 heimbach 1.109
596 jmc 1.113 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
597    
598 jmc 1.122 #ifdef ALLOW_NONHYDROSTATIC
599     C-- Step forward W field in N-H algorithm
600 jmc 1.136 IF ( nonHydrostatic ) THEN
601 jmc 1.122 #ifdef ALLOW_DEBUG
602 jmc 1.123 IF ( debugLevel .GE. debLevB )
603     & CALL DEBUG_CALL('CALC_GW', myThid )
604 jmc 1.122 #endif
605     CALL TIMER_START('CALC_GW [DYNAMICS]',myThid)
606 baylor 1.135 CALL CALC_GW(
607 jmc 1.136 I bi,bj, KappaRU, KappaRV,
608     I myTime, myIter, myThid )
609     ENDIF
610     IF ( nonHydrostatic.OR.implicitIntGravWave )
611     & CALL TIMESTEP_WVEL( bi,bj, myTime, myIter, myThid )
612     IF ( nonHydrostatic )
613 jmc 1.128 & CALL TIMER_STOP ('CALC_GW [DYNAMICS]',myThid)
614 jmc 1.122 #endif
615    
616     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
617    
618 jmc 1.136 C- end of bi,bj loops
619     ENDDO
620     ENDDO
621    
622     #ifdef ALLOW_OBCS
623     IF (useOBCS) THEN
624     CALL OBCS_PRESCRIBE_EXCHANGES(myThid)
625     ENDIF
626     #endif
627    
628 mlosch 1.90 Cml(
629     C In order to compare the variance of phiHydLow of a p/z-coordinate
630     C run with etaH of a z/p-coordinate run the drift of phiHydLow
631     C has to be removed by something like the following subroutine:
632     C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
633     C & 'phiHydLow', myThid )
634     Cml)
635 adcroft 1.69
636 jmc 1.113 #ifdef ALLOW_DIAGNOSTICS
637 jmc 1.130 IF ( useDiagnostics ) THEN
638 jmc 1.113
639     CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
640 jmc 1.120 CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0, 1,0,1,1,myThid)
641 molod 1.116
642 jmc 1.120 tmpFac = 1. _d 0
643     CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
644     & 'PHIHYDSQ',0,Nr,0,1,1,myThid)
645 molod 1.116
646 jmc 1.120 CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
647     & 'PHIBOTSQ',0, 1,0,1,1,myThid)
648 jmc 1.113
649     ENDIF
650     #endif /* ALLOW_DIAGNOSTICS */
651 jmc 1.136
652 edhill 1.104 #ifdef ALLOW_DEBUG
653 heimbach 1.98 If ( debugLevel .GE. debLevB ) THEN
654 adcroft 1.69 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
655 adcroft 1.73 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
656 adcroft 1.69 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
657     CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
658     CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
659     CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
660 jmc 1.115 CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
661     CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
662     CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
663     CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
664     #ifndef ALLOW_ADAMSBASHFORTH_3
665     CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
666     CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
667     CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
668     CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
669     #endif
670 adcroft 1.70 ENDIF
671 adcroft 1.69 #endif
672 cnh 1.1
673 jmc 1.125 #ifdef DYNAMICS_GUGV_EXCH_CHECK
674     C- jmc: For safety checking only: This Exchange here should not change
675     C the solution. If solution changes, it means something is wrong,
676     C but it does not mean that it is less wrong with this exchange.
677     IF ( debugLevel .GT. debLevB ) THEN
678     CALL EXCH_UV_XYZ_RL(gU,gV,.TRUE.,myThid)
679     ENDIF
680     #endif
681    
682 jmc 1.123 #ifdef ALLOW_DEBUG
683     IF ( debugLevel .GE. debLevB )
684     & CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
685     #endif
686    
687 cnh 1.1 RETURN
688     END

  ViewVC Help
Powered by ViewVC 1.1.22