1 |
C $Header: /u/gcmpack/MITgcm/model/src/do_oceanic_phys.F,v 1.76 2008/10/23 18:14:59 dfer 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 "GRID.h" |
40 |
#include "DYNVARS.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 |
# ifdef ALLOW_SALT_PLUME |
74 |
# include "SALT_PLUME.h" |
75 |
# endif |
76 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
77 |
|
78 |
C !INPUT/OUTPUT PARAMETERS: |
79 |
C == Routine arguments == |
80 |
C myTime :: Current time in simulation |
81 |
C myIter :: Current iteration number in simulation |
82 |
C myThid :: Thread number for this instance of the routine. |
83 |
_RL myTime |
84 |
INTEGER myIter |
85 |
INTEGER myThid |
86 |
|
87 |
C !LOCAL VARIABLES: |
88 |
C == Local variables |
89 |
C rhoK, rhoKm1 :: Density at current level, and level above |
90 |
C iMin, iMax :: Ranges and sub-block indices on which calculations |
91 |
C jMin, jMax are applied. |
92 |
C bi, bj :: tile indices |
93 |
C i,j,k :: loop indices |
94 |
_RL rhoKp1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
95 |
_RL rhoKm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
96 |
_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 |
INTEGER i, j, k |
103 |
INTEGER doDiagsRho |
104 |
#ifdef ALLOW_DIAGNOSTICS |
105 |
LOGICAL DIAGNOSTICS_IS_ON |
106 |
EXTERNAL DIAGNOSTICS_IS_ON |
107 |
#endif /* ALLOW_DIAGNOSTICS */ |
108 |
|
109 |
CEOP |
110 |
|
111 |
#ifdef ALLOW_AUTODIFF_TAMC |
112 |
C-- dummy statement to end declaration part |
113 |
itdkey = 1 |
114 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
115 |
|
116 |
#ifdef ALLOW_DEBUG |
117 |
IF ( debugLevel .GE. debLevB ) |
118 |
& CALL DEBUG_ENTER('DO_OCEANIC_PHYS',myThid) |
119 |
#endif |
120 |
|
121 |
doDiagsRho = 0 |
122 |
#ifdef ALLOW_DIAGNOSTICS |
123 |
IF ( useDiagnostics .AND. fluidIsWater ) THEN |
124 |
IF ( DIAGNOSTICS_IS_ON('WRHOMASS',myThid) ) |
125 |
& doDiagsRho = doDiagsRho + 1 |
126 |
IF ( DIAGNOSTICS_IS_ON('DRHODR ',myThid) ) |
127 |
& doDiagsRho = doDiagsRho + 2 |
128 |
IF ( DIAGNOSTICS_IS_ON('MXLDEPTH',myThid) ) |
129 |
& doDiagsRho = doDiagsRho + 4 |
130 |
ENDIF |
131 |
#endif /* ALLOW_DIAGNOSTICS */ |
132 |
|
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 |
CADJ & kind = isbyte |
140 |
cph-adj-test) |
141 |
CADJ STORE atemp,aqh,precip = comlev1, key = ikey_dynamics, |
142 |
CADJ & kind = isbyte |
143 |
CADJ STORE swdown,lwdown = comlev1, key = ikey_dynamics, |
144 |
CADJ & kind = isbyte |
145 |
cph# ifdef EXF_READ_EVAP |
146 |
CADJ STORE evap = comlev1, key = ikey_dynamics, |
147 |
CADJ & kind = isbyte |
148 |
cph# endif |
149 |
CADJ STORE uvel,vvel = comlev1, key = ikey_dynamics, |
150 |
CADJ & kind = isbyte |
151 |
# ifdef SEAICE_ALLOW_DYNAMICS |
152 |
CADJ STORE uice = comlev1, key = ikey_dynamics, |
153 |
CADJ & kind = isbyte |
154 |
CADJ STORE vice = comlev1, key = ikey_dynamics, |
155 |
CADJ & kind = isbyte |
156 |
# ifdef SEAICE_ALLOW_EVP |
157 |
CADJ STORE seaice_sigma1 = comlev1, key = ikey_dynamics, |
158 |
CADJ & kind = isbyte |
159 |
CADJ STORE seaice_sigma2 = comlev1, key = ikey_dynamics, |
160 |
CADJ & kind = isbyte |
161 |
CADJ STORE seaice_sigma12 = comlev1, key = ikey_dynamics, |
162 |
CADJ & kind = isbyte |
163 |
# endif |
164 |
# endif |
165 |
# ifdef SEAICE_SALINITY |
166 |
CADJ STORE salt = comlev1, key = ikey_dynamics, |
167 |
CADJ & kind = isbyte |
168 |
# endif |
169 |
# ifdef ATMOSPHERIC_LOADING |
170 |
CADJ STORE pload = comlev1, key = ikey_dynamics, |
171 |
CADJ & kind = isbyte |
172 |
CADJ STORE siceload = comlev1, key = ikey_dynamics, |
173 |
CADJ & kind = isbyte |
174 |
# endif |
175 |
# ifdef NONLIN_FRSURF |
176 |
CADJ STORE recip_hfacc = comlev1, key = ikey_dynamics, |
177 |
CADJ & kind = isbyte |
178 |
# endif |
179 |
# endif |
180 |
# ifdef ALLOW_DEBUG |
181 |
IF ( debugLevel .GE. debLevB ) |
182 |
& CALL DEBUG_CALL('SEAICE_MODEL',myThid) |
183 |
# endif |
184 |
CALL TIMER_START('SEAICE_MODEL [DO_OCEANIC_PHYS]', myThid) |
185 |
CALL SEAICE_MODEL( myTime, myIter, myThid ) |
186 |
CALL TIMER_STOP ('SEAICE_MODEL [DO_OCEANIC_PHYS]', myThid) |
187 |
# ifdef ALLOW_COST |
188 |
CALL SEAICE_COST_SENSI ( myTime, myIter, myThid ) |
189 |
# endif |
190 |
ENDIF |
191 |
#endif /* ALLOW_SEAICE */ |
192 |
|
193 |
#ifdef ALLOW_AUTODIFF_TAMC |
194 |
CADJ STORE sst, sss = comlev1, key = ikey_dynamics, |
195 |
CADJ & kind = isbyte |
196 |
CADJ STORE qsw = comlev1, key = ikey_dynamics, |
197 |
CADJ & kind = isbyte |
198 |
# ifdef ALLOW_SEAICE |
199 |
CADJ STORE area = comlev1, key = ikey_dynamics, |
200 |
CADJ & kind = isbyte |
201 |
# endif |
202 |
#endif |
203 |
|
204 |
#if (defined ALLOW_THSICE) && !(defined ALLOW_ATM2D) |
205 |
IF ( useThSIce .AND. fluidIsWater ) THEN |
206 |
#ifdef ALLOW_DEBUG |
207 |
IF ( debugLevel .GE. debLevB ) |
208 |
& CALL DEBUG_CALL('THSICE_MAIN',myThid) |
209 |
#endif |
210 |
C-- Step forward Therm.Sea-Ice variables |
211 |
C and modify forcing terms including effects from ice |
212 |
CALL TIMER_START('THSICE_MAIN [DO_OCEANIC_PHYS]', myThid) |
213 |
CALL THSICE_MAIN( myTime, myIter, myThid ) |
214 |
CALL TIMER_STOP( 'THSICE_MAIN [DO_OCEANIC_PHYS]', myThid) |
215 |
ENDIF |
216 |
#endif /* ALLOW_THSICE */ |
217 |
|
218 |
#ifdef ALLOW_SHELFICE |
219 |
IF ( useShelfIce .AND. fluidIsWater ) THEN |
220 |
#ifdef ALLOW_DEBUG |
221 |
IF ( debugLevel .GE. debLevB ) |
222 |
& CALL DEBUG_CALL('SHELFICE_THERMODYNAMICS',myThid) |
223 |
#endif |
224 |
C compute temperature and (virtual) salt flux at the |
225 |
C shelf-ice ocean interface |
226 |
CALL TIMER_START('SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]', |
227 |
& myThid) |
228 |
CALL SHELFICE_THERMODYNAMICS( myTime, myIter, myThid ) |
229 |
CALL TIMER_STOP( 'SHELFICE_THERMODYNAMICS [DO_OCEANIC_PHYS]', |
230 |
& myThid) |
231 |
ENDIF |
232 |
#endif /* ALLOW_SHELFICE */ |
233 |
|
234 |
C-- Freeze water at the surface |
235 |
#ifdef ALLOW_AUTODIFF_TAMC |
236 |
CADJ STORE theta = comlev1, key = ikey_dynamics, |
237 |
CADJ & kind = isbyte |
238 |
#endif |
239 |
IF ( allowFreezing ) THEN |
240 |
CALL FREEZE_SURFACE( myTime, myIter, myThid ) |
241 |
ENDIF |
242 |
|
243 |
#ifdef ALLOW_OCN_COMPON_INTERF |
244 |
C-- Apply imported data (from coupled interface) to forcing fields |
245 |
C jmc: do not know precisely where to put this call (bf or af thSIce ?) |
246 |
IF ( useCoupler ) THEN |
247 |
CALL OCN_APPLY_IMPORT( .TRUE., myTime, myIter, myThid ) |
248 |
ENDIF |
249 |
#endif /* ALLOW_OCN_COMPON_INTERF */ |
250 |
|
251 |
#ifdef ALLOW_BALANCE_FLUXES |
252 |
C balance fluxes |
253 |
IF ( balanceEmPmR ) |
254 |
& CALL REMOVE_MEAN_RS( 1, EmPmR, maskH, maskH, rA, drF, |
255 |
& 'EmPmR', myTime, myThid ) |
256 |
IF ( balanceQnet ) |
257 |
& CALL REMOVE_MEAN_RS( 1, Qnet, maskH, maskH, rA, drF, |
258 |
& 'Qnet ', myTime, myThid ) |
259 |
#endif /* ALLOW_BALANCE_FLUXES */ |
260 |
|
261 |
#ifdef ALLOW_AUTODIFF_TAMC |
262 |
C-- HPF directive to help TAMC |
263 |
CHPF$ INDEPENDENT |
264 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
265 |
DO bj=myByLo(myThid),myByHi(myThid) |
266 |
#ifdef ALLOW_AUTODIFF_TAMC |
267 |
C-- HPF directive to help TAMC |
268 |
CHPF$ INDEPENDENT |
269 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
270 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
271 |
|
272 |
#ifdef ALLOW_AUTODIFF_TAMC |
273 |
act1 = bi - myBxLo(myThid) |
274 |
max1 = myBxHi(myThid) - myBxLo(myThid) + 1 |
275 |
act2 = bj - myByLo(myThid) |
276 |
max2 = myByHi(myThid) - myByLo(myThid) + 1 |
277 |
act3 = myThid - 1 |
278 |
max3 = nTx*nTy |
279 |
act4 = ikey_dynamics - 1 |
280 |
itdkey = (act1 + 1) + act2*max1 |
281 |
& + act3*max1*max2 |
282 |
& + act4*max1*max2*max3 |
283 |
#else /* ALLOW_AUTODIFF_TAMC */ |
284 |
C if fluid is not water, by-pass find_rho, gmredi, surfaceForcing |
285 |
C and all vertical mixing schemes, but keep OBCS_CALC |
286 |
IF ( fluidIsWater ) THEN |
287 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
288 |
|
289 |
C-- Set up work arrays with valid (i.e. not NaN) values |
290 |
C These inital values do not alter the numerical results. They |
291 |
C just ensure that all memory references are to valid floating |
292 |
C point numbers. This prevents spurious hardware signals due to |
293 |
C uninitialised but inert locations. |
294 |
|
295 |
#ifdef ALLOW_AUTODIFF_TAMC |
296 |
DO j=1-OLy,sNy+OLy |
297 |
DO i=1-OLx,sNx+OLx |
298 |
rhoKm1 (i,j) = 0. _d 0 |
299 |
rhoKp1 (i,j) = 0. _d 0 |
300 |
ENDDO |
301 |
ENDDO |
302 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
303 |
|
304 |
DO k=1,Nr |
305 |
DO j=1-OLy,sNy+OLy |
306 |
DO i=1-OLx,sNx+OLx |
307 |
C This is currently used by GMRedi, IVDC, MXL-depth and Diagnostics |
308 |
sigmaX(i,j,k) = 0. _d 0 |
309 |
sigmaY(i,j,k) = 0. _d 0 |
310 |
sigmaR(i,j,k) = 0. _d 0 |
311 |
#ifdef ALLOW_AUTODIFF_TAMC |
312 |
cph all the following init. are necessary for TAF |
313 |
cph although some of these are re-initialised later. |
314 |
c rhoInSitu(i,j,k,bi,bj) = 0. |
315 |
IVDConvCount(i,j,k,bi,bj) = 0. |
316 |
# ifdef ALLOW_GMREDI |
317 |
Kwx(i,j,k,bi,bj) = 0. _d 0 |
318 |
Kwy(i,j,k,bi,bj) = 0. _d 0 |
319 |
Kwz(i,j,k,bi,bj) = 0. _d 0 |
320 |
# ifdef GM_NON_UNITY_DIAGONAL |
321 |
Kux(i,j,k,bi,bj) = 0. _d 0 |
322 |
Kvy(i,j,k,bi,bj) = 0. _d 0 |
323 |
# endif |
324 |
# ifdef GM_EXTRA_DIAGONAL |
325 |
Kuz(i,j,k,bi,bj) = 0. _d 0 |
326 |
Kvz(i,j,k,bi,bj) = 0. _d 0 |
327 |
# endif |
328 |
# ifdef GM_BOLUS_ADVEC |
329 |
GM_PsiX(i,j,k,bi,bj) = 0. _d 0 |
330 |
GM_PsiY(i,j,k,bi,bj) = 0. _d 0 |
331 |
# endif |
332 |
# ifdef GM_VISBECK_VARIABLE_K |
333 |
VisbeckK(i,j,bi,bj) = 0. _d 0 |
334 |
# endif |
335 |
# endif /* ALLOW_GMREDI */ |
336 |
# ifdef ALLOW_KPP |
337 |
KPPdiffKzS(i,j,k,bi,bj) = 0. _d 0 |
338 |
KPPdiffKzT(i,j,k,bi,bj) = 0. _d 0 |
339 |
# endif /* ALLOW_KPP */ |
340 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
341 |
ENDDO |
342 |
ENDDO |
343 |
ENDDO |
344 |
DO j=1-OLy,sNy+OLy |
345 |
DO i=1-OLx,sNx+OLx |
346 |
#ifdef ALLOW_AUTODIFF_TAMC |
347 |
# ifdef ALLOW_SALT_PLUME |
348 |
saltPlumeDepth(i,j,bi,bj) = 0. _d 0 |
349 |
saltPlumeFlux(i,j,bi,bj) = 0. _d 0 |
350 |
# endif |
351 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
352 |
ENDDO |
353 |
ENDDO |
354 |
|
355 |
iMin = 1-OLx |
356 |
iMax = sNx+OLx |
357 |
jMin = 1-OLy |
358 |
jMax = sNy+OLy |
359 |
|
360 |
#ifdef ALLOW_AUTODIFF_TAMC |
361 |
CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, |
362 |
CADJ & kind = isbyte |
363 |
CADJ STORE salt (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, |
364 |
CADJ & kind = isbyte |
365 |
CADJ STORE totphihyd(:,:,:,bi,bj) |
366 |
CADJ & = comlev1_bibj, key=itdkey, |
367 |
CADJ & kind = isbyte |
368 |
# ifdef ALLOW_KPP |
369 |
CADJ STORE uvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, |
370 |
CADJ & kind = isbyte |
371 |
CADJ STORE vvel (:,:,:,bi,bj) = comlev1_bibj, key=itdkey, |
372 |
CADJ & kind = isbyte |
373 |
# endif |
374 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
375 |
|
376 |
#ifdef ALLOW_DEBUG |
377 |
IF ( debugLevel .GE. debLevB ) |
378 |
& CALL DEBUG_MSG('ENTERING UPWARD K LOOP',myThid) |
379 |
#endif |
380 |
|
381 |
C-- Start of diagnostic loop |
382 |
DO k=Nr,1,-1 |
383 |
|
384 |
#ifdef ALLOW_AUTODIFF_TAMC |
385 |
C? Patrick, is this formula correct now that we change the loop range? |
386 |
C? Do we still need this? |
387 |
cph kkey formula corrected. |
388 |
cph Needed for rhoK, rhoKm1, in the case useGMREDI. |
389 |
kkey = (itdkey-1)*Nr + k |
390 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
391 |
|
392 |
C-- Always compute density (stored in common block) here; even when it is not |
393 |
C needed here, will be used anyway in calc_phi_hyd (data flow easier this way) |
394 |
#ifdef ALLOW_DEBUG |
395 |
IF ( debugLevel .GE. debLevB ) |
396 |
& CALL DEBUG_CALL('FIND_RHO_2D',myThid) |
397 |
#endif |
398 |
#ifdef ALLOW_AUTODIFF_TAMC |
399 |
IF ( fluidIsWater ) THEN |
400 |
CADJ STORE theta(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, |
401 |
CADJ & kind = isbyte |
402 |
CADJ STORE salt(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, |
403 |
CADJ & kind = isbyte |
404 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
405 |
#ifdef ALLOW_DOWN_SLOPE |
406 |
IF ( useDOWN_SLOPE ) THEN |
407 |
CALL DWNSLP_CALC_RHO( |
408 |
I theta, salt, |
409 |
O rhoInSitu(1-OLx,1-OLy,k,bi,bj), |
410 |
I k, bi, bj, myTime, myIter, myThid ) |
411 |
ELSE |
412 |
#endif /* ALLOW_DOWN_SLOPE */ |
413 |
CALL FIND_RHO_2D( |
414 |
I iMin, iMax, jMin, jMax, k, |
415 |
I theta(1-OLx,1-OLy,k,bi,bj), |
416 |
I salt (1-OLx,1-OLy,k,bi,bj), |
417 |
O rhoInSitu(1-OLx,1-OLy,k,bi,bj), |
418 |
I k, bi, bj, myThid ) |
419 |
#ifdef ALLOW_DOWN_SLOPE |
420 |
ENDIF |
421 |
#endif /* ALLOW_DOWN_SLOPE */ |
422 |
#ifdef ALLOW_AUTODIFF_TAMC |
423 |
ELSE |
424 |
C- fluid is not water: |
425 |
DO j=1-OLy,sNy+OLy |
426 |
DO i=1-OLx,sNx+OLx |
427 |
rhoInSitu(i,j,k,bi,bj) = 0. |
428 |
ENDDO |
429 |
ENDDO |
430 |
ENDIF |
431 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
432 |
|
433 |
C-- Calculate gradients of potential density for isoneutral |
434 |
C slope terms (e.g. GM/Redi tensor or IVDC diffusivity) |
435 |
IF ( useGMRedi .OR. (k.GT.1 .AND. ivdc_kappa.NE.0.) |
436 |
& .OR. useSALT_PLUME .OR. doDiagsRho.GE.1 ) THEN |
437 |
IF (k.GT.1) THEN |
438 |
#ifdef ALLOW_AUTODIFF_TAMC |
439 |
CADJ STORE theta(:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, |
440 |
CADJ & kind = isbyte |
441 |
CADJ STORE salt (:,:,k-1,bi,bj) = comlev1_bibj_k, key=kkey, |
442 |
CADJ & kind = isbyte |
443 |
CADJ STORE rhokm1 (bi,bj) = comlev1_bibj_k, key=kkey, |
444 |
CADJ & kind = isbyte |
445 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
446 |
CALL FIND_RHO_2D( |
447 |
I iMin, iMax, jMin, jMax, k, |
448 |
I theta(1-OLx,1-OLy,k-1,bi,bj), |
449 |
I salt (1-OLx,1-OLy,k-1,bi,bj), |
450 |
O rhoKm1, |
451 |
I k-1, bi, bj, myThid ) |
452 |
ENDIF |
453 |
#ifdef ALLOW_DEBUG |
454 |
IF ( debugLevel .GE. debLevB ) |
455 |
& CALL DEBUG_CALL('GRAD_SIGMA',myThid) |
456 |
#endif |
457 |
cph Avoid variable aliasing for adjoint !!! |
458 |
DO j=jMin,jMax |
459 |
DO i=iMin,iMax |
460 |
rhoKp1(i,j) = rhoInSitu(i,j,k,bi,bj) |
461 |
ENDDO |
462 |
ENDDO |
463 |
CALL GRAD_SIGMA( |
464 |
I bi, bj, iMin, iMax, jMin, jMax, k, |
465 |
I rhoInSitu(1-OLx,1-OLy,k,bi,bj), rhoKm1, rhoKp1, |
466 |
O sigmaX, sigmaY, sigmaR, |
467 |
I myThid ) |
468 |
#ifdef ALLOW_AUTODIFF_TAMC |
469 |
#ifdef GMREDI_WITH_STABLE_ADJOINT |
470 |
cgf zero out adjoint fields to stabilize pkg/gmredi adjoint |
471 |
cgf -> cuts adjoint dependency from slope to state |
472 |
CALL ZERO_ADJ_LOC( Nr, sigmaX, myThid) |
473 |
CALL ZERO_ADJ_LOC( Nr, sigmaY, myThid) |
474 |
CALL ZERO_ADJ_LOC( Nr, sigmaR, myThid) |
475 |
#endif |
476 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
477 |
ENDIF |
478 |
|
479 |
C-- Implicit Vertical Diffusion for Convection |
480 |
c ==> should use sigmaR !!! |
481 |
IF (k.GT.1 .AND. ivdc_kappa.NE.0.) THEN |
482 |
#ifdef ALLOW_DEBUG |
483 |
IF ( debugLevel .GE. debLevB ) |
484 |
& CALL DEBUG_CALL('CALC_IVDC',myThid) |
485 |
#endif |
486 |
CALL CALC_IVDC( |
487 |
I bi, bj, iMin, iMax, jMin, jMax, k, |
488 |
I rhoKm1, rhoInSitu(1-OLx,1-OLy,k,bi,bj), |
489 |
I myTime, myIter, myThid) |
490 |
ENDIF |
491 |
|
492 |
#ifdef ALLOW_DIAGNOSTICS |
493 |
IF ( MOD(doDiagsRho,2).EQ.1 ) THEN |
494 |
CALL DIAGS_RHO_L( k, bi, bj, |
495 |
I rhoInSitu(1-OLx,1-OLy,k,bi,bj), |
496 |
I rhoKm1, wVel, |
497 |
I myTime, myIter, myThid ) |
498 |
ENDIF |
499 |
#endif |
500 |
|
501 |
C-- end of diagnostic k loop (Nr:1) |
502 |
ENDDO |
503 |
|
504 |
#ifdef ALLOW_AUTODIFF_TAMC |
505 |
CADJ STORE IVDConvCount(:,:,:,bi,bj) |
506 |
CADJ & = comlev1_bibj, key=itdkey, |
507 |
CADJ & kind = isbyte |
508 |
#endif |
509 |
|
510 |
C-- Diagnose Mixed Layer Depth: |
511 |
IF ( useGMRedi .OR. doDiagsRho.GE.4 ) THEN |
512 |
CALL CALC_OCE_MXLAYER( |
513 |
I rhoInSitu(1-OLx,1-OLy,1,bi,bj), sigmaR, |
514 |
I bi, bj, myTime, myIter, myThid ) |
515 |
ENDIF |
516 |
|
517 |
#ifdef ALLOW_SALT_PLUME |
518 |
IF ( useSALT_PLUME ) THEN |
519 |
CALL SALT_PLUME_CALC_DEPTH( |
520 |
I rhoInSitu(1-OLx,1-OLy,1,bi,bj), sigmaR, |
521 |
I bi, bj, myTime, myIter, myThid ) |
522 |
ENDIF |
523 |
#endif /* ALLOW_SALT_PLUME */ |
524 |
|
525 |
#ifdef ALLOW_DIAGNOSTICS |
526 |
IF ( MOD(doDiagsRho,4).GE.2 ) THEN |
527 |
CALL DIAGNOSTICS_FILL (sigmaR, 'DRHODR ', 0, Nr, |
528 |
& 2, bi, bj, myThid) |
529 |
ENDIF |
530 |
#endif /* ALLOW_DIAGNOSTICS */ |
531 |
|
532 |
C-- Determines forcing terms based on external fields |
533 |
C relaxation terms, etc. |
534 |
#ifdef ALLOW_DEBUG |
535 |
IF ( debugLevel .GE. debLevB ) |
536 |
& CALL DEBUG_CALL('EXTERNAL_FORCING_SURF',myThid) |
537 |
#endif |
538 |
#ifdef ALLOW_AUTODIFF_TAMC |
539 |
CADJ STORE EmPmR(:,:,bi,bj) |
540 |
CADJ & = comlev1_bibj, key=itdkey, |
541 |
CADJ & kind = isbyte |
542 |
# ifdef EXACT_CONSERV |
543 |
CADJ STORE PmEpR(:,:,bi,bj) |
544 |
CADJ & = comlev1_bibj, key=itdkey, |
545 |
CADJ & kind = isbyte |
546 |
# endif |
547 |
# ifdef NONLIN_FRSURF |
548 |
CADJ STORE hFac_surfC(:,:,bi,bj) |
549 |
CADJ & = comlev1_bibj, key=itdkey, |
550 |
CADJ & kind = isbyte |
551 |
CADJ STORE recip_hFacC(:,:,:,bi,bj) |
552 |
CADJ & = comlev1_bibj, key=itdkey, |
553 |
CADJ & kind = isbyte |
554 |
# endif |
555 |
#endif |
556 |
CALL EXTERNAL_FORCING_SURF( |
557 |
I bi, bj, iMin, iMax, jMin, jMax, |
558 |
I myTime, myIter, myThid ) |
559 |
#ifdef ALLOW_AUTODIFF_TAMC |
560 |
# ifdef EXACT_CONSERV |
561 |
cph-test |
562 |
cphCADJ STORE PmEpR(:,:,bi,bj) |
563 |
cphCADJ & = comlev1_bibj, key=itdkey, |
564 |
cphCADJ & kind = isbyte |
565 |
# endif |
566 |
#endif |
567 |
|
568 |
#ifdef ALLOW_AUTODIFF_TAMC |
569 |
cph needed for KPP |
570 |
CADJ STORE surfaceForcingU(:,:,bi,bj) |
571 |
CADJ & = comlev1_bibj, key=itdkey, |
572 |
CADJ & kind = isbyte |
573 |
CADJ STORE surfaceForcingV(:,:,bi,bj) |
574 |
CADJ & = comlev1_bibj, key=itdkey, |
575 |
CADJ & kind = isbyte |
576 |
CADJ STORE surfaceForcingS(:,:,bi,bj) |
577 |
CADJ & = comlev1_bibj, key=itdkey, |
578 |
CADJ & kind = isbyte |
579 |
CADJ STORE surfaceForcingT(:,:,bi,bj) |
580 |
CADJ & = comlev1_bibj, key=itdkey, |
581 |
CADJ & kind = isbyte |
582 |
CADJ STORE surfaceForcingTice(:,:,bi,bj) |
583 |
CADJ & = comlev1_bibj, key=itdkey, |
584 |
CADJ & kind = isbyte |
585 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
586 |
|
587 |
#ifdef ALLOW_KPP |
588 |
C-- Compute KPP mixing coefficients |
589 |
IF (useKPP) THEN |
590 |
#ifdef ALLOW_DEBUG |
591 |
IF ( debugLevel .GE. debLevB ) |
592 |
& CALL DEBUG_CALL('KPP_CALC',myThid) |
593 |
#endif |
594 |
CALL TIMER_START('KPP_CALC [DO_OCEANIC_PHYS]', myThid) |
595 |
CALL KPP_CALC( |
596 |
I bi, bj, myTime, myIter, myThid ) |
597 |
CALL TIMER_STOP ('KPP_CALC [DO_OCEANIC_PHYS]', myThid) |
598 |
#ifdef ALLOW_AUTODIFF_TAMC |
599 |
ELSE |
600 |
CALL KPP_CALC_DUMMY( |
601 |
I bi, bj, myTime, myIter, myThid ) |
602 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
603 |
ENDIF |
604 |
|
605 |
#endif /* ALLOW_KPP */ |
606 |
|
607 |
#ifdef ALLOW_PP81 |
608 |
C-- Compute PP81 mixing coefficients |
609 |
IF (usePP81) THEN |
610 |
#ifdef ALLOW_DEBUG |
611 |
IF ( debugLevel .GE. debLevB ) |
612 |
& CALL DEBUG_CALL('PP81_CALC',myThid) |
613 |
#endif |
614 |
CALL PP81_CALC( |
615 |
I bi, bj, myTime, myThid ) |
616 |
ENDIF |
617 |
#endif /* ALLOW_PP81 */ |
618 |
|
619 |
#ifdef ALLOW_MY82 |
620 |
C-- Compute MY82 mixing coefficients |
621 |
IF (useMY82) THEN |
622 |
#ifdef ALLOW_DEBUG |
623 |
IF ( debugLevel .GE. debLevB ) |
624 |
& CALL DEBUG_CALL('MY82_CALC',myThid) |
625 |
#endif |
626 |
CALL MY82_CALC( |
627 |
I bi, bj, myTime, myThid ) |
628 |
ENDIF |
629 |
#endif /* ALLOW_MY82 */ |
630 |
|
631 |
#ifdef ALLOW_GGL90 |
632 |
C-- Compute GGL90 mixing coefficients |
633 |
IF (useGGL90) THEN |
634 |
#ifdef ALLOW_DEBUG |
635 |
IF ( debugLevel .GE. debLevB ) |
636 |
& CALL DEBUG_CALL('GGL90_CALC',myThid) |
637 |
#endif |
638 |
CALL TIMER_START('GGL90_CALC [DO_OCEANIC_PHYS]', myThid) |
639 |
CALL GGL90_CALC( |
640 |
I bi, bj, myTime, myThid ) |
641 |
CALL TIMER_STOP ('GGL90_CALC [DO_OCEANIC_PHYS]', myThid) |
642 |
ENDIF |
643 |
#endif /* ALLOW_GGL90 */ |
644 |
|
645 |
#ifdef ALLOW_TIMEAVE |
646 |
IF ( taveFreq.GT. 0. _d 0 ) THEN |
647 |
CALL TIMEAVE_SURF_FLUX( bi, bj, myTime, myIter, myThid) |
648 |
ENDIF |
649 |
IF (taveFreq.GT.0. .AND. ivdc_kappa.NE.0.) THEN |
650 |
CALL TIMEAVE_CUMULATE(ConvectCountTave, IVDConvCount, |
651 |
I Nr, deltaTclock, bi, bj, myThid) |
652 |
ENDIF |
653 |
#endif /* ALLOW_TIMEAVE */ |
654 |
|
655 |
#ifdef ALLOW_GMREDI |
656 |
#ifdef ALLOW_AUTODIFF_TAMC |
657 |
# ifndef GM_EXCLUDE_CLIPPING |
658 |
cph storing here is needed only for one GMREDI_OPTIONS: |
659 |
cph define GM_BOLUS_ADVEC |
660 |
cph keep it although TAF says you dont need to. |
661 |
cph but I've avoided the #ifdef for now, in case more things change |
662 |
CADJ STORE sigmaX(:,:,:) = comlev1_bibj, key=itdkey, |
663 |
CADJ & kind = isbyte |
664 |
CADJ STORE sigmaY(:,:,:) = comlev1_bibj, key=itdkey, |
665 |
CADJ & kind = isbyte |
666 |
CADJ STORE sigmaR(:,:,:) = comlev1_bibj, key=itdkey, |
667 |
CADJ & kind = isbyte |
668 |
# endif |
669 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
670 |
|
671 |
C-- Calculate iso-neutral slopes for the GM/Redi parameterisation |
672 |
IF (useGMRedi) THEN |
673 |
#ifdef ALLOW_DEBUG |
674 |
IF ( debugLevel .GE. debLevB ) |
675 |
& CALL DEBUG_CALL('GMREDI_CALC_TENSOR',myThid) |
676 |
#endif |
677 |
CALL GMREDI_CALC_TENSOR( |
678 |
I iMin, iMax, jMin, jMax, |
679 |
I sigmaX, sigmaY, sigmaR, |
680 |
I bi, bj, myTime, myIter, myThid ) |
681 |
#ifdef ALLOW_AUTODIFF_TAMC |
682 |
ELSE |
683 |
CALL GMREDI_CALC_TENSOR_DUMMY( |
684 |
I iMin, iMax, jMin, jMax, |
685 |
I sigmaX, sigmaY, sigmaR, |
686 |
I bi, bj, myTime, myIter, myThid ) |
687 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
688 |
ENDIF |
689 |
#endif /* ALLOW_GMREDI */ |
690 |
|
691 |
#ifdef ALLOW_DOWN_SLOPE |
692 |
IF ( useDOWN_SLOPE ) THEN |
693 |
C-- Calculate Downsloping Flow for Down_Slope parameterization |
694 |
IF ( usingPCoords ) THEN |
695 |
CALL DWNSLP_CALC_FLOW( |
696 |
I bi, bj, kSurfC, rhoInSitu, |
697 |
I myTime, myIter, myThid ) |
698 |
ELSE |
699 |
CALL DWNSLP_CALC_FLOW( |
700 |
I bi, bj, kLowC, rhoInSitu, |
701 |
I myTime, myIter, myThid ) |
702 |
ENDIF |
703 |
ENDIF |
704 |
#endif /* ALLOW_DOWN_SLOPE */ |
705 |
|
706 |
#ifndef ALLOW_AUTODIFF_TAMC |
707 |
C--- if fluid Is Water: end |
708 |
ENDIF |
709 |
#endif |
710 |
|
711 |
#ifdef ALLOW_OBCS |
712 |
C-- Calculate future values on open boundaries |
713 |
IF (useOBCS) THEN |
714 |
#ifdef ALLOW_DEBUG |
715 |
IF ( debugLevel .GE. debLevB ) |
716 |
& CALL DEBUG_CALL('OBCS_CALC',myThid) |
717 |
#endif |
718 |
CALL OBCS_CALC( bi, bj, myTime+deltaTclock, myIter+1, |
719 |
I uVel, vVel, wVel, theta, salt, |
720 |
I myThid ) |
721 |
ENDIF |
722 |
#endif /* ALLOW_OBCS */ |
723 |
|
724 |
C-- end bi,bj loops. |
725 |
ENDDO |
726 |
ENDDO |
727 |
|
728 |
#ifdef ALLOW_KPP |
729 |
IF (useKPP) THEN |
730 |
CALL KPP_DO_EXCH( myThid ) |
731 |
ENDIF |
732 |
#endif /* ALLOW_KPP */ |
733 |
|
734 |
#ifdef ALLOW_DIAGNOSTICS |
735 |
IF ( fluidIsWater .AND. useDiagnostics ) THEN |
736 |
CALL DIAGS_RHO_G( |
737 |
I rhoInSitu, uVel, vVel, |
738 |
I myTime, myIter, myThid ) |
739 |
CALL DIAGS_OCEANIC_SURF_FLUX( myTime, myIter, myThid ) |
740 |
ENDIF |
741 |
IF ( ivdc_kappa.NE.0 .AND. useDiagnostics ) THEN |
742 |
CALL DIAGNOSTICS_FILL( IVDConvCount, 'CONVADJ ', |
743 |
& 0, Nr, 0, 1, 1, myThid ) |
744 |
ENDIF |
745 |
#endif |
746 |
|
747 |
#ifdef ALLOW_DEBUG |
748 |
IF ( debugLevel .GE. debLevB ) |
749 |
& CALL DEBUG_LEAVE('DO_OCEANIC_PHYS',myThid) |
750 |
#endif |
751 |
|
752 |
RETURN |
753 |
END |