23 #include "../config/config.hpp"
28 #ifdef MFEM_USE_GINKGO
29 #include <ginkgo/ginkgo.hpp>
64 gko::size_type size,
Vector *mfem_vec,
65 bool ownership =
false)
66 : gko::matrix::Dense<double>(
68 gko::
dim<2> {size, 1},
69 gko::Array<double>::view(exec,
72 exec != exec->get_master() ?
true :
false)),
84 std::function<void(Vector *)>>(
89 using deleter = gko::null_deleter<Vector>;
91 std::function<void(Vector *)>>(
96 static std::unique_ptr<VectorWrapper>
create(
97 std::shared_ptr<const gko::Executor> exec,
100 bool ownership =
false)
102 return std::unique_ptr<VectorWrapper>(
114 virtual std::unique_ptr<gko::matrix::Dense<double>>
119 this->
wrapped_vec.get()->GetMemory().GetMemoryType());
137 std::shared_ptr<const gko::Executor> exec,
138 const gko::dim<2> &size,
139 gko::size_type stride)
const override
146 "VectorWrapper cannot be created with stride > 1");
149 gko::size_type total_size = size[0]*size[1];
152 this->
wrapped_vec.get()->GetMemory().GetMemoryType());
160 this->get_executor(), total_size, mfem_vec,
167 const gko::span &rows,
168 const gko::span &columns,
169 const gko::size_type stride)
override
172 gko::size_type num_rows = rows.end - rows.begin;
173 gko::size_type num_cols = columns.end - columns.begin;
177 if (num_cols > 1 || stride > 1)
179 throw gko::BadDimension(
180 __FILE__, __LINE__, __func__,
"new_submatrix", num_rows,
182 "VectorWrapper submatrix must have one column and stride = 1");
184 int data_size =
static_cast<int>(num_rows * num_cols);
185 int start =
static_cast<int>(rows.begin);
196 this->get_executor(), data_size, mfem_vec,
201 std::unique_ptr<Vector, std::function<void(Vector *)>>
wrapped_vec;
213 :
public gko::EnableLinOp<OperatorWrapper>,
214 public gko::EnableCreateMethod<OperatorWrapper>
218 gko::size_type size = 0,
221 gko::EnableCreateMethod<OperatorWrapper>()
223 this->wrapped_oper = oper;
227 void apply_impl(
const gko::LinOp *
b, gko::LinOp *x)
const override;
229 const gko::LinOp *beta, gko::LinOp *x)
const override;
237 template <
typename ValueType=
double>
238 double get_norm(
const gko::matrix::Dense<ValueType> *norm)
241 auto cpu_norm = clone(norm->get_executor()->get_master(), norm);
243 return cpu_norm->at(0, 0);
248 template <
typename ValueType=
double>
252 auto exec = b->get_executor();
254 auto b_norm = gko::initialize<gko::matrix::Dense<ValueType>>({0.0}, exec);
256 b->compute_norm2(lend(b_norm));
258 return std::pow(
get_norm(lend(b_norm)),2);
272 template <
typename ValueType=
double>
279 if (compute_real_residual)
281 mfem::out <<
"Iteration log with real residual norms:" << std::endl;
285 mfem::out <<
"Iteration log with residual norms:" << std::endl;
287 mfem::out <<
'|' << std::setw(10) <<
"Iteration" <<
'|' << std::setw(25)
288 <<
"Residual Norm" <<
'|' << std::endl;
291 mfem::out <<
'|' << std::setfill(
'-') << std::setw(11) <<
'|' <<
292 std::setw(26) <<
'|' << std::setfill(
' ') << std::endl;
295 for (std::size_t i = 0; i < iterations.size(); i++)
297 mfem::out <<
'|' << std::setw(10) << iterations[i] <<
'|'
298 << std::setw(25) << residual_norms[i] <<
'|' << std::endl;
302 mfem::out.unsetf(std::ios_base::floatfield);
304 mfem::out <<
'|' << std::setfill(
'-') << std::setw(11) <<
'|' <<
305 std::setw(26) <<
'|' << std::setfill(
' ') << std::endl;
313 const gko::size_type &iteration,
314 const gko::LinOp *residual,
315 const gko::LinOp *solution,
316 const gko::LinOp *residual_norm,
317 const gko::LinOp *implicit_sq_residual_norm)
const override
321 if (solution && compute_real_residual)
324 auto exec = matrix->get_executor();
327 matrix->apply(gko::lend(solution), gko::lend(res));
330 if (dynamic_cast<const VectorWrapper*>(b))
332 const VectorWrapper *b_cast = gko::as<const VectorWrapper>(b);
340 auto neg_one = gko::initialize<gko_dense>({-1.0}, exec);
341 res->add_scaled(gko::lend(neg_one), gko::lend(b));
351 if (implicit_sq_residual_norm)
353 auto dense_norm = gko::as<gko_dense>(implicit_sq_residual_norm);
355 residual_norms.push_back(
get_norm(dense_norm));
358 else if (residual_norm)
360 auto dense_norm = gko::as<gko_dense>(residual_norm);
362 residual_norms.push_back(
get_norm(dense_norm));
367 auto dense_residual = gko::as<gko_dense>(residual);
371 residual_norms.push_back(norm);
375 iterations.push_back(iteration);
380 const gko::size_type &iteration,
381 const gko::LinOp *residual,
382 const gko::LinOp *solution,
383 const gko::LinOp *residual_norm)
const override
391 const gko::LinOp *matrix,
const gko_dense *b,
392 bool compute_real_residual=
false)
393 : gko::log::Logger(exec,
394 gko::log::Logger::iteration_complete_mask),
397 compute_real_residual{compute_real_residual}
399 if (compute_real_residual ==
true)
401 if (dynamic_cast<const VectorWrapper*>(b))
403 const VectorWrapper *b_cast = gko::as<const VectorWrapper>(b);
404 res = std::move(gko_dense::create_with_config_of(b_cast).release());
408 res = std::move(gko::clone(b).release());
415 const gko::LinOp *matrix;
421 mutable std::vector<ValueType> residual_norms{};
423 mutable std::vector<std::size_t> iterations{};
426 const bool compute_real_residual;
474 return this->executor;
478 std::shared_ptr<gko::Executor> executor;
681 std::shared_ptr<gko::stop::ResidualNorm<>::Factory>
688 std::shared_ptr<gko::stop::ResidualNorm<>::Factory>
695 std::shared_ptr<gko::stop::ImplicitResidualNorm<>::Factory>
702 std::shared_ptr<gko::stop::ImplicitResidualNorm<>::Factory>
753 initialize_ginkgo_log(gko::matrix::Dense<double>*
b)
const;
759 std::shared_ptr<gko::LinOp> system_oper;
769 template<
typename SolverType>
780 gko::as<typename SolverType::Factory>(
solver_gen)->get_parameters().criteria =
792 gko::as<typename SolverType::Factory>(
solver_gen)->get_parameters().criteria =
804 gko::as<typename SolverType::Factory>(
solver_gen)->get_parameters().criteria =
968 using gko::solver::cb_gmres::storage_precision;
998 storage_precision prec = storage_precision::reduce1);
1017 storage_precision prec = storage_precision::reduce1);
1078 const std::string &storage_opt =
"none",
1079 const double accuracy = 1.e-1,
1080 const int max_block_size = 32
1112 const std::string &factorization_type =
"exact",
1113 const int sweeps = 0,
1114 const bool skip_sort =
false
1154 const std::string &factorization_type =
"exact",
1155 const int sweeps = 0,
1156 const int sparsity_power = 1,
1157 const bool skip_sort =
false
1189 const std::string &factorization_type =
"exact",
1190 const int sweeps = 0,
1191 const bool skip_sort =
false
1231 const std::string &factorization_type =
"exact",
1232 const int sweeps = 0,
1233 const int sparsity_power = 1,
1234 const bool skip_sort =
false
1254 const Solver &mfem_precond
1264 MFEM_ABORT(
"Ginkgo::MFEMPreconditioner must be constructed "
1265 "with the MFEM Operator that it will wrap as an argument;\n"
1266 "calling SetOperator() is not allowed.");
1273 #endif // MFEM_GINKGO
bool use_implicit_res_norm
std::shared_ptr< gko::log::Convergence<> > convergence_logger
std::shared_ptr< gko::Executor > GetExecutor() const
void SetRelTol(double rtol)
std::shared_ptr< ResidualLogger<> > residual_logger
virtual ~GinkgoExecutor()=default
virtual void SetOperator(const Operator &op)
std::shared_ptr< gko::stop::ResidualNorm<>::Factory > abs_criterion
GinkgoIterativeSolver(GinkgoExecutor &exec, bool use_implicit_res_norm)
CBGMRESSolver(GinkgoExecutor &exec, int dim=0, storage_precision prec=storage_precision::reduce1)
std::shared_ptr< gko::LinOpFactory > precond_gen
CGSSolver(GinkgoExecutor &exec)
virtual std::unique_ptr< gko::matrix::Dense< double > > create_with_same_config() const override
MFEMPreconditioner(GinkgoExecutor &exec, const Solver &mfem_precond)
static std::unique_ptr< VectorWrapper > create(std::shared_ptr< const gko::Executor > exec, gko::size_type size, Vector *mfem_vec, bool ownership=false)
JacobiPreconditioner(GinkgoExecutor &exec, const std::string &storage_opt="none", const double accuracy=1.e-1, const int max_block_size=32)
std::shared_ptr< gko::Executor > executor
GinkgoExecutor(ExecType exec_type)
virtual void SetOperator(const Operator &op)
int GetNumIterations() const
void SetMaxIter(int max_it)
gko::matrix::Dense< ValueType > gko_dense
virtual void UseDevice(bool use_dev) const
Enable execution of Vector operations using the mfem::Device.
BICGSTABSolver(GinkgoExecutor &exec)
virtual void Mult(const Vector &x, Vector &y) const
std::shared_ptr< gko::stop::ImplicitResidualNorm<>::Factory > imp_rel_criterion
bool HasGeneratedPreconditioner() const
IluPreconditioner(GinkgoExecutor &exec, const std::string &factorization_type="exact", const int sweeps=0, const bool skip_sort=false)
virtual std::unique_ptr< gko::matrix::Dense< double > > create_submatrix_impl(const gko::span &rows, const gko::span &columns, const gko::size_type stride) override
const std::shared_ptr< gko::LinOpFactory > GetFactory() const
Vector & get_mfem_vec_ref()
bool UsesVectorWrappers() const
std::shared_ptr< gko::Executor > executor
void on_iteration_complete(const gko::LinOp *op, const gko::size_type &iteration, const gko::LinOp *residual, const gko::LinOp *solution, const gko::LinOp *residual_norm) const override
double GetFinalNorm() const
bool sub_op_needs_wrapped_vecs
void update_stop_factory()
CGSolver(GinkgoExecutor &exec)
void SetPrintLevel(int print_lvl)
virtual std::unique_ptr< gko::matrix::Dense< double > > create_with_type_of_impl(std::shared_ptr< const gko::Executor > exec, const gko::dim< 2 > &size, gko::size_type stride) const override
ResidualLogger(std::shared_ptr< const gko::Executor > exec, const gko::LinOp *matrix, const gko_dense *b, bool compute_real_residual=false)
VectorWrapper(std::shared_ptr< const gko::Executor > exec, gko::size_type size, Vector *mfem_vec, bool ownership=false)
GMRESSolver(GinkgoExecutor &exec, int dim=0)
virtual ~GinkgoPreconditioner()=default
const std::shared_ptr< gko::LinOpFactory > GetFactory() const
double compute_norm(const gko::matrix::Dense< ValueType > *b)
std::shared_ptr< gko::stop::ResidualNorm<>::Factory > rel_criterion
double get_norm(const gko::matrix::Dense< ValueType > *norm)
IRSolver(GinkgoExecutor &exec)
const Vector & get_mfem_vec_const_ref() const
IcPreconditioner(GinkgoExecutor &exec, const std::string &factorization_type="exact", const int sweeps=0, const bool skip_sort=false)
std::shared_ptr< gko::stop::Combined::Factory > combined_factory
std::shared_ptr< gko::LinOpFactory > solver_gen
bool has_generated_precond
FCGSolver(GinkgoExecutor &exec)
IcIsaiPreconditioner(GinkgoExecutor &exec, const std::string &factorization_type="exact", const int sweeps=0, const int sparsity_power=1, const bool skip_sort=false)
void on_iteration_complete(const gko::LinOp *, const gko::size_type &iteration, const gko::LinOp *residual, const gko::LinOp *solution, const gko::LinOp *residual_norm, const gko::LinOp *implicit_sq_residual_norm) const override
virtual void Mult(const Vector &x, Vector &y) const
IluIsaiPreconditioner(GinkgoExecutor &exec, const std::string &factorization_type="exact", const int sweeps=0, const int sparsity_power=1, const bool skip_sort=false)
void SetAbsTol(double atol)
GinkgoPreconditioner(GinkgoExecutor &exec)
OperatorWrapper(std::shared_ptr< const gko::Executor > exec, gko::size_type size=0, const Operator *oper=NULL)
virtual ~GinkgoIterativeSolver()=default
std::shared_ptr< gko::LinOp > generated_precond
EnableGinkgoSolver(GinkgoExecutor &exec, bool use_implicit_res_norm)
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
std::shared_ptr< gko::stop::ImplicitResidualNorm<>::Factory > imp_abs_criterion
virtual void SetOperator(const Operator &op)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
const std::shared_ptr< gko::LinOp > GetGeneratedPreconditioner() const
void apply_impl(const gko::LinOp *b, gko::LinOp *x) const override
std::shared_ptr< gko::LinOp > solver
void operator()(pointer ptr) const noexcept