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

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

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


Revision 1.22 - (hide annotations) (download)
Fri Feb 10 07:56:20 2006 UTC (18 years, 3 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint58a_post
Changes since 1.21: +14 -1 lines
o add code to balance EmPmP and Qnet at the end of do_ocean_physics. Useful
  if bulk formulae are used in long integration (especially EmPmR). Turn
  on with balanceEmPmR = .true. or balanceQnet = .true. in data, PARM01
  if balancePrintMean, the imbalance that is substracted is printed to
  STDOUT.

1 mlosch 1.22 C $Header: /u/gcmpack/MITgcm/model/src/do_oceanic_phys.F,v 1.21 2006/02/07 11:47:49 mlosch Exp $
2 jmc 1.1 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 jmc 1.20 #ifdef ALLOW_TIMEAVE
39     #include "TIMEAVE_STATV.h"
40     #endif
41 mlosch 1.22 #if defined (ALLOW_BALANCE_FLUXES) && !(defined ALLOW_AUTODIFF_TAMC)
42     #include "FFIELDS.h"
43     #endif
44 jmc 1.1
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 heimbach 1.10 # ifdef EXACT_CONSERV
60     # include "SURFACE.h"
61     # endif
62 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
63    
64     C !INPUT/OUTPUT PARAMETERS:
65     C == Routine arguments ==
66 jmc 1.18 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 jmc 1.1 _RL myTime
70     INTEGER myIter
71     INTEGER myThid
72    
73     C !LOCAL VARIABLES:
74     C == Local variables
75 jmc 1.18 C rhoK, rhoKM1 :: Density at current level, and level above
76     C iMin, iMax :: Ranges and sub-block indices on which calculations
77 jmc 1.1 C jMin, jMax are applied.
78 jmc 1.18 C bi, bj :: tile indices
79     C i,j,k :: loop indices
80 jmc 1.1 _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 jmc 1.18 INTEGER i, j, k
89 jmc 1.17 INTEGER doDiagsRho
90     #ifdef ALLOW_DIAGNOSTICS
91     LOGICAL DIAGNOSTICS_IS_ON
92     EXTERNAL DIAGNOSTICS_IS_ON
93     #endif /* ALLOW_DIAGNOSTICS */
94 jmc 1.1
95     CEOP
96    
97 heimbach 1.12 #ifdef ALLOW_AUTODIFF_TAMC
98     C-- dummy statement to end declaration part
99     itdkey = 1
100     #endif /* ALLOW_AUTODIFF_TAMC */
101    
102 jmc 1.1 #ifdef ALLOW_DEBUG
103 jmc 1.5 IF ( debugLevel .GE. debLevB )
104 jmc 1.1 & CALL DEBUG_ENTER('DO_OCEANIC_PHYS',myThid)
105     #endif
106    
107 jmc 1.17 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 jmc 1.5 #ifdef ALLOW_THSICE
120 jmc 1.14 IF ( useThSIce .AND. fluidIsWater ) THEN
121 jmc 1.5 #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 mlosch 1.21 #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 jmc 1.5 C-- Freeze water at the surface
150     #ifdef ALLOW_AUTODIFF_TAMC
151     CADJ STORE theta = comlev1, key = ikey_dynamics
152     #endif
153 heimbach 1.12 IF ( allowFreezing
154     & .AND. .NOT. useSEAICE
155 jmc 1.5 & .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 edhill 1.7 C jmc: do not know precisely where to put this call (bf or af thSIce ?)
163 jmc 1.5 IF ( useCoupler ) THEN
164     CALL OCN_APPLY_IMPORT( .TRUE., myTime, myIter, myThid )
165     ENDIF
166     # endif
167     #endif /* COMPONENT_MODULE */
168    
169 jmc 1.1 #ifdef ALLOW_AUTODIFF_TAMC
170     C-- HPF directive to help TAMC
171     CHPF$ INDEPENDENT
172     #endif /* ALLOW_AUTODIFF_TAMC */
173     DO bj=myByLo(myThid),myByHi(myThid)
174     #ifdef ALLOW_AUTODIFF_TAMC
175 heimbach 1.15 C-- HPF directive to help TAMC
176     CHPF$ INDEPENDENT
177 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
178     DO bi=myBxLo(myThid),myBxHi(myThid)
179    
180     #ifdef ALLOW_AUTODIFF_TAMC
181     act1 = bi - myBxLo(myThid)
182     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
183     act2 = bj - myByLo(myThid)
184     max2 = myByHi(myThid) - myByLo(myThid) + 1
185     act3 = myThid - 1
186     max3 = nTx*nTy
187     act4 = ikey_dynamics - 1
188     itdkey = (act1 + 1) + act2*max1
189     & + act3*max1*max2
190     & + act4*max1*max2*max3
191     #endif /* ALLOW_AUTODIFF_TAMC */
192    
193     C-- Set up work arrays with valid (i.e. not NaN) values
194     C These inital values do not alter the numerical results. They
195     C just ensure that all memory references are to valid floating
196     C point numbers. This prevents spurious hardware signals due to
197     C uninitialised but inert locations.
198    
199     DO j=1-OLy,sNy+OLy
200     DO i=1-OLx,sNx+OLx
201     rhok (i,j) = 0. _d 0
202     rhoKM1 (i,j) = 0. _d 0
203     ENDDO
204     ENDDO
205    
206     DO k=1,Nr
207     DO j=1-OLy,sNy+OLy
208     DO i=1-OLx,sNx+OLx
209     C This is currently also used by IVDC and Diagnostics
210     sigmaX(i,j,k) = 0. _d 0
211     sigmaY(i,j,k) = 0. _d 0
212     sigmaR(i,j,k) = 0. _d 0
213     #ifdef ALLOW_AUTODIFF_TAMC
214     cph all the following init. are necessary for TAF
215     cph although some of these are re-initialised later.
216     IVDConvCount(i,j,k,bi,bj) = 0.
217     # ifdef ALLOW_GMREDI
218     Kwx(i,j,k,bi,bj) = 0. _d 0
219     Kwy(i,j,k,bi,bj) = 0. _d 0
220     Kwz(i,j,k,bi,bj) = 0. _d 0
221     # ifdef GM_NON_UNITY_DIAGONAL
222     Kux(i,j,k,bi,bj) = 0. _d 0
223     Kvy(i,j,k,bi,bj) = 0. _d 0
224     # endif
225     # ifdef GM_EXTRA_DIAGONAL
226     Kuz(i,j,k,bi,bj) = 0. _d 0
227     Kvz(i,j,k,bi,bj) = 0. _d 0
228     # endif
229     # ifdef GM_BOLUS_ADVEC
230     GM_PsiX(i,j,k,bi,bj) = 0. _d 0
231     GM_PsiY(i,j,k,bi,bj) = 0. _d 0
232     # endif
233     # ifdef GM_VISBECK_VARIABLE_K
234     VisbeckK(i,j,bi,bj) = 0. _d 0
235     # endif
236     # endif /* ALLOW_GMREDI */
237     #endif /* ALLOW_AUTODIFF_TAMC */
238     ENDDO
239     ENDDO
240     ENDDO
241    
242     iMin = 1-OLx
243     iMax = sNx+OLx
244     jMin = 1-OLy
245     jMax = sNy+OLy
246    
247     #ifdef ALLOW_AUTODIFF_TAMC
248     CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
249     CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
250 heimbach 1.12 CADJ STORE totphihyd(:,:,:,bi,bj)
251 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
252 heimbach 1.10 # ifdef ALLOW_KPP
253 jmc 1.1 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
254     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
255 heimbach 1.10 # endif
256     # ifdef EXACT_CONSERV
257     CADJ STORE pmepr(:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
258     # endif
259 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
260    
261     #ifdef ALLOW_DEBUG
262     IF ( debugLevel .GE. debLevB )
263     & CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
264     #endif
265    
266     C-- Start of diagnostic loop
267     DO k=Nr,1,-1
268    
269     #ifdef ALLOW_AUTODIFF_TAMC
270     C? Patrick, is this formula correct now that we change the loop range?
271     C? Do we still need this?
272     cph kkey formula corrected.
273     cph Needed for rhok, rhokm1, in the case useGMREDI.
274     kkey = (itdkey-1)*Nr + k
275     #endif /* ALLOW_AUTODIFF_TAMC */
276    
277     C-- Calculate gradients of potential density for isoneutral
278     C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
279     c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
280 jmc 1.17 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.)
281     & .OR. doDiagsRho.GE.1 ) THEN
282 jmc 1.1 #ifdef ALLOW_DEBUG
283     IF ( debugLevel .GE. debLevB )
284     & CALL DEBUG_CALL('FIND_RHO',myThid)
285     #endif
286     #ifdef ALLOW_AUTODIFF_TAMC
287     CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
288     CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
289     #endif /* ALLOW_AUTODIFF_TAMC */
290     CALL FIND_RHO(
291     I bi, bj, iMin, iMax, jMin, jMax, k, k,
292     I theta, salt,
293     O rhoK,
294     I myThid )
295    
296     IF (k.GT.1) THEN
297     #ifdef ALLOW_AUTODIFF_TAMC
298     CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
299     CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
300     #endif /* ALLOW_AUTODIFF_TAMC */
301     CALL FIND_RHO(
302     I bi, bj, iMin, iMax, jMin, jMax, k-1, k,
303     I theta, salt,
304     O rhoKm1,
305     I myThid )
306     ENDIF
307     #ifdef ALLOW_DEBUG
308     IF ( debugLevel .GE. debLevB )
309     & CALL DEBUG_CALL('GRAD_SIGMA',myThid)
310     #endif
311     CALL GRAD_SIGMA(
312     I bi, bj, iMin, iMax, jMin, jMax, k,
313     I rhoK, rhoKm1, rhoK,
314     O sigmaX, sigmaY, sigmaR,
315     I myThid )
316     ENDIF
317    
318     #ifdef ALLOW_AUTODIFF_TAMC
319 heimbach 1.12 ctest# ifndef GM_EXCLUDE_CLIPPING
320 jmc 1.1 CADJ STORE rhok (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
321 heimbach 1.12 ctest# endif
322 jmc 1.1 CADJ STORE rhokm1 (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
323     #endif /* ALLOW_AUTODIFF_TAMC */
324     C-- Implicit Vertical Diffusion for Convection
325     c ==> should use sigmaR !!!
326     IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
327     #ifdef ALLOW_DEBUG
328     IF ( debugLevel .GE. debLevB )
329     & CALL DEBUG_CALL('CALC_IVDC',myThid)
330     #endif
331     CALL CALC_IVDC(
332     I bi, bj, iMin, iMax, jMin, jMax, k,
333     I rhoKm1, rhoK,
334     I myTime, myIter, myThid)
335     ENDIF
336    
337 jmc 1.17 #ifdef ALLOW_DIAGNOSTICS
338     IF ( doDiagsRho.GE.2 ) THEN
339     CALL DIAGS_RHO( k, bi, bj,
340     I rhoK, rhoKm1,
341     I myTime, myIter, myThid)
342     ENDIF
343     #endif
344    
345 jmc 1.1 C-- end of diagnostic k loop (Nr:1)
346     ENDDO
347    
348 jmc 1.8 #ifdef ALLOW_DIAGNOSTICS
349 jmc 1.17 c IF ( useDiagnostics .AND.
350     c & (useGMRedi .OR. ivdc_kappa.NE.0.) ) THEN
351     IF ( doDiagsRho.GE.1 ) THEN
352 jmc 1.16 CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR ', 0, Nr,
353     & 2, bi, bj, myThid)
354 jmc 1.8 ENDIF
355     #endif
356    
357 jmc 1.1 #ifdef ALLOW_OBCS
358     C-- Calculate future values on open boundaries
359     IF (useOBCS) THEN
360     #ifdef ALLOW_DEBUG
361     IF ( debugLevel .GE. debLevB )
362     & CALL DEBUG_CALL('OBCS_CALC',myThid)
363     #endif
364 heimbach 1.11 CALL OBCS_CALC( bi, bj, myTime+deltaTclock, myIter+1,
365 jmc 1.1 I uVel, vVel, wVel, theta, salt,
366     I myThid )
367     ENDIF
368     #endif /* ALLOW_OBCS */
369    
370     #ifndef ALLOW_AUTODIFF_TAMC
371 jmc 1.14 IF ( fluidIsWater ) THEN
372 jmc 1.1 #endif
373     C-- Determines forcing terms based on external fields
374     C relaxation terms, etc.
375     #ifdef ALLOW_DEBUG
376     IF ( debugLevel .GE. debLevB )
377     & CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
378     #endif
379     CALL EXTERNAL_FORCING_SURF(
380     I bi, bj, iMin, iMax, jMin, jMax,
381     I myTime, myIter, myThid )
382     #ifndef ALLOW_AUTODIFF_TAMC
383     ENDIF
384     #endif
385    
386     #ifdef ALLOW_AUTODIFF_TAMC
387     cph needed for KPP
388 jmc 1.4 CADJ STORE surfaceForcingU(:,:,bi,bj)
389 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
390 jmc 1.4 CADJ STORE surfaceForcingV(:,:,bi,bj)
391 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
392 jmc 1.4 CADJ STORE surfaceForcingS(:,:,bi,bj)
393 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
394 jmc 1.4 CADJ STORE surfaceForcingT(:,:,bi,bj)
395 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
396 jmc 1.4 CADJ STORE surfaceForcingTice(:,:,bi,bj)
397 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
398 heimbach 1.13
399 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
400    
401     #ifdef ALLOW_GMREDI
402    
403     #ifdef ALLOW_AUTODIFF_TAMC
404 heimbach 1.10 # ifndef GM_EXCLUDE_CLIPPING
405 jmc 1.1 cph storing here is needed only for one GMREDI_OPTIONS:
406     cph define GM_BOLUS_ADVEC
407 heimbach 1.10 cph keep it although TAF says you dont need to.
408 jmc 1.1 cph but I've avoided the #ifdef for now, in case more things change
409     CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
410     CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
411 heimbach 1.12 CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
412 heimbach 1.10 # endif
413 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
414    
415     C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
416     IF (useGMRedi) THEN
417     #ifdef ALLOW_DEBUG
418     IF ( debugLevel .GE. debLevB )
419     & CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
420     #endif
421     CALL GMREDI_CALC_TENSOR(
422     I bi, bj, iMin, iMax, jMin, jMax,
423     I sigmaX, sigmaY, sigmaR,
424     I myThid )
425     #ifdef ALLOW_AUTODIFF_TAMC
426     ELSE
427     CALL GMREDI_CALC_TENSOR_DUMMY(
428     I bi, bj, iMin, iMax, jMin, jMax,
429     I sigmaX, sigmaY, sigmaR,
430     I myThid )
431     #endif /* ALLOW_AUTODIFF_TAMC */
432     ENDIF
433    
434     #endif /* ALLOW_GMREDI */
435    
436     #ifdef ALLOW_KPP
437     C-- Compute KPP mixing coefficients
438     IF (useKPP) THEN
439     #ifdef ALLOW_DEBUG
440     IF ( debugLevel .GE. debLevB )
441     & CALL DEBUG_CALL('KPP_CALC',myThid)
442     #endif
443     CALL KPP_CALC(
444     I bi, bj, myTime, myThid )
445     #ifdef ALLOW_AUTODIFF_TAMC
446     ELSE
447     CALL KPP_CALC_DUMMY(
448     I bi, bj, myTime, myThid )
449     #endif /* ALLOW_AUTODIFF_TAMC */
450     ENDIF
451    
452     #endif /* ALLOW_KPP */
453    
454 mlosch 1.6 #ifdef ALLOW_PP81
455     C-- Compute PP81 mixing coefficients
456     IF (usePP81) THEN
457     #ifdef ALLOW_DEBUG
458     IF ( debugLevel .GE. debLevB )
459     & CALL DEBUG_CALL('PP81_CALC',myThid)
460     #endif
461     CALL PP81_CALC(
462     I bi, bj, myTime, myThid )
463     ENDIF
464     #endif /* ALLOW_PP81 */
465    
466     #ifdef ALLOW_MY82
467     C-- Compute MY82 mixing coefficients
468     IF (useMY82) THEN
469     #ifdef ALLOW_DEBUG
470     IF ( debugLevel .GE. debLevB )
471     & CALL DEBUG_CALL('MY82_CALC',myThid)
472     #endif
473     CALL MY82_CALC(
474     I bi, bj, myTime, myThid )
475     ENDIF
476     #endif /* ALLOW_MY82 */
477    
478 mlosch 1.9 #ifdef ALLOW_GGL90
479     C-- Compute GGL90 mixing coefficients
480     IF (useGGL90) THEN
481     #ifdef ALLOW_DEBUG
482     IF ( debugLevel .GE. debLevB )
483     & CALL DEBUG_CALL('GGL90_CALC',myThid)
484     #endif
485     CALL GGL90_CALC(
486     I bi, bj, myTime, myThid )
487     ENDIF
488     #endif /* ALLOW_GGL90 */
489    
490 jmc 1.20 #ifdef ALLOW_TIMEAVE
491     IF ( taveFreq.GT. 0. _d 0 .AND. fluidIsWater ) THEN
492     CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid)
493     ENDIF
494     IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
495     CALL TIMEAVE_CUMULATE(ConvectCountTave, IVDConvCount,
496     I Nr, deltaTclock, bi, bj, myThid)
497     ENDIF
498     #endif /* ALLOW_TIMEAVE */
499    
500 jmc 1.1 C-- end bi,bj loops.
501     ENDDO
502     ENDDO
503    
504 mlosch 1.22 #ifdef ALLOW_BALANCE_FLUXES
505     C balance fluxes
506     if ( balanceEmPmR )
507     & CALL REMOVE_MEAN_RS( 1, EmPmR, maskH, maskH, rA, drF,
508     & 'EmPmR', myTime, myThid )
509     if ( balanceQnet )
510     & CALL REMOVE_MEAN_RS( 1, Qnet, maskH, maskH, rA, drF,
511     & 'Qnet ', myTime, myThid )
512     #endif /* ALLOW_BALANCE_FLUXES */
513    
514 jmc 1.18 #ifdef ALLOW_DIAGNOSTICS
515     IF ( fluidIsWater .AND. useDiagnostics ) THEN
516     CALL DIAGS_OCEANIC_SURF_FLUX( myTime, myIter, myThid )
517     ENDIF
518 jmc 1.19 IF ( ivdc_kappa.NE.0 .AND. useDiagnostics ) THEN
519     CALL DIAGNOSTICS_FILL( IVDConvCount,'CONVADJ ',
520     & 0, Nr, 0, 1, 1, myThid )
521     ENDIF
522 jmc 1.18 #endif
523    
524 jmc 1.1 #ifdef ALLOW_DEBUG
525     IF ( debugLevel .GE. debLevB )
526     & CALL DEBUG_LEAVE('DO_OCEANIC_PHYS',myThid)
527     #endif
528    
529     RETURN
530     END

  ViewVC Help
Powered by ViewVC 1.1.22