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