1 |
C $Header: /u/gcmpack/MITgcm/pkg/layers/layers_calc.F,v 1.5 2010/12/04 23:50:32 dfer Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "LAYERS_OPTIONS.h" |
5 |
#ifdef ALLOW_GMREDI |
6 |
#include "GMREDI_OPTIONS.h" |
7 |
#endif |
8 |
|
9 |
CBOP 0 |
10 |
C !ROUTINE: LAYERS_CALC |
11 |
|
12 |
C !INTERFACE: |
13 |
SUBROUTINE LAYERS_CALC( |
14 |
I myTime, myIter, myThid ) |
15 |
|
16 |
C !DESCRIPTION: |
17 |
C =================================================================== |
18 |
C Calculate the transport in isopycnal layers. |
19 |
C This is the meat of the LAYERS package. |
20 |
C =================================================================== |
21 |
|
22 |
C !USES: |
23 |
IMPLICIT NONE |
24 |
#include "SIZE.h" |
25 |
#include "EEPARAMS.h" |
26 |
#include "PARAMS.h" |
27 |
#include "GRID.h" |
28 |
#include "DYNVARS.h" |
29 |
#include "LAYERS_SIZE.h" |
30 |
#include "LAYERS.h" |
31 |
#ifdef ALLOW_GMREDI |
32 |
# include "GMREDI.h" |
33 |
#endif |
34 |
|
35 |
C !INPUT PARAMETERS: |
36 |
C myTime :: Current time in simulation |
37 |
C myIter :: Current iteration number |
38 |
C myThid :: my Thread Id number |
39 |
_RL myTime |
40 |
INTEGER myIter |
41 |
INTEGER myThid |
42 |
CEOP |
43 |
|
44 |
#ifdef ALLOW_LAYERS |
45 |
|
46 |
C !LOCAL VARIABLES: |
47 |
C bi, bj :: tile indices |
48 |
C i,j :: horizontal indices |
49 |
C k :: vertical index for model grid |
50 |
C kci :: index from CellIndex |
51 |
C kg :: index for looping though layers_G |
52 |
C kk :: vertical index for ZZ (fine) grid |
53 |
C kgu,kgv :: vertical index for isopycnal grid |
54 |
C TatV :: temperature at U point |
55 |
C TatV :: temperature at V point |
56 |
|
57 |
INTEGER bi, bj |
58 |
INTEGER i,j,k,kk,kg,kci |
59 |
INTEGER kgu(sNx+1,sNy+1), kgv(sNx+1,sNy+1) |
60 |
_RL TatU, TatV |
61 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
62 |
#if (defined ALLOW_GMREDI) && (defined GM_BOLUS_ADVEC) |
63 |
INTEGER kcip1 |
64 |
_RL delPsi, maskp1 |
65 |
#endif |
66 |
|
67 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
68 |
|
69 |
C --- The tile loops |
70 |
DO bj=myByLo(myThid),myByHi(myThid) |
71 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
72 |
|
73 |
C Initialize the search indices |
74 |
DO j = 1,sNy+1 |
75 |
DO i = 1,sNx+1 |
76 |
C The temperature index (layer_G) goes from cold to warm. |
77 |
C The water column goes from warm (k=1) to cold (k=Nr). |
78 |
C So initialize the search with the warmest value. |
79 |
kgu(i,j) = Nlayers |
80 |
kgv(i,j) = Nlayers |
81 |
ENDDO |
82 |
ENDDO |
83 |
|
84 |
C Reset the arrays |
85 |
DO kg=1,Nlayers |
86 |
DO j = 1,sNy+1 |
87 |
DO i = 1,sNx+1 |
88 |
#ifdef LAYERS_UFLUX |
89 |
layers_UFlux(i,j,kg,bi,bj) = 0. _d 0 |
90 |
#ifdef LAYERS_THICKNESS |
91 |
layers_HU(i,j,kg,bi,bj) = 0. _d 0 |
92 |
#endif /* LAYERS_THICKNESS */ |
93 |
#endif /* LAYERS_UFLUX */ |
94 |
#ifdef LAYERS_VFLUX |
95 |
layers_VFlux(i,j,kg,bi,bj) = 0. _d 0 |
96 |
#ifdef LAYERS_THICKNESS |
97 |
layers_HV(i,j,kg,bi,bj) = 0. _d 0 |
98 |
#endif /* LAYERS_THICKNESS */ |
99 |
#endif /* LAYERS_VFLUX */ |
100 |
ENDDO |
101 |
ENDDO |
102 |
ENDDO |
103 |
|
104 |
C _RL theta(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
105 |
C Sometimes it is done this way |
106 |
C DO j=1-Oly+1,sNy+Oly-1 |
107 |
C DO i=1-Olx+1,sNx+Olx-1 |
108 |
DO kk=1,NZZ |
109 |
k = MapIndex(kk) |
110 |
kci = CellIndex(kk) |
111 |
DO j = 1,sNy+1 |
112 |
DO i = 1,sNx+1 |
113 |
|
114 |
#ifdef LAYERS_UFLUX |
115 |
C ------ Find theta at the U point (west) on the fine Z grid |
116 |
IF (LAYER_nb .EQ. 1) THEN |
117 |
TatU = MapFact(kk) * |
118 |
& 0.5 _d 0 * (theta(i-1,j,k,bi,bj)+theta(i,j,k,bi,bj)) + |
119 |
& (1-MapFact(kk)) * |
120 |
& 0.5 _d 0 * (theta(i-1,j,k+1,bi,bj)+theta(i,j,k+1,bi,bj)) |
121 |
ELSEIF (LAYER_nb .EQ. 2) THEN |
122 |
TatU = MapFact(kk) * |
123 |
& 0.5 _d 0 * (salt(i-1,j,k,bi,bj)+salt(i,j,k,bi,bj)) + |
124 |
& (1-MapFact(kk)) * |
125 |
& 0.5 _d 0 * (salt(i-1,j,k+1,bi,bj)+salt(i,j,k+1,bi,bj)) |
126 |
ENDIF |
127 |
|
128 |
C ------ Now that we know T everywhere, determine the binning. |
129 |
|
130 |
IF (TatU .GE. layers_G(Nlayers)) THEN |
131 |
C the point is in the hottest bin or hotter |
132 |
kgu(i,j) = Nlayers |
133 |
ELSE IF (TatU .LT. layers_G(2)) THEN |
134 |
C the point is in the coldest bin or colder |
135 |
kgu(i,j) = 1 |
136 |
ELSE IF ( (TatU .GE. layers_G(kgu(i,j))) |
137 |
& .AND. (TatU .LT. layers_G(kgu(i,j)+1)) ) THEN |
138 |
C already on the right bin -- do nothing |
139 |
ELSE IF (TatU .GE. layers_G(kgu(i,j))) THEN |
140 |
C have to hunt for the right bin by getting hotter |
141 |
DO WHILE (TatU .GE. layers_G(kgu(i,j)+1)) |
142 |
kgu(i,j) = kgu(i,j) + 1 |
143 |
ENDDO |
144 |
C now layers_G(kgu(i,j)+1) < TatU <= layers_G(kgu(i,j)+1) |
145 |
ELSE IF (TatU .LT. layers_G(kgu(i,j)+1)) THEN |
146 |
C have to hunt for the right bin by getting colder |
147 |
DO WHILE (TatU .LT. layers_G(kgu(i,j))) |
148 |
kgu(i,j) = kgu(i,j) - 1 |
149 |
ENDDO |
150 |
C now layers_G(kgu(i,j)+1) <= TatU < layers_G(kgu(i,j)+1) |
151 |
ELSE |
152 |
C that should have covered all the options |
153 |
WRITE(msgBuf,'(A,1E14.6)') |
154 |
& 'S/R LAYERS_CALC: Couldnt find a bin in layers_G for TatU=', |
155 |
& TatU |
156 |
CALL PRINT_ERROR( msgBuf, myThid ) |
157 |
STOP 'ABNORMAL END: S/R LAYERS_INIT_FIXED' |
158 |
END IF |
159 |
|
160 |
C ------ Augment the bin values |
161 |
layers_UFlux(i,j,kgu(i,j),bi,bj) = |
162 |
& layers_UFlux(i,j,kgu(i,j),bi,bj) + |
163 |
& dZZf(kk) * uVel(i,j,kci,bi,bj) * hFacW(i,j,kci,bi,bj) |
164 |
|
165 |
#if (defined ALLOW_GMREDI) && (defined GM_BOLUS_ADVEC) |
166 |
IF (GM_AdvForm) THEN |
167 |
kcip1 = MIN(kci+1,Nr) |
168 |
maskp1 = 1. |
169 |
IF (kci.GE.Nr) maskp1 = 0. |
170 |
delPsi = GM_PsiX(i,j,kcip1,bi,bj)*maskp1 |
171 |
& - GM_PsiX(i,j, kci, bi,bj) |
172 |
layers_UFlux(i,j,kgu(i,j),bi,bj) = |
173 |
& layers_UFlux(i,j,kgu(i,j),bi,bj) |
174 |
& + delPsi*recip_drF(kci)*_recip_hFacW(i,j,kci,bi,bj) |
175 |
& * dZZf(kk)*hFacW(i,j,kci,bi,bj) |
176 |
ENDIF |
177 |
#endif |
178 |
|
179 |
#ifdef LAYERS_THICKNESS |
180 |
layers_HU(i,j,kgu(i,j),bi,bj) = layers_HU(i,j,kgu(i,j),bi,bj) |
181 |
& + dZZf(kk) * hFacW(i,j,kci,bi,bj) |
182 |
#endif /* LAYERS_THICKNESS */ |
183 |
|
184 |
#endif /* LAYERS_UFLUX */ |
185 |
|
186 |
#ifdef LAYERS_VFLUX |
187 |
C ------ Find theta at the V point (south) on the fine Z grid |
188 |
IF (LAYER_nb .EQ. 1) THEN |
189 |
TatV = MapFact(kk) * |
190 |
& 0.5 _d 0 * (theta(i,j-1,k,bi,bj)+theta(i,j,k,bi,bj)) + |
191 |
& (1-MapFact(kk)) * |
192 |
& 0.5 _d 0 * (theta(i,j-1,k+1,bi,bj)+theta(i,j,k+1,bi,bj)) |
193 |
ELSEIF (LAYER_nb .EQ. 2) THEN |
194 |
TatV = MapFact(kk) * |
195 |
& 0.5 _d 0 * (salt(i,j-1,k,bi,bj)+salt(i,j,k,bi,bj)) + |
196 |
& (1-MapFact(kk)) * |
197 |
& 0.5 _d 0 * (salt(i,j-1,k+1,bi,bj)+salt(i,j,k+1,bi,bj)) |
198 |
ENDIF |
199 |
|
200 |
C ------ Now that we know T everywhere, determine the binning |
201 |
IF (TatV .GE. layers_G(Nlayers)) THEN |
202 |
C the point is in the hottest bin or hotter |
203 |
kgv(i,j) = Nlayers |
204 |
ELSE IF (TatV .LT. layers_G(2)) THEN |
205 |
C the point is in the coldest bin or colder |
206 |
kgv(i,j) = 1 |
207 |
ELSE IF ( (TatV .GE. layers_G(kgv(i,j))) |
208 |
& .AND. (TatV .LT. layers_G(kgv(i,j)+1)) ) THEN |
209 |
C already on the right bin -- do nothing |
210 |
ELSE IF (TatV .GE. layers_G(kgv(i,j))) THEN |
211 |
C have to hunt for the right bin by getting hotter |
212 |
DO WHILE (TatV .GE. layers_G(kgv(i,j)+1)) |
213 |
kgv(i,j) = kgv(i,j) + 1 |
214 |
ENDDO |
215 |
C now layers_G(kgv(i,j)+1) < TatV <= layers_G(kgv(i,j)+1) |
216 |
ELSE IF (TatV .LT. layers_G(kgv(i,j)+1)) THEN |
217 |
C have to hunt for the right bin by getting colder |
218 |
DO WHILE (TatV .LT. layers_G(kgv(i,j))) |
219 |
kgv(i,j) = kgv(i,j) - 1 |
220 |
ENDDO |
221 |
C now layers_G(kgv(i,j)+1) <= TatV < layers_G(kgv(i,j)+1) |
222 |
ELSE |
223 |
C that should have covered all the options |
224 |
WRITE(msgBuf,'(A,1E14.6)') |
225 |
& 'S/R LAYERS_CALC: Couldnt find a bin in layers_G for TatV=', |
226 |
& TatV |
227 |
CALL PRINT_ERROR( msgBuf, myThid ) |
228 |
STOP 'ABNORMAL END: S/R LAYERS_INIT_FIXED' |
229 |
END IF |
230 |
|
231 |
C ------ Augment the bin values |
232 |
layers_VFlux(i,j,kgv(i,j),bi,bj) = |
233 |
& layers_VFlux(i,j,kgv(i,j),bi,bj) |
234 |
& + dZZf(kk) * vVel(i,j,kci,bi,bj) * hFacS(i,j,kci,bi,bj) |
235 |
|
236 |
#if (defined ALLOW_GMREDI) && (defined GM_BOLUS_ADVEC) |
237 |
IF (GM_AdvForm) THEN |
238 |
kcip1 = MIN(kci+1,Nr) |
239 |
maskp1 = 1. |
240 |
IF (kci.GE.Nr) maskp1 = 0. |
241 |
delPsi = GM_PsiX(i,j,kcip1,bi,bj)*maskp1 |
242 |
& - GM_PsiX(i,j, kci, bi,bj) |
243 |
layers_UFlux(i,j,kgu(i,j),bi,bj) = |
244 |
& layers_UFlux(i,j,kgu(i,j),bi,bj) |
245 |
& + delPsi*recip_drF(kci)*_recip_hFacW(i,j,kci,bi,bj) |
246 |
& * dZZf(kk)*hFacW(i,j,kci,bi,bj) |
247 |
ENDIF |
248 |
#endif |
249 |
|
250 |
#ifdef LAYERS_THICKNESS |
251 |
layers_HV(i,j,kgv(i,j),bi,bj) = layers_HV(i,j,kgv(i,j),bi,bj) |
252 |
& + dZZf(kk) * hFacS(i,j,kci,bi,bj) |
253 |
#endif /* LAYERS_THICKNESS */ |
254 |
|
255 |
#endif /* LAYERS_VFLUX */ |
256 |
|
257 |
ENDDO |
258 |
ENDDO |
259 |
ENDDO |
260 |
|
261 |
#ifdef ALLOW_TIMEAVE |
262 |
C-- Time-average |
263 |
IF ( layers_taveFreq.GT.0. ) THEN |
264 |
|
265 |
#ifdef LAYERS_UFLUX |
266 |
CALL TIMEAVE_CUMULATE( layers_UFlux_T, layers_UFlux, Nlayers, |
267 |
& deltaTclock, bi, bj, myThid ) |
268 |
#ifdef LAYERS_THICKNESS |
269 |
CALL TIMEAVE_CUMULATE( layers_HU_T, layers_HU, Nlayers, |
270 |
& deltaTclock, bi, bj, myThid ) |
271 |
#endif /* LAYERS_THICKNESS */ |
272 |
#endif /* LAYERS_UFLUX */ |
273 |
#ifdef LAYERS_VFLUX |
274 |
CALL TIMEAVE_CUMULATE( layers_VFlux_T, layers_VFlux, Nlayers, |
275 |
& deltaTclock, bi, bj, myThid ) |
276 |
#ifdef LAYERS_THICKNESS |
277 |
CALL TIMEAVE_CUMULATE( layers_HV_T, layers_HV, Nlayers, |
278 |
& deltaTclock, bi, bj, myThid ) |
279 |
#endif /* LAYERS_THICKNESS */ |
280 |
#endif /* LAYERS_VFLUX */ |
281 |
|
282 |
layers_TimeAve(bi,bj)=layers_TimeAve(bi,bj)+deltaTclock |
283 |
|
284 |
ENDIF |
285 |
#endif /* ALLOW_TIMEAVE */ |
286 |
|
287 |
C --- End bi,bj loop |
288 |
ENDDO |
289 |
ENDDO |
290 |
|
291 |
#endif /* ALLOW_LAYERS */ |
292 |
|
293 |
RETURN |
294 |
END |