/[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.4 - (hide annotations) (download)
Thu May 1 21:30:47 2014 UTC (11 years, 3 months ago) by atn
Branch: MAIN
Changes since 1.3: +4 -5 lines
bug fixing

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

  ViewVC Help
Powered by ViewVC 1.1.22