/[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.47 - (show annotations) (download)
Sun Jun 3 22:53:50 2007 UTC (17 years ago) by jmc
Branch: MAIN
Changes since 1.46: +56 -50 lines
- move GMRedi after KPP and all vertical mixing scheme.
- add call for computing mixed-layer depth (but hidden from TAF for now)

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