/[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.58 - (show annotations) (download)
Fri Sep 28 19:42:06 2007 UTC (16 years, 8 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint59i, checkpoint59h, checkpoint59j
Changes since 1.57: +2 -5 lines
Put back IF(useSEAICE) for adjoint (seems benign).

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

  ViewVC Help
Powered by ViewVC 1.1.22