/[MITgcm]/MITgcm/pkg/diagnostics/diag_cg2d.F
ViewVC logotype

Contents of /MITgcm/pkg/diagnostics/diag_cg2d.F

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


Revision 1.6 - (show annotations) (download)
Tue Nov 18 23:53:04 2014 UTC (9 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65g, HEAD
Changes since 1.5: +9 -8 lines
add preconditioner off-diagonal factor (diagCG_pcOffDFac) as run-time params

1 C $Header: /u/gcmpack/MITgcm/pkg/diagnostics/diag_cg2d.F,v 1.5 2012/03/19 18:22:09 jmc Exp $
2 C $Name: $
3
4 #include "DIAG_OPTIONS.h"
5 #undef DEBUG_DIAG_CG2D
6
7 CBOP
8 C !ROUTINE: DIAG_CG2D
9 C !INTERFACE:
10 SUBROUTINE DIAG_CG2D(
11 I aW2d, aS2d, b2d,
12 I offDiagFactor, residCriter,
13 O firstResidual, minResidual, lastResidual,
14 U x2d, numIters,
15 O nIterMin,
16 I printResidFrq, myThid )
17 C !DESCRIPTION: \bv
18 C *==========================================================*
19 C | SUBROUTINE CG2D
20 C | o Two-dimensional grid problem conjugate-gradient
21 C | inverter (with preconditioner).
22 C *==========================================================*
23 C | Con. grad is an iterative procedure for solving Ax = b.
24 C | It requires the A be symmetric.
25 C | This implementation assumes A is a five-diagonal
26 C | matrix of the form that arises in the discrete
27 C | representation of the del^2 operator in a
28 C | two-dimensional space.
29 C *==========================================================*
30 C \ev
31
32 C !USES:
33 IMPLICIT NONE
34 C === Global data ===
35 #include "SIZE.h"
36 #include "EEPARAMS.h"
37 #include "PARAMS.h"
38
39 C !INPUT/OUTPUT PARAMETERS:
40 C === Routine arguments ===
41 C b2d :: The source term or "right hand side"
42 C x2d :: The solution
43 C offDiagFactor :: preconditioner off-diagonal factor
44 C residCriter :: residual target for convergence
45 C firstResidual :: the initial residual before any iterations
46 C minResidual :: the lowest residual reached
47 C lastResidual :: the actual residual reached
48 C numIters :: Entry: the maximum number of iterations allowed
49 C Exit: the actual number of iterations used
50 C nIterMin :: iteration number corresponding to lowest residual
51 C printResidFrq :: Frequency for printing residual in CG iterations
52 C myThid :: my Thread Id number
53 _RS aW2d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
54 _RS aS2d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
55 _RL b2d (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
56 _RL residCriter
57 _RL offDiagFactor
58 _RL firstResidual
59 _RL minResidual
60 _RL lastResidual
61 _RL x2d (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
62 INTEGER numIters
63 INTEGER nIterMin
64 INTEGER printResidFrq
65 INTEGER myThid
66
67 C !LOCAL VARIABLES:
68 C === Local variables ====
69 C bi, bj :: tile indices
70 C eta_qrN :: Used in computing search directions
71 C eta_qrNM1 suffix N and NM1 denote current and
72 C cgBeta previous iterations respectively.
73 C alpha
74 C sumRHS :: Sum of right-hand-side. Sometimes this is a
75 C useful debuggin/trouble shooting diagnostic.
76 C For neumann problems sumRHS needs to be ~0.
77 C or they converge at a non-zero residual.
78 C err :: Measure of current residual of Ax - b, usually the norm.
79 C i, j, it2d :: Loop counters ( it2d counts CG iterations )
80 INTEGER bi, bj
81 INTEGER i, j, it2d
82 _RL err, errTile(nSx,nSy)
83 _RL eta_qrN, eta_qrNtile(nSx,nSy)
84 _RL eta_qrNM1
85 _RL cgBeta
86 _RL alpha, alphaTile(nSx,nSy)
87 _RL sumRHS, sumRHStile(nSx,nSy)
88 _RL pW_tmp, pS_tmp
89 CHARACTER*(MAX_LEN_MBUF) msgBuf
90 LOGICAL printResidual
91 CEOP
92 _RS aC2d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
93 _RS pW (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
94 _RS pS (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
95 _RS pC (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
96 _RL q2d(0:sNx+1,0:sNy+1,nSx,nSy)
97 #ifdef DEBUG_DIAG_CG2D
98 CHARACTER*(10) sufx
99 _RL r2d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
100 #else
101 _RL r2d(0:sNx+1,0:sNy+1,nSx,nSy)
102 #endif
103 _RL s2d(0:sNx+1,0:sNy+1,nSx,nSy)
104 _RL x2dm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
105
106 #ifdef ALLOW_DEBUG
107 IF (debugMode) CALL DEBUG_ENTER('DIAG_CG2D',myThid)
108 #endif
109
110 C-- Set matrice main diagnonal:
111 DO bj=myByLo(myThid),myByHi(myThid)
112 DO bi=myBxLo(myThid),myBxHi(myThid)
113 C- Initialise overlap regions (in case EXCH call do not fill them all)
114 DO j = 1-OLy,sNy+OLy
115 DO i = 1-OLx,sNx+OLx
116 aC2d(i,j,bi,bj) = 0.
117 ENDDO
118 ENDDO
119 DO j=1,sNy
120 DO i=1,sNx
121 aC2d(i,j,bi,bj) = -( ( aW2d(i,j,bi,bj)+aW2d(i+1,j,bi,bj) )
122 & +( aS2d(i,j,bi,bj)+aS2d(i,j+1,bi,bj) )
123 & )
124 ENDDO
125 ENDDO
126 ENDDO
127 ENDDO
128 CALL EXCH_XY_RS(aC2d, myThid)
129
130 C-- Initialise preconditioner
131 DO bj=myByLo(myThid),myByHi(myThid)
132 DO bi=myBxLo(myThid),myBxHi(myThid)
133 DO j=1,sNy+1
134 DO i=1,sNx+1
135 IF ( aC2d(i,j,bi,bj) .EQ. 0. ) THEN
136 pC(i,j,bi,bj) = 1. _d 0
137 ELSE
138 pC(i,j,bi,bj) = 1. _d 0 / aC2d(i,j,bi,bj)
139 ENDIF
140 pW_tmp = aC2d(i,j,bi,bj)+aC2d(i-1,j,bi,bj)
141 IF ( pW_tmp .EQ. 0. ) THEN
142 pW(i,j,bi,bj) = 0.
143 ELSE
144 pW(i,j,bi,bj) = -aW2d(i,j,bi,bj)*offDiagFactor
145 & *4. _d 0/( pW_tmp*pW_tmp )
146 ENDIF
147 pS_tmp = aC2d(i,j,bi,bj)+aC2d(i,j-1,bi,bj)
148 IF ( pS_tmp .EQ. 0. ) THEN
149 pS(i,j,bi,bj) = 0.
150 ELSE
151 pS(i,j,bi,bj) = -aS2d(i,j,bi,bj)*offDiagFactor
152 & *4. _d 0/( pS_tmp*pS_tmp )
153 ENDIF
154 c pC(i,j,bi,bj) = 1.
155 c pW(i,j,bi,bj) = 0.
156 c pS(i,j,bi,bj) = 0.
157 ENDDO
158 ENDDO
159 ENDDO
160 ENDDO
161
162 C-- Initialise inverter
163 eta_qrNM1 = 1. _d 0
164
165 CALL EXCH_XY_RL( x2d, myThid )
166
167 C-- Initial residual calculation
168 DO bj=myByLo(myThid),myByHi(myThid)
169 DO bi=myBxLo(myThid),myBxHi(myThid)
170 DO j=0,sNy+1
171 DO i=0,sNx+1
172 r2d(i,j,bi,bj) = 0.
173 s2d(i,j,bi,bj) = 0.
174 x2dm(i,j,bi,bj) = x2d(i,j,bi,bj)
175 ENDDO
176 ENDDO
177 sumRHStile(bi,bj) = 0. _d 0
178 errTile(bi,bj) = 0. _d 0
179 DO j=1,sNy
180 DO i=1,sNx
181 r2d(i,j,bi,bj) = b2d(i,j,bi,bj) -
182 & (aW2d(i ,j ,bi,bj)*x2d(i-1,j ,bi,bj)
183 & +aW2d(i+1,j ,bi,bj)*x2d(i+1,j ,bi,bj)
184 & +aS2d(i ,j ,bi,bj)*x2d(i ,j-1,bi,bj)
185 & +aS2d(i ,j+1,bi,bj)*x2d(i ,j+1,bi,bj)
186 & +aC2d(i ,j ,bi,bj)*x2d(i ,j ,bi,bj)
187 & )
188 errTile(bi,bj) = errTile(bi,bj)
189 & + r2d(i,j,bi,bj)*r2d(i,j,bi,bj)
190 sumRHStile(bi,bj) = sumRHStile(bi,bj) + b2d(i,j,bi,bj)
191 ENDDO
192 ENDDO
193 ENDDO
194 ENDDO
195 #ifdef DEBUG_DIAG_CG2D
196 CALL EXCH_XY_RL ( r2d, myThid )
197 #else
198 CALL EXCH_S3D_RL( r2d, 1, myThid )
199 #endif
200 CALL GLOBAL_SUM_TILE_RL( errTile, err, myThid )
201 IF ( printResidFrq.GE.1 )
202 & CALL GLOBAL_SUM_TILE_RL( sumRHStile, sumRHS, myThid )
203 err = SQRT(err)
204 it2d = 0
205 firstResidual = err
206 minResidual = err
207 nIterMin = it2d
208
209 printResidual = .FALSE.
210 IF ( debugLevel .GE. debLevZero ) THEN
211 _BEGIN_MASTER( myThid )
212 printResidual = printResidFrq.GE.1
213 IF ( printResidual ) THEN
214 WRITE(msgBuf,'(2A,I6,A,1PE17.9,A,1PE14.6)')' diag_cg2d:',
215 & ' iter=', it2d, ' ; resid.=', err, ' ; sumRHS=', sumRHS
216 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
217 & SQUEEZE_RIGHT, myThid )
218 ENDIF
219 _END_MASTER( myThid )
220 ENDIF
221 #ifdef DEBUG_DIAG_CG2D
222 IF ( printResidFrq.GE.1 ) THEN
223 WRITE(sufx,'(I10.10)') 0
224 CALL WRITE_FLD_XY_RL( 'r2d.',sufx, r2d, 1, myThid )
225 ENDIF
226 #endif
227
228 IF ( err .LT. residCriter ) GOTO 11
229
230 C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
231 DO 10 it2d=1, numIters
232
233 C-- Solve preconditioning equation and update
234 C-- conjugate direction vector "s".
235 DO bj=myByLo(myThid),myByHi(myThid)
236 DO bi=myBxLo(myThid),myBxHi(myThid)
237 eta_qrNtile(bi,bj) = 0. _d 0
238 DO j=1,sNy
239 DO i=1,sNx
240 q2d(i,j,bi,bj) =
241 & pC(i ,j ,bi,bj)*r2d(i ,j ,bi,bj)
242 & +pW(i ,j ,bi,bj)*r2d(i-1,j ,bi,bj)
243 & +pW(i+1,j ,bi,bj)*r2d(i+1,j ,bi,bj)
244 & +pS(i ,j ,bi,bj)*r2d(i ,j-1,bi,bj)
245 & +pS(i ,j+1,bi,bj)*r2d(i ,j+1,bi,bj)
246 eta_qrNtile(bi,bj) = eta_qrNtile(bi,bj)
247 & +q2d(i,j,bi,bj)*r2d(i,j,bi,bj)
248 ENDDO
249 ENDDO
250 ENDDO
251 ENDDO
252
253 CALL GLOBAL_SUM_TILE_RL( eta_qrNtile,eta_qrN,myThid )
254 cgBeta = eta_qrN/eta_qrNM1
255 eta_qrNM1 = eta_qrN
256
257 DO bj=myByLo(myThid),myByHi(myThid)
258 DO bi=myBxLo(myThid),myBxHi(myThid)
259 DO j=1,sNy
260 DO i=1,sNx
261 s2d(i,j,bi,bj) = q2d(i,j,bi,bj)
262 & + cgBeta*s2d(i,j,bi,bj)
263 ENDDO
264 ENDDO
265 ENDDO
266 ENDDO
267
268 C-- Do exchanges that require messages i.e. between processes.
269 CALL EXCH_S3D_RL( s2d, 1, myThid )
270
271 C== Evaluate laplace operator on conjugate gradient vector
272 C== q = A.s
273 DO bj=myByLo(myThid),myByHi(myThid)
274 DO bi=myBxLo(myThid),myBxHi(myThid)
275 alphaTile(bi,bj) = 0. _d 0
276 DO j=1,sNy
277 DO i=1,sNx
278 q2d(i,j,bi,bj) =
279 & aW2d(i ,j ,bi,bj)*s2d(i-1,j ,bi,bj)
280 & +aW2d(i+1,j ,bi,bj)*s2d(i+1,j ,bi,bj)
281 & +aS2d(i ,j ,bi,bj)*s2d(i ,j-1,bi,bj)
282 & +aS2d(i ,j+1,bi,bj)*s2d(i ,j+1,bi,bj)
283 & +aC2d(i ,j ,bi,bj)*s2d(i ,j ,bi,bj)
284 alphaTile(bi,bj) = alphaTile(bi,bj)
285 & + s2d(i,j,bi,bj)*q2d(i,j,bi,bj)
286 ENDDO
287 ENDDO
288 ENDDO
289 ENDDO
290 CALL GLOBAL_SUM_TILE_RL( alphaTile, alpha, myThid )
291 alpha = eta_qrN/alpha
292
293 C== Update solution and residual vectors
294 C Now compute "interior" points.
295 DO bj=myByLo(myThid),myByHi(myThid)
296 DO bi=myBxLo(myThid),myBxHi(myThid)
297 errTile(bi,bj) = 0. _d 0
298 DO j=1,sNy
299 DO i=1,sNx
300 x2d(i,j,bi,bj)=x2d(i,j,bi,bj)+alpha*s2d(i,j,bi,bj)
301 r2d(i,j,bi,bj)=r2d(i,j,bi,bj)-alpha*q2d(i,j,bi,bj)
302 errTile(bi,bj) = errTile(bi,bj)
303 & + r2d(i,j,bi,bj)*r2d(i,j,bi,bj)
304 ENDDO
305 ENDDO
306 ENDDO
307 ENDDO
308
309 CALL GLOBAL_SUM_TILE_RL( errTile, err, myThid )
310 err = SQRT(err)
311 IF ( printResidual ) THEN
312 IF ( MOD( it2d-1, printResidFrq ).EQ.0 ) THEN
313 WRITE(msgBuf,'(A,I6,A,1PE17.9)')
314 & ' diag_cg2d: iter=', it2d, ' ; resid.=', err
315 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
316 & SQUEEZE_RIGHT, myThid )
317 ENDIF
318 ENDIF
319 IF ( err .LT. residCriter ) GOTO 11
320 IF ( err .LT. minResidual ) THEN
321 C- Store lowest residual solution
322 minResidual = err
323 nIterMin = it2d
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 x2dm(i,j,bi,bj) = x2d(i,j,bi,bj)
329 ENDDO
330 ENDDO
331 ENDDO
332 ENDDO
333 ENDIF
334
335 #ifdef DEBUG_DIAG_CG2D
336 CALL EXCH_XY_RL( r2d, myThid )
337 IF ( printResidFrq.GE.1 ) THEN
338 WRITE(sufx,'(I10.10)') it2d
339 CALL WRITE_FLD_XY_RL( 'r2d.',sufx, r2d, 1, myThid )
340 CALL WRITE_FLD_XY_RL( 'x2d.',sufx, x2d, 1, myThid )
341 ENDIF
342 #else
343 CALL EXCH_S3D_RL( r2d, 1, myThid )
344 #endif
345
346 10 CONTINUE
347 it2d = numIters
348 11 CONTINUE
349
350 C-- Return parameters to caller
351 lastResidual = err
352 numIters = it2d
353
354 IF ( err .GT. minResidual ) THEN
355 C- use the lowest residual solution (instead of current one <-> last residual)
356 DO bj=myByLo(myThid),myByHi(myThid)
357 DO bi=myBxLo(myThid),myBxHi(myThid)
358 DO j=1,sNy
359 DO i=1,sNx
360 x2d(i,j,bi,bj) = x2dm(i,j,bi,bj)
361 ENDDO
362 ENDDO
363 ENDDO
364 ENDDO
365 ENDIF
366 c CALL EXCH_XY_RL( x2d, myThid )
367
368 #ifdef ALLOW_DEBUG
369 IF (debugMode) CALL DEBUG_LEAVE('DIAG_CG2D',myThid)
370 #endif
371
372 RETURN
373 END

  ViewVC Help
Powered by ViewVC 1.1.22