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.6.0
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
31namespace petsc
32{
34void error(int error_code, std::string filename, std::string petsc_function);
35
43std::vector<Vec>
44create_vectors(MPI_Comm comm,
45 const std::vector<std::span<const PetscScalar>>& x);
46
52Vec create_vector(const common::IndexMap& map, int bs);
53
62Vec create_vector(MPI_Comm comm, std::array<std::int64_t, 2> range,
63 std::span<const std::int64_t> ghosts, int bs);
64
74Vec create_vector_wrap(const common::IndexMap& map, int bs,
75 std::span<const PetscScalar> x);
76
80template <typename Allocator>
81Vec create_vector_wrap(const la::Vector<PetscScalar, Allocator>& x)
82{
83 assert(x.map());
84 return create_vector_wrap(*x.map(), x.bs(), x.array());
85}
86
98std::vector<IS> create_index_sets(
99 const std::vector<
100 std::pair<std::reference_wrapper<const common::IndexMap>, int>>& maps);
101
103std::vector<std::vector<PetscScalar>> get_local_vectors(
104 const Vec x,
105 const std::vector<
106 std::pair<std::reference_wrapper<const common::IndexMap>, int>>& maps);
107
109void scatter_local_vectors(
110 Vec x, const std::vector<std::span<const PetscScalar>>& x_b,
111 const std::vector<
112 std::pair<std::reference_wrapper<const common::IndexMap>, int>>& maps);
113
116Mat create_matrix(MPI_Comm comm, const SparsityPattern& sp,
117 const std::string& type = std::string());
118
124MatNullSpace create_nullspace(MPI_Comm comm, std::span<const Vec> basis);
125
131namespace options
132{
134void set(std::string option);
135
137template <typename T>
138void set(std::string option, const T value)
139{
140 if (option[0] != '-')
141 option = '-' + option;
142
143 PetscErrorCode ierr;
144 ierr = PetscOptionsSetValue(nullptr, option.c_str(),
145 boost::lexical_cast<std::string>(value).c_str());
146 if (ierr != 0)
147 petsc::error(ierr, __FILE__, "PetscOptionsSetValue");
148}
149
151void clear(std::string option);
152
154void clear();
155} // namespace options
156
163{
164public:
169 Vector(const common::IndexMap& map, int bs);
170
171 // Delete copy constructor to avoid accidental copying of 'heavy' data
172 Vector(const Vector& x) = delete;
173
175 Vector(Vector&& x);
176
188 Vector(Vec x, bool inc_ref_count);
189
191 virtual ~Vector();
192
193 // Assignment operator (disabled)
194 Vector& operator=(const Vector& x) = delete;
195
197 Vector& operator=(Vector&& x);
198
201 Vector copy() const;
202
204 std::int64_t size() const;
205
207 std::int32_t local_size() const;
208
210 std::array<std::int64_t, 2> local_range() const;
211
213 MPI_Comm comm() const;
214
216 void set_options_prefix(std::string options_prefix);
217
220 std::string get_options_prefix() const;
221
223 void set_from_options();
224
226 Vec vec() const;
227
228private:
229 // PETSc Vec pointer
230 Vec _x;
231};
232
236{
237public:
239 Operator(Mat A, bool inc_ref_count);
240
241 // Copy constructor (deleted)
242 Operator(const Operator& A) = delete;
243
245 Operator(Operator&& A);
246
248 virtual ~Operator();
249
251 Operator& operator=(const Operator& A) = delete;
252
255
258 std::array<std::int64_t, 2> size() const;
259
265 Vec create_vector(std::size_t dim) const;
266
268 Mat mat() const;
269
270protected:
271 // PETSc Mat pointer
272 Mat _matA;
273};
274
280class Matrix : public Operator
281{
282public:
287 static auto set_fn(Mat A, InsertMode mode)
288 {
289 return [A, mode, cache = std::vector<PetscInt>()](
290 const std::span<const std::int32_t>& rows,
291 const std::span<const std::int32_t>& cols,
292 const std::span<const PetscScalar>& vals) mutable -> int
293 {
294 PetscErrorCode ierr;
295#ifdef PETSC_USE_64BIT_INDICES
296 cache.resize(rows.size() + cols.size());
297 std::copy(rows.begin(), rows.end(), cache.begin());
298 std::copy(cols.begin(), cols.end(),
299 std::next(cache.begin(), rows.size()));
300 const PetscInt* _rows = cache.data();
301 const PetscInt* _cols = cache.data() + rows.size();
302 ierr = MatSetValuesLocal(A, rows.size(), _rows, cols.size(), _cols,
303 vals.data(), mode);
304#else
305 ierr = MatSetValuesLocal(A, rows.size(), rows.data(), cols.size(),
306 cols.data(), vals.data(), mode);
307#endif
308
309#ifndef NDEBUG
310 if (ierr != 0)
311 petsc::error(ierr, __FILE__, "MatSetValuesLocal");
312#endif
313 return ierr;
314 };
315 }
316
322 static auto set_block_fn(Mat A, InsertMode mode)
323 {
324 return [A, mode, cache = std::vector<PetscInt>()](
325 const std::span<const std::int32_t>& rows,
326 const std::span<const std::int32_t>& cols,
327 const std::span<const PetscScalar>& vals) mutable -> int
328 {
329 PetscErrorCode ierr;
330#ifdef PETSC_USE_64BIT_INDICES
331 cache.resize(rows.size() + cols.size());
332 std::copy(rows.begin(), rows.end(), cache.begin());
333 std::copy(cols.begin(), cols.end(),
334 std::next(cache.begin(), rows.size()));
335 const PetscInt* _rows = cache.data();
336 const PetscInt* _cols = cache.data() + rows.size();
337 ierr = MatSetValuesBlockedLocal(A, rows.size(), _rows, cols.size(), _cols,
338 vals.data(), mode);
339#else
340 ierr = MatSetValuesBlockedLocal(A, rows.size(), rows.data(), cols.size(),
341 cols.data(), vals.data(), mode);
342#endif
343
344#ifndef NDEBUG
345 if (ierr != 0)
346 petsc::error(ierr, __FILE__, "MatSetValuesBlockedLocal");
347#endif
348 return ierr;
349 };
350 }
351
360 static auto set_block_expand_fn(Mat A, int bs0, int bs1, InsertMode mode)
361 {
362 return [A, bs0, bs1, mode, cache0 = std::vector<PetscInt>(),
363 cache1 = std::vector<PetscInt>()](
364 const std::span<const std::int32_t>& rows,
365 const std::span<const std::int32_t>& cols,
366 const std::span<const PetscScalar>& vals) mutable -> int
367 {
368 PetscErrorCode ierr;
369 cache0.resize(bs0 * rows.size());
370 cache1.resize(bs1 * cols.size());
371 for (std::size_t i = 0; i < rows.size(); ++i)
372 for (int k = 0; k < bs0; ++k)
373 cache0[bs0 * i + k] = bs0 * rows[i] + k;
374
375 for (std::size_t i = 0; i < cols.size(); ++i)
376 for (int k = 0; k < bs1; ++k)
377 cache1[bs1 * i + k] = bs1 * cols[i] + k;
378
379 ierr = MatSetValuesLocal(A, cache0.size(), cache0.data(), cache1.size(),
380 cache1.data(), vals.data(), mode);
381#ifndef NDEBUG
382 if (ierr != 0)
383 petsc::error(ierr, __FILE__, "MatSetValuesLocal");
384#endif
385 return ierr;
386 };
387 }
388
390 Matrix(MPI_Comm comm, const SparsityPattern& sp,
391 const std::string& type = std::string());
392
397 Matrix(Mat A, bool inc_ref_count);
398
399 // Copy constructor (deleted)
400 Matrix(const Matrix& A) = delete;
401
403 Matrix(Matrix&& A) = default;
404
406 ~Matrix() = default;
407
409 Matrix& operator=(const Matrix& A) = delete;
410
412 Matrix& operator=(Matrix&& A) = default;
413
417 enum class AssemblyType : std::int32_t
418 {
419 FINAL,
420 FLUSH
421 };
422
428 void apply(AssemblyType type);
429
431 double norm(Norm norm_type) const;
432
433 //--- Special PETSc Functions ---
434
437 void set_options_prefix(std::string options_prefix);
438
441 std::string get_options_prefix() const;
442
444 void set_from_options();
445};
446
450{
451public:
454 explicit KrylovSolver(MPI_Comm comm);
455
459 KrylovSolver(KSP ksp, bool inc_ref_count);
460
461 // Copy constructor (deleted)
462 KrylovSolver(const KrylovSolver& solver) = delete;
463
465 KrylovSolver(KrylovSolver&& solver);
466
469
470 // Assignment operator (deleted)
471 KrylovSolver& operator=(const KrylovSolver&) = delete;
472
474 KrylovSolver& operator=(KrylovSolver&& solver);
475
477 void set_operator(const Mat A);
478
480 void set_operators(const Mat A, const Mat P);
481
484 int solve(Vec x, const Vec b, bool transpose = false) const;
485
488 void set_options_prefix(std::string options_prefix);
489
492 std::string get_options_prefix() const;
493
495 void set_from_options() const;
496
498 KSP ksp() const;
499
500private:
501 // PETSc solver pointer
502 KSP _ksp;
503};
504} // namespace petsc
505} // namespace dolfinx::la
This class represents the distribution index arrays across processes. An index array is a contiguous ...
Definition: IndexMap.h:64
This class provides a sparsity pattern data structure that can be used to initialize sparse matrices....
Definition: SparsityPattern.h:34
This class implements Krylov methods for linear systems of the form Ax = b. It is a wrapper for the K...
Definition: petsc.h:450
void set_operator(const Mat A)
Set operator (Mat)
Definition: petsc.cpp:650
KSP ksp() const
Return PETSc KSP pointer.
Definition: petsc.cpp:778
std::string get_options_prefix() const
Returns the prefix used by PETSc when searching the PETSc options database.
Definition: petsc.cpp:760
void set_operators(const Mat A, const Mat P)
Set operator and preconditioner matrix (Mat)
Definition: petsc.cpp:652
void set_from_options() const
Set options from PETSc options database.
Definition: petsc.cpp:770
~KrylovSolver()
Destructor.
Definition: petsc.cpp:638
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:661
void set_options_prefix(std::string options_prefix)
Sets the prefix used by PETSc when searching the PETSc options database.
Definition: petsc.cpp:751
It is a simple wrapper for a PETSc matrix pointer (Mat). Its main purpose is to assist memory managem...
Definition: petsc.h:281
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:322
void set_from_options()
Call PETSc function MatSetFromOptions on the PETSc Mat object.
Definition: petsc.cpp:606
double norm(Norm norm_type) const
Return norm of matrix.
Definition: petsc.cpp:550
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:287
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:598
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:360
~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:418
void apply(AssemblyType type)
Finalize assembly of tensor. The following values are recognized for the mode parameter:
Definition: petsc.cpp:575
void set_options_prefix(std::string options_prefix)
Sets the prefix used by PETSc when searching the options database.
Definition: petsc.cpp:592
This class is a base class for matrices that can be used in petsc::KrylovSolver.
Definition: petsc.h:236
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:497
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:507
Mat mat() const
Return PETSc Mat pointer.
Definition: petsc.cpp:535
virtual ~Operator()
Destructor.
Definition: petsc.cpp:483
A simple wrapper for a PETSc vector pointer (Vec). Its main purpose is to assist with memory/lifetime...
Definition: petsc.h:163
void set_from_options()
Call PETSc function VecSetFromOptions on the underlying Vec object.
Definition: petsc.cpp:462
std::int64_t size() const
Return global size of the vector.
Definition: petsc.cpp:409
std::string get_options_prefix() const
Returns the prefix used by PETSc when searching the options database.
Definition: petsc.cpp:453
Vector copy() const
Create a copy of the vector.
Definition: petsc.cpp:399
virtual ~Vector()
Destructor.
Definition: petsc.cpp:387
MPI_Comm comm() const
Return MPI communicator.
Definition: petsc.cpp:437
std::array< std::int64_t, 2 > local_range() const
Return ownership range for calling rank.
Definition: petsc.cpp:427
std::int32_t local_size() const
Return local size of vector (belonging to the call rank)
Definition: petsc.cpp:418
Vec vec() const
Return pointer to PETSc Vec object.
Definition: petsc.cpp:469
void set_options_prefix(std::string options_prefix)
Sets the prefix used by PETSc when searching the options database.
Definition: petsc.cpp:446
Miscellaneous classes, functions and types.
Mat create_matrix(const Form< PetscScalar > &a, const std::string &type=std::string())
Create a matrix.
Definition: petsc.cpp:21
void clear()
Clear PETSc global options database.
Definition: petsc.cpp:363
void set(std::string option)
Set PETSc option that takes no value.
Definition: petsc.cpp:347
Linear algebra interface.
Definition: sparsitybuild.h:15
Norm
Norm types.
Definition: utils.h:17