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