1 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/ini_cg3d.F,v 1.7.2.2 2001/01/30 21:02:59 adcroft Exp $ |
2 |
|
3 |
#include "CPP_OPTIONS.h" |
4 |
|
5 |
CStartOfInterface |
6 |
SUBROUTINE INI_CG3D( myThid ) |
7 |
C /==========================================================\ |
8 |
C | SUBROUTINE INI_CG3D | |
9 |
C | o Initialise 3d conjugate gradient solver operators. | |
10 |
C |==========================================================| |
11 |
C | These arrays are purely a function of the basin geom. | |
12 |
C | We set then here once and them use then repeatedly. | |
13 |
C \==========================================================/ |
14 |
IMPLICIT NONE |
15 |
|
16 |
C === Global variables === |
17 |
#include "SIZE.h" |
18 |
#include "EEPARAMS.h" |
19 |
#include "PARAMS.h" |
20 |
#include "GRID.h" |
21 |
#include "DYNVARS.h" |
22 |
#include "CG3D.h" |
23 |
#ifdef ALLOW_OBCS |
24 |
#include "OBCS.h" |
25 |
#endif |
26 |
|
27 |
C === Routine arguments === |
28 |
C myThid - Thread no. that called this routine. |
29 |
INTEGER myThid |
30 |
CEndOfInterface |
31 |
|
32 |
#ifdef ALLOW_NONHYDROSTATIC |
33 |
|
34 |
C === Local variables === |
35 |
C xG, yG - Global coordinate location. |
36 |
C zG |
37 |
C iG, jG - Global coordinate index |
38 |
C bi,bj - Loop counters |
39 |
C faceArea - Temporary used to hold cell face areas. |
40 |
C I,J,K |
41 |
C myNorm - Work variable used in clculating normalisation factor |
42 |
CHARACTER*(MAX_LEN_MBUF) msgBuf |
43 |
INTEGER bi, bj |
44 |
INTEGER I, J, K |
45 |
_RL faceArea |
46 |
_RS myNorm |
47 |
_RL aCw, aCs |
48 |
_RL theRecip_Dr |
49 |
_RL aU, aL, aW, aE, aN, aS, aC |
50 |
|
51 |
CcnhDebugStarts |
52 |
_RL phi(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
53 |
CcnhDebugEnds |
54 |
|
55 |
C-- Initialise laplace operator |
56 |
C aW3d: Ax/dX |
57 |
C aS3d: Ay/dY |
58 |
C aV3d: Ar/dR |
59 |
myNorm = 0. _d 0 |
60 |
DO bj=myByLo(myThid),myByHi(myThid) |
61 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
62 |
DO K=1,Nr |
63 |
DO J=1,sNy |
64 |
DO I=1,sNx |
65 |
faceArea = _dyG(I,J,bi,bj)*drF(K) |
66 |
& *_hFacW(I,J,K,bi,bj) |
67 |
aW3d(I,J,K,bi,bj) = Gravity*faceArea*recip_dxC(I,J,bi,bj) |
68 |
faceArea = _dxG(I,J,bi,bj)*drF(K) |
69 |
& *_hFacS(I,J,K,bi,bj) |
70 |
aS3d(I,J,K,bi,bj) = Gravity*faceArea*recip_dyC(I,J,bi,bj) |
71 |
myNorm = MAX(ABS(aW3d(I,J,K,bi,bj)),myNorm) |
72 |
myNorm = MAX(ABS(aS3d(I,J,K,bi,bj)),myNorm) |
73 |
ENDDO |
74 |
ENDDO |
75 |
ENDDO |
76 |
DO K=1,1 |
77 |
DO J=1,sNy |
78 |
DO I=1,sNx |
79 |
aV3d(I,J,K,bi,bj) = 0. |
80 |
myNorm = MAX(ABS(aV3d(I,J,K,bi,bj)),myNorm) |
81 |
ENDDO |
82 |
ENDDO |
83 |
ENDDO |
84 |
DO K=2,Nr |
85 |
DO J=1,sNy |
86 |
DO I=1,sNx |
87 |
IF ( _hFacC(i,j,k,bi,bj) .EQ. 0. ) THEN |
88 |
faceArea = 0. |
89 |
ELSE |
90 |
faceArea = _rA(I,J,bi,bj) |
91 |
ENDIF |
92 |
theRecip_Dr = |
93 |
& drC(K) |
94 |
caja & drF(K )*_hFacC(i,j,k ,bi,bj)*0.5 |
95 |
caja & +drF(K-1)*_hFacC(i,j,k-1,bi,bj)*0.5 |
96 |
IF ( theRecip_Dr .NE. 0. ) |
97 |
& theRecip_Dr = 1. _d 0/theRecip_Dr |
98 |
aV3d(I,J,K,bi,bj) = |
99 |
& Gravity*faceArea*theRecip_Dr |
100 |
myNorm = MAX(ABS(aV3d(I,J,K,bi,bj)),myNorm) |
101 |
ENDDO |
102 |
ENDDO |
103 |
ENDDO |
104 |
#ifdef ALLOW_OBCS |
105 |
DO K=1,Nr |
106 |
IF (useOBCS) THEN |
107 |
DO I=1,sNx |
108 |
IF (OB_Jn(I,bi,bj).NE.0) THEN |
109 |
aS3d(I,OB_Jn(I,bi,bj),K,bi,bj)=0. |
110 |
aS3d(I,OB_Jn(I,bi,bj)+1,K,bi,bj)=0. |
111 |
aW3d(I,OB_Jn(I,bi,bj),K,bi,bj)=0. |
112 |
aW3d(I+1,OB_Jn(I,bi,bj),K,bi,bj)=0. |
113 |
aV3d(I,OB_Jn(I,bi,bj),K,bi,bj)=0. |
114 |
ENDIF |
115 |
IF (OB_Js(I,bi,bj).NE.0) THEN |
116 |
aS3d(I,OB_Js(I,bi,bj)+1,K,bi,bj)=0. |
117 |
aS3d(I,OB_Js(I,bi,bj),K,bi,bj)=0. |
118 |
aW3d(I,OB_Js(I,bi,bj),K,bi,bj)=0. |
119 |
aW3d(I+1,OB_Js(I,bi,bj),K,bi,bj)=0. |
120 |
aV3d(I,OB_Js(I,bi,bj),K,bi,bj)=0. |
121 |
ENDIF |
122 |
ENDDO |
123 |
DO J=1,sNy |
124 |
IF (OB_Ie(J,bi,bj).NE.0) THEN |
125 |
aW3d(OB_Ie(J,bi,bj),J,K,bi,bj)=0. |
126 |
aW3d(OB_Ie(J,bi,bj)+1,J,K,bi,bj)=0. |
127 |
aS3d(OB_Ie(J,bi,bj),J,K,bi,bj)=0. |
128 |
aS3d(OB_Ie(J,bi,bj),J+1,K,bi,bj)=0. |
129 |
aV3d(OB_Ie(J,bi,bj),J,K,bi,bj)=0. |
130 |
ENDIF |
131 |
IF (OB_Iw(J,bi,bj).NE.0) THEN |
132 |
aW3d(OB_Iw(J,bi,bj)+1,J,K,bi,bj)=0. |
133 |
aW3d(OB_Iw(J,bi,bj),J,K,bi,bj)=0. |
134 |
aS3d(OB_Iw(J,bi,bj),J,K,bi,bj)=0. |
135 |
aS3d(OB_Iw(J,bi,bj),J+1,K,bi,bj)=0. |
136 |
aV3d(OB_Iw(J,bi,bj),J,K,bi,bj)=0. |
137 |
ENDIF |
138 |
ENDDO |
139 |
ENDIF |
140 |
ENDDO |
141 |
#endif |
142 |
ENDDO |
143 |
ENDDO |
144 |
_GLOBAL_MAX_R4( myNorm, myThid ) |
145 |
IF ( myNorm .NE. 0. _d 0 ) THEN |
146 |
myNorm = 1. _d 0/myNorm |
147 |
ELSE |
148 |
myNorm = 1. _d 0 |
149 |
ENDIF |
150 |
cg3dNorm = myNorm |
151 |
_BEGIN_MASTER( myThid ) |
152 |
CcnhDebugStarts |
153 |
WRITE(msgBuf,'(A,E40.25)') '// CG3D normalisation factor = ' |
154 |
& , cg3dNorm |
155 |
CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1) |
156 |
WRITE(msgBuf,*) ' ' |
157 |
CALL PRINT_MESSAGE( msgBuf,standardMessageUnit,SQUEEZE_RIGHT,1) |
158 |
CcnhDebugEnds |
159 |
_END_MASTER( myThid ) |
160 |
DO bj=myByLo(myThid),myByHi(myThid) |
161 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
162 |
DO K=1,Nr |
163 |
DO J=1,sNy |
164 |
DO I=1,sNx |
165 |
aW3d(I,J,K,bi,bj) = aW3d(I,J,K,bi,bj)*myNorm |
166 |
aS3d(I,J,K,bi,bj) = aS3d(I,J,K,bi,bj)*myNorm |
167 |
aV3d(I,J,K,bi,bj) = aV3d(I,J,K,bi,bj)*myNorm |
168 |
ENDDO |
169 |
ENDDO |
170 |
ENDDO |
171 |
ENDDO |
172 |
ENDDO |
173 |
|
174 |
C-- Update overlap regions |
175 |
_EXCH_XYZ_R4(aW3d, myThid) |
176 |
_EXCH_XYZ_R4(aS3d, myThid) |
177 |
_EXCH_XYZ_R4(aV3d, myThid) |
178 |
CcnhDebugStarts |
179 |
C CALL PLOT_FIELD_XYZRS( aW3d, 'AW3D INI_CG3D.1' , Nr, 1, myThid ) |
180 |
C CALL PLOT_FIELD_XYZRS( aS3d, 'AS3D INI_CG3D.1' , Nr, 1, myThid ) |
181 |
CcnhDebugEnds |
182 |
|
183 |
C-- Initialise preconditioner |
184 |
C For now PC is just the identity. Change to |
185 |
C be LU factorization of d2/dz2 later. Note |
186 |
C check for consistency with S/R CG3D before |
187 |
C assuming zML is lower and zMU is upper! |
188 |
DO bj=myByLo(myThid),myByHi(myThid) |
189 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
190 |
DO K=1,Nr |
191 |
DO J=1,sNy |
192 |
DO I=1,sNx |
193 |
aW = aW3d(I ,J,K,bi,bj) |
194 |
aE = aW3d(I+1,J,K,bi,bj) |
195 |
aN = aS3d(I,J+1,K,bi,bj) |
196 |
aS = aS3d(I,J ,K,bi,bj) |
197 |
IF ( K .NE. 1 ) THEN |
198 |
aU = aV3d(I,J,K,bi,bj) |
199 |
ELSE |
200 |
aU = 0 |
201 |
ENDIF |
202 |
IF ( K .NE. Nr ) THEN |
203 |
aL = aV3d(I,J,K+1,bi,bj) |
204 |
ELSE |
205 |
aL = 0 |
206 |
ENDIF |
207 |
aC = -aW-aE-aN-aS-aU-aL |
208 |
IF ( K .EQ. 1 .AND. aC .NE. 0. ) THEN |
209 |
aC = aC |
210 |
& -freeSurfFac*myNorm* horiVertRatio* |
211 |
& rA(I,J,bi,bj)/deltaTMom/deltaTMom |
212 |
ENDIF |
213 |
IF ( aC .NE. 0. ) THEN |
214 |
zMC(i,j,k,bi,bj) = aC |
215 |
zMU(i,j,k,bi,bj) = aL |
216 |
zML(i,j,k,bi,bj) = aU |
217 |
CcnhDebugStarts |
218 |
C zMC(i,j,k,bi,bj) = 1. |
219 |
C zMU(i,j,k,bi,bj) = 0. |
220 |
C zML(i,j,k,bi,bj) = 0. |
221 |
CcnhDebugEnds |
222 |
ELSE |
223 |
zMC(i,j,k,bi,bj) = 1. _d 0 |
224 |
zMU(i,j,k,bi,bj) = 0. |
225 |
zML(i,j,k,bi,bj) = 0. |
226 |
ENDIF |
227 |
ENDDO |
228 |
ENDDO |
229 |
ENDDO |
230 |
DO J=1,sNy |
231 |
DO I=1,sNx |
232 |
zMC(i,j,1,bi,bj)= |
233 |
& 1. _d 0 / zMC(i,j,1,bi,bj) |
234 |
zMU(i,j,1,bi,bj)= |
235 |
& zMU(i,j,1,bi,bj)*zMC(i,j,1,bi,bj) |
236 |
ENDDO |
237 |
ENDDO |
238 |
DO K=2,Nr |
239 |
DO J=1,sNy |
240 |
DO I=1,sNx |
241 |
zMC(i,j,k,bi,bj) = 1. _d 0 / |
242 |
& (zMC(i,j,k,bi,bj)-zML(i,j,k,bi,bj)*zMU(i,j,k-1,bi,bj)) |
243 |
zMU(i,j,k,bi,bj)=zMU(i,j,k,bi,bj)*zMC(i,j,k,bi,bj) |
244 |
ENDDO |
245 |
ENDDO |
246 |
ENDDO |
247 |
DO K=1,Nr |
248 |
DO J=1,sNy |
249 |
DO I=1,sNx |
250 |
aW = aW3d(I ,J,K,bi,bj) |
251 |
aE = aW3d(I+1,J,K,bi,bj) |
252 |
aN = aS3d(I,J+1,K,bi,bj) |
253 |
aS = aS3d(I,J ,K,bi,bj) |
254 |
IF ( K .NE. 1 ) THEN |
255 |
aU = aV3d(I,J,K-1,bi,bj) |
256 |
ELSE |
257 |
aU = 0 |
258 |
ENDIF |
259 |
IF ( K .NE. Nr ) THEN |
260 |
aL = aV3d(I,J,K+1,bi,bj) |
261 |
ELSE |
262 |
aL = 0 |
263 |
ENDIF |
264 |
aC = -aW-aE-aN-aS-aU-aL |
265 |
IF ( aC .EQ. 0. ) THEN |
266 |
zMC(i,j,k,bi,bj) = 1. |
267 |
zML(i,j,k,bi,bj) = 0. |
268 |
zMU(i,j,k,bi,bj) = 0. |
269 |
CcnhDebugStarts |
270 |
C ELSE |
271 |
C zMC(i,j,k,bi,bj) = 1. |
272 |
C zML(i,j,k,bi,bj) = 0. |
273 |
C zMU(i,j,k,bi,bj) = 0. |
274 |
CcnhDEbugEnds |
275 |
ENDIF |
276 |
ENDDO |
277 |
ENDDO |
278 |
ENDDO |
279 |
ENDDO |
280 |
ENDDO |
281 |
C-- Update overlap regions |
282 |
_EXCH_XYZ_R4(zMC, myThid) |
283 |
_EXCH_XYZ_R4(zML, myThid) |
284 |
_EXCH_XYZ_R4(zMU, myThid) |
285 |
|
286 |
CcnhDebugStarts |
287 |
DO k=1,nr |
288 |
DO j=1-OLy,sNy+OLy |
289 |
DO i=1-OLx,sNx+OLx |
290 |
phi(i,j,1,1) = zMc(i,j,k,1,1) |
291 |
ENDDO |
292 |
ENDDO |
293 |
C CALL PLOT_FIELD_XYRS( phi, 'zMC INI_CG3D.1' , 1, myThid ) |
294 |
ENDDO |
295 |
C CALL PLOT_FIELD_XYRS( zMU, 'zMU INI_CG3D.1' , Nr, 1, myThid ) |
296 |
C CALL PLOT_FIELD_XYRS( zML, 'zML INI_CG3D.1' , Nr, 1, myThid ) |
297 |
CcnhDebugEnds |
298 |
|
299 |
C-- Set default values for initial guess and RHS |
300 |
IF ( startTime .EQ. 0 ) THEN |
301 |
DO bj=myByLo(myThid),myByHi(myThid) |
302 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
303 |
DO K=1,Nr |
304 |
DO J=1-OLy,sNy+OLy |
305 |
DO I=1-OLx,sNx+OLx |
306 |
cg3d_x(I,J,K,bi,bj) = 0. _d 0 |
307 |
cg3d_b(I,J,K,bi,bj) = 0. _d 0 |
308 |
#ifdef INCLUDE_CD_CODE |
309 |
cg3d_xNM1(I,J,K,bi,bj) = 0. _d 0 |
310 |
#endif |
311 |
ENDDO |
312 |
ENDDO |
313 |
ENDDO |
314 |
ENDDO |
315 |
ENDDO |
316 |
C-- Update overlap regions |
317 |
_EXCH_XYZ_R8(cg3d_x, myThid) |
318 |
_EXCH_XYZ_R8(cg3d_b, myThid) |
319 |
#ifdef INCLUDE_CD_CODE |
320 |
_EXCH_XYZ_R8(cg3d_xNM1, myThid) |
321 |
#endif |
322 |
ENDIF |
323 |
|
324 |
#endif /* ALLOW_NONHYDROSTATIC */ |
325 |
|
326 |
RETURN |
327 |
END |