/[escript]/branches/windows_from_1456_trunk_1574_merged_in/paso/src/BiCGStab.c
ViewVC logotype

Contents of /branches/windows_from_1456_trunk_1574_merged_in/paso/src/BiCGStab.c

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1577 - (show annotations)
Wed May 28 11:13:51 2008 UTC (11 years, 4 months ago) by phornby
File MIME type: text/plain
File size: 8264 byte(s)
After checking in the merged files, rename 1544 merge branch to a 1574 merged branch.


1
2 /* $Id$ */
3
4 /*******************************************************
5 *
6 * Copyright 2003-2007 by ACceSS MNRF
7 * Copyright 2007 by University of Queensland
8 *
9 * http://esscc.uq.edu.au
10 * Primary Business: Queensland, Australia
11 * Licensed under the Open Software License version 3.0
12 * http://www.opensource.org/licenses/osl-3.0.php
13 *
14 *******************************************************/
15
16 /*
17 Crude modifications and translations for Paso by Matt Davies and Lutz Gross
18 */
19
20 #include "Paso.h"
21 #include "SystemMatrix.h"
22 #include "Solver.h"
23 #ifdef _OPENMP
24 #include <omp.h>
25 #endif
26
27 /* -- Iterative template routine --
28 * Univ. of Tennessee and Oak Ridge National Laboratory
29 * October 1, 1993
30 * Details of this algorithm are described in "Templates for the
31 * Solution of Linear Systems: Building Blocks for Iterative
32 * Methods", Barrett, Berry, Chan, Demmel, Donato, Dongarra,
33 * Eijkhout, Pozo, Romine, and van der Vorst, SIAM Publications,
34 * 1993. (ftp netlib2.cs.utk.edu; cd linalg; get templates.ps).
35 *
36 * Purpose
37 * =======
38 *
39 * BICGSTAB solves the linear system A*x = b using the
40 * BiConjugate Gradient Stabilized iterative method with
41 * preconditioning.
42 *
43 * Convergence test: norm( b - A*x )< TOL.
44 * For other measures, see the above reference.
45 *
46 * Arguments
47 * =========
48 *
49 * A (input)
50 *
51 * R (input) DOUBLE PRECISION array, dimension N.
52 * On entry, residual of inital guess X
53 *
54 * X (input/output) DOUBLE PRECISION array, dimension N.
55 * On input, the initial guess.
56 *
57 * ITER (input/output) INT
58 * On input, the maximum iterations to be performed.
59 * On output, actual number of iterations performed.
60 *
61 * RESID (input/output) DOUBLE PRECISION
62 * On input, the allowable convergence measure for
63 * norm( b - A*x )
64 * On output, the final value of this measure.
65 *
66 * return value
67 *
68 * = SOLVER_NO_ERROR: Successful exit. Iterated approximate solution returned.
69 * = SOLVER_MAXITER_REACHED
70 * = SOLVER_INPUT_ERROR Illegal parameter:
71 * = SOLVER_BREAKDOWN: If parameters RHO or OMEGA become smaller
72 * = SOLVER_MEMORY_ERROR : If parameters RHO or OMEGA become smaller
73 *
74 * ==============================================================
75 */
76
77 err_t Paso_Solver_BiCGStab(
78 Paso_SystemMatrix * A,
79 double * r,
80 double * x,
81 dim_t *iter,
82 double * tolerance,
83 Paso_Performance* pp) {
84
85
86 /* Local variables */
87 double *rtld=NULL,*p=NULL,*v=NULL,*t=NULL,*phat=NULL,*shat=NULL,*s=NULL, *buf1=NULL, *buf0=NULL;
88 double beta,norm_of_residual,sum_1,sum_2,sum_3,sum_4,norm_of_residual_global;
89 double alpha, omega, omegaNumtr, omegaDenumtr, rho, tol, rho1, loc_sum[2], sum[2];
90 dim_t num_iter=0,maxit,num_iter_global;
91 dim_t i0;
92 bool_t breakFlag=FALSE, maxIterFlag=FALSE, convergeFlag=FALSE;
93 dim_t status = SOLVER_NO_ERROR;
94 double *resid = tolerance;
95 dim_t n = Paso_SystemMatrix_getTotalNumRows(A);
96
97 /* Executable Statements */
98
99 /* allocate memory: */
100 rtld=TMPMEMALLOC(n,double);
101 p=TMPMEMALLOC(n,double);
102 v=TMPMEMALLOC(n,double);
103 t=TMPMEMALLOC(n,double);
104 phat=TMPMEMALLOC(n,double);
105 shat=TMPMEMALLOC(n,double);
106 s=TMPMEMALLOC(n,double);
107 /* Test the input parameters. */
108
109 if (n < 0) {
110 status = SOLVER_INPUT_ERROR;
111 } else if (rtld==NULL || p==NULL || v==NULL || t==NULL || phat==NULL || shat==NULL || s==NULL) {
112 status = SOLVER_MEMORY_ERROR;
113 } else {
114
115 /* now bicgstab starts : */
116 maxit = *iter;
117 tol = *resid;
118
119 num_iter =0;
120 convergeFlag=FALSE;
121 maxIterFlag=FALSE;
122 breakFlag=FALSE;
123
124 /* initialize arrays */
125
126 #pragma omp parallel for private(i0) schedule(static)
127 for (i0 = 0; i0 < n; i0++) {
128 rtld[i0]=0;
129 p[i0]=0;
130 v[i0]=0;
131 t[i0]=0;
132 phat[i0]=0;
133 shat[i0]=0;
134 rtld[i0] = r[i0];
135 }
136
137 /* Perform BiConjugate Gradient Stabilized iteration. */
138
139 L10:
140 ++(num_iter);
141 sum_1 = 0;
142 sum_2 = 0;
143 sum_3 = 0;
144 sum_4 = 0;
145 omegaNumtr = 0.0;
146 omegaDenumtr = 0.0;
147 #pragma omp parallel for private(i0) reduction(+:sum_1) schedule(static)
148 for (i0 = 0; i0 < n; i0++) sum_1 += rtld[i0] * r[i0];
149 #ifdef PASO_MPI
150 loc_sum[0] = sum_1;
151 MPI_Allreduce(loc_sum, &sum_1, 1, MPI_DOUBLE, MPI_SUM, A->mpi_info->comm);
152 #endif
153 rho = sum_1;
154
155 if (! (breakFlag = (ABS(rho) <= TOLERANCE_FOR_SCALARS))) {
156 /* Compute vector P. */
157
158 if (num_iter > 1) {
159 beta = rho / rho1 * (alpha / omega);
160 #pragma omp parallel for private(i0) schedule(static)
161 for (i0 = 0; i0 < n; i0++) p[i0] = r[i0] + beta * (p[i0] - omega * v[i0]);
162 } else {
163 #pragma omp parallel for private(i0) schedule(static)
164 for (i0 = 0; i0 < n; i0++) p[i0] = r[i0];
165 }
166
167 /* Compute direction adjusting vector PHAT and scalar ALPHA. */
168
169 Paso_Solver_solvePreconditioner(A,&phat[0], &p[0]);
170 Paso_SystemMatrix_MatrixVector_CSR_OFFSET0(ONE, A, &phat[0],ZERO, &v[0]);
171
172 #pragma omp parallel for private(i0) reduction(+:sum_2) schedule(static)
173 for (i0 = 0; i0 < n; i0++) sum_2 += rtld[i0] * v[i0];
174 #ifdef PASO_MPI
175 loc_sum[0] = sum_2;
176 MPI_Allreduce(loc_sum, &sum_2, 1, MPI_DOUBLE, MPI_SUM, A->mpi_info->comm);
177 #endif
178 if (! (breakFlag = (ABS(sum_2) <= TOLERANCE_FOR_SCALARS))) {
179 alpha = rho / sum_2;
180
181 #pragma omp parallel for private(i0) reduction(+:sum_3) schedule(static)
182 for (i0 = 0; i0 < n; i0++) {
183 r[i0] -= alpha * v[i0];
184 s[i0] = r[i0];
185 sum_3 += s[i0] * s[i0];
186 }
187 #ifdef PASO_MPI
188 loc_sum[0] = sum_3;
189 MPI_Allreduce(loc_sum, &sum_3, 1, MPI_DOUBLE, MPI_SUM, A->mpi_info->comm);
190 #endif
191 norm_of_residual = sqrt(sum_3);
192
193 /* Early check for tolerance. */
194 if ( (convergeFlag = (norm_of_residual <= tol)) ) {
195 #pragma omp parallel for private(i0) schedule(static)
196 for (i0 = 0; i0 < n; i0++) x[i0] += alpha * phat[i0];
197 maxIterFlag = FALSE;
198 breakFlag = FALSE;
199 } else {
200 /* Compute stabilizer vector SHAT and scalar OMEGA. */
201 Paso_Solver_solvePreconditioner(A,&shat[0], &s[0]);
202 Paso_SystemMatrix_MatrixVector_CSR_OFFSET0(ONE, A, &shat[0],ZERO,&t[0]);
203
204 #pragma omp parallel for private(i0) reduction(+:omegaNumtr,omegaDenumtr) schedule(static)
205 for (i0 = 0; i0 < n; i0++) {
206 omegaNumtr +=t[i0] * s[i0];
207 omegaDenumtr += t[i0] * t[i0];
208 }
209 #ifdef PASO_MPI
210 loc_sum[0] = omegaNumtr;
211 loc_sum[1] = omegaDenumtr;
212 MPI_Allreduce(loc_sum, sum, 2, MPI_DOUBLE, MPI_SUM, A->mpi_info->comm);
213 omegaNumtr=sum[0];
214 omegaDenumtr=sum[1];
215 #endif
216 if (! (breakFlag = (ABS(omegaDenumtr) <= TOLERANCE_FOR_SCALARS))) {
217 omega = omegaNumtr / omegaDenumtr;
218
219 #pragma omp parallel for private(i0) reduction(+:sum_4) schedule(static)
220 for (i0 = 0; i0 < n; i0++) {
221 x[i0] += alpha * phat[i0] + omega * shat[i0];
222 r[i0] = s[i0]-omega * t[i0];
223 sum_4 += r[i0] * r[i0];
224 }
225 #ifdef PASO_MPI
226 loc_sum[0] = sum_4;
227 MPI_Allreduce(loc_sum, &sum_4, 1, MPI_DOUBLE, MPI_SUM, A->mpi_info->comm);
228 #endif
229 norm_of_residual = sqrt(sum_4);
230 convergeFlag = norm_of_residual <= tol;
231 maxIterFlag = num_iter == maxit;
232 breakFlag = (ABS(omega) <= TOLERANCE_FOR_SCALARS);
233 }
234 }
235 }
236 if (!(convergeFlag || maxIterFlag || breakFlag)) {
237 rho1 = rho;
238 goto L10;
239 }
240 }
241 /* end of iteration */
242 num_iter_global=num_iter;
243 norm_of_residual_global=norm_of_residual;
244 if (maxIterFlag) {
245 status = SOLVER_MAXITER_REACHED;
246 } else if (breakFlag) {
247 status = SOLVER_BREAKDOWN;
248 }
249 }
250 TMPMEMFREE(rtld);
251 TMPMEMFREE(p);
252 TMPMEMFREE(v);
253 TMPMEMFREE(t);
254 TMPMEMFREE(phat);
255 TMPMEMFREE(shat);
256 TMPMEMFREE(s);
257 *iter=num_iter_global;
258 *resid=norm_of_residual_global;
259
260 /* End of BICGSTAB */
261 return status;
262 }

Properties

Name Value
svn:eol-style native
svn:keywords Author Date Id Revision

  ViewVC Help
Powered by ViewVC 1.1.26