hypre/parcsr_ls/F90_HYPRE_parcsr_bicgstab.c
2000-07-13 22:59:17 +00:00

232 lines
9.1 KiB
C

/*BHEADER**********************************************************************
* (c) 1998 The Regents of the University of California
*
* See the file COPYRIGHT_and_DISCLAIMER for a complete copyright
* notice, contact person, and disclaimer.
*
* $Revision$
*********************************************************************EHEADER*/
/******************************************************************************
*
* HYPRE_ParCSRBiCGSTAB Fortran interface
*
*****************************************************************************/
#include "headers.h"
#include "fortran.h"
/*--------------------------------------------------------------------------
* HYPRE_ParCSRBiCGSTABCreate
*--------------------------------------------------------------------------*/
void
hypre_F90_IFACE(hypre_parcsrbicgstabcreate)( int *comm,
long int *solver,
int *ierr )
{
*ierr = (int) ( HYPRE_ParCSRBiCGSTABCreate( (MPI_Comm) *comm,
(HYPRE_Solver *) solver ) );
}
/*--------------------------------------------------------------------------
* HYPRE_ParCSRBiCGSTABDestroy
*--------------------------------------------------------------------------*/
void
hypre_F90_IFACE(hypre_parcsrbicgstabdestroy)( long int *solver,
int *ierr )
{
*ierr = (int) ( HYPRE_ParCSRBiCGSTABDestroy( (HYPRE_Solver) *solver ) );
}
/*--------------------------------------------------------------------------
* HYPRE_ParCSRBiCGSTABSetup
*--------------------------------------------------------------------------*/
void
hypre_F90_IFACE(hypre_parcsrbicgstabsetup)( long int *solver,
long int *A,
long int *b,
long int *x,
int *ierr )
{
*ierr = (int) ( HYPRE_ParCSRBiCGSTABSetup( (HYPRE_Solver) *solver,
(HYPRE_ParCSRMatrix) *A,
(HYPRE_ParVector) *b,
(HYPRE_ParVector) *x ) );
}
/*--------------------------------------------------------------------------
* HYPRE_ParCSRBiCGSTABSolve
*--------------------------------------------------------------------------*/
void
hypre_F90_IFACE(hypre_parcsrbicgstabsolve)( long int *solver,
long int *A,
long int *b,
long int *x,
int *ierr )
{
*ierr = (int) ( HYPRE_ParCSRBiCGSTABSolve( (HYPRE_Solver) *solver,
(HYPRE_ParCSRMatrix) *A,
(HYPRE_ParVector) *b,
(HYPRE_ParVector) *x ) );
}
/*--------------------------------------------------------------------------
* HYPRE_ParCSRBiCGSTABSetTol
*--------------------------------------------------------------------------*/
void
hypre_F90_IFACE(hypre_parcsrbicgstabsettol)( long int *solver,
double *tol,
int *ierr )
{
*ierr = (int) ( HYPRE_ParCSRBiCGSTABSetTol( (HYPRE_Solver) *solver,
(double) *tol ) );
}
/*--------------------------------------------------------------------------
* HYPRE_ParCSRBiCGSTABSetMinIter
*--------------------------------------------------------------------------*/
void
hypre_F90_IFACE(hypre_parcsrbicgstabsetminiter)( long int *solver,
int *min_iter,
int *ierr )
{
*ierr = (int) ( HYPRE_ParCSRBiCGSTABSetMinIter( (HYPRE_Solver) *solver,
(int) *min_iter ) );
}
/*--------------------------------------------------------------------------
* HYPRE_ParCSRBiCGSTABSetMaxIter
*--------------------------------------------------------------------------*/
void
hypre_F90_IFACE(hypre_parcsrbicgstabsetmaxiter)( long int *solver,
int *max_iter,
int *ierr )
{
*ierr = (int) ( HYPRE_ParCSRBiCGSTABSetMaxIter( (HYPRE_Solver) *solver,
(int) *max_iter ) );
}
/*--------------------------------------------------------------------------
* HYPRE_ParCSRBiCGSTABSetPrecond
*--------------------------------------------------------------------------*/
void
hypre_F90_IFACE(hypre_parcsrbicgstabsetprecond)( long int *solver,
int *precond_id,
long int *precond_solver,
int *ierr )
{
/*------------------------------------------------------------
* The precond_id flags mean :
* 0 - no preconditioner
* 1 - set up a ds preconditioner
* 2 - set up an amg preconditioner
* 3 - set up a pilut preconditioner
* 4 - set up a parasails preconditioner
*------------------------------------------------------------*/
if (*precond_id == 0)
{
*ierr = 0;
}
if (*precond_id == 1)
{
*ierr = (int)
( HYPRE_ParCSRBiCGSTABSetPrecond( (HYPRE_Solver) *solver,
HYPRE_ParCSRDiagScale,
HYPRE_ParCSRDiagScaleSetup,
NULL ) );
}
else if (*precond_id == 2)
{
*ierr = (int) ( HYPRE_ParCSRBiCGSTABSetPrecond( (HYPRE_Solver) *solver,
HYPRE_BoomerAMGSolve,
HYPRE_BoomerAMGSetup,
(void *) *precond_solver ) );
}
else if (*precond_id == 3)
{
*ierr = (int)
( HYPRE_ParCSRBiCGSTABSetPrecond( (HYPRE_Solver) *solver,
HYPRE_ParCSRPilutSolve,
HYPRE_ParCSRPilutSetup,
(void *) *precond_solver ) );
}
else if (*precond_id == 4)
{
*ierr = (int)
( HYPRE_ParCSRBiCGSTABSetPrecond( (HYPRE_Solver) *solver,
HYPRE_ParCSRParaSailsSolve,
HYPRE_ParCSRParaSailsSetup,
(void *) *precond_solver ) );
}
else
{
*ierr = -1;
}
}
/*--------------------------------------------------------------------------
* HYPRE_ParCSRBiCGSTABGetPrecond
*--------------------------------------------------------------------------*/
void
hypre_F90_IFACE(hypre_parcsrbicgstabgetprecond)( long int *solver,
long int *precond_solver_ptr,
int *ierr )
{
*ierr = (int)
( HYPRE_ParCSRBiCGSTABGetPrecond( (HYPRE_Solver) *solver,
(HYPRE_Solver *) precond_solver_ptr ) );
}
/*--------------------------------------------------------------------------
* HYPRE_ParCSRBiCGSTABSetLogging
*--------------------------------------------------------------------------*/
void
hypre_F90_IFACE(hypre_parcsrbicgstabsetlogging)( long int *solver,
int *logging,
int *ierr )
{
*ierr = (int) ( HYPRE_ParCSRBiCGSTABSetLogging( (HYPRE_Solver) *solver,
(int) *logging ) );
}
/*--------------------------------------------------------------------------
* HYPRE_ParCSRBiCGSTABGetNumIter
*--------------------------------------------------------------------------*/
void
hypre_F90_IFACE(hypre_parcsrbicgstabgetnumiter)( long int *solver,
int *num_iterations,
int *ierr )
{
*ierr = (int) ( HYPRE_ParCSRBiCGSTABGetNumIterations(
(HYPRE_Solver) *solver,
(int *) num_iterations ) );
}
/*--------------------------------------------------------------------------
* HYPRE_ParCSRBiCGSTABGetFinalRel
*--------------------------------------------------------------------------*/
void
hypre_F90_IFACE(hypre_parcsrbicgstabgetfinalrel)( long int *solver,
double *norm,
int *ierr )
{
*ierr = (int) ( HYPRE_ParCSRBiCGSTABGetFinalRelativeResidualNorm(
(HYPRE_Solver) *solver,
(double *) norm ) );
}