/[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.6 - (show annotations) (download)
Fri Jun 12 19:33:33 1998 UTC (25 years, 11 months ago) by cnh
Branch: MAIN
Changes since 1.5: +5 -5 lines
Chages to make default setup correct for 4 degreee global comparison

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

  ViewVC Help
Powered by ViewVC 1.1.22