13#include <dolfinx/la/petsc.h>
29template <dolfinx::scalar T, std::
floating_po
int U>
41template <std::
floating_po
int T>
43 const std::string& type = std::string())
58template <std::
floating_po
int T>
61 const std::string& type = std::string())
64 std::array<std::vector<std::shared_ptr<const FunctionSpace<T>>>, 2> V
66 std::array<std::vector<int>, 2> bs_dofs;
67 for (std::size_t i = 0; i < 2; ++i)
70 bs_dofs[i].push_back(_V->dofmap()->bs());
74 std::shared_ptr<const mesh::Mesh<T>> mesh;
75 std::vector<std::vector<std::unique_ptr<la::SparsityPattern>>> patterns(
77 for (std::size_t row = 0; row < V[0].size(); ++row)
79 for (std::size_t col = 0; col < V[1].size(); ++col)
83 patterns[row].push_back(std::make_unique<la::SparsityPattern>(
89 patterns[row].push_back(
nullptr);
94 throw std::runtime_error(
"Could not find a Mesh.");
97 std::array<std::vector<std::pair<
98 std::reference_wrapper<const common::IndexMap>,
int>>,
101 for (std::size_t d = 0; d < 2; ++d)
103 for (
auto space : V[d])
105 maps[d].emplace_back(*space->dofmap()->index_map,
106 space->dofmap()->index_map_bs());
111 std::vector<std::vector<const la::SparsityPattern*>> p(V[0].size());
112 for (std::size_t row = 0; row < V[0].size(); ++row)
113 for (std::size_t col = 0; col < V[1].size(); ++col)
114 p[row].push_back(patterns[row][col].get());
127 std::array<std::vector<PetscInt>, 2> _maps;
128 for (
int d = 0; d < 2; ++d)
139 auto [rank_offset, local_offset, ghosts, _]
141 for (std::size_t f = 0; f < maps[d].size(); ++f)
144 const int bs = maps[d][f].second;
145 const std::int32_t size_local = bs * map.
size_local();
147 for (std::int32_t i = 0; i < size_local; ++i)
148 _maps[d].push_back(i + rank_offset + local_offset[f]);
149 for (std::size_t i = size_local; i < bs * global.size(); ++i)
150 _maps[d].push_back(ghosts[f][i - size_local]);
155 ISLocalToGlobalMapping petsc_local_to_global0;
156 ISLocalToGlobalMappingCreate(MPI_COMM_SELF, 1, _maps[0].size(),
157 _maps[0].data(), PETSC_COPY_VALUES,
158 &petsc_local_to_global0);
161 MatSetLocalToGlobalMapping(A, petsc_local_to_global0,
162 petsc_local_to_global0);
163 ISLocalToGlobalMappingDestroy(&petsc_local_to_global0);
168 ISLocalToGlobalMapping petsc_local_to_global1;
169 ISLocalToGlobalMappingCreate(MPI_COMM_SELF, 1, _maps[1].size(),
170 _maps[1].data(), PETSC_COPY_VALUES,
171 &petsc_local_to_global1);
172 MatSetLocalToGlobalMapping(A, petsc_local_to_global0,
173 petsc_local_to_global1);
174 ISLocalToGlobalMappingDestroy(&petsc_local_to_global0);
175 ISLocalToGlobalMappingDestroy(&petsc_local_to_global1);
184template <std::
floating_po
int T>
187 const std::vector<std::vector<std::string>>& types)
192 std::vector<std::vector<std::string>> _types(
193 a.size(), std::vector<std::string>(a.front().size()));
199 int cols = a.front().size();
200 std::vector<Mat> mats(rows * cols,
nullptr);
201 std::shared_ptr<const mesh::Mesh<T>> mesh;
202 for (
int i = 0; i < rows; ++i)
204 for (
int j = 0; j < cols; ++j)
215 throw std::runtime_error(
"Could not find a Mesh.");
219 MatCreate(mesh->comm(), &A);
220 MatSetType(A, MATNEST);
221 MatNestSetSubMats(A, rows,
nullptr, cols,
nullptr, mats.data());
225 for (std::size_t i = 0; i < mats.size(); ++i)
228 MatDestroy(&mats[i]);
239 std::pair<std::reference_wrapper<const common::IndexMap>,
int>>& maps);
244 std::pair<std::reference_wrapper<const common::IndexMap>,
int>>& maps);
260template <std::
floating_po
int T>
263 std::span<const PetscScalar> constants,
264 const std::map<std::pair<IntegralType, int>,
265 std::pair<std::span<const PetscScalar>,
int>>& coeffs)
268 VecGhostGetLocalForm(b, &b_local);
270 VecGetSize(b_local, &n);
271 PetscScalar* array =
nullptr;
272 VecGetArray(b_local, &array);
273 std::span<PetscScalar> _b(array, n);
275 VecRestoreArray(b_local, &array);
276 VecGhostRestoreLocalForm(b, &b_local);
288template <std::
floating_po
int T>
292 VecGhostGetLocalForm(b, &b_local);
294 VecGetSize(b_local, &n);
295 PetscScalar* array =
nullptr;
296 VecGetArray(b_local, &array);
297 std::span<PetscScalar> _b(array, n);
299 VecRestoreArray(b_local, &array);
300 VecGhostRestoreLocalForm(b, &b_local);
324template <std::
floating_po
int T>
327 const std::vector<std::span<const PetscScalar>>& constants,
328 const std::vector<std::map<std::pair<IntegralType, int>,
329 std::pair<std::span<const PetscScalar>,
int>>>&
333 const std::vector<Vec>& x0, PetscScalar scale)
336 VecGhostGetLocalForm(b, &b_local);
338 VecGetSize(b_local, &n);
339 PetscScalar* array =
nullptr;
340 VecGetArray(b_local, &array);
341 std::span<PetscScalar> _b(array, n);
347 std::vector<std::span<const PetscScalar>> x0_ref;
348 std::vector<Vec> x0_local(a.size());
349 std::vector<const PetscScalar*> x0_array(a.size());
350 for (std::size_t i = 0; i < a.size(); ++i)
353 VecGhostGetLocalForm(x0[i], &x0_local[i]);
355 VecGetSize(x0_local[i], &n);
356 VecGetArrayRead(x0_local[i], &x0_array[i]);
357 x0_ref.emplace_back(x0_array[i], n);
360 std::vector x0_tmp(x0_ref.begin(), x0_ref.end());
363 for (std::size_t i = 0; i < x0_local.size(); ++i)
365 VecRestoreArrayRead(x0_local[i], &x0_array[i]);
366 VecGhostRestoreLocalForm(x0[i], &x0_local[i]);
370 VecRestoreArray(b_local, &array);
371 VecGhostRestoreLocalForm(b, &b_local);
392template <std::
floating_po
int T>
399 const std::vector<Vec>& x0, PetscScalar scale)
402 VecGhostGetLocalForm(b, &b_local);
404 VecGetSize(b_local, &n);
405 PetscScalar* array =
nullptr;
406 VecGetArray(b_local, &array);
407 std::span<PetscScalar> _b(array, n);
410 fem::apply_lifting<PetscScalar>(_b, a, bcs1, {}, scale);
413 std::vector<std::span<const PetscScalar>> x0_ref;
414 std::vector<Vec> x0_local(a.size());
415 std::vector<const PetscScalar*> x0_array(a.size());
416 for (std::size_t i = 0; i < a.size(); ++i)
419 VecGhostGetLocalForm(x0[i], &x0_local[i]);
421 VecGetSize(x0_local[i], &n);
422 VecGetArrayRead(x0_local[i], &x0_array[i]);
423 x0_ref.emplace_back(x0_array[i], n);
426 std::vector x0_tmp(x0_ref.begin(), x0_ref.end());
427 fem::apply_lifting<PetscScalar>(_b, a, bcs1, x0_tmp, scale);
429 for (std::size_t i = 0; i < x0_local.size(); ++i)
431 VecRestoreArrayRead(x0_local[i], &x0_array[i]);
432 VecGhostRestoreLocalForm(x0[i], &x0_local[i]);
436 VecRestoreArray(b_local, &array);
437 VecGhostRestoreLocalForm(b, &b_local);
450template <std::
floating_po
int T>
454 const Vec x0, PetscScalar scale = 1)
457 VecGetLocalSize(b, &n);
458 PetscScalar* array =
nullptr;
459 VecGetArray(b, &array);
460 std::span<PetscScalar> _b(array, n);
464 VecGhostGetLocalForm(x0, &x0_local);
466 VecGetSize(x0_local, &n);
467 const PetscScalar* array =
nullptr;
468 VecGetArrayRead(x0_local, &array);
469 std::span<const PetscScalar> _x0(array, n);
471 VecRestoreArrayRead(x0_local, &array);
472 VecGhostRestoreLocalForm(x0, &x0_local);
477 VecRestoreArray(b, &array);
This class represents the distribution index arrays across processes. An index array is a contiguous ...
Definition IndexMap.h:64
std::int32_t size_local() const noexcept
Number of indices owned by this process.
Definition IndexMap.cpp:391
std::vector< std::int64_t > global_indices() const
Build list of indices with global indexing.
Definition IndexMap.cpp:448
Object for setting (strong) Dirichlet boundary conditions.
Definition DirichletBC.h:251
Sparsity pattern data structure that can be used to initialize sparse matrices. After assembly,...
Definition SparsityPattern.h:26
void finalize()
Finalize sparsity pattern and communicate off-process entries.
Definition SparsityPattern.cpp:226
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 common::IndexMap >, int > > &maps)
Compute layout data and ghost indices for a stacked (concatenated) index map, i.e....
Definition IndexMap.cpp:145
void set_bc(Vec b, const std::vector< std::shared_ptr< const DirichletBC< PetscScalar, T > > > &bcs, const Vec x0, PetscScalar scale=1)
Set bc values in owned (local) part of the PETSc vector, multiplied by 'scale'. The vectors b and x0 ...
Definition petsc.h:451
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:325
Mat create_matrix_block(const std::vector< std::vector< const Form< PetscScalar, T > * > > &a, const std::string &type=std::string())
Initialise a monolithic matrix for an array of bilinear forms.
Definition petsc.h:59
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:261
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:185
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:15
Mat create_matrix(const Form< PetscScalar, T > &a, const std::string &type=std::string())
Create a matrix.
Definition petsc.h:42
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:60
Finite element method functionality.
Definition assemble_matrix_impl.h:25
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:109
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)
Modify b such that:
Definition assembler.h:150
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:135
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)
Set bc values in owned (local) part of the vector, multiplied by 'scale'. The vectors b and x0 must h...
Definition assembler.h:438
la::SparsityPattern create_sparsity_pattern(const Form< T, U > &a)
Create a sparsity pattern for a given form.
Definition utils.h:160
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:343
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