15#include <dolfinx/la/petsc.h>
32template <dolfinx::scalar T, std::
floating_po
int U>
44template <std::
floating_po
int T>
45Mat create_matrix(
const Form<PetscScalar, T>& a,
46 std::string type = std::string())
48 la::SparsityPattern pattern = fem::create_sparsity_pattern(a);
50 return la::petsc::create_matrix(a.mesh()->comm(), pattern, type);
62template <std::
floating_po
int T>
63Mat create_matrix_block(
64 const std::vector<std::vector<
const Form<PetscScalar, T>*>>& a,
65 std::string type = std::string())
68 std::array<std::vector<std::shared_ptr<const FunctionSpace<T>>>, 2> V
70 std::array<std::vector<int>, 2> bs_dofs;
71 for (std::size_t i = 0; i < 2; ++i)
74 bs_dofs[i].push_back(_V->dofmap()->bs());
78 std::shared_ptr<const mesh::Mesh<T>> mesh;
79 std::vector<std::vector<std::unique_ptr<la::SparsityPattern>>> patterns(
81 for (std::size_t row = 0; row < V[0].size(); ++row)
83 for (std::size_t col = 0; col < V[1].size(); ++col)
85 if (
const Form<PetscScalar, T>* form = a[row][col]; form)
87 patterns[row].push_back(std::make_unique<la::SparsityPattern>(
93 patterns[row].push_back(
nullptr);
98 throw std::runtime_error(
"Could not find a Mesh.");
101 std::array<std::vector<std::pair<
102 std::reference_wrapper<const common::IndexMap>,
int>>,
105 for (std::size_t d = 0; d < 2; ++d)
107 for (
auto& space : V[d])
109 maps[d].emplace_back(*space->dofmap()->index_map,
110 space->dofmap()->index_map_bs());
115 std::vector<std::vector<const la::SparsityPattern*>> p(V[0].
size());
116 for (std::size_t row = 0; row < V[0].size(); ++row)
117 for (std::size_t col = 0; col < V[1].size(); ++col)
118 p[row].push_back(patterns[row][col].get());
120 la::SparsityPattern pattern(mesh->comm(), p, maps, bs_dofs);
127 Mat A = la::petsc::create_matrix(mesh->comm(), pattern, type);
131 std::array<std::vector<PetscInt>, 2> _maps;
132 for (
int d = 0; d < 2; ++d)
142 std::pair<std::reference_wrapper<const common::IndexMap>,
int>>& map
144 std::vector<PetscInt>& _map = _maps[d];
147 const auto [rank_offset, local_offset, ghosts, _]
148 = common::stack_index_maps(map);
149 for (std::size_t f = 0; f < map.size(); ++f)
151 auto offset = local_offset[f];
152 const common::IndexMap& imap = map[f].first.get();
153 int bs = map[f].second;
154 for (std::int32_t i = 0; i < bs * imap.size_local(); ++i)
155 _map.push_back(i + rank_offset + offset);
156 for (std::int32_t i = 0; i < bs * imap.num_ghosts(); ++i)
157 _map.push_back(ghosts[f][i]);
162 ISLocalToGlobalMapping petsc_local_to_global0;
163 ISLocalToGlobalMappingCreate(MPI_COMM_SELF, 1, _maps[0].
size(),
164 _maps[0].data(), PETSC_COPY_VALUES,
165 &petsc_local_to_global0);
168 MatSetLocalToGlobalMapping(A, petsc_local_to_global0,
169 petsc_local_to_global0);
170 ISLocalToGlobalMappingDestroy(&petsc_local_to_global0);
175 ISLocalToGlobalMapping petsc_local_to_global1;
176 ISLocalToGlobalMappingCreate(MPI_COMM_SELF, 1, _maps[1].
size(),
177 _maps[1].data(), PETSC_COPY_VALUES,
178 &petsc_local_to_global1);
179 MatSetLocalToGlobalMapping(A, petsc_local_to_global0,
180 petsc_local_to_global1);
181 ISLocalToGlobalMappingDestroy(&petsc_local_to_global0);
182 ISLocalToGlobalMappingDestroy(&petsc_local_to_global1);
191template <std::
floating_po
int T>
192Mat create_matrix_nest(
193 const std::vector<std::vector<
const Form<PetscScalar, T>*>>& a,
194 const std::vector<std::vector<std::string>>& types)
199 std::vector<std::vector<std::string>> _types(
200 a.size(), std::vector<std::string>(a.front().size()));
206 int cols = a.front().size();
207 std::vector<Mat> mats(rows * cols,
nullptr);
208 std::shared_ptr<const mesh::Mesh<T>> mesh;
209 for (
int i = 0; i < rows; ++i)
211 for (
int j = 0; j < cols; ++j)
213 if (
const Form<PetscScalar, T>* form = a[i][j]; form)
215 mats[i * cols + j] = create_matrix(*form, _types[i][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]);
244Vec create_vector_block(
246 std::pair<std::reference_wrapper<const common::IndexMap>,
int>>& maps);
249Vec create_vector_nest(
251 std::pair<std::reference_wrapper<const common::IndexMap>,
int>>& maps);
267template <std::
floating_po
int T>
269 Vec b,
const Form<PetscScalar, T>& L,
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);
281 fem::assemble_vector(_b, L, constants, coeffs);
282 VecRestoreArray(b_local, &array);
283 VecGhostRestoreLocalForm(b, &b_local);
295template <std::
floating_po
int T>
296void assemble_vector(Vec b,
const Form<PetscScalar, T>& L)
299 VecGhostGetLocalForm(b, &b_local);
301 VecGetSize(b_local, &n);
302 PetscScalar* array =
nullptr;
303 VecGetArray(b_local, &array);
304 std::span<PetscScalar> _b(array, n);
305 fem::assemble_vector(_b, L);
306 VecRestoreArray(b_local, &array);
307 VecGhostRestoreLocalForm(b, &b_local);
331template <std::
floating_po
int T>
335 std::optional<std::reference_wrapper<
const Form<PetscScalar, T>>>>
337 const std::vector<std::span<const PetscScalar>>& constants,
338 const std::vector<std::map<std::pair<IntegralType, int>,
339 std::pair<std::span<const PetscScalar>,
int>>>&
342 std::vector<std::reference_wrapper<
const DirichletBC<PetscScalar, T>>>>&
344 const std::vector<Vec>& x0, PetscScalar alpha)
347 VecGhostGetLocalForm(b, &b_local);
349 VecGetSize(b_local, &n);
350 PetscScalar* array =
nullptr;
351 VecGetArray(b_local, &array);
352 std::span<PetscScalar> _b(array, n);
355 fem::apply_lifting(_b, a, constants, coeffs, bcs1, {}, alpha);
358 std::vector<std::span<const PetscScalar>> x0_ref;
359 std::vector<Vec> x0_local(a.size());
360 std::vector<const PetscScalar*> x0_array(a.size());
361 for (std::size_t i = 0; i < a.size(); ++i)
364 VecGhostGetLocalForm(x0[i], &x0_local[i]);
366 VecGetSize(x0_local[i], &n);
367 VecGetArrayRead(x0_local[i], &x0_array[i]);
368 x0_ref.emplace_back(x0_array[i], n);
371 std::vector x0_tmp(x0_ref.begin(), x0_ref.end());
372 fem::apply_lifting(_b, a, constants, coeffs, bcs1, x0_tmp, alpha);
374 for (std::size_t i = 0; i < x0_local.size(); ++i)
376 VecRestoreArrayRead(x0_local[i], &x0_array[i]);
377 VecGhostRestoreLocalForm(x0[i], &x0_local[i]);
381 VecRestoreArray(b_local, &array);
382 VecGhostRestoreLocalForm(b, &b_local);
403template <std::
floating_po
int T>
407 std::optional<std::reference_wrapper<
const Form<PetscScalar, double>>>>
409 const std::vector<std::vector<
410 std::reference_wrapper<
const DirichletBC<PetscScalar, double>>>>& bcs1,
411 const std::vector<Vec>& x0, PetscScalar alpha)
414 VecGhostGetLocalForm(b, &b_local);
416 VecGetSize(b_local, &n);
417 PetscScalar* array =
nullptr;
418 VecGetArray(b_local, &array);
419 std::span<PetscScalar> _b(array, n);
422 fem::apply_lifting<PetscScalar>(_b, a, bcs1, {}, alpha);
425 std::vector<std::span<const PetscScalar>> x0_ref;
426 std::vector<Vec> x0_local(a.size());
427 std::vector<const PetscScalar*> x0_array(a.size());
428 for (std::size_t i = 0; i < a.size(); ++i)
431 VecGhostGetLocalForm(x0[i], &x0_local[i]);
433 VecGetSize(x0_local[i], &n);
434 VecGetArrayRead(x0_local[i], &x0_array[i]);
435 x0_ref.emplace_back(x0_array[i], n);
438 std::vector x0_tmp(x0_ref.begin(), x0_ref.end());
439 fem::apply_lifting<PetscScalar>(_b, a, bcs1, x0_tmp, alpha);
441 for (std::size_t i = 0; i < x0_local.size(); ++i)
443 VecRestoreArrayRead(x0_local[i], &x0_array[i]);
444 VecGhostRestoreLocalForm(x0[i], &x0_local[i]);
448 VecRestoreArray(b_local, &array);
449 VecGhostRestoreLocalForm(b, &b_local);
462template <std::
floating_po
int T>
465 const std::vector<std::reference_wrapper<
const DirichletBC<PetscScalar, T>>> bcs,
466 const Vec x0, PetscScalar alpha = 1)
469 VecGetLocalSize(b, &n);
470 PetscScalar* array =
nullptr;
471 VecGetArray(b, &array);
472 std::span<PetscScalar> _b(array, n);
476 VecGhostGetLocalForm(x0, &x0_local);
478 VecGetSize(x0_local, &n);
479 const PetscScalar* array =
nullptr;
480 VecGetArrayRead(x0_local, &array);
481 std::span<const PetscScalar> _x0(array, n);
483 bc.set(_b, _x0, alpha);
484 VecRestoreArrayRead(x0_local, &array);
485 VecGhostRestoreLocalForm(x0, &x0_local);
490 bc->set(_b, std::nullopt, alpha);
492 VecRestoreArray(b, &array);
Functions supporting finite element method operations.
int size(MPI_Comm comm)
Definition MPI.cpp:72
Miscellaneous classes, functions and types.
Definition dolfinx_common.h:8
Finite element method functionality.
Definition assemble_matrix_impl.h:26
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:127
la::SparsityPattern create_sparsity_pattern(const Form< T, U > &a)
Create a sparsity pattern for a given form.
Definition utils.h:152