1 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/cg2d.F,v 1.14 1998/10/28 03:11:36 cnh Exp $ |
2 |
|
3 |
#include "CPP_EEOPTIONS.h" |
4 |
|
5 |
SUBROUTINE CG2D( |
6 |
I cg2d_b, |
7 |
U cg2d_x, |
8 |
I myThid ) |
9 |
C /==========================================================\ |
10 |
C | SUBROUTINE CG2D | |
11 |
C | o Two-dimensional grid problem conjugate-gradient | |
12 |
C | inverter (with preconditioner). | |
13 |
C |==========================================================| |
14 |
C | Con. grad is an iterative procedure for solving Ax = b. | |
15 |
C | It requires the A be symmetric. | |
16 |
C | This implementation assumes A is a five-diagonal | |
17 |
C | matrix of the form that arises in the discrete | |
18 |
C | representation of the del^2 operator in a | |
19 |
C | two-dimensional space. | |
20 |
C | Notes: | |
21 |
C | ====== | |
22 |
C | This implementation can support shared-memory | |
23 |
C | multi-threaded execution. In order to do this COMMON | |
24 |
C | blocks are used for many of the arrays - even ones that | |
25 |
C | are only used for intermedaite results. This design is | |
26 |
C | OK if you want to all the threads to collaborate on | |
27 |
C | solving the same problem. On the other hand if you want | |
28 |
C | the threads to solve several different problems | |
29 |
C | concurrently this implementation will not work. | |
30 |
C \==========================================================/ |
31 |
|
32 |
C === Global data === |
33 |
#include "SIZE.h" |
34 |
#include "EEPARAMS.h" |
35 |
#include "PARAMS.h" |
36 |
#include "GRID.h" |
37 |
#include "CG2D_INTERNAL.h" |
38 |
|
39 |
C === Routine arguments === |
40 |
C myThid - Thread on which I am working. |
41 |
INTEGER myThid |
42 |
_RL cg2d_x(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
43 |
_RL cg2d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
44 |
|
45 |
|
46 |
C === Local variables ==== |
47 |
C actualIts - Number of iterations taken |
48 |
C actualResidual - residual |
49 |
C bi - Block index in X and Y. |
50 |
C bj |
51 |
C etaN - Used in computing search directions |
52 |
C etaNM1 suffix N and NM1 denote current and |
53 |
C cgBeta previous iterations respectively. |
54 |
C alpha |
55 |
C sumRHS - Sum of right-hand-side. Sometimes this is a |
56 |
C useful debuggin/trouble shooting diagnostic. |
57 |
C For neumann problems sumRHS needs to be ~0. |
58 |
C or they converge at a non-zero residual. |
59 |
C err - Measure of residual of Ax - b, usually the norm. |
60 |
C I, J, N - Loop counters ( N counts CG iterations ) |
61 |
INTEGER actualIts |
62 |
_RL actualResidual |
63 |
INTEGER bi, bj |
64 |
INTEGER I, J, it2d |
65 |
_RL err |
66 |
_RL etaN |
67 |
_RL etaNM1 |
68 |
_RL cgBeta |
69 |
_RL alpha |
70 |
_RL sumRHS |
71 |
_RL rhsMax |
72 |
_RL rhsNorm |
73 |
|
74 |
INTEGER OLw |
75 |
INTEGER OLe |
76 |
INTEGER OLn |
77 |
INTEGER OLs |
78 |
INTEGER exchWidthX |
79 |
INTEGER exchWidthY |
80 |
INTEGER myNz |
81 |
|
82 |
|
83 |
CcnhDebugStarts |
84 |
CHARACTER*(MAX_LEN_FNAM) suff |
85 |
CcnhDebugEnds |
86 |
|
87 |
|
88 |
C-- Initialise inverter |
89 |
etaNBuf(1,myThid) = 0. _d 0 |
90 |
errBuf(1,myThid) = 0. _d 0 |
91 |
sumRHSBuf(1,myThid) = 0. _d 0 |
92 |
etaNM1 = 1. _d 0 |
93 |
|
94 |
CcnhDebugStarts |
95 |
C _EXCH_XY_R8( cg2d_b, myThid ) |
96 |
C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.0 CG2D_B' , 1, myThid ) |
97 |
C suff = 'unnormalised' |
98 |
C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid) |
99 |
C STOP |
100 |
CcnhDebugEnds |
101 |
|
102 |
C-- Normalise RHS |
103 |
rhsMax = 0. _d 0 |
104 |
DO bj=myByLo(myThid),myByHi(myThid) |
105 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
106 |
DO J=1,sNy |
107 |
DO I=1,sNx |
108 |
cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*cg2dNorm |
109 |
rhsMax = MAX(ABS(cg2d_b(I,J,bi,bj)),rhsMax) |
110 |
ENDDO |
111 |
ENDDO |
112 |
ENDDO |
113 |
ENDDO |
114 |
rhsMaxBuf(1,myThid) = rhsMax |
115 |
_GLOBAL_MAX_R8( rhsMaxbuf, rhsMax, myThid ) |
116 |
rhsMax = rhsMaxBuf(1,1) |
117 |
rhsNorm = 1. _d 0 |
118 |
IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax |
119 |
DO bj=myByLo(myThid),myByHi(myThid) |
120 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
121 |
DO J=1,sNy |
122 |
DO I=1,sNx |
123 |
cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*rhsNorm |
124 |
cg2d_x(I,J,bi,bj) = cg2d_x(I,J,bi,bj)*rhsNorm |
125 |
ENDDO |
126 |
ENDDO |
127 |
ENDDO |
128 |
ENDDO |
129 |
|
130 |
C-- Update overlaps |
131 |
_EXCH_XY_R8( cg2d_b, myThid ) |
132 |
_EXCH_XY_R8( cg2d_x, myThid ) |
133 |
CcnhDebugStarts |
134 |
C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.1 CG2D_B' , 1, myThid ) |
135 |
C suff = 'normalised' |
136 |
C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid) |
137 |
CcnhDebugEnds |
138 |
|
139 |
C-- Initial residual calculation |
140 |
err = 0. _d 0 |
141 |
sumRHS = 0. _d 0 |
142 |
DO bj=myByLo(myThid),myByHi(myThid) |
143 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
144 |
DO J=1,sNy |
145 |
DO I=1,sNx |
146 |
cg2d_s(I,J,bi,bj) = 0. |
147 |
cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) - |
148 |
& (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj) |
149 |
& +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj) |
150 |
& +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj) |
151 |
& +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj) |
152 |
& -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
153 |
& -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
154 |
& -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
155 |
& -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj) |
156 |
& -freeSurfFac*_rA(i,j,bi,bj)* horiVertRatio* |
157 |
& cg2d_x(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm |
158 |
& ) |
159 |
err = err + |
160 |
& cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj) |
161 |
sumRHS = sumRHS + |
162 |
& cg2d_b(I,J,bi,bj) |
163 |
ENDDO |
164 |
ENDDO |
165 |
ENDDO |
166 |
ENDDO |
167 |
C _EXCH_XY_R8( cg2d_r, myThid ) |
168 |
OLw = 1 |
169 |
OLe = 1 |
170 |
OLn = 1 |
171 |
OLs = 1 |
172 |
exchWidthX = 1 |
173 |
exchWidthY = 1 |
174 |
myNz = 1 |
175 |
CALL EXCH_RL( cg2d_r, |
176 |
I OLw, OLe, OLs, OLn, myNz, |
177 |
I exchWidthX, exchWidthY, |
178 |
I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid ) |
179 |
C _EXCH_XY_R8( cg2d_s, myThid ) |
180 |
OLw = 1 |
181 |
OLe = 1 |
182 |
OLn = 1 |
183 |
OLs = 1 |
184 |
exchWidthX = 1 |
185 |
exchWidthY = 1 |
186 |
myNz = 1 |
187 |
CALL EXCH_RL( cg2d_s, |
188 |
I OLw, OLe, OLs, OLn, myNz, |
189 |
I exchWidthX, exchWidthY, |
190 |
I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid ) |
191 |
sumRHSBuf(1,myThid) = sumRHS |
192 |
_GLOBAL_SUM_R8( sumRHSBuf , sumRHS, myThid ) |
193 |
sumRHS = sumRHSBuf(1,1) |
194 |
errBuf(1,myThid) = err |
195 |
C WRITE(6,*) ' mythid, err = ', mythid, SQRT(err) |
196 |
_GLOBAL_SUM_R8( errBuf , err , myThid ) |
197 |
err = errBuf(1,1) |
198 |
|
199 |
_BEGIN_MASTER( myThid ) |
200 |
write(0,'(A,1PE30.20)') 'cg2d: Sum(rhs) = ',sumRHS |
201 |
_END_MASTER( ) |
202 |
|
203 |
actualIts = 0 |
204 |
actualResidual = SQRT(err) |
205 |
C _BARRIER |
206 |
_BEGIN_MASTER( myThid ) |
207 |
WRITE(0,'(A,I6,1PE30.20)') ' CG2D iters, err = ', |
208 |
& actualIts, actualResidual |
209 |
_END_MASTER( ) |
210 |
|
211 |
C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
212 |
DO 10 it2d=1, cg2dMaxIters |
213 |
|
214 |
CcnhDebugStarts |
215 |
C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' residual = ', |
216 |
C & actualResidual |
217 |
CcnhDebugEnds |
218 |
IF ( err .LT. cg2dTargetResidual ) GOTO 11 |
219 |
C-- Solve preconditioning equation and update |
220 |
C-- conjugate direction vector "s". |
221 |
etaN = 0. _d 0 |
222 |
DO bj=myByLo(myThid),myByHi(myThid) |
223 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
224 |
DO J=1,sNy |
225 |
DO I=1,sNx |
226 |
cg2d_q(I,J,bi,bj) = |
227 |
& pC(I ,J ,bi,bj)*cg2d_r(I ,J ,bi,bj) |
228 |
& +pW(I ,J ,bi,bj)*cg2d_r(I-1,J ,bi,bj) |
229 |
& +pW(I+1,J ,bi,bj)*cg2d_r(I+1,J ,bi,bj) |
230 |
& +pS(I ,J ,bi,bj)*cg2d_r(I ,J-1,bi,bj) |
231 |
& +pS(I ,J+1,bi,bj)*cg2d_r(I ,J+1,bi,bj) |
232 |
CcnhDebugStarts |
233 |
C cg2d_q(I,J,bi,bj) = cg2d_r(I ,J ,bi,bj) |
234 |
CcnhDebugEnds |
235 |
etaN = etaN |
236 |
& +cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj) |
237 |
ENDDO |
238 |
ENDDO |
239 |
ENDDO |
240 |
ENDDO |
241 |
|
242 |
etanBuf(1,myThid) = etaN |
243 |
_GLOBAL_SUM_R8(etaNbuf,etaN, myThid) |
244 |
etaN = etaNBuf(1,1) |
245 |
CcnhDebugStarts |
246 |
C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' etaN = ',etaN |
247 |
CcnhDebugEnds |
248 |
cgBeta = etaN/etaNM1 |
249 |
CcnhDebugStarts |
250 |
C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' beta = ',cgBeta |
251 |
CcnhDebugEnds |
252 |
etaNM1 = etaN |
253 |
|
254 |
DO bj=myByLo(myThid),myByHi(myThid) |
255 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
256 |
DO J=1,sNy |
257 |
DO I=1,sNx |
258 |
cg2d_s(I,J,bi,bj) = cg2d_q(I,J,bi,bj) |
259 |
& + cgBeta*cg2d_s(I,J,bi,bj) |
260 |
ENDDO |
261 |
ENDDO |
262 |
ENDDO |
263 |
ENDDO |
264 |
|
265 |
C-- Do exchanges that require messages i.e. between |
266 |
C-- processes. |
267 |
C _EXCH_XY_R8( cg2d_s, myThid ) |
268 |
OLw = 1 |
269 |
OLe = 1 |
270 |
OLn = 1 |
271 |
OLs = 1 |
272 |
exchWidthX = 1 |
273 |
exchWidthY = 1 |
274 |
myNz = 1 |
275 |
CALL EXCH_RL( cg2d_s, |
276 |
I OLw, OLe, OLs, OLn, myNz, |
277 |
I exchWidthX, exchWidthY, |
278 |
I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid ) |
279 |
|
280 |
|
281 |
C== Evaluate laplace operator on conjugate gradient vector |
282 |
C== q = A.s |
283 |
alpha = 0. _d 0 |
284 |
DO bj=myByLo(myThid),myByHi(myThid) |
285 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
286 |
DO J=1,sNy |
287 |
DO I=1,sNx |
288 |
cg2d_q(I,J,bi,bj) = |
289 |
& aW2d(I ,J ,bi,bj)*cg2d_s(I-1,J ,bi,bj) |
290 |
& +aW2d(I+1,J ,bi,bj)*cg2d_s(I+1,J ,bi,bj) |
291 |
& +aS2d(I ,J ,bi,bj)*cg2d_s(I ,J-1,bi,bj) |
292 |
& +aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J+1,bi,bj) |
293 |
& -aW2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj) |
294 |
& -aW2d(I+1,J ,bi,bj)*cg2d_s(I ,J ,bi,bj) |
295 |
& -aS2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj) |
296 |
& -aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J ,bi,bj) |
297 |
& -freeSurfFac*_rA(i,j,bi,bj)* horiVertRatio* |
298 |
& cg2d_s(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm |
299 |
alpha = alpha+cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj) |
300 |
ENDDO |
301 |
ENDDO |
302 |
ENDDO |
303 |
ENDDO |
304 |
alphaBuf(1,myThid) = alpha |
305 |
_GLOBAL_SUM_R8(alphaBuf,alpha,myThid) |
306 |
alpha = alphaBuf(1,1) |
307 |
CcnhDebugStarts |
308 |
C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' SUM(s*q)= ',alpha |
309 |
CcnhDebugEnds |
310 |
alpha = etaN/alpha |
311 |
CcnhDebugStarts |
312 |
C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' alpha= ',alpha |
313 |
CcnhDebugEnds |
314 |
|
315 |
C== Update solution and residual vectors |
316 |
C Now compute "interior" points. |
317 |
err = 0. _d 0 |
318 |
DO bj=myByLo(myThid),myByHi(myThid) |
319 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
320 |
DO J=1,sNy |
321 |
DO I=1,sNx |
322 |
cg2d_x(I,J,bi,bj)=cg2d_x(I,J,bi,bj)+alpha*cg2d_s(I,J,bi,bj) |
323 |
cg2d_r(I,J,bi,bj)=cg2d_r(I,J,bi,bj)-alpha*cg2d_q(I,J,bi,bj) |
324 |
err = err+cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj) |
325 |
ENDDO |
326 |
ENDDO |
327 |
ENDDO |
328 |
ENDDO |
329 |
|
330 |
errBuf(1,myThid) = err |
331 |
_GLOBAL_SUM_R8( errBuf , err , myThid ) |
332 |
err = errBuf(1,1) |
333 |
err = SQRT(err) |
334 |
actualIts = it2d |
335 |
actualResidual = err |
336 |
IF ( err .LT. cg2dTargetResidual ) GOTO 11 |
337 |
C _EXCH_XY_R8(cg2d_r, myThid ) |
338 |
OLw = 1 |
339 |
OLe = 1 |
340 |
OLn = 1 |
341 |
OLs = 1 |
342 |
exchWidthX = 1 |
343 |
exchWidthY = 1 |
344 |
myNz = 1 |
345 |
CALL EXCH_RL( cg2d_r, |
346 |
I OLw, OLe, OLs, OLn, myNz, |
347 |
I exchWidthX, exchWidthY, |
348 |
I FORWARD_SIMULATION, EXCH_IGNORE_CORNERS, myThid ) |
349 |
|
350 |
10 CONTINUE |
351 |
11 CONTINUE |
352 |
|
353 |
C-- Un-normalise the answer |
354 |
DO bj=myByLo(myThid),myByHi(myThid) |
355 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
356 |
DO J=1,sNy |
357 |
DO I=1,sNx |
358 |
cg2d_x(I ,J ,bi,bj) = cg2d_x(I ,J ,bi,bj)/rhsNorm |
359 |
ENDDO |
360 |
ENDDO |
361 |
ENDDO |
362 |
ENDDO |
363 |
|
364 |
_EXCH_XY_R8(cg2d_x, myThid ) |
365 |
_BEGIN_MASTER( myThid ) |
366 |
WRITE(0,'(A,I6,1PE30.20)') ' CG2D iters, err = ', |
367 |
& actualIts, actualResidual |
368 |
_END_MASTER( ) |
369 |
|
370 |
CcnhDebugStarts |
371 |
C CALL PLOT_FIELD_XYRL( cg2d_x, 'CALC_MOM_RHS CG2D_X' , 1, myThid ) |
372 |
C err = 0. _d 0 |
373 |
C DO bj=myByLo(myThid),myByHi(myThid) |
374 |
C DO bi=myBxLo(myThid),myBxHi(myThid) |
375 |
C DO J=1,sNy |
376 |
C DO I=1,sNx |
377 |
C cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) - |
378 |
C & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj) |
379 |
C & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj) |
380 |
C & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj) |
381 |
C & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj) |
382 |
C & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
383 |
C & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
384 |
C & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
385 |
C & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj)) |
386 |
C err = err + |
387 |
C & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj) |
388 |
C ENDDO |
389 |
C ENDDO |
390 |
C ENDDO |
391 |
C ENDDO |
392 |
C errBuf(1,myThid) = err |
393 |
C _GLOBAL_SUM_R8( errBuf , err , myThid ) |
394 |
C err = errBuf(1,1) |
395 |
C write(0,*) 'cg2d: Ax - b = ',SQRT(err) |
396 |
CcnhDebugEnds |
397 |
|
398 |
|
399 |
END |