DOLFINx 0.9.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>
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
61template <std::floating_point T>
63 const std::vector<std::vector<const Form<PetscScalar, T>*>>& a,
64 std::string type = std::string())
65{
66 // Extract and check row/column ranges
67 std::array<std::vector<std::shared_ptr<const FunctionSpace<T>>>, 2> V
69 std::array<std::vector<int>, 2> bs_dofs;
70 for (std::size_t i = 0; i < 2; ++i)
71 {
72 for (auto& _V : V[i])
73 bs_dofs[i].push_back(_V->dofmap()->bs());
74 }
75
76 // Build sparsity pattern for each block
77 std::shared_ptr<const mesh::Mesh<T>> mesh;
78 std::vector<std::vector<std::unique_ptr<la::SparsityPattern>>> patterns(
79 V[0].size());
80 for (std::size_t row = 0; row < V[0].size(); ++row)
81 {
82 for (std::size_t col = 0; col < V[1].size(); ++col)
83 {
84 if (const Form<PetscScalar, T>* form = a[row][col]; form)
85 {
86 patterns[row].push_back(std::make_unique<la::SparsityPattern>(
88 if (!mesh)
89 mesh = form->mesh();
90 }
91 else
92 patterns[row].push_back(nullptr);
93 }
94 }
95
96 if (!mesh)
97 throw std::runtime_error("Could not find a Mesh.");
98
99 // Compute offsets for the fields
100 std::array<std::vector<std::pair<
101 std::reference_wrapper<const common::IndexMap>, int>>,
102 2>
103 maps;
104 for (std::size_t d = 0; d < 2; ++d)
105 {
106 for (auto& space : V[d])
107 {
108 maps[d].emplace_back(*space->dofmap()->index_map,
109 space->dofmap()->index_map_bs());
110 }
111 }
112
113 // Create merged sparsity pattern
114 std::vector<std::vector<const la::SparsityPattern*>> p(V[0].size());
115 for (std::size_t row = 0; row < V[0].size(); ++row)
116 for (std::size_t col = 0; col < V[1].size(); ++col)
117 p[row].push_back(patterns[row][col].get());
118
119 la::SparsityPattern pattern(mesh->comm(), p, maps, bs_dofs);
120 pattern.finalize();
121
122 // FIXME: Add option to pass customised local-to-global map to PETSc
123 // Mat constructor
124
125 // Initialise matrix
126 Mat A = la::petsc::create_matrix(mesh->comm(), pattern, type);
127
128 // Create row and column local-to-global maps (field0, field1, field2,
129 // etc), i.e. ghosts of field0 appear before owned indices of field1
130 std::array<std::vector<PetscInt>, 2> _maps;
131 for (int d = 0; d < 2; ++d)
132 {
133 // TODO: Index map concatenation has already been computed inside
134 // the SparsityPattern constructor, but we also need it here to
135 // build the PETSc local-to-global map. Compute outside and pass
136 // into SparsityPattern constructor.
137 // TODO: avoid concatenating the same maps twice in case that V[0]
138 // == V[1].
139
140 const std::vector<
141 std::pair<std::reference_wrapper<const common::IndexMap>, int>>& map
142 = maps[d];
143 std::vector<PetscInt>& _map = _maps[d];
144
145 // Concatenate the block index map in the row and column directions
146 const auto [rank_offset, local_offset, ghosts, _]
148 for (std::size_t f = 0; f < map.size(); ++f)
149 {
150 auto offset = local_offset[f];
151 const common::IndexMap& imap = map[f].first.get();
152 int bs = map[f].second;
153 for (std::int32_t i = 0; i < bs * imap.size_local(); ++i)
154 _map.push_back(i + rank_offset + offset);
155 for (std::int32_t i = 0; i < bs * imap.num_ghosts(); ++i)
156 _map.push_back(ghosts[f][i]);
157 }
158 }
159
160 // Create PETSc local-to-global map/index sets and attach to matrix
161 ISLocalToGlobalMapping petsc_local_to_global0;
162 ISLocalToGlobalMappingCreate(MPI_COMM_SELF, 1, _maps[0].size(),
163 _maps[0].data(), PETSC_COPY_VALUES,
164 &petsc_local_to_global0);
165 if (V[0] == V[1])
166 {
167 MatSetLocalToGlobalMapping(A, petsc_local_to_global0,
168 petsc_local_to_global0);
169 ISLocalToGlobalMappingDestroy(&petsc_local_to_global0);
170 }
171 else
172 {
173
174 ISLocalToGlobalMapping petsc_local_to_global1;
175 ISLocalToGlobalMappingCreate(MPI_COMM_SELF, 1, _maps[1].size(),
176 _maps[1].data(), PETSC_COPY_VALUES,
177 &petsc_local_to_global1);
178 MatSetLocalToGlobalMapping(A, petsc_local_to_global0,
179 petsc_local_to_global1);
180 ISLocalToGlobalMappingDestroy(&petsc_local_to_global0);
181 ISLocalToGlobalMappingDestroy(&petsc_local_to_global1);
182 }
183
184 return A;
185}
186
190template <std::floating_point T>
192 const std::vector<std::vector<const Form<PetscScalar, T>*>>& a,
193 const std::vector<std::vector<std::string>>& types)
194{
195 // Extract and check row/column ranges
197
198 std::vector<std::vector<std::string>> _types(
199 a.size(), std::vector<std::string>(a.front().size()));
200 if (!types.empty())
201 _types = types;
202
203 // Loop over each form and create matrix
204 int rows = a.size();
205 int cols = a.front().size();
206 std::vector<Mat> mats(rows * cols, nullptr);
207 std::shared_ptr<const mesh::Mesh<T>> mesh;
208 for (int i = 0; i < rows; ++i)
209 {
210 for (int j = 0; j < cols; ++j)
211 {
212 if (const Form<PetscScalar, T>* form = a[i][j]; form)
213 {
214 mats[i * cols + j] = create_matrix(*form, _types[i][j]);
215 mesh = form->mesh();
216 }
217 }
218 }
219
220 if (!mesh)
221 throw std::runtime_error("Could not find a Mesh.");
222
223 // Initialise block (MatNest) matrix
224 Mat A;
225 MatCreate(mesh->comm(), &A);
226 MatSetType(A, MATNEST);
227 MatNestSetSubMats(A, rows, nullptr, cols, nullptr, mats.data());
228 MatSetUp(A);
229
230 // De-reference Mat objects
231 for (std::size_t i = 0; i < mats.size(); ++i)
232 {
233 if (mats[i])
234 MatDestroy(&mats[i]);
235 }
236
237 return A;
238}
239
244 const std::vector<
245 std::pair<std::reference_wrapper<const common::IndexMap>, int>>& maps);
246
249 const std::vector<
250 std::pair<std::reference_wrapper<const common::IndexMap>, int>>& maps);
251
252// -- Vectors ----------------------------------------------------------------
253
266template <std::floating_point T>
268 Vec b, const Form<PetscScalar, T>& L,
269 std::span<const PetscScalar> constants,
270 const std::map<std::pair<IntegralType, int>,
271 std::pair<std::span<const PetscScalar>, int>>& coeffs)
272{
273 Vec b_local;
274 VecGhostGetLocalForm(b, &b_local);
275 PetscInt n = 0;
276 VecGetSize(b_local, &n);
277 PetscScalar* array = nullptr;
278 VecGetArray(b_local, &array);
279 std::span<PetscScalar> _b(array, n);
280 fem::assemble_vector(_b, L, constants, coeffs);
281 VecRestoreArray(b_local, &array);
282 VecGhostRestoreLocalForm(b, &b_local);
283}
284
294template <std::floating_point T>
296{
297 Vec b_local;
298 VecGhostGetLocalForm(b, &b_local);
299 PetscInt n = 0;
300 VecGetSize(b_local, &n);
301 PetscScalar* array = nullptr;
302 VecGetArray(b_local, &array);
303 std::span<PetscScalar> _b(array, n);
305 VecRestoreArray(b_local, &array);
306 VecGhostRestoreLocalForm(b, &b_local);
307}
308
309// FIXME: clarify how x0 is used
310// FIXME: if bcs entries are set
311
312// FIXME: need to pass an array of Vec for x0?
313// FIXME: clarify zeroing of vector
314
330template <std::floating_point T>
332 Vec b, const std::vector<std::shared_ptr<const Form<PetscScalar, T>>>& a,
333 const std::vector<std::span<const PetscScalar>>& constants,
334 const std::vector<std::map<std::pair<IntegralType, int>,
335 std::pair<std::span<const PetscScalar>, int>>>&
336 coeffs,
337 const std::vector<
338 std::vector<std::shared_ptr<const DirichletBC<PetscScalar, T>>>>& bcs1,
339 const std::vector<Vec>& x0, PetscScalar alpha)
340{
341 Vec b_local;
342 VecGhostGetLocalForm(b, &b_local);
343 PetscInt n = 0;
344 VecGetSize(b_local, &n);
345 PetscScalar* array = nullptr;
346 VecGetArray(b_local, &array);
347 std::span<PetscScalar> _b(array, n);
348
349 if (x0.empty())
350 fem::apply_lifting(_b, a, constants, coeffs, bcs1, {}, alpha);
351 else
352 {
353 std::vector<std::span<const PetscScalar>> x0_ref;
354 std::vector<Vec> x0_local(a.size());
355 std::vector<const PetscScalar*> x0_array(a.size());
356 for (std::size_t i = 0; i < a.size(); ++i)
357 {
358 assert(x0[i]);
359 VecGhostGetLocalForm(x0[i], &x0_local[i]);
360 PetscInt n = 0;
361 VecGetSize(x0_local[i], &n);
362 VecGetArrayRead(x0_local[i], &x0_array[i]);
363 x0_ref.emplace_back(x0_array[i], n);
364 }
365
366 std::vector x0_tmp(x0_ref.begin(), x0_ref.end());
367 fem::apply_lifting(_b, a, constants, coeffs, bcs1, x0_tmp, alpha);
368
369 for (std::size_t i = 0; i < x0_local.size(); ++i)
370 {
371 VecRestoreArrayRead(x0_local[i], &x0_array[i]);
372 VecGhostRestoreLocalForm(x0[i], &x0_local[i]);
373 }
374 }
375
376 VecRestoreArray(b_local, &array);
377 VecGhostRestoreLocalForm(b, &b_local);
378}
379
380// FIXME: clarify how x0 is used
381// FIXME: if bcs entries are set
382
383// FIXME: need to pass an array of Vec for x0?
384// FIXME: clarify zeroing of vector
385
398template <std::floating_point T>
400 Vec b,
401 const std::vector<std::shared_ptr<const Form<PetscScalar, double>>>& a,
402 const std::vector<
403 std::vector<std::shared_ptr<const DirichletBC<PetscScalar, double>>>>&
404 bcs1,
405 const std::vector<Vec>& x0, PetscScalar alpha)
406{
407 Vec b_local;
408 VecGhostGetLocalForm(b, &b_local);
409 PetscInt n = 0;
410 VecGetSize(b_local, &n);
411 PetscScalar* array = nullptr;
412 VecGetArray(b_local, &array);
413 std::span<PetscScalar> _b(array, n);
414
415 if (x0.empty())
416 fem::apply_lifting<PetscScalar>(_b, a, bcs1, {}, alpha);
417 else
418 {
419 std::vector<std::span<const PetscScalar>> x0_ref;
420 std::vector<Vec> x0_local(a.size());
421 std::vector<const PetscScalar*> x0_array(a.size());
422 for (std::size_t i = 0; i < a.size(); ++i)
423 {
424 assert(x0[i]);
425 VecGhostGetLocalForm(x0[i], &x0_local[i]);
426 PetscInt n = 0;
427 VecGetSize(x0_local[i], &n);
428 VecGetArrayRead(x0_local[i], &x0_array[i]);
429 x0_ref.emplace_back(x0_array[i], n);
430 }
431
432 std::vector x0_tmp(x0_ref.begin(), x0_ref.end());
433 fem::apply_lifting<PetscScalar>(_b, a, bcs1, x0_tmp, alpha);
434
435 for (std::size_t i = 0; i < x0_local.size(); ++i)
436 {
437 VecRestoreArrayRead(x0_local[i], &x0_array[i]);
438 VecGhostRestoreLocalForm(x0[i], &x0_local[i]);
439 }
440 }
441
442 VecRestoreArray(b_local, &array);
443 VecGhostRestoreLocalForm(b, &b_local);
444}
445
446// -- Setting bcs ------------------------------------------------------------
447
448// FIXME: Move these function elsewhere?
449
450// FIXME: clarify x0
451// FIXME: clarify what happens with ghosts
452
456template <std::floating_point T>
458 Vec b,
459 const std::vector<std::shared_ptr<const DirichletBC<PetscScalar, T>>>& bcs,
460 const Vec x0, PetscScalar alpha = 1)
461{
462 PetscInt n = 0;
463 VecGetLocalSize(b, &n);
464 PetscScalar* array = nullptr;
465 VecGetArray(b, &array);
466 std::span<PetscScalar> _b(array, n);
467 if (x0)
468 {
469 Vec x0_local;
470 VecGhostGetLocalForm(x0, &x0_local);
471 PetscInt n = 0;
472 VecGetSize(x0_local, &n);
473 const PetscScalar* array = nullptr;
474 VecGetArrayRead(x0_local, &array);
475 std::span<const PetscScalar> _x0(array, n);
476 for (auto& bc : bcs)
477 bc->set(_b, _x0, alpha);
478 VecRestoreArrayRead(x0_local, &array);
479 VecGhostRestoreLocalForm(x0, &x0_local);
480 }
481 else
482 {
483 for (auto& bc : bcs)
484 bc->set(_b, std::nullopt, alpha);
485 }
486 VecRestoreArray(b, &array);
487}
488
489} // namespace petsc
490} // namespace dolfinx::fem
491
492#endif
Definition IndexMap.h:94
std::int32_t num_ghosts() const noexcept
Number of ghost indices on this process.
Definition IndexMap.cpp:929
std::int32_t size_local() const noexcept
Number of indices owned by this process.
Definition IndexMap.cpp:931
Definition petsc.h:32
A representation of finite element variational forms.
Definition Form.h:139
Definition SparsityPattern.h:26
void finalize()
Finalize sparsity pattern and communicate off-process entries.
Definition SparsityPattern.cpp:245
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:648
Mat create_matrix(const Form< PetscScalar, T > &a, std::string type=std::string())
Create a matrix.
Definition petsc.h:44
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 alpha)
Modify RHS vector to account for Dirichlet boundary conditions.
Definition petsc.h:331
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:267
void set_bc(Vec b, const std::vector< std::shared_ptr< const DirichletBC< PetscScalar, T > > > &bcs, const Vec x0, PetscScalar alpha=1)
Definition petsc.h:457
Mat create_matrix_block(const std::vector< std::vector< const Form< PetscScalar, T > * > > &a, std::string type=std::string())
Initialise a monolithic matrix for an array of bilinear forms.
Definition petsc.h:62
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:191
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:111
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:122
la::SparsityPattern create_sparsity_pattern(const Form< T, U > &a)
Create a sparsity pattern for a given form.
Definition utils.h:147
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 alpha)
Definition assembler.h:152
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:235