Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

[SofaPreconditioner] FIX ShewchukPCGLinearSolver #737

Merged
merged 3 commits into from
Jul 31, 2018
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
Original file line number Diff line number Diff line change
Expand Up @@ -61,6 +61,7 @@ ShewchukPCGLinearSolver<TMatrix,TVector>::ShewchukPCGLinearSolver()
, f_build_precond( initData(&f_build_precond,true,"build_precond","Build the preconditioners, if false build the preconditioner only at the initial step") )
, f_preconditioners( initData(&f_preconditioners, "preconditioners", "If not empty: path to the solvers to use as preconditioners") )
, f_graph( initData(&f_graph,"graph","Graph of residuals at each iteration") )
, m_preconditioners(0)
{
f_graph.setWidget("graph");
// f_graph.setReadOnly(true);
Expand All @@ -74,10 +75,10 @@ void ShewchukPCGLinearSolver<TMatrix,TVector>::init()
if (! f_preconditioners.getValue().empty()) {
BaseContext * c = this->getContext();

c->get(preconditioners, f_preconditioners.getValue());
c->get(m_preconditioners, f_preconditioners.getValue());

if (preconditioners) sout << "Found " << f_preconditioners.getValue() << sendl;
else serr << "Solver \"" << f_preconditioners.getValue() << "\" not found." << sendl;
if (m_preconditioners) msg_info() << "Found " << f_preconditioners.getValue();
else msg_error() << "Solver \"" << f_preconditioners.getValue() << "\" not found.";
}

first = true;
Expand All @@ -93,10 +94,10 @@ void ShewchukPCGLinearSolver<TMatrix,TVector>::setSystemMBKMatrix(const core::Me

sofa::helper::AdvancedTimer::stepEnd("PCG::setSystemMBKMatrix");

if (preconditioners==NULL) return;
if (m_preconditioners==NULL) return;

if (first) { //We initialize all the preconditioners for the first step
preconditioners->setSystemMBKMatrix(mparams);
m_preconditioners->setSystemMBKMatrix(mparams);
first = false;
next_refresh_step = 1;

Expand All @@ -106,7 +107,7 @@ void ShewchukPCGLinearSolver<TMatrix,TVector>::setSystemMBKMatrix(const core::Me

if (f_update_step.getValue()>0) {
if (next_refresh_step>=f_update_step.getValue()) {
preconditioners->setSystemMBKMatrix(mparams);
m_preconditioners->setSystemMBKMatrix(mparams);
next_refresh_step=1;
} else {
next_refresh_step++;
Expand All @@ -116,7 +117,7 @@ void ShewchukPCGLinearSolver<TMatrix,TVector>::setSystemMBKMatrix(const core::Me
sofa::helper::AdvancedTimer::stepEnd("PCG::PrecondSetSystemMBKMatrix");
}

preconditioners->updateSystemMatrix();
m_preconditioners->updateSystemMatrix();
}

template<>
Expand Down Expand Up @@ -161,7 +162,7 @@ void ShewchukPCGLinearSolver<TMatrix,TVector>::solve (Matrix& M, Vector& x, Vect
Vector& w = *vtmp.createTempVector();
Vector& s = *vtmp.createTempVector();

bool apply_precond = preconditioners!=NULL && f_use_precond.getValue();
bool apply_precond = m_preconditioners!=NULL && f_use_precond.getValue();

double b_norm = b.dot(b);
double tol = f_tolerance.getValue() * b_norm;
Expand All @@ -171,9 +172,9 @@ void ShewchukPCGLinearSolver<TMatrix,TVector>::solve (Matrix& M, Vector& x, Vect

if (apply_precond) {
sofa::helper::AdvancedTimer::stepBegin("PCGLinearSolver::apply Precond");
preconditioners->setSystemLHVector(w);
preconditioners->setSystemRHVector(r);
preconditioners->solveSystem();
m_preconditioners->setSystemLHVector(w);
m_preconditioners->setSystemRHVector(r);
m_preconditioners->solveSystem();

sofa::helper::AdvancedTimer::stepEnd("PCGLinearSolver::apply Precond");
} else {
Expand All @@ -194,9 +195,9 @@ void ShewchukPCGLinearSolver<TMatrix,TVector>::solve (Matrix& M, Vector& x, Vect

if (apply_precond) {
sofa::helper::AdvancedTimer::stepBegin("PCGLinearSolver::apply Precond");
preconditioners->setSystemLHVector(s);
preconditioners->setSystemRHVector(r);
preconditioners->solveSystem();
m_preconditioners->setSystemLHVector(s);
m_preconditioners->setSystemRHVector(r);
m_preconditioners->solveSystem();

sofa::helper::AdvancedTimer::stepEnd("PCGLinearSolver::apply Precond");
} else {
Expand Down Expand Up @@ -229,7 +230,6 @@ SOFA_DECL_CLASS(ShewchukPCGLinearSolver)
int ShewchukPCGLinearSolverClass = core::RegisterObject("Linear system solver using the conjugate gradient iterative algorithm")
.add< ShewchukPCGLinearSolver<GraphScatteredMatrix,GraphScatteredVector> >(true)
.addAlias("PCGLinearSolver");
;

} // namespace linearsolver

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -70,7 +70,7 @@ class ShewchukPCGLinearSolver : public sofa::component::linearsolver::MatrixLine

private :
unsigned next_refresh_step;
sofa::core::behavior::LinearSolver* preconditioners;
sofa::core::behavior::LinearSolver* m_preconditioners;
bool first;
int newton_iter;

Expand Down