/[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.67 - (show annotations) (download)
Wed Jun 11 18:29:08 2008 UTC (16 years ago) by gforget
Branch: MAIN
CVS Tags: checkpoint60, checkpoint61, checkpoint61b, checkpoint61a
Changes since 1.66: +4 -4 lines
More flexibility in zeroadj.

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

  ViewVC Help
Powered by ViewVC 1.1.22