/[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.124 - (show annotations) (download)
Sun Sep 11 20:52:09 2005 UTC (18 years, 8 months ago) by jmc
Branch: MAIN
Changes since 1.123: +3 -3 lines
call impldiff with tracerId=-1,-2 for gU,gV resp. (to do diagnostics inside)

1 C $Header: /u/gcmpack/MITgcm/model/src/dynamics.F,v 1.123 2005/08/16 22:52:06 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 #ifdef ALLOW_DEBUG
223 IF ( debugLevel .GE. debLevB )
224 & CALL DEBUG_ENTER( 'DYNAMICS', myThid )
225 #endif
226
227 C-- Call to routine for calculation of
228 C Eliassen-Palm-flux-forced U-tendency,
229 C if desired:
230 #ifdef INCLUDE_EP_FORCING_CODE
231 CALL CALC_EP_FORCING(myThid)
232 #endif
233
234 #ifdef ALLOW_AUTODIFF_TAMC
235 C-- HPF directive to help TAMC
236 CHPF$ INDEPENDENT
237 #endif /* ALLOW_AUTODIFF_TAMC */
238
239 DO bj=myByLo(myThid),myByHi(myThid)
240
241 #ifdef ALLOW_AUTODIFF_TAMC
242 C-- HPF directive to help TAMC
243 CHPF$ INDEPENDENT, NEW (fVerU,fVerV
244 CHPF$& ,phiHydF
245 CHPF$& ,KappaRU,KappaRV
246 CHPF$& )
247 #endif /* ALLOW_AUTODIFF_TAMC */
248
249 DO bi=myBxLo(myThid),myBxHi(myThid)
250
251 #ifdef ALLOW_AUTODIFF_TAMC
252 act1 = bi - myBxLo(myThid)
253 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
254 act2 = bj - myByLo(myThid)
255 max2 = myByHi(myThid) - myByLo(myThid) + 1
256 act3 = myThid - 1
257 max3 = nTx*nTy
258 act4 = ikey_dynamics - 1
259 idynkey = (act1 + 1) + act2*max1
260 & + act3*max1*max2
261 & + act4*max1*max2*max3
262 #endif /* ALLOW_AUTODIFF_TAMC */
263
264 C-- Set up work arrays with valid (i.e. not NaN) values
265 C These inital values do not alter the numerical results. They
266 C just ensure that all memory references are to valid floating
267 C point numbers. This prevents spurious hardware signals due to
268 C uninitialised but inert locations.
269
270 DO k=1,Nr
271 DO j=1-OLy,sNy+OLy
272 DO i=1-OLx,sNx+OLx
273 KappaRU(i,j,k) = 0. _d 0
274 KappaRV(i,j,k) = 0. _d 0
275 #ifdef ALLOW_AUTODIFF_TAMC
276 cph(
277 c-- need some re-initialisation here to break dependencies
278 cph)
279 gU(i,j,k,bi,bj) = 0. _d 0
280 gV(i,j,k,bi,bj) = 0. _d 0
281 #endif
282 ENDDO
283 ENDDO
284 ENDDO
285 DO j=1-OLy,sNy+OLy
286 DO i=1-OLx,sNx+OLx
287 fVerU (i,j,1) = 0. _d 0
288 fVerU (i,j,2) = 0. _d 0
289 fVerV (i,j,1) = 0. _d 0
290 fVerV (i,j,2) = 0. _d 0
291 phiHydF (i,j) = 0. _d 0
292 phiHydC (i,j) = 0. _d 0
293 dPhiHydX(i,j) = 0. _d 0
294 dPhiHydY(i,j) = 0. _d 0
295 phiSurfX(i,j) = 0. _d 0
296 phiSurfY(i,j) = 0. _d 0
297 guDissip(i,j) = 0. _d 0
298 gvDissip(i,j) = 0. _d 0
299 ENDDO
300 ENDDO
301
302 C-- Start computation of dynamics
303 iMin = 0
304 iMax = sNx+1
305 jMin = 0
306 jMax = sNy+1
307
308 #ifdef ALLOW_AUTODIFF_TAMC
309 CADJ STORE wvel (:,:,:,bi,bj) =
310 CADJ & comlev1_bibj, key = idynkey, byte = isbyte
311 #endif /* ALLOW_AUTODIFF_TAMC */
312
313 C-- Explicit part of the Surface Potentiel Gradient (add in TIMESTEP)
314 C (note: this loop will be replaced by CALL CALC_GRAD_ETA)
315 IF (implicSurfPress.NE.1.) THEN
316 CALL CALC_GRAD_PHI_SURF(
317 I bi,bj,iMin,iMax,jMin,jMax,
318 I etaN,
319 O phiSurfX,phiSurfY,
320 I myThid )
321 ENDIF
322
323 #ifdef ALLOW_AUTODIFF_TAMC
324 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
325 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=idynkey, byte=isbyte
326 #ifdef ALLOW_KPP
327 CADJ STORE KPPviscAz (:,:,:,bi,bj)
328 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
329 #endif /* ALLOW_KPP */
330 #endif /* ALLOW_AUTODIFF_TAMC */
331
332 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
333 C-- Calculate the total vertical diffusivity
334 DO k=1,Nr
335 CALL CALC_VISCOSITY(
336 I bi,bj,iMin,iMax,jMin,jMax,k,
337 O KappaRU,KappaRV,
338 I myThid)
339 ENDDO
340 #endif
341
342 #ifdef ALLOW_AUTODIFF_TAMC
343 CADJ STORE KappaRU(:,:,:)
344 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
345 CADJ STORE KappaRV(:,:,:)
346 CADJ & = comlev1_bibj, key=idynkey, byte=isbyte
347 #endif /* ALLOW_AUTODIFF_TAMC */
348
349 C-- Start of dynamics loop
350 DO k=1,Nr
351
352 C-- km1 Points to level above k (=k-1)
353 C-- kup Cycles through 1,2 to point to layer above
354 C-- kDown Cycles through 2,1 to point to current layer
355
356 km1 = MAX(1,k-1)
357 kp1 = MIN(k+1,Nr)
358 kup = 1+MOD(k+1,2)
359 kDown= 1+MOD(k,2)
360
361 #ifdef ALLOW_AUTODIFF_TAMC
362 kkey = (idynkey-1)*Nr + k
363 c
364 CADJ STORE totphihyd (:,:,k,bi,bj)
365 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
366 CADJ STORE theta (:,:,k,bi,bj)
367 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
368 CADJ STORE salt (:,:,k,bi,bj)
369 CADJ & = comlev1_bibj_k, key=kkey, byte=isbyte
370 #endif /* ALLOW_AUTODIFF_TAMC */
371
372 C-- Integrate hydrostatic balance for phiHyd with BC of
373 C phiHyd(z=0)=0
374 CALL CALC_PHI_HYD(
375 I bi,bj,iMin,iMax,jMin,jMax,k,
376 I theta, salt,
377 U phiHydF,
378 O phiHydC, dPhiHydX, dPhiHydY,
379 I myTime, myIter, myThid )
380
381 C-- Calculate accelerations in the momentum equations (gU, gV, ...)
382 C and step forward storing the result in gU, gV, etc...
383 IF ( momStepping ) THEN
384 #ifdef ALLOW_MOM_FLUXFORM
385 IF (.NOT. vectorInvariantMomentum) CALL MOM_FLUXFORM(
386 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
387 I KappaRU, KappaRV,
388 U fVerU, fVerV,
389 O guDissip, gvDissip,
390 I myTime, myIter, myThid)
391 #endif
392 #ifdef ALLOW_MOM_VECINV
393 IF (vectorInvariantMomentum) CALL MOM_VECINV(
394 I bi,bj,iMin,iMax,jMin,jMax,k,kup,kDown,
395 I KappaRU, KappaRV,
396 U fVerU, fVerV,
397 O guDissip, gvDissip,
398 I myTime, myIter, myThid)
399 #endif
400 CALL TIMESTEP(
401 I bi,bj,iMin,iMax,jMin,jMax,k,
402 I dPhiHydX,dPhiHydY, phiSurfX, phiSurfY,
403 I guDissip, gvDissip,
404 I myTime, myIter, myThid)
405
406 #ifdef ALLOW_OBCS
407 C-- Apply open boundary conditions
408 IF (useOBCS) THEN
409 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
410 ENDIF
411 #endif /* ALLOW_OBCS */
412
413 ENDIF
414
415
416 C-- end of dynamics k loop (1:Nr)
417 ENDDO
418
419 C-- Implicit Vertical advection & viscosity
420 #ifdef INCLUDE_IMPLVERTADV_CODE
421 IF ( momImplVertAdv ) THEN
422 CALL MOM_U_IMPLICIT_R( kappaRU,
423 I bi, bj, myTime, myIter, myThid )
424 CALL MOM_V_IMPLICIT_R( kappaRV,
425 I bi, bj, myTime, myIter, myThid )
426 ELSEIF ( implicitViscosity ) THEN
427 #else /* INCLUDE_IMPLVERTADV_CODE */
428 IF ( implicitViscosity ) THEN
429 #endif /* INCLUDE_IMPLVERTADV_CODE */
430 #ifdef ALLOW_AUTODIFF_TAMC
431 CADJ STORE KappaRU(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
432 CADJ STORE gU(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
433 #endif /* ALLOW_AUTODIFF_TAMC */
434 CALL IMPLDIFF(
435 I bi, bj, iMin, iMax, jMin, jMax,
436 I -1, KappaRU,recip_HFacW,
437 U gU,
438 I myThid )
439 #ifdef ALLOW_AUTODIFF_TAMC
440 CADJ STORE KappaRV(:,:,:) = comlev1_bibj , key=idynkey, byte=isbyte
441 CADJ STORE gV(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
442 #endif /* ALLOW_AUTODIFF_TAMC */
443 CALL IMPLDIFF(
444 I bi, bj, iMin, iMax, jMin, jMax,
445 I -2, KappaRV,recip_HFacS,
446 U gV,
447 I myThid )
448 ENDIF
449
450 #ifdef ALLOW_OBCS
451 C-- Apply open boundary conditions
452 IF ( useOBCS .AND.(implicitViscosity.OR.momImplVertAdv) ) THEN
453 DO K=1,Nr
454 CALL OBCS_APPLY_UV( bi, bj, k, gU, gV, myThid )
455 ENDDO
456 ENDIF
457 #endif /* ALLOW_OBCS */
458
459 #ifdef ALLOW_CD_CODE
460 IF (implicitViscosity.AND.useCDscheme) THEN
461 #ifdef ALLOW_AUTODIFF_TAMC
462 CADJ STORE vVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
463 #endif /* ALLOW_AUTODIFF_TAMC */
464 CALL IMPLDIFF(
465 I bi, bj, iMin, iMax, jMin, jMax,
466 I 0, KappaRU,recip_HFacW,
467 U vVelD,
468 I myThid )
469 #ifdef ALLOW_AUTODIFF_TAMC
470 CADJ STORE uVelD(:,:,:,bi,bj) = comlev1_bibj , key=idynkey, byte=isbyte
471 #endif /* ALLOW_AUTODIFF_TAMC */
472 CALL IMPLDIFF(
473 I bi, bj, iMin, iMax, jMin, jMax,
474 I 0, KappaRV,recip_HFacS,
475 U uVelD,
476 I myThid )
477 ENDIF
478 #endif /* ALLOW_CD_CODE */
479 C-- End implicit Vertical advection & viscosity
480
481 ENDDO
482 ENDDO
483
484 #ifdef ALLOW_OBCS
485 IF (useOBCS) THEN
486 CALL OBCS_PRESCRIBE_EXCHANGES(myThid)
487 ENDIF
488 #endif
489
490 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
491
492 #ifdef ALLOW_NONHYDROSTATIC
493 C-- Step forward W field in N-H algorithm
494 IF ( momStepping .AND. nonHydrostatic ) THEN
495 #ifdef ALLOW_DEBUG
496 IF ( debugLevel .GE. debLevB )
497 & CALL DEBUG_CALL('CALC_GW', myThid )
498 #endif
499 CALL TIMER_START('CALC_GW [DYNAMICS]',myThid)
500 CALL CALC_GW( myTime, myIter, myThid )
501 CALL TIMER_STOP ('CALC_GW [DYNAMICS]',myThid)
502 ENDIF
503 #endif
504
505 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
506
507 Cml(
508 C In order to compare the variance of phiHydLow of a p/z-coordinate
509 C run with etaH of a z/p-coordinate run the drift of phiHydLow
510 C has to be removed by something like the following subroutine:
511 C CALL REMOVE_MEAN_RL( 1, phiHydLow, maskH, maskH, rA, drF,
512 C & 'phiHydLow', myThid )
513 Cml)
514
515 #ifdef ALLOW_DIAGNOSTICS
516 IF ( usediagnostics ) THEN
517
518 CALL DIAGNOSTICS_FILL(totPhihyd,'PHIHYD ',0,Nr,0,1,1,myThid)
519 CALL DIAGNOSTICS_FILL(phiHydLow,'PHIBOT ',0, 1,0,1,1,myThid)
520
521 tmpFac = 1. _d 0
522 CALL DIAGNOSTICS_SCALE_FILL(totPhihyd,tmpFac,2,
523 & 'PHIHYDSQ',0,Nr,0,1,1,myThid)
524
525 CALL DIAGNOSTICS_SCALE_FILL(phiHydLow,tmpFac,2,
526 & 'PHIBOTSQ',0, 1,0,1,1,myThid)
527
528 ENDIF
529 #endif /* ALLOW_DIAGNOSTICS */
530
531 #ifdef ALLOW_DEBUG
532 If ( debugLevel .GE. debLevB ) THEN
533 CALL DEBUG_STATS_RL(1,EtaN,'EtaN (DYNAMICS)',myThid)
534 CALL DEBUG_STATS_RL(Nr,uVel,'Uvel (DYNAMICS)',myThid)
535 CALL DEBUG_STATS_RL(Nr,vVel,'Vvel (DYNAMICS)',myThid)
536 CALL DEBUG_STATS_RL(Nr,wVel,'Wvel (DYNAMICS)',myThid)
537 CALL DEBUG_STATS_RL(Nr,theta,'Theta (DYNAMICS)',myThid)
538 CALL DEBUG_STATS_RL(Nr,salt,'Salt (DYNAMICS)',myThid)
539 CALL DEBUG_STATS_RL(Nr,gU,'Gu (DYNAMICS)',myThid)
540 CALL DEBUG_STATS_RL(Nr,gV,'Gv (DYNAMICS)',myThid)
541 CALL DEBUG_STATS_RL(Nr,gT,'Gt (DYNAMICS)',myThid)
542 CALL DEBUG_STATS_RL(Nr,gS,'Gs (DYNAMICS)',myThid)
543 #ifndef ALLOW_ADAMSBASHFORTH_3
544 CALL DEBUG_STATS_RL(Nr,guNm1,'GuNm1 (DYNAMICS)',myThid)
545 CALL DEBUG_STATS_RL(Nr,gvNm1,'GvNm1 (DYNAMICS)',myThid)
546 CALL DEBUG_STATS_RL(Nr,gtNm1,'GtNm1 (DYNAMICS)',myThid)
547 CALL DEBUG_STATS_RL(Nr,gsNm1,'GsNm1 (DYNAMICS)',myThid)
548 #endif
549 ENDIF
550 #endif
551
552 #ifdef ALLOW_DEBUG
553 IF ( debugLevel .GE. debLevB )
554 & CALL DEBUG_LEAVE( 'DYNAMICS', myThid )
555 #endif
556
557 RETURN
558 END

  ViewVC Help
Powered by ViewVC 1.1.22