/[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.1 - (show annotations) (download)
Tue Jun 14 00:18:36 2011 UTC (12 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint63, checkpoint62z
first attempt to compute Velocity Potential (from UVELMASS & VVELMASS diags).

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

  ViewVC Help
Powered by ViewVC 1.1.22