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

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

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


Revision 1.27 - (show annotations) (download)
Wed Mar 8 06:36:39 2006 UTC (18 years, 3 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58c_post, checkpoint58b_post
Changes since 1.26: +14 -1 lines
o Another overhaul of store dirs. for NLFS to eliminate
  "hidden" recomputations.
o TBD: "hidden" mom_vecinv recomp. in dynamics

1 C $Header: /u/gcmpack/MITgcm/model/src/do_oceanic_phys.F,v 1.26 2006/03/06 18:25:49 heimbach Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 #ifdef ALLOW_AUTODIFF_TAMC
8 # ifdef ALLOW_GMREDI
9 # include "GMREDI_OPTIONS.h"
10 # endif
11 # ifdef ALLOW_KPP
12 # include "KPP_OPTIONS.h"
13 # endif
14 #endif /* ALLOW_AUTODIFF_TAMC */
15
16 CBOP
17 C !ROUTINE: DO_OCEANIC_PHYS
18 C !INTERFACE:
19 SUBROUTINE DO_OCEANIC_PHYS(myTime, myIter, myThid)
20 C !DESCRIPTION: \bv
21 C *==========================================================*
22 C | SUBROUTINE DO_OCEANIC_PHYS
23 C | o Controlling routine for oceanic physics and
24 C | parameterization
25 C *==========================================================*
26 C | o originally, part of S/R thermodynamics
27 C *==========================================================*
28 C \ev
29
30 C !USES:
31 IMPLICIT NONE
32 C == Global variables ===
33 #include "SIZE.h"
34 #include "EEPARAMS.h"
35 #include "PARAMS.h"
36 #include "DYNVARS.h"
37 #include "GRID.h"
38 #ifdef ALLOW_TIMEAVE
39 #include "TIMEAVE_STATV.h"
40 #endif
41 #if defined (ALLOW_BALANCE_FLUXES) && !(defined ALLOW_AUTODIFF_TAMC)
42 #include "FFIELDS.h"
43 #endif
44
45 #ifdef ALLOW_AUTODIFF_TAMC
46 # include "tamc.h"
47 # include "tamc_keys.h"
48 # include "FFIELDS.h"
49 # include "EOS.h"
50 # ifdef ALLOW_KPP
51 # include "KPP.h"
52 # endif
53 # ifdef ALLOW_GMREDI
54 # include "GMREDI.h"
55 # endif
56 # ifdef ALLOW_EBM
57 # include "EBM.h"
58 # endif
59 # ifdef EXACT_CONSERV
60 # include "SURFACE.h"
61 # endif
62 #endif /* ALLOW_AUTODIFF_TAMC */
63
64 C !INPUT/OUTPUT PARAMETERS:
65 C == Routine arguments ==
66 C myTime :: Current time in simulation
67 C myIter :: Current iteration number in simulation
68 C myThid :: Thread number for this instance of the routine.
69 _RL myTime
70 INTEGER myIter
71 INTEGER myThid
72
73 C !LOCAL VARIABLES:
74 C == Local variables
75 C rhoK, rhoKM1 :: Density at current level, and level above
76 C iMin, iMax :: Ranges and sub-block indices on which calculations
77 C jMin, jMax are applied.
78 C bi, bj :: tile indices
79 C i,j,k :: loop indices
80 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81 _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
82 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
83 _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
84 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
85 INTEGER iMin, iMax
86 INTEGER jMin, jMax
87 INTEGER bi, bj
88 INTEGER i, j, k
89 INTEGER doDiagsRho
90 #ifdef ALLOW_DIAGNOSTICS
91 LOGICAL DIAGNOSTICS_IS_ON
92 EXTERNAL DIAGNOSTICS_IS_ON
93 #endif /* ALLOW_DIAGNOSTICS */
94
95 CEOP
96
97 #ifdef ALLOW_AUTODIFF_TAMC
98 C-- dummy statement to end declaration part
99 itdkey = 1
100 #endif /* ALLOW_AUTODIFF_TAMC */
101
102 #ifdef ALLOW_DEBUG
103 IF ( debugLevel .GE. debLevB )
104 & CALL DEBUG_ENTER('DO_OCEANIC_PHYS',myThid)
105 #endif
106
107 doDiagsRho = 0
108 #ifdef ALLOW_DIAGNOSTICS
109 IF ( useDiagnostics .AND. fluidIsWater ) THEN
110 IF ( DIAGNOSTICS_IS_ON('DRHODR ',myThid) ) doDiagsRho = 1
111 IF ( DIAGNOSTICS_IS_ON('RHOANOSQ',myThid) .OR.
112 & DIAGNOSTICS_IS_ON('URHOMASS',myThid) .OR.
113 & DIAGNOSTICS_IS_ON('VRHOMASS',myThid) .OR.
114 & DIAGNOSTICS_IS_ON('WRHOMASS',myThid) .OR.
115 & DIAGNOSTICS_IS_ON('WRHOMASS',myThid) ) doDiagsRho = 2
116 ENDIF
117 #endif /* ALLOW_DIAGNOSTICS */
118
119 #ifdef ALLOW_THSICE
120 IF ( useThSIce .AND. fluidIsWater ) THEN
121 #ifdef ALLOW_DEBUG
122 IF ( debugLevel .GE. debLevB )
123 & CALL DEBUG_CALL('THSICE_MAIN',myThid)
124 #endif
125 C-- Step forward Therm.Sea-Ice variables
126 C and modify forcing terms including effects from ice
127 CALL TIMER_START('THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
128 CALL THSICE_MAIN( myTime, myIter, myThid )
129 CALL TIMER_STOP( 'THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
130 ENDIF
131 #endif /* ALLOW_THSICE */
132
133 #ifdef ALLOW_SHELFICE
134 IF ( useShelfIce .AND. fluidIsWater ) THEN
135 #ifdef ALLOW_DEBUG
136 IF ( debugLevel .GE. debLevB )
137 & CALL DEBUG_CALL('SHELFICE_THERMODYNAMICS',myThid)
138 #endif
139 C compute temperature and (virtual) salt flux at the
140 C shelf-ice ocean interface
141 CALL TIMER_START('SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',
142 & myThid)
143 CALL SHELFICE_THERMODYNAMICS( myTime, myIter, myThid )
144 CALL TIMER_STOP( 'SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',
145 & myThid)
146 ENDIF
147 #endif /* ALLOW_SHELFICE */
148
149 C-- Freeze water at the surface
150 #ifdef ALLOW_AUTODIFF_TAMC
151 CADJ STORE theta = comlev1, key = ikey_dynamics
152 #endif
153 IF ( allowFreezing
154 & .AND. .NOT. useSEAICE
155 & .AND. .NOT. useThSIce ) THEN
156 CALL FREEZE_SURFACE( myTime, myIter, myThid )
157 ENDIF
158
159 #ifdef COMPONENT_MODULE
160 # ifndef ALLOW_AIM
161 C-- Apply imported data (from coupled interface) to forcing fields
162 C jmc: do not know precisely where to put this call (bf or af thSIce ?)
163 IF ( useCoupler ) THEN
164 CALL OCN_APPLY_IMPORT( .TRUE., myTime, myIter, myThid )
165 ENDIF
166 # endif
167 #endif /* COMPONENT_MODULE */
168
169 #ifdef ALLOW_BALANCE_FLUXES
170 C balance fluxes
171 IF ( balanceEmPmR )
172 & CALL REMOVE_MEAN_RS( 1, EmPmR, maskH, maskH, rA, drF,
173 & 'EmPmR', myTime, myThid )
174 IF ( balanceQnet )
175 & CALL REMOVE_MEAN_RS( 1, Qnet, maskH, maskH, rA, drF,
176 & 'Qnet ', myTime, myThid )
177 #endif /* ALLOW_BALANCE_FLUXES */
178
179 #ifdef ALLOW_AUTODIFF_TAMC
180 C-- HPF directive to help TAMC
181 CHPF$ INDEPENDENT
182 #endif /* ALLOW_AUTODIFF_TAMC */
183 DO bj=myByLo(myThid),myByHi(myThid)
184 #ifdef ALLOW_AUTODIFF_TAMC
185 C-- HPF directive to help TAMC
186 CHPF$ INDEPENDENT
187 #endif /* ALLOW_AUTODIFF_TAMC */
188 DO bi=myBxLo(myThid),myBxHi(myThid)
189
190 #ifdef ALLOW_AUTODIFF_TAMC
191 act1 = bi - myBxLo(myThid)
192 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
193 act2 = bj - myByLo(myThid)
194 max2 = myByHi(myThid) - myByLo(myThid) + 1
195 act3 = myThid - 1
196 max3 = nTx*nTy
197 act4 = ikey_dynamics - 1
198 itdkey = (act1 + 1) + act2*max1
199 & + act3*max1*max2
200 & + act4*max1*max2*max3
201 #endif /* ALLOW_AUTODIFF_TAMC */
202
203 C-- Set up work arrays with valid (i.e. not NaN) values
204 C These inital values do not alter the numerical results. They
205 C just ensure that all memory references are to valid floating
206 C point numbers. This prevents spurious hardware signals due to
207 C uninitialised but inert locations.
208
209 DO j=1-OLy,sNy+OLy
210 DO i=1-OLx,sNx+OLx
211 rhok (i,j) = 0. _d 0
212 rhoKM1 (i,j) = 0. _d 0
213 ENDDO
214 ENDDO
215
216 DO k=1,Nr
217 DO j=1-OLy,sNy+OLy
218 DO i=1-OLx,sNx+OLx
219 C This is currently also used by IVDC and Diagnostics
220 sigmaX(i,j,k) = 0. _d 0
221 sigmaY(i,j,k) = 0. _d 0
222 sigmaR(i,j,k) = 0. _d 0
223 #ifdef ALLOW_AUTODIFF_TAMC
224 cph all the following init. are necessary for TAF
225 cph although some of these are re-initialised later.
226 IVDConvCount(i,j,k,bi,bj) = 0.
227 # ifdef ALLOW_GMREDI
228 Kwx(i,j,k,bi,bj) = 0. _d 0
229 Kwy(i,j,k,bi,bj) = 0. _d 0
230 Kwz(i,j,k,bi,bj) = 0. _d 0
231 # ifdef GM_NON_UNITY_DIAGONAL
232 Kux(i,j,k,bi,bj) = 0. _d 0
233 Kvy(i,j,k,bi,bj) = 0. _d 0
234 # endif
235 # ifdef GM_EXTRA_DIAGONAL
236 Kuz(i,j,k,bi,bj) = 0. _d 0
237 Kvz(i,j,k,bi,bj) = 0. _d 0
238 # endif
239 # ifdef GM_BOLUS_ADVEC
240 GM_PsiX(i,j,k,bi,bj) = 0. _d 0
241 GM_PsiY(i,j,k,bi,bj) = 0. _d 0
242 # endif
243 # ifdef GM_VISBECK_VARIABLE_K
244 VisbeckK(i,j,bi,bj) = 0. _d 0
245 # endif
246 # endif /* ALLOW_GMREDI */
247 #endif /* ALLOW_AUTODIFF_TAMC */
248 ENDDO
249 ENDDO
250 ENDDO
251
252 iMin = 1-OLx
253 iMax = sNx+OLx
254 jMin = 1-OLy
255 jMax = sNy+OLy
256
257 #ifdef ALLOW_AUTODIFF_TAMC
258 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
259 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
260 CADJ STORE totphihyd(:,:,:,bi,bj)
261 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
262 # ifdef ALLOW_KPP
263 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
264 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
265 # endif
266 #endif /* ALLOW_AUTODIFF_TAMC */
267
268 #ifdef ALLOW_DEBUG
269 IF ( debugLevel .GE. debLevB )
270 & CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
271 #endif
272
273 C-- Start of diagnostic loop
274 DO k=Nr,1,-1
275
276 #ifdef ALLOW_AUTODIFF_TAMC
277 C? Patrick, is this formula correct now that we change the loop range?
278 C? Do we still need this?
279 cph kkey formula corrected.
280 cph Needed for rhok, rhokm1, in the case useGMREDI.
281 kkey = (itdkey-1)*Nr + k
282 #endif /* ALLOW_AUTODIFF_TAMC */
283
284 C-- Calculate gradients of potential density for isoneutral
285 C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
286 c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
287 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.)
288 & .OR. doDiagsRho.GE.1 ) THEN
289 #ifdef ALLOW_DEBUG
290 IF ( debugLevel .GE. debLevB )
291 & CALL DEBUG_CALL('FIND_RHO',myThid)
292 #endif
293 #ifdef ALLOW_AUTODIFF_TAMC
294 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
295 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
296 #endif /* ALLOW_AUTODIFF_TAMC */
297 CALL FIND_RHO(
298 I bi, bj, iMin, iMax, jMin, jMax, k, k,
299 I theta, salt,
300 O rhoK,
301 I myThid )
302
303 IF (k.GT.1) THEN
304 #ifdef ALLOW_AUTODIFF_TAMC
305 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
306 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
307 #endif /* ALLOW_AUTODIFF_TAMC */
308 CALL FIND_RHO(
309 I bi, bj, iMin, iMax, jMin, jMax, k-1, k,
310 I theta, salt,
311 O rhoKm1,
312 I myThid )
313 ENDIF
314 #ifdef ALLOW_DEBUG
315 IF ( debugLevel .GE. debLevB )
316 & CALL DEBUG_CALL('GRAD_SIGMA',myThid)
317 #endif
318 CALL GRAD_SIGMA(
319 I bi, bj, iMin, iMax, jMin, jMax, k,
320 I rhoK, rhoKm1, rhoK,
321 O sigmaX, sigmaY, sigmaR,
322 I myThid )
323 ENDIF
324
325 #ifdef ALLOW_AUTODIFF_TAMC
326 ctest# ifndef GM_EXCLUDE_CLIPPING
327 CADJ STORE rhok (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
328 ctest# endif
329 CADJ STORE rhokm1 (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
330 #endif /* ALLOW_AUTODIFF_TAMC */
331 C-- Implicit Vertical Diffusion for Convection
332 c ==> should use sigmaR !!!
333 IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
334 #ifdef ALLOW_DEBUG
335 IF ( debugLevel .GE. debLevB )
336 & CALL DEBUG_CALL('CALC_IVDC',myThid)
337 #endif
338 CALL CALC_IVDC(
339 I bi, bj, iMin, iMax, jMin, jMax, k,
340 I rhoKm1, rhoK,
341 I myTime, myIter, myThid)
342 ENDIF
343
344 #ifdef ALLOW_DIAGNOSTICS
345 IF ( doDiagsRho.GE.2 ) THEN
346 CALL DIAGS_RHO( k, bi, bj,
347 I rhoK, rhoKm1,
348 I myTime, myIter, myThid)
349 ENDIF
350 #endif
351
352 C-- end of diagnostic k loop (Nr:1)
353 ENDDO
354
355 #ifdef ALLOW_DIAGNOSTICS
356 c IF ( useDiagnostics .AND.
357 c & (useGMRedi .OR. ivdc_kappa.NE.0.) ) THEN
358 IF ( doDiagsRho.GE.1 ) THEN
359 CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR ', 0, Nr,
360 & 2, bi, bj, myThid)
361 ENDIF
362 #endif
363
364 #ifdef ALLOW_OBCS
365 C-- Calculate future values on open boundaries
366 IF (useOBCS) THEN
367 #ifdef ALLOW_DEBUG
368 IF ( debugLevel .GE. debLevB )
369 & CALL DEBUG_CALL('OBCS_CALC',myThid)
370 #endif
371 CALL OBCS_CALC( bi, bj, myTime+deltaTclock, myIter+1,
372 I uVel, vVel, wVel, theta, salt,
373 I myThid )
374 ENDIF
375 #endif /* ALLOW_OBCS */
376
377 #ifndef ALLOW_AUTODIFF_TAMC
378 IF ( fluidIsWater ) THEN
379 #endif
380 C-- Determines forcing terms based on external fields
381 C relaxation terms, etc.
382 #ifdef ALLOW_DEBUG
383 IF ( debugLevel .GE. debLevB )
384 & CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
385 #endif
386 #ifdef ALLOW_AUTODIFF_TAMC
387 CADJ STORE EmPmR(:,:,bi,bj)
388 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
389 # ifdef EXACT_CONSERV
390 CADJ STORE PmEpR(:,:,bi,bj)
391 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
392 # endif
393 # ifdef NONLIN_FRSURF
394 CADJ STORE hFac_surfC(:,:,bi,bj)
395 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
396 CADJ STORE recip_hFacC(:,:,:,bi,bj)
397 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
398 # endif
399 #endif
400 CALL EXTERNAL_FORCING_SURF(
401 I bi, bj, iMin, iMax, jMin, jMax,
402 I myTime, myIter, myThid )
403 #ifndef ALLOW_AUTODIFF_TAMC
404 ENDIF
405 #endif
406 #ifdef ALLOW_AUTODIFF_TAMC
407 # ifdef EXACT_CONSERV
408 cph-test
409 cphCADJ STORE PmEpR(:,:,bi,bj)
410 cphCADJ & = comlev1_bibj, key=itdkey, byte=isbyte
411 # endif
412 #endif
413
414 #ifdef ALLOW_AUTODIFF_TAMC
415 cph needed for KPP
416 CADJ STORE surfaceForcingU(:,:,bi,bj)
417 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
418 CADJ STORE surfaceForcingV(:,:,bi,bj)
419 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
420 CADJ STORE surfaceForcingS(:,:,bi,bj)
421 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
422 CADJ STORE surfaceForcingT(:,:,bi,bj)
423 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
424 CADJ STORE surfaceForcingTice(:,:,bi,bj)
425 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
426 #endif /* ALLOW_AUTODIFF_TAMC */
427
428 #ifdef ALLOW_GMREDI
429
430 #ifdef ALLOW_AUTODIFF_TAMC
431 # ifndef GM_EXCLUDE_CLIPPING
432 cph storing here is needed only for one GMREDI_OPTIONS:
433 cph define GM_BOLUS_ADVEC
434 cph keep it although TAF says you dont need to.
435 cph but I've avoided the #ifdef for now, in case more things change
436 CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
437 CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
438 CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
439 # endif
440 #endif /* ALLOW_AUTODIFF_TAMC */
441
442 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
443 IF (useGMRedi) THEN
444 #ifdef ALLOW_DEBUG
445 IF ( debugLevel .GE. debLevB )
446 & CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
447 #endif
448 CALL GMREDI_CALC_TENSOR(
449 I bi, bj, iMin, iMax, jMin, jMax,
450 I sigmaX, sigmaY, sigmaR,
451 I myThid )
452 #ifdef ALLOW_AUTODIFF_TAMC
453 ELSE
454 CALL GMREDI_CALC_TENSOR_DUMMY(
455 I bi, bj, iMin, iMax, jMin, jMax,
456 I sigmaX, sigmaY, sigmaR,
457 I myThid )
458 #endif /* ALLOW_AUTODIFF_TAMC */
459 ENDIF
460
461 #endif /* ALLOW_GMREDI */
462
463 #ifdef ALLOW_KPP
464 C-- Compute KPP mixing coefficients
465 IF (useKPP) THEN
466 #ifdef ALLOW_DEBUG
467 IF ( debugLevel .GE. debLevB )
468 & CALL DEBUG_CALL('KPP_CALC',myThid)
469 #endif
470 CALL KPP_CALC(
471 I bi, bj, myTime, myThid )
472 #ifdef ALLOW_AUTODIFF_TAMC
473 ELSE
474 CALL KPP_CALC_DUMMY(
475 I bi, bj, myTime, myThid )
476 #endif /* ALLOW_AUTODIFF_TAMC */
477 ENDIF
478
479 #endif /* ALLOW_KPP */
480
481 #ifdef ALLOW_PP81
482 C-- Compute PP81 mixing coefficients
483 IF (usePP81) THEN
484 #ifdef ALLOW_DEBUG
485 IF ( debugLevel .GE. debLevB )
486 & CALL DEBUG_CALL('PP81_CALC',myThid)
487 #endif
488 CALL PP81_CALC(
489 I bi, bj, myTime, myThid )
490 ENDIF
491 #endif /* ALLOW_PP81 */
492
493 #ifdef ALLOW_MY82
494 C-- Compute MY82 mixing coefficients
495 IF (useMY82) THEN
496 #ifdef ALLOW_DEBUG
497 IF ( debugLevel .GE. debLevB )
498 & CALL DEBUG_CALL('MY82_CALC',myThid)
499 #endif
500 CALL MY82_CALC(
501 I bi, bj, myTime, myThid )
502 ENDIF
503 #endif /* ALLOW_MY82 */
504
505 #ifdef ALLOW_GGL90
506 C-- Compute GGL90 mixing coefficients
507 IF (useGGL90) THEN
508 #ifdef ALLOW_DEBUG
509 IF ( debugLevel .GE. debLevB )
510 & CALL DEBUG_CALL('GGL90_CALC',myThid)
511 #endif
512 CALL GGL90_CALC(
513 I bi, bj, myTime, myThid )
514 ENDIF
515 #endif /* ALLOW_GGL90 */
516
517 #ifdef ALLOW_TIMEAVE
518 IF ( taveFreq.GT. 0. _d 0 .AND. fluidIsWater ) THEN
519 CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid)
520 ENDIF
521 IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
522 CALL TIMEAVE_CUMULATE(ConvectCountTave, IVDConvCount,
523 I Nr, deltaTclock, bi, bj, myThid)
524 ENDIF
525 #endif /* ALLOW_TIMEAVE */
526
527 C-- end bi,bj loops.
528 ENDDO
529 ENDDO
530
531 #ifdef ALLOW_DIAGNOSTICS
532 IF ( fluidIsWater .AND. useDiagnostics ) THEN
533 CALL DIAGS_OCEANIC_SURF_FLUX( myTime, myIter, myThid )
534 ENDIF
535 IF ( ivdc_kappa.NE.0 .AND. useDiagnostics ) THEN
536 CALL DIAGNOSTICS_FILL( IVDConvCount,'CONVADJ ',
537 & 0, Nr, 0, 1, 1, myThid )
538 ENDIF
539 #endif
540
541 #ifdef ALLOW_DEBUG
542 IF ( debugLevel .GE. debLevB )
543 & CALL DEBUG_LEAVE('DO_OCEANIC_PHYS',myThid)
544 #endif
545
546 RETURN
547 END

  ViewVC Help
Powered by ViewVC 1.1.22