24#include "fasp_functs.h"
30#include "BlaSpmvMatFree.inl"
70 REAL solve_start, solve_end, solve_time;
76 printf(
"### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
77 printf(
"### DEBUG: rhs/sol size: %d %d\n", b->
row, x->
row);
81 ITS_CHECK(MaxIt, tol);
84 switch (itsolver_type) {
87 iter = fasp_solver_pcg(mf, b, x, pc, tol, abstol, MaxIt, stop_type, prtlvl);
121 printf(
"### ERROR: Unknown iterative solver type %d! [%s]\n", itsolver_type,
128 solve_time = solve_end - solve_start;
133 printf(
"### DEBUG: [--End--] %s ...\n", __FUNCTION__);
163 REAL solve_start, solve_end, solve_time;
166 printf(
"### DEBUG: [-Begin-] %s ...\n", __FUNCTION__);
167 printf(
"### DEBUG: rhs/sol size: %d %d\n", b->
row, x->
row);
176 solve_time = solve_end - solve_start;
181 printf(
"### DEBUG: [--End--] %s ...\n", __FUNCTION__);
203 switch (matrix_format) {
206 mf->
fct = fasp_blas_mxv_csr;
210 mf->
fct = fasp_blas_mxv_bsr;
214 mf->
fct = fasp_blas_mxv_str;
218 mf->
fct = fasp_blas_mxv_blc;
222 mf->
fct = fasp_blas_mxv_csrl;
226 printf(
"### ERROR: Unknown matrix format %d!\n", matrix_format);
void fasp_cputime(const char *message, const REAL cputime)
Print CPU walltime.
void fasp_gettime(REAL *time)
Get system time.
INT fasp_solver_pbcgs(mxv_matfree *mf, dvector *b, dvector *u, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT StopType, const SHORT PrtLvl)
Preconditioned BiCGstab method for solving Au=b.
INT fasp_solver_pgcg(mxv_matfree *mf, dvector *b, dvector *u, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT StopType, const SHORT PrtLvl)
Preconditioned generilzed conjugate gradient (GCG) method for solving Au=b.
INT fasp_solver_pgmres(mxv_matfree *mf, dvector *b, dvector *x, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT restart, const SHORT StopType, const SHORT PrtLvl)
Solve "Ax=b" using PGMRES (right preconditioned) iterative method.
INT fasp_solver_pminres(mxv_matfree *mf, dvector *b, dvector *u, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT StopType, const SHORT PrtLvl)
A preconditioned minimal residual (Minres) method for solving Au=b.
INT fasp_solver_pvfgmres(mxv_matfree *mf, dvector *b, dvector *x, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, const SHORT restart, const SHORT StopType, const SHORT PrtLvl)
Solve "Ax=b" using PFGMRES(right preconditioned) iterative method in which the restart parameter can ...
INT fasp_solver_pvgmres(mxv_matfree *mf, dvector *b, dvector *x, precond *pc, const REAL tol, const REAL abstol, const INT MaxIt, SHORT restart, const SHORT StopType, const SHORT PrtLvl)
Solve "Ax=b" using PGMRES(right preconditioned) iterative method in which the restart parameter can b...
void fasp_solver_matfree_init(INT matrix_format, mxv_matfree *mf, void *A)
Initialize MatFree (or non-specified format) itsovlers.
INT fasp_solver_krylov(mxv_matfree *mf, dvector *b, dvector *x, ITS_param *itparam)
Solve Ax=b by standard Krylov methods – without preconditioner.
INT fasp_solver_itsolver(mxv_matfree *mf, dvector *b, dvector *x, precond *pc, ITS_param *itparam)
Solve Ax=b by preconditioned Krylov methods for CSR matrices.
Main header file for the FASP project.
#define SHORT
FASP integer and floating point numbers.
Header file for FASP block matrices.
#define FASP_SUCCESS
Definition of return status and error messages.
#define ERROR_SOLVER_TYPE
#define ERROR_DATA_STRUCTURE
Parameters for iterative solvers.
Vector with n entries of REAL type.
Matrix-vector multiplication, replace the actual matrix.
void * data
data for MxV, can be a Matrix or something else
void(* fct)(const void *, const REAL *, REAL *)
action for MxV, void function pointer
Preconditioner data and action.