/[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.36 - (show annotations) (download)
Fri Dec 29 00:20:12 2006 UTC (17 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58t_post
Changes since 1.35: +41 -41 lines
if fluid is not water, by-pass all the calls, except OBCS_CALC
 which was moved down to the end of bi,bj loops ;

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

  ViewVC Help
Powered by ViewVC 1.1.22