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