53 struct d_solver solver;
61 double *x = (
double*)(xp->
b);
62 double *b = (
double*)(bp->
b);
75 double normb, normx, normr, alpha, beta;
78 memset( &solver, 0,
sizeof(
struct d_solver) );
81 if ( !(pastix_data->
steps & STEP_NUMFACT) ) {
85 n = pastix_data->
bcsc->n;
89 gradr = (
double *)solver.malloc(n *
sizeof(
double));
90 gradr2 = (
double *)solver.malloc(n *
sizeof(
double));
91 gradp = (
double *)solver.malloc(n *
sizeof(
double));
92 grady = (
double *)solver.malloc(n *
sizeof(
double));
93 gradv = (
double *)solver.malloc(n *
sizeof(
double));
94 grads = (
double *)solver.malloc(n *
sizeof(
double));
95 gradz = (
double *)solver.malloc(n *
sizeof(
double));
96 gradt = (
double *)solver.malloc(n *
sizeof(
double));
97 grad2 = (
double *)solver.malloc(n *
sizeof(
double));
98 grad3 = (
double *)solver.malloc(n *
sizeof(
double));
103 sgrad = solver.malloc( n *
sizeof(
float) );
106 clockInit(refine_clk);clockStart(refine_clk);
108 normb = solver.norm( pastix_data, n, b );
112 normx = solver.norm( pastix_data, n, x );
115 solver.copy( pastix_data, n, b, gradr );
117 solver.spmv( pastix_data,
PastixNoTrans, -1., x, 1., gradr );
119 normr = solver.norm( pastix_data, n, gradr );
122 solver.copy( pastix_data, n, gradr, gradr2 );
124 solver.copy( pastix_data, n, gradr, gradp );
127 resid_b = normr / normb;
129 while ((resid_b > eps) && (nb_iter < itermax))
131 clockStop((refine_clk));
136 solver.copy( pastix_data, n, gradp, grady );
138 solver.spsv( pastix_data, grady, sgrad );
142 solver.spmv( pastix_data,
PastixNoTrans, 1.0, grady, 0., gradv );
145 alpha = solver.dot( pastix_data, n, gradv, gradr2 );
146 beta = solver.dot( pastix_data, n, gradr, gradr2 );
147 alpha = beta / alpha;
150 solver.copy( pastix_data, n, gradr, grads );
151 solver.axpy( pastix_data, n, -alpha, gradv, grads );
154 solver.copy( pastix_data, n, grads, gradz );
156 solver.spsv( pastix_data, gradz, sgrad );
160 solver.spmv( pastix_data,
PastixNoTrans, 1.0, gradz, 0., gradt );
164 solver.copy( pastix_data, n, gradt, grad2 );
166 solver.spsv( pastix_data, grad2, sgrad );
171 v1 = solver.dot( pastix_data, n, grad2, gradz );
172 v2 = solver.dot( pastix_data, n, grad2, grad2 );
177 solver.axpy( pastix_data, n, alpha, grady, x );
180 solver.axpy( pastix_data, n, w, gradz, x );
183 solver.copy( pastix_data, n, grads, gradr );
184 solver.axpy( pastix_data, n, -w, gradt, gradr );
188 v1 = solver.dot( pastix_data, n, gradr, gradr2 );
196 solver.axpy( pastix_data, n, -w, gradv, gradp );
199 solver.scal( pastix_data, n, beta, gradp );
200 solver.axpy( pastix_data, n, 1., gradr, gradp );
202 normr = solver.norm( pastix_data, n, gradr );
203 resid_b = normr / normb;
205 clockStop((refine_clk));
208 ( pastix_data->
procnum == 0 ) ) {
209 solver.output_oneiter( t0, t3, resid_b, nb_iter );
213 solver.output_final(pastix_data, resid_b, nb_iter, t3, x, x);
215 solver.free((
void*) gradr);
216 solver.free((
void*) gradr2);
217 solver.free((
void*) gradp);
218 solver.free((
void*) grady);
219 solver.free((
void*) gradv);
220 solver.free((
void*) grads);
221 solver.free((
void*) gradz);
222 solver.free((
void*) gradt);
223 solver.free((
void*) grad2);
224 solver.free((
void*) grad3);
225 solver.free((
void*) sgrad);
BEGIN_C_DECLS typedef int pastix_int_t
@ DPARM_EPSILON_REFINEMENT
void d_refine_init(struct d_solver *, pastix_data_t *)
Initiate functions pointers to define basic operations.
pastix_int_t d_bicgstab_smp(pastix_data_t *pastix_data, pastix_rhs_t xp, pastix_rhs_t bp)
Main PaStiX data structure.
Main PaStiX RHS structure.