/[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.54 - (show annotations) (download)
Sat Aug 18 21:34:01 2007 UTC (16 years, 9 months ago) by heimbach
Branch: MAIN
Changes since 1.53: +2 -4 lines
Update NLFS adjoint.

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

  ViewVC Help
Powered by ViewVC 1.1.22