/[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.30 - (hide annotations) (download)
Wed Sep 6 15:30:25 2006 UTC (17 years, 9 months ago) by jscott
Branch: MAIN
CVS Tags: checkpoint58q_post, checkpoint58p_post
Changes since 1.29: +2 -2 lines
add hooks for atm2d package

1 jscott 1.30 C $Header: /u/gcmpack/MITgcm/model/src/do_oceanic_phys.F,v 1.29 2006/06/15 16:27:54 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     #include "DYNVARS.h"
40     #include "GRID.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     # include "EOS.h"
53     # ifdef ALLOW_KPP
54     # include "KPP.h"
55     # endif
56     # ifdef ALLOW_GMREDI
57     # include "GMREDI.h"
58     # endif
59     # ifdef ALLOW_EBM
60     # include "EBM.h"
61     # endif
62 heimbach 1.10 # ifdef EXACT_CONSERV
63     # include "SURFACE.h"
64     # endif
65 jmc 1.29 # ifdef ALLOW_EXF
66     # include "ctrl.h"
67     # include "exf_fields.h"
68     # include "exf_clim_fields.h"
69     # ifdef ALLOW_BULKFORMULAE
70     # include "exf_constants.h"
71     # endif
72     # endif
73     # ifdef ALLOW_SEAICE
74     # include "SEAICE.h"
75     # endif
76 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
77    
78     C !INPUT/OUTPUT PARAMETERS:
79     C == Routine arguments ==
80 jmc 1.18 C myTime :: Current time in simulation
81     C myIter :: Current iteration number in simulation
82     C myThid :: Thread number for this instance of the routine.
83 jmc 1.1 _RL myTime
84     INTEGER myIter
85     INTEGER myThid
86    
87     C !LOCAL VARIABLES:
88     C == Local variables
89 jmc 1.18 C rhoK, rhoKM1 :: Density at current level, and level above
90     C iMin, iMax :: Ranges and sub-block indices on which calculations
91 jmc 1.1 C jMin, jMax are applied.
92 jmc 1.18 C bi, bj :: tile indices
93     C i,j,k :: loop indices
94 jmc 1.1 _RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95     _RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96     _RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
97     _RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
98     _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
99     INTEGER iMin, iMax
100     INTEGER jMin, jMax
101     INTEGER bi, bj
102 jmc 1.18 INTEGER i, j, k
103 jmc 1.17 INTEGER doDiagsRho
104     #ifdef ALLOW_DIAGNOSTICS
105     LOGICAL DIAGNOSTICS_IS_ON
106     EXTERNAL DIAGNOSTICS_IS_ON
107     #endif /* ALLOW_DIAGNOSTICS */
108 jmc 1.1
109     CEOP
110    
111 heimbach 1.12 #ifdef ALLOW_AUTODIFF_TAMC
112     C-- dummy statement to end declaration part
113     itdkey = 1
114     #endif /* ALLOW_AUTODIFF_TAMC */
115    
116 jmc 1.1 #ifdef ALLOW_DEBUG
117 jmc 1.29 IF ( debugLevel .GE. debLevB )
118     & CALL DEBUG_ENTER('DO_OCEANIC_PHYS',myThid)
119 jmc 1.1 #endif
120    
121 jmc 1.17 doDiagsRho = 0
122     #ifdef ALLOW_DIAGNOSTICS
123     IF ( useDiagnostics .AND. fluidIsWater ) THEN
124     IF ( DIAGNOSTICS_IS_ON('DRHODR ',myThid) ) doDiagsRho = 1
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     ENDIF
131     #endif /* ALLOW_DIAGNOSTICS */
132    
133 jmc 1.29 #ifdef ALLOW_SEAICE
134     C-- Call sea ice model to compute forcing/external data fields. In
135     C addition to computing prognostic sea-ice variables and diagnosing the
136     C forcing/external data fields that drive the ocean model, SEAICE_MODEL
137     C also sets theta to the freezing point under sea-ice. The implied
138     C surface heat flux is then stored in variable surfaceTendencyTice,
139     C which is needed by KPP package (kpp_calc.F and kpp_transport_t.F)
140     C to diagnose surface buoyancy fluxes and for the non-local transport
141     C term. Because this call precedes model thermodynamics, temperature
142     C under sea-ice may not be "exactly" at the freezing point by the time
143     C theta is dumped or time-averaged.
144     IF ( useSEAICE ) THEN
145     #ifdef ALLOW_AUTODIFF_TAMC
146     CADJ STORE aqh,precip,swdown = comlev1, key = ikey_dynamics
147     CADJ STORE theta = comlev1, key = ikey_dynamics
148     # ifdef SEAICE_ALLOW_DYNAMICS
149     CADJ STORE uvel,vvel = comlev1, key = ikey_dynamics
150     # endif
151     #endif
152     #ifdef ALLOW_DEBUG
153     IF ( debugLevel .GE. debLevB )
154     & CALL DEBUG_CALL('SEAICE_MODEL',myThid)
155     #endif
156     CALL TIMER_START('SEAICE_MODEL [DO_OCEANIC_PHYS]', myThid)
157     CALL SEAICE_MODEL( myTime, myIter, myThid )
158     CALL TIMER_STOP ('SEAICE_MODEL [DO_OCEANIC_PHYS]', myThid)
159     #ifdef ALLOW_COST_ICE
160     CALL COST_ICE_TEST ( myTime, myIter, myThid )
161     #endif
162     ENDIF
163     #endif /* ALLOW_SEAICE */
164    
165 jscott 1.30 #if (defined ALLOW_THSICE) && !(defined ALLOW_ATM2D)
166 jmc 1.14 IF ( useThSIce .AND. fluidIsWater ) THEN
167 jmc 1.5 #ifdef ALLOW_DEBUG
168     IF ( debugLevel .GE. debLevB )
169     & CALL DEBUG_CALL('THSICE_MAIN',myThid)
170     #endif
171     C-- Step forward Therm.Sea-Ice variables
172     C and modify forcing terms including effects from ice
173     CALL TIMER_START('THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
174     CALL THSICE_MAIN( myTime, myIter, myThid )
175     CALL TIMER_STOP( 'THSICE_MAIN [DO_OCEANIC_PHYS]', myThid)
176     ENDIF
177     #endif /* ALLOW_THSICE */
178    
179 mlosch 1.21 #ifdef ALLOW_SHELFICE
180     IF ( useShelfIce .AND. fluidIsWater ) THEN
181     #ifdef ALLOW_DEBUG
182     IF ( debugLevel .GE. debLevB )
183     & CALL DEBUG_CALL('SHELFICE_THERMODYNAMICS',myThid)
184     #endif
185     C compute temperature and (virtual) salt flux at the
186     C shelf-ice ocean interface
187     CALL TIMER_START('SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',
188     & myThid)
189     CALL SHELFICE_THERMODYNAMICS( myTime, myIter, myThid )
190     CALL TIMER_STOP( 'SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]',
191     & myThid)
192     ENDIF
193     #endif /* ALLOW_SHELFICE */
194    
195 jmc 1.5 C-- Freeze water at the surface
196     #ifdef ALLOW_AUTODIFF_TAMC
197     CADJ STORE theta = comlev1, key = ikey_dynamics
198     #endif
199 heimbach 1.12 IF ( allowFreezing
200     & .AND. .NOT. useSEAICE
201 jmc 1.5 & .AND. .NOT. useThSIce ) THEN
202     CALL FREEZE_SURFACE( myTime, myIter, myThid )
203     ENDIF
204    
205 jmc 1.28 #ifdef ALLOW_OCN_COMPON_INTERF
206 jmc 1.5 C-- Apply imported data (from coupled interface) to forcing fields
207 jmc 1.28 C jmc: do not know precisely where to put this call (bf or af thSIce ?)
208 jmc 1.5 IF ( useCoupler ) THEN
209     CALL OCN_APPLY_IMPORT( .TRUE., myTime, myIter, myThid )
210     ENDIF
211 jmc 1.28 #endif /* ALLOW_OCN_COMPON_INTERF */
212 jmc 1.5
213 jmc 1.25 #ifdef ALLOW_BALANCE_FLUXES
214     C balance fluxes
215     IF ( balanceEmPmR )
216     & CALL REMOVE_MEAN_RS( 1, EmPmR, maskH, maskH, rA, drF,
217     & 'EmPmR', myTime, myThid )
218     IF ( balanceQnet )
219     & CALL REMOVE_MEAN_RS( 1, Qnet, maskH, maskH, rA, drF,
220     & 'Qnet ', myTime, myThid )
221     #endif /* ALLOW_BALANCE_FLUXES */
222    
223 jmc 1.1 #ifdef ALLOW_AUTODIFF_TAMC
224     C-- HPF directive to help TAMC
225     CHPF$ INDEPENDENT
226     #endif /* ALLOW_AUTODIFF_TAMC */
227     DO bj=myByLo(myThid),myByHi(myThid)
228     #ifdef ALLOW_AUTODIFF_TAMC
229 heimbach 1.15 C-- HPF directive to help TAMC
230     CHPF$ INDEPENDENT
231 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
232     DO bi=myBxLo(myThid),myBxHi(myThid)
233    
234     #ifdef ALLOW_AUTODIFF_TAMC
235     act1 = bi - myBxLo(myThid)
236     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
237     act2 = bj - myByLo(myThid)
238     max2 = myByHi(myThid) - myByLo(myThid) + 1
239     act3 = myThid - 1
240     max3 = nTx*nTy
241     act4 = ikey_dynamics - 1
242     itdkey = (act1 + 1) + act2*max1
243     & + act3*max1*max2
244     & + act4*max1*max2*max3
245     #endif /* ALLOW_AUTODIFF_TAMC */
246    
247     C-- Set up work arrays with valid (i.e. not NaN) values
248     C These inital values do not alter the numerical results. They
249     C just ensure that all memory references are to valid floating
250     C point numbers. This prevents spurious hardware signals due to
251     C uninitialised but inert locations.
252    
253     DO j=1-OLy,sNy+OLy
254     DO i=1-OLx,sNx+OLx
255     rhok (i,j) = 0. _d 0
256     rhoKM1 (i,j) = 0. _d 0
257     ENDDO
258     ENDDO
259    
260     DO k=1,Nr
261     DO j=1-OLy,sNy+OLy
262     DO i=1-OLx,sNx+OLx
263     C This is currently also used by IVDC and Diagnostics
264     sigmaX(i,j,k) = 0. _d 0
265     sigmaY(i,j,k) = 0. _d 0
266     sigmaR(i,j,k) = 0. _d 0
267     #ifdef ALLOW_AUTODIFF_TAMC
268     cph all the following init. are necessary for TAF
269     cph although some of these are re-initialised later.
270     IVDConvCount(i,j,k,bi,bj) = 0.
271     # ifdef ALLOW_GMREDI
272     Kwx(i,j,k,bi,bj) = 0. _d 0
273     Kwy(i,j,k,bi,bj) = 0. _d 0
274     Kwz(i,j,k,bi,bj) = 0. _d 0
275     # ifdef GM_NON_UNITY_DIAGONAL
276     Kux(i,j,k,bi,bj) = 0. _d 0
277     Kvy(i,j,k,bi,bj) = 0. _d 0
278     # endif
279     # ifdef GM_EXTRA_DIAGONAL
280     Kuz(i,j,k,bi,bj) = 0. _d 0
281     Kvz(i,j,k,bi,bj) = 0. _d 0
282     # endif
283     # ifdef GM_BOLUS_ADVEC
284     GM_PsiX(i,j,k,bi,bj) = 0. _d 0
285     GM_PsiY(i,j,k,bi,bj) = 0. _d 0
286     # endif
287     # ifdef GM_VISBECK_VARIABLE_K
288     VisbeckK(i,j,bi,bj) = 0. _d 0
289     # endif
290     # endif /* ALLOW_GMREDI */
291     #endif /* ALLOW_AUTODIFF_TAMC */
292     ENDDO
293     ENDDO
294     ENDDO
295    
296     iMin = 1-OLx
297     iMax = sNx+OLx
298     jMin = 1-OLy
299     jMax = sNy+OLy
300    
301     #ifdef ALLOW_AUTODIFF_TAMC
302     CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
303     CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
304 heimbach 1.12 CADJ STORE totphihyd(:,:,:,bi,bj)
305 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
306 heimbach 1.10 # ifdef ALLOW_KPP
307 jmc 1.1 CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
308     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
309 heimbach 1.10 # endif
310 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
311    
312     #ifdef ALLOW_DEBUG
313     IF ( debugLevel .GE. debLevB )
314     & CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
315     #endif
316    
317     C-- Start of diagnostic loop
318     DO k=Nr,1,-1
319    
320     #ifdef ALLOW_AUTODIFF_TAMC
321     C? Patrick, is this formula correct now that we change the loop range?
322     C? Do we still need this?
323     cph kkey formula corrected.
324     cph Needed for rhok, rhokm1, in the case useGMREDI.
325     kkey = (itdkey-1)*Nr + k
326     #endif /* ALLOW_AUTODIFF_TAMC */
327    
328     C-- Calculate gradients of potential density for isoneutral
329     C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
330     c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN
331 jmc 1.17 IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.)
332     & .OR. doDiagsRho.GE.1 ) THEN
333 jmc 1.1 #ifdef ALLOW_DEBUG
334     IF ( debugLevel .GE. debLevB )
335     & CALL DEBUG_CALL('FIND_RHO',myThid)
336     #endif
337     #ifdef ALLOW_AUTODIFF_TAMC
338     CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
339     CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
340     #endif /* ALLOW_AUTODIFF_TAMC */
341     CALL FIND_RHO(
342     I bi, bj, iMin, iMax, jMin, jMax, k, k,
343     I theta, salt,
344     O rhoK,
345     I myThid )
346    
347     IF (k.GT.1) THEN
348     #ifdef ALLOW_AUTODIFF_TAMC
349     CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
350     CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
351     #endif /* ALLOW_AUTODIFF_TAMC */
352     CALL FIND_RHO(
353     I bi, bj, iMin, iMax, jMin, jMax, k-1, k,
354     I theta, salt,
355     O rhoKm1,
356     I myThid )
357     ENDIF
358     #ifdef ALLOW_DEBUG
359     IF ( debugLevel .GE. debLevB )
360     & CALL DEBUG_CALL('GRAD_SIGMA',myThid)
361     #endif
362     CALL GRAD_SIGMA(
363     I bi, bj, iMin, iMax, jMin, jMax, k,
364     I rhoK, rhoKm1, rhoK,
365     O sigmaX, sigmaY, sigmaR,
366     I myThid )
367     ENDIF
368    
369     #ifdef ALLOW_AUTODIFF_TAMC
370 heimbach 1.12 ctest# ifndef GM_EXCLUDE_CLIPPING
371 jmc 1.1 CADJ STORE rhok (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
372 heimbach 1.12 ctest# endif
373 jmc 1.1 CADJ STORE rhokm1 (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte
374     #endif /* ALLOW_AUTODIFF_TAMC */
375     C-- Implicit Vertical Diffusion for Convection
376     c ==> should use sigmaR !!!
377     IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
378     #ifdef ALLOW_DEBUG
379     IF ( debugLevel .GE. debLevB )
380     & CALL DEBUG_CALL('CALC_IVDC',myThid)
381     #endif
382     CALL CALC_IVDC(
383     I bi, bj, iMin, iMax, jMin, jMax, k,
384     I rhoKm1, rhoK,
385     I myTime, myIter, myThid)
386     ENDIF
387    
388 jmc 1.17 #ifdef ALLOW_DIAGNOSTICS
389     IF ( doDiagsRho.GE.2 ) THEN
390     CALL DIAGS_RHO( k, bi, bj,
391     I rhoK, rhoKm1,
392     I myTime, myIter, myThid)
393     ENDIF
394     #endif
395    
396 jmc 1.1 C-- end of diagnostic k loop (Nr:1)
397     ENDDO
398    
399 jmc 1.8 #ifdef ALLOW_DIAGNOSTICS
400 jmc 1.17 c IF ( useDiagnostics .AND.
401     c & (useGMRedi .OR. ivdc_kappa.NE.0.) ) THEN
402     IF ( doDiagsRho.GE.1 ) THEN
403 jmc 1.16 CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR ', 0, Nr,
404     & 2, bi, bj, myThid)
405 jmc 1.8 ENDIF
406     #endif
407    
408 jmc 1.1 #ifdef ALLOW_OBCS
409     C-- Calculate future values on open boundaries
410     IF (useOBCS) THEN
411     #ifdef ALLOW_DEBUG
412     IF ( debugLevel .GE. debLevB )
413     & CALL DEBUG_CALL('OBCS_CALC',myThid)
414     #endif
415 heimbach 1.11 CALL OBCS_CALC( bi, bj, myTime+deltaTclock, myIter+1,
416 jmc 1.1 I uVel, vVel, wVel, theta, salt,
417     I myThid )
418     ENDIF
419     #endif /* ALLOW_OBCS */
420    
421     #ifndef ALLOW_AUTODIFF_TAMC
422 jmc 1.14 IF ( fluidIsWater ) THEN
423 jmc 1.1 #endif
424     C-- Determines forcing terms based on external fields
425     C relaxation terms, etc.
426     #ifdef ALLOW_DEBUG
427     IF ( debugLevel .GE. debLevB )
428     & CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
429     #endif
430 heimbach 1.23 #ifdef ALLOW_AUTODIFF_TAMC
431     CADJ STORE EmPmR(:,:,bi,bj)
432     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
433 heimbach 1.26 # ifdef EXACT_CONSERV
434 heimbach 1.23 CADJ STORE PmEpR(:,:,bi,bj)
435     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
436 heimbach 1.26 # endif
437 heimbach 1.27 # ifdef NONLIN_FRSURF
438     CADJ STORE hFac_surfC(:,:,bi,bj)
439     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
440     CADJ STORE recip_hFacC(:,:,:,bi,bj)
441     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
442     # endif
443 heimbach 1.23 #endif
444 jmc 1.1 CALL EXTERNAL_FORCING_SURF(
445     I bi, bj, iMin, iMax, jMin, jMax,
446     I myTime, myIter, myThid )
447     #ifndef ALLOW_AUTODIFF_TAMC
448     ENDIF
449     #endif
450 heimbach 1.27 #ifdef ALLOW_AUTODIFF_TAMC
451     # ifdef EXACT_CONSERV
452     cph-test
453     cphCADJ STORE PmEpR(:,:,bi,bj)
454     cphCADJ & = comlev1_bibj, key=itdkey, byte=isbyte
455     # endif
456     #endif
457 jmc 1.1
458     #ifdef ALLOW_AUTODIFF_TAMC
459     cph needed for KPP
460 jmc 1.4 CADJ STORE surfaceForcingU(:,:,bi,bj)
461 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
462 jmc 1.4 CADJ STORE surfaceForcingV(:,:,bi,bj)
463 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
464 jmc 1.4 CADJ STORE surfaceForcingS(:,:,bi,bj)
465 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
466 jmc 1.4 CADJ STORE surfaceForcingT(:,:,bi,bj)
467 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
468 jmc 1.4 CADJ STORE surfaceForcingTice(:,:,bi,bj)
469 jmc 1.1 CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
470     #endif /* ALLOW_AUTODIFF_TAMC */
471    
472     #ifdef ALLOW_GMREDI
473    
474     #ifdef ALLOW_AUTODIFF_TAMC
475 heimbach 1.10 # ifndef GM_EXCLUDE_CLIPPING
476 jmc 1.1 cph storing here is needed only for one GMREDI_OPTIONS:
477     cph define GM_BOLUS_ADVEC
478 heimbach 1.10 cph keep it although TAF says you dont need to.
479 jmc 1.1 cph but I've avoided the #ifdef for now, in case more things change
480     CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
481     CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
482 heimbach 1.12 CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
483 heimbach 1.10 # endif
484 jmc 1.1 #endif /* ALLOW_AUTODIFF_TAMC */
485    
486     C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
487     IF (useGMRedi) THEN
488     #ifdef ALLOW_DEBUG
489     IF ( debugLevel .GE. debLevB )
490     & CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
491     #endif
492     CALL GMREDI_CALC_TENSOR(
493     I bi, bj, iMin, iMax, jMin, jMax,
494     I sigmaX, sigmaY, sigmaR,
495     I myThid )
496     #ifdef ALLOW_AUTODIFF_TAMC
497     ELSE
498     CALL GMREDI_CALC_TENSOR_DUMMY(
499     I bi, bj, iMin, iMax, jMin, jMax,
500     I sigmaX, sigmaY, sigmaR,
501     I myThid )
502     #endif /* ALLOW_AUTODIFF_TAMC */
503     ENDIF
504    
505     #endif /* ALLOW_GMREDI */
506    
507     #ifdef ALLOW_KPP
508     C-- Compute KPP mixing coefficients
509     IF (useKPP) THEN
510     #ifdef ALLOW_DEBUG
511     IF ( debugLevel .GE. debLevB )
512     & CALL DEBUG_CALL('KPP_CALC',myThid)
513     #endif
514     CALL KPP_CALC(
515     I bi, bj, myTime, myThid )
516     #ifdef ALLOW_AUTODIFF_TAMC
517     ELSE
518     CALL KPP_CALC_DUMMY(
519     I bi, bj, myTime, myThid )
520     #endif /* ALLOW_AUTODIFF_TAMC */
521     ENDIF
522    
523     #endif /* ALLOW_KPP */
524    
525 mlosch 1.6 #ifdef ALLOW_PP81
526     C-- Compute PP81 mixing coefficients
527     IF (usePP81) THEN
528     #ifdef ALLOW_DEBUG
529     IF ( debugLevel .GE. debLevB )
530     & CALL DEBUG_CALL('PP81_CALC',myThid)
531     #endif
532     CALL PP81_CALC(
533     I bi, bj, myTime, myThid )
534     ENDIF
535     #endif /* ALLOW_PP81 */
536    
537     #ifdef ALLOW_MY82
538     C-- Compute MY82 mixing coefficients
539     IF (useMY82) THEN
540     #ifdef ALLOW_DEBUG
541     IF ( debugLevel .GE. debLevB )
542     & CALL DEBUG_CALL('MY82_CALC',myThid)
543     #endif
544     CALL MY82_CALC(
545     I bi, bj, myTime, myThid )
546     ENDIF
547     #endif /* ALLOW_MY82 */
548    
549 mlosch 1.9 #ifdef ALLOW_GGL90
550     C-- Compute GGL90 mixing coefficients
551     IF (useGGL90) THEN
552     #ifdef ALLOW_DEBUG
553     IF ( debugLevel .GE. debLevB )
554     & CALL DEBUG_CALL('GGL90_CALC',myThid)
555     #endif
556     CALL GGL90_CALC(
557     I bi, bj, myTime, myThid )
558     ENDIF
559     #endif /* ALLOW_GGL90 */
560    
561 jmc 1.20 #ifdef ALLOW_TIMEAVE
562     IF ( taveFreq.GT. 0. _d 0 .AND. fluidIsWater ) THEN
563     CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid)
564     ENDIF
565     IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
566     CALL TIMEAVE_CUMULATE(ConvectCountTave, IVDConvCount,
567     I Nr, deltaTclock, bi, bj, myThid)
568     ENDIF
569     #endif /* ALLOW_TIMEAVE */
570    
571 jmc 1.1 C-- end bi,bj loops.
572     ENDDO
573     ENDDO
574    
575 jmc 1.18 #ifdef ALLOW_DIAGNOSTICS
576     IF ( fluidIsWater .AND. useDiagnostics ) THEN
577     CALL DIAGS_OCEANIC_SURF_FLUX( myTime, myIter, myThid )
578     ENDIF
579 jmc 1.19 IF ( ivdc_kappa.NE.0 .AND. useDiagnostics ) THEN
580     CALL DIAGNOSTICS_FILL( IVDConvCount,'CONVADJ ',
581     & 0, Nr, 0, 1, 1, myThid )
582     ENDIF
583 jmc 1.18 #endif
584    
585 jmc 1.1 #ifdef ALLOW_DEBUG
586 jmc 1.29 IF ( debugLevel .GE. debLevB )
587     & CALL DEBUG_LEAVE('DO_OCEANIC_PHYS',myThid)
588 jmc 1.1 #endif
589    
590     RETURN
591     END

  ViewVC Help
Powered by ViewVC 1.1.22