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())
60template <std::
floating_po
int T>
63 std::string type = std::string())
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)
72 bs_dofs[i].push_back(_V->dofmap()->bs());
76 std::shared_ptr<const mesh::Mesh<T>> mesh;
77 std::vector<std::vector<std::unique_ptr<la::SparsityPattern>>> patterns(
79 for (std::size_t row = 0; row < V[0].size(); ++row)
81 for (std::size_t col = 0; col < V[1].size(); ++col)
85 patterns[row].push_back(std::make_unique<la::SparsityPattern>(
91 patterns[row].push_back(
nullptr);
96 throw std::runtime_error(
"Could not find a Mesh.");
99 std::array<std::vector<std::pair<
100 std::reference_wrapper<const common::IndexMap>,
int>>,
103 for (std::size_t d = 0; d < 2; ++d)
105 for (
auto& space : V[d])
107 maps[d].emplace_back(*space->dofmap()->index_map,
108 space->dofmap()->index_map_bs());
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());
129 std::array<std::vector<PetscInt>, 2> _maps;
130 for (
int d = 0; d < 2; ++d)
140 std::pair<std::reference_wrapper<const common::IndexMap>,
int>>& map
142 std::vector<PetscInt>& _map = _maps[d];
145 const auto [rank_offset, local_offset, ghosts, _]
147 for (std::size_t f = 0; f < map.size(); ++f)
149 auto offset = local_offset[f];
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]);
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);
166 MatSetLocalToGlobalMapping(A, petsc_local_to_global0,
167 petsc_local_to_global0);
168 ISLocalToGlobalMappingDestroy(&petsc_local_to_global0);
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);
189template <std::
floating_po
int T>
192 const std::vector<std::vector<std::string>>& types)
197 std::vector<std::vector<std::string>> _types(
198 a.size(), std::vector<std::string>(a.front().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)
209 for (
int j = 0; j < cols; ++j)
220 throw std::runtime_error(
"Could not find a Mesh.");
224 MatCreate(mesh->comm(), &A);
225 MatSetType(A, MATNEST);
226 MatNestSetSubMats(A, rows,
nullptr, cols,
nullptr, mats.data());
230 for (std::size_t i = 0; i < mats.size(); ++i)
233 MatDestroy(&mats[i]);
244 std::pair<std::reference_wrapper<const common::IndexMap>,
int>>& maps);
249 std::pair<std::reference_wrapper<const common::IndexMap>,
int>>& maps);
265template <std::
floating_po
int T>
268 std::span<const PetscScalar> constants,
269 const std::map<std::pair<IntegralType, int>,
270 std::pair<std::span<const PetscScalar>,
int>>& coeffs)
273 VecGhostGetLocalForm(b, &b_local);
275 VecGetSize(b_local, &n);
276 PetscScalar* array =
nullptr;
277 VecGetArray(b_local, &array);
278 std::span<PetscScalar> _b(array, n);
280 VecRestoreArray(b_local, &array);
281 VecGhostRestoreLocalForm(b, &b_local);
293template <std::
floating_po
int T>
297 VecGhostGetLocalForm(b, &b_local);
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);
329template <std::
floating_po
int T>
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>>>&
338 const std::vector<Vec>& x0, PetscScalar scale)
341 VecGhostGetLocalForm(b, &b_local);
343 VecGetSize(b_local, &n);
344 PetscScalar* array =
nullptr;
345 VecGetArray(b_local, &array);
346 std::span<PetscScalar> _b(array, n);
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)
358 VecGhostGetLocalForm(x0[i], &x0_local[i]);
360 VecGetSize(x0_local[i], &n);
361 VecGetArrayRead(x0_local[i], &x0_array[i]);
362 x0_ref.emplace_back(x0_array[i], n);
365 std::vector x0_tmp(x0_ref.begin(), x0_ref.end());
368 for (std::size_t i = 0; i < x0_local.size(); ++i)
370 VecRestoreArrayRead(x0_local[i], &x0_array[i]);
371 VecGhostRestoreLocalForm(x0[i], &x0_local[i]);
375 VecRestoreArray(b_local, &array);
376 VecGhostRestoreLocalForm(b, &b_local);
397template <std::
floating_po
int T>
404 const std::vector<Vec>& x0, PetscScalar scale)
407 VecGhostGetLocalForm(b, &b_local);
409 VecGetSize(b_local, &n);
410 PetscScalar* array =
nullptr;
411 VecGetArray(b_local, &array);
412 std::span<PetscScalar> _b(array, n);
415 fem::apply_lifting<PetscScalar>(_b, a, bcs1, {}, scale);
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)
424 VecGhostGetLocalForm(x0[i], &x0_local[i]);
426 VecGetSize(x0_local[i], &n);
427 VecGetArrayRead(x0_local[i], &x0_array[i]);
428 x0_ref.emplace_back(x0_array[i], n);
431 std::vector x0_tmp(x0_ref.begin(), x0_ref.end());
432 fem::apply_lifting<PetscScalar>(_b, a, bcs1, x0_tmp, scale);
434 for (std::size_t i = 0; i < x0_local.size(); ++i)
436 VecRestoreArrayRead(x0_local[i], &x0_array[i]);
437 VecGhostRestoreLocalForm(x0[i], &x0_local[i]);
441 VecRestoreArray(b_local, &array);
442 VecGhostRestoreLocalForm(b, &b_local);
455template <std::
floating_po
int T>
459 const Vec x0, PetscScalar scale = 1)
462 VecGetLocalSize(b, &n);
463 PetscScalar* array =
nullptr;
464 VecGetArray(b, &array);
465 std::span<PetscScalar> _b(array, n);
469 VecGhostGetLocalForm(x0, &x0_local);
471 VecGetSize(x0_local, &n);
472 const PetscScalar* array =
nullptr;
473 VecGetArrayRead(x0_local, &array);
474 std::span<const PetscScalar> _x0(array, n);
476 VecRestoreArrayRead(x0_local, &array);
477 VecGhostRestoreLocalForm(x0, &x0_local);
482 VecRestoreArray(b, &array);
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
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