/[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.158 - (hide annotations) (download)
Wed Jun 8 01:21:14 2011 UTC (12 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62z, checkpoint63, checkpoint63d, checkpoint63a, checkpoint63b, checkpoint63c
Changes since 1.157: +3 -3 lines
refine debugLevel criteria when printing messages

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

  ViewVC Help
Powered by ViewVC 1.1.22