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