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