File petsc_interface.hxx#

Classes to wrap PETSc matrices and vectors, providing a convenient interface to them. In particular, they will internally convert between BOUT++ indices and PETSc ones, making it far easier to set up a linear system.

Copyright 2019 C. MacMackin

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/.

Functions

template<class T>
void swap(PetscVector<T> &first, PetscVector<T> &second)#

Move reference to one Vec from first to second and vice versa.

template<class T>
inline MPI_Comm getComm([[maybe_unused]] const T &field)#
template<>
inline MPI_Comm getComm([[maybe_unused]] const FieldPerp &field)#
template<class T>
void swap(PetscMatrix<T> &first, PetscMatrix<T> &second)#

Move reference to one Mat from first to second and vice versa.

template<class T>
PetscVector<T> operator*(const PetscMatrix<T> &mat, const PetscVector<T> &vec)#

Performs matrix-multiplication on the supplied vector

template<class T>
class PetscVector#

A class which wraps PETSc vector objects, allowing them to be indexed using the BOUT++ scheme. Note that boundaries are only included in the vector to a depth of 1.

Public Types

using ind_type = typename T::ind_type#

Public Functions

inline PetscVector()#
~PetscVector() = default#
inline PetscVector(const PetscVector<T> &vec)#
inline PetscVector(PetscVector<T> &&vec) noexcept#
inline PetscVector(const T &field, IndexerPtr<T> indConverter)#

Construct from a field, copying over the field values.

inline PetscVector(const PetscVector<T> &other, Vec *vec)#

Construct a vector like v, but using data from a raw PETSc Vec. That Vec (not a copy) will then be owned by the new object.

inline PetscVector<T> &operator=(const PetscVector<T> &rhs)#

Copy assignment.

inline PetscVector<T> &operator=(PetscVector<T> &&rhs) noexcept#

Move assignment.

inline PetscVector<T> &operator=(const T &field)#
inline BoutReal &operator()(const ind_type &index)#
inline BoutReal operator()(const ind_type &index) const#
inline void assemble()#
inline void destroy()#
inline T toField() const#

Returns a field constructed from the contents of this vector.

inline Vec *get()#

Provides a reference to the raw PETSc Vec object.

inline const Vec *get() const#

Private Members

PetscLib lib = {}#
std::unique_ptr<Vec, VectorDeleter> vector = nullptr#
IndexerPtr<T> indexConverter = {}#
CELL_LOC location = CELL_LOC::deflt#
bool initialised = false#
Array<BoutReal> vector_values = {}#

Friends

friend void swap(PetscVector<T> &first, PetscVector<T> &second)#

Move reference to one Vec from first to second and vice versa.

struct VectorDeleter#

Public Functions

inline void operator()(Vec *vec) const#
template<class T>
class PetscMatrix#

A class which wraps PETSc vector objects, allowing them to be indexed using the BOUT++ scheme. It provides the option of setting a y-offset that interpolates onto field lines.

Public Types

using ind_type = typename T::ind_type#

Public Functions

inline PetscMatrix()#

Default constructor does nothing.

~PetscMatrix() = default#
inline PetscMatrix(const PetscMatrix<T> &mat)#

Copy constructor.

inline PetscMatrix(PetscMatrix<T> &&mat) noexcept#

Move constrcutor.

inline PetscMatrix(IndexerPtr<T> indConverter, bool preallocate = true)#
inline PetscMatrix<T> &operator=(PetscMatrix<T> rhs)#

Copy assignment.

inline PetscMatrix<T> &operator=(PetscMatrix<T> &&rhs) noexcept#

Move assignment.

inline Element operator()(const ind_type &index1, const ind_type &index2)#
inline BoutReal operator()(const ind_type &index1, const ind_type &index2) const#
inline void assemble()#
inline void partialAssemble()#
inline void destroy()#
inline PetscMatrix<T> yup(int index = 0)#
inline PetscMatrix<T> ydown(int index = 0)#
inline PetscMatrix<T> ynext(int dir)#
inline Mat *get()#

Provides a reference to the raw PETSc Mat object.

inline const Mat *get() const#

Private Members

PetscLib lib#
std::shared_ptr<Mat> matrix = nullptr#
IndexerPtr<T> indexConverter = {}#
ParallelTransform *pt = {}#
int yoffset = 0#
bool initialised = false#

Friends

friend void swap(PetscMatrix<T> &first, PetscMatrix<T> &second)#

Move reference to one Mat from first to second and vice versa.

class Element#
#include <petsc_interface.hxx>

A class which is used to assign to a particular element of a PETSc matrix, potentially with a y-offset. It is meant to be transient and will be destroyed immediately after use. In general you should not try to assign an instance to a variable.

The Element object will store a copy of the value which has been assigned to it. This is because mixing calls to the PETSc function MatGetValues() and MatSetValues() without intermediate matrix assembly will cause errors. Thus, if the user wishes to get the value of the matrix, it must be stored here.

Public Functions

Element() = delete#
~Element() = default#
Element(Element&&) noexcept = default#
Element &operator=(Element&&) noexcept = default#
Element(const Element &other) = default#
inline Element(Mat *matrix, PetscInt row, PetscInt col, std::vector<PetscInt> position = {}, std::vector<BoutReal> weight = {})#
inline Element &operator=(const Element &other)#
inline Element &operator=(BoutReal val)#
inline Element &operator+=(BoutReal val)#
inline operator BoutReal() const#

Private Functions

inline void setValues(BoutReal val, InsertMode mode)#

Private Members

Mat *petscMatrix#
PetscInt petscRow#
PetscInt petscCol#
PetscScalar value#
std::vector<PetscInt> positions#
std::vector<BoutReal> weights#
struct MatrixDeleter#

Public Functions

inline void operator()(Mat *mat) const#