/[MITgcm]/MITgcm/model/src/dynamics.F
ViewVC logotype

Contents of /MITgcm/model/src/dynamics.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.122 - (show annotations) (download)
Sat Jul 30 23:39:48 2005 UTC (18 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57p_post
Changes since 1.121: +28 -4 lines
call CALC_GW from DYNAMICS (instead of from FORWARD_STEP)

1 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.121 2005/07/30 22:09:38 jmc Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 CBOP
8 C !ROUTINE: DYNAMICS
9 C !INTERFACE:
10 SUBROUTINE DYNAMICS(myTime, myIter, myThid)
11 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 | W[n] = W* + dt x d/dz P (NH mode)
34 C |
35 C | "Calculation of Gs"
36 C | ===================
37 C | This is where all the accelerations and tendencies (ie.
38 C | physics, parameterizations etc...) are calculated
39 C | rho = rho ( theta[n], salt[n] )
40 C | b = b(rho, theta)
41 C | K31 = K31 ( rho )
42 C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
43 C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
44 C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
45 C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
46 C |
47 C | "Time-stepping" or "Prediction"
48 C | ================================
49 C | The models variables are stepped forward with the appropriate
50 C | time-stepping scheme (currently we use Adams-Bashforth II)
51 C | - For momentum, the result is always *only* a "prediction"
52 C | in that the flow may be divergent and will be "corrected"
53 C | later with a surface pressure gradient.
54 C | - Normally for tracers the result is the new field at time
55 C | level [n+1} *BUT* in the case of implicit diffusion the result
56 C | is also *only* a prediction.
57 C | - We denote "predictors" with an asterisk (*).
58 C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
59 C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
60 C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
61 C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
62 C | With implicit diffusion:
63 C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
64 C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
65 C | (1 + dt * K * d_zz) theta[n] = theta*
66 C | (1 + dt * K * d_zz) salt[n] = salt*
67 C |
68 C *==========================================================*
69 C \ev
70 C !USES:
71 IMPLICIT NONE
72 C == Global variables ===
73 #include "SIZE.h"
74 #include "EEPARAMS.h"
75 #include "PARAMS.h"
76 #include "DYNVARS.h"
77 #ifdef ALLOW_CD_CODE
78 #include "CD_CODE_VARS.h"
79 #endif
80 #include "GRID.h"
81 #ifdef ALLOW_AUTODIFF_TAMC
82 # include "tamc.h"
83 # include "tamc_keys.h"
84 # include "FFIELDS.h"
85 # include "EOS.h"
86 # ifdef ALLOW_KPP
87 # include "KPP.h"
88 # endif
89 #endif /* ALLOW_AUTODIFF_TAMC */
90
91 C !CALLING SEQUENCE:
92 C DYNAMICS()
93 C |
94 C |-- CALC_EP_FORCING
95 C |
96 C |-- CALC_GRAD_PHI_SURF
97 C |
98 C |-- CALC_VISCOSITY
99 C |
100 C |-- CALC_PHI_HYD
101 C |
102 C |-- MOM_FLUXFORM
103 C |
104 C |-- MOM_VECINV
105 C |
106 C |-- TIMESTEP
107 C |
108 C |-- OBCS_APPLY_UV
109 C |
110 C |-- MOM_U_IMPLICIT_R
111 C |-- MOM_V_IMPLICIT_R
112 C |
113 C |-- IMPLDIFF
114 C |
115 C |-- OBCS_APPLY_UV
116 C |
117 C |-- CALC_GW
118 C |
119 C |-- DIAGNOSTICS_FILL
120 C |-- DEBUG_STATS_RL
121
122 C !INPUT/OUTPUT PARAMETERS:
123 C == Routine arguments ==
124 C myTime - Current time in simulation
125 C myIter - Current iteration number in simulation
126 C myThid - Thread number for this instance of the routine.
127 _RL myTime
128 INTEGER myIter
129 INTEGER myThid
130
131 C !LOCAL VARIABLES:
132 C == Local variables
133 C fVer[UV] o fVer: Vertical flux term - note fVer
134 C is "pipelined" in the vertical
135 C so we need an fVer for each
136 C variable.
137 C phiHydC :: hydrostatic potential anomaly at cell center
138 C In z coords phiHyd is the hydrostatic potential
139 C (=pressure/rho0) anomaly
140 C In p coords phiHyd is the geopotential height anomaly.
141 C phiHydF :: hydrostatic potential anomaly at middle between 2 centers
142 C dPhiHydX,Y :: Gradient (X & Y directions) of hydrostatic potential anom.
143 C phiSurfX, :: gradient of Surface potential (Pressure/rho, ocean)
144 C phiSurfY or geopotential (atmos) in X and Y direction
145 C guDissip :: dissipation tendency (all explicit terms), u component
146 C gvDissip :: dissipation tendency (all explicit terms), v component
147 C iMin, iMax - Ranges and sub-block indices on which calculations
148 C jMin, jMax are applied.
149 C bi, bj
150 C k, kup, - Index for layer above and below. kup and kDown
151 C kDown, km1 are switched with layer to be the appropriate
152 C index into fVerTerm.
153 _RL fVerU (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
154 _RL fVerV (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
155 _RL phiHydF (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
156 _RL phiHydC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
157 _RL dPhiHydX(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
158 _RL dPhiHydY(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
159 _RL phiSurfX(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
160 _RL phiSurfY(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
161 _RL guDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
162 _RL gvDissip(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
163 _RL KappaRU (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
164 _RL KappaRV (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr)
165
166 INTEGER iMin, iMax
167 INTEGER jMin, jMax
168 INTEGER bi, bj
169 INTEGER i, j
170 INTEGER k, km1, kp1, kup, kDown
171
172 #ifdef ALLOW_DIAGNOSTICS
173 _RL tmpFac
174 #endif /* ALLOW_DIAGNOSTICS */
175
176
177 C--- The algorithm...
178 C
179 C "Correction Step"
180 C =================
181 C Here we update the horizontal velocities with the surface
182 C pressure such that the resulting flow is either consistent
183 C with the free-surface evolution or the rigid-lid:
184 C U[n] = U* + dt x d/dx P
185 C V[n] = V* + dt x d/dy P
186 C
187 C "Calculation of Gs"
188 C ===================
189 C This is where all the accelerations and tendencies (ie.
190 C physics, parameterizations etc...) are calculated
191 C rho = rho ( theta[n], salt[n] )
192 C b = b(rho, theta)
193 C K31 = K31 ( rho )
194 C Gu[n] = Gu( u[n], v[n], wVel, b, ... )
195 C Gv[n] = Gv( u[n], v[n], wVel, b, ... )
196 C Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
197 C Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
198 C
199 C "Time-stepping" or "Prediction"
200 C ================================
201 C The models variables are stepped forward with the appropriate
202 C time-stepping scheme (currently we use Adams-Bashforth II)
203 C - For momentum, the result is always *only* a "prediction"
204 C in that the flow may be divergent and will be "corrected"
205 C later with a surface pressure gradient.
206 C - Normally for tracers the result is the new field at time
207 C level [n+1} *BUT* in the case of implicit diffusion the result
208 C is also *only* a prediction.
209 C - We denote "predictors" with an asterisk (*).
210 C U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
211 C V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
212 C theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
213 C salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
214 C With implicit diffusion:
215 C theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
216 C salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
217 C (1 + dt * K * d_zz) theta[n] = theta*
218 C (1 + dt * K * d_zz) salt[n] = salt*
219 C---
220 CEOP
221
222 C-- Call to routine for calculation of
223 C Eliassen-Palm-flux-forced U-tendency,
224 C if desired:
225 #ifdef INCLUDE_EP_FORCING_CODE
226 CALL CALC_EP_FORCING(myThid)
227 #endif
228
229 #ifdef ALLOW_AUTODIFF_TAMC
230 C-- HPF directive to help TAMC
231 CHPF$ INDEPENDENT
232 #endif /* ALLOW_AUTODIFF_TAMC */
233
234 DO bj=myByLo(myThid),myByHi(myThid)
235
236 #ifdef ALLOW_AUTODIFF_TAMC
237 C-- HPF directive to help TAMC
238 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
239 CHPF$& ,phiHydF
240 CHPF$& ,KappaRU,KappaRV
241 CHPF$& )
242 #endif /* ALLOW_AUTODIFF_TAMC */
243
244 DO bi=myBxLo(myThid),myBxHi(myThid)
245
246 #ifdef ALLOW_AUTODIFF_TAMC
247 act1 = bi - myBxLo(myThid)
248 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
249 act2 = bj - myByLo(myThid)
250 max2 = myByHi(myThid) - myByLo(myThid) + 1
251 act3 = myThid - 1
252 max3 = nTx*nTy
253 act4 = ikey_dynamics - 1
254 idynkey = (act1 + 1) + act2*max1
255 & + act3*max1*max2
256 & + act4*max1*max2*max3
257 #endif /* ALLOW_AUTODIFF_TAMC */
258
259 C-- Set up work arrays with valid (i.e. not NaN) values
260 C These inital values do not alter the numerical results. They
261 C just ensure that all memory references are to valid floating
262 C point numbers. This prevents spurious hardware signals due to
263 C uninitialised but inert locations.
264
265 DO k=1,Nr
266 DO j=1-OLy,sNy+OLy
267 DO i=1-OLx,sNx+OLx
268 KappaRU(i,j,k) = 0. _d 0
269 KappaRV(i,j,k) = 0. _d 0
270 #ifdef ALLOW_AUTODIFF_TAMC
271 cph(
272 c-- need some re-initialisation here to break dependencies
273 cph)
274 gU(i,j,k,bi,bj) = 0. _d 0
275 gV(i,j,k,bi,bj) = 0. _d 0
276 #endif
277 ENDDO
278 ENDDO
279 ENDDO
280 DO j=1-OLy,sNy+OLy
281 DO i=1-OLx,sNx+OLx
282 fVerU (i,j,1) = 0. _d 0
283 fVerU (i,j,2) = 0. _d 0
284 fVerV (i,j,1) = 0. _d 0
285 fVerV (i,j,2) = 0. _d 0
286 phiHydF (i,j) = 0. _d 0
287 phiHydC (i,j) = 0. _d 0
288 dPhiHydX(i,j) = 0. _d 0
289 dPhiHydY(i,j) = 0. _d 0
290 phiSurfX(i,j) = 0. _d 0
291 phiSurfY(i,j) = 0. _d 0
292 guDissip(i,j) = 0. _d 0
293 gvDissip(i,j) = 0. _d 0
294 ENDDO
295 ENDDO
296
297 C-- Start computation of dynamics
298 iMin = 0
299 iMax = sNx+1
300 jMin = 0
301 jMax = sNy+1
302
303 #ifdef ALLOW_AUTODIFF_TAMC
304 CADJ STORE wvel (:,:,:,bi,bj) =
305 CADJ & comlev1_bibj, key = idynkey, byte = isbyte
306 #endif /* ALLOW_AUTODIFF_TAMC */
307
308 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
309 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
310 IF (implicSurfPress.NE.1.) THEN
311 CALL CALC_GRAD_PHI_SURF(
312 I bi,bj,iMin,iMax,jMin,jMax,
313 I etaN,
314 O phiSurfX,phiSurfY,
315 I myThid )
316 ENDIF
317
318 #ifdef ALLOW_AUTODIFF_TAMC
319 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
320 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
321 #ifdef ALLOW_KPP
322 CADJ STORE KPPviscAz (:,:,:,bi,bj)
323 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
324 #endif /* ALLOW_KPP */
325 #endif /* ALLOW_AUTODIFF_TAMC */
326
327 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
328 C-- Calculate the total vertical diffusivity
329 DO k=1,Nr
330 CALL CALC_VISCOSITY(
331 I bi,bj,iMin,iMax,jMin,jMax,k,
332 O KappaRU,KappaRV,
333 I myThid)
334 ENDDO
335 #endif
336
337 #ifdef ALLOW_AUTODIFF_TAMC
338 CADJ STORE KappaRU(:,:,:)
339 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
340 CADJ STORE KappaRV(:,:,:)
341 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
342 #endif /* ALLOW_AUTODIFF_TAMC */
343
344 C-- Start of dynamics loop
345 DO k=1,Nr
346
347 C-- km1 Points to level above k (=k-1)
348 C-- kup Cycles through 1,2 to point to layer above
349 C-- kDown Cycles through 2,1 to point to current layer
350
351 km1 = MAX(1,k-1)
352 kp1 = MIN(k+1,Nr)
353 kup = 1+MOD(k+1,2)
354 kDown= 1+MOD(k,2)
355
356 #ifdef ALLOW_AUTODIFF_TAMC
357 kkey = (idynkey-1)*Nr + k
358 c
359 CADJ STORE totphihyd (:,:,k,bi,bj)
360 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
361 CADJ STORE theta (:,:,k,bi,bj)
362 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
363 CADJ STORE salt (:,:,k,bi,bj)
364 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
365 #endif /* ALLOW_AUTODIFF_TAMC */
366
367 C-- Integrate hydrostatic balance for phiHyd with BC of
368 C phiHyd(z=0)=0
369 CALL CALC_PHI_HYD(
370 I bi,bj,iMin,iMax,jMin,jMax,k,
371 I theta, salt,
372 U phiHydF,
373 O phiHydC, dPhiHydX, dPhiHydY,
374 I myTime, myIter, myThid )
375
376 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
377 C and step forward storing the result in gU, gV, etc...
378 IF ( momStepping ) THEN
379 #ifdef ALLOW_MOM_FLUXFORM
380 IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
381 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
382 I KappaRU, KappaRV,
383 U fVerU, fVerV,
384 O guDissip, gvDissip,
385 I myTime, myIter, myThid)
386 #endif
387 #ifdef ALLOW_MOM_VECINV
388 IF (vectorInvariantMomentum) CALL MOM_VECINV(
389 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
390 I KappaRU, KappaRV,
391 U fVerU, fVerV,
392 O guDissip, gvDissip,
393 I myTime, myIter, myThid)
394 #endif
395 CALL TIMESTEP(
396 I bi,bj,iMin,iMax,jMin,jMax,k,
397 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
398 I guDissip, gvDissip,
399 I myTime, myIter, myThid)
400
401 #ifdef ALLOW_OBCS
402 C-- Apply open boundary conditions
403 IF (useOBCS) THEN
404 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
405 ENDIF
406 #endif /* ALLOW_OBCS */
407
408 ENDIF
409
410
411 C-- end of dynamics k loop (1:Nr)
412 ENDDO
413
414 C-- Implicit Vertical advection & viscosity
415 #ifdef INCLUDE_IMPLVERTADV_CODE
416 IF ( momImplVertAdv ) THEN
417 CALL MOM_U_IMPLICIT_R( kappaRU,
418 I bi, bj, myTime, myIter, myThid )
419 CALL MOM_V_IMPLICIT_R( kappaRV,
420 I bi, bj, myTime, myIter, myThid )
421 ELSEIF ( implicitViscosity ) THEN
422 #else /* INCLUDE_IMPLVERTADV_CODE */
423 IF ( implicitViscosity ) THEN
424 #endif /* INCLUDE_IMPLVERTADV_CODE */
425 #ifdef ALLOW_AUTODIFF_TAMC
426 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
427 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
428 #endif /* ALLOW_AUTODIFF_TAMC */
429 CALL IMPLDIFF(
430 I bi, bj, iMin, iMax, jMin, jMax,
431 I 0, KappaRU,recip_HFacW,
432 U gU,
433 I myThid )
434 #ifdef ALLOW_AUTODIFF_TAMC
435 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
436 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
437 #endif /* ALLOW_AUTODIFF_TAMC */
438 CALL IMPLDIFF(
439 I bi, bj, iMin, iMax, jMin, jMax,
440 I 0, KappaRV,recip_HFacS,
441 U gV,
442 I myThid )
443 ENDIF
444
445 #ifdef ALLOW_OBCS
446 C-- Apply open boundary conditions
447 IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
448 DO K=1,Nr
449 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
450 ENDDO
451 ENDIF
452 #endif /* ALLOW_OBCS */
453
454 #ifdef ALLOW_CD_CODE
455 IF (implicitViscosity.AND.useCDscheme) THEN
456 #ifdef ALLOW_AUTODIFF_TAMC
457 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
458 #endif /* ALLOW_AUTODIFF_TAMC */
459 CALL IMPLDIFF(
460 I bi, bj, iMin, iMax, jMin, jMax,
461 I 0, KappaRU,recip_HFacW,
462 U vVelD,
463 I myThid )
464 #ifdef ALLOW_AUTODIFF_TAMC
465 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
466 #endif /* ALLOW_AUTODIFF_TAMC */
467 CALL IMPLDIFF(
468 I bi, bj, iMin, iMax, jMin, jMax,
469 I 0, KappaRV,recip_HFacS,
470 U uVelD,
471 I myThid )
472 ENDIF
473 #endif /* ALLOW_CD_CODE */
474 C-- End implicit Vertical advection & viscosity
475
476 ENDDO
477 ENDDO
478
479 #ifdef ALLOW_OBCS
480 IF (useOBCS) THEN
481 CALL OBCS_PRESCRIBE_EXCHANGES(myThid)
482 ENDIF
483 #endif
484
485 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
486
487 #ifdef ALLOW_NONHYDROSTATIC
488 C-- Step forward W field in N-H algorithm
489 IF ( momStepping .AND. nonHydrostatic ) THEN
490 #ifdef ALLOW_DEBUG
491 IF ( debugLevel .GE. debLevB )
492 & CALL DEBUG_CALL('CALC_GW',myThid)
493 #endif
494 CALL TIMER_START('CALC_GW [DYNAMICS]',myThid)
495 CALL CALC_GW( myTime, myIter, myThid )
496 CALL TIMER_STOP ('CALC_GW [DYNAMICS]',myThid)
497 ENDIF
498 #endif
499
500 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
501
502 Cml(
503 C In order to compare the variance of phiHydLow of a p/z-coordinate
504 C run with etaH of a z/p-coordinate run the drift of phiHydLow
505 C has to be removed by something like the following subroutine:
506 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
507 C & 'phiHydLow', myThid )
508 Cml)
509
510 #ifdef ALLOW_DIAGNOSTICS
511 IF ( usediagnostics ) THEN
512
513 CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
514 CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0, 1,0,1,1,myThid)
515
516 tmpFac = 1. _d 0
517 CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
518 & 'PHIHYDSQ',0,Nr,0,1,1,myThid)
519
520 CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
521 & 'PHIBOTSQ',0, 1,0,1,1,myThid)
522
523 ENDIF
524 #endif /* ALLOW_DIAGNOSTICS */
525
526 #ifdef ALLOW_DEBUG
527 If ( debugLevel .GE. debLevB ) THEN
528 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
529 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
530 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
531 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
532 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
533 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
534 CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
535 CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
536 CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
537 CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
538 #ifndef ALLOW_ADAMSBASHFORTH_3
539 CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
540 CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
541 CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
542 CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
543 #endif
544 ENDIF
545 #endif
546
547 RETURN
548 END

  ViewVC Help
Powered by ViewVC 1.1.22