/[MITgcm]/MITgcm_contrib/atnguyen/code_21Dec2012_saltplume/do_oceanic_phys.F
ViewVC logotype

Contents of /MITgcm_contrib/atnguyen/code_21Dec2012_saltplume/do_oceanic_phys.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.6 - (show annotations) (download)
Sun May 4 09:26:06 2014 UTC (11 years, 3 months ago) by atn
Branch: MAIN
CVS Tags: HEAD
Changes since 1.5: +15 -6 lines
minor fix. add initialization.

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

  ViewVC Help
Powered by ViewVC 1.1.22