/[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.144 - (show annotations) (download)
Sun Nov 2 21:21:59 2014 UTC (9 years, 7 months ago) by gforget
Branch: MAIN
CVS Tags: checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65l, checkpoint65m, checkpoint65g
Changes since 1.143: +2 -2 lines
- packages_boot.F :
    - add useCTRL, useECCO to run time parameters list.
    - default is true if ALLOW_AUTODIFF, false otherwise.
    - IF (useECCO) useCAL = .TRUE.
    - IF (useGrdchk) useCTRL = .TRUE.
- do_oceanic_phys.F, packages_check.F, packages_init_fixed.F,
  packages_init_variables.F, the_main_loop.F : add useECCO switch
- forward_step.F, load_fields_driver.F, packages_check.F, packages_init_fixed.F,
  packages_init_variables.F, the_model_main.F : add useCTRL switch
- initialise_varia.F : add ALLOW_CTRL bracket

1 C $Header: /u/gcmpack/MITgcm/model/src/do_oceanic_phys.F,v 1.143 2014/09/04 22:30:46 jmc Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6 #ifdef ALLOW_AUTODIFF
7 # include "AUTODIFF_OPTIONS.h"
8 #endif
9 #ifdef ALLOW_CTRL
10 # include "CTRL_OPTIONS.h"
11 #endif
12 #ifdef ALLOW_SALT_PLUME
13 # include "SALT_PLUME_OPTIONS.h"
14 #endif
15 #ifdef ALLOW_ECCO
16 # include "ECCO_OPTIONS.h"
17 #endif
18
19 #ifdef ALLOW_AUTODIFF
20 # ifdef ALLOW_GMREDI
21 # include "GMREDI_OPTIONS.h"
22 # endif
23 # ifdef ALLOW_KPP
24 # include "KPP_OPTIONS.h"
25 # endif
26 # ifdef ALLOW_SEAICE
27 # include "SEAICE_OPTIONS.h"
28 # endif
29 # ifdef ALLOW_EXF
30 # include "EXF_OPTIONS.h"
31 # endif
32 #endif /* ALLOW_AUTODIFF */
33
34 CBOP
35 C !ROUTINE: DO_OCEANIC_PHYS
36 C !INTERFACE:
37 SUBROUTINE DO_OCEANIC_PHYS(myTime, myIter, myThid)
38 C !DESCRIPTION: \bv
39 C *==========================================================*
40 C | SUBROUTINE DO_OCEANIC_PHYS
41 C | o Controlling routine for oceanic physics and
42 C | parameterization
43 C *==========================================================*
44 C | o originally, part of S/R thermodynamics
45 C *==========================================================*
46 C \ev
47
48 C !CALLING SEQUENCE:
49 C DO_OCEANIC_PHYS
50 C |
51 C |-- OBCS_CALC
52 C |
53 C |-- FRAZIL_CALC_RHS
54 C |
55 C |-- THSICE_MAIN
56 C |
57 C |-- SEAICE_FAKE
58 C |-- SEAICE_MODEL
59 C |-- SEAICE_COST_SENSI
60 C |
61 C |-- SHELFICE_THERMODYNAMICS
62 C |
63 C |-- ICEFRONT_THERMODYNAMICS
64 C |
65 C |-- SALT_PLUME_DO_EXCH
66 C |
67 C |-- FREEZE_SURFACE
68 C |
69 C |-- OCN_APPLY_IMPORT
70 C |
71 C |-- EXTERNAL_FORCING_SURF
72 C |
73 C |- k loop (Nr:1):
74 C | - DWNSLP_CALC_RHO
75 C | - BBL_CALC_RHO
76 C | - FIND_RHO_2D @ p(k)
77 C | - FIND_RHO_2D @ p(k-1)
78 C | - GRAD_SIGMA
79 C | - CALC_IVDC
80 C | - DIAGS_RHO_L
81 C |- end k loop.
82 C |
83 C |-- CALC_OCE_MXLAYER
84 C |
85 C |-- SALT_PLUME_CALC_DEPTH
86 C |-- SALT_PLUME_VOLFRAC
87 C |-- SALT_PLUME_APPLY
88 C |-- SALT_PLUME_APPLY
89 C |-- SALT_PLUME_FORCING_SURF
90 C |
91 C |-- KPP_CALC
92 C |-- KPP_CALC_DUMMY
93 C |
94 C |-- PP81_CALC
95 C |
96 C |-- KL10_CALC
97 C |
98 C |-- MY82_CALC
99 C |
100 C |-- GGL90_CALC
101 C |
102 C |-- TIMEAVE_SURF_FLUX
103 C |
104 C |-- GMREDI_CALC_TENSOR
105 C |-- GMREDI_CALC_TENSOR_DUMMY
106 C |
107 C |-- DWNSLP_CALC_FLOW
108 C |-- DWNSLP_CALC_FLOW
109 C |
110 C |-- BBL_CALC_RHS
111 C |
112 C |-- MYPACKAGE_CALC_RHS
113 C |
114 C |-- GMREDI_DO_EXCH
115 C |
116 C |-- KPP_DO_EXCH
117 C |
118 C |-- DIAGS_RHO_G
119 C |-- DIAGS_OCEANIC_SURF_FLUX
120 C |-- SALT_PLUME_DIAGNOSTICS_FILL
121 C |
122 C |-- ECCO_PHYS
123
124 C !USES:
125 IMPLICIT NONE
126 C == Global variables ===
127 #include "SIZE.h"
128 #include "EEPARAMS.h"
129 #include "PARAMS.h"
130 #include "GRID.h"
131 #include "DYNVARS.h"
132 #ifdef ALLOW_TIMEAVE
133 # include "TIMEAVE_STATV.h"
134 #endif
135 #ifdef ALLOW_OFFLINE
136 # include "OFFLINE_SWITCH.h"
137 #endif
138
139 #ifdef ALLOW_AUTODIFF
140 # include "AUTODIFF_MYFIELDS.h"
141 # include "tamc.h"
142 # include "tamc_keys.h"
143 # include "FFIELDS.h"
144 # include "SURFACE.h"
145 # include "EOS.h"
146 # ifdef ALLOW_KPP
147 # include "KPP.h"
148 # endif
149 # ifdef ALLOW_GGL90
150 # include "GGL90.h"
151 # endif
152 # ifdef ALLOW_GMREDI
153 # include "GMREDI.h"
154 # endif
155 # ifdef ALLOW_EBM
156 # include "EBM.h"
157 # endif
158 # ifdef ALLOW_EXF
159 # include "ctrl.h"
160 # include "EXF_FIELDS.h"
161 # ifdef ALLOW_BULKFORMULAE
162 # include "EXF_CONSTANTS.h"
163 # endif
164 # endif
165 # ifdef ALLOW_SEAICE
166 # include "SEAICE_SIZE.h"
167 # include "SEAICE.h"
168 # include "SEAICE_PARAMS.h"
169 # endif
170 # ifdef ALLOW_THSICE
171 # include "THSICE_VARS.h"
172 # endif
173 # ifdef ALLOW_SALT_PLUME
174 # include "SALT_PLUME.h"
175 # endif
176 # ifdef ALLOW_ECCO
177 # ifdef ALLOW_SIGMAR_COST_CONTRIBUTION
178 # include "ecco_cost.h"
179 # endif
180 # endif
181 #endif /* ALLOW_AUTODIFF */
182
183 C !INPUT/OUTPUT PARAMETERS:
184 C == Routine arguments ==
185 C myTime :: Current time in simulation
186 C myIter :: Current iteration number in simulation
187 C myThid :: Thread number for this instance of the routine.
188 _RL myTime
189 INTEGER myIter
190 INTEGER myThid
191
192 C !LOCAL VARIABLES:
193 C == Local variables
194 C rhoK, rhoKm1 :: Density at current level, and level above
195 C iMin, iMax :: Ranges and sub-block indices on which calculations
196 C jMin, jMax are applied.
197 C bi, bj :: tile indices
198 C i,j,k :: loop indices
199 _RL rhoKp1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
200 _RL rhoKm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
201 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
202 _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
203 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
204 INTEGER iMin, iMax
205 INTEGER jMin, jMax
206 INTEGER bi, bj
207 INTEGER i, j, k
208 INTEGER doDiagsRho
209 LOGICAL calcGMRedi, calcKPP, calcConvect
210 #ifdef ALLOW_DIAGNOSTICS
211 LOGICAL DIAGNOSTICS_IS_ON
212 EXTERNAL DIAGNOSTICS_IS_ON
213 #endif /* ALLOW_DIAGNOSTICS */
214 #ifdef ALLOW_AUTODIFF
215 _RL thetaRef
216 #endif /* ALLOW_AUTODIFF */
217 CEOP
218
219 #ifdef ALLOW_AUTODIFF_TAMC
220 C-- dummy statement to end declaration part
221 itdkey = 1
222 #endif /* ALLOW_AUTODIFF_TAMC */
223
224 #ifdef ALLOW_DEBUG
225 IF (debugMode) CALL DEBUG_ENTER('DO_OCEANIC_PHYS',myThid)
226 #endif
227
228 doDiagsRho = 0
229 #ifdef ALLOW_DIAGNOSTICS
230 IF ( useDiagnostics .AND. fluidIsWater ) THEN
231 IF ( DIAGNOSTICS_IS_ON('MXLDEPTH',myThid) )
232 & doDiagsRho = doDiagsRho + 1
233 IF ( DIAGNOSTICS_IS_ON('DRHODR ',myThid) )
234 & doDiagsRho = doDiagsRho + 2
235 IF ( DIAGNOSTICS_IS_ON('WdRHO_P ',myThid) )
236 & doDiagsRho = doDiagsRho + 4
237 IF ( DIAGNOSTICS_IS_ON('WdRHOdP ',myThid) )
238 & doDiagsRho = doDiagsRho + 8
239 ENDIF
240 #endif /* ALLOW_DIAGNOSTICS */
241
242 calcGMRedi = useGMRedi
243 calcKPP = useKPP
244 calcConvect = ivdc_kappa.NE.0.
245 #ifdef ALLOW_OFFLINE
246 IF ( useOffLine ) THEN
247 calcGMRedi = useGMRedi .AND. .NOT.offlineLoadGMRedi
248 calcKPP = useKPP .AND. .NOT.offlineLoadKPP
249 calcConvect=calcConvect.AND. .NOT.offlineLoadConvec
250 ENDIF
251 #endif /* ALLOW_OFFLINE */
252
253 #ifdef ALLOW_OBCS
254 IF (useOBCS) THEN
255 C-- Calculate future values on open boundaries
256 C-- moved before SEAICE_MODEL call since SEAICE_MODEL needs seaice-obcs fields
257 # ifdef ALLOW_AUTODIFF_TAMC
258 CADJ STORE theta = comlev1, key=ikey_dynamics, kind=isbyte
259 CADJ STORE salt = comlev1, key=ikey_dynamics, kind=isbyte
260 # endif
261 # ifdef ALLOW_DEBUG
262 IF (debugMode) CALL DEBUG_CALL('OBCS_CALC',myThid)
263 # endif
264 CALL OBCS_CALC( myTime+deltaTClock, myIter+1,
265 I uVel, vVel, wVel, theta, salt, myThid )
266 ENDIF
267 #endif /* ALLOW_OBCS */
268
269 #ifdef ALLOW_AUTODIFF
270 # ifdef ALLOW_SALT_PLUME
271 DO bj=myByLo(myThid),myByHi(myThid)
272 DO bi=myBxLo(myThid),myBxHi(myThid)
273 DO j=1-OLy,sNy+OLy
274 DO i=1-OLx,sNx+OLx
275 saltPlumeDepth(i,j,bi,bj) = 0. _d 0
276 saltPlumeFlux(i,j,bi,bj) = 0. _d 0
277 ENDDO
278 ENDDO
279 ENDDO
280 ENDDO
281 # endif
282 # ifdef ALLOW_ECCO
283 # ifdef ALLOW_SIGMAR_COST_CONTRIBUTION
284 DO bj=myByLo(myThid),myByHi(myThid)
285 DO bi=myBxLo(myThid),myBxHi(myThid)
286 DO k=1,Nr
287 DO j=1-OLy,sNy+OLy
288 DO i=1-OLx,sNx+OLx
289 sigmaRfield(i,j,k,bi,bj) = 0. _d 0
290 ENDDO
291 ENDDO
292 ENDDO
293 ENDDO
294 ENDDO
295 # endif
296 # endif
297 #endif /* ALLOW_AUTODIFF */
298
299 #ifdef ALLOW_FRAZIL
300 IF ( useFRAZIL ) THEN
301 C-- Freeze water in the ocean interior and let it rise to the surface
302 CALL FRAZIL_CALC_RHS( myTime, myIter, myThid )
303 ENDIF
304 #endif /* ALLOW_FRAZIL */
305
306 #ifndef OLD_THSICE_CALL_SEQUENCE
307 #if (defined ALLOW_THSICE) && !(defined ALLOW_ATM2D)
308 IF ( useThSIce .AND. fluidIsWater ) THEN
309 # ifdef ALLOW_AUTODIFF_TAMC
310 CADJ STORE uice,vice = comlev1, key = ikey_dynamics,
311 CADJ & kind = isbyte
312 CADJ STORE iceMask,iceHeight = comlev1, key = ikey_dynamics,
313 CADJ & kind = isbyte
314 CADJ STORE snowHeight, Tsrf = comlev1, key = ikey_dynamics,
315 CADJ & kind = isbyte
316 CADJ STORE Qice1, Qice2 = comlev1, key = ikey_dynamics,
317 CADJ & kind = isbyte
318 CADJ STORE sHeating, snowAge = comlev1, key = ikey_dynamics,
319 CADJ & kind = isbyte
320 CADJ STORE hocemxl = comlev1, key = ikey_dynamics,
321 CADJ & kind = isbyte
322 CADJ STORE icflxsw = comlev1, key = ikey_dynamics,
323 CADJ & kind = isbyte
324 CADJ STORE salt,theta = comlev1, key = ikey_dynamics,
325 CADJ & kind = isbyte
326 CADJ STORE uvel,vvel = comlev1, key = ikey_dynamics,
327 CADJ & kind = isbyte
328 CADJ STORE qnet,qsw, empmr = comlev1, key = ikey_dynamics,
329 CADJ & kind = isbyte
330 CADJ STORE atemp,aqh,precip = comlev1, key = ikey_dynamics,
331 CADJ & kind = isbyte
332 CADJ STORE swdown,lwdown = comlev1, key = ikey_dynamics,
333 CADJ & kind = isbyte
334 # ifdef NONLIN_FRSURF
335 CADJ STORE hFac_surfC = comlev1, key = ikey_dynamics,
336 CADJ & kind = isbyte
337 # endif
338 # endif /* ALLOW_AUTODIFF_TAMC */
339 # ifdef ALLOW_DEBUG
340 IF (debugMode) CALL DEBUG_CALL('THSICE_MAIN',myThid)
341 # endif
342 C-- Step forward Therm.Sea-Ice variables
343 C and modify forcing terms including effects from ice
344 CALL TIMER_START('THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
345 CALL THSICE_MAIN( myTime, myIter, myThid )
346 CALL TIMER_STOP( 'THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
347 ENDIF
348 #endif /* ALLOW_THSICE */
349 #endif /* ndef OLD_THSICE_CALL_SEQUENCE */
350
351 #ifdef ALLOW_SEAICE
352 # ifdef ALLOW_AUTODIFF
353 CADJ STORE area = comlev1, key=ikey_dynamics, kind=isbyte
354 CADJ STORE fu,fv = comlev1, key=ikey_dynamics, kind=isbyte
355 CADJ STORE qnet = comlev1, key=ikey_dynamics, kind=isbyte
356 CADJ STORE qsw = comlev1, key=ikey_dynamics, kind=isbyte
357 CADJ STORE theta = comlev1, key=ikey_dynamics, kind=isbyte
358 CADJ STORE salt = comlev1, key=ikey_dynamics, kind=isbyte
359 #if (defined ALLOW_EXF) && (defined ALLOW_ATM_TEMP)
360 CADJ STORE evap = comlev1, key=ikey_dynamics, kind=isbyte
361 #endif
362 IF ( .NOT.useSEAICE .AND. SEAICEadjMODE .EQ. -1 ) THEN
363 CALL SEAICE_FAKE( myTime, myIter, myThid )
364 ENDIF
365 CADJ STORE area = comlev1, key=ikey_dynamics, kind=isbyte
366 CADJ STORE fu,fv = comlev1, key=ikey_dynamics, kind=isbyte
367 CADJ STORE qnet = comlev1, key=ikey_dynamics, kind=isbyte
368 CADJ STORE qsw = comlev1, key=ikey_dynamics, kind=isbyte
369 CADJ STORE theta = comlev1, key=ikey_dynamics, kind=isbyte
370 CADJ STORE salt = comlev1, key=ikey_dynamics, kind=isbyte
371 #if (defined ALLOW_EXF) && (defined ALLOW_ATM_TEMP)
372 CADJ STORE evap = comlev1, key=ikey_dynamics, kind=isbyte
373 #endif
374 # endif /* ALLOW_AUTODIFF */
375 #endif /* ALLOW_SEAICE */
376
377 #ifdef ALLOW_SEAICE
378 IF ( useSEAICE ) THEN
379 # ifdef ALLOW_AUTODIFF_TAMC
380 cph-adj-test(
381 CADJ STORE area = comlev1, key=ikey_dynamics, kind=isbyte
382 CADJ STORE hsnow = comlev1, key=ikey_dynamics, kind=isbyte
383 CADJ STORE heff = comlev1, key=ikey_dynamics, kind=isbyte
384 CADJ STORE tices = comlev1, key=ikey_dynamics, kind=isbyte
385 CADJ STORE empmr, qnet = comlev1, key=ikey_dynamics, kind=isbyte
386 CADJ STORE qsw,saltflux = comlev1, key=ikey_dynamics, kind=isbyte
387 CADJ STORE fu, fv = comlev1, key=ikey_dynamics, kind=isbyte
388 cCADJ STORE theta = comlev1, key=ikey_dynamics, kind=isbyte
389 cCADJ STORE salt = comlev1, key=ikey_dynamics, kind=isbyte
390 cph-adj-test)
391 c#ifdef ALLOW_EXF
392 CADJ STORE atemp,aqh,precip = comlev1, key = ikey_dynamics,
393 CADJ & kind = isbyte
394 CADJ STORE swdown,lwdown = comlev1, key = ikey_dynamics,
395 CADJ & kind = isbyte
396 CADJ STORE evap = comlev1, key = ikey_dynamics,
397 CADJ & kind = isbyte
398 CADJ STORE uwind,vwind = comlev1, key = ikey_dynamics,
399 CADJ & kind = isbyte
400 c#endif
401 CADJ STORE uvel,vvel = comlev1, key = ikey_dynamics,
402 CADJ & kind = isbyte
403 # ifdef SEAICE_CGRID
404 CADJ STORE stressdivergencex = comlev1, key = ikey_dynamics,
405 CADJ & kind = isbyte
406 CADJ STORE stressdivergencey = comlev1, key = ikey_dynamics,
407 CADJ & kind = isbyte
408 # endif
409 # ifdef SEAICE_ALLOW_DYNAMICS
410 CADJ STORE uice = comlev1, key = ikey_dynamics,
411 CADJ & kind = isbyte
412 CADJ STORE vice = comlev1, key = ikey_dynamics,
413 CADJ & kind = isbyte
414 CADJ STORE dwatn = comlev1, key = ikey_dynamics,
415 CADJ & kind = isbyte
416 # ifdef SEAICE_ALLOW_EVP
417 CADJ STORE seaice_sigma1 = comlev1, key = ikey_dynamics,
418 CADJ & kind = isbyte
419 CADJ STORE seaice_sigma2 = comlev1, key = ikey_dynamics,
420 CADJ & kind = isbyte
421 CADJ STORE seaice_sigma12 = comlev1, key = ikey_dynamics,
422 CADJ & kind = isbyte
423 # endif
424 # endif
425 # ifdef SEAICE_VARIABLE_SALINITY
426 CADJ STORE hsalt = comlev1, key = ikey_dynamics,
427 CADJ & kind = isbyte
428 # endif
429 # ifdef ATMOSPHERIC_LOADING
430 CADJ STORE pload = comlev1, key = ikey_dynamics,
431 CADJ & kind = isbyte
432 CADJ STORE siceload = comlev1, key = ikey_dynamics,
433 CADJ & kind = isbyte
434 # endif
435 # ifdef NONLIN_FRSURF
436 CADJ STORE recip_hfacc = comlev1, key = ikey_dynamics,
437 CADJ & kind = isbyte
438 # endif
439 # ifdef ANNUAL_BALANCE
440 CADJ STORE balance_itcount = comlev1, key = ikey_dynamics,
441 CADJ & kind = isbyte
442 # endif /* ANNUAL_BALANCE */
443 # ifdef ALLOW_THSICE
444 C-- store thSIce vars before advection (called from SEAICE_MODEL) update them:
445 CADJ STORE iceMask,iceHeight = comlev1, key = ikey_dynamics,
446 CADJ & kind = isbyte
447 CADJ STORE snowHeight, hOceMxL = comlev1, key = ikey_dynamics,
448 CADJ & kind = isbyte
449 CADJ STORE Qice1, Qice2 = comlev1, key = ikey_dynamics,
450 CADJ & kind = isbyte
451 # endif /* ALLOW_THSICE */
452 # endif /* ALLOW_AUTODIFF_TAMC */
453 # ifdef ALLOW_DEBUG
454 IF (debugMode) CALL DEBUG_CALL('SEAICE_MODEL',myThid)
455 # endif
456 CALL TIMER_START('SEAICE_MODEL [DO_OCEANIC_PHYS]', myThid)
457 CALL SEAICE_MODEL( myTime, myIter, myThid )
458 CALL TIMER_STOP ('SEAICE_MODEL [DO_OCEANIC_PHYS]', myThid)
459 # ifdef ALLOW_COST
460 CALL SEAICE_COST_SENSI ( myTime, myIter, myThid )
461 # endif
462 ENDIF
463 #endif /* ALLOW_SEAICE */
464
465 #ifdef ALLOW_AUTODIFF_TAMC
466 CADJ STORE sst, sss = comlev1, key = ikey_dynamics,
467 CADJ & kind = isbyte
468 CADJ STORE qsw = comlev1, key = ikey_dynamics,
469 CADJ & kind = isbyte
470 # ifdef ALLOW_SEAICE
471 CADJ STORE area = comlev1, key = ikey_dynamics,
472 CADJ & kind = isbyte
473 # endif
474 #endif
475
476 #ifdef OLD_THSICE_CALL_SEQUENCE
477 #if (defined ALLOW_THSICE) && !(defined ALLOW_ATM2D)
478 IF ( useThSIce .AND. fluidIsWater ) THEN
479 # ifdef ALLOW_AUTODIFF_TAMC
480 cph(
481 # ifdef NONLIN_FRSURF
482 CADJ STORE uice,vice = comlev1, key = ikey_dynamics,
483 CADJ & kind = isbyte
484 CADJ STORE salt,theta = comlev1, key = ikey_dynamics,
485 CADJ & kind = isbyte
486 CADJ STORE qnet,qsw, empmr = comlev1, key = ikey_dynamics,
487 CADJ & kind = isbyte
488 CADJ STORE hFac_surfC = comlev1, key = ikey_dynamics,
489 CADJ & kind = isbyte
490 # endif
491 # endif
492 # ifdef ALLOW_DEBUG
493 IF (debugMode) CALL DEBUG_CALL('THSICE_MAIN',myThid)
494 # endif
495 C-- Step forward Therm.Sea-Ice variables
496 C and modify forcing terms including effects from ice
497 CALL TIMER_START('THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
498 CALL THSICE_MAIN( myTime, myIter, myThid )
499 CALL TIMER_STOP( 'THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
500 ENDIF
501 #endif /* ALLOW_THSICE */
502 #endif /* OLD_THSICE_CALL_SEQUENCE */
503
504 #ifdef ALLOW_SHELFICE
505 # ifdef ALLOW_AUTODIFF_TAMC
506 CADJ STORE salt, theta = comlev1, key = ikey_dynamics,
507 CADJ & kind = isbyte
508 # endif
509 IF ( useShelfIce .AND. fluidIsWater ) THEN
510 #ifdef ALLOW_DEBUG
511 IF (debugMode) CALL DEBUG_CALL('SHELFICE_THERMODYNAMICS',myThid)
512 #endif
513 C compute temperature and (virtual) salt flux at the
514 C shelf-ice ocean interface
515 CALL TIMER_START('SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',
516 & myThid)
517 CALL SHELFICE_THERMODYNAMICS( myTime, myIter, myThid )
518 CALL TIMER_STOP( 'SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',
519 & myThid)
520 ENDIF
521 #endif /* ALLOW_SHELFICE */
522
523 #ifdef ALLOW_ICEFRONT
524 IF ( useICEFRONT .AND. fluidIsWater ) THEN
525 #ifdef ALLOW_DEBUG
526 IF (debugMode) CALL DEBUG_CALL('ICEFRONT_THERMODYNAMICS',myThid)
527 #endif
528 C compute temperature and (virtual) salt flux at the
529 C ice-front ocean interface
530 CALL TIMER_START('ICEFRONT_THERMODYNAMICS [DO_OCEANIC_PHYS]',
531 & myThid)
532 CALL ICEFRONT_THERMODYNAMICS( myTime, myIter, myThid )
533 CALL TIMER_STOP( 'ICEFRONT_THERMODYNAMICS [DO_OCEANIC_PHYS]',
534 & myThid)
535 ENDIF
536 #endif /* ALLOW_ICEFRONT */
537
538 #ifdef ALLOW_SALT_PLUME
539 IF ( useSALT_PLUME ) THEN
540 Catn: exchanging saltPlumeFlux:
541 CALL SALT_PLUME_DO_EXCH( myTime, myIter, myThid )
542 ENDIF
543 #endif /* ALLOW_SALT_PLUME */
544
545 C-- Freeze water at the surface
546 IF ( allowFreezing ) THEN
547 #ifdef ALLOW_AUTODIFF_TAMC
548 CADJ STORE theta = comlev1, key = ikey_dynamics,
549 CADJ & kind = isbyte
550 #endif
551 CALL FREEZE_SURFACE( myTime, myIter, myThid )
552 ENDIF
553
554 #ifdef ALLOW_OCN_COMPON_INTERF
555 C-- Apply imported data (from coupled interface) to forcing fields
556 C jmc: do not know precisely where to put this call (bf or af thSIce ?)
557 IF ( useCoupler ) THEN
558 CALL OCN_APPLY_IMPORT( .TRUE., myTime, myIter, myThid )
559 ENDIF
560 #endif /* ALLOW_OCN_COMPON_INTERF */
561
562 iMin = 1-OLx
563 iMax = sNx+OLx
564 jMin = 1-OLy
565 jMax = sNy+OLy
566
567 C--- Determines forcing terms based on external fields
568 C relaxation terms, etc.
569 #ifdef ALLOW_DEBUG
570 IF (debugMode) CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
571 #endif
572 #ifdef ALLOW_AUTODIFF
573 CADJ STORE salt, theta = comlev1, key = ikey_dynamics,
574 CADJ & kind = isbyte
575 #else /* ALLOW_AUTODIFF */
576 C-- if fluid is not water, by-pass surfaceForcing, find_rho, gmredi
577 C and all vertical mixing schemes, but keep OBCS_CALC
578 IF ( fluidIsWater ) THEN
579 #endif /* ALLOW_AUTODIFF */
580 CALL EXTERNAL_FORCING_SURF(
581 I iMin, iMax, jMin, jMax,
582 I myTime, myIter, myThid )
583
584 #ifdef ALLOW_AUTODIFF_TAMC
585 C-- HPF directive to help TAMC
586 CHPF$ INDEPENDENT
587 #endif /* ALLOW_AUTODIFF_TAMC */
588 DO bj=myByLo(myThid),myByHi(myThid)
589 #ifdef ALLOW_AUTODIFF_TAMC
590 C-- HPF directive to help TAMC
591 CHPF$ INDEPENDENT
592 #endif /* ALLOW_AUTODIFF_TAMC */
593 DO bi=myBxLo(myThid),myBxHi(myThid)
594
595 #ifdef ALLOW_AUTODIFF_TAMC
596 act1 = bi - myBxLo(myThid)
597 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
598 act2 = bj - myByLo(myThid)
599 max2 = myByHi(myThid) - myByLo(myThid) + 1
600 act3 = myThid - 1
601 max3 = nTx*nTy
602 act4 = ikey_dynamics - 1
603 itdkey = (act1 + 1) + act2*max1
604 & + act3*max1*max2
605 & + act4*max1*max2*max3
606 #endif /* ALLOW_AUTODIFF_TAMC */
607
608 C-- Set up work arrays with valid (i.e. not NaN) values
609 C These inital values do not alter the numerical results. They
610 C just ensure that all memory references are to valid floating
611 C point numbers. This prevents spurious hardware signals due to
612 C uninitialised but inert locations.
613 DO k=1,Nr
614 DO j=1-OLy,sNy+OLy
615 DO i=1-OLx,sNx+OLx
616 C This is currently used by GMRedi, IVDC, MXL-depth and Diagnostics
617 sigmaX(i,j,k) = 0. _d 0
618 sigmaY(i,j,k) = 0. _d 0
619 sigmaR(i,j,k) = 0. _d 0
620 ENDDO
621 ENDDO
622 ENDDO
623
624 #ifdef ALLOW_AUTODIFF
625 DO j=1-OLy,sNy+OLy
626 DO i=1-OLx,sNx+OLx
627 rhoKm1 (i,j) = 0. _d 0
628 rhoKp1 (i,j) = 0. _d 0
629 ENDDO
630 ENDDO
631 cph all the following init. are necessary for TAF
632 cph although some of these are re-initialised later.
633 DO k=1,Nr
634 DO j=1-OLy,sNy+OLy
635 DO i=1-OLx,sNx+OLx
636 rhoInSitu(i,j,k,bi,bj) = 0.
637 # ifdef ALLOW_GGL90
638 GGL90viscArU(i,j,k,bi,bj) = 0. _d 0
639 GGL90viscArV(i,j,k,bi,bj) = 0. _d 0
640 GGL90diffKr(i,j,k,bi,bj) = 0. _d 0
641 # endif /* ALLOW_GGL90 */
642 # ifdef ALLOW_SALT_PLUME
643 # ifdef SALT_PLUME_VOLUME
644 SPforcingS(i,j,k,bi,bj) = 0. _d 0
645 SPforcingT(i,j,k,bi,bj) = 0. _d 0
646 # endif
647 # endif /* ALLOW_SALT_PLUME */
648 ENDDO
649 ENDDO
650 ENDDO
651 #ifdef ALLOW_OFFLINE
652 IF ( calcConvect ) THEN
653 #endif
654 DO k=1,Nr
655 DO j=1-OLy,sNy+OLy
656 DO i=1-OLx,sNx+OLx
657 IVDConvCount(i,j,k,bi,bj) = 0.
658 ENDDO
659 ENDDO
660 ENDDO
661 #ifdef ALLOW_OFFLINE
662 ENDIF
663 IF ( calcGMRedi ) THEN
664 #endif
665 # ifdef ALLOW_GMREDI
666 DO k=1,Nr
667 DO j=1-OLy,sNy+OLy
668 DO i=1-OLx,sNx+OLx
669 Kwx(i,j,k,bi,bj) = 0. _d 0
670 Kwy(i,j,k,bi,bj) = 0. _d 0
671 Kwz(i,j,k,bi,bj) = 0. _d 0
672 # ifdef GM_NON_UNITY_DIAGONAL
673 Kux(i,j,k,bi,bj) = 0. _d 0
674 Kvy(i,j,k,bi,bj) = 0. _d 0
675 # endif
676 # ifdef GM_EXTRA_DIAGONAL
677 Kuz(i,j,k,bi,bj) = 0. _d 0
678 Kvz(i,j,k,bi,bj) = 0. _d 0
679 # endif
680 # ifdef GM_BOLUS_ADVEC
681 GM_PsiX(i,j,k,bi,bj) = 0. _d 0
682 GM_PsiY(i,j,k,bi,bj) = 0. _d 0
683 # endif
684 # ifdef GM_VISBECK_VARIABLE_K
685 VisbeckK(i,j,bi,bj) = 0. _d 0
686 # endif
687 ENDDO
688 ENDDO
689 ENDDO
690 # endif /* ALLOW_GMREDI */
691 #ifdef ALLOW_OFFLINE
692 ENDIF
693 IF ( calcKPP ) THEN
694 #endif
695 # ifdef ALLOW_KPP
696 DO k=1,Nr
697 DO j=1-OLy,sNy+OLy
698 DO i=1-OLx,sNx+OLx
699 KPPdiffKzS(i,j,k,bi,bj) = 0. _d 0
700 KPPdiffKzT(i,j,k,bi,bj) = 0. _d 0
701 ENDDO
702 ENDDO
703 ENDDO
704 # endif /* ALLOW_KPP */
705 #ifdef ALLOW_OFFLINE
706 ENDIF
707 #endif
708 #endif /* ALLOW_AUTODIFF */
709
710 #ifdef ALLOW_AUTODIFF_TAMC
711 CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey,
712 CADJ & kind = isbyte
713 CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey,
714 CADJ & kind = isbyte
715 CADJ STORE totphihyd(:,:,:,bi,bj)
716 CADJ & = comlev1_bibj, key=itdkey,
717 CADJ & kind = isbyte
718 # ifdef ALLOW_KPP
719 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey,
720 CADJ & kind = isbyte
721 CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey,
722 CADJ & kind = isbyte
723 # endif
724 # ifdef ALLOW_SALT_PLUME
725 CADJ STORE saltplumedepth(:,:,bi,bj) = comlev1_bibj, key=itdkey,
726 CADJ & kind = isbyte
727 CADJ STORE saltplumeflux(:,:,bi,bj) = comlev1_bibj, key=itdkey,
728 CADJ & kind = isbyte
729 # endif
730 #endif /* ALLOW_AUTODIFF_TAMC */
731
732 C-- Always compute density (stored in common block) here; even when it is not
733 C needed here, will be used anyway in calc_phi_hyd (data flow easier this way)
734 #ifdef ALLOW_DEBUG
735 IF (debugMode) CALL DEBUG_CALL('FIND_RHO_2D (xNr)',myThid)
736 #endif
737 #ifdef ALLOW_AUTODIFF
738 IF ( fluidIsWater ) THEN
739 #endif /* ALLOW_AUTODIFF */
740 #ifdef ALLOW_DOWN_SLOPE
741 IF ( useDOWN_SLOPE ) THEN
742 DO k=1,Nr
743 CALL DWNSLP_CALC_RHO(
744 I theta, salt,
745 O rhoInSitu(1-OLx,1-OLy,k,bi,bj),
746 I k, bi, bj, myTime, myIter, myThid )
747 ENDDO
748 ENDIF
749 #endif /* ALLOW_DOWN_SLOPE */
750 #ifdef ALLOW_BBL
751 IF ( useBBL ) THEN
752 C pkg/bbl requires in-situ bbl density for depths equal to and deeper than the bbl.
753 C To reduce computation and storage requirement, these densities are stored in the
754 C dry grid boxes of rhoInSitu. See BBL_CALC_RHO for details.
755 DO k=Nr,1,-1
756 CALL BBL_CALC_RHO(
757 I theta, salt,
758 O rhoInSitu,
759 I k, bi, bj, myTime, myIter, myThid )
760
761 ENDDO
762 ENDIF
763 #endif /* ALLOW_BBL */
764 IF ( .NOT. ( useDOWN_SLOPE .OR. useBBL ) ) THEN
765 DO k=1,Nr
766 CALL FIND_RHO_2D(
767 I iMin, iMax, jMin, jMax, k,
768 I theta(1-OLx,1-OLy,k,bi,bj),
769 I salt (1-OLx,1-OLy,k,bi,bj),
770 O rhoInSitu(1-OLx,1-OLy,k,bi,bj),
771 I k, bi, bj, myThid )
772 ENDDO
773 ENDIF
774 #ifdef ALLOW_AUTODIFF
775 ELSE
776 C- fluid is not water:
777 DO k=1,Nr
778 IF ( select_rStar.GE.1 .OR. selectSigmaCoord.GE.1 ) THEN
779 C- isothermal (theta=const) reference state
780 thetaRef = thetaConst
781 ELSE
782 C- horizontally uniform (tRef) reference state
783 thetaRef = tRef(k)
784 ENDIF
785 DO j=1-OLy,sNy+OLy
786 DO i=1-OLx,sNx+OLx
787 rhoInSitu(i,j,k,bi,bj) =
788 & ( theta(i,j,k,bi,bj)
789 & *( salt(i,j,k,bi,bj)*atm_Rq + oneRL )
790 & - thetaRef )*maskC(i,j,k,bi,bj)
791 ENDDO
792 ENDDO
793 ENDDO
794 ENDIF
795 #endif /* ALLOW_AUTODIFF */
796
797 #ifdef ALLOW_DEBUG
798 IF (debugMode) CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
799 #endif
800
801 C-- Start of diagnostic loop
802 DO k=Nr,1,-1
803
804 #ifdef ALLOW_AUTODIFF_TAMC
805 C? Patrick, is this formula correct now that we change the loop range?
806 C? Do we still need this?
807 cph kkey formula corrected.
808 cph Needed for rhoK, rhoKm1, in the case useGMREDI.
809 kkey = (itdkey-1)*Nr + k
810 #endif /* ALLOW_AUTODIFF_TAMC */
811
812 c#ifdef ALLOW_AUTODIFF_TAMC
813 cCADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey,
814 cCADJ & kind = isbyte
815 cCADJ STORE salt(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey,
816 cCADJ & kind = isbyte
817 c#endif /* ALLOW_AUTODIFF_TAMC */
818
819 C-- Calculate gradients of potential density for isoneutral
820 C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
821 IF ( calcGMRedi .OR. (k.GT.1 .AND. calcConvect)
822 & .OR. usePP81 .OR. useKL10
823 & .OR. useMY82 .OR. useGGL90
824 & .OR. useSALT_PLUME .OR. doDiagsRho.GE.1 ) THEN
825 IF (k.GT.1) THEN
826 #ifdef ALLOW_AUTODIFF_TAMC
827 CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey,
828 CADJ & kind = isbyte
829 CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey,
830 CADJ & kind = isbyte
831 CADJ STORE rhokm1 (bi,bj) = comlev1_bibj_k, key=kkey,
832 CADJ & kind = isbyte
833 #endif /* ALLOW_AUTODIFF_TAMC */
834 CALL FIND_RHO_2D(
835 I iMin, iMax, jMin, jMax, k,
836 I theta(1-OLx,1-OLy,k-1,bi,bj),
837 I salt (1-OLx,1-OLy,k-1,bi,bj),
838 O rhoKm1,
839 I k-1, bi, bj, myThid )
840 ENDIF
841 #ifdef ALLOW_DEBUG
842 IF (debugMode) CALL DEBUG_CALL('GRAD_SIGMA',myThid)
843 #endif
844 cph Avoid variable aliasing for adjoint !!!
845 DO j=jMin,jMax
846 DO i=iMin,iMax
847 rhoKp1(i,j) = rhoInSitu(i,j,k,bi,bj)
848 ENDDO
849 ENDDO
850 CALL GRAD_SIGMA(
851 I bi, bj, iMin, iMax, jMin, jMax, k,
852 I rhoInSitu(1-OLx,1-OLy,k,bi,bj), rhoKm1, rhoKp1,
853 O sigmaX, sigmaY, sigmaR,
854 I myThid )
855 #ifdef ALLOW_ECCO
856 # ifdef ALLOW_SIGMAR_COST_CONTRIBUTION
857 DO j=jMin,jMax
858 DO i=iMin,iMax
859 sigmaRfield(i,j,k,bi,bj)=sigmaR(i,j,k)
860 ENDDO
861 ENDDO
862 # endif
863 #endif /* ALLOW_ECCO */
864 #ifdef ALLOW_AUTODIFF
865 #ifdef GMREDI_WITH_STABLE_ADJOINT
866 cgf zero out adjoint fields to stabilize pkg/gmredi adjoint
867 cgf -> cuts adjoint dependency from slope to state
868 CALL ZERO_ADJ_LOC( Nr, sigmaX, myThid)
869 CALL ZERO_ADJ_LOC( Nr, sigmaY, myThid)
870 CALL ZERO_ADJ_LOC( Nr, sigmaR, myThid)
871 #endif
872 #endif /* ALLOW_AUTODIFF */
873 ENDIF
874
875 C-- Implicit Vertical Diffusion for Convection
876 IF (k.GT.1 .AND. calcConvect) THEN
877 #ifdef ALLOW_DEBUG
878 IF (debugMode) CALL DEBUG_CALL('CALC_IVDC',myThid)
879 #endif
880 CALL CALC_IVDC(
881 I bi, bj, iMin, iMax, jMin, jMax, k,
882 I sigmaR,
883 I myTime, myIter, myThid)
884 ENDIF
885
886 #ifdef ALLOW_DIAGNOSTICS
887 IF ( doDiagsRho.GE.4 ) THEN
888 CALL DIAGS_RHO_L( doDiagsRho, k, bi, bj,
889 I rhoInSitu(1-OLx,1-OLy,1,bi,bj),
890 I rhoKm1, wVel,
891 I myTime, myIter, myThid )
892 ENDIF
893 #endif
894
895 C-- end of diagnostic k loop (Nr:1)
896 ENDDO
897
898 #ifdef ALLOW_AUTODIFF_TAMC
899 CADJ STORE IVDConvCount(:,:,:,bi,bj)
900 CADJ & = comlev1_bibj, key=itdkey,
901 CADJ & kind = isbyte
902 #endif
903
904 C-- Diagnose Mixed Layer Depth:
905 IF ( calcGMRedi .OR. MOD(doDiagsRho,2).EQ.1 ) THEN
906 CALL CALC_OCE_MXLAYER(
907 I rhoInSitu(1-OLx,1-OLy,1,bi,bj), sigmaR,
908 I bi, bj, myTime, myIter, myThid )
909 ENDIF
910
911 #ifdef ALLOW_SALT_PLUME
912 IF ( useSALT_PLUME ) THEN
913 CALL SALT_PLUME_CALC_DEPTH(
914 I rhoInSitu(1-OLx,1-OLy,1,bi,bj), sigmaR,
915 I bi, bj, myTime, myIter, myThid )
916 #ifdef SALT_PLUME_VOLUME
917 CALL SALT_PLUME_VOLFRAC(
918 I bi, bj, myTime, myIter, myThid )
919 C-- get forcings for kpp
920 CALL SALT_PLUME_APPLY(
921 I 1, bi, bj, recip_hFacC(1-OLx,1-OLy,1,bi,bj),
922 I theta, 0,
923 I myTime, myIter, myThid )
924 CALL SALT_PLUME_APPLY(
925 I 2, bi, bj, recip_hFacC(1-OLx,1-OLy,1,bi,bj),
926 I salt, 0,
927 I myTime, myIter, myThid )
928 C-- need to call this S/R from here to apply just before kpp
929 CALL SALT_PLUME_FORCING_SURF(
930 I bi, bj, iMin, iMax, jMin, jMax,
931 I myTime, myIter, myThid )
932 #endif /* SALT_PLUME_VOLUME */
933 ENDIF
934 #endif /* ALLOW_SALT_PLUME */
935
936 #ifdef ALLOW_DIAGNOSTICS
937 IF ( MOD(doDiagsRho,4).GE.2 ) THEN
938 CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR ', 0, Nr,
939 & 2, bi, bj, myThid)
940 ENDIF
941 #endif /* ALLOW_DIAGNOSTICS */
942
943 C-- This is where EXTERNAL_FORCING_SURF(bi,bj) used to be called;
944 C now called earlier, before bi,bj loop.
945
946 #ifdef ALLOW_AUTODIFF_TAMC
947 cph needed for KPP
948 CADJ STORE surfaceForcingU(:,:,bi,bj)
949 CADJ & = comlev1_bibj, key=itdkey,
950 CADJ & kind = isbyte
951 CADJ STORE surfaceForcingV(:,:,bi,bj)
952 CADJ & = comlev1_bibj, key=itdkey,
953 CADJ & kind = isbyte
954 CADJ STORE surfaceForcingS(:,:,bi,bj)
955 CADJ & = comlev1_bibj, key=itdkey,
956 CADJ & kind = isbyte
957 CADJ STORE surfaceForcingT(:,:,bi,bj)
958 CADJ & = comlev1_bibj, key=itdkey,
959 CADJ & kind = isbyte
960 CADJ STORE surfaceForcingTice(:,:,bi,bj)
961 CADJ & = comlev1_bibj, key=itdkey,
962 CADJ & kind = isbyte
963 #endif /* ALLOW_AUTODIFF_TAMC */
964
965 #ifdef ALLOW_KPP
966 C-- Compute KPP mixing coefficients
967 IF ( calcKPP ) THEN
968 #ifdef ALLOW_DEBUG
969 IF (debugMode) CALL DEBUG_CALL('KPP_CALC',myThid)
970 #endif
971 CALL TIMER_START('KPP_CALC [DO_OCEANIC_PHYS]', myThid)
972 CALL KPP_CALC(
973 I bi, bj, myTime, myIter, myThid )
974 CALL TIMER_STOP ('KPP_CALC [DO_OCEANIC_PHYS]', myThid)
975 #if (defined ALLOW_AUTODIFF) && !(defined ALLOW_OFFLINE)
976 ELSE
977 CALL KPP_CALC_DUMMY(
978 I bi, bj, myTime, myIter, myThid )
979 #endif /* ALLOW_AUTODIFF and not ALLOW_OFFLINE */
980 ENDIF
981 #endif /* ALLOW_KPP */
982
983 #ifdef ALLOW_PP81
984 C-- Compute PP81 mixing coefficients
985 IF (usePP81) THEN
986 #ifdef ALLOW_DEBUG
987 IF (debugMode) CALL DEBUG_CALL('PP81_CALC',myThid)
988 #endif
989 CALL PP81_CALC(
990 I bi, bj, sigmaR, myTime, myIter, myThid )
991 ENDIF
992 #endif /* ALLOW_PP81 */
993
994 #ifdef ALLOW_KL10
995 C-- Compute KL10 mixing coefficients
996 IF (useKL10) THEN
997 #ifdef ALLOW_DEBUG
998 IF (debugMode) CALL DEBUG_CALL('KL10_CALC',myThid)
999 #endif
1000 CALL KL10_CALC(
1001 I bi, bj, sigmaR, myTime, myIter, myThid )
1002 ENDIF
1003 #endif /* ALLOW_KL10 */
1004
1005 #ifdef ALLOW_MY82
1006 C-- Compute MY82 mixing coefficients
1007 IF (useMY82) THEN
1008 #ifdef ALLOW_DEBUG
1009 IF (debugMode) CALL DEBUG_CALL('MY82_CALC',myThid)
1010 #endif
1011 CALL MY82_CALC(
1012 I bi, bj, sigmaR, myTime, myIter, myThid )
1013 ENDIF
1014 #endif /* ALLOW_MY82 */
1015
1016 #ifdef ALLOW_GGL90
1017 #ifdef ALLOW_AUTODIFF_TAMC
1018 CADJ STORE GGL90TKE (:,:,:,bi,bj) = comlev1_bibj, key=itdkey,
1019 CADJ & kind = isbyte
1020 #endif /* ALLOW_AUTODIFF_TAMC */
1021 C-- Compute GGL90 mixing coefficients
1022 IF (useGGL90) THEN
1023 #ifdef ALLOW_DEBUG
1024 IF (debugMode) CALL DEBUG_CALL('GGL90_CALC',myThid)
1025 #endif
1026 CALL TIMER_START('GGL90_CALC [DO_OCEANIC_PHYS]', myThid)
1027 CALL GGL90_CALC(
1028 I bi, bj, sigmaR, myTime, myIter, myThid )
1029 CALL TIMER_STOP ('GGL90_CALC [DO_OCEANIC_PHYS]', myThid)
1030 ENDIF
1031 #endif /* ALLOW_GGL90 */
1032
1033 #ifdef ALLOW_TIMEAVE
1034 IF ( taveFreq.GT. 0. _d 0 ) THEN
1035 CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid)
1036 ENDIF
1037 IF ( taveFreq.GT.0. .AND. calcConvect ) THEN
1038 CALL TIMEAVE_CUMULATE(ConvectCountTave, IVDConvCount,
1039 I Nr, deltaTClock, bi, bj, myThid)
1040 ENDIF
1041 #endif /* ALLOW_TIMEAVE */
1042
1043 #ifdef ALLOW_GMREDI
1044 #ifdef ALLOW_AUTODIFF_TAMC
1045 # ifndef GM_EXCLUDE_CLIPPING
1046 cph storing here is needed only for one GMREDI_OPTIONS:
1047 cph define GM_BOLUS_ADVEC
1048 cph keep it although TAF says you dont need to.
1049 cph but I have avoided the #ifdef for now, in case more things change
1050 CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey,
1051 CADJ & kind = isbyte
1052 CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey,
1053 CADJ & kind = isbyte
1054 CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey,
1055 CADJ & kind = isbyte
1056 # endif
1057 #endif /* ALLOW_AUTODIFF_TAMC */
1058
1059 C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
1060 IF ( calcGMRedi ) THEN
1061 #ifdef ALLOW_DEBUG
1062 IF (debugMode) CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
1063 #endif
1064 CALL GMREDI_CALC_TENSOR(
1065 I iMin, iMax, jMin, jMax,
1066 I sigmaX, sigmaY, sigmaR,
1067 I bi, bj, myTime, myIter, myThid )
1068 #if (defined ALLOW_AUTODIFF) && !(defined ALLOW_OFFLINE)
1069 ELSE
1070 CALL GMREDI_CALC_TENSOR_DUMMY(
1071 I iMin, iMax, jMin, jMax,
1072 I sigmaX, sigmaY, sigmaR,
1073 I bi, bj, myTime, myIter, myThid )
1074 #endif /* ALLOW_AUTODIFF and not ALLOW_OFFLINE */
1075 ENDIF
1076 #endif /* ALLOW_GMREDI */
1077
1078 #ifdef ALLOW_DOWN_SLOPE
1079 IF ( useDOWN_SLOPE ) THEN
1080 C-- Calculate Downsloping Flow for Down_Slope parameterization
1081 IF ( usingPCoords ) THEN
1082 CALL DWNSLP_CALC_FLOW(
1083 I bi, bj, kSurfC, rhoInSitu,
1084 I myTime, myIter, myThid )
1085 ELSE
1086 CALL DWNSLP_CALC_FLOW(
1087 I bi, bj, kLowC, rhoInSitu,
1088 I myTime, myIter, myThid )
1089 ENDIF
1090 ENDIF
1091 #endif /* ALLOW_DOWN_SLOPE */
1092
1093 C-- end bi,bj loops.
1094 ENDDO
1095 ENDDO
1096
1097 #ifndef ALLOW_AUTODIFF
1098 C--- if fluid Is Water: end
1099 ENDIF
1100 #endif
1101
1102 #ifdef ALLOW_BBL
1103 IF ( useBBL ) THEN
1104 CALL BBL_CALC_RHS(
1105 I myTime, myIter, myThid )
1106 ENDIF
1107 #endif /* ALLOW_BBL */
1108
1109 #ifdef ALLOW_MYPACKAGE
1110 IF ( useMYPACKAGE ) THEN
1111 CALL MYPACKAGE_CALC_RHS(
1112 I myTime, myIter, myThid )
1113 ENDIF
1114 #endif /* ALLOW_MYPACKAGE */
1115
1116 #ifdef ALLOW_GMREDI
1117 IF ( calcGMRedi ) THEN
1118 CALL GMREDI_DO_EXCH( myTime, myIter, myThid )
1119 ENDIF
1120 #endif /* ALLOW_GMREDI */
1121
1122 #ifdef ALLOW_KPP
1123 IF ( calcKPP ) THEN
1124 CALL KPP_DO_EXCH( myThid )
1125 ENDIF
1126 #endif /* ALLOW_KPP */
1127
1128 #ifdef ALLOW_DIAGNOSTICS
1129 IF ( fluidIsWater .AND. useDiagnostics ) THEN
1130 CALL DIAGS_RHO_G(
1131 I rhoInSitu, uVel, vVel, wVel,
1132 I myTime, myIter, myThid )
1133 ENDIF
1134 IF ( useDiagnostics ) THEN
1135 CALL DIAGS_OCEANIC_SURF_FLUX( myTime, myIter, myThid )
1136 ENDIF
1137 IF ( calcConvect .AND. useDiagnostics ) THEN
1138 CALL DIAGNOSTICS_FILL( IVDConvCount, 'CONVADJ ',
1139 & 0, Nr, 0, 1, 1, myThid )
1140 ENDIF
1141 #ifdef ALLOW_SALT_PLUME
1142 IF ( useDiagnostics )
1143 & CALL SALT_PLUME_DIAGNOSTICS_FILL(bi,bj,myThid)
1144 #endif
1145 #endif
1146
1147 #ifdef ALLOW_ECCO
1148 IF ( useECCO ) CALL ECCO_PHYS( myThid )
1149 #endif
1150
1151 #ifdef ALLOW_DEBUG
1152 IF (debugMode) CALL DEBUG_LEAVE('DO_OCEANIC_PHYS',myThid)
1153 #endif
1154
1155 RETURN
1156 END

  ViewVC Help
Powered by ViewVC 1.1.22