/[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.49 - (show annotations) (download)
Mon Jun 4 18:48:17 2007 UTC (17 years ago) by heimbach
Branch: MAIN
Changes since 1.48: +2 -18 lines
Re-arrange store directives.

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

  ViewVC Help
Powered by ViewVC 1.1.22