15#include <dolfinx/la/petsc.h>
31template <dolfinx::scalar T, std::
floating_po
int U>
43template <std::
floating_po
int T>
45 std::string type = std::string())
61template <std::
floating_po
int T>
64 std::string type = std::string())
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)
73 bs_dofs[i].push_back(_V->dofmap()->bs());
77 std::shared_ptr<const mesh::Mesh<T>> mesh;
78 std::vector<std::vector<std::unique_ptr<la::SparsityPattern>>> patterns(
80 for (std::size_t row = 0; row < V[0].size(); ++row)
82 for (std::size_t col = 0; col < V[1].size(); ++col)
86 patterns[row].push_back(std::make_unique<la::SparsityPattern>(
92 patterns[row].push_back(
nullptr);
97 throw std::runtime_error(
"Could not find a Mesh.");
100 std::array<std::vector<std::pair<
101 std::reference_wrapper<const common::IndexMap>,
int>>,
104 for (std::size_t d = 0; d < 2; ++d)
106 for (
auto& space : V[d])
108 maps[d].emplace_back(*space->dofmap()->index_map,
109 space->dofmap()->index_map_bs());
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());
130 std::array<std::vector<PetscInt>, 2> _maps;
131 for (
int d = 0; d < 2; ++d)
141 std::pair<std::reference_wrapper<const common::IndexMap>,
int>>& map
143 std::vector<PetscInt>& _map = _maps[d];
146 const auto [rank_offset, local_offset, ghosts, _]
148 for (std::size_t f = 0; f < map.size(); ++f)
150 auto offset = local_offset[f];
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]);
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);
167 MatSetLocalToGlobalMapping(A, petsc_local_to_global0,
168 petsc_local_to_global0);
169 ISLocalToGlobalMappingDestroy(&petsc_local_to_global0);
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);
190template <std::
floating_po
int T>
193 const std::vector<std::vector<std::string>>& types)
198 std::vector<std::vector<std::string>> _types(
199 a.size(), std::vector<std::string>(a.front().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)
210 for (
int j = 0; j < cols; ++j)
221 throw std::runtime_error(
"Could not find a Mesh.");
225 MatCreate(mesh->comm(), &A);
226 MatSetType(A, MATNEST);
227 MatNestSetSubMats(A, rows,
nullptr, cols,
nullptr, mats.data());
231 for (std::size_t i = 0; i < mats.size(); ++i)
234 MatDestroy(&mats[i]);
245 std::pair<std::reference_wrapper<const common::IndexMap>,
int>>& maps);
250 std::pair<std::reference_wrapper<const common::IndexMap>,
int>>& maps);
266template <std::
floating_po
int T>
269 std::span<const PetscScalar> constants,
270 const std::map<std::pair<IntegralType, int>,
271 std::pair<std::span<const PetscScalar>,
int>>& coeffs)
274 VecGhostGetLocalForm(b, &b_local);
276 VecGetSize(b_local, &n);
277 PetscScalar* array =
nullptr;
278 VecGetArray(b_local, &array);
279 std::span<PetscScalar> _b(array, n);
281 VecRestoreArray(b_local, &array);
282 VecGhostRestoreLocalForm(b, &b_local);
294template <std::
floating_po
int T>
298 VecGhostGetLocalForm(b, &b_local);
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);
330template <std::
floating_po
int T>
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>>>&
339 const std::vector<Vec>& x0, PetscScalar alpha)
342 VecGhostGetLocalForm(b, &b_local);
344 VecGetSize(b_local, &n);
345 PetscScalar* array =
nullptr;
346 VecGetArray(b_local, &array);
347 std::span<PetscScalar> _b(array, n);
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)
359 VecGhostGetLocalForm(x0[i], &x0_local[i]);
361 VecGetSize(x0_local[i], &n);
362 VecGetArrayRead(x0_local[i], &x0_array[i]);
363 x0_ref.emplace_back(x0_array[i], n);
366 std::vector x0_tmp(x0_ref.begin(), x0_ref.end());
369 for (std::size_t i = 0; i < x0_local.size(); ++i)
371 VecRestoreArrayRead(x0_local[i], &x0_array[i]);
372 VecGhostRestoreLocalForm(x0[i], &x0_local[i]);
376 VecRestoreArray(b_local, &array);
377 VecGhostRestoreLocalForm(b, &b_local);
398template <std::
floating_po
int T>
405 const std::vector<Vec>& x0, PetscScalar alpha)
408 VecGhostGetLocalForm(b, &b_local);
410 VecGetSize(b_local, &n);
411 PetscScalar* array =
nullptr;
412 VecGetArray(b_local, &array);
413 std::span<PetscScalar> _b(array, n);
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)
425 VecGhostGetLocalForm(x0[i], &x0_local[i]);
427 VecGetSize(x0_local[i], &n);
428 VecGetArrayRead(x0_local[i], &x0_array[i]);
429 x0_ref.emplace_back(x0_array[i], n);
432 std::vector x0_tmp(x0_ref.begin(), x0_ref.end());
435 for (std::size_t i = 0; i < x0_local.size(); ++i)
437 VecRestoreArrayRead(x0_local[i], &x0_array[i]);
438 VecGhostRestoreLocalForm(x0[i], &x0_local[i]);
442 VecRestoreArray(b_local, &array);
443 VecGhostRestoreLocalForm(b, &b_local);
456template <std::
floating_po
int T>
460 const Vec x0, PetscScalar alpha = 1)
463 VecGetLocalSize(b, &n);
464 PetscScalar* array =
nullptr;
465 VecGetArray(b, &array);
466 std::span<PetscScalar> _b(array, n);
470 VecGhostGetLocalForm(x0, &x0_local);
472 VecGetSize(x0_local, &n);
473 const PetscScalar* array =
nullptr;
474 VecGetArrayRead(x0_local, &array);
475 std::span<const PetscScalar> _x0(array, n);
477 bc->set(_b, _x0, alpha);
478 VecRestoreArrayRead(x0_local, &array);
479 VecGhostRestoreLocalForm(x0, &x0_local);
484 bc->set(_b, std::nullopt, alpha);
486 VecRestoreArray(b, &array);
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 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