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