/[MITgcm]/MITgcm_contrib/MPMice/beaufort/code/do_oceanic_phys.F
ViewVC logotype

Annotation of /MITgcm_contrib/MPMice/beaufort/code/do_oceanic_phys.F

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


Revision 1.2 - (hide annotations) (download)
Sun May 31 03:43:51 2009 UTC (16 years, 2 months ago) by dimitri
Branch: MAIN
Changes since 1.1: +30 -24 lines
Updating subroutine for consistency with checkpoint59q

1 dimitri 1.2 C $Header: /u/gcmpack/MITgcm/model/src/do_oceanic_phys.F,v 1.65 2008/05/03 00:07:08 heimbach Exp $
2 dimitri 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     # ifdef ALLOW_SEAICE
15     # include "SEAICE_OPTIONS.h"
16     # endif
17     #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     C | SUBROUTINE DO_OCEANIC_PHYS
26     C | o Controlling routine for oceanic physics and
27     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     #ifdef ALLOW_TIMEAVE
42     #include "TIMEAVE_STATV.h"
43     #endif
44     #if defined (ALLOW_BALANCE_FLUXES) && !(defined ALLOW_AUTODIFF_TAMC)
45     #include "FFIELDS.h"
46     #endif
47    
48     #ifdef ALLOW_AUTODIFF_TAMC
49     # include "tamc.h"
50     # include "tamc_keys.h"
51     # include "FFIELDS.h"
52     # include "SURFACE.h"
53     # 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     # ifdef ALLOW_EXF
64     # include "ctrl.h"
65     # include "EXF_FIELDS.h"
66     # ifdef ALLOW_BULKFORMULAE
67     # include "EXF_CONSTANTS.h"
68     # endif
69     # endif
70     # ifdef ALLOW_SEAICE
71     # include "SEAICE.h"
72     # endif
73     #endif /* ALLOW_AUTODIFF_TAMC */
74    
75     C !INPUT/OUTPUT PARAMETERS:
76     C == Routine arguments ==
77     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     _RL myTime
81     INTEGER myIter
82     INTEGER myThid
83    
84     C !LOCAL VARIABLES:
85     C == Local variables
86     C rhoK, rhoKm1 :: Density at current level, and level above
87     C iMin, iMax :: Ranges and sub-block indices on which calculations
88     C jMin, jMax are applied.
89     C bi, bj :: tile indices
90     C i,j,k :: loop indices
91     _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     _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     INTEGER i, j, k
101     INTEGER doDiagsRho
102     #ifdef ALLOW_DIAGNOSTICS
103     LOGICAL DIAGNOSTICS_IS_ON
104     EXTERNAL DIAGNOSTICS_IS_ON
105     #endif /* ALLOW_DIAGNOSTICS */
106    
107     CEOP
108    
109     #ifdef ALLOW_AUTODIFF_TAMC
110     C-- dummy statement to end declaration part
111     itdkey = 1
112     #endif /* ALLOW_AUTODIFF_TAMC */
113    
114     #ifdef ALLOW_DEBUG
115     IF ( debugLevel .GE. debLevB )
116     & CALL DEBUG_ENTER('DO_OCEANIC_PHYS',myThid)
117     #endif
118    
119     doDiagsRho = 0
120     #ifdef ALLOW_DIAGNOSTICS
121     IF ( useDiagnostics .AND. fluidIsWater ) THEN
122     IF ( DIAGNOSTICS_IS_ON('RHOANOSQ',myThid) .OR.
123     & DIAGNOSTICS_IS_ON('URHOMASS',myThid) .OR.
124     & DIAGNOSTICS_IS_ON('VRHOMASS',myThid) .OR.
125     & DIAGNOSTICS_IS_ON('WRHOMASS',myThid) .OR.
126     & DIAGNOSTICS_IS_ON('WRHOMASS',myThid) ) doDiagsRho = 2
127     IF ( doDiagsRho.EQ.0 .AND.
128     & DIAGNOSTICS_IS_ON('MXLDEPTH',myThid) ) doDiagsRho = 1
129     IF ( doDiagsRho.EQ.0 .AND.
130     & DIAGNOSTICS_IS_ON('DRHODR ',myThid) ) doDiagsRho = 1
131     ENDIF
132     #endif /* ALLOW_DIAGNOSTICS */
133    
134     #ifdef ALLOW_CPL_MPMICE
135     CALL CPL_MPMICE( myTime, myIter, myThid )
136     #endif /* ALLOW_CPL_MPMICE */
137    
138     #ifdef ALLOW_SEAICE
139     IF ( useSEAICE ) THEN
140 dimitri 1.2 # ifdef ALLOW_AUTODIFF_TAMC
141     cph-adj-test(
142     CADJ STORE area,empmr,qsw,theta = comlev1, key = ikey_dynamics
143     cph-adj-test)
144 dimitri 1.1 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     CADJ STORE evap = comlev1, key = ikey_dynamics
148     cph# endif
149     CADJ STORE uvel,vvel = comlev1, key = ikey_dynamics
150 dimitri 1.2 # 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 dimitri 1.1 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 dimitri 1.2 # endif
158     # endif
159     # ifdef SEAICE_SALINITY
160 dimitri 1.1 CADJ STORE salt = comlev1, key = ikey_dynamics
161 dimitri 1.2 # endif
162     # ifdef ATMOSPHERIC_LOADING
163     CADJ STORE pload = comlev1, key = ikey_dynamics
164 dimitri 1.1 CADJ STORE siceload = comlev1, key = ikey_dynamics
165 dimitri 1.2 # endif
166     # ifdef NONLIN_FRSURF
167 dimitri 1.1 CADJ STORE recip_hfacc = comlev1, key = ikey_dynamics
168 dimitri 1.2 # endif
169 dimitri 1.1 # endif
170 dimitri 1.2 # ifdef ALLOW_DEBUG
171 dimitri 1.1 IF ( debugLevel .GE. debLevB )
172     & CALL DEBUG_CALL('SEAICE_MODEL',myThid)
173 dimitri 1.2 # endif
174 dimitri 1.1 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 dimitri 1.2 # ifdef ALLOW_COST
178 dimitri 1.1 CALL SEAICE_COST_SENSI ( myTime, myIter, myThid )
179 dimitri 1.2 # endif
180 dimitri 1.1 ENDIF
181     #endif /* ALLOW_SEAICE */
182    
183 dimitri 1.2 #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 dimitri 1.1 #if (defined ALLOW_THSICE) && !(defined ALLOW_ATM2D)
192     IF ( useThSIce .AND. fluidIsWater ) THEN
193     #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     #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     C compute temperature and (virtual) salt flux at the
212     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     C-- Freeze water at the surface
222     #ifdef ALLOW_AUTODIFF_TAMC
223     CADJ STORE theta = comlev1, key = ikey_dynamics
224     #endif
225     IF ( allowFreezing
226     & .AND. .NOT. useSEAICE
227     & .AND. .NOT. useThSIce ) THEN
228     CALL FREEZE_SURFACE( myTime, myIter, myThid )
229     ENDIF
230    
231     #ifdef ALLOW_OCN_COMPON_INTERF
232     C-- Apply imported data (from coupled interface) to forcing fields
233     C jmc: do not know precisely where to put this call (bf or af thSIce ?)
234     IF ( useCoupler ) THEN
235     CALL OCN_APPLY_IMPORT( .TRUE., myTime, myIter, myThid )
236     ENDIF
237     #endif /* ALLOW_OCN_COMPON_INTERF */
238    
239     #ifdef ALLOW_BALANCE_FLUXES
240     C balance fluxes
241     IF ( balanceEmPmR )
242     & CALL REMOVE_MEAN_RS( 1, EmPmR, maskH, maskH, rA, drF,
243     & 'EmPmR', myTime, myThid )
244     IF ( balanceQnet )
245     & CALL REMOVE_MEAN_RS( 1, Qnet, maskH, maskH, rA, drF,
246     & 'Qnet ', myTime, myThid )
247     #endif /* ALLOW_BALANCE_FLUXES */
248    
249     #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     C-- HPF directive to help TAMC
256     CHPF$ INDEPENDENT
257     #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     #else /* ALLOW_AUTODIFF_TAMC */
272     C if fluid is not water, by-pass find_rho, gmredi, surfaceForcing
273     C and all vertical mixing schemes, but keep OBCS_CALC
274     IF ( fluidIsWater ) THEN
275     #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     DO j=1-OLy,sNy+OLy
284     DO i=1-OLx,sNx+OLx
285     rhoK (i,j) = 0. _d 0
286     rhoKm1 (i,j) = 0. _d 0
287     rhoKp1 (i,j) = 0. _d 0
288     ENDDO
289     ENDDO
290    
291     DO k=1,Nr
292     DO j=1-OLy,sNy+OLy
293     DO i=1-OLx,sNx+OLx
294     C This is currently also used by IVDC and Diagnostics
295     sigmaX(i,j,k) = 0. _d 0
296     sigmaY(i,j,k) = 0. _d 0
297     sigmaR(i,j,k) = 0. _d 0
298     #ifdef ALLOW_AUTODIFF_TAMC
299     cph all the following init. are necessary for TAF
300     cph although some of these are re-initialised later.
301     IVDConvCount(i,j,k,bi,bj) = 0.
302     # ifdef ALLOW_GMREDI
303     Kwx(i,j,k,bi,bj) = 0. _d 0
304     Kwy(i,j,k,bi,bj) = 0. _d 0
305     Kwz(i,j,k,bi,bj) = 0. _d 0
306     # ifdef GM_NON_UNITY_DIAGONAL
307     Kux(i,j,k,bi,bj) = 0. _d 0
308     Kvy(i,j,k,bi,bj) = 0. _d 0
309     # endif
310     # ifdef GM_EXTRA_DIAGONAL
311     Kuz(i,j,k,bi,bj) = 0. _d 0
312     Kvz(i,j,k,bi,bj) = 0. _d 0
313     # endif
314     # ifdef GM_BOLUS_ADVEC
315     GM_PsiX(i,j,k,bi,bj) = 0. _d 0
316     GM_PsiY(i,j,k,bi,bj) = 0. _d 0
317     # endif
318     # ifdef GM_VISBECK_VARIABLE_K
319     VisbeckK(i,j,bi,bj) = 0. _d 0
320     # endif
321     # endif /* ALLOW_GMREDI */
322     # ifdef ALLOW_KPP
323     KPPdiffKzS(i,j,k,bi,bj) = 0. _d 0
324     KPPdiffKzT(i,j,k,bi,bj) = 0. _d 0
325     # endif /* ALLOW_KPP */
326     #endif /* ALLOW_AUTODIFF_TAMC */
327     ENDDO
328     ENDDO
329     ENDDO
330    
331     iMin = 1-OLx
332     iMax = sNx+OLx
333     jMin = 1-OLy
334     jMax = sNy+OLy
335    
336     #ifdef ALLOW_AUTODIFF_TAMC
337     CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
338     CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
339     CADJ STORE totphihyd(:,:,:,bi,bj)
340     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
341     # ifdef ALLOW_KPP
342     CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
343     CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
344     # endif
345     #endif /* ALLOW_AUTODIFF_TAMC */
346    
347     #ifdef ALLOW_DEBUG
348     IF ( debugLevel .GE. debLevB )
349     & CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid)
350     #endif
351    
352     C-- Start of diagnostic loop
353     DO k=Nr,1,-1
354    
355     #ifdef ALLOW_AUTODIFF_TAMC
356     C? Patrick, is this formula correct now that we change the loop range?
357     C? Do we still need this?
358     cph kkey formula corrected.
359     cph Needed for rhoK, rhoKm1, in the case useGMREDI.
360     kkey = (itdkey-1)*Nr + k
361     #endif /* ALLOW_AUTODIFF_TAMC */
362    
363     C-- Calculate gradients of potential density for isoneutral
364     C slope terms (e.g. GM/Redi tensor or IVDC diffusivity)
365     IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.)
366     & .OR. useSALT_PLUME .OR. doDiagsRho.GE.1 ) THEN
367     #ifdef ALLOW_DEBUG
368     IF ( debugLevel .GE. debLevB )
369     & CALL DEBUG_CALL('FIND_RHO',myThid)
370     #endif
371     #ifdef ALLOW_AUTODIFF_TAMC
372     CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
373     CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
374     #endif /* ALLOW_AUTODIFF_TAMC */
375     CALL FIND_RHO(
376     I bi, bj, iMin, iMax, jMin, jMax, k, k,
377     I theta, salt,
378     O rhoK,
379     I myThid )
380    
381     IF (k.GT.1) THEN
382     #ifdef ALLOW_AUTODIFF_TAMC
383     CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
384     CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte
385     #endif /* ALLOW_AUTODIFF_TAMC */
386     CALL FIND_RHO(
387     I bi, bj, iMin, iMax, jMin, jMax, k-1, k,
388     I theta, salt,
389     O rhoKm1,
390     I myThid )
391     ENDIF
392     #ifdef ALLOW_DEBUG
393     IF ( debugLevel .GE. debLevB )
394     & CALL DEBUG_CALL('GRAD_SIGMA',myThid)
395     #endif
396     cph Avoid variable aliasing for adjoint !!!
397     DO j=jMin,jMax
398     DO i=iMin,iMax
399     rhoKp1(i,j) = rhoK(i,j)
400     ENDDO
401     ENDDO
402     CALL GRAD_SIGMA(
403     I bi, bj, iMin, iMax, jMin, jMax, k,
404     I rhoK, rhoKm1, rhoKp1,
405     O sigmaX, sigmaY, sigmaR,
406     I myThid )
407     ENDIF
408    
409     C-- Implicit Vertical Diffusion for Convection
410     c ==> should use sigmaR !!!
411     IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN
412     #ifdef ALLOW_DEBUG
413     IF ( debugLevel .GE. debLevB )
414     & CALL DEBUG_CALL('CALC_IVDC',myThid)
415     #endif
416     CALL CALC_IVDC(
417     I bi, bj, iMin, iMax, jMin, jMax, k,
418     I rhoKm1, rhoK,
419     I myTime, myIter, myThid)
420     ENDIF
421    
422     #ifdef ALLOW_DIAGNOSTICS
423     IF ( doDiagsRho.GE.2 ) THEN
424     CALL DIAGS_RHO( k, bi, bj,
425     I rhoK, rhoKm1,
426     I myTime, myIter, myThid)
427     ENDIF
428     #endif
429    
430     C-- end of diagnostic k loop (Nr:1)
431     ENDDO
432    
433     #ifdef ALLOW_AUTODIFF_TAMC
434     CADJ STORE IVDConvCount(:,:,:,bi,bj)
435     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
436     #endif
437    
438     C-- Diagnose Mixed Layer Depth:
439     IF ( useGMRedi .OR. doDiagsRho.GE.1 ) THEN
440     CALL CALC_OCE_MXLAYER( rhoK, sigmaR,
441     & bi, bj, myTime, myIter, myThid )
442     ENDIF
443    
444     #ifdef ALLOW_SALT_PLUME
445     IF ( useSALT_PLUME ) THEN
446     CALL SALT_PLUME_CALC_DEPTH( rhoK, sigmaR,
447     & bi, bj, myTime, myIter, myThid )
448     ENDIF
449     #endif /* ALLOW_SALT_PLUME */
450    
451     #ifdef ALLOW_DIAGNOSTICS
452     IF ( doDiagsRho.GE.1 ) THEN
453     CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR ', 0, Nr,
454     & 2, bi, bj, myThid)
455     ENDIF
456     #endif /* ALLOW_DIAGNOSTICS */
457    
458     C-- Determines forcing terms based on external fields
459     C relaxation terms, etc.
460     #ifdef ALLOW_DEBUG
461     IF ( debugLevel .GE. debLevB )
462     & CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid)
463     #endif
464     #ifdef ALLOW_AUTODIFF_TAMC
465     CADJ STORE EmPmR(:,:,bi,bj)
466     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
467     # ifdef EXACT_CONSERV
468     CADJ STORE PmEpR(:,:,bi,bj)
469     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
470     # endif
471     # ifdef NONLIN_FRSURF
472     CADJ STORE hFac_surfC(:,:,bi,bj)
473     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
474     CADJ STORE recip_hFacC(:,:,:,bi,bj)
475     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
476     # endif
477     #endif
478     CALL EXTERNAL_FORCING_SURF(
479     I bi, bj, iMin, iMax, jMin, jMax,
480     I myTime, myIter, myThid )
481     #ifdef ALLOW_AUTODIFF_TAMC
482     # ifdef EXACT_CONSERV
483     cph-test
484     cphCADJ STORE PmEpR(:,:,bi,bj)
485     cphCADJ & = comlev1_bibj, key=itdkey, byte=isbyte
486     # endif
487     #endif
488    
489     #ifdef ALLOW_AUTODIFF_TAMC
490     cph needed for KPP
491     CADJ STORE surfaceForcingU(:,:,bi,bj)
492     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
493     CADJ STORE surfaceForcingV(:,:,bi,bj)
494     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
495     CADJ STORE surfaceForcingS(:,:,bi,bj)
496     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
497     CADJ STORE surfaceForcingT(:,:,bi,bj)
498     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
499     CADJ STORE surfaceForcingTice(:,:,bi,bj)
500     CADJ & = comlev1_bibj, key=itdkey, byte=isbyte
501     #endif /* ALLOW_AUTODIFF_TAMC */
502    
503     #ifdef ALLOW_KPP
504     C-- Compute KPP mixing coefficients
505     IF (useKPP) THEN
506     #ifdef ALLOW_DEBUG
507     IF ( debugLevel .GE. debLevB )
508     & CALL DEBUG_CALL('KPP_CALC',myThid)
509     #endif
510     CALL KPP_CALC(
511     I bi, bj, myTime, myIter, myThid )
512     #ifdef ALLOW_AUTODIFF_TAMC
513     ELSE
514     CALL KPP_CALC_DUMMY(
515     I bi, bj, myTime, myIter, myThid )
516     #endif /* ALLOW_AUTODIFF_TAMC */
517     ENDIF
518    
519     #endif /* ALLOW_KPP */
520    
521     #ifdef ALLOW_PP81
522     C-- Compute PP81 mixing coefficients
523     IF (usePP81) THEN
524     #ifdef ALLOW_DEBUG
525     IF ( debugLevel .GE. debLevB )
526     & CALL DEBUG_CALL('PP81_CALC',myThid)
527     #endif
528     CALL PP81_CALC(
529     I bi, bj, myTime, myThid )
530     ENDIF
531     #endif /* ALLOW_PP81 */
532    
533     #ifdef ALLOW_MY82
534     C-- Compute MY82 mixing coefficients
535     IF (useMY82) THEN
536     #ifdef ALLOW_DEBUG
537     IF ( debugLevel .GE. debLevB )
538     & CALL DEBUG_CALL('MY82_CALC',myThid)
539     #endif
540     CALL MY82_CALC(
541     I bi, bj, myTime, myThid )
542     ENDIF
543     #endif /* ALLOW_MY82 */
544    
545     #ifdef ALLOW_GGL90
546     C-- Compute GGL90 mixing coefficients
547     IF (useGGL90) THEN
548     #ifdef ALLOW_DEBUG
549     IF ( debugLevel .GE. debLevB )
550     & CALL DEBUG_CALL('GGL90_CALC',myThid)
551     #endif
552     CALL GGL90_CALC(
553     I bi, bj, myTime, myThid )
554     ENDIF
555     #endif /* ALLOW_GGL90 */
556    
557     #ifdef ALLOW_TIMEAVE
558     IF ( taveFreq.GT. 0. _d 0 ) THEN
559     CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid)
560     ENDIF
561     IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN
562     CALL TIMEAVE_CUMULATE(ConvectCountTave, IVDConvCount,
563     I Nr, deltaTclock, bi, bj, myThid)
564     ENDIF
565     #endif /* ALLOW_TIMEAVE */
566    
567     #ifdef ALLOW_GMREDI
568     #ifdef ALLOW_AUTODIFF_TAMC
569     # ifndef GM_EXCLUDE_CLIPPING
570     cph storing here is needed only for one GMREDI_OPTIONS:
571     cph define GM_BOLUS_ADVEC
572     cph keep it although TAF says you dont need to.
573     cph but I've avoided the #ifdef for now, in case more things change
574     CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
575     CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
576     CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte
577     # endif
578     #endif /* ALLOW_AUTODIFF_TAMC */
579    
580     C-- Calculate iso-neutral slopes for the GM/Redi parameterisation
581     IF (useGMRedi) THEN
582     #ifdef ALLOW_DEBUG
583     IF ( debugLevel .GE. debLevB )
584     & CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid)
585     #endif
586     CALL GMREDI_CALC_TENSOR(
587     c I bi, bj, iMin, iMax, jMin, jMax,
588     c I sigmaX, sigmaY, sigmaR,
589     c I myThid )
590     I iMin, iMax, jMin, jMax,
591     I sigmaX, sigmaY, sigmaR,
592     I bi, bj, myTime, myIter, myThid )
593     #ifdef ALLOW_AUTODIFF_TAMC
594     ELSE
595     CALL GMREDI_CALC_TENSOR_DUMMY(
596     c I bi, bj, iMin, iMax, jMin, jMax,
597     c I sigmaX, sigmaY, sigmaR,
598     c I myThid )
599     I iMin, iMax, jMin, jMax,
600     I sigmaX, sigmaY, sigmaR,
601     I bi, bj, myTime, myIter, myThid )
602     #endif /* ALLOW_AUTODIFF_TAMC */
603     ENDIF
604     #endif /* ALLOW_GMREDI */
605    
606     #ifndef ALLOW_AUTODIFF_TAMC
607     C--- if fluid Is Water: end
608     ENDIF
609     #endif
610    
611     #ifdef ALLOW_OBCS
612     C-- Calculate future values on open boundaries
613     IF (useOBCS) THEN
614     #ifdef ALLOW_DEBUG
615     IF ( debugLevel .GE. debLevB )
616     & CALL DEBUG_CALL('OBCS_CALC',myThid)
617     #endif
618     CALL OBCS_CALC( bi, bj, myTime+deltaTclock, myIter+1,
619     I uVel, vVel, wVel, theta, salt,
620     I myThid )
621     ENDIF
622     #endif /* ALLOW_OBCS */
623    
624     C-- end bi,bj loops.
625     ENDDO
626     ENDDO
627    
628     #ifdef ALLOW_KPP
629     IF (useKPP) THEN
630     CALL KPP_DO_EXCH( myThid )
631     ENDIF
632     #endif /* ALLOW_KPP */
633    
634     #ifdef ALLOW_DIAGNOSTICS
635     IF ( fluidIsWater .AND. useDiagnostics ) THEN
636     CALL DIAGS_OCEANIC_SURF_FLUX( myTime, myIter, myThid )
637     ENDIF
638     IF ( ivdc_kappa.NE.0 .AND. useDiagnostics ) THEN
639     CALL DIAGNOSTICS_FILL( IVDConvCount,'CONVADJ ',
640     & 0, Nr, 0, 1, 1, myThid )
641     ENDIF
642     #endif
643    
644     #ifdef ALLOW_DEBUG
645     IF ( debugLevel .GE. debLevB )
646     & CALL DEBUG_LEAVE('DO_OCEANIC_PHYS',myThid)
647     #endif
648    
649     RETURN
650     END

  ViewVC Help
Powered by ViewVC 1.1.22