1 /******************************************************************************
2  * Copyright 1998-2019 Lawrence Livermore National Security, LLC and other
3  * HYPRE Project Developers. See the top-level COPYRIGHT file for details.
4  *
5  * SPDX-License-Identifier: (Apache-2.0 OR MIT)
6  ******************************************************************************/
7 
8 /******************************************************************************
9  *
10  * BiCGSTAB bicgstab
11  *
12  *****************************************************************************/
13 
14 #include "krylov.h"
15 #include "_hypre_utilities.h"
16 
17 /*--------------------------------------------------------------------------
18  * hypre_BiCGSTABFunctionsCreate
19  *--------------------------------------------------------------------------*/
20 
21 hypre_BiCGSTABFunctions *
hypre_BiCGSTABFunctionsCreate(void * (* CreateVector)(void * vvector),HYPRE_Int (* DestroyVector)(void * vvector),void * (* MatvecCreate)(void * A,void * x),HYPRE_Int (* Matvec)(void * matvec_data,HYPRE_Complex alpha,void * A,void * x,HYPRE_Complex beta,void * y),HYPRE_Int (* MatvecDestroy)(void * matvec_data),HYPRE_Real (* InnerProd)(void * x,void * y),HYPRE_Int (* CopyVector)(void * x,void * y),HYPRE_Int (* ClearVector)(void * x),HYPRE_Int (* ScaleVector)(HYPRE_Complex alpha,void * x),HYPRE_Int (* Axpy)(HYPRE_Complex alpha,void * x,void * y),HYPRE_Int (* CommInfo)(void * A,HYPRE_Int * my_id,HYPRE_Int * num_procs),HYPRE_Int (* PrecondSetup)(void * vdata,void * A,void * b,void * x),HYPRE_Int (* Precond)(void * vdata,void * A,void * b,void * x))22 hypre_BiCGSTABFunctionsCreate(
23    void *     (*CreateVector)  ( void *vvector ),
24    HYPRE_Int  (*DestroyVector) ( void *vvector ),
25    void *     (*MatvecCreate)  ( void *A , void *x ),
26    HYPRE_Int  (*Matvec)        ( void *matvec_data , HYPRE_Complex alpha , void *A ,
27                                  void *x , HYPRE_Complex beta , void *y ),
28    HYPRE_Int  (*MatvecDestroy) ( void *matvec_data ),
29    HYPRE_Real (*InnerProd)     ( void *x , void *y ),
30    HYPRE_Int  (*CopyVector)    ( void *x , void *y ),
31    HYPRE_Int  (*ClearVector)   ( void *x ),
32    HYPRE_Int  (*ScaleVector)   ( HYPRE_Complex alpha , void *x ),
33    HYPRE_Int  (*Axpy)          ( HYPRE_Complex alpha , void *x , void *y ),
34    HYPRE_Int  (*CommInfo)      ( void *A , HYPRE_Int *my_id ,
35                                  HYPRE_Int *num_procs ),
36    HYPRE_Int  (*PrecondSetup)  ( void *vdata, void *A, void *b, void *x ),
37    HYPRE_Int  (*Precond)       ( void *vdata, void *A, void *b, void *x )
38    )
39 {
40    hypre_BiCGSTABFunctions * bicgstab_functions;
41    bicgstab_functions = (hypre_BiCGSTABFunctions *)
42       hypre_CTAlloc( hypre_BiCGSTABFunctions,  1 , HYPRE_MEMORY_HOST);
43 
44    bicgstab_functions->CreateVector = CreateVector;
45    bicgstab_functions->DestroyVector = DestroyVector;
46    bicgstab_functions->MatvecCreate = MatvecCreate;
47    bicgstab_functions->Matvec = Matvec;
48    bicgstab_functions->MatvecDestroy = MatvecDestroy;
49    bicgstab_functions->InnerProd = InnerProd;
50    bicgstab_functions->CopyVector = CopyVector;
51    bicgstab_functions->ClearVector = ClearVector;
52    bicgstab_functions->ScaleVector = ScaleVector;
53    bicgstab_functions->Axpy = Axpy;
54    bicgstab_functions->CommInfo = CommInfo;
55    bicgstab_functions->precond_setup = PrecondSetup;
56    bicgstab_functions->precond = Precond;
57 
58    return bicgstab_functions;
59 }
60 
61 /*--------------------------------------------------------------------------
62  * hypre_BiCGSTABCreate
63  *--------------------------------------------------------------------------*/
64 
65 void *
hypre_BiCGSTABCreate(hypre_BiCGSTABFunctions * bicgstab_functions)66 hypre_BiCGSTABCreate( hypre_BiCGSTABFunctions * bicgstab_functions )
67 {
68    hypre_BiCGSTABData *bicgstab_data;
69 
70    HYPRE_ANNOTATE_FUNC_BEGIN;
71 
72    bicgstab_data = hypre_CTAlloc( hypre_BiCGSTABData,  1, HYPRE_MEMORY_HOST);
73    bicgstab_data->functions = bicgstab_functions;
74 
75    /* set defaults */
76    (bicgstab_data -> tol)            = 1.0e-06;
77    (bicgstab_data -> min_iter)       = 0;
78    (bicgstab_data -> max_iter)       = 1000;
79    (bicgstab_data -> stop_crit)      = 0; /* rel. residual norm */
80    (bicgstab_data -> a_tol)          = 0.0;
81    (bicgstab_data -> precond_data)   = NULL;
82    (bicgstab_data -> logging)        = 0;
83    (bicgstab_data -> print_level)    = 0;
84    (bicgstab_data -> hybrid)         = 0;
85    (bicgstab_data -> p)              = NULL;
86    (bicgstab_data -> q)              = NULL;
87    (bicgstab_data -> r)              = NULL;
88    (bicgstab_data -> r0)             = NULL;
89    (bicgstab_data -> s)              = NULL;
90    (bicgstab_data -> v)             = NULL;
91    (bicgstab_data -> matvec_data)    = NULL;
92    (bicgstab_data -> norms)          = NULL;
93    (bicgstab_data -> log_file_name)  = NULL;
94 
95    HYPRE_ANNOTATE_FUNC_END;
96 
97    return (void *) bicgstab_data;
98 }
99 
100 /*--------------------------------------------------------------------------
101  * hypre_BiCGSTABDestroy
102  *--------------------------------------------------------------------------*/
103 
104 HYPRE_Int
hypre_BiCGSTABDestroy(void * bicgstab_vdata)105 hypre_BiCGSTABDestroy( void *bicgstab_vdata )
106 {
107    hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData *)bicgstab_vdata;
108 
109    HYPRE_ANNOTATE_FUNC_BEGIN;
110 
111    if (bicgstab_data)
112    {
113       hypre_BiCGSTABFunctions *bicgstab_functions = bicgstab_data->functions;
114       if ( (bicgstab_data -> norms) != NULL )
115             hypre_TFree(bicgstab_data -> norms, HYPRE_MEMORY_HOST);
116 
117       (*(bicgstab_functions->MatvecDestroy))(bicgstab_data -> matvec_data);
118 
119       (*(bicgstab_functions->DestroyVector))(bicgstab_data -> r);
120       (*(bicgstab_functions->DestroyVector))(bicgstab_data -> r0);
121       (*(bicgstab_functions->DestroyVector))(bicgstab_data -> s);
122       (*(bicgstab_functions->DestroyVector))(bicgstab_data -> v);
123       (*(bicgstab_functions->DestroyVector))(bicgstab_data -> p);
124       (*(bicgstab_functions->DestroyVector))(bicgstab_data -> q);
125 
126       hypre_TFree(bicgstab_data, HYPRE_MEMORY_HOST);
127       hypre_TFree(bicgstab_functions, HYPRE_MEMORY_HOST);
128    }
129 
130    HYPRE_ANNOTATE_FUNC_END;
131 
132    return(hypre_error_flag);
133 }
134 
135 /*--------------------------------------------------------------------------
136  * hypre_BiCGSTABSetup
137  *--------------------------------------------------------------------------*/
138 
139 HYPRE_Int
hypre_BiCGSTABSetup(void * bicgstab_vdata,void * A,void * b,void * x)140 hypre_BiCGSTABSetup( void *bicgstab_vdata,
141                   void *A,
142                   void *b,
143                   void *x         )
144 {
145    hypre_BiCGSTABData      *bicgstab_data      = (hypre_BiCGSTABData *)bicgstab_vdata;
146    hypre_BiCGSTABFunctions *bicgstab_functions = bicgstab_data->functions;
147 
148    HYPRE_Int            max_iter         = (bicgstab_data -> max_iter);
149    HYPRE_Int          (*precond_setup)(void*,void*,void*,void*) = (bicgstab_functions -> precond_setup);
150    void          *precond_data     = (bicgstab_data -> precond_data);
151 
152    HYPRE_ANNOTATE_FUNC_BEGIN;
153 
154    (bicgstab_data -> A) = A;
155 
156    /*--------------------------------------------------
157     * The arguments for NewVector are important to
158     * maintain consistency between the setup and
159     * compute phases of matvec and the preconditioner.
160     *--------------------------------------------------*/
161 
162    if ((bicgstab_data -> p) == NULL)
163       (bicgstab_data -> p) = (*(bicgstab_functions->CreateVector))(b);
164    if ((bicgstab_data -> q) == NULL)
165       (bicgstab_data -> q) = (*(bicgstab_functions->CreateVector))(b);
166    if ((bicgstab_data -> r) == NULL)
167       (bicgstab_data -> r) = (*(bicgstab_functions->CreateVector))(b);
168    if ((bicgstab_data -> r0) == NULL)
169       (bicgstab_data -> r0) = (*(bicgstab_functions->CreateVector))(b);
170    if ((bicgstab_data -> s) == NULL)
171       (bicgstab_data -> s) = (*(bicgstab_functions->CreateVector))(b);
172    if ((bicgstab_data -> v) == NULL)
173       (bicgstab_data -> v) = (*(bicgstab_functions->CreateVector))(b);
174 
175    if ((bicgstab_data -> matvec_data) == NULL)
176       (bicgstab_data -> matvec_data) =
177          (*(bicgstab_functions->MatvecCreate))(A, x);
178 
179    precond_setup(precond_data, A, b, x);
180 
181    /*-----------------------------------------------------
182     * Allocate space for log info
183     *-----------------------------------------------------*/
184 
185    if ((bicgstab_data->logging)>0 || (bicgstab_data->print_level) > 0)
186    {
187       if ((bicgstab_data -> norms) != NULL)
188          hypre_TFree (bicgstab_data -> norms, HYPRE_MEMORY_HOST);
189       (bicgstab_data -> norms) = hypre_CTAlloc(HYPRE_Real,  max_iter + 1, HYPRE_MEMORY_HOST);
190    }
191    if ((bicgstab_data -> print_level) > 0)
192    {
193       if ((bicgstab_data -> log_file_name) == NULL)
194 		  (bicgstab_data -> log_file_name) = (char*)"bicgstab.out.log";
195    }
196 
197    HYPRE_ANNOTATE_FUNC_END;
198 
199    return hypre_error_flag;
200 }
201 
202 /*--------------------------------------------------------------------------
203  * hypre_BiCGSTABSolve
204  *-------------------------------------------------------------------------*/
205 
206 HYPRE_Int
hypre_BiCGSTABSolve(void * bicgstab_vdata,void * A,void * b,void * x)207 hypre_BiCGSTABSolve(void  *bicgstab_vdata,
208                  void  *A,
209                  void  *b,
210 		 void  *x)
211 {
212    hypre_BiCGSTABData      *bicgstab_data      = (hypre_BiCGSTABData*)bicgstab_vdata;
213    hypre_BiCGSTABFunctions *bicgstab_functions = bicgstab_data->functions;
214 
215    HYPRE_Int               min_iter     = (bicgstab_data -> min_iter);
216    HYPRE_Int 		     max_iter     = (bicgstab_data -> max_iter);
217    HYPRE_Int 		     stop_crit    = (bicgstab_data -> stop_crit);
218    HYPRE_Int 		     hybrid    = (bicgstab_data -> hybrid);
219    HYPRE_Real 	     r_tol     = (bicgstab_data -> tol);
220    HYPRE_Real 	     cf_tol       = (bicgstab_data -> cf_tol);
221    void             *matvec_data  = (bicgstab_data -> matvec_data);
222    HYPRE_Real        a_tol        = (bicgstab_data -> a_tol);
223 
224 
225 
226    void             *r            = (bicgstab_data -> r);
227    void             *r0           = (bicgstab_data -> r0);
228    void             *s            = (bicgstab_data -> s);
229    void             *v           = (bicgstab_data -> v);
230    void             *p            = (bicgstab_data -> p);
231    void             *q            = (bicgstab_data -> q);
232 
233    HYPRE_Int 	           (*precond)(void*,void*,void*,void*)   = (bicgstab_functions -> precond);
234    HYPRE_Int 	            *precond_data = (HYPRE_Int*)(bicgstab_data -> precond_data);
235 
236    /* logging variables */
237    HYPRE_Int             logging        = (bicgstab_data -> logging);
238    HYPRE_Int             print_level    = (bicgstab_data -> print_level);
239    HYPRE_Real     *norms          = (bicgstab_data -> norms);
240    /*   char           *log_file_name  = (bicgstab_data -> log_file_name);
241      FILE           *fp; */
242 
243    HYPRE_Int        iter;
244    HYPRE_Int        my_id, num_procs;
245    HYPRE_Real alpha, beta, gamma, epsilon, temp, res, r_norm, b_norm;
246    HYPRE_Real epsmac = HYPRE_REAL_MIN;
247    HYPRE_Real ieee_check = 0.;
248    HYPRE_Real cf_ave_0 = 0.0;
249    HYPRE_Real cf_ave_1 = 0.0;
250    HYPRE_Real weight;
251    HYPRE_Real r_norm_0;
252    HYPRE_Real den_norm;
253    HYPRE_Real gamma_numer;
254    HYPRE_Real gamma_denom;
255 
256    HYPRE_ANNOTATE_FUNC_BEGIN;
257 
258    (bicgstab_data -> converged) = 0;
259 
260    (*(bicgstab_functions->CommInfo))(A,&my_id,&num_procs);
261    if (logging > 0 || print_level > 0)
262    {
263       norms          = (bicgstab_data -> norms);
264       /* log_file_name  = (bicgstab_data -> log_file_name);
265          fp = fopen(log_file_name,"w"); */
266    }
267 
268    /* initialize work arrays */
269    (*(bicgstab_functions->CopyVector))(b,r0);
270 
271    /* compute initial residual */
272 
273    (*(bicgstab_functions->Matvec))(matvec_data,-1.0, A, x, 1.0, r0);
274    (*(bicgstab_functions->CopyVector))(r0,r);
275    (*(bicgstab_functions->CopyVector))(r0,p);
276 
277    b_norm = sqrt((*(bicgstab_functions->InnerProd))(b,b));
278 
279    /* Since it is does not diminish performance, attempt to return an error flag
280       and notify users when they supply bad input. */
281    if (b_norm != 0.) ieee_check = b_norm/b_norm; /* INF -> NaN conversion */
282    if (ieee_check != ieee_check)
283    {
284       /* ...INFs or NaNs in input can make ieee_check a NaN.  This test
285          for ieee_check self-equality works on all IEEE-compliant compilers/
286          machines, c.f. page 8 of "Lecture Notes on the Status of IEEE 754"
287          by W. Kahan, May 31, 1996.  Currently (July 2002) this paper may be
288          found at http://HTTP.CS.Berkeley.EDU/~wkahan/ieee754status/IEEE754.PDF */
289       if (logging > 0 || print_level > 0)
290       {
291         hypre_printf("\n\nERROR detected by Hypre ...  BEGIN\n");
292         hypre_printf("ERROR -- hypre_BiCGSTABSolve: INFs and/or NaNs detected in input.\n");
293         hypre_printf("User probably placed non-numerics in supplied b.\n");
294         hypre_printf("Returning error flag += 101.  Program not terminated.\n");
295         hypre_printf("ERROR detected by Hypre ...  END\n\n\n");
296       }
297       hypre_error(HYPRE_ERROR_GENERIC);
298       HYPRE_ANNOTATE_FUNC_END;
299 
300       return hypre_error_flag;
301    }
302 
303    res = (*(bicgstab_functions->InnerProd))(r0,r0);
304    r_norm = sqrt(res);
305    r_norm_0 = r_norm;
306 
307    /* Since it is does not diminish performance, attempt to return an error flag
308       and notify users when they supply bad input. */
309    if (r_norm != 0.) ieee_check = r_norm/r_norm; /* INF -> NaN conversion */
310    if (ieee_check != ieee_check)
311    {
312       /* ...INFs or NaNs in input can make ieee_check a NaN.  This test
313          for ieee_check self-equality works on all IEEE-compliant compilers/
314          machines, c.f. page 8 of "Lecture Notes on the Status of IEEE 754"
315          by W. Kahan, May 31, 1996.  Currently (July 2002) this paper may be
316          found at http://HTTP.CS.Berkeley.EDU/~wkahan/ieee754status/IEEE754.PDF */
317       if (logging > 0 || print_level > 0)
318       {
319         hypre_printf("\n\nERROR detected by Hypre ...  BEGIN\n");
320         hypre_printf("ERROR -- hypre_BiCGSTABSolve: INFs and/or NaNs detected in input.\n");
321         hypre_printf("User probably placed non-numerics in supplied A or x_0.\n");
322         hypre_printf("Returning error flag += 101.  Program not terminated.\n");
323         hypre_printf("ERROR detected by Hypre ...  END\n\n\n");
324       }
325 
326       hypre_error(HYPRE_ERROR_GENERIC);
327       HYPRE_ANNOTATE_FUNC_END;
328 
329       return hypre_error_flag;
330    }
331 
332    if (logging > 0 || print_level > 0)
333    {
334       norms[0] = r_norm;
335       if (print_level > 0 && my_id == 0)
336       {
337    	     hypre_printf("L2 norm of b: %e\n", b_norm);
338          if (b_norm == 0.0)
339             hypre_printf("Rel_resid_norm actually contains the residual norm\n");
340          hypre_printf("Initial L2 norm of residual: %e\n", r_norm);
341       }
342    }
343    iter = 0;
344 
345    if (b_norm > 0.0)
346    {
347       /* convergence criterion |r_i| <= r_tol*|b| if |b| > 0 */
348       den_norm = b_norm;
349    }
350    else
351    {
352       /* convergence criterion |r_i| <= r_tol*|r0| if |b| = 0 */
353       den_norm = r_norm;
354    };
355 
356    /* convergence criterion |r_i| <= r_tol/a_tol , absolute residual norm*/
357    if (stop_crit)
358    {
359       if (a_tol == 0.0) /* this is for backwards compatibility
360                            (accomodating setting stop_crit to 1, but not setting a_tol) -
361                            eventually we will get rid of the stop_crit flag as with GMRES */
362          epsilon = r_tol;
363       else
364          epsilon = a_tol; /* this means new interface fcn called */
365 
366    }
367    else /* default convergence test (stop_crit = 0)*/
368    {
369 
370       /* convergence criteria: |r_i| <= max( a_tol, r_tol * den_norm)
371       den_norm = |r_0| or |b|
372       note: default for a_tol is 0.0, so relative residual criteria is used unless
373             user also specifies a_tol or sets r_tol = 0.0, which means absolute
374             tol only is checked  */
375 
376       epsilon = hypre_max(a_tol, r_tol*den_norm);
377 
378    }
379 
380 
381    if (print_level > 0 && my_id == 0)
382    {
383       if (b_norm > 0.0)
384          {hypre_printf("=============================================\n\n");
385           hypre_printf("Iters     resid.norm     conv.rate  rel.res.norm\n");
386           hypre_printf("-----    ------------    ---------- ------------\n");
387       }
388       else
389          {hypre_printf("=============================================\n\n");
390           hypre_printf("Iters     resid.norm     conv.rate\n");
391           hypre_printf("-----    ------------    ----------\n");
392 
393       }
394    }
395 
396    (bicgstab_data -> num_iterations) = iter;
397    if (b_norm > 0.0)
398       (bicgstab_data -> rel_residual_norm) = r_norm/b_norm;
399    /* check for convergence before starting */
400    if (r_norm == 0.0)
401    {
402        HYPRE_ANNOTATE_FUNC_END;
403        return hypre_error_flag;
404    }
405    else if (r_norm <= epsilon && iter >= min_iter)
406    {
407        if (print_level > 0 && my_id == 0)
408        {
409           hypre_printf("\n\n");
410           hypre_printf("Tolerance and min_iter requirements satisfied by initial data.\n");
411           hypre_printf("Final L2 norm of residual: %e\n\n", r_norm);
412        }
413        (bicgstab_data -> converged) = 1;
414        HYPRE_ANNOTATE_FUNC_END;
415 
416        return hypre_error_flag;
417    }
418    /* Start BiCGStab iterations */
419    while (iter < max_iter)
420    {
421         iter++;
422 
423 	(*(bicgstab_functions->ClearVector))(v);
424         precond(precond_data, A, p, v);
425         (*(bicgstab_functions->Matvec))(matvec_data,1.0,A,v,0.0,q);
426       	temp = (*(bicgstab_functions->InnerProd))(r0,q);
427       	if (fabs(temp) >= epsmac)
428 	   alpha = res/temp;
429 	else
430 	{
431 	   hypre_error_w_msg(HYPRE_ERROR_GENERIC,"BiCGSTAB broke down!! divide by near zero\n");
432            HYPRE_ANNOTATE_FUNC_END;
433 
434 	   return hypre_error_flag;
435 	}
436 	(*(bicgstab_functions->Axpy))(alpha,v,x);
437 	(*(bicgstab_functions->Axpy))(-alpha,q,r);
438 	(*(bicgstab_functions->ClearVector))(v);
439         precond(precond_data, A, r, v);
440         (*(bicgstab_functions->Matvec))(matvec_data,1.0,A,v,0.0,s);
441       	/* Handle case when gamma = 0.0/0.0 as 0.0 and not NAN */
442         gamma_numer = (*(bicgstab_functions->InnerProd))(r,s);
443         gamma_denom = (*(bicgstab_functions->InnerProd))(s,s);
444         if ((gamma_numer == 0.0) && (gamma_denom == 0.0))
445             gamma = 0.0;
446         else
447             gamma= gamma_numer/gamma_denom;
448 	(*(bicgstab_functions->Axpy))(gamma,v,x);
449 	(*(bicgstab_functions->Axpy))(-gamma,s,r);
450     /* residual is now updated, must immediately check for convergence */
451 	r_norm = sqrt((*(bicgstab_functions->InnerProd))(r,r));
452 	if (logging > 0 || print_level > 0)
453 	{
454 	   norms[iter] = r_norm;
455 	}
456         if (print_level > 0 && my_id == 0)
457         {
458            if (b_norm > 0.0)
459               hypre_printf("% 5d    %e    %f   %e\n", iter, norms[iter],
460                       norms[iter]/norms[iter-1], norms[iter]/b_norm);
461            else
462               hypre_printf("% 5d    %e    %f\n", iter, norms[iter],
463 		                             norms[iter]/norms[iter-1]);
464 	}
465     /* check for convergence, evaluate actual residual */
466 	if (r_norm <= epsilon && iter >= min_iter)
467         {
468 	   (*(bicgstab_functions->CopyVector))(b,r);
469            (*(bicgstab_functions->Matvec))(matvec_data,-1.0,A,x,1.0,r);
470 	   r_norm = sqrt((*(bicgstab_functions->InnerProd))(r,r));
471 	   if (r_norm <= epsilon)
472            {
473               if (print_level > 0 && my_id == 0)
474               {
475                  hypre_printf("\n\n");
476                  hypre_printf("Final L2 norm of residual: %e\n\n", r_norm);
477               }
478               (bicgstab_data -> converged) = 1;
479               break;
480           }
481        }
482     /*--------------------------------------------------------------------
483      * Optional test to see if adequate progress is being made.
484      * The average convergence factor is recorded and compared
485      * against the tolerance 'cf_tol'. The weighting factor is
486      * intended to pay more attention to the test when an accurate
487      * estimate for average convergence factor is available.
488      *--------------------------------------------------------------------*/
489        if (cf_tol > 0.0)
490        {
491           cf_ave_0 = cf_ave_1;
492           cf_ave_1 = pow( r_norm / r_norm_0, 1.0/(2.0*iter));
493 
494           weight   = fabs(cf_ave_1 - cf_ave_0);
495           weight   = weight / hypre_max(cf_ave_1, cf_ave_0);
496           weight   = 1.0 - weight;
497           if (weight * cf_ave_1 > cf_tol) break;
498        }
499 
500        if (fabs(res) >= epsmac)
501           beta = 1.0/res;
502        else
503        {
504 	  hypre_error_w_msg(HYPRE_ERROR_GENERIC,"BiCGSTAB broke down!! res=0 \n");
505           HYPRE_ANNOTATE_FUNC_END;
506 
507 	  return hypre_error_flag;
508        }
509        res = (*(bicgstab_functions->InnerProd))(r0,r);
510        beta *= res;
511        (*(bicgstab_functions->Axpy))(-gamma,q,p);
512        if (fabs(gamma) >= epsmac)
513            (*(bicgstab_functions->ScaleVector))((beta*alpha/gamma),p);
514        else
515        {
516 	  hypre_error_w_msg(HYPRE_ERROR_GENERIC,"BiCGSTAB broke down!! gamma=0 \n");
517           HYPRE_ANNOTATE_FUNC_END;
518 
519 	  return hypre_error_flag;
520        }
521        (*(bicgstab_functions->Axpy))(1.0,r,p);
522    } /* end while loop */
523 
524    (bicgstab_data -> num_iterations) = iter;
525    if (b_norm > 0.0)
526       (bicgstab_data -> rel_residual_norm) = r_norm/b_norm;
527    if (b_norm == 0.0)
528       (bicgstab_data -> rel_residual_norm) = r_norm;
529 
530    if (iter >= max_iter && r_norm > epsilon && epsilon > 0 && hybrid != -1) hypre_error(HYPRE_ERROR_CONV);
531 
532    HYPRE_ANNOTATE_FUNC_END;
533 
534    return hypre_error_flag;
535 }
536 
537 /*--------------------------------------------------------------------------
538  * hypre_BiCGSTABSetTol
539  *--------------------------------------------------------------------------*/
540 
541 HYPRE_Int
hypre_BiCGSTABSetTol(void * bicgstab_vdata,HYPRE_Real tol)542 hypre_BiCGSTABSetTol( void   *bicgstab_vdata,
543                    HYPRE_Real  tol       )
544 {
545 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
546 
547    (bicgstab_data -> tol) = tol;
548 
549    return hypre_error_flag;
550 }
551 /*--------------------------------------------------------------------------
552  * hypre_BiCGSTABSetAbsoluteTol
553  *--------------------------------------------------------------------------*/
554 
555 HYPRE_Int
hypre_BiCGSTABSetAbsoluteTol(void * bicgstab_vdata,HYPRE_Real a_tol)556 hypre_BiCGSTABSetAbsoluteTol( void   *bicgstab_vdata,
557                    HYPRE_Real  a_tol       )
558 {
559 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
560 
561    (bicgstab_data -> a_tol) = a_tol;
562 
563    return hypre_error_flag;
564 }
565 
566 /*--------------------------------------------------------------------------
567  * hypre_BiCGSTABSetConvergenceFactorTol
568  *--------------------------------------------------------------------------*/
569 
570 HYPRE_Int
hypre_BiCGSTABSetConvergenceFactorTol(void * bicgstab_vdata,HYPRE_Real cf_tol)571 hypre_BiCGSTABSetConvergenceFactorTol( void   *bicgstab_vdata,
572                    HYPRE_Real  cf_tol       )
573 {
574 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
575 
576    (bicgstab_data -> cf_tol) = cf_tol;
577 
578    return hypre_error_flag;
579 }
580 
581 /*--------------------------------------------------------------------------
582  * hypre_BiCGSTABSetMinIter
583  *--------------------------------------------------------------------------*/
584 
585 HYPRE_Int
hypre_BiCGSTABSetMinIter(void * bicgstab_vdata,HYPRE_Int min_iter)586 hypre_BiCGSTABSetMinIter( void *bicgstab_vdata,
587                        HYPRE_Int   min_iter  )
588 {
589 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
590 
591    (bicgstab_data -> min_iter) = min_iter;
592 
593    return hypre_error_flag;
594 }
595 
596 /*--------------------------------------------------------------------------
597  * hypre_BiCGSTABSetMaxIter
598  *--------------------------------------------------------------------------*/
599 
600 HYPRE_Int
hypre_BiCGSTABSetMaxIter(void * bicgstab_vdata,HYPRE_Int max_iter)601 hypre_BiCGSTABSetMaxIter( void *bicgstab_vdata,
602                        HYPRE_Int   max_iter  )
603 {
604 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
605 
606    (bicgstab_data -> max_iter) = max_iter;
607 
608    return hypre_error_flag;
609 }
610 
611 /*--------------------------------------------------------------------------
612  * hypre_BiCGSTABSetStopCrit
613  *--------------------------------------------------------------------------*/
614 
615 HYPRE_Int
hypre_BiCGSTABSetStopCrit(void * bicgstab_vdata,HYPRE_Int stop_crit)616 hypre_BiCGSTABSetStopCrit( void   *bicgstab_vdata,
617                         HYPRE_Int  stop_crit       )
618 {
619 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
620 
621    (bicgstab_data -> stop_crit) = stop_crit;
622 
623    return hypre_error_flag;
624 }
625 
626 /*--------------------------------------------------------------------------
627  * hypre_BiCGSTABSetPrecond
628  *--------------------------------------------------------------------------*/
629 
630 HYPRE_Int
hypre_BiCGSTABSetPrecond(void * bicgstab_vdata,HYPRE_Int (* precond)(void *,void *,void *,void *),HYPRE_Int (* precond_setup)(void *,void *,void *,void *),void * precond_data)631 hypre_BiCGSTABSetPrecond( void  *bicgstab_vdata,
632 						  HYPRE_Int  (*precond)(void*,void*,void*,void*),
633 						  HYPRE_Int  (*precond_setup)(void*,void*,void*,void*),
634                        void  *precond_data )
635 {
636 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
637    hypre_BiCGSTABFunctions *bicgstab_functions = bicgstab_data->functions;
638 
639 
640    (bicgstab_functions -> precond)        = precond;
641    (bicgstab_functions -> precond_setup)  = precond_setup;
642    (bicgstab_data -> precond_data)   = precond_data;
643 
644    return hypre_error_flag;
645 }
646 
647 /*--------------------------------------------------------------------------
648  * hypre_BiCGSTABGetPrecond
649  *--------------------------------------------------------------------------*/
650 
651 HYPRE_Int
hypre_BiCGSTABGetPrecond(void * bicgstab_vdata,HYPRE_Solver * precond_data_ptr)652 hypre_BiCGSTABGetPrecond( void         *bicgstab_vdata,
653                        HYPRE_Solver *precond_data_ptr )
654 {
655 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
656 
657    *precond_data_ptr = (HYPRE_Solver)(bicgstab_data -> precond_data);
658 
659    return hypre_error_flag;
660 }
661 
662 /*--------------------------------------------------------------------------
663  * hypre_BiCGSTABSetLogging
664  *--------------------------------------------------------------------------*/
665 
666 HYPRE_Int
hypre_BiCGSTABSetLogging(void * bicgstab_vdata,HYPRE_Int logging)667 hypre_BiCGSTABSetLogging( void *bicgstab_vdata,
668                        HYPRE_Int   logging)
669 {
670 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
671 
672    (bicgstab_data -> logging) = logging;
673 
674    return hypre_error_flag;
675 }
676 
677 HYPRE_Int
hypre_BiCGSTABSetHybrid(void * bicgstab_vdata,HYPRE_Int logging)678 hypre_BiCGSTABSetHybrid( void *bicgstab_vdata,
679                        HYPRE_Int   logging)
680 {
681 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
682 
683    (bicgstab_data -> hybrid) = logging;
684 
685    return hypre_error_flag;
686 }
687 
688 /*--------------------------------------------------------------------------
689  * hypre_BiCGSTABSetPrintLevel
690  *--------------------------------------------------------------------------*/
691 
692 HYPRE_Int
hypre_BiCGSTABSetPrintLevel(void * bicgstab_vdata,HYPRE_Int print_level)693 hypre_BiCGSTABSetPrintLevel( void *bicgstab_vdata,
694                        HYPRE_Int   print_level)
695 {
696 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
697 
698    (bicgstab_data -> print_level) = print_level;
699 
700    return hypre_error_flag;
701 }
702 
703 /*--------------------------------------------------------------------------
704  * hypre_BiCGSTABGetConverged
705  *--------------------------------------------------------------------------*/
706 
707 HYPRE_Int
hypre_BiCGSTABGetConverged(void * bicgstab_vdata,HYPRE_Int * converged)708 hypre_BiCGSTABGetConverged( void *bicgstab_vdata,
709                              HYPRE_Int  *converged )
710 {
711 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
712 
713    *converged = (bicgstab_data -> converged);
714 
715    return hypre_error_flag;
716 }
717 
718 /*--------------------------------------------------------------------------
719  * hypre_BiCGSTABGetNumIterations
720  *--------------------------------------------------------------------------*/
721 
722 HYPRE_Int
hypre_BiCGSTABGetNumIterations(void * bicgstab_vdata,HYPRE_Int * num_iterations)723 hypre_BiCGSTABGetNumIterations( void *bicgstab_vdata,
724                              HYPRE_Int  *num_iterations )
725 {
726 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
727 
728    *num_iterations = (bicgstab_data -> num_iterations);
729 
730    return hypre_error_flag;
731 }
732 
733 /*--------------------------------------------------------------------------
734  * hypre_BiCGSTABGetFinalRelativeResidualNorm
735  *--------------------------------------------------------------------------*/
736 
737 HYPRE_Int
hypre_BiCGSTABGetFinalRelativeResidualNorm(void * bicgstab_vdata,HYPRE_Real * relative_residual_norm)738 hypre_BiCGSTABGetFinalRelativeResidualNorm( void   *bicgstab_vdata,
739                                          HYPRE_Real *relative_residual_norm )
740 {
741 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
742 
743    *relative_residual_norm = (bicgstab_data -> rel_residual_norm);
744 
745    return hypre_error_flag;
746 }
747 
748 /*--------------------------------------------------------------------------
749  * hypre_BiCGSTABGetResidual
750  *--------------------------------------------------------------------------*/
751 
752 HYPRE_Int
hypre_BiCGSTABGetResidual(void * bicgstab_vdata,void ** residual)753 hypre_BiCGSTABGetResidual( void   *bicgstab_vdata,
754                            void **residual )
755 {
756 	hypre_BiCGSTABData *bicgstab_data = (hypre_BiCGSTABData  *)bicgstab_vdata;
757 
758    *residual = (bicgstab_data -> r);
759 
760    return hypre_error_flag;
761 }
762