File hypre_interface.hxx#

Defines

HypreMalloc(P, SIZE)#
HypreFree(P)#

Functions

BOUT_ENUM_CLASS(HYPRE_SOLVER_TYPE, gmres, bicgstab, pcg)#
namespace bout

Provides access to the Hypre library, handling initialisation and finalisation.

Usage

#include <bout/hyprelib.hxx>

class MyClass { public:

private: HypreLib lib; };

This will then automatically initialise Hypre the first time an object is created, and finalise it when the last object is destroyed.

Copyright 2012 B.D.Dudson, S.Farley, M.V.Umansky, X.Q.Xu

Contact: Ben Dudson, bd512@york.ac.uk

This file is part of BOUT++.

BOUT++ is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.

BOUT++ is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.

You should have received a copy of the GNU Lesser General Public License along with BOUT++. If not, see http://www.gnu.org/licenses/.

Explicit inversion of a 3x3 matrix a

If the matrix is singular (ill conditioned), the determinant is return. Otherwise, an empty std::optional is return

SNB model

template<class T>
class HypreVector#

Public Types

using ind_type = typename T::ind_type#

Public Functions

HypreVector() = default#
inline ~HypreVector()#
HypreVector(const HypreVector<T>&) = delete#
auto operator=(const HypreVector<T>&) = delete#
inline HypreVector(HypreVector<T> &&other)#
inline HypreVector<T> &operator=(HypreVector<T> &&other)#
inline HypreVector<T> &operator=(const T &f)#
inline explicit HypreVector(const T &f, IndexerPtr<T> indConverter)#

Construct from a field, copying over the field values.

inline explicit HypreVector(IndexerPtr<T> indConverter)#

Construct a vector with given index set, but don’t set any values.

inline void assemble()#
inline void writeCacheToHypre()#
inline void readCacheFromHypre()#
inline T toField()#
inline void importValuesFromField(const T &f)#
inline HYPRE_IJVector get()#
inline const HYPRE_IJVector &get() const#
inline HYPRE_ParVector getParallel()#
inline const HYPRE_ParVector &getParallel() const#
inline Element operator()(ind_type &index)#

Private Members

MPI_Comm comm#
HYPRE_BigInt jlower#
HYPRE_BigInt jupper#
HYPRE_BigInt vsize#
HYPRE_IJVector hypre_vector = {nullptr}#
HYPRE_ParVector parallel_vector = {nullptr}#
IndexerPtr<T> indexConverter#
CELL_LOC location#
bool initialised = {false}#
bool have_indices = {false}#
HYPRE_BigInt *I = {nullptr}#
HYPRE_Complex *V = {nullptr}#
HypreLib hyprelib = {}#

Friends

inline friend void swap(HypreVector<T> &lhs, HypreVector<T> &rhs)#
class Element#

Public Functions

Element() = delete#
Element(const Element&) = default#
Element(Element&&) = default#
inline Element(HypreVector<T> &vector_, int index_)#
inline Element &operator=(const Element &other)#
inline Element &operator=(BoutReal value_)#
inline Element &operator+=(BoutReal value_)#
inline operator BoutReal() const#

Private Members

HypreVector<T> *vector#
HYPRE_IJVector hypre_vector = {nullptr}#
HYPRE_BigInt index#
int vec_i#
HYPRE_Complex value#
template<class T>
class HypreMatrix#

Public Types

using ind_type = typename T::ind_type#

Public Functions

HypreMatrix() = default#
HypreMatrix(const HypreMatrix<T>&) = delete#
inline HypreMatrix(HypreMatrix<T> &&other)#
HypreMatrix<T> &operator=(const HypreMatrix<T>&) = delete#
inline HypreMatrix<T> &operator=(HypreMatrix<T> &&other)#
inline explicit HypreMatrix(IndexerPtr<T> indConverter, bool preallocate = true)#

Construct a matrix capable of operating on the specified field, preallocating memory if requeted and possible.

note: preallocate not currently used, but here to match PetscMatrix interface

inline BoutReal getVal(const ind_type &row, const ind_type &column) const#
inline BoutReal getVal(const HYPRE_BigInt row, const HYPRE_BigInt column) const#
inline void setVal(const ind_type &row, const ind_type &column, BoutReal value)#
inline void setVal(const HYPRE_BigInt row, const HYPRE_BigInt column, BoutReal value)#
inline void addVal(const ind_type &row, const ind_type &column, BoutReal value)#
inline void addVal(const HYPRE_BigInt row, const HYPRE_BigInt column, BoutReal value)#
inline BoutReal operator()(const ind_type &row, const ind_type &column) const#
inline BoutReal operator()(const HYPRE_BigInt &row, const HYPRE_BigInt &column) const#
inline Element operator()(const ind_type &row, const ind_type &column)#
inline void assemble()#
inline bool isAssembled() const#
inline void print() const#

Call HYPRE_IJMatrixPrint.

