/[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.62 - (show annotations) (download)
Tue Apr 22 15:18:00 2008 UTC (16 years, 1 month ago) by heimbach
Branch: MAIN
Changes since 1.61: +22 -14 lines
Some stuff needed for latest v3 production.

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

  ViewVC Help
Powered by ViewVC 1.1.22