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

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

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


Revision 1.3 - (hide annotations) (download)
Mon Aug 1 20:40:43 2011 UTC (12 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint63a
Changes since 1.2: +43 -11 lines
Add some debugging pieces of code

1 jmc 1.3 C $Header: /u/gcmpack/MITgcm/pkg/diagnostics/diag_cg2d.F,v 1.2 2011/07/22 19:53:39 jmc Exp $
2 jmc 1.1 C $Name: $
3    
4     #include "DIAG_OPTIONS.h"
5 jmc 1.3 #undef DEBUG_DIAG_CG2D
6 jmc 1.1
7     CBOP
8     C !ROUTINE: DIAG_CG2D
9     C !INTERFACE:
10     SUBROUTINE DIAG_CG2D(
11     I aW2d, aS2d, b2d,
12     I residCriter,
13 jmc 1.2 O firstResidual, minResidual, lastResidual,
14 jmc 1.1 U x2d, numIters,
15 jmc 1.2 O nIterMin,
16     I printResidFrq, myThid )
17 jmc 1.1 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 jmc 1.2 C b2d :: The source term or "right hand side"
42     C x2d :: The solution
43 jmc 1.1 C firstResidual :: the initial residual before any iterations
44 jmc 1.2 C minResidual :: the lowest residual reached
45 jmc 1.1 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 jmc 1.2 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 jmc 1.1 _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 jmc 1.2 _RL minResidual
57 jmc 1.1 _RL lastResidual
58     _RL x2d (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
59     INTEGER numIters
60 jmc 1.2 INTEGER nIterMin
61     INTEGER printResidFrq
62 jmc 1.1 INTEGER myThid
63    
64     C !LOCAL VARIABLES:
65     C === Local variables ====
66 jmc 1.2 C bi, bj :: tile indices
67 jmc 1.1 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 jmc 1.2 C err :: Measure of current residual of Ax - b, usually the norm.
76 jmc 1.1 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 jmc 1.3 _RL diagCG_pcOffDFac
87 jmc 1.1 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(1-1:sNx+1,1-1:sNy+1,nSx,nSy)
95 jmc 1.3 #ifdef DEBUG_DIAG_CG2D
96     CHARACTER*(10) sufx
97     _RL r2d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
98     #else
99 jmc 1.1 _RL r2d(1-1:sNx+1,1-1:sNy+1,nSx,nSy)
100 jmc 1.3 #endif
101 jmc 1.1 _RL s2d(1-1:sNx+1,1-1:sNy+1,nSx,nSy)
102 jmc 1.2 _RL x2dm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
103 jmc 1.1
104     C-- Set matrice main diagnonal:
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     aC2d(i,j,bi,bj) = -( ( aW2d(i,j,bi,bj)+aW2d(i+1,j,bi,bj) )
110     & +( aS2d(i,j,bi,bj)+aS2d(i,j+1,bi,bj) )
111     & )
112     ENDDO
113     ENDDO
114     ENDDO
115     ENDDO
116     CALL EXCH_XY_RS(aC2d, myThid)
117    
118     C-- Initialise preconditioner
119 jmc 1.3 diagCG_pcOffDFac = cg2dpcOffDFac
120 jmc 1.1 DO bj=myByLo(myThid),myByHi(myThid)
121     DO bi=myBxLo(myThid),myBxHi(myThid)
122     DO j=1,sNy+1
123     DO i=1,sNx+1
124     IF ( aC2d(i,j,bi,bj) .EQ. 0. ) THEN
125     pC(i,j,bi,bj) = 1. _d 0
126     ELSE
127     pC(i,j,bi,bj) = 1. _d 0 / aC2d(i,j,bi,bj)
128     ENDIF
129     pW_tmp = aC2d(i,j,bi,bj)+aC2d(i-1,j,bi,bj)
130     IF ( pW_tmp .EQ. 0. ) THEN
131     pW(i,j,bi,bj) = 0.
132     ELSE
133     pW(i,j,bi,bj) = -aW2d(i,j,bi,bj)
134 jmc 1.3 & /( (diagCG_pcOffDFac*pW_tmp)**2 )
135 jmc 1.1 ENDIF
136     pS_tmp = aC2d(i,j,bi,bj)+aC2d(i,j-1,bi,bj)
137     IF ( pS_tmp .EQ. 0. ) THEN
138     pS(i,j,bi,bj) = 0.
139     ELSE
140     pS(i,j,bi,bj) = -aS2d(i,j,bi,bj)
141 jmc 1.3 & /( (diagCG_pcOffDFac*pS_tmp)**2 )
142 jmc 1.1 ENDIF
143 jmc 1.3 c pC(i,j,bi,bj) = 1.
144     c pW(i,j,bi,bj) = 0.
145     c pS(i,j,bi,bj) = 0.
146 jmc 1.1 ENDDO
147     ENDDO
148     ENDDO
149     ENDDO
150    
151     C-- Initialise inverter
152     eta_qrNM1 = 1. _d 0
153    
154     CALL EXCH_XY_RL( x2d, myThid )
155    
156     C-- Initial residual calculation
157     DO bj=myByLo(myThid),myByHi(myThid)
158     DO bi=myBxLo(myThid),myBxHi(myThid)
159     DO j=1-1,sNy+1
160     DO i=1-1,sNx+1
161     s2d(i,j,bi,bj) = 0.
162 jmc 1.2 x2dm(i,j,bi,bj) = x2d(i,j,bi,bj)
163 jmc 1.1 ENDDO
164     ENDDO
165     sumRHStile(bi,bj) = 0. _d 0
166     errTile(bi,bj) = 0. _d 0
167     DO j=1,sNy
168     DO i=1,sNx
169     r2d(i,j,bi,bj) = b2d(i,j,bi,bj) -
170     & (aW2d(i ,j ,bi,bj)*x2d(i-1,j ,bi,bj)
171     & +aW2d(i+1,j ,bi,bj)*x2d(i+1,j ,bi,bj)
172     & +aS2d(i ,j ,bi,bj)*x2d(i ,j-1,bi,bj)
173     & +aS2d(i ,j+1,bi,bj)*x2d(i ,j+1,bi,bj)
174     & +aC2d(i ,j ,bi,bj)*x2d(i ,j ,bi,bj)
175     & )
176     errTile(bi,bj) = errTile(bi,bj)
177     & + r2d(i,j,bi,bj)*r2d(i,j,bi,bj)
178     sumRHStile(bi,bj) = sumRHStile(bi,bj) + b2d(i,j,bi,bj)
179     ENDDO
180     ENDDO
181     ENDDO
182     ENDDO
183 jmc 1.3 #ifdef DEBUG_DIAG_CG2D
184     CALL EXCH_XY_RL ( r2d, myThid )
185     #else
186 jmc 1.1 CALL EXCH_S3D_RL( r2d, 1, myThid )
187 jmc 1.3 #endif
188 jmc 1.1 CALL GLOBAL_SUM_TILE_RL( errTile, err, myThid )
189 jmc 1.3 IF ( printResidFrq.GE.1 )
190     & CALL GLOBAL_SUM_TILE_RL( sumRHStile, sumRHS, myThid )
191 jmc 1.1 err = SQRT(err)
192 jmc 1.2 it2d = 0
193     firstResidual = err
194     minResidual = err
195     nIterMin = it2d
196 jmc 1.1
197     printResidual = .FALSE.
198     IF ( debugLevel .GE. debLevZero ) THEN
199     _BEGIN_MASTER( myThid )
200 jmc 1.2 printResidual = printResidFrq.GE.1
201 jmc 1.3 IF ( printResidual ) THEN
202     WRITE(msgBuf,'(2A,I6,A,1PE17.9,A,1PE14.6)')' diag_cg2d:',
203     & ' iter=', it2d, ' ; resid.=', err, ' ; sumRHS=', sumRHS
204     CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
205     & SQUEEZE_RIGHT, myThid )
206     ENDIF
207 jmc 1.1 _END_MASTER( myThid )
208     ENDIF
209 jmc 1.3 #ifdef DEBUG_DIAG_CG2D
210     IF ( printResidFrq.GE.1 ) THEN
211     WRITE(sufx,'(I10.10)') 0
212     CALL WRITE_FLD_XY_RL( 'r2d.',sufx, r2d, 1, myThid )
213     ENDIF
214     #endif
215 jmc 1.1
216     IF ( err .LT. residCriter ) GOTO 11
217    
218     C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
219     DO 10 it2d=1, numIters
220    
221     C-- Solve preconditioning equation and update
222     C-- conjugate direction vector "s".
223     DO bj=myByLo(myThid),myByHi(myThid)
224     DO bi=myBxLo(myThid),myBxHi(myThid)
225     eta_qrNtile(bi,bj) = 0. _d 0
226     DO j=1,sNy
227     DO i=1,sNx
228     q2d(i,j,bi,bj) =
229     & pC(i ,j ,bi,bj)*r2d(i ,j ,bi,bj)
230     & +pW(i ,j ,bi,bj)*r2d(i-1,j ,bi,bj)
231     & +pW(i+1,j ,bi,bj)*r2d(i+1,j ,bi,bj)
232     & +pS(i ,j ,bi,bj)*r2d(i ,j-1,bi,bj)
233     & +pS(i ,j+1,bi,bj)*r2d(i ,j+1,bi,bj)
234     eta_qrNtile(bi,bj) = eta_qrNtile(bi,bj)
235     & +q2d(i,j,bi,bj)*r2d(i,j,bi,bj)
236     ENDDO
237     ENDDO
238     ENDDO
239     ENDDO
240    
241     CALL GLOBAL_SUM_TILE_RL( eta_qrNtile,eta_qrN,myThid )
242     cgBeta = eta_qrN/eta_qrNM1
243     eta_qrNM1 = eta_qrN
244    
245     DO bj=myByLo(myThid),myByHi(myThid)
246     DO bi=myBxLo(myThid),myBxHi(myThid)
247     DO j=1,sNy
248     DO i=1,sNx
249     s2d(i,j,bi,bj) = q2d(i,j,bi,bj)
250     & + cgBeta*s2d(i,j,bi,bj)
251     ENDDO
252     ENDDO
253     ENDDO
254     ENDDO
255    
256     C-- Do exchanges that require messages i.e. between processes.
257     CALL EXCH_S3D_RL( s2d, 1, myThid )
258    
259     C== Evaluate laplace operator on conjugate gradient vector
260     C== q = A.s
261     DO bj=myByLo(myThid),myByHi(myThid)
262     DO bi=myBxLo(myThid),myBxHi(myThid)
263     alphaTile(bi,bj) = 0. _d 0
264     DO j=1,sNy
265     DO i=1,sNx
266     q2d(i,j,bi,bj) =
267     & aW2d(i ,j ,bi,bj)*s2d(i-1,j ,bi,bj)
268     & +aW2d(i+1,j ,bi,bj)*s2d(i+1,j ,bi,bj)
269     & +aS2d(i ,j ,bi,bj)*s2d(i ,j-1,bi,bj)
270     & +aS2d(i ,j+1,bi,bj)*s2d(i ,j+1,bi,bj)
271     & +aC2d(i ,j ,bi,bj)*s2d(i ,j ,bi,bj)
272     alphaTile(bi,bj) = alphaTile(bi,bj)
273     & + s2d(i,j,bi,bj)*q2d(i,j,bi,bj)
274     ENDDO
275     ENDDO
276     ENDDO
277     ENDDO
278     CALL GLOBAL_SUM_TILE_RL( alphaTile, alpha, myThid )
279     alpha = eta_qrN/alpha
280    
281     C== Update solution and residual vectors
282     C Now compute "interior" points.
283     DO bj=myByLo(myThid),myByHi(myThid)
284     DO bi=myBxLo(myThid),myBxHi(myThid)
285     errTile(bi,bj) = 0. _d 0
286     DO j=1,sNy
287     DO i=1,sNx
288     x2d(i,j,bi,bj)=x2d(i,j,bi,bj)+alpha*s2d(i,j,bi,bj)
289     r2d(i,j,bi,bj)=r2d(i,j,bi,bj)-alpha*q2d(i,j,bi,bj)
290     errTile(bi,bj) = errTile(bi,bj)
291     & + r2d(i,j,bi,bj)*r2d(i,j,bi,bj)
292     ENDDO
293     ENDDO
294     ENDDO
295     ENDDO
296    
297     CALL GLOBAL_SUM_TILE_RL( errTile, err, myThid )
298     err = SQRT(err)
299     IF ( printResidual ) THEN
300 jmc 1.2 IF ( MOD( it2d-1, printResidFrq ).EQ.0 ) THEN
301 jmc 1.3 WRITE(msgBuf,'(A,I6,A,1PE17.9)')
302     & ' diag_cg2d: iter=', it2d, ' ; resid.=', err
303 jmc 1.1 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
304     & SQUEEZE_RIGHT, myThid )
305     ENDIF
306     ENDIF
307     IF ( err .LT. residCriter ) GOTO 11
308 jmc 1.2 IF ( err .LT. minResidual ) THEN
309     C- Store lowest residual solution
310     minResidual = err
311     nIterMin = it2d
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     x2dm(i,j,bi,bj) = x2d(i,j,bi,bj)
317     ENDDO
318     ENDDO
319     ENDDO
320     ENDDO
321     ENDIF
322 jmc 1.1
323 jmc 1.3 #ifdef DEBUG_DIAG_CG2D
324     CALL EXCH_XY_RL( r2d, myThid )
325     IF ( printResidFrq.GE.1 ) THEN
326     WRITE(sufx,'(I10.10)') it2d
327     CALL WRITE_FLD_XY_RL( 'r2d.',sufx, r2d, 1, myThid )
328     CALL WRITE_FLD_XY_RL( 'x2d.',sufx, x2d, 1, myThid )
329     ENDIF
330     #else
331 jmc 1.1 CALL EXCH_S3D_RL( r2d, 1, myThid )
332 jmc 1.3 #endif
333 jmc 1.1
334     10 CONTINUE
335 jmc 1.2 it2d = numIters
336 jmc 1.1 11 CONTINUE
337    
338 jmc 1.2 C-- Return parameters to caller
339     lastResidual = err
340     numIters = it2d
341    
342     IF ( err .GT. minResidual ) THEN
343     C- use the lowest residual solution (instead of current one <-> last residual)
344     DO bj=myByLo(myThid),myByHi(myThid)
345     DO bi=myBxLo(myThid),myBxHi(myThid)
346     DO j=1,sNy
347     DO i=1,sNx
348     x2d(i,j,bi,bj) = x2dm(i,j,bi,bj)
349     ENDDO
350     ENDDO
351     ENDDO
352     ENDDO
353     ENDIF
354 jmc 1.1 c CALL EXCH_XY_RL( x2d, myThid )
355    
356     RETURN
357     END

  ViewVC Help
Powered by ViewVC 1.1.22