1 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/cg2d.F,v 1.10 1998/09/06 17:35:19 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 |
C-- Initialise inverter |
70 |
etaNBuf(1,myThid) = 0. D0 |
71 |
errBuf(1,myThid) = 0. D0 |
72 |
sumRHSBuf(1,myThid) = 0. D0 |
73 |
etaNM1 = 1. D0 |
74 |
|
75 |
CcnhDebugStarts |
76 |
C _EXCH_XY_R8( cg2d_b, myThid ) |
77 |
C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.0 CG2D_B' , 1, myThid ) |
78 |
CcnhDebugEnds |
79 |
|
80 |
C-- Normalise RHS |
81 |
rhsMax = 0. _d 0 |
82 |
DO bj=myByLo(myThid),myByHi(myThid) |
83 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
84 |
DO J=1,sNy |
85 |
DO I=1,sNx |
86 |
cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*cg2dNorm |
87 |
rhsMax = MAX(ABS(cg2d_b(I,J,bi,bj)),rhsMax) |
88 |
ENDDO |
89 |
ENDDO |
90 |
ENDDO |
91 |
ENDDO |
92 |
rhsMaxBuf(1,myThid) = rhsMax |
93 |
_GLOBAL_MAX_R8( rhsMaxbuf, rhsMax, myThid ) |
94 |
rhsMax = rhsMaxBuf(1,1) |
95 |
rhsNorm = 1. _d 0 |
96 |
IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax |
97 |
DO bj=myByLo(myThid),myByHi(myThid) |
98 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
99 |
DO J=1,sNy |
100 |
DO I=1,sNx |
101 |
cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*rhsNorm |
102 |
cg2d_x(I,J,bi,bj) = cg2d_x(I,J,bi,bj)*rhsNorm |
103 |
ENDDO |
104 |
ENDDO |
105 |
ENDDO |
106 |
ENDDO |
107 |
|
108 |
C-- Update overlaps |
109 |
_EXCH_XY_R8( cg2d_b, myThid ) |
110 |
_EXCH_XY_R8( cg2d_x, myThid ) |
111 |
CcnhDebugStarts |
112 |
C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.1 CG2D_B' , 1, myThid ) |
113 |
CcnhDebugEnds |
114 |
|
115 |
C-- Initial residual calculation |
116 |
err = 0. _d 0 |
117 |
sumRHS = 0. _d 0 |
118 |
DO bj=myByLo(myThid),myByHi(myThid) |
119 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
120 |
DO J=1,sNy |
121 |
DO I=1,sNx |
122 |
cg2d_s(I,J,bi,bj) = 0. |
123 |
cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) - |
124 |
& (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj) |
125 |
& +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj) |
126 |
& +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj) |
127 |
& +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj) |
128 |
& -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
129 |
& -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
130 |
& -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
131 |
& -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj) |
132 |
& -freeSurfFac*_rA(i,j,bi,bj)* horiVertRatio* |
133 |
& cg2d_x(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm |
134 |
& ) |
135 |
err = err + |
136 |
& cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj) |
137 |
sumRHS = sumRHS + |
138 |
& cg2d_b(I,J,bi,bj) |
139 |
ENDDO |
140 |
ENDDO |
141 |
ENDDO |
142 |
ENDDO |
143 |
_EXCH_XY_R8( cg2d_r, myThid ) |
144 |
_EXCH_XY_R8( cg2d_s, myThid ) |
145 |
sumRHSBuf(1,myThid) = sumRHS |
146 |
_GLOBAL_SUM_R8( sumRHSBuf , sumRHS, myThid ) |
147 |
sumRHS = sumRHSBuf(1,1) |
148 |
errBuf(1,myThid) = err |
149 |
C WRITE(6,*) ' mythid, err = ', mythid, SQRT(err) |
150 |
_GLOBAL_SUM_R8( errBuf , err , myThid ) |
151 |
err = errBuf(1,1) |
152 |
write(0,*) 'cg2d: Sum(rhs) = ',sumRHS |
153 |
|
154 |
actualIts = 0 |
155 |
actualResidual = SQRT(err) |
156 |
C _BARRIER |
157 |
_BEGIN_MASTER( myThid ) |
158 |
WRITE(0,*) ' CG2D iters, err = ', actualIts, actualResidual |
159 |
_END_MASTER( ) |
160 |
|
161 |
C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
162 |
DO 10 it2d=1, cg2dMaxIters |
163 |
|
164 |
CcnhDebugStarts |
165 |
C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' residual = ',actualResidual |
166 |
CcnhDebugEnds |
167 |
IF ( err .LT. cg2dTargetResidual ) GOTO 11 |
168 |
C-- Solve preconditioning equation and update |
169 |
C-- conjugate direction vector "s". |
170 |
etaN = 0. _d 0 |
171 |
DO bj=myByLo(myThid),myByHi(myThid) |
172 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
173 |
DO J=1,sNy |
174 |
DO I=1,sNx |
175 |
cg2d_q(I,J,bi,bj) = |
176 |
& pC(I ,J ,bi,bj)*cg2d_r(I ,J ,bi,bj) |
177 |
& +pW(I ,J ,bi,bj)*cg2d_r(I-1,J ,bi,bj) |
178 |
& +pW(I+1,J ,bi,bj)*cg2d_r(I+1,J ,bi,bj) |
179 |
& +pS(I ,J ,bi,bj)*cg2d_r(I ,J-1,bi,bj) |
180 |
& +pS(I ,J+1,bi,bj)*cg2d_r(I ,J+1,bi,bj) |
181 |
CcnhDebugStarts |
182 |
C cg2d_q(I,J,bi,bj) = cg2d_r(I ,J ,bi,bj) |
183 |
CcnhDebugEnds |
184 |
etaN = etaN |
185 |
& +cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj) |
186 |
ENDDO |
187 |
ENDDO |
188 |
ENDDO |
189 |
ENDDO |
190 |
|
191 |
etanBuf(1,myThid) = etaN |
192 |
_GLOBAL_SUM_R8(etaNbuf,etaN, myThid) |
193 |
etaN = etaNBuf(1,1) |
194 |
CcnhDebugStarts |
195 |
C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' etaN = ',etaN |
196 |
CcnhDebugEnds |
197 |
cgBeta = etaN/etaNM1 |
198 |
CcnhDebugStarts |
199 |
C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' beta = ',cgBeta |
200 |
CcnhDebugEnds |
201 |
etaNM1 = etaN |
202 |
|
203 |
DO bj=myByLo(myThid),myByHi(myThid) |
204 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
205 |
DO J=1,sNy |
206 |
DO I=1,sNx |
207 |
cg2d_s(I,J,bi,bj) = cg2d_q(I,J,bi,bj) + cgBeta*cg2d_s(I,J,bi,bj) |
208 |
ENDDO |
209 |
ENDDO |
210 |
ENDDO |
211 |
ENDDO |
212 |
|
213 |
C-- Do exchanges that require messages i.e. between |
214 |
C-- processes. |
215 |
_EXCH_XY_R8( cg2d_s, myThid ) |
216 |
|
217 |
C== Evaluate laplace operator on conjugate gradient vector |
218 |
C== q = A.s |
219 |
alpha = 0. _d 0 |
220 |
DO bj=myByLo(myThid),myByHi(myThid) |
221 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
222 |
DO J=1,sNy |
223 |
DO I=1,sNx |
224 |
cg2d_q(I,J,bi,bj) = |
225 |
& aW2d(I ,J ,bi,bj)*cg2d_s(I-1,J ,bi,bj) |
226 |
& +aW2d(I+1,J ,bi,bj)*cg2d_s(I+1,J ,bi,bj) |
227 |
& +aS2d(I ,J ,bi,bj)*cg2d_s(I ,J-1,bi,bj) |
228 |
& +aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J+1,bi,bj) |
229 |
& -aW2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj) |
230 |
& -aW2d(I+1,J ,bi,bj)*cg2d_s(I ,J ,bi,bj) |
231 |
& -aS2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj) |
232 |
& -aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J ,bi,bj) |
233 |
& -freeSurfFac*_rA(i,j,bi,bj)* horiVertRatio* |
234 |
& cg2d_s(I ,J ,bi,bj)/deltaTMom/deltaTMom*cg2dNorm |
235 |
alpha = alpha+cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj) |
236 |
ENDDO |
237 |
ENDDO |
238 |
ENDDO |
239 |
ENDDO |
240 |
alphaBuf(1,myThid) = alpha |
241 |
_GLOBAL_SUM_R8(alphaBuf,alpha,myThid) |
242 |
alpha = alphaBuf(1,1) |
243 |
CcnhDebugStarts |
244 |
C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' SUM(s*q)= ',alpha |
245 |
CcnhDebugEnds |
246 |
alpha = etaN/alpha |
247 |
CcnhDebugStarts |
248 |
C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' alpha= ',alpha |
249 |
CcnhDebugEnds |
250 |
|
251 |
C== Update solution and residual vectors |
252 |
C Now compute "interior" points. |
253 |
err = 0. _d 0 |
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_x(I,J,bi,bj)=cg2d_x(I,J,bi,bj)+alpha*cg2d_s(I,J,bi,bj) |
259 |
cg2d_r(I,J,bi,bj)=cg2d_r(I,J,bi,bj)-alpha*cg2d_q(I,J,bi,bj) |
260 |
err = err+cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj) |
261 |
ENDDO |
262 |
ENDDO |
263 |
ENDDO |
264 |
ENDDO |
265 |
|
266 |
errBuf(1,myThid) = err |
267 |
_GLOBAL_SUM_R8( errBuf , err , myThid ) |
268 |
err = errBuf(1,1) |
269 |
err = SQRT(err) |
270 |
actualIts = it2d |
271 |
actualResidual = err |
272 |
IF ( err .LT. cg2dTargetResidual ) GOTO 11 |
273 |
_EXCH_XY_R8(cg2d_r, myThid ) |
274 |
10 CONTINUE |
275 |
11 CONTINUE |
276 |
|
277 |
C-- Un-normalise the answer |
278 |
DO bj=myByLo(myThid),myByHi(myThid) |
279 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
280 |
DO J=1,sNy |
281 |
DO I=1,sNx |
282 |
cg2d_x(I ,J ,bi,bj) = cg2d_x(I ,J ,bi,bj)/rhsNorm |
283 |
ENDDO |
284 |
ENDDO |
285 |
ENDDO |
286 |
ENDDO |
287 |
|
288 |
_EXCH_XY_R8(cg2d_x, myThid ) |
289 |
_BEGIN_MASTER( myThid ) |
290 |
WRITE(0,*) ' CG2D iters, err = ', actualIts, actualResidual |
291 |
_END_MASTER( ) |
292 |
|
293 |
CcnhDebugStarts |
294 |
C CALL PLOT_FIELD_XYRL( cg2d_x, 'CALC_MOM_RHS CG2D_X' , 1, myThid ) |
295 |
C err = 0. _d 0 |
296 |
C DO bj=myByLo(myThid),myByHi(myThid) |
297 |
C DO bi=myBxLo(myThid),myBxHi(myThid) |
298 |
C DO J=1,sNy |
299 |
C DO I=1,sNx |
300 |
C cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) - |
301 |
C & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj) |
302 |
C & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj) |
303 |
C & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj) |
304 |
C & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj) |
305 |
C & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
306 |
C & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
307 |
C & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
308 |
C & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj)) |
309 |
C err = err + |
310 |
C & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj) |
311 |
C ENDDO |
312 |
C ENDDO |
313 |
C ENDDO |
314 |
C ENDDO |
315 |
C errBuf(1,myThid) = err |
316 |
C _GLOBAL_SUM_R8( errBuf , err , myThid ) |
317 |
C err = errBuf(1,1) |
318 |
C write(0,*) 'cg2d: Ax - b = ',SQRT(err) |
319 |
CcnhDebugEnds |
320 |
|
321 |
|
322 |
END |