/[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.33 - (hide annotations) (download)
Wed Nov 8 16:01:04 2006 UTC (17 years, 6 months ago) by heimbach
Branch: MAIN
CVS Tags: mitgcm_mapl_00, checkpoint58r_post, checkpoint58s_post
Changes since 1.32: +4 -1 lines
Adding store.

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

  ViewVC Help
Powered by ViewVC 1.1.22