/[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.2 - (show annotations) (download)
Fri Jul 22 19:53:39 2011 UTC (12 years, 9 months ago) by jmc
Branch: MAIN
Changes since 1.1: +55 -24 lines
- fix mask for OBCS (still problems in stream-function with OBCS);
- add specific parameter (default = main code CG2D params) for solver;
- in case of poor convergence, use solution corresponding to lowest residual.

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

  ViewVC Help
Powered by ViewVC 1.1.22