/[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.44 - (hide annotations) (download)
Thu May 3 21:40:48 2007 UTC (17 years, 1 month ago) by jmc
Branch: MAIN
Changes since 1.43: +3 -3 lines
add myIter to KPP_CALC call.

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

  ViewVC Help
Powered by ViewVC 1.1.22