/[MITgcm]/MITgcm/pkg/layers/layers_fluxcalc.F
ViewVC logotype

Contents of /MITgcm/pkg/layers/layers_fluxcalc.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.9 - (show annotations) (download)
Sat Feb 9 14:52:29 2013 UTC (11 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64x, checkpoint64o, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64n, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l
Changes since 1.8: +41 -36 lines
- fix so that it also compiles with g77 & pgf77
- avoid unused variables

1 C $Header: /u/gcmpack/MITgcm/pkg/layers/layers_fluxcalc.F,v 1.8 2013/02/08 07:56:14 mlosch 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_LOCATE
13
14 CBOP 0
15 C !ROUTINE: LAYERS_FLUXCALC
16 C !INTERFACE:
17 SUBROUTINE LAYERS_FLUXCALC(
18 I uVel,vVel,tracer,iLa,
19 O UH,VH,Hw,Hs,PIw,PIs,U,V,
20 I myThid )
21
22 C !DESCRIPTION: \bv
23 C *==========================================================*
24 C | SUBROUTINE LAYERS_FLUXCALC
25 C | Calculate the transport in isotracer layers, for a chosen
26 C | tracer. This is the meat of the LAYERS package.
27 C *==========================================================*
28 C \ev
29
30 C !USES:
31 IMPLICIT NONE
32 C == Global variables ===
33 #include "SIZE.h"
34 #include "EEPARAMS.h"
35 #include "PARAMS.h"
36 #include "GRID.h"
37 #include "LAYERS_SIZE.h"
38 #include "LAYERS.h"
39 #ifdef ALLOW_GMREDI
40 # include "GMREDI.h"
41 #endif
42
43 C !INPUT PARAMETERS:
44 C myThid :: my Thread Id number
45 C uVel :: zonal velocity (m/s, i=1 held at western face)
46 C vVel :: meridional velocity (m/s, j=1 held at southern face)
47 C tracer :: potential temperature, salt or potential density prho
48 C UH :: U integrated over layer (m^2/s)
49 C VH :: V integrated over layer (m^2/s)
50 C Hw :: Layer thickness at the U point (m)
51 C Hs :: Layer thickness at the V point (m)
52 C PIw :: 1 if layer exists, 0 otherwise (at U point)
53 C PIs :: 1 if layer exists, 0 otherwise (at V point)
54 C U :: average U over layer (m/s)
55 C V :: average V over layer (m/s)
56 C iLa :: layer coordinate index
57 INTEGER myThid
58 _RL uVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr, nSx,nSy)
59 _RL vVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr, nSx,nSy)
60 _RL tracer (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr, nSx,nSy)
61 _RL UH (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
62 _RL VH (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
63 _RL Hw (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
64 _RL Hs (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
65 _RL PIw (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
66 _RL PIs (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
67 _RL U (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
68 _RL V (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nlayers,nSx,nSy)
69 INTEGER iLa
70 CEOP
71
72 #ifdef ALLOW_LAYERS
73
74 C !LOCAL VARIABLES:
75 C bi, bj :: tile indices
76 C i,j :: horizontal indices
77 C k :: vertical index for model grid
78 C kci :: index from CellIndex
79 C kg :: index for looping though layers_bounds
80 C kk :: vertical index for ZZ (fine) grid
81 C kgu,kgv :: vertical index for isopycnal grid
82 C kloc :: local copy of kgu/v to reduce accesses to index arrays
83 C mSteps :: maximum number of steps for bisection method
84 C prho :: potential density referenced to layers_krho pressure
85 C TatU :: temperature at U point
86 C TatV :: temperature at V point
87
88 INTEGER bi, bj
89 INTEGER i,j,k,kk,kg,kci,kp1,kloc
90 INTEGER mSteps
91 INTEGER kgu(sNx+1,sNy+1), kgv(sNx+1,sNy+1)
92 _RL TatU(sNx+1,sNy+1), TatV(sNx+1,sNy+1)
93 #ifdef ALLOW_GMREDI
94 INTEGER kcip1
95 _RL delPsi, maskp1
96 #endif
97 LOGICAL errorFlag
98 CHARACTER*(MAX_LEN_MBUF) msgBuf
99
100 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
101
102 C compute maximum number of steps for bisection method (approx.
103 C log2(Nlayers)) as log2(Nlayers) + 1 for safety
104 mSteps = int(log10(dble(Nlayers))/log10(2. _d 0))+1
105
106 C --- The tile loops
107 DO bj=myByLo(myThid),myByHi(myThid)
108 DO bi=myBxLo(myThid),myBxHi(myThid)
109
110 C Initialize the search indices
111 DO j = 1,sNy+1
112 DO i = 1,sNx+1
113 C The temperature index (layer_G) goes from cold to warm.
114 C The water column goes from warm (k=1) to cold (k=Nr).
115 C So initialize the search with the warmest value.
116 kgu(i,j) = Nlayers
117 kgv(i,j) = Nlayers
118 ENDDO
119 ENDDO
120
121 C Reset the arrays
122 DO kg=1,Nlayers
123 DO j = 1-OLy,sNy+OLy
124 DO i = 1-OLx,sNx+OLx
125 #ifdef LAYERS_UFLUX
126 UH(i,j,kg,bi,bj) = 0. _d 0
127 #ifdef LAYERS_THICKNESS
128 Hw(i,j,kg,bi,bj) = 0. _d 0
129 #endif /* LAYERS_THICKNESS */
130 #endif /* UH */
131 #ifdef LAYERS_VFLUX
132 VH(i,j,kg,bi,bj) = 0. _d 0
133 #ifdef LAYERS_THICKNESS
134 Hs(i,j,kg,bi,bj) = 0. _d 0
135 #endif /* LAYERS_THICKNESS */
136 #endif /* VH */
137 ENDDO
138 ENDDO
139 ENDDO
140
141 DO kk=1,NZZ
142 k = MapIndex(kk)
143 kci = CellIndex(kk)
144 #ifdef ALLOW_GMREDI
145 kcip1 = MIN(kci+1,Nr)
146 maskp1 = 1.
147 IF (kci.GE.Nr) maskp1 = 0.
148 #endif /* ALLOW_GMREDI */
149 #ifdef LAYERS_UFLUX
150 DO j = 1,sNy+1
151 DO i = 1,sNx+1
152
153 C ------ Find theta at the U point (west) on the fine Z grid
154 kp1=k+1
155 IF (hFacW(i,j,kp1,bi,bj) .EQ. 0.) kp1=k
156 TatU(i,j) = MapFact(kk) *
157 & 0.5 _d 0 * (tracer(i-1,j,k,bi,bj)+tracer(i,j,k,bi,bj)) +
158 & (1-MapFact(kk)) *
159 & 0.5 _d 0 * (tracer(i-1,j,kp1,bi,bj)+tracer(i,j,kp1,bi,bj))
160
161 ENDDO
162 ENDDO
163 C ------ Now that we know T everywhere, determine the binning.
164 C find the layer indices kgu
165 CALL LAYERS_LOCATE(
166 I layers_bounds(1,iLa),Nlayers,mSteps,sNx,sNy,TatU,
167 O kgu,
168 I myThid )
169 #ifndef TARGET_NEC_SX
170 C check for failures
171 IF ( debugLevel .GE. debLevC ) THEN
172 errorFlag = .FALSE.
173 DO j = 1,sNy+1
174 DO i = 1,sNx+1
175 IF ( kgu(i,j) .LE. 0 ) THEN
176 WRITE(msgBuf,'(2A,I3,A,I3,A,1E14.6)')
177 & 'S/R LAYERS_LOCATE: Could not find a bin in ',
178 & 'layers_bounds for TatU(',i,',',j,',)=',TatU(i,j)
179 CALL PRINT_ERROR( msgBuf, myThid )
180 errorFlag = .TRUE.
181 ENDIF
182 ENDDO
183 ENDDO
184 IF ( errorFlag ) STOP 'ABNORMAL END: S/R LAYERS_FLUXCALC'
185 ENDIF
186 #endif /* ndef TARGET_NEC_SX */
187 C
188 DO j = 1,sNy+1
189 DO i = 1,sNx+1
190
191 kloc = kgu(i,j)
192 C ------ Augment the bin values
193 UH(i,j,kloc,bi,bj) =
194 & UH(i,j,kloc,bi,bj) +
195 & dZZf(kk) * uVel(i,j,kci,bi,bj) * hFacW(i,j,kci,bi,bj)
196
197 #ifdef ALLOW_GMREDI
198 IF ( layers_bolus(iLa) ) THEN
199 IF ( .NOT.GM_AdvForm ) THEN
200 delPsi = 0.25 _d 0 *(
201 & ( rA(i-1,j,bi,bj)*Kwx(i-1,j,kcip1,bi,bj)
202 & +rA( i ,j,bi,bj)*Kwx( i ,j,kcip1,bi,bj)
203 & ) * maskW(i,j,kcip1,bi,bj) * maskp1
204 & - ( rA(i-1,j,bi,bj)*Kwx(i-1,j, kci ,bi,bj)
205 & +rA( i ,j,bi,bj)*Kwx( i ,j, kci ,bi,bj)
206 & ) * maskW(i,j, kci ,bi,bj)
207 & ) * recip_rAw(i,j,bi,bj)
208 #ifdef GM_BOLUS_ADVEC
209 ELSE
210 delPsi = GM_PsiX(i,j,kcip1,bi,bj)*maskp1
211 & - GM_PsiX(i,j, kci, bi,bj)
212 #endif
213 ENDIF
214 UH(i,j,kloc,bi,bj) = UH(i,j,kloc,bi,bj)
215 & + delPsi*recip_drF(kci)*_recip_hFacW(i,j,kci,bi,bj)
216 & * dZZf(kk)*hFacW(i,j,kci,bi,bj)
217 ENDIF
218 #endif /* ALLOW_GMREDI */
219
220 #ifdef LAYERS_THICKNESS
221 Hw(i,j,kloc,bi,bj) = Hw(i,j,kloc,bi,bj)
222 & + dZZf(kk) * hFacW(i,j,kci,bi,bj)
223 #endif /* LAYERS_THICKNESS */
224
225 ENDDO
226 ENDDO
227 #endif /* LAYERS_UFLUX */
228
229 #ifdef LAYERS_VFLUX
230 DO j = 1,sNy+1
231 DO i = 1,sNx+1
232 C ------ Find theta at the V point (south) on the fine Z grid
233 kp1=k+1
234 IF (hFacS(i,j,kp1,bi,bj) .EQ. 0.) kp1=k
235 TatV(i,j) = MapFact(kk) *
236 & 0.5 _d 0 * (tracer(i,j-1,k,bi,bj)+tracer(i,j,k,bi,bj)) +
237 & (1-MapFact(kk)) *
238 & 0.5 _d 0 * (tracer(i,j-1,kp1,bi,bj)+tracer(i,j,kp1,bi,bj))
239
240 ENDDO
241 ENDDO
242 C ------ Now that we know T everywhere, determine the binning.
243 C find the layer indices kgv
244 CALL LAYERS_LOCATE(
245 I layers_bounds(1,iLa),Nlayers,mSteps,sNx,sNy,TatV,
246 O kgv,
247 I myThid )
248 #ifndef TARGET_NEC_SX
249 IF ( debugLevel .GE. debLevC ) THEN
250 C check for failures
251 errorFlag = .FALSE.
252 DO j = 1,sNy+1
253 DO i = 1,sNx+1
254 IF ( kgv(i,j) .LE. 0 ) THEN
255 WRITE(msgBuf,'(2A,I3,A,I3,A,1E14.6)')
256 & 'S/R LAYERS_LOCATE: Could not find a bin in ',
257 & 'layers_bounds for TatV(',i,',',j,',)=',TatV(i,j)
258 CALL PRINT_ERROR( msgBuf, myThid )
259 errorFlag = .TRUE.
260 ENDIF
261 ENDDO
262 ENDDO
263 IF ( errorFlag ) STOP 'ABNORMAL END: S/R LAYERS_FLUXCALC'
264 ENDIF
265 #endif /* ndef TARGET_NEC_SX */
266 C
267 DO j = 1,sNy+1
268 DO i = 1,sNx+1
269
270 kloc = kgv(i,j)
271 C ------ Augment the bin values
272 VH(i,j,kloc,bi,bj) =
273 & VH(i,j,kloc,bi,bj)
274 & + dZZf(kk) * vVel(i,j,kci,bi,bj) * hFacS(i,j,kci,bi,bj)
275
276 #ifdef ALLOW_GMREDI
277 IF ( layers_bolus(iLa) ) THEN
278 IF ( .NOT.GM_AdvForm ) THEN
279 delPsi = 0.25 _d 0 *(
280 & ( rA(i,j-1,bi,bj)*Kwy(i,j-1,kcip1,bi,bj)
281 & +rA(i, j ,bi,bj)*Kwy(i, j ,kcip1,bi,bj)
282 & ) * maskS(i,j,kcip1,bi,bj) * maskp1
283 & - ( rA(i,j-1,bi,bj)*Kwy(i,j-1, kci ,bi,bj)
284 & +rA(i, j ,bi,bj)*Kwy(i, j , kci ,bi,bj)
285 & ) * maskS(i,j, kci ,bi,bj)
286 & ) * recip_rAs(i,j,bi,bj)
287 #ifdef GM_BOLUS_ADVEC
288 ELSE
289 delPsi = GM_PsiY(i,j,kcip1,bi,bj)*maskp1
290 & - GM_PsiY(i,j, kci, bi,bj)
291 #endif
292 ENDIF
293 VH(i,j,kloc,bi,bj) = VH(i,j,kloc,bi,bj)
294 & + delPsi*recip_drF(kci)*_recip_hFacS(i,j,kci,bi,bj)
295 & * dZZf(kk)*hFacS(i,j,kci,bi,bj)
296 ENDIF
297 #endif /* ALLOW_GMREDI */
298
299 #ifdef LAYERS_THICKNESS
300 Hs(i,j,kloc,bi,bj) = Hs(i,j,kloc,bi,bj)
301 & + dZZf(kk) * hFacS(i,j,kci,bi,bj)
302 #endif /* LAYERS_THICKNESS */
303
304 ENDDO
305 ENDDO
306 #endif /* LAYERS_VFLUX */
307 ENDDO
308
309 C-- Now that we know the thicknesses, compute the heaviside function
310 C-- (Needs another loop through Ng)
311 #ifdef LAYERS_THICKNESS
312 DO kg=1,Nlayers
313 DO j = 1,sNy+1
314 DO i = 1,sNx+1
315 #ifdef LAYERS_UFLUX
316 IF (Hw(i,j,kg,bi,bj) .GT. 0.) THEN
317 PIw(i,j,kg,bi,bj) = 1. _d 0
318 U(i,j,kg,bi,bj) =
319 & UH(i,j,kg,bi,bj) / Hw(i,j,kg,bi,bj)
320 ENDIF
321 #endif /* LAYERS_UFLUX */
322 #ifdef LAYERS_VFLUX
323 IF (Hs(i,j,kg,bi,bj) .GT. 0.) THEN
324 PIs(i,j,kg,bi,bj) = 1. _d 0
325 V(i,j,kg,bi,bj) =
326 & VH(i,j,kg,bi,bj) / Hs(i,j,kg,bi,bj)
327 ENDIF
328 #endif /* LAYERS_VFLUX */
329 ENDDO
330 ENDDO
331 ENDDO
332 #endif /* LAYERS_THICKNESS */
333
334 C --- End bi,bj loop
335 ENDDO
336 ENDDO
337
338 RETURN
339 END
340
341 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
342
343 CBOP
344 SUBROUTINE LAYERS_LOCATE(
345 I xx,n,m,sNx,sNy,x,
346 O k,
347 I myThid )
348
349 C !DESCRIPTION: \bv
350 C *==========================================================*
351 C | Find the index(-array) k such that x is between xx(k)
352 C | and xx(k+1) by bisection, following Press et al.,
353 C | Numerical Recipes in Fortran. xx must be monotonic.
354 C *==========================================================*
355 C \ev
356
357 C !USES:
358 IMPLICIT NONE
359 C !INPUT PARAMETERS:
360 C xx :: array of bin-boundaries (layers_boundaries)
361 C n :: length of xx
362 C m :: int(log2(n)) + 1 = length of bisection loop
363 C sNx,sNy :: size of index array and input x
364 C x :: input array of values
365 C k :: index array (output)
366 C myThid :: my Thread Id number
367 INTEGER n,m,sNx,sNy
368 _RL xx(1:n+1)
369 _RL x(snx+1,sny+1)
370 INTEGER k(snx+1,sny+1)
371 INTEGER myThid
372
373 C !LOCAL VARIABLES:
374 C i,j :: horizontal indices
375 C l :: bisection loop index
376 C kl,ku,km :: work arrays and variables
377 INTEGER i,j
378 CEOP
379 #ifdef TARGET_NEC_SX
380 INTEGER l, km
381 INTEGER kl(sNx+1,sNy+1), ku(sNx+1,sNy+1)
382
383 C bisection, following Press et al., Numerical Recipes in Fortran,
384 C mostly, because it can be vectorized
385 DO j = 1,sNy+1
386 DO i = 1,sNx+1
387 kl(i,j)=1
388 ku(i,j)=n+1
389 END DO
390 END DO
391 DO l = 1,m
392 DO j = 1,sNy+1
393 DO i = 1,sNx+1
394 IF (ku(i,j)-kl(i,j).GT.1) THEN
395 km=(ku(i,j)+kl(i,j))/2
396 CML IF ((xx(n).GE.xx(1)).EQV.(x(i,j).GE.xx(km))) THEN
397 IF ( ((xx(n).GE.xx(1)).AND.(x(i,j).GE.xx(km))).OR.
398 & ((xx(n).GE.xx(1)).AND.(x(i,j).GE.xx(km))) ) THEN
399 kl(i,j)=km
400 ELSE
401 ku(i,j)=km
402 END IF
403 END IF
404 END DO
405 END DO
406 END DO
407 DO j = 1,sNy+1
408 DO i = 1,sNx+1
409 IF ( x(i,j).LT.xx(2) ) THEN
410 k(i,j)=1
411 ELSE IF ( x(i,j).GE.xx(n) ) THEN
412 k(i,j)=n
413 ELSE
414 k(i,j)=kl(i,j)
415 END IF
416 END DO
417 END DO
418 #else
419 C the old way
420 DO j = 1,sNy+1
421 DO i = 1,sNx+1
422 IF (x(i,j) .GE. xx(n)) THEN
423 C the point is in the hottest bin or hotter
424 k(i,j) = n
425 ELSE IF (x(i,j) .LT. xx(2)) THEN
426 C the point is in the coldest bin or colder
427 k(i,j) = 1
428 ELSE IF ( (x(i,j) .GE. xx(k(i,j)))
429 & .AND. (x(i,j) .LT. xx(k(i,j)+1)) ) THEN
430 C already on the right bin -- do nothing
431 ELSE IF (x(i,j) .GE. xx(k(i,j))) THEN
432 C have to hunt for the right bin by getting hotter
433 DO WHILE (x(i,j) .GE. xx(k(i,j)+1))
434 k(i,j) = k(i,j) + 1
435 ENDDO
436 C now xx(k) < x <= xx(k+1)
437 ELSE IF (x(i,j) .LT. xx(k(i,j)+1)) THEN
438 C have to hunt for the right bin by getting colder
439 DO WHILE (x(i,j) .LT. xx(k(i,j)))
440 k(i,j) = k(i,j) - 1
441 ENDDO
442 C now xx(k) <= x < xx(k+1)
443 ELSE
444 C that should have covered all the options
445 k(i,j) = -1
446 ENDIF
447
448 ENDDO
449 ENDDO
450 #endif /* TARGET_NEC_SX */
451
452 #endif /* ALLOW_LAYERS */
453
454 RETURN
455 END

  ViewVC Help
Powered by ViewVC 1.1.22