1 |
C $Header: /u/gcmpack/MITgcm/model/src/do_oceanic_phys.F,v 1.27 2006/03/08 06:36:39 heimbach Exp $ |
2 |
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 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
15 |
|
16 |
CBOP |
17 |
C !ROUTINE: DO_OCEANIC_PHYS |
18 |
C !INTERFACE: |
19 |
SUBROUTINE DO_OCEANIC_PHYS(myTime, myIter, myThid) |
20 |
C !DESCRIPTION: \bv |
21 |
C *==========================================================* |
22 |
C | SUBROUTINE DO_OCEANIC_PHYS |
23 |
C | o Controlling routine for oceanic physics and |
24 |
C | parameterization |
25 |
C *==========================================================* |
26 |
C | o originally, part of S/R thermodynamics |
27 |
C *==========================================================* |
28 |
C \ev |
29 |
|
30 |
C !USES: |
31 |
IMPLICIT NONE |
32 |
C == Global variables === |
33 |
#include "SIZE.h" |
34 |
#include "EEPARAMS.h" |
35 |
#include "PARAMS.h" |
36 |
#include "DYNVARS.h" |
37 |
#include "GRID.h" |
38 |
#ifdef ALLOW_TIMEAVE |
39 |
#include "TIMEAVE_STATV.h" |
40 |
#endif |
41 |
#if defined (ALLOW_BALANCE_FLUXES) && !(defined ALLOW_AUTODIFF_TAMC) |
42 |
#include "FFIELDS.h" |
43 |
#endif |
44 |
|
45 |
#ifdef ALLOW_AUTODIFF_TAMC |
46 |
# include "tamc.h" |
47 |
# include "tamc_keys.h" |
48 |
# include "FFIELDS.h" |
49 |
# include "EOS.h" |
50 |
# ifdef ALLOW_KPP |
51 |
# include "KPP.h" |
52 |
# endif |
53 |
# ifdef ALLOW_GMREDI |
54 |
# include "GMREDI.h" |
55 |
# endif |
56 |
# ifdef ALLOW_EBM |
57 |
# include "EBM.h" |
58 |
# endif |
59 |
# ifdef EXACT_CONSERV |
60 |
# include "SURFACE.h" |
61 |
# endif |
62 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
63 |
|
64 |
C !INPUT/OUTPUT PARAMETERS: |
65 |
C == Routine arguments == |
66 |
C myTime :: Current time in simulation |
67 |
C myIter :: Current iteration number in simulation |
68 |
C myThid :: Thread number for this instance of the routine. |
69 |
_RL myTime |
70 |
INTEGER myIter |
71 |
INTEGER myThid |
72 |
|
73 |
C !LOCAL VARIABLES: |
74 |
C == Local variables |
75 |
C rhoK, rhoKM1 :: Density at current level, and level above |
76 |
C iMin, iMax :: Ranges and sub-block indices on which calculations |
77 |
C jMin, jMax are applied. |
78 |
C bi, bj :: tile indices |
79 |
C i,j,k :: loop indices |
80 |
_RL rhokm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
81 |
_RL rhok (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
82 |
_RL sigmaX (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
83 |
_RL sigmaY (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
84 |
_RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
85 |
INTEGER iMin, iMax |
86 |
INTEGER jMin, jMax |
87 |
INTEGER bi, bj |
88 |
INTEGER i, j, k |
89 |
INTEGER doDiagsRho |
90 |
#ifdef ALLOW_DIAGNOSTICS |
91 |
LOGICAL DIAGNOSTICS_IS_ON |
92 |
EXTERNAL DIAGNOSTICS_IS_ON |
93 |
#endif /* ALLOW_DIAGNOSTICS */ |
94 |
|
95 |
CEOP |
96 |
|
97 |
#ifdef ALLOW_AUTODIFF_TAMC |
98 |
C-- dummy statement to end declaration part |
99 |
itdkey = 1 |
100 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
101 |
|
102 |
#ifdef ALLOW_DEBUG |
103 |
IF ( debugLevel .GE. debLevB ) |
104 |
& CALL DEBUG_ENTER('DO_OCEANIC_PHYS',myThid) |
105 |
#endif |
106 |
|
107 |
doDiagsRho = 0 |
108 |
#ifdef ALLOW_DIAGNOSTICS |
109 |
IF ( useDiagnostics .AND. fluidIsWater ) THEN |
110 |
IF ( DIAGNOSTICS_IS_ON('DRHODR ',myThid) ) doDiagsRho = 1 |
111 |
IF ( DIAGNOSTICS_IS_ON('RHOANOSQ',myThid) .OR. |
112 |
& DIAGNOSTICS_IS_ON('URHOMASS',myThid) .OR. |
113 |
& DIAGNOSTICS_IS_ON('VRHOMASS',myThid) .OR. |
114 |
& DIAGNOSTICS_IS_ON('WRHOMASS',myThid) .OR. |
115 |
& DIAGNOSTICS_IS_ON('WRHOMASS',myThid) ) doDiagsRho = 2 |
116 |
ENDIF |
117 |
#endif /* ALLOW_DIAGNOSTICS */ |
118 |
|
119 |
#ifdef ALLOW_THSICE |
120 |
IF ( useThSIce .AND. fluidIsWater ) THEN |
121 |
#ifdef ALLOW_DEBUG |
122 |
IF ( debugLevel .GE. debLevB ) |
123 |
& CALL DEBUG_CALL('THSICE_MAIN',myThid) |
124 |
#endif |
125 |
C-- Step forward Therm.Sea-Ice variables |
126 |
C and modify forcing terms including effects from ice |
127 |
CALL TIMER_START('THSICE_MAIN [DO_OCEANIC_PHYS]', myThid) |
128 |
CALL THSICE_MAIN( myTime, myIter, myThid ) |
129 |
CALL TIMER_STOP( 'THSICE_MAIN [DO_OCEANIC_PHYS]', myThid) |
130 |
ENDIF |
131 |
#endif /* ALLOW_THSICE */ |
132 |
|
133 |
#ifdef ALLOW_SHELFICE |
134 |
IF ( useShelfIce .AND. fluidIsWater ) THEN |
135 |
#ifdef ALLOW_DEBUG |
136 |
IF ( debugLevel .GE. debLevB ) |
137 |
& CALL DEBUG_CALL('SHELFICE_THERMODYNAMICS',myThid) |
138 |
#endif |
139 |
C compute temperature and (virtual) salt flux at the |
140 |
C shelf-ice ocean interface |
141 |
CALL TIMER_START('SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]', |
142 |
& myThid) |
143 |
CALL SHELFICE_THERMODYNAMICS( myTime, myIter, myThid ) |
144 |
CALL TIMER_STOP( 'SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]', |
145 |
& myThid) |
146 |
ENDIF |
147 |
#endif /* ALLOW_SHELFICE */ |
148 |
|
149 |
C-- Freeze water at the surface |
150 |
#ifdef ALLOW_AUTODIFF_TAMC |
151 |
CADJ STORE theta = comlev1, key = ikey_dynamics |
152 |
#endif |
153 |
IF ( allowFreezing |
154 |
& .AND. .NOT. useSEAICE |
155 |
& .AND. .NOT. useThSIce ) THEN |
156 |
CALL FREEZE_SURFACE( myTime, myIter, myThid ) |
157 |
ENDIF |
158 |
|
159 |
#ifdef ALLOW_OCN_COMPON_INTERF |
160 |
C-- Apply imported data (from coupled interface) to forcing fields |
161 |
C jmc: do not know precisely where to put this call (bf or af thSIce ?) |
162 |
IF ( useCoupler ) THEN |
163 |
CALL OCN_APPLY_IMPORT( .TRUE., myTime, myIter, myThid ) |
164 |
ENDIF |
165 |
#endif /* ALLOW_OCN_COMPON_INTERF */ |
166 |
|
167 |
#ifdef ALLOW_BALANCE_FLUXES |
168 |
C balance fluxes |
169 |
IF ( balanceEmPmR ) |
170 |
& CALL REMOVE_MEAN_RS( 1, EmPmR, maskH, maskH, rA, drF, |
171 |
& 'EmPmR', myTime, myThid ) |
172 |
IF ( balanceQnet ) |
173 |
& CALL REMOVE_MEAN_RS( 1, Qnet, maskH, maskH, rA, drF, |
174 |
& 'Qnet ', myTime, myThid ) |
175 |
#endif /* ALLOW_BALANCE_FLUXES */ |
176 |
|
177 |
#ifdef ALLOW_AUTODIFF_TAMC |
178 |
C-- HPF directive to help TAMC |
179 |
CHPF$ INDEPENDENT |
180 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
181 |
DO bj=myByLo(myThid),myByHi(myThid) |
182 |
#ifdef ALLOW_AUTODIFF_TAMC |
183 |
C-- HPF directive to help TAMC |
184 |
CHPF$ INDEPENDENT |
185 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
186 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
187 |
|
188 |
#ifdef ALLOW_AUTODIFF_TAMC |
189 |
act1 = bi - myBxLo(myThid) |
190 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
191 |
act2 = bj - myByLo(myThid) |
192 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
193 |
act3 = myThid - 1 |
194 |
max3 = nTx*nTy |
195 |
act4 = ikey_dynamics - 1 |
196 |
itdkey = (act1 + 1) + act2*max1 |
197 |
& + act3*max1*max2 |
198 |
& + act4*max1*max2*max3 |
199 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
200 |
|
201 |
C-- Set up work arrays with valid (i.e. not NaN) values |
202 |
C These inital values do not alter the numerical results. They |
203 |
C just ensure that all memory references are to valid floating |
204 |
C point numbers. This prevents spurious hardware signals due to |
205 |
C uninitialised but inert locations. |
206 |
|
207 |
DO j=1-OLy,sNy+OLy |
208 |
DO i=1-OLx,sNx+OLx |
209 |
rhok (i,j) = 0. _d 0 |
210 |
rhoKM1 (i,j) = 0. _d 0 |
211 |
ENDDO |
212 |
ENDDO |
213 |
|
214 |
DO k=1,Nr |
215 |
DO j=1-OLy,sNy+OLy |
216 |
DO i=1-OLx,sNx+OLx |
217 |
C This is currently also used by IVDC and Diagnostics |
218 |
sigmaX(i,j,k) = 0. _d 0 |
219 |
sigmaY(i,j,k) = 0. _d 0 |
220 |
sigmaR(i,j,k) = 0. _d 0 |
221 |
#ifdef ALLOW_AUTODIFF_TAMC |
222 |
cph all the following init. are necessary for TAF |
223 |
cph although some of these are re-initialised later. |
224 |
IVDConvCount(i,j,k,bi,bj) = 0. |
225 |
# ifdef ALLOW_GMREDI |
226 |
Kwx(i,j,k,bi,bj) = 0. _d 0 |
227 |
Kwy(i,j,k,bi,bj) = 0. _d 0 |
228 |
Kwz(i,j,k,bi,bj) = 0. _d 0 |
229 |
# ifdef GM_NON_UNITY_DIAGONAL |
230 |
Kux(i,j,k,bi,bj) = 0. _d 0 |
231 |
Kvy(i,j,k,bi,bj) = 0. _d 0 |
232 |
# endif |
233 |
# ifdef GM_EXTRA_DIAGONAL |
234 |
Kuz(i,j,k,bi,bj) = 0. _d 0 |
235 |
Kvz(i,j,k,bi,bj) = 0. _d 0 |
236 |
# endif |
237 |
# ifdef GM_BOLUS_ADVEC |
238 |
GM_PsiX(i,j,k,bi,bj) = 0. _d 0 |
239 |
GM_PsiY(i,j,k,bi,bj) = 0. _d 0 |
240 |
# endif |
241 |
# ifdef GM_VISBECK_VARIABLE_K |
242 |
VisbeckK(i,j,bi,bj) = 0. _d 0 |
243 |
# endif |
244 |
# endif /* ALLOW_GMREDI */ |
245 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
246 |
ENDDO |
247 |
ENDDO |
248 |
ENDDO |
249 |
|
250 |
iMin = 1-OLx |
251 |
iMax = sNx+OLx |
252 |
jMin = 1-OLy |
253 |
jMax = sNy+OLy |
254 |
|
255 |
#ifdef ALLOW_AUTODIFF_TAMC |
256 |
CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte |
257 |
CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte |
258 |
CADJ STORE totphihyd(:,:,:,bi,bj) |
259 |
CADJ & = comlev1_bibj, key=itdkey, byte=isbyte |
260 |
# ifdef ALLOW_KPP |
261 |
CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte |
262 |
CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte |
263 |
# endif |
264 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
265 |
|
266 |
#ifdef ALLOW_DEBUG |
267 |
IF ( debugLevel .GE. debLevB ) |
268 |
& CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid) |
269 |
#endif |
270 |
|
271 |
C-- Start of diagnostic loop |
272 |
DO k=Nr,1,-1 |
273 |
|
274 |
#ifdef ALLOW_AUTODIFF_TAMC |
275 |
C? Patrick, is this formula correct now that we change the loop range? |
276 |
C? Do we still need this? |
277 |
cph kkey formula corrected. |
278 |
cph Needed for rhok, rhokm1, in the case useGMREDI. |
279 |
kkey = (itdkey-1)*Nr + k |
280 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
281 |
|
282 |
C-- Calculate gradients of potential density for isoneutral |
283 |
C slope terms (e.g. GM/Redi tensor or IVDC diffusivity) |
284 |
c IF ( k.GT.1 .AND. (useGMRedi.OR.ivdc_kappa.NE.0.) ) THEN |
285 |
IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) |
286 |
& .OR. doDiagsRho.GE.1 ) THEN |
287 |
#ifdef ALLOW_DEBUG |
288 |
IF ( debugLevel .GE. debLevB ) |
289 |
& CALL DEBUG_CALL('FIND_RHO',myThid) |
290 |
#endif |
291 |
#ifdef ALLOW_AUTODIFF_TAMC |
292 |
CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
293 |
CADJ STORE salt (:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
294 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
295 |
CALL FIND_RHO( |
296 |
I bi, bj, iMin, iMax, jMin, jMax, k, k, |
297 |
I theta, salt, |
298 |
O rhoK, |
299 |
I myThid ) |
300 |
|
301 |
IF (k.GT.1) THEN |
302 |
#ifdef ALLOW_AUTODIFF_TAMC |
303 |
CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
304 |
CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, byte=isbyte |
305 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
306 |
CALL FIND_RHO( |
307 |
I bi, bj, iMin, iMax, jMin, jMax, k-1, k, |
308 |
I theta, salt, |
309 |
O rhoKm1, |
310 |
I myThid ) |
311 |
ENDIF |
312 |
#ifdef ALLOW_DEBUG |
313 |
IF ( debugLevel .GE. debLevB ) |
314 |
& CALL DEBUG_CALL('GRAD_SIGMA',myThid) |
315 |
#endif |
316 |
CALL GRAD_SIGMA( |
317 |
I bi, bj, iMin, iMax, jMin, jMax, k, |
318 |
I rhoK, rhoKm1, rhoK, |
319 |
O sigmaX, sigmaY, sigmaR, |
320 |
I myThid ) |
321 |
ENDIF |
322 |
|
323 |
#ifdef ALLOW_AUTODIFF_TAMC |
324 |
ctest# ifndef GM_EXCLUDE_CLIPPING |
325 |
CADJ STORE rhok (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte |
326 |
ctest# endif |
327 |
CADJ STORE rhokm1 (:,:) = comlev1_bibj_k , key=kkey, byte=isbyte |
328 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
329 |
C-- Implicit Vertical Diffusion for Convection |
330 |
c ==> should use sigmaR !!! |
331 |
IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN |
332 |
#ifdef ALLOW_DEBUG |
333 |
IF ( debugLevel .GE. debLevB ) |
334 |
& CALL DEBUG_CALL('CALC_IVDC',myThid) |
335 |
#endif |
336 |
CALL CALC_IVDC( |
337 |
I bi, bj, iMin, iMax, jMin, jMax, k, |
338 |
I rhoKm1, rhoK, |
339 |
I myTime, myIter, myThid) |
340 |
ENDIF |
341 |
|
342 |
#ifdef ALLOW_DIAGNOSTICS |
343 |
IF ( doDiagsRho.GE.2 ) THEN |
344 |
CALL DIAGS_RHO( k, bi, bj, |
345 |
I rhoK, rhoKm1, |
346 |
I myTime, myIter, myThid) |
347 |
ENDIF |
348 |
#endif |
349 |
|
350 |
C-- end of diagnostic k loop (Nr:1) |
351 |
ENDDO |
352 |
|
353 |
#ifdef ALLOW_DIAGNOSTICS |
354 |
c IF ( useDiagnostics .AND. |
355 |
c & (useGMRedi .OR. ivdc_kappa.NE.0.) ) THEN |
356 |
IF ( doDiagsRho.GE.1 ) THEN |
357 |
CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR ', 0, Nr, |
358 |
& 2, bi, bj, myThid) |
359 |
ENDIF |
360 |
#endif |
361 |
|
362 |
#ifdef ALLOW_OBCS |
363 |
C-- Calculate future values on open boundaries |
364 |
IF (useOBCS) THEN |
365 |
#ifdef ALLOW_DEBUG |
366 |
IF ( debugLevel .GE. debLevB ) |
367 |
& CALL DEBUG_CALL('OBCS_CALC',myThid) |
368 |
#endif |
369 |
CALL OBCS_CALC( bi, bj, myTime+deltaTclock, myIter+1, |
370 |
I uVel, vVel, wVel, theta, salt, |
371 |
I myThid ) |
372 |
ENDIF |
373 |
#endif /* ALLOW_OBCS */ |
374 |
|
375 |
#ifndef ALLOW_AUTODIFF_TAMC |
376 |
IF ( fluidIsWater ) THEN |
377 |
#endif |
378 |
C-- Determines forcing terms based on external fields |
379 |
C relaxation terms, etc. |
380 |
#ifdef ALLOW_DEBUG |
381 |
IF ( debugLevel .GE. debLevB ) |
382 |
& CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid) |
383 |
#endif |
384 |
#ifdef ALLOW_AUTODIFF_TAMC |
385 |
CADJ STORE EmPmR(:,:,bi,bj) |
386 |
CADJ & = comlev1_bibj, key=itdkey, byte=isbyte |
387 |
# ifdef EXACT_CONSERV |
388 |
CADJ STORE PmEpR(:,:,bi,bj) |
389 |
CADJ & = comlev1_bibj, key=itdkey, byte=isbyte |
390 |
# endif |
391 |
# ifdef NONLIN_FRSURF |
392 |
CADJ STORE hFac_surfC(:,:,bi,bj) |
393 |
CADJ & = comlev1_bibj, key=itdkey, byte=isbyte |
394 |
CADJ STORE recip_hFacC(:,:,:,bi,bj) |
395 |
CADJ & = comlev1_bibj, key=itdkey, byte=isbyte |
396 |
# endif |
397 |
#endif |
398 |
CALL EXTERNAL_FORCING_SURF( |
399 |
I bi, bj, iMin, iMax, jMin, jMax, |
400 |
I myTime, myIter, myThid ) |
401 |
#ifndef ALLOW_AUTODIFF_TAMC |
402 |
ENDIF |
403 |
#endif |
404 |
#ifdef ALLOW_AUTODIFF_TAMC |
405 |
# ifdef EXACT_CONSERV |
406 |
cph-test |
407 |
cphCADJ STORE PmEpR(:,:,bi,bj) |
408 |
cphCADJ & = comlev1_bibj, key=itdkey, byte=isbyte |
409 |
# endif |
410 |
#endif |
411 |
|
412 |
#ifdef ALLOW_AUTODIFF_TAMC |
413 |
cph needed for KPP |
414 |
CADJ STORE surfaceForcingU(:,:,bi,bj) |
415 |
CADJ & = comlev1_bibj, key=itdkey, byte=isbyte |
416 |
CADJ STORE surfaceForcingV(:,:,bi,bj) |
417 |
CADJ & = comlev1_bibj, key=itdkey, byte=isbyte |
418 |
CADJ STORE surfaceForcingS(:,:,bi,bj) |
419 |
CADJ & = comlev1_bibj, key=itdkey, byte=isbyte |
420 |
CADJ STORE surfaceForcingT(:,:,bi,bj) |
421 |
CADJ & = comlev1_bibj, key=itdkey, byte=isbyte |
422 |
CADJ STORE surfaceForcingTice(:,:,bi,bj) |
423 |
CADJ & = comlev1_bibj, key=itdkey, byte=isbyte |
424 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
425 |
|
426 |
#ifdef ALLOW_GMREDI |
427 |
|
428 |
#ifdef ALLOW_AUTODIFF_TAMC |
429 |
# ifndef GM_EXCLUDE_CLIPPING |
430 |
cph storing here is needed only for one GMREDI_OPTIONS: |
431 |
cph define GM_BOLUS_ADVEC |
432 |
cph keep it although TAF says you dont need to. |
433 |
cph but I've avoided the #ifdef for now, in case more things change |
434 |
CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte |
435 |
CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte |
436 |
CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, byte=isbyte |
437 |
# endif |
438 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
439 |
|
440 |
C-- Calculate iso-neutral slopes for the GM/Redi parameterisation |
441 |
IF (useGMRedi) THEN |
442 |
#ifdef ALLOW_DEBUG |
443 |
IF ( debugLevel .GE. debLevB ) |
444 |
& CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid) |
445 |
#endif |
446 |
CALL GMREDI_CALC_TENSOR( |
447 |
I bi, bj, iMin, iMax, jMin, jMax, |
448 |
I sigmaX, sigmaY, sigmaR, |
449 |
I myThid ) |
450 |
#ifdef ALLOW_AUTODIFF_TAMC |
451 |
ELSE |
452 |
CALL GMREDI_CALC_TENSOR_DUMMY( |
453 |
I bi, bj, iMin, iMax, jMin, jMax, |
454 |
I sigmaX, sigmaY, sigmaR, |
455 |
I myThid ) |
456 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
457 |
ENDIF |
458 |
|
459 |
#endif /* ALLOW_GMREDI */ |
460 |
|
461 |
#ifdef ALLOW_KPP |
462 |
C-- Compute KPP mixing coefficients |
463 |
IF (useKPP) THEN |
464 |
#ifdef ALLOW_DEBUG |
465 |
IF ( debugLevel .GE. debLevB ) |
466 |
& CALL DEBUG_CALL('KPP_CALC',myThid) |
467 |
#endif |
468 |
CALL KPP_CALC( |
469 |
I bi, bj, myTime, myThid ) |
470 |
#ifdef ALLOW_AUTODIFF_TAMC |
471 |
ELSE |
472 |
CALL KPP_CALC_DUMMY( |
473 |
I bi, bj, myTime, myThid ) |
474 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
475 |
ENDIF |
476 |
|
477 |
#endif /* ALLOW_KPP */ |
478 |
|
479 |
#ifdef ALLOW_PP81 |
480 |
C-- Compute PP81 mixing coefficients |
481 |
IF (usePP81) THEN |
482 |
#ifdef ALLOW_DEBUG |
483 |
IF ( debugLevel .GE. debLevB ) |
484 |
& CALL DEBUG_CALL('PP81_CALC',myThid) |
485 |
#endif |
486 |
CALL PP81_CALC( |
487 |
I bi, bj, myTime, myThid ) |
488 |
ENDIF |
489 |
#endif /* ALLOW_PP81 */ |
490 |
|
491 |
#ifdef ALLOW_MY82 |
492 |
C-- Compute MY82 mixing coefficients |
493 |
IF (useMY82) THEN |
494 |
#ifdef ALLOW_DEBUG |
495 |
IF ( debugLevel .GE. debLevB ) |
496 |
& CALL DEBUG_CALL('MY82_CALC',myThid) |
497 |
#endif |
498 |
CALL MY82_CALC( |
499 |
I bi, bj, myTime, myThid ) |
500 |
ENDIF |
501 |
#endif /* ALLOW_MY82 */ |
502 |
|
503 |
#ifdef ALLOW_GGL90 |
504 |
C-- Compute GGL90 mixing coefficients |
505 |
IF (useGGL90) THEN |
506 |
#ifdef ALLOW_DEBUG |
507 |
IF ( debugLevel .GE. debLevB ) |
508 |
& CALL DEBUG_CALL('GGL90_CALC',myThid) |
509 |
#endif |
510 |
CALL GGL90_CALC( |
511 |
I bi, bj, myTime, myThid ) |
512 |
ENDIF |
513 |
#endif /* ALLOW_GGL90 */ |
514 |
|
515 |
#ifdef ALLOW_TIMEAVE |
516 |
IF ( taveFreq.GT. 0. _d 0 .AND. fluidIsWater ) THEN |
517 |
CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid) |
518 |
ENDIF |
519 |
IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN |
520 |
CALL TIMEAVE_CUMULATE(ConvectCountTave, IVDConvCount, |
521 |
I Nr, deltaTclock, bi, bj, myThid) |
522 |
ENDIF |
523 |
#endif /* ALLOW_TIMEAVE */ |
524 |
|
525 |
C-- end bi,bj loops. |
526 |
ENDDO |
527 |
ENDDO |
528 |
|
529 |
#ifdef ALLOW_DIAGNOSTICS |
530 |
IF ( fluidIsWater .AND. useDiagnostics ) THEN |
531 |
CALL DIAGS_OCEANIC_SURF_FLUX( myTime, myIter, myThid ) |
532 |
ENDIF |
533 |
IF ( ivdc_kappa.NE.0 .AND. useDiagnostics ) THEN |
534 |
CALL DIAGNOSTICS_FILL( IVDConvCount,'CONVADJ ', |
535 |
& 0, Nr, 0, 1, 1, myThid ) |
536 |
ENDIF |
537 |
#endif |
538 |
|
539 |
#ifdef ALLOW_DEBUG |
540 |
IF ( debugLevel .GE. debLevB ) |
541 |
& CALL DEBUG_LEAVE('DO_OCEANIC_PHYS',myThid) |
542 |
#endif |
543 |
|
544 |
RETURN |
545 |
END |