/[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.28 - (show annotations) (download)
Tue Mar 6 17:02:57 2001 UTC (23 years, 2 months ago) by jmc
Branch: MAIN
Changes since 1.27: +14 -14 lines
change local variable names etaN & etaNM1 (now global state variables)

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

  ViewVC Help
Powered by ViewVC 1.1.22