#include <iostream>
 
 
{
    rhs.setZero(N,1);
 
    mat.resize(N,N);
    mat.setZero();
    real_t meshSize = 1./(N+1);
    real_t pi = EIGEN_PI;
 
    
    mat.reservePerColumn( 3 ); 
 
    mat(0,0) = 2;
    mat(0,1) = -1;
    mat(N-1, N-1) = 2;
    mat(N-1, N-2) = -1;
    {
        mat(k,k) = 2;
        mat(k,k-1) = -1;
        mat(k,k+1) = -1;
    }
        rhs(k,0) = pi*pi*meshSize*meshSize*math::cos(meshSize*(1+k)*pi);
 
    
    mat.makeCompressed();
}
 
template<typename SolverType>
void gsIterativeSolverInfo(const SolverType &method,
                           real_t error, double time, bool& succeeded )
{
    gsInfo << 
" Computed res. error  : " << error << 
"\n";
 
    gsInfo << 
" Time to solve:       : " << time << 
"\n";
 
    if ( method.error() <= method.tolerance() && error <= method.tolerance() )
    {
    }
    else
    {
        succeeded = false;
    }
}
 
int main(int argc, char *argv[])
{
 
    bool succeeded = true;
 
    
    real_t tol = std::pow(10.0, - REAL_DIG * 0.75);
 
    gsCmdLine cmd(
"Solves a 1D PDE with a Courant discretization with several solvers.");
 
    cmd.addInt ("n", "number", "Number of unknowns",                  N  );
    cmd.addReal("",  "tol",    "Tolerance for the iterative solvers", tol);
 
    try { cmd.getValues(argc,argv); } catch (int rv) { return rv; }
 
 
    
    poissonDiscretization(mat, rhs, N);
 
    
    
 
 
    
 
#ifndef gsGmp_ENABLED
 
    
 
    opt.
setInt (
"MaxIterations", 3*N);
 
 
    gsInfo << 
"Testing G+Smo's linear solvers:\n";
 
 
 
    
    MinRes.setOptions(opt);
 
    
    x0.setZero(N,1);
 
    
    gsInfo << 
"\nMinRes: Started solving... ";
 
    MinRes.solve(rhs, x0);
    gsIterativeSolverInfo(MinRes, (mat*x0-rhs).norm()/rhs.norm(), clock.
stop(), succeeded);
 
 
 
    
    MinResIR.setOptions(opt);
    MinResIR.setInexactResidual(true);
 
    
    x0.setZero(N,1);
 
    
    gsInfo << 
"\nMinResIR: Started solving... ";
 
    MinResIR.solve(rhs, x0);
    gsIterativeSolverInfo(MinResIR, (mat*x0-rhs).norm()/rhs.norm(), clock.
stop(), succeeded);
 
 
    
    GMResSolver.setOptions(opt);
 
    
    x0.setZero(N,1);
 
    if (N < 200)
    {
        
        gsInfo << 
"\nGMRes: Started solving... ";
 
        GMResSolver.solve(rhs,x0);
        gsIterativeSolverInfo(GMResSolver, (mat*x0-rhs).norm()/rhs.norm(), clock.
stop(), succeeded);
 
    }
    else
        gsInfo << 
"\nSkipping GMRes due to high number of iterations...\n";
 
 
 
    
    CGSolver.setOptions(opt);
 
    
    x0.setZero(N,1);
 
    
    gsInfo << 
"\nCG: Started solving... ";
 
    CGSolver.solve(rhs,x0);
    gsIterativeSolverInfo(CGSolver, (mat*x0-rhs).norm()/rhs.norm(), clock.
stop(), succeeded);
 
 
    
    MRQLPSolver.setOptions(opt);
 
    
    x0.setZero(N,1);
 
    
    gsInfo << 
"\nMRQLPSolver: Started solving... ";
 
    MRQLPSolver.solve(rhs,x0);
    gsIterativeSolverInfo(MRQLPSolver, (mat*x0-rhs).norm()/rhs.norm(), clock.
stop(), succeeded);
 
 
 
    gsInfo << 
"\nTesting Eigen's interative solvers:\n";
 
 
    gsSparseSolver<>::CGIdentity EigenCGIsolver;
    EigenCGIsolver.setMaxIterations(maxIters);
    EigenCGIsolver.setTolerance(tol);
    gsInfo << 
"\nEigen's CG + identity prec.: Started solving... ";
 
    EigenCGIsolver.compute(mat);
    x0 = EigenCGIsolver.solve(rhs);
    gsIterativeSolverInfo(EigenCGIsolver, (mat*x0-rhs).norm()/rhs.norm(), clock.
stop(), succeeded);
 
 
    gsSparseSolver<>::CGDiagonal EigenCGDsolver;
    EigenCGDsolver.setMaxIterations(maxIters);
    EigenCGDsolver.setTolerance(tol);
    gsInfo << 
"\nEigen's CG + diagonal prec.: Started solving... ";
 
    EigenCGDsolver.compute(mat);
    x0 = EigenCGDsolver.solve(rhs);
    gsIterativeSolverInfo(EigenCGDsolver, (mat*x0-rhs).norm()/rhs.norm(), clock.
stop(), succeeded);
 
 
    gsSparseSolver<>::BiCGSTABIdentity EigenBCGIsolver;
    EigenBCGIsolver.setMaxIterations(maxIters);
    EigenBCGIsolver.setTolerance(tol);
    gsInfo << 
"\nEigen's bi conjugate gradient stabilized + identity prec.: Started solving... ";
 
    EigenBCGIsolver.compute(mat);
    x0 = EigenBCGIsolver.solve(rhs);
    gsIterativeSolverInfo(EigenBCGIsolver, (mat*x0-rhs).norm()/rhs.norm(), clock.
stop(), succeeded);
 
 
    gsSparseSolver<>::BiCGSTABDiagonal EigenBCGDsolver;
    EigenBCGDsolver.setMaxIterations(maxIters);
    EigenBCGDsolver.setTolerance(tol);
    gsInfo << 
"\nEigen's bi conjugate gradient stabilized solver + diagonal prec.: Started solving... ";
 
    EigenBCGDsolver.compute(mat);
    x0 = EigenBCGDsolver.solve(rhs);
    gsIterativeSolverInfo(EigenBCGDsolver, (mat*x0-rhs).norm()/rhs.norm(), clock.
stop(), succeeded);
 
 
    gsSparseSolver<>::BiCGSTABILUT EigenBCGILUsolver;
    
    EigenBCGILUsolver.setMaxIterations(maxIters);
    EigenBCGILUsolver.setTolerance(tol);
    gsInfo << 
"\nEigen's bi conjugate gradient stabilized solver + ILU prec.: Started solving... ";
 
    EigenBCGILUsolver.compute(mat);
    x0 = EigenBCGILUsolver.solve(rhs);
    gsIterativeSolverInfo(EigenBCGILUsolver, (mat*x0-rhs).norm()/rhs.norm(), clock.
stop(), succeeded);
 
 
    gsSparseSolver<>::SimplicialLDLT EigenSLDLTsolver;
    gsInfo << 
"\nEigen's Simplicial LDLT: Started solving... ";
 
    EigenSLDLTsolver.compute(mat);
    x0 = EigenSLDLTsolver.solve(rhs);
    gsInfo << 
"Eigen's Simplicial LDLT: Time to solve       : " << clock.
stop() << 
"\n";
 
 
    gsSparseSolver<>::QR solverQR;
    gsInfo << 
"\nEigen's QR: Started solving... ";
 
    solverQR.compute(mat);
    x0 = solverQR.solve(rhs);
    gsInfo << 
"Eigen's QR: Time to solve       : " << clock.
stop() << 
"\n";
 
 
#endif
 
    gsSparseSolver<>::LU solverLU;
    gsInfo << 
"\nEigen's LU: Started solving... ";
 
    solverLU.compute(mat);
    x0 = solverLU.solve(rhs);
    gsInfo << 
"Eigen's LU: Time to solve       : " << clock.
stop() << 
"\n";
 
 
 
    return succeeded ? EXIT_SUCCESS : EXIT_FAILURE;
}
Class for command-line argument parsing.
Definition gsCmdLine.h:57
 
