/[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.48 - (hide annotations) (download)
Sun Jun 3 23:56:53 2007 UTC (17 years ago) by jmc
Branch: MAIN
CVS Tags: checkpoint59c
Changes since 1.47: +3 -1 lines
forgot the #ifndef ALLOW_AUTODIFF_TAMC in the previous check-in

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

  ViewVC Help
Powered by ViewVC 1.1.22