/[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.55 - (hide annotations) (download)
Sat Aug 18 21:51:10 2007 UTC (16 years, 9 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint59g, checkpoint59f
Changes since 1.54: +4 -1 lines
Update NLFS adjoint.

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

  ViewVC Help
Powered by ViewVC 1.1.22