The conjugate gradient method.
Definition gsConjugateGradient.h:30
 
The generalized minimal residual (GMRES) method.
Definition gsGMRes.h:25
 
double stop()
Return elapsed time in seconds.
Definition gsStopwatch.h:83
 
void restart()
Start taking the time.
Definition gsStopwatch.h:80
 
static uPtr make(index_t dim)
Make function returning a smart pointer.
Definition gsLinearOperator.h:132
 
static gsOptionList defaultOptions()
Returns a list of default options.
Definition gsIterativeSolver.h:92
 
memory::shared_ptr< gsLinearOperator > Ptr
Shared pointer for gsLinearOperator.
Definition gsLinearOperator.h:33
 
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
 
The minimal residual (MinRes-QLP) method.
Definition gsMinResQLP.h:34
 
The minimal residual (MinRes) method.
Definition gsMinimalResidual.h:26
 
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
 
void setReal(const std::string &label, const Real &value)
Sets an existing option label to be equal to value.
Definition gsOptionList.cpp:166
 
void setInt(const std::string &label, const index_t &value)
Sets an existing option label to be equal to value.
Definition gsOptionList.cpp:158
 
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
 
Main header to be included by clients using the G+Smo library.
 
#define index_t
Definition gsConfig.h:32
 
#define gsInfo
Definition gsDebug.h:43
 
The G+Smo namespace, containing all definitions for the library.