/[MITgcm]/MITgcm/model/src/cg2d.F
ViewVC logotype

Contents of /MITgcm/model/src/cg2d.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.1 - (show annotations) (download)
Wed Apr 22 19:15:30 1998 UTC (26 years, 1 month ago) by cnh
Branch: MAIN
Branch point for: cnh
Initial revision

1 C $Id$
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 "CG2D.h"
35
36 C === Routine arguments ===
37 C myThid - Thread on which I am working.
38 INTEGER myThid
39
40 C === Local variables ====
41 C actualIts - Number of iterations taken
42 C actualResidual - residual
43 C bi - Block index in X and Y.
44 C bj
45 C etaN - Used in computing search directions
46 C etaNM1 suffix N and NM1 denote current and
47 C cgBeta previous iterations respectively.
48 C alpha
49 C sumRHS - Sum of right-hand-side. Sometimes this is a
50 C useful debuggin/trouble shooting diagnostic.
51 C For neumann problems sumRHS needs to be ~0.
52 C or they converge at a non-zero residual.
53 C err - Measure of residual of Ax - b, usually the norm.
54 C I, J, N - Loop counters ( N counts CG iterations )
55 INTEGER actualIts
56 REAL actualResidual
57 INTEGER bi, bj
58 INTEGER I, J, it2d
59 REAL err
60 REAL etaN
61 REAL etaNM1
62 REAL cgBeta
63 REAL alpha
64 REAL sumRHS
65 REAL rhsMax
66 REAL rhsNorm
67
68 C-- Initialise inverter
69 etaNBuf(1,myThid) = 0. D0
70 errBuf(1,myThid) = 0. D0
71 sumRHSBuf(1,myThid) = 0. D0
72 etaNM1 = 1. D0
73
74 C-- Normalise RHS
75 rhsMax = 0. _d 0
76 DO bj=myByLo(myThid),myByHi(myThid)
77 DO bi=myBxLo(myThid),myBxHi(myThid)
78 DO J=1,sNy
79 DO I=1,sNx
80 cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*cg2dNorm
81 rhsMax = MAX(ABS(cg2d_b(I,J,bi,bj)),rhsMax)
82 ENDDO
83 ENDDO
84 ENDDO
85 ENDDO
86 rhsMaxBuf(1,myThid) = rhsMax
87 _GLOBAL_MAX_R8( rhsMaxbuf, rhsMax, myThid )
88 rhsMax = rhsMaxBuf(1,1)
89 rhsNorm = 1. _d 0
90 IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax
91 DO bj=myByLo(myThid),myByHi(myThid)
92 DO bi=myBxLo(myThid),myBxHi(myThid)
93 DO J=1,sNy
94 DO I=1,sNx
95 cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*rhsNorm
96 cg2d_x(I,J,bi,bj) = cg2d_x(I,J,bi,bj)*rhsNorm
97 ENDDO
98 ENDDO
99 ENDDO
100 ENDDO
101
102 C-- Update overlaps
103 _EXCH_XY_R8( cg2d_b, myThid )
104 _EXCH_XY_R8( cg2d_x, myThid )
105 CcnhDebugStarts
106 C CALL PLOT_FIELD_XYR8( cg2d_b, 'CG2D.1 CG2D_B' , 1, myThid )
107 CcnhDebugEnds
108
109 C-- Initial residual calculation
110 err = 0. _d 0
111 sumRHS = 0. _d 0
112 DO bj=myByLo(myThid),myByHi(myThid)
113 DO bi=myBxLo(myThid),myBxHi(myThid)
114 DO J=1,sNy
115 DO I=1,sNx
116 cg2d_s(I,J,bi,bj) = 0.
117 cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
118 & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
119 & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
120 & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
121 & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
122 & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
123 & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
124 & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
125 & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj))
126 err = err +
127 & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
128 sumRHS = sumRHS +
129 & cg2d_b(I,J,bi,bj)
130 ENDDO
131 ENDDO
132 ENDDO
133 ENDDO
134 _EXCH_XY_R8( cg2d_r, myThid )
135 _EXCH_XY_R8( cg2d_s, myThid )
136 sumRHSBuf(1,myThid) = sumRHS
137 _GLOBAL_SUM_R8( sumRHSBuf , sumRHS, myThid )
138 sumRHS = sumRHSBuf(1,1)
139 errBuf(1,myThid) = err
140 C WRITE(6,*) ' mythid, err = ', mythid, SQRT(err)
141 _GLOBAL_SUM_R8( errBuf , err , myThid )
142 err = errBuf(1,1)
143 C write(0,*) 'cg2d: Sum(rhs) = ',sumRHS
144
145 actualIts = 0
146 actualResidual = SQRT(err)
147 C _BARRIER
148 _BEGIN_MASTER( myThid )
149 WRITE(0,*) ' CG2D iters, err = ', actualIts, actualResidual
150 _END_MASTER( )
151
152 C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
153 DO 10 it2d=1, cg2dMaxIters
154
155 CcnhDebugStarts
156 C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' residual = ',actualResidual
157 CcnhDebugEnds
158 IF ( err .LT. cg2dTargetResidual ) GOTO 11
159 C-- Solve preconditioning equation and update
160 C-- conjugate direction vector "s".
161 etaN = 0. _d 0
162 DO bj=myByLo(myThid),myByHi(myThid)
163 DO bi=myBxLo(myThid),myBxHi(myThid)
164 DO J=1,sNy
165 DO I=1,sNx
166 cg2d_q(I,J,bi,bj) =
167 C & pW(I ,J ,bi,bj)*cg2d_r(I-1,J ,bi,bj)
168 C & +pW(I+1,J ,bi,bj)*cg2d_r(I+1,J ,bi,bj)
169 C & +pS(I ,J ,bi,bj)*cg2d_r(I ,J-1,bi,bj)
170 C & +pS(I ,J+1,bi,bj)*cg2d_r(I ,J+1,bi,bj)
171 C & +pC(I ,J ,bi,bj)*cg2d_r(I ,J ,bi,bj)
172 & cg2d_r(I ,J ,bi,bj)
173 etaN = etaN
174 & +cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
175 ENDDO
176 ENDDO
177 ENDDO
178 ENDDO
179
180 etanBuf(1,myThid) = etaN
181 _GLOBAL_SUM_R8(etaNbuf,etaN, myThid)
182 etaN = etaNBuf(1,1)
183 CcnhDebugStarts
184 C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' etaN = ',etaN
185 CcnhDebugEnds
186 cgBeta = etaN/etaNM1
187 CcnhDebugStarts
188 C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' beta = ',cgBeta
189 CcnhDebugEnds
190 etaNM1 = etaN
191
192 DO bj=myByLo(myThid),myByHi(myThid)
193 DO bi=myBxLo(myThid),myBxHi(myThid)
194 DO J=1,sNy
195 DO I=1,sNx
196 cg2d_s(I,J,bi,bj) = cg2d_q(I,J,bi,bj) + cgBeta*cg2d_s(I,J,bi,bj)
197 ENDDO
198 ENDDO
199 ENDDO
200 ENDDO
201
202 C-- Do exchanges that require messages i.e. between
203 C-- processes.
204 _EXCH_XY_R8( cg2d_s, myThid )
205
206 C== Evaluate laplace operator on conjugate gradient vector
207 C== q = A.s
208 alpha = 0. _d 0
209 DO bj=myByLo(myThid),myByHi(myThid)
210 DO bi=myBxLo(myThid),myBxHi(myThid)
211 DO J=1,sNy
212 DO I=1,sNx
213 cg2d_q(I,J,bi,bj) =
214 & aW2d(I ,J ,bi,bj)*cg2d_s(I-1,J ,bi,bj)
215 & +aW2d(I+1,J ,bi,bj)*cg2d_s(I+1,J ,bi,bj)
216 & +aS2d(I ,J ,bi,bj)*cg2d_s(I ,J-1,bi,bj)
217 & +aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J+1,bi,bj)
218 & -aW2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
219 & -aW2d(I+1,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
220 & -aS2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj)
221 & -aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J ,bi,bj)
222 alpha = alpha+cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj)
223 ENDDO
224 ENDDO
225 ENDDO
226 ENDDO
227 alphaBuf(1,myThid) = alpha
228 _GLOBAL_SUM_R8(alphaBuf,alpha,myThid)
229 alpha = alphaBuf(1,1)
230 CcnhDebugStarts
231 C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' SUM(s*q)= ',alpha
232 CcnhDebugEnds
233 alpha = etaN/alpha
234 CcnhDebugStarts
235 C WRITE(0,*) ' CG2D: Iteration ',it2d-1,' alpha= ',alpha
236 CcnhDebugEnds
237
238 C== Update solution and residual vectors
239 C Now compute "interior" points.
240 err = 0. _d 0
241 DO bj=myByLo(myThid),myByHi(myThid)
242 DO bi=myBxLo(myThid),myBxHi(myThid)
243 DO J=1,sNy
244 DO I=1,sNx
245 cg2d_x(I,J,bi,bj)=cg2d_x(I,J,bi,bj)+alpha*cg2d_s(I,J,bi,bj)
246 cg2d_r(I,J,bi,bj)=cg2d_r(I,J,bi,bj)-alpha*cg2d_q(I,J,bi,bj)
247 err = err+cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
248 ENDDO
249 ENDDO
250 ENDDO
251 ENDDO
252
253 errBuf(1,myThid) = err
254 _GLOBAL_SUM_R8( errBuf , err , myThid )
255 err = errBuf(1,1)
256 err = SQRT(err)
257 actualIts = it2d
258 actualResidual = err
259 IF ( err .LT. cg2dTargetResidual ) GOTO 11
260 _EXCH_XY_R8(cg2d_r, myThid )
261 10 CONTINUE
262 11 CONTINUE
263
264 C-- Un-normalise the answer
265 DO bj=myByLo(myThid),myByHi(myThid)
266 DO bi=myBxLo(myThid),myBxHi(myThid)
267 DO J=1,sNy
268 DO I=1,sNx
269 cg2d_x(I ,J ,bi,bj) = cg2d_x(I ,J ,bi,bj)/rhsNorm
270 ENDDO
271 ENDDO
272 ENDDO
273 ENDDO
274
275 _EXCH_XY_R8(cg2d_x, myThid )
276 _BEGIN_MASTER( myThid )
277 WRITE(6,*) ' CG2D iters, err = ', actualIts, actualResidual
278 _END_MASTER( )
279
280 CcnhDebugStarts
281 C CALL PLOT_FIELD_XYR8( cg2d_x, 'CALC_MOM_RHS CG2D_X' , 1, myThid )
282 C err = 0. _d 0
283 C DO bj=myByLo(myThid),myByHi(myThid)
284 C DO bi=myBxLo(myThid),myBxHi(myThid)
285 C DO J=1,sNy
286 C DO I=1,sNx
287 C cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) -
288 C & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj)
289 C & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj)
290 C & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj)
291 C & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj)
292 C & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
293 C & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
294 C & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj)
295 C & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj))
296 C err = err +
297 C & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj)
298 C ENDDO
299 C ENDDO
300 C ENDDO
301 C ENDDO
302 C errBuf(1,myThid) = err
303 C _GLOBAL_SUM_R8( errBuf , err , myThid )
304 C err = errBuf(1,1)
305 C write(0,*) 'cg2d: Ax - b = ',SQRT(err)
306 CcnhDebugEnds
307
308
309 END

  ViewVC Help
Powered by ViewVC 1.1.22