Skip to content

Commit

Permalink
added is dense flag
Browse files Browse the repository at this point in the history
  • Loading branch information
teseoch committed Oct 19, 2023
1 parent 66a559f commit 33c61b0
Show file tree
Hide file tree
Showing 8 changed files with 25 additions and 3 deletions.
2 changes: 2 additions & 0 deletions src/polysolve/linear/CuSolverDN.cuh
Original file line number Diff line number Diff line change
Expand Up @@ -42,6 +42,8 @@ namespace polysolve::linear
// Factorize system matrix (dense, preferred)
virtual void factorize_dense(const Eigen::MatrixXd &A) override;

bool is_dense() const override { return true; }

// Solve the linear system Ax = b
virtual void solve(const Ref<const VectorXd> b, Ref<VectorXd> x) override;

Expand Down
2 changes: 2 additions & 0 deletions src/polysolve/linear/EigenSolver.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -97,6 +97,8 @@ namespace polysolve::linear
// Constructor requires a solver name used for finding parameters in the json file passed to set_parameters
EigenDenseSolver(const std::string &name) { m_Name = name; }

bool is_dense() const override { return true; }

public:
// Get info on the last solve step
virtual void get_info(json &params) const override;
Expand Down
1 change: 1 addition & 0 deletions src/polysolve/linear/FEMSolver.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -103,6 +103,7 @@ namespace polysolve::linear
const bool remove_zero_cols,
const bool skip_last_cols)
{
assert(!solver.is_dense());
// Let Γ be the set of Dirichlet dofs.
// To implement nonzero Dirichlet boundary conditions, we seek to replace
// the linear system Au = f with a new system Ãx = g, where
Expand Down
3 changes: 3 additions & 0 deletions src/polysolve/linear/Solver.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -91,6 +91,9 @@ namespace polysolve::linear
// Factorize system matrix of a dense matrix
virtual void factorize_dense(const Eigen::MatrixXd &A) {}

// If solver uses dense matrices
virtual bool is_dense() const { return false; }

//
// @brief { Solve the linear system Ax = b }
//
Expand Down
1 change: 1 addition & 0 deletions src/polysolve/nonlinear/BFGS.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -12,6 +12,7 @@ namespace polysolve::nonlinear
: Superclass(solver_params, characteristic_length, logger)
{
linear_solver = polysolve::linear::Solver::create(linear_solver_params, logger);
assert(linear_solver->is_dense());
}

std::string BFGS::descent_strategy_name(int descent_strategy) const
Expand Down
1 change: 1 addition & 0 deletions src/polysolve/nonlinear/DenseNewton.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -14,6 +14,7 @@ namespace polysolve::nonlinear
logger)
{
linear_solver = polysolve::linear::Solver::create(linear_solver_params, logger);
assert(linear_solver->is_dense());
}

double DenseNewton::solve_linear_system(Problem &objFunc,
Expand Down
1 change: 1 addition & 0 deletions src/polysolve/nonlinear/SparseNewton.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -16,6 +16,7 @@ namespace polysolve::nonlinear
logger)
{
linear_solver = polysolve::linear::Solver::create(linear_solver_params, logger);
assert(!linear_solver->is_dense());
}

double SparseNewton::solve_linear_system(Problem &objFunc,
Expand Down
17 changes: 14 additions & 3 deletions tests/test_linear_solver.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -106,6 +106,8 @@ TEST_CASE("all", "[solver]")
Eigen::SparseMatrix<double> A;
const bool ok = loadMarket(A, path + "/A_2.mat");
REQUIRE(ok);
json solver_info;
Eigen::MatrixXd A_dense(A);

auto solvers = Solver::available_solvers();
for (const auto &s : solvers)
Expand All @@ -130,13 +132,22 @@ TEST_CASE("all", "[solver]")
Eigen::VectorXd x(b.size());
x.setZero();

solver->analyze_pattern(A, A.rows());
solver->factorize(A);
if (solver->is_dense())
{
solver->analyze_pattern_dense(A, A.rows());
solver->factorize_dense(A);
}
else
{
solver->analyze_pattern(A, A.rows());
solver->factorize(A);
}

solver->solve(b, x);

REQUIRE(solver->name() == s);

// solver->get_info(solver_info);
solver->get_info(solver_info);

// std::cout<<"Solver error: "<<x<<std::endl;
const double err = (A * x - b).norm();
Expand Down

0 comments on commit 33c61b0

Please sign in to comment.