/[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.40 - (show annotations) (download)
Mon Apr 16 23:31:59 2007 UTC (17 years, 1 month ago) by jmc
Branch: MAIN
Changes since 1.39: +4 -4 lines
move EXF header files from lower_case.h to UPPER_CASE.h ;

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

  ViewVC Help
Powered by ViewVC 1.1.22