/[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.61 - (hide annotations) (download)
Wed Nov 28 09:26:16 2007 UTC (16 years, 6 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint59p, checkpoint59m, checkpoint59l, checkpoint59o, checkpoint59n, checkpoint59k
Changes since 1.60: +6 -5 lines
added check, readparms, and more diagnostice to pkg/salt_plume
also changed package flag from useSaltPlume to useSALT_PLUME

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

  ViewVC Help
Powered by ViewVC 1.1.22