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