/[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.18 - (show annotations) (download)
Wed Dec 9 16:11:51 1998 UTC (25 years, 5 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint19
Changes since 1.17: +2 -1 lines
Added IMPLICIT NONE in a lot of subroutines.
Also corrected the recip_Rhonil bug: we didn't set it in ini_parms.F

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

  ViewVC Help
Powered by ViewVC 1.1.22