Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/d0/d75/fem_2petsc_8h_source.html
DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
petsc.h
1// Copyright (C) 2018-2021 Garth N. Wells
2//
3// This file is part of DOLFINx (https://www.fenicsproject.org)
4//
5// SPDX-License-Identifier: LGPL-3.0-or-later
6
7#pragma once
8
9#ifdef HAS_PETSC
10
11#include "Form.h"
12#include "assembler.h"
13#include "utils.h"
14#include <concepts>
15#include <dolfinx/la/petsc.h>
16#include <map>
17#include <memory>
18#include <petscmat.h>
19#include <petscvec.h>
20#include <span>
21#include <utility>
22#include <vector>
23
24namespace dolfinx::common
25{
26class IndexMap;
27}
28
29namespace dolfinx::fem
30{
31template <dolfinx::scalar T, std::floating_point U>
32class DirichletBC;
33
35namespace petsc
36{
43template <std::floating_point T>
45 std::string type = std::string())
46{
48 pattern.finalize();
49 return la::petsc::create_matrix(a.mesh()->comm(), pattern, type);
50}
51
60template <std::floating_point T>
62 const std::vector<std::vector<const Form<PetscScalar, T>*>>& a,
63 std::string type = std::string())
64{
65 // Extract and check row/column ranges
66 std::array<std::vector<std::shared_ptr<const FunctionSpace<T>>>, 2> V
68 std::array<std::vector<int>, 2> bs_dofs;
69 for (std::size_t i = 0; i < 2; ++i)
70 {
71 for (auto& _V : V[i])
72 bs_dofs[i].push_back(_V->dofmap()->bs());
73 }
74
75 // Build sparsity pattern for each block
76 std::shared_ptr<const mesh::Mesh<T>> mesh;
77 std::vector<std::vector<std::unique_ptr<la::SparsityPattern>>> patterns(
78 V[0].size());
79 for (std::size_t row = 0; row < V[0].size(); ++row)
80 {
81 for (std::size_t col = 0; col < V[1].size(); ++col)
82 {
83 if (const Form<PetscScalar, T>* form = a[row][col]; form)
84 {
85 patterns[row].push_back(std::make_unique<la::SparsityPattern>(
87 if (!mesh)
88 mesh = form->mesh();
89 }
90 else
91 patterns[row].push_back(nullptr);
92 }
93 }
94
95 if (!mesh)
96 throw std::runtime_error("Could not find a Mesh.");
97
98 // Compute offsets for the fields
99 std::array<std::vector<std::pair<
100 std::reference_wrapper<const common::IndexMap>, int>>,
101 2>
102 maps;
103 for (std::size_t d = 0; d < 2; ++d)
104 {
105 for (auto& space : V[d])
106 {
107 maps[d].emplace_back(*space->dofmap()->index_map,
108 space->dofmap()->index_map_bs());
109 }
110 }
111
112 // Create merged sparsity pattern
113 std::vector<std::vector<const la::SparsityPattern*>> p(V[0].size());
114 for (std::size_t row = 0; row < V[0].size(); ++row)
115 for (std::size_t col = 0; col < V[1].size(); ++col)
116 p[row].push_back(patterns[row][col].get());
117
118 la::SparsityPattern pattern(mesh->comm(), p, maps, bs_dofs);
119 pattern.finalize();
120
121 // FIXME: Add option to pass customised local-to-global map to PETSc
122 // Mat constructor
123
124 // Initialise matrix
125 Mat A = la::petsc::create_matrix(mesh->comm(), pattern, type);
126
127 // Create row and column local-to-global maps (field0, field1, field2,
128 // etc), i.e. ghosts of field0 appear before owned indices of field1
129 std::array<std::vector<PetscInt>, 2> _maps;
130 for (int d = 0; d < 2; ++d)
131 {
132 // TODO: Index map concatenation has already been computed inside
133 // the SparsityPattern constructor, but we also need it here to
134 // build the PETSc local-to-global map. Compute outside and pass
135 // into SparsityPattern constructor.
136 // TODO: avoid concatenating the same maps twice in case that V[0]
137 // == V[1].
138
139 const std::vector<
140 std::pair<std::reference_wrapper<const common::IndexMap>, int>>& map
141 = maps[d];
142 std::vector<PetscInt>& _map = _maps[d];
143
144 // Concatenate the block index map in the row and column directions
145 const auto [rank_offset, local_offset, ghosts, _]
147 for (std::size_t f = 0; f < map.size(); ++f)
148 {
149 auto offset = local_offset[f];
150 const common::IndexMap& imap = map[f].first.get();
151 int bs = map[f].second;
152 for (std::int32_t i = 0; i < bs * imap.size_local(); ++i)
153 _map.push_back(i + rank_offset + offset);
154 for (std::int32_t i = 0; i < bs * imap.num_ghosts(); ++i)
155 _map.push_back(ghosts[f][i]);
156 }
157 }
158
159 // Create PETSc local-to-global map/index sets and attach to matrix
160 ISLocalToGlobalMapping petsc_local_to_global0;
161 ISLocalToGlobalMappingCreate(MPI_COMM_SELF, 1, _maps[0].size(),
162 _maps[0].data(), PETSC_COPY_VALUES,
163 &petsc_local_to_global0);
164 if (V[0] == V[1])
165 {
166 MatSetLocalToGlobalMapping(A, petsc_local_to_global0,
167 petsc_local_to_global0);
168 ISLocalToGlobalMappingDestroy(&petsc_local_to_global0);
169 }
170 else
171 {
172
173 ISLocalToGlobalMapping petsc_local_to_global1;
174 ISLocalToGlobalMappingCreate(MPI_COMM_SELF, 1, _maps[1].size(),
175 _maps[1].data(), PETSC_COPY_VALUES,
176 &petsc_local_to_global1);
177 MatSetLocalToGlobalMapping(A, petsc_local_to_global0,
178 petsc_local_to_global1);
179 ISLocalToGlobalMappingDestroy(&petsc_local_to_global0);
180 ISLocalToGlobalMappingDestroy(&petsc_local_to_global1);
181 }
182
183 return A;
184}
185
189template <std::floating_point T>
191 const std::vector<std::vector<const Form<PetscScalar, T>*>>& a,
192 const std::vector<std::vector<std::string>>& types)
193{
194 // Extract and check row/column ranges
196
197 std::vector<std::vector<std::string>> _types(
198 a.size(), std::vector<std::string>(a.front().size()));
199 if (!types.empty())
200 _types = types;
201
202 // Loop over each form and create matrix
203 int rows = a.size();
204 int cols = a.front().size();
205 std::vector<Mat> mats(rows * cols, nullptr);
206 std::shared_ptr<const mesh::Mesh<T>> mesh;
207 for (int i = 0; i < rows; ++i)
208 {
209 for (int j = 0; j < cols; ++j)
210 {
211 if (const Form<PetscScalar, T>* form = a[i][j]; form)
212 {
213 mats[i * cols + j] = create_matrix(*form, _types[i][j]);
214 mesh = form->mesh();
215 }
216 }
217 }
218
219 if (!mesh)
220 throw std::runtime_error("Could not find a Mesh.");
221
222 // Initialise block (MatNest) matrix
223 Mat A;
224 MatCreate(mesh->comm(), &A);
225 MatSetType(A, MATNEST);
226 MatNestSetSubMats(A, rows, nullptr, cols, nullptr, mats.data());
227 MatSetUp(A);
228
229 // De-reference Mat objects
230 for (std::size_t i = 0; i < mats.size(); ++i)
231 {
232 if (mats[i])
233 MatDestroy(&mats[i]);
234 }
235
236 return A;
237}
238
243 const std::vector<
244 std::pair<std::reference_wrapper<const common::IndexMap>, int>>& maps);
245
248 const std::vector<
249 std::pair<std::reference_wrapper<const common::IndexMap>, int>>& maps);
250
251// -- Vectors ----------------------------------------------------------------
252
265template <std::floating_point T>
267 Vec b, const Form<PetscScalar, T>& L,
268 std::span<const PetscScalar> constants,
269 const std::map<std::pair<IntegralType, int>,
270 std::pair<std::span<const PetscScalar>, int>>& coeffs)
271{
272 Vec b_local;
273 VecGhostGetLocalForm(b, &b_local);
274 PetscInt n = 0;
275 VecGetSize(b_local, &n);
276 PetscScalar* array = nullptr;
277 VecGetArray(b_local, &array);
278 std::span<PetscScalar> _b(array, n);
279 fem::assemble_vector(_b, L, constants, coeffs);
280 VecRestoreArray(b_local, &array);
281 VecGhostRestoreLocalForm(b, &b_local);
282}
283
293template <std::floating_point T>
295{
296 Vec b_local;
297 VecGhostGetLocalForm(b, &b_local);
298 PetscInt n = 0;
299 VecGetSize(b_local, &n);
300 PetscScalar* array = nullptr;
301 VecGetArray(b_local, &array);
302 std::span<PetscScalar> _b(array, n);
304 VecRestoreArray(b_local, &array);
305 VecGhostRestoreLocalForm(b, &b_local);
306}
307
308// FIXME: clarify how x0 is used
309// FIXME: if bcs entries are set
310
311// FIXME: need to pass an array of Vec for x0?
312// FIXME: clarify zeroing of vector
313
329template <std::floating_point T>
331 Vec b, const std::vector<std::shared_ptr<const Form<PetscScalar, T>>>& a,
332 const std::vector<std::span<const PetscScalar>>& constants,
333 const std::vector<std::map<std::pair<IntegralType, int>,
334 std::pair<std::span<const PetscScalar>, int>>>&
335 coeffs,
336 const std::vector<
337 std::vector<std::shared_ptr<const DirichletBC<PetscScalar, T>>>>& bcs1,
338 const std::vector<Vec>& x0, PetscScalar scale)
339{
340 Vec b_local;
341 VecGhostGetLocalForm(b, &b_local);
342 PetscInt n = 0;
343 VecGetSize(b_local, &n);
344 PetscScalar* array = nullptr;
345 VecGetArray(b_local, &array);
346 std::span<PetscScalar> _b(array, n);
347
348 if (x0.empty())
349 fem::apply_lifting(_b, a, constants, coeffs, bcs1, {}, scale);
350 else
351 {
352 std::vector<std::span<const PetscScalar>> x0_ref;
353 std::vector<Vec> x0_local(a.size());
354 std::vector<const PetscScalar*> x0_array(a.size());
355 for (std::size_t i = 0; i < a.size(); ++i)
356 {
357 assert(x0[i]);
358 VecGhostGetLocalForm(x0[i], &x0_local[i]);
359 PetscInt n = 0;
360 VecGetSize(x0_local[i], &n);
361 VecGetArrayRead(x0_local[i], &x0_array[i]);
362 x0_ref.emplace_back(x0_array[i], n);
363 }
364
365 std::vector x0_tmp(x0_ref.begin(), x0_ref.end());
366 fem::apply_lifting(_b, a, constants, coeffs, bcs1, x0_tmp, scale);
367
368 for (std::size_t i = 0; i < x0_local.size(); ++i)
369 {
370 VecRestoreArrayRead(x0_local[i], &x0_array[i]);
371 VecGhostRestoreLocalForm(x0[i], &x0_local[i]);
372 }
373 }
374
375 VecRestoreArray(b_local, &array);
376 VecGhostRestoreLocalForm(b, &b_local);
377}
378
379// FIXME: clarify how x0 is used
380// FIXME: if bcs entries are set
381
382// FIXME: need to pass an array of Vec for x0?
383// FIXME: clarify zeroing of vector
384
397template <std::floating_point T>
399 Vec b,
400 const std::vector<std::shared_ptr<const Form<PetscScalar, double>>>& a,
401 const std::vector<
402 std::vector<std::shared_ptr<const DirichletBC<PetscScalar, double>>>>&
403 bcs1,
404 const std::vector<Vec>& x0, PetscScalar scale)
405{
406 Vec b_local;
407 VecGhostGetLocalForm(b, &b_local);
408 PetscInt n = 0;
409 VecGetSize(b_local, &n);
410 PetscScalar* array = nullptr;
411 VecGetArray(b_local, &array);
412 std::span<PetscScalar> _b(array, n);
413
414 if (x0.empty())
415 fem::apply_lifting<PetscScalar>(_b, a, bcs1, {}, scale);
416 else
417 {
418 std::vector<std::span<const PetscScalar>> x0_ref;
419 std::vector<Vec> x0_local(a.size());
420 std::vector<const PetscScalar*> x0_array(a.size());
421 for (std::size_t i = 0; i < a.size(); ++i)
422 {
423 assert(x0[i]);
424 VecGhostGetLocalForm(x0[i], &x0_local[i]);
425 PetscInt n = 0;
426 VecGetSize(x0_local[i], &n);
427 VecGetArrayRead(x0_local[i], &x0_array[i]);
428 x0_ref.emplace_back(x0_array[i], n);
429 }
430
431 std::vector x0_tmp(x0_ref.begin(), x0_ref.end());
432 fem::apply_lifting<PetscScalar>(_b, a, bcs1, x0_tmp, scale);
433
434 for (std::size_t i = 0; i < x0_local.size(); ++i)
435 {
436 VecRestoreArrayRead(x0_local[i], &x0_array[i]);
437 VecGhostRestoreLocalForm(x0[i], &x0_local[i]);
438 }
439 }
440
441 VecRestoreArray(b_local, &array);
442 VecGhostRestoreLocalForm(b, &b_local);
443}
444
445// -- Setting bcs ------------------------------------------------------------
446
447// FIXME: Move these function elsewhere?
448
449// FIXME: clarify x0
450// FIXME: clarify what happens with ghosts
451
455template <std::floating_point T>
457 Vec b,
458 const std::vector<std::shared_ptr<const DirichletBC<PetscScalar, T>>>& bcs,
459 const Vec x0, PetscScalar scale = 1)
460{
461 PetscInt n = 0;
462 VecGetLocalSize(b, &n);
463 PetscScalar* array = nullptr;
464 VecGetArray(b, &array);
465 std::span<PetscScalar> _b(array, n);
466 if (x0)
467 {
468 Vec x0_local;
469 VecGhostGetLocalForm(x0, &x0_local);
470 PetscInt n = 0;
471 VecGetSize(x0_local, &n);
472 const PetscScalar* array = nullptr;
473 VecGetArrayRead(x0_local, &array);
474 std::span<const PetscScalar> _x0(array, n);
475 fem::set_bc(_b, bcs, _x0, scale);
476 VecRestoreArrayRead(x0_local, &array);
477 VecGhostRestoreLocalForm(x0, &x0_local);
478 }
479 else
480 fem::set_bc(_b, bcs, scale);
481
482 VecRestoreArray(b, &array);
483}
484
485} // namespace petsc
486} // namespace dolfinx::fem
487
488#endif
Definition IndexMap.h:94
std::int32_t num_ghosts() const noexcept
Number of ghost indices on this process.
Definition IndexMap.cpp:866
std::int32_t size_local() const noexcept
Number of indices owned by this process.
Definition IndexMap.cpp:868
Definition DirichletBC.h:259
A representation of finite element variational forms.
Definition Form.h:120
Definition SparsityPattern.h:26
void finalize()
Finalize sparsity pattern and communicate off-process entries.
Definition SparsityPattern.cpp:226
Functions supporting finite element method operations.
Miscellaneous classes, functions and types.
Definition dolfinx_common.h:8
std::tuple< std::int64_t, std::vector< std::int32_t >, std::vector< std::vector< std::int64_t > >, std::vector< std::vector< int > > > stack_index_maps(const std::vector< std::pair< std::reference_wrapper< const IndexMap >, int > > &maps)
Compute layout data and ghost indices for a stacked (concatenated) index map, i.e....
Definition IndexMap.cpp:583
void set_bc(Vec b, const std::vector< std::shared_ptr< const DirichletBC< PetscScalar, T > > > &bcs, const Vec x0, PetscScalar scale=1)
Definition petsc.h:456
void apply_lifting(Vec b, const std::vector< std::shared_ptr< const Form< PetscScalar, T > > > &a, const std::vector< std::span< const PetscScalar > > &constants, const std::vector< std::map< std::pair< IntegralType, int >, std::pair< std::span< const PetscScalar >, int > > > &coeffs, const std::vector< std::vector< std::shared_ptr< const DirichletBC< PetscScalar, T > > > > &bcs1, const std::vector< Vec > &x0, PetscScalar scale)
Modify RHS vector to account for Dirichlet boundary conditions.
Definition petsc.h:330
Mat create_matrix(const Form< PetscScalar, T > &a, std::string type=std::string())
Definition petsc.h:44
void assemble_vector(Vec b, const Form< PetscScalar, T > &L, std::span< const PetscScalar > constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const PetscScalar >, int > > &coeffs)
Assemble linear form into an already allocated PETSc vector.
Definition petsc.h:266
Mat create_matrix_block(const std::vector< std::vector< const Form< PetscScalar, T > * > > &a, std::string type=std::string())
Definition petsc.h:61
Mat create_matrix_nest(const std::vector< std::vector< const Form< PetscScalar, T > * > > &a, const std::vector< std::vector< std::string > > &types)
Create nested (MatNest) matrix.
Definition petsc.h:190
Vec create_vector_block(const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int > > &maps)
Definition petsc.cpp:17
Vec create_vector_nest(const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int > > &maps)
Create nested (VecNest) vector. Vector is not zeroed.
Definition petsc.cpp:62
Finite element method functionality.
Definition assemble_matrix_impl.h:26
void assemble_vector(std::span< T > b, const Form< T, U > &L, std::span< const T > constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > &coefficients)
Assemble linear form into a vector.
Definition assembler.h:108
void apply_lifting(std::span< T > b, const std::vector< std::shared_ptr< const Form< T, U > > > &a, const std::vector< std::span< const T > > &constants, const std::vector< std::map< std::pair< IntegralType, int >, std::pair< std::span< const T >, int > > > &coeffs, const std::vector< std::vector< std::shared_ptr< const DirichletBC< T, U > > > > &bcs1, const std::vector< std::span< const T > > &x0, T scale)
Definition assembler.h:149
std::vector< std::vector< std::array< std::shared_ptr< const FunctionSpace< U > >, 2 > > > extract_function_spaces(const std::vector< std::vector< const Form< T, U > * > > &a)
Extract test (0) and trial (1) function spaces pairs for each bilinear form for a rectangular array o...
Definition utils.h:125
void set_bc(std::span< T > b, const std::vector< std::shared_ptr< const DirichletBC< T, U > > > &bcs, std::span< const T > x0, T scale=1)
Definition assembler.h:437
la::SparsityPattern create_sparsity_pattern(const Form< T, U > &a)
Create a sparsity pattern for a given form.
Definition utils.h:150
std::array< std::vector< std::shared_ptr< const FunctionSpace< T > > >, 2 > common_function_spaces(const std::vector< std::vector< std::array< std::shared_ptr< const FunctionSpace< T > >, 2 > > > &V)
Definition FunctionSpace.h:385
Mat create_matrix(MPI_Comm comm, const SparsityPattern &sp, std::string type=std::string())
Definition petsc.cpp:234