15#include <dolfinx/la/petsc.h>
33template <dolfinx::scalar T, std::
floating_po
int U>
45template <std::
floating_po
int T>
47 std::optional<std::string> type = std::nullopt)
64template <std::
floating_po
int T>
67 std::optional<std::string> type = std::nullopt)
70 std::array<std::vector<std::shared_ptr<const FunctionSpace<T>>>, 2> V
72 std::array<std::vector<int>, 2> bs_dofs;
73 for (std::size_t i = 0; i < 2; ++i)
76 bs_dofs[i].push_back(_V->dofmap()->bs());
80 std::shared_ptr<const mesh::Mesh<T>>
mesh;
81 std::vector<std::vector<std::unique_ptr<la::SparsityPattern>>> patterns(
83 for (std::size_t row = 0; row < V[0].size(); ++row)
85 for (std::size_t col = 0; col < V[1].size(); ++col)
89 patterns[row].push_back(std::make_unique<la::SparsityPattern>(
95 patterns[row].push_back(
nullptr);
100 throw std::runtime_error(
"Could not find a Mesh.");
103 std::array<std::vector<std::pair<
104 std::reference_wrapper<const common::IndexMap>,
int>>,
107 for (std::size_t d = 0; d < 2; ++d)
109 for (
auto& space : V[d])
111 maps[d].emplace_back(*space->dofmap()->index_map,
112 space->dofmap()->index_map_bs());
117 std::vector<std::vector<const la::SparsityPattern*>> p(V[0].size());
118 for (std::size_t row = 0; row < V[0].size(); ++row)
119 for (std::size_t col = 0; col < V[1].size(); ++col)
120 p[row].push_back(patterns[row][col].get());
133 std::array<std::vector<PetscInt>, 2> _maps;
134 for (
int d = 0; d < 2; ++d)
144 std::pair<std::reference_wrapper<const common::IndexMap>,
int>>& map
146 std::vector<PetscInt>& _map = _maps[d];
149 const auto [rank_offset, local_offset, ghosts, _]
151 for (std::size_t f = 0; f < map.size(); ++f)
153 auto offset = local_offset[f];
155 int bs = map[f].second;
156 for (std::int32_t i = 0; i < bs * imap.
size_local(); ++i)
157 _map.push_back(i + rank_offset + offset);
158 for (std::int32_t i = 0; i < bs * imap.
num_ghosts(); ++i)
159 _map.push_back(ghosts[f][i]);
164 ISLocalToGlobalMapping petsc_local_to_global0;
165 ISLocalToGlobalMappingCreate(MPI_COMM_SELF, 1, _maps[0].size(),
166 _maps[0].data(), PETSC_COPY_VALUES,
167 &petsc_local_to_global0);
170 MatSetLocalToGlobalMapping(A, petsc_local_to_global0,
171 petsc_local_to_global0);
172 ISLocalToGlobalMappingDestroy(&petsc_local_to_global0);
177 ISLocalToGlobalMapping petsc_local_to_global1;
178 ISLocalToGlobalMappingCreate(MPI_COMM_SELF, 1, _maps[1].size(),
179 _maps[1].data(), PETSC_COPY_VALUES,
180 &petsc_local_to_global1);
181 MatSetLocalToGlobalMapping(A, petsc_local_to_global0,
182 petsc_local_to_global1);
183 ISLocalToGlobalMappingDestroy(&petsc_local_to_global0);
184 ISLocalToGlobalMappingDestroy(&petsc_local_to_global1);
193template <std::
floating_po
int T>
196 std::optional<std::vector<std::vector<std::optional<std::string>>>> types)
203 int cols = a.front().size();
204 std::vector<Mat> mats(rows * cols,
nullptr);
205 std::shared_ptr<const mesh::Mesh<T>>
mesh;
206 for (
int i = 0; i < rows; ++i)
208 for (
int j = 0; j < cols; ++j)
213 mats[i * cols + j] =
create_matrix(*form, types->at(i).at(j));
222 throw std::runtime_error(
"Could not find a Mesh.");
226 MatCreate(
mesh->comm(), &A);
227 MatSetType(A, MATNEST);
228 MatNestSetSubMats(A, rows,
nullptr, cols,
nullptr, mats.data());
232 for (std::size_t i = 0; i < mats.size(); ++i)
235 MatDestroy(&mats[i]);
246 std::pair<std::reference_wrapper<const common::IndexMap>,
int>>& maps);
251 std::pair<std::reference_wrapper<const common::IndexMap>,
int>>& maps);
267template <std::
floating_po
int T>
270 std::span<const PetscScalar> constants,
271 const std::map<std::pair<IntegralType, int>,
272 std::pair<std::span<const PetscScalar>,
int>>& coeffs)
275 VecGhostGetLocalForm(b, &b_local);
277 VecGetSize(b_local, &n);
278 PetscScalar* array =
nullptr;
279 VecGetArray(b_local, &array);
280 std::span<PetscScalar> _b(array, n);
282 VecRestoreArray(b_local, &array);
283 VecGhostRestoreLocalForm(b, &b_local);
296template <std::
floating_po
int T>
300 VecGhostGetLocalForm(b, &b_local);
302 VecGetSize(b_local, &n);
303 PetscScalar* array =
nullptr;
304 VecGetArray(b_local, &array);
305 std::span<PetscScalar> _b(array, n);
307 VecRestoreArray(b_local, &array);
308 VecGhostRestoreLocalForm(b, &b_local);
332template <std::
floating_po
int T>
338 const std::vector<std::span<const PetscScalar>>& constants,
339 const std::vector<std::map<std::pair<IntegralType, int>,
340 std::pair<std::span<const PetscScalar>,
int>>>&
345 const std::vector<Vec>& x0, PetscScalar alpha)
348 VecGhostGetLocalForm(b, &b_local);
350 VecGetSize(b_local, &n);
351 PetscScalar* array =
nullptr;
352 VecGetArray(b_local, &array);
353 std::span<PetscScalar> _b(array, n);
359 std::vector<std::span<const PetscScalar>> x0_ref;
360 std::vector<Vec> x0_local(a.size());
361 std::vector<const PetscScalar*> x0_array(a.size());
362 for (std::size_t i = 0; i < a.size(); ++i)
365 VecGhostGetLocalForm(x0[i], &x0_local[i]);
367 VecGetSize(x0_local[i], &n);
368 VecGetArrayRead(x0_local[i], &x0_array[i]);
369 x0_ref.emplace_back(x0_array[i], n);
372 std::vector x0_tmp(x0_ref.begin(), x0_ref.end());
375 for (std::size_t i = 0; i < x0_local.size(); ++i)
377 VecRestoreArrayRead(x0_local[i], &x0_array[i]);
378 VecGhostRestoreLocalForm(x0[i], &x0_local[i]);
382 VecRestoreArray(b_local, &array);
383 VecGhostRestoreLocalForm(b, &b_local);
404template <std::
floating_po
int T>
410 const std::vector<std::vector<
412 const std::vector<Vec>& x0, PetscScalar alpha)
415 VecGhostGetLocalForm(b, &b_local);
417 VecGetSize(b_local, &n);
418 PetscScalar* array =
nullptr;
419 VecGetArray(b_local, &array);
420 std::span<PetscScalar> _b(array, n);
426 std::vector<std::span<const PetscScalar>> x0_ref;
427 std::vector<Vec> x0_local(a.size());
428 std::vector<const PetscScalar*> x0_array(a.size());
429 for (std::size_t i = 0; i < a.size(); ++i)
432 VecGhostGetLocalForm(x0[i], &x0_local[i]);
434 VecGetSize(x0_local[i], &n);
435 VecGetArrayRead(x0_local[i], &x0_array[i]);
436 x0_ref.emplace_back(x0_array[i], n);
439 std::vector x0_tmp(x0_ref.begin(), x0_ref.end());
442 for (std::size_t i = 0; i < x0_local.size(); ++i)
444 VecRestoreArrayRead(x0_local[i], &x0_array[i]);
445 VecGhostRestoreLocalForm(x0[i], &x0_local[i]);
449 VecRestoreArray(b_local, &array);
450 VecGhostRestoreLocalForm(b, &b_local);
463template <std::
floating_po
int T>
468 const Vec x0, PetscScalar alpha = 1)
471 VecGetLocalSize(b, &n);
472 PetscScalar* array =
nullptr;
473 VecGetArray(b, &array);
474 std::span<PetscScalar> _b(array, n);
478 VecGhostGetLocalForm(x0, &x0_local);
480 VecGetSize(x0_local, &n);
481 const PetscScalar* array =
nullptr;
482 VecGetArrayRead(x0_local, &array);
483 std::span<const PetscScalar> _x0(array, n);
485 bc.set(_b, _x0, alpha);
486 VecRestoreArrayRead(x0_local, &array);
487 VecGhostRestoreLocalForm(x0, &x0_local);
492 bc->set(_b, std::nullopt, alpha);
494 VecRestoreArray(b, &array);
Functions supporting assembly of finite element fem::Form and fem::Expression.
std::int32_t num_ghosts() const noexcept
Number of ghost indices on this process.
Definition IndexMap.cpp:931
std::int32_t size_local() const noexcept
Number of indices owned by this process.
Definition IndexMap.cpp:933
Definition DirichletBC.h:258
Definition SparsityPattern.h:26
void finalize()
Finalize sparsity pattern and communicate off-process entries.
Definition SparsityPattern.cpp:244
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:650
Helper functions for assembly into PETSc data structures.
Definition petsc.h:38
Mat create_matrix(const Form< PetscScalar, T > &a, std::optional< std::string > type=std::nullopt)
Create a matrix.
Definition petsc.h:46
Mat create_matrix_block(const std::vector< std::vector< const Form< PetscScalar, T > * > > &a, std::optional< std::string > type=std::nullopt)
Initialise a monolithic matrix for an array of bilinear forms.
Definition petsc.h:65
void set_bc(Vec b, const std::vector< std::reference_wrapper< const DirichletBC< PetscScalar, T > > > bcs, const Vec x0, PetscScalar alpha=1)
Definition petsc.h:464
Mat create_matrix_nest(const std::vector< std::vector< const Form< PetscScalar, T > * > > &a, std::optional< std::vector< std::vector< std::optional< std::string > > > > types)
Create nested (MatNest) matrix.
Definition petsc.h:194
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:268
Vec create_vector_block(const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int > > &maps)
Initialise monolithic vector. Vector is not zeroed.
Definition petsc.cpp:17
void apply_lifting(Vec b, std::vector< std::optional< std::reference_wrapper< 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::reference_wrapper< const DirichletBC< PetscScalar, T > > > > &bcs1, const std::vector< Vec > &x0, PetscScalar alpha)
Modify RHS vector to account for Dirichlet boundary conditions.
Definition petsc.h:333
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_expression_impl.h:23
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:172
void assemble_vector(V &&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:218
void apply_lifting(V &&b, std::vector< std::optional< std::reference_wrapper< 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::reference_wrapper< const DirichletBC< T, U > > > > &bcs1, const std::vector< std::span< const T > > &x0, T alpha)
Modify the right-hand side vector to account for constraints (Dirichlet boundary condition constraint...
Definition assembler.h:326
la::SparsityPattern create_sparsity_pattern(const Form< T, U > &a)
Create a sparsity pattern for a given form.
Definition utils.h:197
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)
Extract FunctionSpaces for (0) rows blocks and (1) columns blocks from a rectangular array of (test,...
Definition FunctionSpace.h:433
Mat create_matrix(MPI_Comm comm, const SparsityPattern &sp, std::optional< std::string > type=std::nullopt)
Definition petsc.cpp:235
Mesh data structures and algorithms on meshes.
Definition DofMap.h:32