Skip to content
Snippets Groups Projects
Commit 46eb2e5f authored by Adrien Crovato's avatar Adrien Crovato
Browse files

Add type() for LinearSolver and Writer. Refactor Gmres. Make linear sovler...

Add type() for LinearSolver and Writer. Refactor Gmres. Make linear sovler public in MshDeform. Minor adjustments
parent f6971d23
No related branches found
No related tags found
1 merge request!6amfe v1.0.5
Pipeline #5578 passed
......@@ -97,9 +97,7 @@ namespace std {
%include "wNode.h"
%nodefault tbox::Element; // Element is pure virtual but swig doesn't know it.
%warnfilter(509); // Shadowed overloaded method MATTYPE. Code compiling and running.
%include "wElement.h"
%warnfilter(+509);
%include "wTag.h"
%include "wQuad4.h"
%include "wTetra4.h"
......
......@@ -94,7 +94,7 @@ void Dss::print(int cod)
void Dss::write(std::ostream &out) const
{
out << "DSS (MKL)" << std::endl;
out << "Intel MKL DSS (PARDISO)" << std::endl;
}
#endif // USE_MKL
......@@ -40,6 +40,7 @@ class TBOX_API Dss : public LinearSolver
public:
Dss();
~Dss();
virtual LSOLTYPE type() const { return LSOLTYPE::DSS; }
#ifndef SWIG
virtual void analyze(Eigen::SparseMatrix<double> const &A) override;
......
......@@ -17,14 +17,14 @@
#include "wGmres.h"
using namespace tbox;
Gmres::Gmres() : LinearSolver()
Gmres::Gmres(int _fill, double _tsh, int _rst, double _tol, int _mit) : LinearSolver(), fill(_fill), tsh(_tsh), rst(_rst), tol(_tol), mit(_mit)
{
// set the preconditioner and solver parameters to "good" starting values
solver.preconditioner().setFillfactor(1);
solver.preconditioner().setDroptol(1e-6);
solver.set_restart(50);
solver.setTolerance(1e-8);
solver.setMaxIterations(1000);
// set the preconditioner and solver parameters
solver.preconditioner().setFillfactor(fill);
solver.preconditioner().setDroptol(tsh);
solver.set_restart(rst);
solver.setTolerance(tol);
solver.setMaxIterations(mit);
}
void Gmres::analyze(Eigen::SparseMatrix<double> const &A)
......@@ -45,7 +45,13 @@ void Gmres::compute(Eigen::SparseMatrix<double> const &A, Eigen::Map<Eigen::Vect
x = solver.solve(b);
}
void Gmres::setTolerance(double t)
{
tol = t;
solver.setTolerance(t);
}
void Gmres::write(std::ostream &out) const
{
out << "GMRES (Eigen)" << std::endl;
out << "Eigen GMRES(" << rst << "," << static_cast<int>(log10(tol)) << ")/ILUT(" << fill << "," << static_cast<int>(log10(tsh)) << ")" << std::endl;
}
......@@ -30,11 +30,18 @@ namespace tbox
*/
class TBOX_API Gmres : public LinearSolver
{
int fill; // fill-in factor for ILUT
double tsh; // drop tolerance threshold for ILUT
int rst; // number of iterations before restart for GMRES
double tol; // relative tolerance for GMRES
int mit; // maximum number of iterations for GMRES
Eigen::GMRES<Eigen::SparseMatrix<double>, Eigen::IncompleteLUT<double>> solver;
public:
Gmres();
Gmres(int _fill = 1, double _tsh = 1e-6, int _rst = 30, double _tol = 1e-8, int _mit = 1000);
Gmres(const Gmres &gmres) : Gmres(gmres.fill, gmres.tsh, gmres.rst, gmres.tol, gmres.mit) {}
virtual ~Gmres() {}
virtual LSOLTYPE type() const override { return LSOLTYPE::GMRES_ILUT; }
#ifndef SWIG
virtual void analyze(Eigen::SparseMatrix<double> const &A) override;
......@@ -45,17 +52,10 @@ public:
virtual double getError() override { return solver.error(); }
virtual int getIterations() override { return static_cast<int>(solver.iterations()); }
void setTolerance(double t);
virtual void write(std::ostream &out) const override;
#endif
void setFillFactor(int f)
{
solver.preconditioner().setFillfactor(f);
}
void setDropTol(double t) { solver.preconditioner().setDroptol(t); }
void setTolerance(double t) { solver.setTolerance(t); }
void setRestart(int r) { solver.set_restart(r); }
void setMaxIterations(int n) { solver.setMaxIterations(n); }
};
} // namespace tbox
......
......@@ -34,6 +34,7 @@ class TBOX_API GmshExport : public MshExport
public:
GmshExport(std::shared_ptr<MshData> _msh);
virtual ~GmshExport();
virtual WRTTYPE type() const override { return WRTTYPE::GMSH; }
virtual void save(std::string const &fname) const override;
#ifndef SWIG
......@@ -45,4 +46,4 @@ public:
} // namespace tbox
#endif //WGMSHUTILS_H
#endif // WGMSHEXPORT_H
......@@ -24,6 +24,20 @@
namespace tbox
{
/**
* @brief Solver type
* @authors Adrien Crovato
*/
enum class LSOLTYPE
{
UNDEFINED = 0,
SPARSE_LU = 1,
GMRES_ILUT = 2,
PARDISO = 3,
DSS = 4,
MUMPS = 5
};
/**
* @brief Base interface class for linear solvers
* @authors Adrien Crovato
......@@ -33,6 +47,7 @@ class TBOX_API LinearSolver : public fwk::wSharedObject
public:
LinearSolver() {}
virtual ~LinearSolver() {}
virtual LSOLTYPE type() const { return LSOLTYPE::UNDEFINED; }
#ifndef SWIG
virtual void analyze(Eigen::SparseMatrix<double> const &A);
......
......@@ -32,10 +32,10 @@ using namespace tbox;
MshDeform::MshDeform(std::shared_ptr<MshData> _msh,
std::shared_ptr<tbox::LinearSolver> _linsol,
int _nDim, int nthrds) : wSharedObject(), linsol(_linsol),
int _nDim, int nthrds) : wSharedObject(),
nDim(_nDim), nthreads(nthrds),
field(false), fixed(false),
moving(false), msh(_msh)
field(false), fixed(false), moving(false),
msh(_msh), linsol(_linsol)
{
// Check problem dimension
if (nDim != 2 && nDim != 3)
......
......@@ -34,7 +34,6 @@ namespace tbox
class TBOX_API MshDeform : public fwk::wSharedObject
{
private:
std::shared_ptr<tbox::LinearSolver> linsol; ///< linear (inner) solver
size_t mshSize; ///< number of nodes in the mesh
int nDim; ///< dimension of the problem (2 or 3)
int nthreads; ///< number of threads
......@@ -62,7 +61,8 @@ private:
Eigen::MatrixXd buildK(tbox::Element const &e, Eigen::MatrixXd const &H);
public:
std::shared_ptr<MshData> msh; ///< mesh
std::shared_ptr<MshData> msh; ///< mesh
std::shared_ptr<tbox::LinearSolver> linsol; ///< linear (inner) solver
MshDeform(std::shared_ptr<MshData> _msh, std::shared_ptr<tbox::LinearSolver> _linsol, int _nDim, int nthrds = 1);
virtual ~MshDeform() { std::cout << "~MshDeform()\n"; }
......@@ -85,4 +85,4 @@ public:
} // namespace tbox
#endif //WMSHDEFORM_H
\ No newline at end of file
#endif // WMSHDEFORM_H
\ No newline at end of file
......@@ -25,6 +25,17 @@
namespace tbox
{
/**
* @brief Writer type
* @authors Adrien Crovato
*/
enum class WRTTYPE
{
UNDEFINED = 0,
GMSH = 1,
VTK = 2
};
/**
* @brief Base class to write mesh
* @authors Adrien Crovato, Romain Boman
......@@ -38,6 +49,7 @@ public:
MshExport(std::shared_ptr<MshData> _msh);
virtual ~MshExport() {}
virtual WRTTYPE type() const { return WRTTYPE::UNDEFINED; }
virtual void save(std::string const &fname) const;
#ifndef SWIG
......@@ -47,4 +59,4 @@ public:
} // namespace tbox
#endif //WMSHEXPORT_H
#endif // WMSHEXPORT_H
......@@ -53,6 +53,7 @@ private:
public:
Mumps();
~Mumps();
virtual LSOLTYPE type() const override { return LSOLTYPE::MUMPS; }
#ifndef SWIG
virtual void analyze(Eigen::SparseMatrix<double> const &A) override;
......
......@@ -50,7 +50,7 @@ void Pardiso::setOption(int k, int v)
void Pardiso::write(std::ostream &out) const
{
out << "Pardiso (MKL/Eigen)" << std::endl;
out << "Eigen Intel MKL PARDISO" << std::endl;
}
#endif // USE_MKL
......@@ -38,6 +38,7 @@ class TBOX_API Pardiso : public LinearSolver
public:
Pardiso() : LinearSolver() {}
virtual ~Pardiso() {}
virtual LSOLTYPE type() const override { return LSOLTYPE::PARDISO; }
#ifndef SWIG
virtual void analyze(Eigen::SparseMatrix<double> const &A) override;
......
......@@ -37,5 +37,5 @@ void SparseLu::compute(Eigen::SparseMatrix<double> const &A, Eigen::Map<Eigen::V
void SparseLu::write(std::ostream &out) const
{
out << "SparseLU (Eigen)" << std::endl;
out << "Eigen SparseLU" << std::endl;
}
......@@ -38,6 +38,7 @@ class TBOX_API SparseLu : public LinearSolver
public:
SparseLu() : LinearSolver() {}
virtual ~SparseLu() {}
virtual LSOLTYPE type() const override { return LSOLTYPE::SPARSE_LU; }
#ifndef SWIG
virtual void analyze(Eigen::SparseMatrix<double> const &A) override;
......
......@@ -44,6 +44,7 @@ class TBOXVTK_API VtkExport : public tbox::MshExport
public:
VtkExport(std::shared_ptr<tbox::MshData> _msh);
virtual ~VtkExport();
virtual tbox::WRTTYPE type() const { return tbox::WRTTYPE::VTK; }
virtual void save(std::string const &fname) const override;
#ifndef SWIG
......@@ -55,4 +56,4 @@ public:
} // namespace tboxVtk
#endif //WGMSHUTILS_H
#endif // WVTKEXPORT_H
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment