/[MITgcm]/MITgcm/model/src/do_oceanic_phys.F
ViewVC logotype

Annotation of /MITgcm/model/src/do_oceanic_phys.F

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


Revision 1.69 - (hide annotations) (download)
Sun Aug 17 02:08:24 2008 UTC (15 years, 9 months ago) by jmc
Branch: MAIN
Changes since 1.68: +47 -19 lines
new pkg "down_slope".

1 jmc 1.69 C $Header: /u/gcmpack/MITgcm/model/src/do_oceanic_phys.F,v 1.68 2008/08/11 22:25:52 jmc Exp $
2 jmc 1.1 C $Name: $
3    
4     #include "PACKAGES_CONFIG.h"
5     #include "CPP_OPTIONS.h"
6    
7     #ifdef ALLOW_AUTODIFF_TAMC
8     # ifdef ALLOW_GMREDI
9     # include "GMREDI_OPTIONS.h"
10     # endif
11     # ifdef ALLOW_KPP
12     # include "KPP_OPTIONS.h"
13     # endif
14 jmc 1.29 # ifdef ALLOW_SEAICE
15     # include "SEAICE_OPTIONS.h"
16     # endif
17 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
18    
19     CBOP
20     C !ROUTINE: DO_OCEANIC_PHYS
21     C !INTERFACE:
22     SUBROUTINE DO_OCEANIC_PHYS(myTime, myIter, myThid)
23     C !DESCRIPTION: \bv
24     C *==========================================================*
25 jmc 1.28 C | SUBROUTINE DO_OCEANIC_PHYS
26     C | o Controlling routine for oceanic physics and
27 jmc 1.1 C | parameterization
28     C *==========================================================*
29     C | o originally, part of S/R thermodynamics
30     C *==========================================================*
31     C \ev
32    
33     C !USES:
34     IMPLICIT NONE
35     C == Global variables ===
36     #include "SIZE.h"
37     #include "EEPARAMS.h"
38     #include "PARAMS.h"
39 jmc 1.69 #include "GRID.h"
40 jmc 1.1 #include "DYNVARS.h"
41 jmc 1.20 #ifdef ALLOW_TIMEAVE
42     #include "TIMEAVE_STATV.h"
43     #endif
44 mlosch 1.22 #if defined (ALLOW_BALANCE_FLUXES) && !(defined ALLOW_AUTODIFF_TAMC)
45     #include "FFIELDS.h"
46     #endif
47 jmc 1.1
48     #ifdef ALLOW_AUTODIFF_TAMC
49     # include "tamc.h"
50     # include "tamc_keys.h"
51     # include "FFIELDS.h"
52 heimbach 1.54 # include "SURFACE.h"
53 jmc 1.1 # include "EOS.h"
54     # ifdef ALLOW_KPP
55     # include "KPP.h"
56     # endif
57     # ifdef ALLOW_GMREDI
58     # include "GMREDI.h"
59     # endif
60     # ifdef ALLOW_EBM
61     # include "EBM.h"
62     # endif
63 jmc 1.29 # ifdef ALLOW_EXF
64     # include "ctrl.h"
65 jmc 1.40 # include "EXF_FIELDS.h"
66 jmc 1.29 # ifdef ALLOW_BULKFORMULAE
67 jmc 1.40 # include "EXF_CONSTANTS.h"
68 jmc 1.29 # endif
69     # endif
70     # ifdef ALLOW_SEAICE
71     # include "SEAICE.h"
72     # endif
73 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
74    
75     C !INPUT/OUTPUT PARAMETERS:
76     C == Routine arguments ==
77 jmc 1.18 C myTime :: Current time in simulation
78     C myIter :: Current iteration number in simulation
79     C myThid :: Thread number for this instance of the routine.
80 jmc 1.1 _RL myTime
81     INTEGER myIter
82     INTEGER myThid
83    
84     C !LOCAL VARIABLES:
85     C == Local variables
86 jmc 1.47 C rhoK, rhoKm1 :: Density at current level, and level above
87 jmc 1.18 C iMin, iMax :: Ranges and sub-block indices on which calculations
88 jmc 1.1 C jMin, jMax are applied.
89 jmc 1.18 C bi, bj :: tile indices
90     C i,j,k :: loop indices
91 jmc 1.47 _RL rhoKp1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92     _RL rhoKm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93     _RL rhoK (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94 jmc 1.1 _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
95     _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
96     _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
97     INTEGER iMin, iMax
98     INTEGER jMin, jMax
99     INTEGER bi, bj
100 jmc 1.18 INTEGER i, j, k
101 jmc 1.17 INTEGER doDiagsRho
102     #ifdef ALLOW_DIAGNOSTICS
103     LOGICAL DIAGNOSTICS_IS_ON
104     EXTERNAL DIAGNOSTICS_IS_ON
105     #endif /* ALLOW_DIAGNOSTICS */
106 jmc 1.69 #ifdef ALLOW_DOWN_SLOPE
107     _RL rhoInSitu(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
108     #endif /* ALLOW_DOWN_SLOPE */
109 jmc 1.1
110     CEOP
111    
112 heimbach 1.12 #ifdef ALLOW_AUTODIFF_TAMC
113     C-- dummy statement to end declaration part
114     itdkey = 1
115     #endif /* ALLOW_AUTODIFF_TAMC */
116    
117 jmc 1.1 #ifdef ALLOW_DEBUG
118 jmc 1.29 IF ( debugLevel .GE. debLevB )
119     & CALL DEBUG_ENTER('DO_OCEANIC_PHYS',myThid)
120 jmc 1.1 #endif
121 jmc 1.36
122 jmc 1.17 doDiagsRho = 0
123     #ifdef ALLOW_DIAGNOSTICS
124     IF ( useDiagnostics .AND. fluidIsWater ) THEN
125     IF ( DIAGNOSTICS_IS_ON('RHOANOSQ',myThid) .OR.
126     & DIAGNOSTICS_IS_ON('URHOMASS',myThid) .OR.
127     & DIAGNOSTICS_IS_ON('VRHOMASS',myThid) .OR.
128     & DIAGNOSTICS_IS_ON('WRHOMASS',myThid) .OR.
129     & DIAGNOSTICS_IS_ON('WRHOMASS',myThid) ) doDiagsRho = 2
130 jmc 1.47 IF ( doDiagsRho.EQ.0 .AND.
131     & DIAGNOSTICS_IS_ON('MXLDEPTH',myThid) ) doDiagsRho = 1
132     IF ( doDiagsRho.EQ.0 .AND.
133     & DIAGNOSTICS_IS_ON('DRHODR ',myThid) ) doDiagsRho = 1
134 jmc 1.17 ENDIF
135     #endif /* ALLOW_DIAGNOSTICS */
136    
137 jmc 1.69
138 jmc 1.29 #ifdef ALLOW_SEAICE
139     IF ( useSEAICE ) THEN
140 heimbach 1.62 # ifdef ALLOW_AUTODIFF_TAMC
141 heimbach 1.65 cph-adj-test(
142     CADJ STORE area,empmr,qsw,theta = comlev1, key = ikey_dynamics
143     cph-adj-test)
144 heimbach 1.34 CADJ STORE atemp,aqh,precip = comlev1, key = ikey_dynamics
145     CADJ STORE swdown,lwdown = comlev1, key = ikey_dynamics
146     cph# ifdef EXF_READ_EVAP
147 heimbach 1.33 CADJ STORE evap = comlev1, key = ikey_dynamics
148 heimbach 1.34 cph# endif
149 jmc 1.29 CADJ STORE uvel,vvel = comlev1, key = ikey_dynamics
150 heimbach 1.62 # ifdef SEAICE_ALLOW_DYNAMICS
151     CADJ STORE uice = comlev1, key = ikey_dynamics
152     CADJ STORE vice = comlev1, key = ikey_dynamics
153     # ifdef SEAICE_ALLOW_EVP
154 heimbach 1.50 CADJ STORE seaice_sigma1 = comlev1, key = ikey_dynamics
155     CADJ STORE seaice_sigma2 = comlev1, key = ikey_dynamics
156     CADJ STORE seaice_sigma12 = comlev1, key = ikey_dynamics
157 heimbach 1.62 # endif
158     # endif
159     # ifdef SEAICE_SALINITY
160 heimbach 1.56 CADJ STORE salt = comlev1, key = ikey_dynamics
161 heimbach 1.62 # endif
162     # ifdef ATMOSPHERIC_LOADING
163 heimbach 1.64 CADJ STORE pload = comlev1, key = ikey_dynamics
164 heimbach 1.37 CADJ STORE siceload = comlev1, key = ikey_dynamics
165 heimbach 1.62 # endif
166     # ifdef NONLIN_FRSURF
167 heimbach 1.55 CADJ STORE recip_hfacc = comlev1, key = ikey_dynamics
168 heimbach 1.62 # endif
169 heimbach 1.55 # endif
170 heimbach 1.62 # ifdef ALLOW_DEBUG
171 jmc 1.29 IF ( debugLevel .GE. debLevB )
172     & CALL DEBUG_CALL('SEAICE_MODEL',myThid)
173 heimbach 1.62 # endif
174 jmc 1.29 CALL TIMER_START('SEAICE_MODEL [DO_OCEANIC_PHYS]', myThid)
175     CALL SEAICE_MODEL( myTime, myIter, myThid )
176     CALL TIMER_STOP ('SEAICE_MODEL [DO_OCEANIC_PHYS]', myThid)
177 heimbach 1.62 # ifdef ALLOW_COST
178 heimbach 1.57 CALL SEAICE_COST_SENSI ( myTime, myIter, myThid )
179 heimbach 1.62 # endif
180 heimbach 1.35 ENDIF
181 jmc 1.29 #endif /* ALLOW_SEAICE */
182    
183 heimbach 1.64 #ifdef ALLOW_AUTODIFF_TAMC
184     CADJ STORE sst, sss = comlev1, key = ikey_dynamics
185     CADJ STORE qsw = comlev1, key = ikey_dynamics
186     # ifdef ALLOW_SEAICE
187     CADJ STORE area = comlev1, key = ikey_dynamics
188     # endif
189     #endif
190    
191 jscott 1.30 #if (defined ALLOW_THSICE) && !(defined ALLOW_ATM2D)
192 jmc 1.14 IF ( useThSIce .AND. fluidIsWater ) THEN
193 jmc 1.5 #ifdef ALLOW_DEBUG
194     IF ( debugLevel .GE. debLevB )
195     & CALL DEBUG_CALL('THSICE_MAIN',myThid)
196     #endif
197     C-- Step forward Therm.Sea-Ice variables
198     C and modify forcing terms including effects from ice
199     CALL TIMER_START('THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
200     CALL THSICE_MAIN( myTime, myIter, myThid )
201     CALL TIMER_STOP( 'THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
202     ENDIF
203     #endif /* ALLOW_THSICE */
204    
205 mlosch 1.21 #ifdef ALLOW_SHELFICE
206     IF ( useShelfIce .AND. fluidIsWater ) THEN
207     #ifdef ALLOW_DEBUG
208     IF ( debugLevel .GE. debLevB )
209     & CALL DEBUG_CALL('SHELFICE_THERMODYNAMICS',myThid)
210     #endif
211 jmc 1.47 C compute temperature and (virtual) salt flux at the
212 mlosch 1.21 C shelf-ice ocean interface
213     CALL TIMER_START('SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',
214     & myThid)
215     CALL SHELFICE_THERMODYNAMICS( myTime, myIter, myThid )
216     CALL TIMER_STOP( 'SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',
217     & myThid)
218     ENDIF
219     #endif /* ALLOW_SHELFICE */
220    
221 jmc 1.5 C-- Freeze water at the surface
222     #ifdef ALLOW_AUTODIFF_TAMC
223     CADJ STORE theta = comlev1, key = ikey_dynamics
224     #endif
225 heimbach 1.12 IF ( allowFreezing
226     & .AND. .NOT. useSEAICE
227 jmc 1.5 & .AND. .NOT. useThSIce ) THEN
228     CALL FREEZE_SURFACE( myTime, myIter, myThid )
229     ENDIF
230    
231 jmc 1.28 #ifdef ALLOW_OCN_COMPON_INTERF
232 jmc 1.5 C-- Apply imported data (from coupled interface) to forcing fields
233 jmc 1.28 C jmc: do not know precisely where to put this call (bf or af thSIce ?)
234 jmc 1.36 IF ( useCoupler ) THEN
235 jmc 1.5 CALL OCN_APPLY_IMPORT( .TRUE., myTime, myIter, myThid )
236 jmc 1.36 ENDIF
237 jmc 1.28 #endif /* ALLOW_OCN_COMPON_INTERF */
238 jmc 1.5
239 jmc 1.25 #ifdef ALLOW_BALANCE_FLUXES
240 jmc 1.36 C balance fluxes
241     IF ( balanceEmPmR )
242 jmc 1.25 & CALL REMOVE_MEAN_RS( 1, EmPmR, maskH, maskH, rA, drF,
243     & 'EmPmR', myTime, myThid )
244 jmc 1.36 IF ( balanceQnet )
245 jmc 1.25 & CALL REMOVE_MEAN_RS( 1, Qnet, maskH, maskH, rA, drF,
246     & 'Qnet ', myTime, myThid )
247     #endif /* ALLOW_BALANCE_FLUXES */
248    
249 jmc 1.1 #ifdef ALLOW_AUTODIFF_TAMC
250     C-- HPF directive to help TAMC
251     CHPF$ INDEPENDENT
252     #endif /* ALLOW_AUTODIFF_TAMC */
253     DO bj=myByLo(myThid),myByHi(myThid)
254     #ifdef ALLOW_AUTODIFF_TAMC
255 heimbach 1.15 C-- HPF directive to help TAMC
256     CHPF$ INDEPENDENT
257 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
258     DO bi=myBxLo(myThid),myBxHi(myThid)
259    
260     #ifdef ALLOW_AUTODIFF_TAMC
261     act1 = bi - myBxLo(myThid)
262     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
263     act2 = bj - myByLo(myThid)
264     max2 = myByHi(myThid) - myByLo(myThid) + 1
265     act3 = myThid - 1
266     max3 = nTx*nTy
267     act4 = ikey_dynamics - 1
268     itdkey = (act1 + 1) + act2*max1
269     & + act3*max1*max2
270     & + act4*max1*max2*max3
271 jmc 1.36 #else /* ALLOW_AUTODIFF_TAMC */
272 jmc 1.47 C if fluid is not water, by-pass find_rho, gmredi, surfaceForcing
273 jmc 1.36 C and all vertical mixing schemes, but keep OBCS_CALC
274     IF ( fluidIsWater ) THEN
275 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
276    
277     C-- Set up work arrays with valid (i.e. not NaN) values
278     C These inital values do not alter the numerical results. They
279     C just ensure that all memory references are to valid floating
280     C point numbers. This prevents spurious hardware signals due to
281     C uninitialised but inert locations.
282    
283 jmc 1.69 #ifdef ALLOW_AUTODIFF_TAMC
284 jmc 1.1 DO j=1-OLy,sNy+OLy
285     DO i=1-OLx,sNx+OLx
286 jmc 1.69 rhoKm1 (i,j) = 0. _d 0
287 jmc 1.47 rhoK (i,j) = 0. _d 0
288     rhoKp1 (i,j) = 0. _d 0
289 jmc 1.1 ENDDO
290     ENDDO
291 jmc 1.69 #endif /* ALLOW_AUTODIFF_TAMC */
292 jmc 1.1
293     DO k=1,Nr
294     DO j=1-OLy,sNy+OLy
295     DO i=1-OLx,sNx+OLx
296 jmc 1.69 C This is currently used by GMRedi, IVDC, MXL-depth and Diagnostics
297 jmc 1.1 sigmaX(i,j,k) = 0. _d 0
298     sigmaY(i,j,k) = 0. _d 0
299     sigmaR(i,j,k) = 0. _d 0
300     #ifdef ALLOW_AUTODIFF_TAMC
301     cph all the following init. are necessary for TAF
302     cph although some of these are re-initialised later.
303     IVDConvCount(i,j,k,bi,bj) = 0.
304     # ifdef ALLOW_GMREDI
305     Kwx(i,j,k,bi,bj) = 0. _d 0
306     Kwy(i,j,k,bi,bj) = 0. _d 0
307     Kwz(i,j,k,bi,bj) = 0. _d 0
308     # ifdef GM_NON_UNITY_DIAGONAL
309     Kux(i,j,k,bi,bj) = 0. _d 0
310     Kvy(i,j,k,bi,bj) = 0. _d 0
311     # endif
312     # ifdef GM_EXTRA_DIAGONAL
313     Kuz(i,j,k,bi,bj) = 0. _d 0
314     Kvz(i,j,k,bi,bj) = 0. _d 0
315     # endif
316     # ifdef GM_BOLUS_ADVEC
317     GM_PsiX(i,j,k,bi,bj) = 0. _d 0
318     GM_PsiY(i,j,k,bi,bj) = 0. _d 0
319     # endif
320     # ifdef GM_VISBECK_VARIABLE_K
321     VisbeckK(i,j,bi,bj) = 0. _d 0
322     # endif
323     # endif /* ALLOW_GMREDI */
324 heimbach 1.42 # ifdef ALLOW_KPP
325     KPPdiffKzS(i,j,k,bi,bj) = 0. _d 0
326     KPPdiffKzT(i,j,k,bi,bj) = 0. _d 0
327     # endif /* ALLOW_KPP */
328 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
329     ENDDO
330     ENDDO
331     ENDDO
332    
333     iMin = 1-OLx
334     iMax = sNx+OLx
335     jMin = 1-OLy
336     jMax = sNy+OLy
337    
338     #ifdef ALLOW_AUTODIFF_TAMC
339     CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
340     CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
341 heimbach 1.12 CADJ STORE totphihyd(:,:,:,bi,bj)
342 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
343 heimbach 1.10 # ifdef ALLOW_KPP
344 jmc 1.1 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
345     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
346 heimbach 1.10 # endif
347 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
348    
349     #ifdef ALLOW_DEBUG
350 jmc 1.36 IF ( debugLevel .GE. debLevB )
351 jmc 1.1 & CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
352     #endif
353    
354     C-- Start of diagnostic loop
355     DO k=Nr,1,-1
356    
357     #ifdef ALLOW_AUTODIFF_TAMC
358     C? Patrick, is this formula correct now that we change the loop range?
359     C? Do we still need this?
360 jmc 1.69 cph kkey formula corrected.
361 jmc 1.47 cph Needed for rhoK, rhoKm1, in the case useGMREDI.
362 jmc 1.1 kkey = (itdkey-1)*Nr + k
363     #endif /* ALLOW_AUTODIFF_TAMC */
364    
365 jmc 1.69 #ifdef ALLOW_DOWN_SLOPE
366     IF ( useDOWN_SLOPE ) THEN
367     CALL DWNSLP_CALC_RHO(
368     I theta, salt,
369     O rhoK, rhoInSitu,
370     I k, bi, bj, myTime, myIter, myThid )
371     ENDIF
372     #endif /* ALLOW_DOWN_SLOPE */
373    
374 jmc 1.1 C-- Calculate gradients of potential density for isoneutral
375     C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
376 jmc 1.17 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.)
377 dimitri 1.61 & .OR. useSALT_PLUME .OR. doDiagsRho.GE.1 ) THEN
378 jmc 1.1 #ifdef ALLOW_DEBUG
379 jmc 1.36 IF ( debugLevel .GE. debLevB )
380 jmc 1.68 & CALL DEBUG_CALL('FIND_RHO_2D',myThid)
381 jmc 1.1 #endif
382 jmc 1.69 IF ( .NOT.useDOWN_SLOPE ) THEN
383 jmc 1.1 #ifdef ALLOW_AUTODIFF_TAMC
384     CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
385     CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
386     #endif /* ALLOW_AUTODIFF_TAMC */
387 jmc 1.69 CALL FIND_RHO_2D(
388 jmc 1.68 I iMin, iMax, jMin, jMax, k,
389     I theta(1-OLx,1-OLy,k,bi,bj),
390     I salt (1-OLx,1-OLy,k,bi,bj),
391     O rhoK,
392     I k, bi, bj, myThid )
393 jmc 1.1
394 jmc 1.69 ENDIF
395 jmc 1.1 IF (k.GT.1) THEN
396     #ifdef ALLOW_AUTODIFF_TAMC
397     CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
398     CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
399     #endif /* ALLOW_AUTODIFF_TAMC */
400 jmc 1.68 CALL FIND_RHO_2D(
401     I iMin, iMax, jMin, jMax, k,
402     I theta(1-OLx,1-OLy,k-1,bi,bj),
403     I salt (1-OLx,1-OLy,k-1,bi,bj),
404     O rhoKm1,
405     I k-1, bi, bj, myThid )
406 jmc 1.1 ENDIF
407     #ifdef ALLOW_DEBUG
408 jmc 1.36 IF ( debugLevel .GE. debLevB )
409 jmc 1.1 & CALL DEBUG_CALL('GRAD_SIGMA',myThid)
410     #endif
411 heimbach 1.31 cph Avoid variable aliasing for adjoint !!!
412     DO j=jMin,jMax
413     DO i=iMin,iMax
414 jmc 1.47 rhoKp1(i,j) = rhoK(i,j)
415 heimbach 1.31 ENDDO
416     ENDDO
417 jmc 1.1 CALL GRAD_SIGMA(
418     I bi, bj, iMin, iMax, jMin, jMax, k,
419 heimbach 1.31 I rhoK, rhoKm1, rhoKp1,
420 jmc 1.1 O sigmaX, sigmaY, sigmaR,
421     I myThid )
422 gforget 1.66 #ifdef ALLOW_AUTODIFF_TAMC
423 jmc 1.69 #ifdef GMREDI_WITH_STABLE_ADJOINT
424 gforget 1.66 cgf zero out adjoint fields to stabilize pkg/gmredi adjoint
425     cgf -> cuts adjoint dependency from slope to state
426 jmc 1.69 CALL ZERO_ADJ_LOC( Nr, sigmaX, myThid)
427     CALL ZERO_ADJ_LOC( Nr, sigmaY, myThid)
428     CALL ZERO_ADJ_LOC( Nr, sigmaR, myThid)
429 gforget 1.66 #endif
430     #endif /* ALLOW_AUTODIFF_TAMC */
431 jmc 1.1 ENDIF
432    
433     C-- Implicit Vertical Diffusion for Convection
434     c ==> should use sigmaR !!!
435     IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
436     #ifdef ALLOW_DEBUG
437 jmc 1.36 IF ( debugLevel .GE. debLevB )
438 jmc 1.1 & CALL DEBUG_CALL('CALC_IVDC',myThid)
439     #endif
440     CALL CALC_IVDC(
441     I bi, bj, iMin, iMax, jMin, jMax, k,
442     I rhoKm1, rhoK,
443     I myTime, myIter, myThid)
444     ENDIF
445    
446 jmc 1.17 #ifdef ALLOW_DIAGNOSTICS
447     IF ( doDiagsRho.GE.2 ) THEN
448     CALL DIAGS_RHO( k, bi, bj,
449     I rhoK, rhoKm1,
450     I myTime, myIter, myThid)
451     ENDIF
452     #endif
453    
454 jmc 1.1 C-- end of diagnostic k loop (Nr:1)
455     ENDDO
456    
457 heimbach 1.57 #ifdef ALLOW_AUTODIFF_TAMC
458 jmc 1.69 CADJ STORE IVDConvCount(:,:,:,bi,bj)
459 heimbach 1.57 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
460     #endif
461    
462 jmc 1.47 C-- Diagnose Mixed Layer Depth:
463     IF ( useGMRedi .OR. doDiagsRho.GE.1 ) THEN
464     CALL CALC_OCE_MXLAYER( rhoK, sigmaR,
465     & bi, bj, myTime, myIter, myThid )
466     ENDIF
467 heimbach 1.53
468 dimitri 1.52 #ifdef ALLOW_SALT_PLUME
469 dimitri 1.61 IF ( useSALT_PLUME ) THEN
470 dimitri 1.60 CALL SALT_PLUME_CALC_DEPTH( rhoK, sigmaR,
471 dimitri 1.52 & bi, bj, myTime, myIter, myThid )
472 dimitri 1.60 ENDIF
473 dimitri 1.61 #endif /* ALLOW_SALT_PLUME */
474    
475 jmc 1.8 #ifdef ALLOW_DIAGNOSTICS
476 jmc 1.17 IF ( doDiagsRho.GE.1 ) THEN
477 jmc 1.16 CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR ', 0, Nr,
478     & 2, bi, bj, myThid)
479 jmc 1.8 ENDIF
480 dimitri 1.61 #endif /* ALLOW_DIAGNOSTICS */
481 jmc 1.8
482 jmc 1.1 C-- Determines forcing terms based on external fields
483     C relaxation terms, etc.
484     #ifdef ALLOW_DEBUG
485 jmc 1.36 IF ( debugLevel .GE. debLevB )
486 jmc 1.1 & CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
487     #endif
488 heimbach 1.23 #ifdef ALLOW_AUTODIFF_TAMC
489     CADJ STORE EmPmR(:,:,bi,bj)
490     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
491 heimbach 1.26 # ifdef EXACT_CONSERV
492 heimbach 1.23 CADJ STORE PmEpR(:,:,bi,bj)
493     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
494 heimbach 1.26 # endif
495 heimbach 1.27 # ifdef NONLIN_FRSURF
496     CADJ STORE hFac_surfC(:,:,bi,bj)
497     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
498     CADJ STORE recip_hFacC(:,:,:,bi,bj)
499     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
500     # endif
501 heimbach 1.23 #endif
502 jmc 1.36 CALL EXTERNAL_FORCING_SURF(
503 jmc 1.1 I bi, bj, iMin, iMax, jMin, jMax,
504     I myTime, myIter, myThid )
505 heimbach 1.27 #ifdef ALLOW_AUTODIFF_TAMC
506     # ifdef EXACT_CONSERV
507     cph-test
508     cphCADJ STORE PmEpR(:,:,bi,bj)
509     cphCADJ & = comlev1_bibj, key=itdkey, byte=isbyte
510     # endif
511     #endif
512 jmc 1.1
513     #ifdef ALLOW_AUTODIFF_TAMC
514     cph needed for KPP
515 jmc 1.4 CADJ STORE surfaceForcingU(:,:,bi,bj)
516 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
517 jmc 1.4 CADJ STORE surfaceForcingV(:,:,bi,bj)
518 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
519 jmc 1.4 CADJ STORE surfaceForcingS(:,:,bi,bj)
520 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
521 jmc 1.4 CADJ STORE surfaceForcingT(:,:,bi,bj)
522 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
523 jmc 1.4 CADJ STORE surfaceForcingTice(:,:,bi,bj)
524 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
525     #endif /* ALLOW_AUTODIFF_TAMC */
526    
527     #ifdef ALLOW_KPP
528     C-- Compute KPP mixing coefficients
529     IF (useKPP) THEN
530     #ifdef ALLOW_DEBUG
531 jmc 1.36 IF ( debugLevel .GE. debLevB )
532 jmc 1.1 & CALL DEBUG_CALL('KPP_CALC',myThid)
533     #endif
534     CALL KPP_CALC(
535 jmc 1.44 I bi, bj, myTime, myIter, myThid )
536 jmc 1.1 #ifdef ALLOW_AUTODIFF_TAMC
537     ELSE
538     CALL KPP_CALC_DUMMY(
539 jmc 1.44 I bi, bj, myTime, myIter, myThid )
540 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
541     ENDIF
542    
543     #endif /* ALLOW_KPP */
544    
545 mlosch 1.6 #ifdef ALLOW_PP81
546     C-- Compute PP81 mixing coefficients
547     IF (usePP81) THEN
548     #ifdef ALLOW_DEBUG
549 jmc 1.36 IF ( debugLevel .GE. debLevB )
550 mlosch 1.6 & CALL DEBUG_CALL('PP81_CALC',myThid)
551     #endif
552     CALL PP81_CALC(
553     I bi, bj, myTime, myThid )
554     ENDIF
555     #endif /* ALLOW_PP81 */
556    
557     #ifdef ALLOW_MY82
558     C-- Compute MY82 mixing coefficients
559     IF (useMY82) THEN
560     #ifdef ALLOW_DEBUG
561 jmc 1.36 IF ( debugLevel .GE. debLevB )
562 mlosch 1.6 & CALL DEBUG_CALL('MY82_CALC',myThid)
563     #endif
564     CALL MY82_CALC(
565     I bi, bj, myTime, myThid )
566     ENDIF
567     #endif /* ALLOW_MY82 */
568    
569 mlosch 1.9 #ifdef ALLOW_GGL90
570     C-- Compute GGL90 mixing coefficients
571     IF (useGGL90) THEN
572     #ifdef ALLOW_DEBUG
573 jmc 1.36 IF ( debugLevel .GE. debLevB )
574 mlosch 1.9 & CALL DEBUG_CALL('GGL90_CALC',myThid)
575     #endif
576     CALL GGL90_CALC(
577     I bi, bj, myTime, myThid )
578     ENDIF
579     #endif /* ALLOW_GGL90 */
580    
581 jmc 1.20 #ifdef ALLOW_TIMEAVE
582 jmc 1.36 IF ( taveFreq.GT. 0. _d 0 ) THEN
583 jmc 1.20 CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid)
584     ENDIF
585     IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
586     CALL TIMEAVE_CUMULATE(ConvectCountTave, IVDConvCount,
587     I Nr, deltaTclock, bi, bj, myThid)
588     ENDIF
589     #endif /* ALLOW_TIMEAVE */
590    
591 jmc 1.69 #ifdef ALLOW_GMREDI
592 jmc 1.47 #ifdef ALLOW_AUTODIFF_TAMC
593     # ifndef GM_EXCLUDE_CLIPPING
594     cph storing here is needed only for one GMREDI_OPTIONS:
595     cph define GM_BOLUS_ADVEC
596     cph keep it although TAF says you dont need to.
597     cph but I've avoided the #ifdef for now, in case more things change
598     CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
599     CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
600     CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
601     # endif
602     #endif /* ALLOW_AUTODIFF_TAMC */
603    
604     C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
605     IF (useGMRedi) THEN
606     #ifdef ALLOW_DEBUG
607     IF ( debugLevel .GE. debLevB )
608     & CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
609     #endif
610     CALL GMREDI_CALC_TENSOR(
611 jmc 1.51 I iMin, iMax, jMin, jMax,
612 jmc 1.47 I sigmaX, sigmaY, sigmaR,
613 jmc 1.51 I bi, bj, myTime, myIter, myThid )
614 jmc 1.47 #ifdef ALLOW_AUTODIFF_TAMC
615     ELSE
616     CALL GMREDI_CALC_TENSOR_DUMMY(
617 jmc 1.51 I iMin, iMax, jMin, jMax,
618 jmc 1.47 I sigmaX, sigmaY, sigmaR,
619 jmc 1.51 I bi, bj, myTime, myIter, myThid )
620 jmc 1.47 #endif /* ALLOW_AUTODIFF_TAMC */
621     ENDIF
622 jmc 1.69 #endif /* ALLOW_GMREDI */
623    
624     #ifdef ALLOW_DOWN_SLOPE
625     IF ( useDOWN_SLOPE ) THEN
626     C-- Calculate Downsloping Flow for Down_Slope parameterization
627     IF ( usingPCoords ) THEN
628     CALL DWNSLP_CALC_FLOW(
629     I bi, bj, kSurfC,
630     I rhoInSitu,
631     I myTime, myIter, myThid )
632     ELSE
633     CALL DWNSLP_CALC_FLOW(
634     I bi, bj, kLowC,
635     I rhoInSitu,
636     I myTime, myIter, myThid )
637     ENDIF
638     ENDIF
639     #endif /* ALLOW_DOWN_SLOPE */
640 jmc 1.47
641 jmc 1.36 #ifndef ALLOW_AUTODIFF_TAMC
642     C--- if fluid Is Water: end
643     ENDIF
644     #endif
645    
646     #ifdef ALLOW_OBCS
647     C-- Calculate future values on open boundaries
648     IF (useOBCS) THEN
649     #ifdef ALLOW_DEBUG
650     IF ( debugLevel .GE. debLevB )
651     & CALL DEBUG_CALL('OBCS_CALC',myThid)
652     #endif
653     CALL OBCS_CALC( bi, bj, myTime+deltaTclock, myIter+1,
654     I uVel, vVel, wVel, theta, salt,
655     I myThid )
656     ENDIF
657     #endif /* ALLOW_OBCS */
658    
659 jmc 1.1 C-- end bi,bj loops.
660     ENDDO
661     ENDDO
662    
663 jmc 1.45 #ifdef ALLOW_KPP
664     IF (useKPP) THEN
665     CALL KPP_DO_EXCH( myThid )
666     ENDIF
667     #endif /* ALLOW_KPP */
668    
669 jmc 1.18 #ifdef ALLOW_DIAGNOSTICS
670     IF ( fluidIsWater .AND. useDiagnostics ) THEN
671     CALL DIAGS_OCEANIC_SURF_FLUX( myTime, myIter, myThid )
672     ENDIF
673 jmc 1.19 IF ( ivdc_kappa.NE.0 .AND. useDiagnostics ) THEN
674     CALL DIAGNOSTICS_FILL( IVDConvCount,'CONVADJ ',
675     & 0, Nr, 0, 1, 1, myThid )
676     ENDIF
677 jmc 1.18 #endif
678    
679 jmc 1.1 #ifdef ALLOW_DEBUG
680 jmc 1.29 IF ( debugLevel .GE. debLevB )
681     & CALL DEBUG_LEAVE('DO_OCEANIC_PHYS',myThid)
682 jmc 1.1 #endif
683    
684     RETURN
685     END

  ViewVC Help
Powered by ViewVC 1.1.22