12#include <boost/lexical_cast.hpp> 
   16#include <petscoptions.h> 
   35void error(
int error_code, std::string filename, std::string petsc_function);
 
   46               const std::vector<std::span<const PetscScalar>>& x);
 
   63Vec 
create_vector(MPI_Comm comm, std::array<std::int64_t, 2> range,
 
   64                  std::span<const std::int64_t> ghosts, 
int bs);
 
   76                       std::span<const PetscScalar> x);
 
  101        std::pair<std::reference_wrapper<const common::IndexMap>, 
int>>& maps);
 
  107        std::pair<std::reference_wrapper<const common::IndexMap>, 
int>>& maps);
 
  111    Vec x, 
const std::vector<std::span<const PetscScalar>>& x_b,
 
  113        std::pair<std::reference_wrapper<const common::IndexMap>, 
int>>& maps);
 
  118                  const std::string& type = std::string());
 
  135void set(std::string option);
 
  139void set(std::string option, 
const T value)
 
  141  if (option[0] != 
'-')
 
  142    option = 
'-' + option;
 
  145  ierr = PetscOptionsSetValue(
nullptr, option.c_str(),
 
  146                              boost::lexical_cast<std::string>(value).c_str());
 
 
  152void clear(std::string option);
 
 
  189  Vector(Vec x, 
bool inc_ref_count);
 
  205  std::int64_t 
size() 
const;
 
  214  MPI_Comm 
comm() 
const;
 
 
  240  Operator(Mat A, 
bool inc_ref_count);
 
  259  std::array<std::int64_t, 2> 
size() 
const;
 
 
  288  static auto set_fn(Mat A, InsertMode mode)
 
  290    return [A, mode, cache = std::vector<PetscInt>()](
 
  291               const std::span<const std::int32_t>& rows,
 
  292               const std::span<const std::int32_t>& cols,
 
  293               const std::span<const PetscScalar>& vals) 
