/[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.53 - (show annotations) (download)
Fri Jul 27 22:18:57 2007 UTC (16 years, 10 months ago) by heimbach
Branch: MAIN
Changes since 1.52: +2 -3 lines
Comment all relevant #ifndef ALLOW_AUTODIFF_TAMC
that used to hide exch2 or cubed-sphere specific code
(commented via 'cph-exch2')

1 C $Header: /u/gcmpack/MITgcm/model/src/do_oceanic_phys.F,v 1.52 2007/07/22 23:51:16 dimitri 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 SEAICE_ALLOW_EVP
156 CADJ STORE seaice_sigma1 = comlev1, key = ikey_dynamics
157 CADJ STORE seaice_sigma2 = comlev1, key = ikey_dynamics
158 CADJ STORE seaice_sigma12 = comlev1, key = ikey_dynamics
159 # endif
160 # ifdef ATMOSPHERIC_LOADING
161 CADJ STORE siceload = comlev1, key = ikey_dynamics
162 # endif
163 #endif
164 #ifdef ALLOW_DEBUG
165 IF ( debugLevel .GE. debLevB )
166 & CALL DEBUG_CALL('SEAICE_MODEL',myThid)
167 #endif
168 CALL TIMER_START('SEAICE_MODEL [DO_OCEANIC_PHYS]', myThid)
169 CALL SEAICE_MODEL( myTime, myIter, myThid )
170 CALL TIMER_STOP ('SEAICE_MODEL [DO_OCEANIC_PHYS]', myThid)
171 #ifdef ALLOW_COST_ICE
172 CALL COST_ICE_TEST ( myTime, myIter, myThid )
173 #endif
174 ENDIF
175 #endif /* ALLOW_SEAICE */
176
177 #if (defined ALLOW_THSICE) && !(defined ALLOW_ATM2D)
178 IF ( useThSIce .AND. fluidIsWater ) THEN
179 #ifdef ALLOW_DEBUG
180 IF ( debugLevel .GE. debLevB )
181 & CALL DEBUG_CALL('THSICE_MAIN',myThid)
182 #endif
183 C-- Step forward Therm.Sea-Ice variables
184 C and modify forcing terms including effects from ice
185 CALL TIMER_START('THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
186 CALL THSICE_MAIN( myTime, myIter, myThid )
187 CALL TIMER_STOP( 'THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
188 ENDIF
189 #endif /* ALLOW_THSICE */
190
191 #ifdef ALLOW_SHELFICE
192 IF ( useShelfIce .AND. fluidIsWater ) THEN
193 #ifdef ALLOW_DEBUG
194 IF ( debugLevel .GE. debLevB )
195 & CALL DEBUG_CALL('SHELFICE_THERMODYNAMICS',myThid)
196 #endif
197 C compute temperature and (virtual) salt flux at the
198 C shelf-ice ocean interface
199 CALL TIMER_START('SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',
200 & myThid)
201 CALL SHELFICE_THERMODYNAMICS( myTime, myIter, myThid )
202 CALL TIMER_STOP( 'SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',
203 & myThid)
204 ENDIF
205 #endif /* ALLOW_SHELFICE */
206
207 C-- Freeze water at the surface
208 #ifdef ALLOW_AUTODIFF_TAMC
209 CADJ STORE theta = comlev1, key = ikey_dynamics
210 #endif
211 IF ( allowFreezing
212 & .AND. .NOT. useSEAICE
213 & .AND. .NOT. useThSIce ) THEN
214 CALL FREEZE_SURFACE( myTime, myIter, myThid )
215 ENDIF
216
217 #ifdef ALLOW_OCN_COMPON_INTERF
218 C-- Apply imported data (from coupled interface) to forcing fields
219 C jmc: do not know precisely where to put this call (bf or af thSIce ?)
220 IF ( useCoupler ) THEN
221 CALL OCN_APPLY_IMPORT( .TRUE., myTime, myIter, myThid )
222 ENDIF
223 #endif /* ALLOW_OCN_COMPON_INTERF */
224
225 #ifdef ALLOW_BALANCE_FLUXES
226 C balance fluxes
227 IF ( balanceEmPmR )
228 & CALL REMOVE_MEAN_RS( 1, EmPmR, maskH, maskH, rA, drF,
229 & 'EmPmR', myTime, myThid )
230 IF ( balanceQnet )
231 & CALL REMOVE_MEAN_RS( 1, Qnet, maskH, maskH, rA, drF,
232 & 'Qnet ', myTime, myThid )
233 #endif /* ALLOW_BALANCE_FLUXES */
234
235 #ifdef ALLOW_AUTODIFF_TAMC
236 C-- HPF directive to help TAMC
237 CHPF$ INDEPENDENT
238 #endif /* ALLOW_AUTODIFF_TAMC */
239 DO bj=myByLo(myThid),myByHi(myThid)
240 #ifdef ALLOW_AUTODIFF_TAMC
241 C-- HPF directive to help TAMC
242 CHPF$ INDEPENDENT
243 #endif /* ALLOW_AUTODIFF_TAMC */
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 itdkey = (act1 + 1) + act2*max1
255 & + act3*max1*max2
256 & + act4*max1*max2*max3
257 #else /* ALLOW_AUTODIFF_TAMC */
258 C if fluid is not water, by-pass find_rho, gmredi, surfaceForcing
259 C and all vertical mixing schemes, but keep OBCS_CALC
260 IF ( fluidIsWater ) THEN
261 #endif /* ALLOW_AUTODIFF_TAMC */
262
263 C-- Set up work arrays with valid (i.e. not NaN) values
264 C These inital values do not alter the numerical results. They
265 C just ensure that all memory references are to valid floating
266 C point numbers. This prevents spurious hardware signals due to
267 C uninitialised but inert locations.
268
269 DO j=1-OLy,sNy+OLy
270 DO i=1-OLx,sNx+OLx
271 rhoK (i,j) = 0. _d 0
272 rhoKm1 (i,j) = 0. _d 0
273 rhoKp1 (i,j) = 0. _d 0
274 ENDDO
275 ENDDO
276
277 DO k=1,Nr
278 DO j=1-OLy,sNy+OLy
279 DO i=1-OLx,sNx+OLx
280 C This is currently also used by IVDC and Diagnostics
281 sigmaX(i,j,k) = 0. _d 0
282 sigmaY(i,j,k) = 0. _d 0
283 sigmaR(i,j,k) = 0. _d 0
284 #ifdef ALLOW_AUTODIFF_TAMC
285 cph all the following init. are necessary for TAF
286 cph although some of these are re-initialised later.
287 IVDConvCount(i,j,k,bi,bj) = 0.
288 # ifdef ALLOW_GMREDI
289 Kwx(i,j,k,bi,bj) = 0. _d 0
290 Kwy(i,j,k,bi,bj) = 0. _d 0
291 Kwz(i,j,k,bi,bj) = 0. _d 0
292 # ifdef GM_NON_UNITY_DIAGONAL
293 Kux(i,j,k,bi,bj) = 0. _d 0
294 Kvy(i,j,k,bi,bj) = 0. _d 0
295 # endif
296 # ifdef GM_EXTRA_DIAGONAL
297 Kuz(i,j,k,bi,bj) = 0. _d 0
298 Kvz(i,j,k,bi,bj) = 0. _d 0
299 # endif
300 # ifdef GM_BOLUS_ADVEC
301 GM_PsiX(i,j,k,bi,bj) = 0. _d 0
302 GM_PsiY(i,j,k,bi,bj) = 0. _d 0
303 # endif
304 # ifdef GM_VISBECK_VARIABLE_K
305 VisbeckK(i,j,bi,bj) = 0. _d 0
306 # endif
307 # endif /* ALLOW_GMREDI */
308 # ifdef ALLOW_KPP
309 KPPdiffKzS(i,j,k,bi,bj) = 0. _d 0
310 KPPdiffKzT(i,j,k,bi,bj) = 0. _d 0
311 # endif /* ALLOW_KPP */
312 #endif /* ALLOW_AUTODIFF_TAMC */
313 ENDDO
314 ENDDO
315 ENDDO
316
317 iMin = 1-OLx
318 iMax = sNx+OLx
319 jMin = 1-OLy
320 jMax = sNy+OLy
321
322 #ifdef ALLOW_AUTODIFF_TAMC
323 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
324 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
325 CADJ STORE totphihyd(:,:,:,bi,bj)
326 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
327 # ifdef ALLOW_KPP
328 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
329 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
330 # endif
331 #endif /* ALLOW_AUTODIFF_TAMC */
332
333 #ifdef ALLOW_DEBUG
334 IF ( debugLevel .GE. debLevB )
335 & CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
336 #endif
337
338 C-- Start of diagnostic loop
339 DO k=Nr,1,-1
340
341 #ifdef ALLOW_AUTODIFF_TAMC
342 C? Patrick, is this formula correct now that we change the loop range?
343 C? Do we still need this?
344 cph kkey formula corrected.
345 cph Needed for rhoK, rhoKm1, in the case useGMREDI.
346 kkey = (itdkey-1)*Nr + k
347 #endif /* ALLOW_AUTODIFF_TAMC */
348
349 C-- Calculate gradients of potential density for isoneutral
350 C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
351 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.)
352 & .OR. doDiagsRho.GE.1 ) THEN
353 #ifdef ALLOW_DEBUG
354 IF ( debugLevel .GE. debLevB )
355 & CALL DEBUG_CALL('FIND_RHO',myThid)
356 #endif
357 #ifdef ALLOW_AUTODIFF_TAMC
358 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
359 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
360 #endif /* ALLOW_AUTODIFF_TAMC */
361 CALL FIND_RHO(
362 I bi, bj, iMin, iMax, jMin, jMax, k, k,
363 I theta, salt,
364 O rhoK,
365 I myThid )
366
367 IF (k.GT.1) THEN
368 #ifdef ALLOW_AUTODIFF_TAMC
369 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
370 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
371 #endif /* ALLOW_AUTODIFF_TAMC */
372 CALL FIND_RHO(
373 I bi, bj, iMin, iMax, jMin, jMax, k-1, k,
374 I theta, salt,
375 O rhoKm1,
376 I myThid )
377 ENDIF
378 #ifdef ALLOW_DEBUG
379 IF ( debugLevel .GE. debLevB )
380 & CALL DEBUG_CALL('GRAD_SIGMA',myThid)
381 #endif
382 cph Avoid variable aliasing for adjoint !!!
383 DO j=jMin,jMax
384 DO i=iMin,iMax
385 rhoKp1(i,j) = rhoK(i,j)
386 ENDDO
387 ENDDO
388 CALL GRAD_SIGMA(
389 I bi, bj, iMin, iMax, jMin, jMax, k,
390 I rhoK, rhoKm1, rhoKp1,
391 O sigmaX, sigmaY, sigmaR,
392 I myThid )
393 ENDIF
394
395 #ifdef ALLOW_AUTODIFF_TAMC
396 ctest# ifndef GM_EXCLUDE_CLIPPING
397 CADJ STORE rhoK (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
398 ctest# endif
399 CADJ STORE rhoKm1 (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
400 #endif /* ALLOW_AUTODIFF_TAMC */
401 C-- Implicit Vertical Diffusion for Convection
402 c ==> should use sigmaR !!!
403 IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
404 #ifdef ALLOW_DEBUG
405 IF ( debugLevel .GE. debLevB )
406 & CALL DEBUG_CALL('CALC_IVDC',myThid)
407 #endif
408 CALL CALC_IVDC(
409 I bi, bj, iMin, iMax, jMin, jMax, k,
410 I rhoKm1, rhoK,
411 I myTime, myIter, myThid)
412 ENDIF
413
414 #ifdef ALLOW_DIAGNOSTICS
415 IF ( doDiagsRho.GE.2 ) THEN
416 CALL DIAGS_RHO( k, bi, bj,
417 I rhoK, rhoKm1,
418 I myTime, myIter, myThid)
419 ENDIF
420 #endif
421
422 C-- end of diagnostic k loop (Nr:1)
423 ENDDO
424
425 C-- Diagnose Mixed Layer Depth:
426 IF ( useGMRedi .OR. doDiagsRho.GE.1 ) THEN
427 CALL CALC_OCE_MXLAYER( rhoK, sigmaR,
428 & bi, bj, myTime, myIter, myThid )
429 ENDIF
430
431 #ifdef ALLOW_SALT_PLUME
432 CALL CALC_SALT_PLUME_DEPTH( rhoK, sigmaR,
433 & bi, bj, myTime, myIter, myThid )
434 #endif
435 #ifdef ALLOW_DIAGNOSTICS
436 IF ( doDiagsRho.GE.1 ) THEN
437 CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR ', 0, Nr,
438 & 2, bi, bj, myThid)
439 ENDIF
440 #endif
441
442 C-- Determines forcing terms based on external fields
443 C relaxation terms, etc.
444 #ifdef ALLOW_DEBUG
445 IF ( debugLevel .GE. debLevB )
446 & CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
447 #endif
448 #ifdef ALLOW_AUTODIFF_TAMC
449 CADJ STORE EmPmR(:,:,bi,bj)
450 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
451 # ifdef EXACT_CONSERV
452 CADJ STORE PmEpR(:,:,bi,bj)
453 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
454 # endif
455 # ifdef NONLIN_FRSURF
456 CADJ STORE hFac_surfC(:,:,bi,bj)
457 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
458 CADJ STORE recip_hFacC(:,:,:,bi,bj)
459 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
460 # endif
461 #endif
462 CALL EXTERNAL_FORCING_SURF(
463 I bi, bj, iMin, iMax, jMin, jMax,
464 I myTime, myIter, myThid )
465 #ifdef ALLOW_AUTODIFF_TAMC
466 # ifdef EXACT_CONSERV
467 cph-test
468 cphCADJ STORE PmEpR(:,:,bi,bj)
469 cphCADJ & = comlev1_bibj, key=itdkey, byte=isbyte
470 # endif
471 #endif
472
473 #ifdef ALLOW_AUTODIFF_TAMC
474 cph needed for KPP
475 CADJ STORE surfaceForcingU(:,:,bi,bj)
476 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
477 CADJ STORE surfaceForcingV(:,:,bi,bj)
478 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
479 CADJ STORE surfaceForcingS(:,:,bi,bj)
480 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
481 CADJ STORE surfaceForcingT(:,:,bi,bj)
482 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
483 CADJ STORE surfaceForcingTice(:,:,bi,bj)
484 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
485 #endif /* ALLOW_AUTODIFF_TAMC */
486
487 #ifdef ALLOW_KPP
488 C-- Compute KPP mixing coefficients
489 IF (useKPP) THEN
490 #ifdef ALLOW_DEBUG
491 IF ( debugLevel .GE. debLevB )
492 & CALL DEBUG_CALL('KPP_CALC',myThid)
493 #endif
494 CALL KPP_CALC(
495 I bi, bj, myTime, myIter, myThid )
496 #ifdef ALLOW_AUTODIFF_TAMC
497 ELSE
498 CALL KPP_CALC_DUMMY(
499 I bi, bj, myTime, myIter, myThid )
500 #endif /* ALLOW_AUTODIFF_TAMC */
501 ENDIF
502
503 #endif /* ALLOW_KPP */
504
505 #ifdef ALLOW_PP81
506 C-- Compute PP81 mixing coefficients
507 IF (usePP81) THEN
508 #ifdef ALLOW_DEBUG
509 IF ( debugLevel .GE. debLevB )
510 & CALL DEBUG_CALL('PP81_CALC',myThid)
511 #endif
512 CALL PP81_CALC(
513 I bi, bj, myTime, myThid )
514 ENDIF
515 #endif /* ALLOW_PP81 */
516
517 #ifdef ALLOW_MY82
518 C-- Compute MY82 mixing coefficients
519 IF (useMY82) THEN
520 #ifdef ALLOW_DEBUG
521 IF ( debugLevel .GE. debLevB )
522 & CALL DEBUG_CALL('MY82_CALC',myThid)
523 #endif
524 CALL MY82_CALC(
525 I bi, bj, myTime, myThid )
526 ENDIF
527 #endif /* ALLOW_MY82 */
528
529 #ifdef ALLOW_GGL90
530 C-- Compute GGL90 mixing coefficients
531 IF (useGGL90) THEN
532 #ifdef ALLOW_DEBUG
533 IF ( debugLevel .GE. debLevB )
534 & CALL DEBUG_CALL('GGL90_CALC',myThid)
535 #endif
536 CALL GGL90_CALC(
537 I bi, bj, myTime, myThid )
538 ENDIF
539 #endif /* ALLOW_GGL90 */
540
541 #ifdef ALLOW_TIMEAVE
542 IF ( taveFreq.GT. 0. _d 0 ) THEN
543 CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid)
544 ENDIF
545 IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
546 CALL TIMEAVE_CUMULATE(ConvectCountTave, IVDConvCount,
547 I Nr, deltaTclock, bi, bj, myThid)
548 ENDIF
549 #endif /* ALLOW_TIMEAVE */
550
551 #ifdef ALLOW_GMREDI
552 #ifdef ALLOW_AUTODIFF_TAMC
553 # ifndef GM_EXCLUDE_CLIPPING
554 cph storing here is needed only for one GMREDI_OPTIONS:
555 cph define GM_BOLUS_ADVEC
556 cph keep it although TAF says you dont need to.
557 cph but I've avoided the #ifdef for now, in case more things change
558 CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
559 CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
560 CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
561 # endif
562 #endif /* ALLOW_AUTODIFF_TAMC */
563
564 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
565 IF (useGMRedi) THEN
566 #ifdef ALLOW_DEBUG
567 IF ( debugLevel .GE. debLevB )
568 & CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
569 #endif
570 CALL GMREDI_CALC_TENSOR(
571 c I bi, bj, iMin, iMax, jMin, jMax,
572 c I sigmaX, sigmaY, sigmaR,
573 c I myThid )
574 I iMin, iMax, jMin, jMax,
575 I sigmaX, sigmaY, sigmaR,
576 I bi, bj, myTime, myIter, myThid )
577 #ifdef ALLOW_AUTODIFF_TAMC
578 ELSE
579 CALL GMREDI_CALC_TENSOR_DUMMY(
580 c I bi, bj, iMin, iMax, jMin, jMax,
581 c I sigmaX, sigmaY, sigmaR,
582 c I myThid )
583 I iMin, iMax, jMin, jMax,
584 I sigmaX, sigmaY, sigmaR,
585 I bi, bj, myTime, myIter, myThid )
586 #endif /* ALLOW_AUTODIFF_TAMC */
587 ENDIF
588 #endif /* ALLOW_GMREDI */
589
590 #ifndef ALLOW_AUTODIFF_TAMC
591 C--- if fluid Is Water: end
592 ENDIF
593 #endif
594
595 #ifdef ALLOW_OBCS
596 C-- Calculate future values on open boundaries
597 IF (useOBCS) THEN
598 #ifdef ALLOW_DEBUG
599 IF ( debugLevel .GE. debLevB )
600 & CALL DEBUG_CALL('OBCS_CALC',myThid)
601 #endif
602 CALL OBCS_CALC( bi, bj, myTime+deltaTclock, myIter+1,
603 I uVel, vVel, wVel, theta, salt,
604 I myThid )
605 ENDIF
606 #endif /* ALLOW_OBCS */
607
608 C-- end bi,bj loops.
609 ENDDO
610 ENDDO
611
612 #ifdef ALLOW_KPP
613 IF (useKPP) THEN
614 CALL KPP_DO_EXCH( myThid )
615 ENDIF
616 #endif /* ALLOW_KPP */
617
618 #ifdef ALLOW_DIAGNOSTICS
619 IF ( fluidIsWater .AND. useDiagnostics ) THEN
620 CALL DIAGS_OCEANIC_SURF_FLUX( myTime, myIter, myThid )
621 ENDIF
622 IF ( ivdc_kappa.NE.0 .AND. useDiagnostics ) THEN
623 CALL DIAGNOSTICS_FILL( IVDConvCount,'CONVADJ ',
624 & 0, Nr, 0, 1, 1, myThid )
625 ENDIF
626 #endif
627
628 #ifdef ALLOW_DEBUG
629 IF ( debugLevel .GE. debLevB )
630 & CALL DEBUG_LEAVE('DO_OCEANIC_PHYS',myThid)
631 #endif
632
633 RETURN
634 END

  ViewVC Help
Powered by ViewVC 1.1.22