1 |
C $Header: /u/gcmpack/MITgcm/pkg/layers/layers_fluxcalc.F,v 1.16 2015/06/08 16:26:12 rpa Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "LAYERS_OPTIONS.h" |
5 |
#ifdef ALLOW_GMREDI |
6 |
#include "GMREDI_OPTIONS.h" |
7 |
#endif |
8 |
|
9 |
C-- File layers_fluxcalc.F: |
10 |
C-- Contents |
11 |
C-- o LAYERS_FLUXCALC |
12 |
C-- o LAYERS_DIAPYCNAL |
13 |
C-- o LAYERS_LOCATE |
14 |
|
15 |
CBOP 0 |
16 |
C !ROUTINE: LAYERS_FLUXCALC |
17 |
C !INTERFACE: |
18 |
SUBROUTINE LAYERS_FLUXCALC( |
19 |
I uVel,vVel,tracer,iLa, |
20 |
O UH,VH,Hw,Hs,PIw,PIs,Uw,Vs, |
21 |
I myThid ) |
22 |
|
23 |
C !DESCRIPTION: \bv |
24 |
C *==========================================================* |
25 |
C | SUBROUTINE LAYERS_FLUXCALC |
26 |
C | Calculate the transport in isotracer layers, for a chosen |
27 |
C | tracer. This is the meat of the LAYERS package. |
28 |
C *==========================================================* |
29 |
C \ev |
30 |
|
31 |
C !USES: |
32 |
IMPLICIT NONE |
33 |
C == Global variables === |
34 |
#include "SIZE.h" |
35 |
#include "EEPARAMS.h" |
36 |
#include "PARAMS.h" |
37 |
#include "GRID.h" |
38 |
#include "LAYERS_SIZE.h" |
39 |
#include "LAYERS.h" |
40 |
#ifdef ALLOW_GMREDI |
41 |
# include "GMREDI.h" |
42 |
#endif |
43 |
|
44 |
C !INPUT PARAMETERS: |
45 |
C myThid :: my Thread Id number |
46 |
C uVel :: zonal velocity (m/s, i=1 held at western face) |
47 |
C vVel :: meridional velocity (m/s, j=1 held at southern face) |
48 |
C tracer :: potential temperature, salt or potential density prho |
49 |
C UH :: U integrated over layer (m^2/s) |
50 |
C VH :: V integrated over layer (m^2/s) |
51 |
C Hw :: Layer thickness at the U point (m) |
52 |
C Hs :: Layer thickness at the V point (m) |
53 |
C PIw :: 1 if layer exists, 0 otherwise (at U point) |
54 |
C PIs :: 1 if layer exists, 0 otherwise (at V point) |
55 |
C Uw :: average U over layer (m/s) |
56 |
C Vs :: average V over layer (m/s) |
57 |
C iLa :: layer coordinate index |
58 |
INTEGER myThid |
59 |
_RL uVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr, nSx,nSy) |
60 |
_RL vVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr, nSx,nSy) |
61 |
_RL tracer (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr, nSx,nSy) |
62 |
#ifdef LAYERS_UFLUX |
63 |
_RL UH (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy) |
64 |
# ifdef LAYERS_THICKNESS |
65 |
_RL Hw (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy) |
66 |
_RL PIw (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy) |
67 |
_RL Uw (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy) |
68 |
# else |
69 |
_RL Hw(1), PIw(1), Uw(1) |
70 |
# endif |
71 |
#else |
72 |
_RL UH(1), Hw(1), PIw(1), Uw(1) |
73 |
#endif |
74 |
#ifdef LAYERS_VFLUX |
75 |
_RL VH (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy) |
76 |
# ifdef LAYERS_THICKNESS |
77 |
_RL Hs (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy) |
78 |
_RL PIs (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy) |
79 |
_RL Vs (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy) |
80 |
# else |
81 |
_RL Hs(1), PIs(1), Vs(1) |
82 |
# endif |
83 |
#else |
84 |
_RL VH(1), Hs(1), PIs(1), Vs(1) |
85 |
#endif |
86 |
INTEGER iLa |
87 |
CEOP |
88 |
|
89 |
#ifdef ALLOW_LAYERS |
90 |
|
91 |
C !LOCAL VARIABLES: |
92 |
C bi, bj :: tile indices |
93 |
C i,j :: horizontal indices |
94 |
C k :: vertical index for model grid |
95 |
C kci :: index from CellIndex |
96 |
C kg :: index for looping though layers_bounds |
97 |
C kk :: vertical index for ZZ (fine) grid |
98 |
C kgu,kgv :: vertical index for isopycnal grid |
99 |
C kloc :: local copy of kgu/v to reduce accesses to index arrays |
100 |
C mSteps :: maximum number of steps for bisection method |
101 |
C prho :: pot. density (less 1000) referenced to layers_krho pressure |
102 |
C TatU :: temperature at U point |
103 |
C TatV :: temperature at V point |
104 |
C dzfac :: temporary sublayer thickness |
105 |
C Tloc,Tp1 :: horizontally interpolated tracer values |
106 |
|
107 |
INTEGER bi, bj |
108 |
INTEGER i,j,k,kk,kg,kci,kp1,kloc |
109 |
INTEGER mSteps |
110 |
INTEGER kgu(sNx+1,sNy+1), kgv(sNx+1,sNy+1) |
111 |
_RL TatU(sNx+1,sNy+1), TatV(sNx+1,sNy+1) |
112 |
_RL dzfac |
113 |
#ifdef ALLOW_GMREDI |
114 |
INTEGER kcip1 |
115 |
_RL delPsi, maskp1 |
116 |
#endif |
117 |
LOGICAL errorFlag |
118 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
119 |
|
120 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
121 |
|
122 |
C compute maximum number of steps for bisection method (approx. |
123 |
C log2(Nlayers)) as log2(Nlayers) + 1 for safety |
124 |
mSteps = int(log10(dble(Nlayers))/log10(2. _d 0))+1 |
125 |
|
126 |
C --- The tile loops |
127 |
DO bj=myByLo(myThid),myByHi(myThid) |
128 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
129 |
|
130 |
C Initialize the search indices |
131 |
DO j = 1,sNy+1 |
132 |
DO i = 1,sNx+1 |
133 |
C The temperature index (layer_G) goes from cold to warm. |
134 |
C The water column goes from warm (k=1) to cold (k=Nr). |
135 |
C So initialize the search with the warmest value. |
136 |
kgu(i,j) = Nlayers |
137 |
kgv(i,j) = Nlayers |
138 |
ENDDO |
139 |
ENDDO |
140 |
|
141 |
C Reset the arrays |
142 |
DO kg=1,Nlayers |
143 |
DO j = 1-OLy,sNy+OLy |
144 |
DO i = 1-OLx,sNx+OLx |
145 |
#ifdef LAYERS_UFLUX |
146 |
UH (i,j,kg,bi,bj) = 0. _d 0 |
147 |
#ifdef LAYERS_THICKNESS |
148 |
Hw(i,j,kg,bi,bj) = 0. _d 0 |
149 |
PIw(i,j,kg,bi,bj) = 0. _d 0 |
150 |
Uw(i,j,kg,bi,bj) = 0. _d 0 |
151 |
#endif /* LAYERS_THICKNESS */ |
152 |
#endif /* UH */ |
153 |
#ifdef LAYERS_VFLUX |
154 |
VH (i,j,kg,bi,bj) = 0. _d 0 |
155 |
#ifdef LAYERS_THICKNESS |
156 |
Hs(i,j,kg,bi,bj) = 0. _d 0 |
157 |
PIs(i,j,kg,bi,bj) = 0. _d 0 |
158 |
Vs(i,j,kg,bi,bj) = 0. _d 0 |
159 |
#endif /* LAYERS_THICKNESS */ |
160 |
#endif /* VH */ |
161 |
ENDDO |
162 |
ENDDO |
163 |
ENDDO |
164 |
|
165 |
DO kk=1,NZZ |
166 |
k = MapIndex(kk) |
167 |
kci = CellIndex(kk) |
168 |
#ifdef ALLOW_GMREDI |
169 |
kcip1 = MIN(kci+1,Nr) |
170 |
maskp1 = 1. |
171 |
IF (kci.GE.Nr) maskp1 = 0. |
172 |
#endif /* ALLOW_GMREDI */ |
173 |
#ifdef LAYERS_UFLUX |
174 |
DO j = 1,sNy+1 |
175 |
DO i = 1,sNx+1 |
176 |
|
177 |
C ------ Find theta at the U point (west) on the fine Z grid |
178 |
kp1=k+1 |
179 |
IF (maskW(i,j,kp1,bi,bj).EQ.zeroRS) kp1=k |
180 |
TatU(i,j) = MapFact(kk) * |
181 |
& 0.5 _d 0 * (tracer(i-1,j,k,bi,bj)+tracer(i,j,k,bi,bj)) + |
182 |
& (1. _d 0 -MapFact(kk)) * |
183 |
& 0.5 _d 0 * (tracer(i-1,j,kp1,bi,bj)+tracer(i,j,kp1,bi,bj)) |
184 |
|
185 |
ENDDO |
186 |
ENDDO |
187 |
C ------ Now that we know T everywhere, determine the binning. |
188 |
C find the layer indices kgu |
189 |
CALL LAYERS_LOCATE( |
190 |
I layers_bounds(1,iLa),Nlayers,mSteps,sNx,sNy,TatU, |
191 |
O kgu, |
192 |
I myThid ) |
193 |
#ifndef TARGET_NEC_SX |
194 |
C check for failures |
195 |
IF ( debugLevel .GE. debLevC ) THEN |
196 |
errorFlag = .FALSE. |
197 |
DO j = 1,sNy+1 |
198 |
DO i = 1,sNx+1 |
199 |
IF ( kgu(i,j) .LE. 0 ) THEN |
200 |
WRITE(msgBuf,'(2A,I3,A,I3,A,1E14.6)') |
201 |
& 'S/R LAYERS_LOCATE: Could not find a bin in ', |
202 |
& 'layers_bounds for TatU(',i,',',j,',)=',TatU(i,j) |
203 |
CALL PRINT_ERROR( msgBuf, myThid ) |
204 |
errorFlag = .TRUE. |
205 |
ENDIF |
206 |
ENDDO |
207 |
ENDDO |
208 |
IF ( errorFlag ) STOP 'ABNORMAL END: S/R LAYERS_FLUXCALC' |
209 |
ENDIF |
210 |
#endif /* ndef TARGET_NEC_SX */ |
211 |
C |
212 |
DO j = 1,sNy+1 |
213 |
DO i = 1,sNx+1 |
214 |
|
215 |
kloc = kgu(i,j) |
216 |
dzfac = dZZf(kk) * hFacW(i,j,kci,bi,bj) |
217 |
|
218 |
C ------ Augment the bin values |
219 |
UH(i,j,kloc,bi,bj) = |
220 |
& UH(i,j,kloc,bi,bj) + |
221 |
& dzfac * uVel(i,j,kci,bi,bj) |
222 |
|
223 |
#ifdef ALLOW_GMREDI |
224 |
IF ( layers_bolus(iLa) ) THEN |
225 |
IF ( .NOT.GM_AdvForm ) THEN |
226 |
delPsi = 0.25 _d 0 *( |
227 |
& ( rA(i-1,j,bi,bj)*Kwx(i-1,j,kcip1,bi,bj) |
228 |
& +rA( i ,j,bi,bj)*Kwx( i ,j,kcip1,bi,bj) |
229 |
& ) * maskW(i,j,kcip1,bi,bj) * maskp1 |
230 |
& - ( rA(i-1,j,bi,bj)*Kwx(i-1,j, kci ,bi,bj) |
231 |
& +rA( i ,j,bi,bj)*Kwx( i ,j, kci ,bi,bj) |
232 |
& ) * maskW(i,j, kci ,bi,bj) |
233 |
& ) * recip_rAw(i,j,bi,bj) |
234 |
#ifdef GM_BOLUS_ADVEC |
235 |
ELSE |
236 |
delPsi = GM_PsiX(i,j,kcip1,bi,bj)*maskp1 |
237 |
& - GM_PsiX(i,j, kci, bi,bj) |
238 |
#endif |
239 |
ENDIF |
240 |
UH(i,j,kloc,bi,bj) = UH(i,j,kloc,bi,bj) |
241 |
& + delPsi*recip_drF(kci)*_recip_hFacW(i,j,kci,bi,bj) |
242 |
& * dzfac |
243 |
ENDIF |
244 |
#endif /* ALLOW_GMREDI */ |
245 |
|
246 |
#ifdef LAYERS_THICKNESS |
247 |
Hw(i,j,kloc,bi,bj) = Hw(i,j,kloc,bi,bj) + dzfac |
248 |
#endif /* LAYERS_THICKNESS */ |
249 |
|
250 |
ENDDO |
251 |
ENDDO |
252 |
#endif /* LAYERS_UFLUX */ |
253 |
|
254 |
#ifdef LAYERS_VFLUX |
255 |
DO j = 1,sNy+1 |
256 |
DO i = 1,sNx+1 |
257 |
C ------ Find theta at the V point (south) on the fine Z grid |
258 |
kp1=k+1 |
259 |
IF (maskS(i,j,kp1,bi,bj).EQ.zeroRS) kp1=k |
260 |
TatV(i,j) = MapFact(kk) * |
261 |
& 0.5 _d 0 * (tracer(i,j-1,k,bi,bj)+tracer(i,j,k,bi,bj)) + |
262 |
& (1. _d 0 -MapFact(kk)) * |
263 |
& 0.5 _d 0 * (tracer(i,j-1,kp1,bi,bj)+tracer(i,j,kp1,bi,bj)) |
264 |
|
265 |
ENDDO |
266 |
ENDDO |
267 |
C ------ Now that we know T everywhere, determine the binning. |
268 |
C find the layer indices kgv |
269 |
CALL LAYERS_LOCATE( |
270 |
I layers_bounds(1,iLa),Nlayers,mSteps,sNx,sNy,TatV, |
271 |
O kgv, |
272 |
I myThid ) |
273 |
#ifndef TARGET_NEC_SX |
274 |
IF ( debugLevel .GE. debLevC ) THEN |
275 |
C check for failures |
276 |
errorFlag = .FALSE. |
277 |
DO j = 1,sNy+1 |
278 |
DO i = 1,sNx+1 |
279 |
IF ( kgv(i,j) .LE. 0 ) THEN |
280 |
WRITE(msgBuf,'(2A,I3,A,I3,A,1E14.6)') |
281 |
& 'S/R LAYERS_LOCATE: Could not find a bin in ', |
282 |
& 'layers_bounds for TatV(',i,',',j,',)=',TatV(i,j) |
283 |
CALL PRINT_ERROR( msgBuf, myThid ) |
284 |
errorFlag = .TRUE. |
285 |
ENDIF |
286 |
ENDDO |
287 |
ENDDO |
288 |
IF ( errorFlag ) STOP 'ABNORMAL END: S/R LAYERS_FLUXCALC' |
289 |
ENDIF |
290 |
#endif /* ndef TARGET_NEC_SX */ |
291 |
C |
292 |
DO j = 1,sNy+1 |
293 |
DO i = 1,sNx+1 |
294 |
|
295 |
kloc = kgv(i,j) |
296 |
dzfac = dZZf(kk) * hFacS(i,j,kci,bi,bj) |
297 |
|
298 |
C ------ debugging stuff |
299 |
C IF (i.EQ.10 .AND. j.EQ.10) THEN |
300 |
C WRITE(msgBuf,'(A,I3,A,F10.2,A,F6.2)') |
301 |
C & ' kloc=', kloc, |
302 |
C & ', TatV=',TatV(i,j), |
303 |
C & ', dzfac=',dzfac |
304 |
C CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
305 |
C & SQUEEZE_RIGHT, myThid ) |
306 |
C ENDIF |
307 |
|
308 |
C ------ Augment the bin values |
309 |
VH(i,j,kloc,bi,bj) = |
310 |
& VH(i,j,kloc,bi,bj) + dzfac * vVel(i,j,kci,bi,bj) |
311 |
|
312 |
#ifdef ALLOW_GMREDI |
313 |
IF ( layers_bolus(iLa) ) THEN |
314 |
IF ( .NOT.GM_AdvForm ) THEN |
315 |
delPsi = 0.25 _d 0 *( |
316 |
& ( rA(i,j-1,bi,bj)*Kwy(i,j-1,kcip1,bi,bj) |
317 |
& +rA(i, j ,bi,bj)*Kwy(i, j ,kcip1,bi,bj) |
318 |
& ) * maskS(i,j,kcip1,bi,bj) * maskp1 |
319 |
& - ( rA(i,j-1,bi,bj)*Kwy(i,j-1, kci ,bi,bj) |
320 |
& +rA(i, j ,bi,bj)*Kwy(i, j , kci ,bi,bj) |
321 |
& ) * maskS(i,j, kci ,bi,bj) |
322 |
& ) * recip_rAs(i,j,bi,bj) |
323 |
#ifdef GM_BOLUS_ADVEC |
324 |
ELSE |
325 |
delPsi = GM_PsiY(i,j,kcip1,bi,bj)*maskp1 |
326 |
& - GM_PsiY(i,j, kci, bi,bj) |
327 |
#endif |
328 |
ENDIF |
329 |
VH(i,j,kloc,bi,bj) = VH(i,j,kloc,bi,bj) |
330 |
& + delPsi*recip_drF(kci)*_recip_hFacS(i,j,kci,bi,bj) |
331 |
& * dzfac |
332 |
ENDIF |
333 |
#endif /* ALLOW_GMREDI */ |
334 |
|
335 |
#ifdef LAYERS_THICKNESS |
336 |
Hs(i,j,kloc,bi,bj) = Hs(i,j,kloc,bi,bj) + dzfac |
337 |
#endif /* LAYERS_THICKNESS */ |
338 |
|
339 |
ENDDO |
340 |
ENDDO |
341 |
#endif /* LAYERS_VFLUX */ |
342 |
ENDDO |
343 |
|
344 |
C-- Now that we know the thicknesses, compute the heaviside function |
345 |
C-- (Needs another loop through Ng) |
346 |
#ifdef LAYERS_THICKNESS |
347 |
DO kg=1,Nlayers |
348 |
DO j = 1,sNy+1 |
349 |
DO i = 1,sNx+1 |
350 |
#ifdef LAYERS_UFLUX |
351 |
IF (Hw(i,j,kg,bi,bj) .GT. 0.) THEN |
352 |
PIw(i,j,kg,bi,bj) = 1. _d 0 |
353 |
Uw(i,j,kg,bi,bj) = |
354 |
& UH(i,j,kg,bi,bj) / Hw(i,j,kg,bi,bj) |
355 |
ENDIF |
356 |
#endif /* LAYERS_UFLUX */ |
357 |
#ifdef LAYERS_VFLUX |
358 |
IF (Hs(i,j,kg,bi,bj) .GT. 0.) THEN |
359 |
PIs(i,j,kg,bi,bj) = 1. _d 0 |
360 |
Vs(i,j,kg,bi,bj) = |
361 |
& VH(i,j,kg,bi,bj) / Hs(i,j,kg,bi,bj) |
362 |
ENDIF |
363 |
#endif /* LAYERS_VFLUX */ |
364 |
ENDDO |
365 |
ENDDO |
366 |
ENDDO |
367 |
#endif /* LAYERS_THICKNESS */ |
368 |
|
369 |
C --- End bi,bj loop |
370 |
ENDDO |
371 |
ENDDO |
372 |
|
373 |
RETURN |
374 |
END |
375 |
|
376 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
377 |
CBOP 0 |
378 |
C !ROUTINE: LAYERS_DIAPYCNAL |
379 |
C !INTERFACE: |
380 |
SUBROUTINE LAYERS_DIAPYCNAL( |
381 |
I tracer,iLa, |
382 |
O TtendSurf, TtendDiffh, TtendDiffr, |
383 |
O TtendAdvh, TtendAdvr, Ttendtot, |
384 |
O StendSurf, StendDiffh, StendDiffr, |
385 |
O StendAdvh, StendAdvr, Stendtot, |
386 |
O Hc, PIc, |
387 |
I myThid ) |
388 |
|
389 |
C !DESCRIPTION: \bv |
390 |
C *==========================================================* |
391 |
C | SUBROUTINE LAYERS_DIAPYCNAL |
392 |
C | Calculate the diapycnal velocity in isotracer layers, for a chosen |
393 |
C | tracer. |
394 |
C *==========================================================* |
395 |
C \ev |
396 |
IMPLICIT NONE |
397 |
#include "SIZE.h" |
398 |
#include "EEPARAMS.h" |
399 |
#include "PARAMS.h" |
400 |
#include "GRID.h" |
401 |
#include "LAYERS_SIZE.h" |
402 |
#include "LAYERS.h" |
403 |
|
404 |
C !INPUT PARAMETERS: |
405 |
C myThid :: my Thread Id number |
406 |
C tracer :: potential temperature, salt or potential density prho |
407 |
C iLa :: layer coordinate index |
408 |
C TtendSurf :: temperature tendency due to surface forcing times thickness |
409 |
C TtendDiffh:: temperature tendency due to horiz. diffusion times thickness |
410 |
C TtendDiffr:: temperature tendency due to vert. diffusion times thickness |
411 |
C TtendAdvh:: salinity tendency due to horiz. advection times thickness |
412 |
C TtendAdvr:: salinity tendency due to vert. advection times thickness |
413 |
C StendSurf :: salinity tendency due to surface forcing times thickness |
414 |
C StendDiffh:: salinity tendency due to horiz. diffusion times thickness |
415 |
C StendDiffr:: salinity tendency due to vert. diffusion times thickness |
416 |
C StendAdvh :: salinity tendency due to horiz. advection times thickness |
417 |
C StendAdvr :: salinity tendency due to vert. advection times thickness |
418 |
C Hc :: Layer thickness at the tracer point (m) |
419 |
C PIw :: 1 if layer exists, 0 otherwise (at tracer point) |
420 |
INTEGER iLa, myThid |
421 |
_RL tracer (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr, nSx,nSy) |
422 |
_RL Hc (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy) |
423 |
_RL PIc (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy) |
424 |
_RL TtendSurf (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers-1,nSx,nSy) |
425 |
_RL TtendDiffh(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers-1,nSx,nSy) |
426 |
_RL TtendDiffr(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers-1,nSx,nSy) |
427 |
_RL TtendAdvh (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers-1,nSx,nSy) |
428 |
_RL TtendAdvr (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers-1,nSx,nSy) |
429 |
_RL Ttendtot (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers-1,nSx,nSy) |
430 |
_RL StendSurf (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers-1,nSx,nSy) |
431 |
_RL StendDiffh(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers-1,nSx,nSy) |
432 |
_RL StendDiffr(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers-1,nSx,nSy) |
433 |
_RL StendAdvh (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers-1,nSx,nSy) |
434 |
_RL StendAdvr (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers-1,nSx,nSy) |
435 |
_RL Stendtot (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers-1,nSx,nSy) |
436 |
CEOP |
437 |
|
438 |
#ifdef LAYERS_THERMODYNAMICS |
439 |
|
440 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
441 |
C !LOCAL VARIABLES: |
442 |
C bi, bj :: tile indices |
443 |
C i,j :: horizontal indices |
444 |
C k :: vertical index for model grid |
445 |
C kp1 :: vertical index for model grid next cell |
446 |
C kci :: index from CellIndex |
447 |
C kg :: index for looping though layers_bounds |
448 |
C kk :: vertical index for ZZ (fine) grid |
449 |
C kloc :: local copy of kgu/v to reduce accesses to index arrays |
450 |
C mSteps :: maximum number of steps for bisection method |
451 |
C TatC :: temperature at C point |
452 |
_RL Hcw (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers-1,nSx,nSy) |
453 |
INTEGER bi, bj |
454 |
INTEGER i,j,k,kk,kg,kci,kloc |
455 |
INTEGER mSteps |
456 |
INTEGER kgc(sNx+1,sNy+1) |
457 |
INTEGER kgcw(sNx+1,sNy+1) |
458 |
_RL TatC(sNx+1,sNy+1), dzfac, Tfac, Sfac |
459 |
LOGICAL errorFlag |
460 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
461 |
#ifdef LAYERS_FINEGRID_DIAPYCNAL |
462 |
INTEGER kp1 |
463 |
#endif |
464 |
|
465 |
C -- constants for T and S forcing, gets reset later for rho |
466 |
Tfac = 1. _d 0 |
467 |
Sfac = 1. _d 0 |
468 |
|
469 |
C compute maximum number of steps for bisection method (approx. |
470 |
C log2(Nlayers)) as log2(Nlayers) + 1 for safety |
471 |
mSteps = int(log10(dble(Nlayers))/log10(2. _d 0))+1 |
472 |
|
473 |
C STOP 'DEBUG END: S/R LAYERS_DIAPYCNAL' |
474 |
|
475 |
C --- The tile loops |
476 |
DO bj=myByLo(myThid),myByHi(myThid) |
477 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
478 |
|
479 |
C Initialize the search indices |
480 |
DO j = 1,sNy+1 |
481 |
DO i = 1,sNx+1 |
482 |
C The temperature index (layer_G) goes from cold to warm. |
483 |
C The water column goes from warm (k=1) to cold (k=Nr). |
484 |
C So initialize the search with the warmest value. |
485 |
kgc(i,j) = Nlayers |
486 |
kgcw(i,j) = Nlayers-1 |
487 |
ENDDO |
488 |
ENDDO |
489 |
|
490 |
C Reset the arrays |
491 |
C --- These are at the w point |
492 |
DO kg=1,Nlayers-1 |
493 |
DO j = 1-OLy,sNy+OLy |
494 |
DO i = 1-OLx,sNx+OLx |
495 |
TtendSurf (i,j,kg,bi,bj) = 0. _d 0 |
496 |
TtendDiffh(i,j,kg,bi,bj) = 0. _d 0 |
497 |
TtendDiffr(i,j,kg,bi,bj) = 0. _d 0 |
498 |
TtendAdvh(i,j,kg,bi,bj) = 0. _d 0 |
499 |
TtendAdvr(i,j,kg,bi,bj) = 0. _d 0 |
500 |
Ttendtot(i,j,kg,bi,bj) = 0. _d 0 |
501 |
StendSurf (i,j,kg,bi,bj) = 0. _d 0 |
502 |
StendDiffh(i,j,kg,bi,bj) = 0. _d 0 |
503 |
StendDiffr(i,j,kg,bi,bj) = 0. _d 0 |
504 |
StendAdvh(i,j,kg,bi,bj) = 0. _d 0 |
505 |
StendAdvr(i,j,kg,bi,bj) = 0. _d 0 |
506 |
Stendtot(i,j,kg,bi,bj) = 0. _d 0 |
507 |
Hcw(i,j,kg,bi,bj) = 0. _d 0 |
508 |
ENDDO |
509 |
ENDDO |
510 |
ENDDO |
511 |
C --- These are at the c point |
512 |
DO kg=1,Nlayers |
513 |
DO j = 1-OLy,sNy+OLy |
514 |
DO i = 1-OLx,sNx+OLx |
515 |
Hc(i,j,kg,bi,bj) = 0. _d 0 |
516 |
PIc(i,j,kg,bi,bj) = 0. _d 0 |
517 |
ENDDO |
518 |
ENDDO |
519 |
ENDDO |
520 |
|
521 |
#ifdef LAYERS_FINEGRID_DIAPYCNAL |
522 |
DO kk=1,NZZ |
523 |
k = MapIndex(kk) |
524 |
kci = CellIndex(kk) |
525 |
DO j = 1,sNy+1 |
526 |
DO i = 1,sNx+1 |
527 |
C ------ Find theta at the V point (south) on the fine Z grid |
528 |
kp1=k+1 |
529 |
IF (maskC(i,j,kp1,bi,bj).EQ.zeroRS) kp1=k |
530 |
TatC(i,j) = MapFact(kk) * tracer(i,j,k,bi,bj) + |
531 |
& (1. _d 0 -MapFact(kk)) * tracer(i,j,kp1,bi,bj) |
532 |
ENDDO |
533 |
ENDDO |
534 |
#else |
535 |
DO kk=1,Nr |
536 |
k = kk |
537 |
kci = kk |
538 |
DO j = 1,sNy+1 |
539 |
DO i = 1,sNx+1 |
540 |
TatC(i,j) = tracer(i,j,k,bi,bj) |
541 |
ENDDO |
542 |
ENDDO |
543 |
#endif /* LAYERS_FINEGRID_DIAPYCNAL */ |
544 |
|
545 |
C ------ debugging stuff |
546 |
c IF (i.EQ.38 .AND. j.EQ.4 .AND. bi.EQ.1 .AND. bj.EQ.1) THEN |
547 |
c i=38 |
548 |
c j=4 |
549 |
c WRITE(msgBuf, |
550 |
c & '(A,I3,A,I3,A,I3,A,F7.2,A,F7.2,A,F7.2,A,F7.2,A,F3.1)') |
551 |
c & 'LAYERS_DEBUG: iLa=', iLa, |
552 |
c & ', kk=', kk, |
553 |
c & ', k=', k, |
554 |
c & ', tracer=', tracer(i,j,k,bi,bj), |
555 |
c & ', TatC=',TatC(i,j), |
556 |
c & ', hFacC=',hFacC(i,j,k,bi,bj) |
557 |
c CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
558 |
c & SQUEEZE_RIGHT, myThid ) |
559 |
c ENDIF |
560 |
C ------ Now that we know T everywhere, determine the binning. |
561 |
C find the layer indices kgc for the center point |
562 |
CALL LAYERS_LOCATE( |
563 |
I layers_bounds(1,iLa),Nlayers,mSteps,sNx,sNy,TatC, |
564 |
O kgc, |
565 |
I myThid ) |
566 |
#ifndef TARGET_NEC_SX |
567 |
C check for failures |
568 |
IF ( debugLevel .GE. debLevC ) THEN |
569 |
errorFlag = .FALSE. |
570 |
DO j = 1,sNy+1 |
571 |
DO i = 1,sNx+1 |
572 |
IF ( kgc(i,j) .LE. 0 ) THEN |
573 |
WRITE(msgBuf,'(2A,I3,A,I3,A,1E14.6)') |
574 |
& 'S/R LAYERS_LOCATE: Could not find a bin in ', |
575 |
& 'layers_bounds for TatC(',i,',',j,',)=',TatC(i,j) |
576 |
CALL PRINT_ERROR( msgBuf, myThid ) |
577 |
errorFlag = .TRUE. |
578 |
ENDIF |
579 |
ENDDO |
580 |
ENDDO |
581 |
IF ( errorFlag ) STOP 'ABNORMAL END: S/R LAYERS_DIAPYCNAL' |
582 |
ENDIF |
583 |
#endif /* ndef TARGET_NEC_SX */ |
584 |
|
585 |
C find the layer indices kgcw for the w point |
586 |
CALL LAYERS_LOCATE( |
587 |
I layers_bounds_w(1,iLa),Nlayers-1,mSteps,sNx,sNy,TatC, |
588 |
O kgcw, |
589 |
I myThid ) |
590 |
#ifndef TARGET_NEC_SX |
591 |
C check for failures |
592 |
IF ( debugLevel .GE. debLevC ) THEN |
593 |
errorFlag = .FALSE. |
594 |
DO j = 1,sNy+1 |
595 |
DO i = 1,sNx+1 |
596 |
IF ( kgcw(i,j) .LE. 0 ) THEN |
597 |
WRITE(msgBuf,'(2A,I3,A,I3,A,1E14.6)') |
598 |
& 'S/R LAYERS_LOCATE: Could not find a bin in ', |
599 |
& 'layers_bounds for TatC(',i,',',j,',)=',TatC(i,j) |
600 |
CALL PRINT_ERROR( msgBuf, myThid ) |
601 |
errorFlag = .TRUE. |
602 |
ENDIF |
603 |
ENDDO |
604 |
ENDDO |
605 |
IF ( errorFlag ) STOP 'ABNORMAL END: S/R LAYERS_DIAPYCNAL' |
606 |
ENDIF |
607 |
#endif /* ndef TARGET_NEC_SX */ |
608 |
|
609 |
C ------ Augment the bin values |
610 |
DO j = 1,sNy+1 |
611 |
DO i = 1,sNx+1 |
612 |
#ifdef LAYERS_FINEGRID_DIAPYCNAL |
613 |
dzfac = dZZf(kk) * hFacC(i,j,kci,bi,bj) |
614 |
#else |
615 |
dzfac = dRf(kci) * hFacC(i,j,kci,bi,bj) |
616 |
#endif /* LAYERS_FINEGRID_DIAPYCNAL */ |
617 |
kloc = kgcw(i,j) |
618 |
|
619 |
C ------- Thickness at w point |
620 |
Hcw(i,j,kloc,bi,bj) = Hcw(i,j,kloc,bi,bj) |
621 |
& + dzfac |
622 |
C ------- Thickness at c point |
623 |
Hc(i,j,kgc(i,j),bi,bj) = Hc(i,j,kgc(i,j),bi,bj) |
624 |
& + dzfac |
625 |
|
626 |
C ------- Now rescale dzfac to include the layer coordinate spacing |
627 |
dzfac = dzfac * layers_recip_delta(kloc,iLa) |
628 |
|
629 |
IF ( layers_num(iLa) .EQ. 3 ) THEN |
630 |
Tfac = layers_alpha(i,j,kci,bi,bj) |
631 |
Sfac = layers_beta(i,j,kci,bi,bj) |
632 |
ENDIF |
633 |
IF (kci.EQ.1) THEN |
634 |
C ------- We are in the surface layer |
635 |
TtendSurf(i,j,kloc,bi,bj) = |
636 |
& TtendSurf(i,j,kloc,bi,bj) + |
637 |
& Tfac * dzfac * layers_surfflux(i,j,1,1,bi,bj) |
638 |
StendSurf(i,j,kloc,bi,bj) = |
639 |
& StendSurf(i,j,kloc,bi,bj) + |
640 |
& Sfac * dzfac * layers_surfflux(i,j,1,2,bi,bj) |
641 |
ENDIF |
642 |
|
643 |
#ifdef SHORTWAVE_HEATING |
644 |
TtendSurf(i,j,kloc,bi,bj) = |
645 |
& TtendSurf(i,j,kloc,bi,bj) + |
646 |
& Tfac * dzfac * layers_sw(i,j,kci,1,bi,bj) |
647 |
#endif /* SHORTWAVE_HEATING */ |
648 |
|
649 |
C ------- Diffusion |
650 |
TtendDiffh(i,j,kloc,bi,bj) = |
651 |
& TtendDiffh(i,j,kloc,bi,bj) + dzfac * Tfac * |
652 |
& (layers_dfx(i,j,kci,1,bi,bj)+ |
653 |
& layers_dfy(i,j,kci,1,bi,bj)) |
654 |
TtendDiffr(i,j,kloc,bi,bj) = |
655 |
& TtendDiffr(i,j,kloc,bi,bj) + |
656 |
& dzfac * Tfac * layers_dfr(i,j,kci,1,bi,bj) |
657 |
StendDiffh(i,j,kloc,bi,bj) = |
658 |
& StendDiffh(i,j,kloc,bi,bj) + dzfac * Sfac * |
659 |
& (layers_dfx(i,j,kci,2,bi,bj)+ |
660 |
& layers_dfy(i,j,kci,2,bi,bj)) |
661 |
StendDiffr(i,j,kloc,bi,bj) = |
662 |
& StendDiffr(i,j,kloc,bi,bj) + |
663 |
& dzfac * Sfac * layers_dfr(i,j,kci,2,bi,bj) |
664 |
C ------- Advection |
665 |
TtendAdvh(i,j,kloc,bi,bj) = |
666 |
& TtendAdvh(i,j,kloc,bi,bj) + dzfac * Tfac * |
667 |
& (layers_afx(i,j,kci,1,bi,bj)+ |
668 |
& layers_afy(i,j,kci,1,bi,bj)) |
669 |
TtendAdvr(i,j,kloc,bi,bj) = |
670 |
& TtendAdvr(i,j,kloc,bi,bj) + |
671 |
& dzfac * Tfac * layers_afr(i,j,kci,1,bi,bj) |
672 |
StendAdvh(i,j,kloc,bi,bj) = |
673 |
& StendAdvh(i,j,kloc,bi,bj) + dzfac * Sfac * |
674 |
& (layers_afx(i,j,kci,2,bi,bj)+ |
675 |
& layers_afy(i,j,kci,2,bi,bj)) |
676 |
StendAdvr(i,j,kloc,bi,bj) = |
677 |
& StendAdvr(i,j,kloc,bi,bj) + |
678 |
& dzfac * Sfac * layers_afr(i,j,kci,2,bi,bj) |
679 |
C -------- Total Tendency |
680 |
Ttendtot(i,j,kloc,bi,bj) = |
681 |
& Ttendtot(i,j,kloc,bi,bj) + |
682 |
& dzfac * Tfac * layers_tottend(i,j,kci,1,bi,bj) |
683 |
Stendtot(i,j,kloc,bi,bj) = |
684 |
& Stendtot(i,j,kloc,bi,bj) + |
685 |
& dzfac * Sfac * layers_tottend(i,j,kci,2,bi,bj) |
686 |
ENDDO |
687 |
ENDDO |
688 |
ENDDO |
689 |
|
690 |
C-- Now that we know the thicknesses, compute the heaviside function |
691 |
C-- (Needs another loop through Ng) |
692 |
DO kg=1,Nlayers |
693 |
DO j = 1,sNy+1 |
694 |
DO i = 1,sNx+1 |
695 |
IF (Hc(i,j,kg,bi,bj) .GT. 0.) THEN |
696 |
PIc(i,j,kg,bi,bj) = 1. _d 0 |
697 |
ENDIF |
698 |
ENDDO |
699 |
ENDDO |
700 |
ENDDO |
701 |
|
702 |
C --- End bi,bj loop |
703 |
ENDDO |
704 |
ENDDO |
705 |
|
706 |
#endif /* LAYERS_THERMODYNAMICS */ |
707 |
|
708 |
RETURN |
709 |
END |
710 |
|
711 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
712 |
|
713 |
CBOP |
714 |
SUBROUTINE LAYERS_LOCATE( |
715 |
I xx,n,m,sNx,sNy,x, |
716 |
O k, |
717 |
I myThid ) |
718 |
|
719 |
C !DESCRIPTION: \bv |
720 |
C *==========================================================* |
721 |
C | Find the index(-array) k such that x is between xx(k) |
722 |
C | and xx(k+1) by bisection, following Press et al., |
723 |
C | Numerical Recipes in Fortran. xx must be monotonic. |
724 |
C *==========================================================* |
725 |
C \ev |
726 |
|
727 |
C !USES: |
728 |
IMPLICIT NONE |
729 |
C !INPUT PARAMETERS: |
730 |
C xx :: array of bin-boundaries (layers_boundaries) |
731 |
C n :: length of xx |
732 |
C m :: int(log2(n)) + 1 = length of bisection loop |
733 |
C sNx,sNy :: size of index array and input x |
734 |
C x :: input array of values |
735 |
C k :: index array (output) |
736 |
C myThid :: my Thread Id number |
737 |
INTEGER n,m,sNx,sNy |
738 |
_RL xx(1:n+1) |
739 |
_RL x(snx+1,sny+1) |
740 |
INTEGER k(snx+1,sny+1) |
741 |
INTEGER myThid |
742 |
|
743 |
C !LOCAL VARIABLES: |
744 |
C i,j :: horizontal indices |
745 |
C l :: bisection loop index |
746 |
C kl,ku,km :: work arrays and variables |
747 |
INTEGER i,j |
748 |
CEOP |
749 |
#ifdef TARGET_NEC_SX |
750 |
INTEGER l, km |
751 |
INTEGER kl(sNx+1,sNy+1), ku(sNx+1,sNy+1) |
752 |
|
753 |
C bisection, following Press et al., Numerical Recipes in Fortran, |
754 |
C mostly, because it can be vectorized |
755 |
DO j = 1,sNy+1 |
756 |
DO i = 1,sNx+1 |
757 |
kl(i,j)=1 |
758 |
ku(i,j)=n+1 |
759 |
END DO |
760 |
END DO |
761 |
DO l = 1,m |
762 |
DO j = 1,sNy+1 |
763 |
DO i = 1,sNx+1 |
764 |
IF (ku(i,j)-kl(i,j).GT.1) THEN |
765 |
km=(ku(i,j)+kl(i,j))/2 |
766 |
CML IF ((xx(n).GE.xx(1)).EQV.(x(i,j).GE.xx(km))) THEN |
767 |
IF ( ((xx(n).GE.xx(1)).AND.(x(i,j).GE.xx(km))).OR. |
768 |
& ((xx(n).GE.xx(1)).AND.(x(i,j).GE.xx(km))) ) THEN |
769 |
kl(i,j)=km |
770 |
ELSE |
771 |
ku(i,j)=km |
772 |
END IF |
773 |
END IF |
774 |
END DO |
775 |
END DO |
776 |
END DO |
777 |
DO j = 1,sNy+1 |
778 |
DO i = 1,sNx+1 |
779 |
IF ( x(i,j).LT.xx(2) ) THEN |
780 |
k(i,j)=1 |
781 |
ELSE IF ( x(i,j).GE.xx(n) ) THEN |
782 |
k(i,j)=n |
783 |
ELSE |
784 |
k(i,j)=kl(i,j) |
785 |
END IF |
786 |
END DO |
787 |
END DO |
788 |
#else |
789 |
C the old way |
790 |
DO j = 1,sNy+1 |
791 |
DO i = 1,sNx+1 |
792 |
IF (x(i,j) .GE. xx(n)) THEN |
793 |
C the point is in the hottest bin or hotter |
794 |
k(i,j) = n |
795 |
ELSE IF (x(i,j) .LT. xx(2)) THEN |
796 |
C the point is in the coldest bin or colder |
797 |
k(i,j) = 1 |
798 |
ELSE IF ( (x(i,j) .GE. xx(k(i,j))) |
799 |
& .AND. (x(i,j) .LT. xx(k(i,j)+1)) ) THEN |
800 |
C already on the right bin -- do nothing |
801 |
ELSE IF (x(i,j) .GE. xx(k(i,j))) THEN |
802 |
C have to hunt for the right bin by getting hotter |
803 |
DO WHILE (x(i,j) .GE. xx(k(i,j)+1)) |
804 |
k(i,j) = k(i,j) + 1 |
805 |
ENDDO |
806 |
C now xx(k) < x <= xx(k+1) |
807 |
ELSE IF (x(i,j) .LT. xx(k(i,j)+1)) THEN |
808 |
C have to hunt for the right bin by getting colder |
809 |
DO WHILE (x(i,j) .LT. xx(k(i,j))) |
810 |
k(i,j) = k(i,j) - 1 |
811 |
ENDDO |
812 |
C now xx(k) <= x < xx(k+1) |
813 |
ELSE |
814 |
C that should have covered all the options |
815 |
k(i,j) = -1 |
816 |
ENDIF |
817 |
|
818 |
ENDDO |
819 |
ENDDO |
820 |
#endif /* TARGET_NEC_SX */ |
821 |
|
822 |
#endif /* ALLOW_LAYERS */ |
823 |
|
824 |
RETURN |
825 |
END |
826 |
|