mutable -> 
int 
  296#ifdef PETSC_USE_64BIT_INDICES 
  297      cache.resize(rows.size() + cols.size());
 
  298      std::copy(rows.begin(), rows.end(), cache.begin());
 
  299      std::copy(cols.begin(), cols.end(),
 
  300                std::next(cache.begin(), rows.size()));
 
  301      const PetscInt* _rows = cache.data();
 
  302      const PetscInt* _cols = cache.data() + rows.size();
 
  303      ierr = MatSetValuesLocal(A, rows.size(), _rows, cols.size(), _cols,
 
  306      ierr = MatSetValuesLocal(A, rows.size(), rows.data(), cols.size(),
 
  307                               cols.data(), vals.data(), mode);
 
 
  325    return [A, mode, cache = std::vector<PetscInt>()](
 
  326               const std::span<const std::int32_t>& rows,
 
  327               const std::span<const std::int32_t>& cols,
 
  328               const std::span<const PetscScalar>& vals) 
mutable -> 
int 
  331#ifdef PETSC_USE_64BIT_INDICES 
  332      cache.resize(rows.size() + cols.size());
 
  333      std::copy(rows.begin(), rows.end(), cache.begin());
 
  334      std::copy(cols.begin(), cols.end(),
 
  335                std::next(cache.begin(), rows.size()));
 
  336      const PetscInt* _rows = cache.data();
 
  337      const PetscInt* _cols = cache.data() + rows.size();
 
  338      ierr = MatSetValuesBlockedLocal(A, rows.size(), _rows, cols.size(), _cols,
 
  341      ierr = MatSetValuesBlockedLocal(A, rows.size(), rows.data(), cols.size(),
 
  342                                      cols.data(), vals.data(), mode);
 
  347        petsc::error(ierr, __FILE__, 
"MatSetValuesBlockedLocal");
 
 
  363    return [A, bs0, bs1, mode, cache0 = std::vector<PetscInt>(),
 
  364            cache1 = std::vector<PetscInt>()](
 
  365               const std::span<const std::int32_t>& rows,
 
  366               const std::span<const std::int32_t>& cols,
 
  367               const std::span<const PetscScalar>& vals) 
mutable -> 
int 
  370      cache0.resize(bs0 * rows.size());
 
  371      cache1.resize(bs1 * cols.size());
 
  372      for (std::size_t i = 0; i < rows.size(); ++i)
 
  373        for (
int k = 0; k < bs0; ++k)
 
  374          cache0[bs0 * i + k] = bs0 * rows[i] + k;
 
  376      for (std::size_t i = 0; i < cols.size(); ++i)
 
  377        for (
int k = 0; k < bs1; ++k)
 
  378          cache1[bs1 * i + k] = bs1 * cols[i] + k;
 
  380      ierr = MatSetValuesLocal(A, cache0.size(), cache0.data(), cache1.size(),
 
  381                               cache1.data(), vals.data(), mode);
 
 
  392         const std::string& type = std::string());
 
  398  Matrix(Mat A, 
bool inc_ref_count);
 
 
  485  int solve(Vec x, 
const Vec b, 
bool transpose = 
false) 
const;
 
 
 
This class represents the distribution index arrays across processes. An index array is a contiguous ...
Definition IndexMap.h:64
 
Sparsity pattern data structure that can be used to initialize sparse matrices. After assembly,...
Definition SparsityPattern.h:26
 
Distributed vector.
Definition Vector.h:31
 
constexpr int bs() const
Get block size.
Definition Vector.h:186
 
std::span< const value_type > array() const
Get local part of the vector (const version)
Definition Vector.h:189
 
std::shared_ptr< const common::IndexMap > index_map() const
Get IndexMap.
Definition Vector.h:183
 
This class implements Krylov methods for linear systems of the form Ax = b. It is a wrapper for the K...
Definition petsc.h:451
 
void set_operator(const Mat A)
Set operator (Mat)
Definition petsc.cpp:673
 
KSP ksp() const
Return PETSc KSP pointer.
Definition petsc.cpp:789
 
std::string get_options_prefix() const
Returns the prefix used by PETSc when searching the PETSc options database.
Definition petsc.cpp:771
 
void set_operators(const Mat A, const Mat P)
Set operator and preconditioner matrix (Mat)
Definition petsc.cpp:675
 
void set_from_options() const
Set options from PETSc options database.
Definition petsc.cpp:781
 
~KrylovSolver()
Destructor.
Definition petsc.cpp:661
 
int solve(Vec x, const Vec b, bool transpose=false) const
Solve linear system Ax = b and return number of iterations (A^t x = b if transpose is true)
Definition petsc.cpp:684
 
void set_options_prefix(std::string options_prefix)
Sets the prefix used by PETSc when searching the PETSc options database.
Definition petsc.cpp:762
 
It is a simple wrapper for a PETSc matrix pointer (Mat). Its main purpose is to assist memory managem...
Definition petsc.h:282
 
static auto set_block_fn(Mat A, InsertMode mode)
Return a function with an interface for adding or inserting values into the matrix A using blocked in...
Definition petsc.h:323
 
void set_from_options()
Call PETSc function MatSetFromOptions on the PETSc Mat object.
Definition petsc.cpp:629
 
double norm(Norm norm_type) const
Return norm of matrix.
Definition petsc.cpp:573
 
static auto set_fn(Mat A, InsertMode mode)
Return a function with an interface for adding or inserting values into the matrix A (calls MatSetVal...
Definition petsc.h:288
 
Matrix(Matrix &&A)=default
Move constructor (falls through to base class move constructor)
 
std::string get_options_prefix() const
Returns the prefix used by PETSc when searching the options database.
Definition petsc.cpp:621
 
Matrix & operator=(Matrix &&A)=default
Move assignment operator.
 
static auto set_block_expand_fn(Mat A, int bs0, int bs1, InsertMode mode)
Return a function with an interface for adding or inserting blocked values to the matrix A using non-...
Definition petsc.h:361
 
~Matrix()=default
Destructor.
 
Matrix & operator=(const Matrix &A)=delete
Assignment operator (deleted)
 
AssemblyType
Assembly type FINAL - corresponds to PETSc MAT_FINAL_ASSEMBLY FLUSH - corresponds to PETSc MAT_FLUSH_...
Definition petsc.h:419
 
void apply(AssemblyType type)
Finalize assembly of tensor. The following values are recognized for the mode parameter:
Definition petsc.cpp:598
 
void set_options_prefix(std::string options_prefix)
Sets the prefix used by PETSc when searching the options database.
Definition petsc.cpp:615
 
This class is a base class for matrices that can be used in petsc::KrylovSolver.
Definition petsc.h:237
 
std::array< std::int64_t, 2 > size() const
Return number of rows and columns (num_rows, num_cols). PETSc returns -1 if size has not been set.
Definition petsc.cpp:520
 
Operator & operator=(const Operator &A)=delete
Assignment operator (deleted)
 
Vec create_vector(std::size_t dim) const
Initialize vector to be compatible with the matrix-vector product y = Ax. In the parallel case,...
Definition petsc.cpp:530
 
Mat mat() const
Return PETSc Mat pointer.
Definition petsc.cpp:558
 
virtual ~Operator()
Destructor.
Definition petsc.cpp:506
 
A simple wrapper for a PETSc vector pointer (Vec). Its main purpose is to assist with memory/lifetime...
Definition petsc.h:164
 
void set_from_options()
Call PETSc function VecSetFromOptions on the underlying Vec object.
Definition petsc.cpp:485
 
std::int64_t size() const
Return global size of the vector.
Definition petsc.cpp:432
 
std::string get_options_prefix() const
Returns the prefix used by PETSc when searching the options database.
Definition petsc.cpp:476
 
Vector copy() const
Create a copy of the vector.
Definition petsc.cpp:422
 
virtual ~Vector()
Destructor.
Definition petsc.cpp:410
 
MPI_Comm comm() const
Return MPI communicator.
Definition petsc.cpp:460
 
std::array< std::int64_t, 2 > local_range() const
Return ownership range for calling rank.
Definition petsc.cpp:450
 
std::int32_t local_size() const
Return local size of vector (belonging to the call rank)
Definition petsc.cpp:441
 
Vec vec() const
Return pointer to PETSc Vec object.
Definition petsc.cpp:492
 
void set_options_prefix(std::string options_prefix)
Sets the prefix used by PETSc when searching the options database.
Definition petsc.cpp:469
 
Miscellaneous classes, functions and types.
Definition dolfinx_common.h:8
 
void clear()
Clear PETSc global options database.
Definition petsc.cpp:387
 
void set(std::string option)
Set PETSc option that takes no value.
Definition petsc.cpp:371
 
Vec create_vector_wrap(const common::IndexMap &map, int bs, std::span< const PetscScalar > x)
Create a PETSc Vec that wraps the data in an array.
Definition petsc.cpp:99
 
MatNullSpace create_nullspace(MPI_Comm comm, std::span< const Vec > basis)
Create PETSc MatNullSpace. Caller is responsible for destruction returned object.
Definition petsc.cpp:359
 
void scatter_local_vectors(Vec x, const std::vector< std::span< const PetscScalar > > &x_b, const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int > > &maps)
Scatter local vectors to Vec.
Definition petsc.cpp:189
 
void error(int error_code, std::string filename, std::string petsc_function)
Print error message for PETSc calls that return an error.
Definition petsc.cpp:31
 
std::vector< IS > create_index_sets(const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int > > &maps)
Definition petsc.cpp:126
 
std::vector< Vec > create_vectors(MPI_Comm comm, const std::vector< std::span< const PetscScalar > > &x)
Create PETsc vectors from the local data. The data is copied into the PETSc vectors and is not shared...
Definition petsc.cpp:49
 
Mat create_matrix(MPI_Comm comm, const SparsityPattern &sp, const std::string &type=std::string())
Create a PETSc Mat. Caller is responsible for destroying the returned object.
Definition petsc.cpp:232
 
std::vector< std::vector< PetscScalar > > get_local_vectors(const Vec x, const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int > > &maps)
Copy blocks from Vec into local vectors.
Definition petsc.cpp:146
 
Vec create_vector(const common::IndexMap &map, int bs)
Create a ghosted PETSc Vec.
Definition petsc.cpp:65
 
Linear algebra interface.
Definition sparsitybuild.h:15
 
Norm
Norm types.
Definition utils.h:17