Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/dd/dd6/la_2petsc_8h_source.html
DOLFINx 0.7.3
DOLFINx C++ interface
Loading...
Searching...
No Matches
petsc.h
1// Copyright (C) 2004-2018 Johan Hoffman, Johan Jansson, Anders Logg and
2// Garth N. Wells
3//
4// This file is part of DOLFINx (https://www.fenicsproject.org)
5//
6// SPDX-License-Identifier: LGPL-3.0-or-later
7
8#pragma once
9
10#include "Vector.h"
11#include "utils.h"
12#include <boost/lexical_cast.hpp>
13#include <functional>
14#include <petscksp.h>
15#include <petscmat.h>
16#include <petscoptions.h>
17#include <petscvec.h>
18#include <span>
19#include <string>
20#include <vector>
21
22namespace dolfinx::common
23{
24class IndexMap;
25} // namespace dolfinx::common
26
27namespace dolfinx::la
28{
29class SparsityPattern;
30
32namespace petsc
33{
35void error(int error_code, std::string filename, std::string petsc_function);
36
44std::vector<Vec>
45create_vectors(MPI_Comm comm,
46 const std::vector<std::span<const PetscScalar>>& x);
47
53Vec create_vector(const common::IndexMap& map, int bs);
54
63Vec create_vector(MPI_Comm comm, std::array<std::int64_t, 2> range,
64 std::span<const std::int64_t> ghosts, int bs);
65
75Vec create_vector_wrap(const common::IndexMap& map, int bs,
76 std::span<const PetscScalar> x);
77
81template <class V>
83{
84 assert(x.index_map());
85 return create_vector_wrap(*x.index_map(), x.bs(), x.array());
86}
87
99std::vector<IS> create_index_sets(
100 const std::vector<
101 std::pair<std::reference_wrapper<const common::IndexMap>, int>>& maps);
102
104std::vector<std::vector<PetscScalar>> get_local_vectors(
105 const Vec x,
106 const std::vector<
107 std::pair<std::reference_wrapper<const common::IndexMap>, int>>& maps);
108
111 Vec x, const std::vector<std::span<const PetscScalar>>& x_b,
112 const std::vector<
113 std::pair<std::reference_wrapper<const common::IndexMap>, int>>& maps);
114
117Mat create_matrix(MPI_Comm comm, const SparsityPattern& sp,
118 const std::string& type = std::string());
119
125MatNullSpace create_nullspace(MPI_Comm comm, std::span<const Vec> basis);
126
132namespace options
133{
135void set(std::string option);
136
138template <typename T>
139void set(std::string option, const T value)
140{
141 if (option[0] != '-')
142 option = '-' + option;
143
144 PetscErrorCode ierr;
145 ierr = PetscOptionsSetValue(nullptr, option.c_str(),
146 boost::lexical_cast<std::string>(value).c_str());
147 if (ierr != 0)
148 petsc::error(ierr, __FILE__, "PetscOptionsSetValue");
149}
150
152void clear(std::string option);
153
155void clear();
156} // namespace options
157
164{
165public:
170 Vector(const common::IndexMap& map, int bs);
171
172 // Delete copy constructor to avoid accidental copying of 'heavy' data
173 Vector(const Vector& x) = delete;
174
176 Vector(Vector&& x);
177
189 Vector(Vec x, bool inc_ref_count);
190
192 virtual ~Vector();
193
194 // Assignment operator (disabled)
195 Vector& operator=(const Vector& x) = delete;
196
198 Vector& operator=(Vector&& x);
199
202 Vector copy() const;
203
205 std::int64_t size() const;
206
208 std::int32_t local_size() const;
209
211 std::array<std::int64_t, 2> local_range() const;
212
214 MPI_Comm comm() const;
215
217 void set_options_prefix(std::string options_prefix);
218
221 std::string get_options_prefix() const;
222
224 void set_from_options();
225
227 Vec vec() const;
228
229private:
230 // PETSc Vec pointer
231 Vec _x;
232};
233
237{
238public:
240 Operator(Mat A, bool inc_ref_count);
241
242 // Copy constructor (deleted)
243 Operator(const Operator& A) = delete;
244
246 Operator(Operator&& A);
247
249 virtual ~Operator();
250
252 Operator& operator=(const Operator& A) = delete;
253
256
259 std::array<std::int64_t, 2> size() const;
260
266 Vec create_vector(std::size_t dim) const;
267
269 Mat mat() const;
270
271protected:
272 // PETSc Mat pointer
273 Mat _matA;
274};
275
281class Matrix : public Operator
282{
283public:
288 static auto set_fn(Mat A, InsertMode mode)
289 {
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
294 {
295 PetscErrorCode ierr;
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,
304 vals.data(), mode);
305#else
306 ierr = MatSetValuesLocal(A, rows.size(), rows.data(), cols.size(),
307 cols.data(), vals.data(), mode);
308#endif
309
310#ifndef NDEBUG
311 if (ierr != 0)
312 petsc::error(ierr, __FILE__, "MatSetValuesLocal");
313#endif
314 return ierr;
315 };
316 }
317
323 static auto set_block_fn(Mat A, InsertMode mode)
324 {
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
329 {
330 PetscErrorCode ierr;
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,
339 vals.data(), mode);
340#else
341 ierr = MatSetValuesBlockedLocal(A, rows.size(), rows.data(), cols.size(),
342 cols.data(), vals.data(), mode);
343#endif
344
345#ifndef NDEBUG
346 if (ierr != 0)
347 petsc::error(ierr, __FILE__, "MatSetValuesBlockedLocal");
348#endif
349 return ierr;
350 };
351 }
352
361 static auto set_block_expand_fn(Mat A, int bs0, int bs1, InsertMode mode)
362 {
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
368 {
369 PetscErrorCode ierr;
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;
375
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;
379
380 ierr = MatSetValuesLocal(A, cache0.size(), cache0.data(), cache1.size(),
381 cache1.data(), vals.data(), mode);
382#ifndef NDEBUG
383 if (ierr != 0)
384 petsc::error(ierr, __FILE__, "MatSetValuesLocal");
385#endif
386 return ierr;
387 };
388 }
389
391 Matrix(MPI_Comm comm, const SparsityPattern& sp,
392 const std::string& type = std::string());
393
398 Matrix(Mat A, bool inc_ref_count);
399
400 // Copy constructor (deleted)
401 Matrix(const Matrix& A) = delete;
402
404 Matrix(Matrix&& A) = default;
405
407 ~Matrix() = default;
408
410 Matrix& operator=(const Matrix& A) = delete;
411
413 Matrix& operator=(Matrix&& A) = default;
414
418 enum class AssemblyType : std::int32_t
419 {
420 FINAL,
421 FLUSH
422 };
423
429 void apply(AssemblyType type);
430
432 double norm(Norm norm_type) const;
433
434 //--- Special PETSc Functions ---
435
438 void set_options_prefix(std::string options_prefix);
439
442 std::string get_options_prefix() const;
443
445 void set_from_options();
446};
447
451{
452public:
455 explicit KrylovSolver(MPI_Comm comm);
456
460 KrylovSolver(KSP ksp, bool inc_ref_count);
461
462 // Copy constructor (deleted)
463 KrylovSolver(const KrylovSolver& solver) = delete;
464
466 KrylovSolver(KrylovSolver&& solver);
467
470
471 // Assignment operator (deleted)
472 KrylovSolver& operator=(const KrylovSolver&) = delete;
473
475 KrylovSolver& operator=(KrylovSolver&& solver);
476
478 void set_operator(const Mat A);
479
481 void set_operators(const Mat A, const Mat P);
482
485 int solve(Vec x, const Vec b, bool transpose = false) const;
486
489 void set_options_prefix(std::string options_prefix);
490
493 std::string get_options_prefix() const;
494
496 void set_from_options() const;
497
499 KSP ksp() const;
500
501private:
502 // PETSc solver pointer
503 KSP _ksp;
504};
505} // namespace petsc
506} // namespace dolfinx::la
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