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