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

Annotation 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.5 - (hide annotations) (download)
Fri May 2 05:35:22 2014 UTC (11 years, 3 months ago) by atn
Branch: MAIN
Changes since 1.4: +10 -6 lines
bug fix. move #include back inside autodiff bracket.

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

  ViewVC Help
Powered by ViewVC 1.1.22