1 |
C $Header: /u/gcmpack/MITgcm/model/src/cg2d.F,v 1.38 2004/02/06 20:21:00 adcroft Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CPP_OPTIONS.h" |
5 |
|
6 |
CBOP |
7 |
C !ROUTINE: CG2D |
8 |
C !INTERFACE: |
9 |
SUBROUTINE CG2D( |
10 |
I cg2d_b, |
11 |
U cg2d_x, |
12 |
O firstResidual, |
13 |
O lastResidual, |
14 |
U numIters, |
15 |
I 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 | Notes: |
29 |
C | ====== |
30 |
C | This implementation can support shared-memory |
31 |
C | multi-threaded execution. In order to do this COMMON |
32 |
C | blocks are used for many of the arrays - even ones that |
33 |
C | are only used for intermedaite results. This design is |
34 |
C | OK if you want to all the threads to collaborate on |
35 |
C | solving the same problem. On the other hand if you want |
36 |
C | the threads to solve several different problems |
37 |
C | concurrently this implementation will not work. |
38 |
C *==========================================================* |
39 |
C \ev |
40 |
|
41 |
C !USES: |
42 |
IMPLICIT NONE |
43 |
C === Global data === |
44 |
#include "SIZE.h" |
45 |
#include "EEPARAMS.h" |
46 |
#include "PARAMS.h" |
47 |
#include "GRID.h" |
48 |
#include "CG2D.h" |
49 |
#include "SURFACE.h" |
50 |
|
51 |
C !INPUT/OUTPUT PARAMETERS: |
52 |
C === Routine arguments === |
53 |
C myThid - Thread on which I am working. |
54 |
C cg2d_b - The source term or "right hand side" |
55 |
C cg2d_x - The solution |
56 |
C firstResidual - the initial residual before any iterations |
57 |
C lastResidual - the actual residual reached |
58 |
C numIters - Entry: the maximum number of iterations allowed |
59 |
C Exit: the actual number of iterations used |
60 |
_RL cg2d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
61 |
_RL cg2d_x(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
62 |
_RL firstResidual |
63 |
_RL lastResidual |
64 |
INTEGER numIters |
65 |
INTEGER myThid |
66 |
|
67 |
C !LOCAL VARIABLES: |
68 |
C === Local variables ==== |
69 |
C actualIts - Number of iterations taken |
70 |
C actualResidual - residual |
71 |
C bi - Block index in X and Y. |
72 |
C bj |
73 |
C eta_qrN - Used in computing search directions |
74 |
C eta_qrNM1 suffix N and NM1 denote current and |
75 |
C cgBeta previous iterations respectively. |
76 |
C alpha |
77 |
C sumRHS - Sum of right-hand-side. Sometimes this is a |
78 |
C useful debuggin/trouble shooting diagnostic. |
79 |
C For neumann problems sumRHS needs to be ~0. |
80 |
C or they converge at a non-zero residual. |
81 |
C err - Measure of residual of Ax - b, usually the norm. |
82 |
C I, J, N - Loop counters ( N counts CG iterations ) |
83 |
INTEGER actualIts |
84 |
_RL actualResidual |
85 |
INTEGER bi, bj |
86 |
INTEGER I, J, it2d |
87 |
_RL err,errTile |
88 |
_RL eta_qrN,eta_qrNtile |
89 |
_RL eta_qrNM1 |
90 |
_RL cgBeta |
91 |
_RL alpha,alphaTile |
92 |
_RL sumRHS,sumRHStile |
93 |
_RL rhsMax |
94 |
_RL rhsNorm |
95 |
|
96 |
INTEGER OLw |
97 |
INTEGER OLe |
98 |
INTEGER OLn |
99 |
INTEGER OLs |
100 |
INTEGER exchWidthX |
101 |
INTEGER exchWidthY |
102 |
INTEGER myNz |
103 |
CEOP |
104 |
|
105 |
|
106 |
CcnhDebugStarts |
107 |
C CHARACTER*(MAX_LEN_FNAM) suff |
108 |
CcnhDebugEnds |
109 |
|
110 |
|
111 |
C-- Initialise inverter |
112 |
eta_qrNM1 = 1. _d 0 |
113 |
|
114 |
CcnhDebugStarts |
115 |
C _EXCH_XY_R8( cg2d_b, myThid ) |
116 |
C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.0 CG2D_B' , 1, myThid ) |
117 |
C suff = 'unnormalised' |
118 |
C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid) |
119 |
C STOP |
120 |
CcnhDebugEnds |
121 |
|
122 |
C-- Normalise RHS |
123 |
rhsMax = 0. _d 0 |
124 |
DO bj=myByLo(myThid),myByHi(myThid) |
125 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
126 |
DO J=1,sNy |
127 |
DO I=1,sNx |
128 |
cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*cg2dNorm |
129 |
rhsMax = MAX(ABS(cg2d_b(I,J,bi,bj)),rhsMax) |
130 |
ENDDO |
131 |
ENDDO |
132 |
ENDDO |
133 |
ENDDO |
134 |
|
135 |
IF (cg2dNormaliseRHS) THEN |
136 |
C- Normalise RHS : |
137 |
#ifdef LETS_MAKE_JAM |
138 |
C _GLOBAL_MAX_R8( rhsMax, myThid ) |
139 |
rhsMax=1. |
140 |
#else |
141 |
_GLOBAL_MAX_R8( rhsMax, myThid ) |
142 |
Catm rhsMax=1. |
143 |
#endif |
144 |
rhsNorm = 1. _d 0 |
145 |
IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax |
146 |
DO bj=myByLo(myThid),myByHi(myThid) |
147 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
148 |
DO J=1,sNy |
149 |
DO I=1,sNx |
150 |
cg2d_b(I,J,bi,bj) = cg2d_b(I,J,bi,bj)*rhsNorm |
151 |
cg2d_x(I,J,bi,bj) = cg2d_x(I,J,bi,bj)*rhsNorm |
152 |
ENDDO |
153 |
ENDDO |
154 |
ENDDO |
155 |
ENDDO |
156 |
C- end Normalise RHS |
157 |
ENDIF |
158 |
|
159 |
C-- Update overlaps |
160 |
_EXCH_XY_R8( cg2d_b, myThid ) |
161 |
_EXCH_XY_R8( cg2d_x, myThid ) |
162 |
CcnhDebugStarts |
163 |
C CALL PLOT_FIELD_XYRL( cg2d_b, 'CG2D.1 CG2D_B' , 1, myThid ) |
164 |
C suff = 'normalised' |
165 |
C CALL WRITE_FLD_XY_RL ( 'cg2d_b.',suff, cg2d_b, 1, myThid) |
166 |
CcnhDebugEnds |
167 |
|
168 |
C-- Initial residual calculation |
169 |
err = 0. _d 0 |
170 |
sumRHS = 0. _d 0 |
171 |
DO bj=myByLo(myThid),myByHi(myThid) |
172 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
173 |
sumRHStile = 0. _d 0 |
174 |
errTile = 0. _d 0 |
175 |
DO J=1,sNy |
176 |
DO I=1,sNx |
177 |
cg2d_s(I,J,bi,bj) = 0. |
178 |
cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) - |
179 |
& (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj) |
180 |
& +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj) |
181 |
& +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj) |
182 |
& +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj) |
183 |
& -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
184 |
& -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
185 |
& -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
186 |
& -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj) |
187 |
& -freeSurfFac*_rA(i,j,bi,bj)*recip_Bo(i,j,bi,bj)* |
188 |
& cg2d_x(I ,J ,bi,bj)/deltaTMom/deltaTfreesurf*cg2dNorm |
189 |
& ) |
190 |
errTile = errTile + |
191 |
& cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj) |
192 |
sumRHStile = sumRHStile + |
193 |
& cg2d_b(I,J,bi,bj) |
194 |
ENDDO |
195 |
ENDDO |
196 |
sumRHS = sumRHS + sumRHStile |
197 |
err = err + errTile |
198 |
ENDDO |
199 |
ENDDO |
200 |
C _EXCH_XY_R8( cg2d_r, myThid ) |
201 |
#ifdef LETS_MAKE_JAM |
202 |
CALL EXCH_XY_O1_R8_JAM( cg2d_r ) |
203 |
#else |
204 |
CALL EXCH_XY_RL( cg2d_r, myThid ) |
205 |
#endif |
206 |
C _EXCH_XY_R8( cg2d_s, myThid ) |
207 |
#ifdef LETS_MAKE_JAM |
208 |
CALL EXCH_XY_O1_R8_JAM( cg2d_s ) |
209 |
#else |
210 |
CALL EXCH_XY_RL( cg2d_s, myThid ) |
211 |
#endif |
212 |
_GLOBAL_SUM_R8( sumRHS, myThid ) |
213 |
_GLOBAL_SUM_R8( err , myThid ) |
214 |
err = SQRT(err) |
215 |
actualIts = 0 |
216 |
actualResidual = err |
217 |
|
218 |
IF ( debugLevel .GE. debLevZero ) THEN |
219 |
_BEGIN_MASTER( myThid ) |
220 |
write(*,'(A,1P2E22.14)')' cg2d: Sum(rhs),rhsMax = ', |
221 |
& sumRHS,rhsMax |
222 |
_END_MASTER( ) |
223 |
ENDIF |
224 |
C _BARRIER |
225 |
c _BEGIN_MASTER( myThid ) |
226 |
c WRITE(*,'(A,I6,1PE30.14)') ' CG2D iters, err = ', |
227 |
c & actualIts, actualResidual |
228 |
c _END_MASTER( ) |
229 |
firstResidual=actualResidual |
230 |
|
231 |
IF ( err .LT. cg2dTolerance ) GOTO 11 |
232 |
|
233 |
C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
234 |
DO 10 it2d=1, numIters |
235 |
|
236 |
CcnhDebugStarts |
237 |
C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' residual = ', |
238 |
C & actualResidual |
239 |
CcnhDebugEnds |
240 |
C-- Solve preconditioning equation and update |
241 |
C-- conjugate direction vector "s". |
242 |
eta_qrN = 0. _d 0 |
243 |
DO bj=myByLo(myThid),myByHi(myThid) |
244 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
245 |
eta_qrNtile = 0. _d 0 |
246 |
DO J=1,sNy |
247 |
DO I=1,sNx |
248 |
cg2d_q(I,J,bi,bj) = |
249 |
& pC(I ,J ,bi,bj)*cg2d_r(I ,J ,bi,bj) |
250 |
& +pW(I ,J ,bi,bj)*cg2d_r(I-1,J ,bi,bj) |
251 |
& +pW(I+1,J ,bi,bj)*cg2d_r(I+1,J ,bi,bj) |
252 |
& +pS(I ,J ,bi,bj)*cg2d_r(I ,J-1,bi,bj) |
253 |
& +pS(I ,J+1,bi,bj)*cg2d_r(I ,J+1,bi,bj) |
254 |
CcnhDebugStarts |
255 |
C cg2d_q(I,J,bi,bj) = cg2d_r(I ,J ,bi,bj) |
256 |
CcnhDebugEnds |
257 |
eta_qrNtile = eta_qrNtile |
258 |
& +cg2d_q(I,J,bi,bj)*cg2d_r(I,J,bi,bj) |
259 |
ENDDO |
260 |
ENDDO |
261 |
eta_qrN = eta_qrN + eta_qrNtile |
262 |
ENDDO |
263 |
ENDDO |
264 |
|
265 |
_GLOBAL_SUM_R8(eta_qrN, myThid) |
266 |
CcnhDebugStarts |
267 |
C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' eta_qrN = ',eta_qrN |
268 |
CcnhDebugEnds |
269 |
cgBeta = eta_qrN/eta_qrNM1 |
270 |
CcnhDebugStarts |
271 |
C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' beta = ',cgBeta |
272 |
CcnhDebugEnds |
273 |
eta_qrNM1 = eta_qrN |
274 |
|
275 |
DO bj=myByLo(myThid),myByHi(myThid) |
276 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
277 |
DO J=1,sNy |
278 |
DO I=1,sNx |
279 |
cg2d_s(I,J,bi,bj) = cg2d_q(I,J,bi,bj) |
280 |
& + cgBeta*cg2d_s(I,J,bi,bj) |
281 |
ENDDO |
282 |
ENDDO |
283 |
ENDDO |
284 |
ENDDO |
285 |
|
286 |
C-- Do exchanges that require messages i.e. between |
287 |
C-- processes. |
288 |
C _EXCH_XY_R8( cg2d_s, myThid ) |
289 |
#ifdef LETS_MAKE_JAM |
290 |
CALL EXCH_XY_O1_R8_JAM( cg2d_s ) |
291 |
#else |
292 |
CALL EXCH_XY_RL( cg2d_s, myThid ) |
293 |
#endif |
294 |
|
295 |
C== Evaluate laplace operator on conjugate gradient vector |
296 |
C== q = A.s |
297 |
alpha = 0. _d 0 |
298 |
DO bj=myByLo(myThid),myByHi(myThid) |
299 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
300 |
alphaTile = 0. _d 0 |
301 |
DO J=1,sNy |
302 |
DO I=1,sNx |
303 |
cg2d_q(I,J,bi,bj) = |
304 |
& aW2d(I ,J ,bi,bj)*cg2d_s(I-1,J ,bi,bj) |
305 |
& +aW2d(I+1,J ,bi,bj)*cg2d_s(I+1,J ,bi,bj) |
306 |
& +aS2d(I ,J ,bi,bj)*cg2d_s(I ,J-1,bi,bj) |
307 |
& +aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J+1,bi,bj) |
308 |
& -aW2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj) |
309 |
& -aW2d(I+1,J ,bi,bj)*cg2d_s(I ,J ,bi,bj) |
310 |
& -aS2d(I ,J ,bi,bj)*cg2d_s(I ,J ,bi,bj) |
311 |
& -aS2d(I ,J+1,bi,bj)*cg2d_s(I ,J ,bi,bj) |
312 |
& -freeSurfFac*_rA(i,j,bi,bj)*recip_Bo(i,j,bi,bj)* |
313 |
& cg2d_s(I ,J ,bi,bj)/deltaTMom/deltaTfreesurf*cg2dNorm |
314 |
alphaTile = alphaTile+cg2d_s(I,J,bi,bj)*cg2d_q(I,J,bi,bj) |
315 |
ENDDO |
316 |
ENDDO |
317 |
alpha = alpha + alphaTile |
318 |
ENDDO |
319 |
ENDDO |
320 |
_GLOBAL_SUM_R8(alpha,myThid) |
321 |
CcnhDebugStarts |
322 |
C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' SUM(s*q)= ',alpha |
323 |
CcnhDebugEnds |
324 |
alpha = eta_qrN/alpha |
325 |
CcnhDebugStarts |
326 |
C WRITE(*,*) ' CG2D: Iteration ',it2d-1,' alpha= ',alpha |
327 |
CcnhDebugEnds |
328 |
|
329 |
C== Update solution and residual vectors |
330 |
C Now compute "interior" points. |
331 |
err = 0. _d 0 |
332 |
DO bj=myByLo(myThid),myByHi(myThid) |
333 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
334 |
errTile = 0. _d 0 |
335 |
DO J=1,sNy |
336 |
DO I=1,sNx |
337 |
cg2d_x(I,J,bi,bj)=cg2d_x(I,J,bi,bj)+alpha*cg2d_s(I,J,bi,bj) |
338 |
cg2d_r(I,J,bi,bj)=cg2d_r(I,J,bi,bj)-alpha*cg2d_q(I,J,bi,bj) |
339 |
errTile = errTile+cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj) |
340 |
ENDDO |
341 |
ENDDO |
342 |
err = err + errTile |
343 |
ENDDO |
344 |
ENDDO |
345 |
|
346 |
_GLOBAL_SUM_R8( err , myThid ) |
347 |
err = SQRT(err) |
348 |
actualIts = it2d |
349 |
actualResidual = err |
350 |
IF ( err .LT. cg2dTolerance ) GOTO 11 |
351 |
C _EXCH_XY_R8(cg2d_r, myThid ) |
352 |
#ifdef LETS_MAKE_JAM |
353 |
CALL EXCH_XY_O1_R8_JAM( cg2d_r ) |
354 |
#else |
355 |
CALL EXCH_XY_RL( cg2d_r, myThid ) |
356 |
#endif |
357 |
|
358 |
10 CONTINUE |
359 |
11 CONTINUE |
360 |
|
361 |
IF (cg2dNormaliseRHS) THEN |
362 |
C-- Un-normalise the answer |
363 |
DO bj=myByLo(myThid),myByHi(myThid) |
364 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
365 |
DO J=1,sNy |
366 |
DO I=1,sNx |
367 |
cg2d_x(I ,J ,bi,bj) = cg2d_x(I ,J ,bi,bj)/rhsNorm |
368 |
ENDDO |
369 |
ENDDO |
370 |
ENDDO |
371 |
ENDDO |
372 |
ENDIF |
373 |
|
374 |
C The following exchange was moved up to solve_for_pressure |
375 |
C for compatibility with TAMC. |
376 |
C _EXCH_XY_R8(cg2d_x, myThid ) |
377 |
c _BEGIN_MASTER( myThid ) |
378 |
c WRITE(*,'(A,I6,1PE30.14)') ' CG2D iters, err = ', |
379 |
c & actualIts, actualResidual |
380 |
c _END_MASTER( ) |
381 |
|
382 |
C-- Return parameters to caller |
383 |
lastResidual=actualResidual |
384 |
numIters=actualIts |
385 |
|
386 |
CcnhDebugStarts |
387 |
C CALL PLOT_FIELD_XYRL( cg2d_x, 'CALC_MOM_RHS CG2D_X' , 1, myThid ) |
388 |
C err = 0. _d 0 |
389 |
C DO bj=myByLo(myThid),myByHi(myThid) |
390 |
C DO bi=myBxLo(myThid),myBxHi(myThid) |
391 |
C DO J=1,sNy |
392 |
C DO I=1,sNx |
393 |
C cg2d_r(I,J,bi,bj) = cg2d_b(I,J,bi,bj) - |
394 |
C & (aW2d(I ,J ,bi,bj)*cg2d_x(I-1,J ,bi,bj) |
395 |
C & +aW2d(I+1,J ,bi,bj)*cg2d_x(I+1,J ,bi,bj) |
396 |
C & +aS2d(I ,J ,bi,bj)*cg2d_x(I ,J-1,bi,bj) |
397 |
C & +aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J+1,bi,bj) |
398 |
C & -aW2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
399 |
C & -aW2d(I+1,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
400 |
C & -aS2d(I ,J ,bi,bj)*cg2d_x(I ,J ,bi,bj) |
401 |
C & -aS2d(I ,J+1,bi,bj)*cg2d_x(I ,J ,bi,bj)) |
402 |
C err = err + |
403 |
C & cg2d_r(I,J,bi,bj)*cg2d_r(I,J,bi,bj) |
404 |
C ENDDO |
405 |
C ENDDO |
406 |
C ENDDO |
407 |
C ENDDO |
408 |
C _GLOBAL_SUM_R8( err , myThid ) |
409 |
C write(*,*) 'cg2d: Ax - b = ',SQRT(err) |
410 |
CcnhDebugEnds |
411 |
|
412 |
RETURN |
413 |
END |