/[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.31 - (hide annotations) (download)
Sun Oct 22 01:11:44 2006 UTC (17 years, 7 months ago) by heimbach
Branch: MAIN
Changes since 1.30: +11 -2 lines
Split seaice_init into _fixed, _varia

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

  ViewVC Help
Powered by ViewVC 1.1.22