/[MITgcm]/MITgcm/model/src/cg2d_nsa.F
ViewVC logotype

Contents of /MITgcm/model/src/cg2d_nsa.F

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


Revision 1.6 - (show annotations) (download)
Sun Aug 12 20:24:23 2012 UTC (11 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint63r, checkpoint63s, checkpoint64
Changes since 1.5: +2 -1 lines
add PACKAGES_CONFIG.h wherever ALLOW_AUTODIFF[_TAMC] is used (in model/src)

1 C $Header: /u/gcmpack/MITgcm/model/src/cg2d_nsa.F,v 1.5 2012/07/17 21:31:52 jmc Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 CML THIS DOES NOT WORK +++++
8 #undef ALLOW_LOOP_DIRECTIVE
9 CBOP
10 C !ROUTINE: CG2D_NSA
11 C !INTERFACE:
12 SUBROUTINE CG2D_NSA(
13 U cg2d_b, cg2d_x,
14 O firstResidual, minResidualSq, lastResidual,
15 U numIters, nIterMin,
16 I myThid )
17 C !DESCRIPTION: \bv
18 C *==========================================================*
19 C | SUBROUTINE CG2D_NSA
20 C | o Two-dimensional grid problem conjugate-gradient
21 C | inverter (with preconditioner).
22 C | o This version is used only in the case when the matrix
23 C | operator is not "self-adjoint" (NSA). Any remaining
24 C | residuals will immediately reported to the department
25 C | of homeland security.
26 C *==========================================================*
27 C | Con. grad is an iterative procedure for solving Ax = b.
28 C | It requires the A be symmetric.
29 C | This implementation assumes A is a five-diagonal
30 C | matrix of the form that arises in the discrete
31 C | representation of the del^2 operator in a
32 C | two-dimensional space.
33 C | Notes:
34 C | ======
35 C | This implementation can support shared-memory
36 C | multi-threaded execution. In order to do this COMMON
37 C | blocks are used for many of the arrays - even ones that
38 C | are only used for intermedaite results. This design is
39 C | OK if you want to all the threads to collaborate on
40 C | solving the same problem. On the other hand if you want
41 C | the threads to solve several different problems
42 C | concurrently this implementation will not work.
43 C *==========================================================*
44 C \ev
45
46 C !USES:
47 IMPLICIT NONE
48 C === Global data ===
49 #include "SIZE.h"
50 #include "EEPARAMS.h"
51 #include "PARAMS.h"
52 #include "CG2D.h"
53 #ifdef ALLOW_AUTODIFF_TAMC
54 # include "tamc.h"
55 # include "tamc_keys.h"
56 #endif
57
58 C !INPUT/OUTPUT PARAMETERS:
59 C === Routine arguments ===
60 C cg2d_b :: The source term or "right hand side" (Output: normalised RHS)
61 C cg2d_x :: The solution (Input: first guess)
62 C firstResidual :: the initial residual before any iterations
63 C minResidualSq :: the lowest residual reached (squared)
64 C lastResidual :: the actual residual reached
65 C numIters :: Inp: the maximum number of iterations allowed
66 C Out: the actual number of iterations used
67 C nIterMin :: Inp: decide to store (if >=0) or not (if <0) lowest res. sol.
68 C Out: iteration number corresponding to lowest residual
69 C myThid :: Thread on which I am working.
70 _RL cg2d_b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
71 _RL cg2d_x(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
72 _RL firstResidual
73 _RL minResidualSq
74 _RL lastResidual
75 INTEGER numIters
76 INTEGER nIterMin
77 INTEGER myThid
78
79 #ifdef ALLOW_CG2D_NSA
80 C !LOCAL VARIABLES:
81 C === Local variables ====
82 C bi, bj :: tile index in X and Y.
83 C i, j, it2d :: Loop counters ( it2d counts CG iterations )
84 C actualIts :: actual CG iteration number
85 C err_sq :: Measure of the square of the residual of Ax - b.
86 C eta_qrN :: Used in computing search directions; suffix N and NM1
87 C eta_qrNM1 denote current and previous iterations respectively.
88 C recip_eta_qrNM1 :: reciprocal of eta_qrNM1
89 C cgBeta :: coeff used to update conjugate direction vector "s".
90 C alpha :: coeff used to update solution & residual
91 C alphaSum :: to avoid the statement: alpha = 1./alpha (for TAMC/TAF)
92 C sumRHS :: Sum of right-hand-side. Sometimes this is a useful
93 C debugging/trouble shooting diagnostic. For neumann problems
94 C sumRHS needs to be ~0 or it converge at a non-zero residual.
95 C cg2d_min :: used to store solution corresponding to lowest residual.
96 C msgBuf :: Informational/error message buffer
97 INTEGER bi, bj
98 INTEGER i, j, it2d
99 INTEGER actualIts
100 _RL cg2dTolerance_sq
101 _RL err_sq, errTile(nSx,nSy)
102 _RL eta_qrN, eta_qrNtile(nSx,nSy)
103 _RL eta_qrNM1, recip_eta_qrNM1
104 _RL cgBeta, alpha
105 _RL alphaSum,alphaTile(nSx,nSy)
106 _RL sumRHS, sumRHStile(nSx,nSy)
107 _RL rhsMax, rhsNorm
108 CHARACTER*(MAX_LEN_MBUF) msgBuf
109 LOGICAL printResidual
110 CEOP
111
112 #ifdef ALLOW_AUTODIFF_TAMC
113 IF ( numIters .GT. numItersMax ) THEN
114 WRITE(msgBuf,'(A,I10)')
115 & 'CG2D_NSA: numIters > numItersMax =', numItersMax
116 CALL PRINT_ERROR( msgBuf, myThid )
117 STOP 'ABNORMAL END: S/R CG2D_NSA'
118 ENDIF
119 IF ( cg2dNormaliseRHS ) THEN
120 WRITE(msgBuf,'(A)') 'CG2D_NSA: cg2dNormaliseRHS is disabled'
121 CALL PRINT_ERROR( msgBuf, myThid )
122 WRITE(msgBuf,'(A)')
123 & 'set cg2dTargetResWunit (instead of cg2dTargetResidual)'
124 CALL PRINT_ERROR( msgBuf, myThid )
125 STOP 'ABNORMAL END: S/R CG2D_NSA'
126 ENDIF
127 #endif /* ALLOW_AUTODIFF_TAMC */
128
129 #ifdef ALLOW_AUTODIFF_TAMC
130 act1 = myThid - 1
131 max1 = nTx*nTy
132 act2 = ikey_dynamics - 1
133 ikey = (act1 + 1) + act2*max1
134 #endif /* ALLOW_AUTODIFF_TAMC */
135
136 C-- Initialise auxiliary constant, some output variable and inverter
137 cg2dTolerance_sq = cg2dTolerance*cg2dTolerance
138 minResidualSq = -1. _d 0
139 eta_qrNM1 = 1. _d 0
140 recip_eta_qrNM1= 1. _d 0
141
142 C-- Normalise RHS
143 rhsMax = 0. _d 0
144 DO bj=myByLo(myThid),myByHi(myThid)
145 DO bi=myBxLo(myThid),myBxHi(myThid)
146 DO j=1,sNy
147 DO i=1,sNx
148 cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj)*cg2dNorm
149 rhsMax = MAX(ABS(cg2d_b(i,j,bi,bj)),rhsMax)
150 ENDDO
151 ENDDO
152 ENDDO
153 ENDDO
154
155 #ifndef ALLOW_AUTODIFF_TAMC
156 IF (cg2dNormaliseRHS) THEN
157 C- Normalise RHS :
158 _GLOBAL_MAX_RL( rhsMax, myThid )
159 rhsNorm = 1. _d 0
160 IF ( rhsMax .NE. 0. ) rhsNorm = 1. _d 0 / rhsMax
161 DO bj=myByLo(myThid),myByHi(myThid)
162 DO bi=myBxLo(myThid),myBxHi(myThid)
163 DO j=1,sNy
164 DO i=1,sNx
165 cg2d_b(i,j,bi,bj) = cg2d_b(i,j,bi,bj)*rhsNorm
166 cg2d_x(i,j,bi,bj) = cg2d_x(i,j,bi,bj)*rhsNorm
167 ENDDO
168 ENDDO
169 ENDDO
170 ENDDO
171 C- end Normalise RHS
172 ENDIF
173 #endif /* ndef ALLOW_AUTODIFF_TAMC */
174
175 C-- Update overlaps
176 CALL EXCH_XY_RL( cg2d_x, myThid )
177
178 C-- Initial residual calculation
179 #ifdef ALLOW_AUTODIFF_TAMC
180 CADJ STORE cg2d_b = comlev1_cg2d, key = ikey, byte = isbyte
181 CADJ STORE cg2d_x = comlev1_cg2d, key = ikey, byte = isbyte
182 #endif /* ALLOW_AUTODIFF_TAMC */
183 DO bj=myByLo(myThid),myByHi(myThid)
184 DO bi=myBxLo(myThid),myBxHi(myThid)
185 errTile(bi,bj) = 0. _d 0
186 sumRHStile(bi,bj) = 0. _d 0
187 DO j=0,sNy+1
188 DO i=0,sNx+1
189 cg2d_s(i,j,bi,bj) = 0.
190 ENDDO
191 ENDDO
192 DO j=1,sNy
193 DO i=1,sNx
194 cg2d_r(i,j,bi,bj) = cg2d_b(i,j,bi,bj) -
195 & (aW2d(i ,j ,bi,bj)*cg2d_x(i-1,j ,bi,bj)
196 & +aW2d(i+1,j ,bi,bj)*cg2d_x(i+1,j ,bi,bj)
197 & +aS2d(i ,j ,bi,bj)*cg2d_x(i ,j-1,bi,bj)
198 & +aS2d(i ,j+1,bi,bj)*cg2d_x(i ,j+1,bi,bj)
199 & +aC2d(i ,j ,bi,bj)*cg2d_x(i ,j ,bi,bj)
200 & )
201 errTile(bi,bj) = errTile(bi,bj)
202 & + cg2d_r(i,j,bi,bj)*cg2d_r(i,j,bi,bj)
203 sumRHStile(bi,bj) = sumRHStile(bi,bj) + cg2d_b(i,j,bi,bj)
204 ENDDO
205 ENDDO
206 ENDDO
207 ENDDO
208
209 c CALL EXCH_S3D_RL( cg2d_r, 1, myThid )
210 CALL EXCH_XY_RL ( cg2d_r, myThid )
211 CALL GLOBAL_SUM_TILE_RL( errTile, err_sq, myThid )
212 CALL GLOBAL_SUM_TILE_RL( sumRHStile, sumRHS, myThid )
213 actualIts = 0
214 IF ( err_sq .NE. 0. ) THEN
215 firstResidual = SQRT(err_sq)
216 ELSE
217 firstResidual = 0.
218 ENDIF
219
220 printResidual = .FALSE.
221 IF ( debugLevel .GE. debLevZero ) THEN
222 _BEGIN_MASTER( myThid )
223 printResidual = printResidualFreq.GE.1
224 WRITE(standardmessageunit,'(A,1P2E22.14)')
225 & ' cg2d: Sum(rhs),rhsMax = ', sumRHS,rhsMax
226 _END_MASTER( myThid )
227 ENDIF
228
229 C >>>>>>>>>>>>>>> BEGIN SOLVER <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
230 Cml begin main solver loop
231 #if ((defined ALLOW_AUTODIFF_TAMC) && (defined ALLOW_LOOP_DIRECTIVE))
232 CADJ LOOP = iteration, cg2d_x = comlev_cg2d_iter
233 #endif /* ALLOW_AUTODIFF_TAMC and ALLOW_LOOP_DIRECTIVE */
234 DO it2d=1, numIters
235 #ifdef ALLOW_LOOP_DIRECTIVE
236 CML it2d = 0
237 CML DO WHILE ( err_sq .GT. cg2dTolerance_sq .and. it2d .LT. numIters )
238 CML it2d = it2d+1
239 #endif /* ALLOW_LOOP_DIRECTIVE */
240
241 #ifdef ALLOW_AUTODIFF_TAMC
242 icg2dkey = (ikey-1)*numItersMax + it2d
243 CADJ STORE err_sq = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
244 #endif /* ALLOW_AUTODIFF_TAMC */
245 IF ( err_sq .GE. cg2dTolerance_sq ) THEN
246
247 C-- Solve preconditioning equation and update
248 C-- conjugate direction vector "s".
249 #ifdef ALLOW_AUTODIFF_TAMC
250 CADJ STORE cg2d_r = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
251 CADJ STORE cg2d_s = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
252 #endif /* ALLOW_AUTODIFF_TAMC */
253 DO bj=myByLo(myThid),myByHi(myThid)
254 DO bi=myBxLo(myThid),myBxHi(myThid)
255 eta_qrNtile(bi,bj) = 0. _d 0
256 DO j=1,sNy
257 DO i=1,sNx
258 cg2d_z(i,j,bi,bj) =
259 & pC(i ,j ,bi,bj)*cg2d_r(i ,j ,bi,bj)
260 & +pW(i ,j ,bi,bj)*cg2d_r(i-1,j ,bi,bj)
261 & +pW(i+1,j ,bi,bj)*cg2d_r(i+1,j ,bi,bj)
262 & +pS(i ,j ,bi,bj)*cg2d_r(i ,j-1,bi,bj)
263 & +pS(i ,j+1,bi,bj)*cg2d_r(i ,j+1,bi,bj)
264 eta_qrNtile(bi,bj) = eta_qrNtile(bi,bj)
265 & +cg2d_z(i,j,bi,bj)*cg2d_r(i,j,bi,bj)
266 ENDDO
267 ENDDO
268 ENDDO
269 ENDDO
270
271 CALL GLOBAL_SUM_TILE_RL( eta_qrNtile,eta_qrN,myThid )
272 #ifdef ALLOW_AUTODIFF_TAMC
273 CMLCADJ STORE eta_qrNM1 = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
274 CADJ STORE recip_eta_qrNM1 = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
275 #endif /* ALLOW_AUTODIFF_TAMC */
276 CML cgBeta = eta_qrN/eta_qrNM1
277 cgBeta = eta_qrN*recip_eta_qrNM1
278 Cml store normalisation factor for the next interation (in case there is one).
279 CML store the inverse of the normalization factor for higher precision
280 CML eta_qrNM1 = eta_qrN
281 recip_eta_qrNM1 = 1. _d 0/eta_qrN
282
283 DO bj=myByLo(myThid),myByHi(myThid)
284 DO bi=myBxLo(myThid),myBxHi(myThid)
285 DO j=1,sNy
286 DO i=1,sNx
287 cg2d_s(i,j,bi,bj) = cg2d_z(i,j,bi,bj)
288 & + cgBeta*cg2d_s(i,j,bi,bj)
289 ENDDO
290 ENDDO
291 ENDDO
292 ENDDO
293
294 C-- Do exchanges that require messages i.e. between
295 C-- processes.
296 c CALL EXCH_S3D_RL( cg2d_s, 1, myThid )
297 CALL EXCH_XY_RL ( cg2d_s, myThid )
298
299 C== Evaluate laplace operator on conjugate gradient vector
300 C== q = A.s
301 #ifdef ALLOW_AUTODIFF_TAMC
302 #ifndef ALLOW_LOOP_DIRECTIVE
303 CADJ STORE cg2d_s = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
304 #endif /* not ALLOW_LOOP_DIRECTIVE */
305 #endif /* ALLOW_AUTODIFF_TAMC */
306 DO bj=myByLo(myThid),myByHi(myThid)
307 DO bi=myBxLo(myThid),myBxHi(myThid)
308 alphaTile(bi,bj) = 0. _d 0
309 DO j=1,sNy
310 DO i=1,sNx
311 cg2d_q(i,j,bi,bj) =
312 & aW2d(i ,j ,bi,bj)*cg2d_s(i-1,j ,bi,bj)
313 & +aW2d(i+1,j ,bi,bj)*cg2d_s(i+1,j ,bi,bj)
314 & +aS2d(i ,j ,bi,bj)*cg2d_s(i ,j-1,bi,bj)
315 & +aS2d(i ,j+1,bi,bj)*cg2d_s(i ,j+1,bi,bj)
316 & +aC2d(i ,j ,bi,bj)*cg2d_s(i ,j ,bi,bj)
317 alphaTile(bi,bj) = alphaTile(bi,bj)
318 & + cg2d_s(i,j,bi,bj)*cg2d_q(i,j,bi,bj)
319 ENDDO
320 ENDDO
321 ENDDO
322 ENDDO
323 CALL GLOBAL_SUM_TILE_RL( alphaTile, alphaSum, myThid )
324 alpha = eta_qrN/alphaSum
325
326 C== Update simultaneously solution and residual vectors (and Iter number)
327 C Now compute "interior" points.
328 #ifdef ALLOW_AUTODIFF_TAMC
329 #ifndef ALLOW_LOOP_DIRECTIVE
330 CADJ STORE cg2d_r = comlev1_cg2d_iter, key = icg2dkey, byte = isbyte
331 #endif /* ALLOW_LOOP_DIRECTIVE */
332 #endif /* ALLOW_AUTODIFF_TAMC */
333 DO bj=myByLo(myThid),myByHi(myThid)
334 DO bi=myBxLo(myThid),myBxHi(myThid)
335 errTile(bi,bj) = 0. _d 0
336 DO j=1,sNy
337 DO i=1,sNx
338 cg2d_x(i,j,bi,bj)=cg2d_x(i,j,bi,bj)+alpha*cg2d_s(i,j,bi,bj)
339 cg2d_r(i,j,bi,bj)=cg2d_r(i,j,bi,bj)-alpha*cg2d_q(i,j,bi,bj)
340 errTile(bi,bj) = errTile(bi,bj)
341 & + cg2d_r(i,j,bi,bj)*cg2d_r(i,j,bi,bj)
342 ENDDO
343 ENDDO
344 ENDDO
345 ENDDO
346 actualIts = it2d
347
348 CALL GLOBAL_SUM_TILE_RL( errTile, err_sq, myThid )
349 IF ( printResidual ) THEN
350 IF ( MOD( it2d-1, printResidualFreq ).EQ.0 ) THEN
351 WRITE(msgBuf,'(A,I6,A,1PE21.14)')
352 & ' cg2d: iter=', it2d, ' ; resid.= ', SQRT(err_sq)
353 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
354 & SQUEEZE_RIGHT, myThid )
355 ENDIF
356 ENDIF
357
358 c CALL EXCH_S3D_RL( cg2d_r, 1, myThid )
359 CALL EXCH_XY_RL ( cg2d_r, myThid )
360
361 Cml end of if "err >= cg2dTolerance" block ; end main solver loop
362 ENDIF
363 ENDDO
364
365 #ifndef ALLOW_AUTODIFF_TAMC
366 IF (cg2dNormaliseRHS) THEN
367 C-- Un-normalise the answer
368 DO bj=myByLo(myThid),myByHi(myThid)
369 DO bi=myBxLo(myThid),myBxHi(myThid)
370 DO j=1,sNy
371 DO i=1,sNx
372 cg2d_x(i,j,bi,bj) = cg2d_x(i,j,bi,bj)/rhsNorm
373 ENDDO
374 ENDDO
375 ENDDO
376 ENDDO
377 ENDIF
378 #endif /* ndef ALLOW_AUTODIFF_TAMC */
379
380 C-- Return parameters to caller
381 IF ( err_sq .NE. 0. ) THEN
382 lastResidual = SQRT(err_sq)
383 ELSE
384 lastResidual = 0.
385 ENDIF
386 numIters = actualIts
387
388 #endif /* ALLOW_CG2D_NSA */
389 RETURN
390 END
391
392 #if ((defined ALLOW_AUTODIFF_TAMC) && (defined ALLOW_LOOP_DIRECTIVE))
393
394 C These routines are routinely part of the TAMC/TAF library that is
395 C not included in the MITcgm, therefore they are mimicked here.
396
397 subroutine adstore(chardum,int1,idow,int2,int3,icount)
398
399 implicit none
400
401 #include "SIZE.h"
402 #include "tamc.h"
403
404 character*(*) chardum
405 integer int1, int2, int3, idow, icount
406
407 C the length of this vector must be greater or equal
408 C twice the number of timesteps
409 integer nidow
410 #ifdef ALLOW_TAMC_CHECKPOINTING
411 parameter ( nidow = 2*nchklev_1*nchklev_2*nchklev_3 )
412 #else
413 parameter ( nidow = 1000000 )
414 #endif /* ALLOW_TAMC_CHECKPOINTING */
415 integer istoreidow(nidow)
416 common /istorecommon/ istoreidow
417
418 print *, 'adstore: ', chardum, int1, idow, int2, int3, icount
419
420 if ( icount .gt. nidow ) then
421 print *, 'adstore: error: icount > nidow = ', nidow
422 stop 'ABNORMAL STOP in adstore'
423 endif
424
425 istoreidow(icount) = idow
426
427 return
428 end
429
430 subroutine adresto(chardum,int1,idow,int2,int3,icount)
431
432 implicit none
433
434 #include "SIZE.h"
435 #include "tamc.h"
436
437 character*(*) chardum
438 integer int1, int2, int3, idow, icount
439
440 C the length of this vector must be greater or equal
441 C twice the number of timesteps
442 integer nidow
443 #ifdef ALLOW_TAMC_CHECKPOINTING
444 parameter ( nidow = 2*nchklev_1*nchklev_2*nchklev_3 )
445 #else
446 parameter ( nidow = 1000000 )
447 #endif /* ALLOW_TAMC_CHECKPOINTING */
448 integer istoreidow(nidow)
449 common /istorecommon/ istoreidow
450
451 print *, 'adresto: ', chardum, int1, idow, int2, int3, icount
452
453 if ( icount .gt. nidow ) then
454 print *, 'adstore: error: icount > nidow = ', nidow
455 stop 'ABNORMAL STOP in adstore'
456 endif
457
458 idow = istoreidow(icount)
459
460 return
461 end
462 #endif /* ALLOW_AUTODIFF_TAMC and ALLOW_LOOP_DIRECTIVE */

  ViewVC Help
Powered by ViewVC 1.1.22