/[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.41 - (hide annotations) (download)
Wed Apr 18 19:54:06 2007 UTC (17 years, 1 month ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint59
Changes since 1.40: +1 -2 lines
Removing exf_clim code.

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

  ViewVC Help
Powered by ViewVC 1.1.22