inline HypreMatrix<T> yup(int index = 0)#
inline HypreMatrix<T> ydown(int index = 0)#
inline HypreMatrix<T> ynext(int dir)#
inline HYPRE_IJMatrix get()#
inline const HYPRE_IJMatrix &get() const#
inline HYPRE_ParCSRMatrix getParallel()#
inline const HYPRE_ParCSRMatrix &getParallel() const#
inline void computeAxpby(double alpha, HypreVector<T> &x, double beta, HypreVector<T> &y)#
inline void computeAx(HypreVector<T> &x, HypreVector<T> &y)#

Private Members

MPI_Comm comm#
HYPRE_BigInt ilower#
HYPRE_BigInt iupper#
std::shared_ptr<HYPRE_IJMatrix> hypre_matrix = {nullptr}#
HYPRE_ParCSRMatrix parallel_matrix = {nullptr}#
IndexerPtr<T> index_converter#
CELL_LOC location#
bool initialised = {false}#
int yoffset = {0}#
ParallelTransform *parallel_transform = {nullptr}#
bool assembled = {false}#
HYPRE_BigInt num_rows#
std::vector<HYPRE_BigInt> *I#
std::vector<std::vector<HYPRE_BigInt>> *J#
std::vector<std::vector<HYPRE_Complex>> *V#
HypreLib hyprelib = {}#
class Element#

Public Functions

Element() = delete#
Element(const Element&) = default#
Element(Element&&) = default#
inline Element(HypreMatrix<T> &matrix_, HYPRE_BigInt row_, HYPRE_BigInt column_, std::vector<HYPRE_BigInt> positions_ = {}, std::vector<BoutReal> weights_ = {})#
inline Element &operator=(const Element &other)#
inline Element &operator=(BoutReal value_)#
inline Element &operator+=(BoutReal value_)#
inline operator BoutReal() const#

Private Functions

inline void setValues(BoutReal value_)#
inline void addValues(BoutReal value_)#

Private Members

HypreMatrix *matrix#
HYPRE_IJMatrix hypre_matrix#
HYPRE_BigInt row#
HYPRE_BigInt column#
HYPRE_Complex value#
std::vector<HYPRE_BigInt> positions#
std::vector<BoutReal> weights#
struct MatrixDeleter#

Public Functions

inline void operator()(HYPRE_IJMatrix *matrix) const#
template<class T>
class HypreSystem#

Public Functions

inline HypreSystem(Mesh &mesh, Options &options)#
~HypreSystem() = default#
inline void setRelTol(double tol)#
inline void setAbsTol(double tol)#
inline void setMaxIter(int max_iter)#
inline double getFinalRelResNorm() const#
inline int getNumItersTaken() const#

Return the number of solver iterations taken.

inline int getNumItersTakenAMG() const#

Return the number of BoomerAMG preconditioner iterations taken.

inline void setMatrix(HypreMatrix<T> *A_)#

Set the HypreMatrix to be used in the solver.

inline int setupAMG(HypreMatrix<T> *P_)#

Enable BoomerAMG preconditioner with the given HypreMatrix.

inline void setSolutionVector(HypreVector<T> *x_)#
inline void setRHSVector(HypreVector<T> *b_)#
inline HypreMatrix<T> *getMatrix()#
inline HypreMatrix<T> *getAMGMatrix()#
inline HypreVector<T> *getRHS()#
inline HypreVector<T> *getSolution()#
inline HYPRE_PtrToParSolverFcn SetupFcn() const#
inline HYPRE_PtrToParSolverFcn SolveFcn() const#
inline int solve()#

Private Members

HypreMatrix<T> *P = {nullptr}#
HypreMatrix<T> *A = {nullptr}#
HypreVector<T> *x = {nullptr}#
HypreVector<T> *b = {nullptr}#
MPI_Comm comm#
HYPRE_Solver solver#
HYPRE_Solver precon#
bool solver_setup = {false}#
HYPRE_SOLVER_TYPE solver_type = HYPRE_SOLVER_TYPE::gmres#
HypreLib hyprelib = {}#
decltype(HYPRE_ParCSRGMRESCreate) *solverCreate = {nullptr}#
decltype(HYPRE_ParCSRGMRESSetPrintLevel) *solverSetPrintLevel = {nullptr}#
decltype(HYPRE_ParCSRGMRESSetTol) *solverSetTol = {nullptr}#
decltype(HYPRE_ParCSRGMRESSetAbsoluteTol) *solverSetAbsoluteTol = {nullptr}#
decltype(HYPRE_ParCSRGMRESSetMaxIter) *solverSetMaxIter = {nullptr}#
decltype(HYPRE_ParCSRGMRESGetFinalRelativeResidualNorm) *solverGetFinalRelativeResidualNorm = {nullptr}#
decltype(HYPRE_ParCSRGMRESGetNumIterations) *solverGetNumIterations = {nullptr}#
decltype(HYPRE_ParCSRGMRESSetPrecond) *solverSetPrecond = {nullptr}#
decltype(HYPRE_ParCSRGMRESSetup) *solverSetup = {nullptr}#
decltype(HYPRE_ParCSRGMRESSolve) *solverSolve = {nullptr}#