/[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.34 - (show annotations) (download)
Thu Dec 14 22:31:18 2006 UTC (17 years, 5 months ago) by heimbach
Branch: MAIN
Changes since 1.33: +12 -7 lines
Updating seaice adjoint, step 1 (everything, except SEAICE_EVP).

1 C $Header: /u/gcmpack/MITgcm/model/src/do_oceanic_phys.F,v 1.33 2006/11/08 16:01:04 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 ENDIF
171 #ifdef ALLOW_COST_ICE
172 CALL COST_ICE_TEST ( myTime, myIter, myThid )
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 #endif /* ALLOW_AUTODIFF_TAMC */
257
258 C-- Set up work arrays with valid (i.e. not NaN) values
259 C These inital values do not alter the numerical results. They
260 C just ensure that all memory references are to valid floating
261 C point numbers. This prevents spurious hardware signals due to
262 C uninitialised but inert locations.
263
264 DO j=1-OLy,sNy+OLy
265 DO i=1-OLx,sNx+OLx
266 rhok (i,j) = 0. _d 0
267 rhoKM1 (i,j) = 0. _d 0
268 rhoKP1 (i,j) = 0. _d 0
269 ENDDO
270 ENDDO
271
272 DO k=1,Nr
273 DO j=1-OLy,sNy+OLy
274 DO i=1-OLx,sNx+OLx
275 C This is currently also used by IVDC and Diagnostics
276 sigmaX(i,j,k) = 0. _d 0
277 sigmaY(i,j,k) = 0. _d 0
278 sigmaR(i,j,k) = 0. _d 0
279 #ifdef ALLOW_AUTODIFF_TAMC
280 cph all the following init. are necessary for TAF
281 cph although some of these are re-initialised later.
282 IVDConvCount(i,j,k,bi,bj) = 0.
283 # ifdef ALLOW_GMREDI
284 Kwx(i,j,k,bi,bj) = 0. _d 0
285 Kwy(i,j,k,bi,bj) = 0. _d 0
286 Kwz(i,j,k,bi,bj) = 0. _d 0
287 # ifdef GM_NON_UNITY_DIAGONAL
288 Kux(i,j,k,bi,bj) = 0. _d 0
289 Kvy(i,j,k,bi,bj) = 0. _d 0
290 # endif
291 # ifdef GM_EXTRA_DIAGONAL
292 Kuz(i,j,k,bi,bj) = 0. _d 0
293 Kvz(i,j,k,bi,bj) = 0. _d 0
294 # endif
295 # ifdef GM_BOLUS_ADVEC
296 GM_PsiX(i,j,k,bi,bj) = 0. _d 0
297 GM_PsiY(i,j,k,bi,bj) = 0. _d 0
298 # endif
299 # ifdef GM_VISBECK_VARIABLE_K
300 VisbeckK(i,j,bi,bj) = 0. _d 0
301 # endif
302 # endif /* ALLOW_GMREDI */
303 #endif /* ALLOW_AUTODIFF_TAMC */
304 ENDDO
305 ENDDO
306 ENDDO
307
308 iMin = 1-OLx
309 iMax = sNx+OLx
310 jMin = 1-OLy
311 jMax = sNy+OLy
312
313 #ifdef ALLOW_AUTODIFF_TAMC
314 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
315 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
316 CADJ STORE totphihyd(:,:,:,bi,bj)
317 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
318 # ifdef ALLOW_KPP
319 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
320 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
321 # endif
322 #endif /* ALLOW_AUTODIFF_TAMC */
323
324 #ifdef ALLOW_DEBUG
325 IF ( debugLevel .GE. debLevB )
326 & CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
327 #endif
328
329 C-- Start of diagnostic loop
330 DO k=Nr,1,-1
331
332 #ifdef ALLOW_AUTODIFF_TAMC
333 C? Patrick, is this formula correct now that we change the loop range?
334 C? Do we still need this?
335 cph kkey formula corrected.
336 cph Needed for rhok, rhokm1, in the case useGMREDI.
337 kkey = (itdkey-1)*Nr + k
338 #endif /* ALLOW_AUTODIFF_TAMC */
339
340 C-- Calculate gradients of potential density for isoneutral
341 C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
342 c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
343 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.)
344 & .OR. doDiagsRho.GE.1 ) THEN
345 #ifdef ALLOW_DEBUG
346 IF ( debugLevel .GE. debLevB )
347 & CALL DEBUG_CALL('FIND_RHO',myThid)
348 #endif
349 #ifdef ALLOW_AUTODIFF_TAMC
350 CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
351 CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
352 #endif /* ALLOW_AUTODIFF_TAMC */
353 CALL FIND_RHO(
354 I bi, bj, iMin, iMax, jMin, jMax, k, k,
355 I theta, salt,
356 O rhoK,
357 I myThid )
358
359 IF (k.GT.1) THEN
360 #ifdef ALLOW_AUTODIFF_TAMC
361 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
362 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
363 #endif /* ALLOW_AUTODIFF_TAMC */
364 CALL FIND_RHO(
365 I bi, bj, iMin, iMax, jMin, jMax, k-1, k,
366 I theta, salt,
367 O rhoKm1,
368 I myThid )
369 ENDIF
370 #ifdef ALLOW_DEBUG
371 IF ( debugLevel .GE. debLevB )
372 & CALL DEBUG_CALL('GRAD_SIGMA',myThid)
373 #endif
374 cph Avoid variable aliasing for adjoint !!!
375 DO j=jMin,jMax
376 DO i=iMin,iMax
377 rhoKP1(i,j) = rhoK(i,j)
378 ENDDO
379 ENDDO
380 CALL GRAD_SIGMA(
381 I bi, bj, iMin, iMax, jMin, jMax, k,
382 I rhoK, rhoKm1, rhoKp1,
383 O sigmaX, sigmaY, sigmaR,
384 I myThid )
385 ENDIF
386
387 #ifdef ALLOW_AUTODIFF_TAMC
388 ctest# ifndef GM_EXCLUDE_CLIPPING
389 CADJ STORE rhok (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
390 ctest# endif
391 CADJ STORE rhokm1 (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
392 #endif /* ALLOW_AUTODIFF_TAMC */
393 C-- Implicit Vertical Diffusion for Convection
394 c ==> should use sigmaR !!!
395 IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
396 #ifdef ALLOW_DEBUG
397 IF ( debugLevel .GE. debLevB )
398 & CALL DEBUG_CALL('CALC_IVDC',myThid)
399 #endif
400 CALL CALC_IVDC(
401 I bi, bj, iMin, iMax, jMin, jMax, k,
402 I rhoKm1, rhoK,
403 I myTime, myIter, myThid)
404 ENDIF
405
406 #ifdef ALLOW_DIAGNOSTICS
407 IF ( doDiagsRho.GE.2 ) THEN
408 CALL DIAGS_RHO( k, bi, bj,
409 I rhoK, rhoKm1,
410 I myTime, myIter, myThid)
411 ENDIF
412 #endif
413
414 C-- end of diagnostic k loop (Nr:1)
415 ENDDO
416
417 #ifdef ALLOW_DIAGNOSTICS
418 c IF ( useDiagnostics .AND.
419 c & (useGMRedi .OR. ivdc_kappa.NE.0.) ) THEN
420 IF ( doDiagsRho.GE.1 ) THEN
421 CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR ', 0, Nr,
422 & 2, bi, bj, myThid)
423 ENDIF
424 #endif
425
426 #ifdef ALLOW_OBCS
427 C-- Calculate future values on open boundaries
428 IF (useOBCS) THEN
429 #ifdef ALLOW_DEBUG
430 IF ( debugLevel .GE. debLevB )
431 & CALL DEBUG_CALL('OBCS_CALC',myThid)
432 #endif
433 CALL OBCS_CALC( bi, bj, myTime+deltaTclock, myIter+1,
434 I uVel, vVel, wVel, theta, salt,
435 I myThid )
436 ENDIF
437 #endif /* ALLOW_OBCS */
438
439 #ifndef ALLOW_AUTODIFF_TAMC
440 IF ( fluidIsWater ) THEN
441 #endif
442 C-- Determines forcing terms based on external fields
443 C relaxation terms, etc.
444 #ifdef ALLOW_DEBUG
445 IF ( debugLevel .GE. debLevB )
446 & CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
447 #endif
448 #ifdef ALLOW_AUTODIFF_TAMC
449 CADJ STORE EmPmR(:,:,bi,bj)
450 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
451 # ifdef EXACT_CONSERV
452 CADJ STORE PmEpR(:,:,bi,bj)
453 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
454 # endif
455 # ifdef NONLIN_FRSURF
456 CADJ STORE hFac_surfC(:,:,bi,bj)
457 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
458 CADJ STORE recip_hFacC(:,:,:,bi,bj)
459 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
460 # endif
461 #endif
462 CALL EXTERNAL_FORCING_SURF(
463 I bi, bj, iMin, iMax, jMin, jMax,
464 I myTime, myIter, myThid )
465 #ifndef ALLOW_AUTODIFF_TAMC
466 ENDIF
467 #endif
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_GMREDI
491
492 #ifdef ALLOW_AUTODIFF_TAMC
493 # ifndef GM_EXCLUDE_CLIPPING
494 cph storing here is needed only for one GMREDI_OPTIONS:
495 cph define GM_BOLUS_ADVEC
496 cph keep it although TAF says you dont need to.
497 cph but I've avoided the #ifdef for now, in case more things change
498 CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
499 CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
500 CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
501 # endif
502 #endif /* ALLOW_AUTODIFF_TAMC */
503
504 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
505 IF (useGMRedi) THEN
506 #ifdef ALLOW_DEBUG
507 IF ( debugLevel .GE. debLevB )
508 & CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
509 #endif
510 CALL GMREDI_CALC_TENSOR(
511 I bi, bj, iMin, iMax, jMin, jMax,
512 I sigmaX, sigmaY, sigmaR,
513 I myThid )
514 #ifdef ALLOW_AUTODIFF_TAMC
515 ELSE
516 CALL GMREDI_CALC_TENSOR_DUMMY(
517 I bi, bj, iMin, iMax, jMin, jMax,
518 I sigmaX, sigmaY, sigmaR,
519 I myThid )
520 #endif /* ALLOW_AUTODIFF_TAMC */
521 ENDIF
522
523 #endif /* ALLOW_GMREDI */
524
525 #ifdef ALLOW_KPP
526 C-- Compute KPP mixing coefficients
527 IF (useKPP) THEN
528 #ifdef ALLOW_DEBUG
529 IF ( debugLevel .GE. debLevB )
530 & CALL DEBUG_CALL('KPP_CALC',myThid)
531 #endif
532 CALL KPP_CALC(
533 I bi, bj, myTime, myThid )
534 #ifdef ALLOW_AUTODIFF_TAMC
535 ELSE
536 CALL KPP_CALC_DUMMY(
537 I bi, bj, myTime, myThid )
538 #endif /* ALLOW_AUTODIFF_TAMC */
539 ENDIF
540
541 #endif /* ALLOW_KPP */
542
543 #ifdef ALLOW_PP81
544 C-- Compute PP81 mixing coefficients
545 IF (usePP81) THEN
546 #ifdef ALLOW_DEBUG
547 IF ( debugLevel .GE. debLevB )
548 & CALL DEBUG_CALL('PP81_CALC',myThid)
549 #endif
550 CALL PP81_CALC(
551 I bi, bj, myTime, myThid )
552 ENDIF
553 #endif /* ALLOW_PP81 */
554
555 #ifdef ALLOW_MY82
556 C-- Compute MY82 mixing coefficients
557 IF (useMY82) THEN
558 #ifdef ALLOW_DEBUG
559 IF ( debugLevel .GE. debLevB )
560 & CALL DEBUG_CALL('MY82_CALC',myThid)
561 #endif
562 CALL MY82_CALC(
563 I bi, bj, myTime, myThid )
564 ENDIF
565 #endif /* ALLOW_MY82 */
566
567 #ifdef ALLOW_GGL90
568 C-- Compute GGL90 mixing coefficients
569 IF (useGGL90) THEN
570 #ifdef ALLOW_DEBUG
571 IF ( debugLevel .GE. debLevB )
572 & CALL DEBUG_CALL('GGL90_CALC',myThid)
573 #endif
574 CALL GGL90_CALC(
575 I bi, bj, myTime, myThid )
576 ENDIF
577 #endif /* ALLOW_GGL90 */
578
579 #ifdef ALLOW_TIMEAVE
580 IF ( taveFreq.GT. 0. _d 0 .AND. fluidIsWater ) THEN
581 CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid)
582 ENDIF
583 IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
584 CALL TIMEAVE_CUMULATE(ConvectCountTave, IVDConvCount,
585 I Nr, deltaTclock, bi, bj, myThid)
586 ENDIF
587 #endif /* ALLOW_TIMEAVE */
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