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