/[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.5 - (show annotations) (download)
Mon Mar 19 18:22:09 2012 UTC (12 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint64, checkpoint65, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65d, checkpoint65e
Changes since 1.4: +9 -8 lines
initialise r2d local array (in case there are blank-tiles)

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

  ViewVC Help
Powered by ViewVC 1.1.22