9 #include "FunctionSpace.h"
12 #include <dolfinx/common/IndexMap.h>
13 #include <dolfinx/mesh/Mesh.h>
14 #include <dolfinx/mesh/MeshTags.h>
64 template <
typename X,
typename =
void>
65 struct scalar_value_type
70 struct scalar_value_type<X, std::void_t<typename X::value_type>>
72 typedef typename X::value_type value_type;
74 using scalar_value_type_t =
typename scalar_value_type<T>::value_type;
97 std::vector<std::pair<
98 int, std::function<
void(T*,
const T*,
const T*,
99 const scalar_value_type_t*,
100 const int*,
const std::uint8_t*)>>>,
105 const std::shared_ptr<const mesh::Mesh>&
mesh =
nullptr)
115 if (_mesh != V->mesh())
116 throw std::runtime_error(
"Incompatible mesh");
119 throw std::runtime_error(
"No mesh could be associated with the Form.");
122 for (
auto& integral_type : integrals)
129 for (
auto& integral : integral_type.second.first)
130 _cell_integrals.insert({integral.first, {integral.second, {}}});
133 for (
auto& integral : integral_type.second.first)
135 _exterior_facet_integrals.insert(
136 {integral.first, {integral.second, {}}});
140 for (
auto& integral : integral_type.second.first)
142 _interior_facet_integrals.insert(
143 {integral.first, {integral.second, {}}});
148 if (integral_type.second.second)
150 assert(_mesh == integral_type.second.second->mesh());
151 set_domains(type, *integral_type.second.second);
157 set_default_domains(*_mesh);
172 int rank()
const {
return _function_spaces.size(); }
176 std::shared_ptr<const mesh::Mesh>
mesh()
const {
return _mesh; }
180 const std::vector<std::shared_ptr<const FunctionSpace>>&
183 return _function_spaces;
191 const std::function<void(T*,
const T*,
const T*,
const scalar_value_type_t*,
192 const int*,
const std::uint8_t*)>&
198 return get_kernel_from_integrals(_cell_integrals, i);
200 return get_kernel_from_integrals(_exterior_facet_integrals, i);
202 return get_kernel_from_integrals(_interior_facet_integrals, i);
204 throw std::runtime_error(
205 "Cannot access kernel. Integral type not supported.");
213 std::set<IntegralType> set;
214 if (!_cell_integrals.empty())
216 if (!_exterior_facet_integrals.empty())
218 if (!_interior_facet_integrals.empty())
232 return _cell_integrals.size();
234 return _exterior_facet_integrals.size();
236 return _interior_facet_integrals.size();
238 throw std::runtime_error(
"Integral type not supported.");
250 std::vector<int> ids;
254 std::transform(_cell_integrals.cbegin(), _cell_integrals.cend(),
255 std::back_inserter(ids),
256 [](
auto& integral) { return integral.first; });
259 std::transform(_exterior_facet_integrals.cbegin(),
260 _exterior_facet_integrals.cend(), std::back_inserter(ids),
261 [](
auto& integral) { return integral.first; });
264 std::transform(_interior_facet_integrals.cbegin(),
265 _interior_facet_integrals.cend(), std::back_inserter(ids),
266 [](
auto& integral) { return integral.first; });
269 throw std::runtime_error(
270 "Cannot return IDs. Integral type not supported.");
282 auto it = _cell_integrals.find(i);
283 if (it == _cell_integrals.end())
284 throw std::runtime_error(
"No mesh entities for requested domain index.");
285 return it->second.second;
295 auto it = _exterior_facet_integrals.find(i);
296 if (it == _exterior_facet_integrals.end())
297 throw std::runtime_error(
"No mesh entities for requested domain index.");
298 return it->second.second;
310 auto it = _interior_facet_integrals.find(i);
311 if (it == _interior_facet_integrals.end())
312 throw std::runtime_error(
"No mesh entities for requested domain index.");
313 return it->second.second;
317 const std::vector<std::shared_ptr<const Function<T>>>&
coefficients()
const
319 return _coefficients;
332 std::vector<int> n = {0};
333 for (
const auto& c : _coefficients)
336 throw std::runtime_error(
"Not all form coefficients have been set.");
337 n.push_back(n.back() + c->function_space()->element()->space_dimension());
343 const std::vector<std::shared_ptr<const Constant<T>>>&
constants()
const
353 = std::function<void(T*,
const T*,
const T*,
const scalar_value_type_t*,
354 const int*,
const std::uint8_t*)>;
361 template <
typename U>
362 const std::function<void(T*,
const T*,
const T*,
const scalar_value_type_t*,
363 const int*,
const std::uint8_t*)>&
364 get_kernel_from_integrals(
const U& integrals,
int i)
const
366 auto it = integrals.find(i);
367 if (it == integrals.end())
368 throw std::runtime_error(
"No kernel for requested domain index.");
369 return it->second.first;
378 template <
int num_cells>
379 static std::array<std::array<std::int32_t, 2>, num_cells>
380 get_cell_local_facet_pairs(
381 std::int32_t f,
const std::span<const std::int32_t>& cells,
385 assert(cells.size() == num_cells);
386 std::array<std::array<std::int32_t, 2>, num_cells> cell_local_facet_pairs;
387 for (
int c = 0; c < num_cells; ++c)
390 std::int32_t
cell = cells[c];
391 auto cell_facets = c_to_f.
links(cell);
392 auto facet_it = std::find(cell_facets.begin(), cell_facets.end(), f);
393 assert(facet_it != cell_facets.end());
394 int local_f = std::distance(cell_facets.begin(), facet_it);
395 cell_local_facet_pairs[c] = {
cell, local_f};
398 return cell_local_facet_pairs;
402 template <
typename iterator>
403 void set_cell_domains(
404 std::map<
int, std::pair<kern, std::vector<std::int32_t>>>& integrals,
405 const iterator& tagged_cells_begin,
const iterator& tagged_cells_end,
406 const std::vector<int>& tags)
409 for (
auto c = tagged_cells_begin; c != tagged_cells_end; ++c)
411 const std::size_t pos = std::distance(tagged_cells_begin, c);
412 if (
auto it = integrals.find(tags[pos]); it != integrals.end())
413 it->second.second.push_back(*c);
418 template <
typename iterator>
419 void set_exterior_facet_domains(
420 const mesh::Topology& topology,
421 std::map<
int, std::pair<kern, std::vector<std::int32_t>>>& integrals,
422 const iterator& tagged_facets_begin,
const iterator& tagged_facets_end,
423 const std::vector<int>& tags)
431 int tdim = topology.dim();
432 assert(topology.index_map(tdim));
433 assert(topology.index_map(tdim - 1));
434 const std::vector<std::int32_t> fwd_shared_facets
435 = topology.index_map(tdim)->overlapped()
436 ? std::vector<std::int32_t>()
437 : topology.index_map(tdim - 1)->shared_indices();
439 auto f_to_c = topology.connectivity(tdim - 1, tdim);
441 auto c_to_f = topology.connectivity(tdim, tdim - 1);
443 for (
auto f = tagged_facets_begin; f != tagged_facets_end; ++f)
449 if (f_to_c->num_links(*f) == 1)
451 if (!std::binary_search(fwd_shared_facets.begin(),
452 fwd_shared_facets.end(), *f))
454 const std::size_t pos = std::distance(tagged_facets_begin, f);
455 if (
auto it = integrals.find(tags[pos]); it != integrals.end())
458 const std::array<std::int32_t, 2> pair
459 = get_cell_local_facet_pairs<1>(*f, f_to_c->links(*f),
461 it->second.second.insert(it->second.second.end(), pair.cbegin(),
470 template <
typename iterator>
471 static void set_interior_facet_domains(
472 const mesh::Topology& topology,
473 std::map<
int, std::pair<kern, std::vector<std::int32_t>>>& integrals,
474 const iterator& tagged_facets_begin,
const iterator& tagged_facets_end,
475 const std::vector<int>& tags)
477 int tdim = topology.dim();
478 auto f_to_c = topology.connectivity(tdim - 1, tdim);
480 auto c_to_f = topology.connectivity(tdim, tdim - 1);
482 for (
auto f = tagged_facets_begin; f != tagged_facets_end; ++f)
484 if (f_to_c->num_links(*f) == 2)
486 const std::size_t pos = std::distance(tagged_facets_begin, f);
487 if (
auto it = integrals.find(tags[pos]); it != integrals.end())
489 const std::array<std::array<std::int32_t, 2>, 2> pairs
490 = get_cell_local_facet_pairs<2>(*f, f_to_c->links(*f), *c_to_f);
491 it->second.second.insert(it->second.second.end(), pairs[0].cbegin(),
493 it->second.second.insert(it->second.second.end(), pairs[1].cbegin(),
506 void set_domains(
IntegralType type,
const mesh::MeshTags<int>& marker)
508 std::shared_ptr<const mesh::Mesh>
mesh = marker.mesh();
509 const mesh::Topology& topology =
mesh->topology();
510 const int tdim = topology.dim();
512 if (dim != marker.dim())
514 throw std::runtime_error(
"Invalid MeshTags dimension: "
515 + std::to_string(marker.dim()));
519 const std::vector<int>& tags = marker.values();
520 const std::vector<std::int32_t>& tagged_entities = marker.indices();
521 assert(topology.index_map(dim));
522 const auto entity_end
523 = std::lower_bound(tagged_entities.begin(), tagged_entities.end(),
524 topology.index_map(dim)->size_local());
528 set_cell_domains(_cell_integrals, tagged_entities.cbegin(), entity_end,
532 mesh->topology_mutable().create_connectivity(dim, tdim);
533 mesh->topology_mutable().create_connectivity(tdim, dim);
537 set_exterior_facet_domains(topology, _exterior_facet_integrals,
538 tagged_entities.cbegin(), entity_end, tags);
541 set_interior_facet_domains(topology, _interior_facet_integrals,
542 tagged_entities.cbegin(), entity_end, tags);
545 throw std::runtime_error(
546 "Cannot set domains. Integral type not supported.");
556 void set_default_domains(
const mesh::Mesh&
mesh)
558 const mesh::Topology& topology =
mesh.topology();
559 const int tdim = topology.dim();
563 for (
auto& [domain_id, kernel_cells] : _cell_integrals)
567 std::vector<std::int32_t>&
cells = kernel_cells.second;
568 const int num_cells = topology.index_map(tdim)->size_local();
569 cells.resize(num_cells);
577 if (!_exterior_facet_integrals.empty())
579 mesh.topology_mutable().create_connectivity(tdim - 1, tdim);
580 mesh.topology_mutable().create_connectivity(tdim, tdim - 1);
582 const std::vector<std::int32_t> boundary_facets
583 = _exterior_facet_integrals.empty()
584 ? std::vector<std::int32_t>()
586 for (
auto& [domain_id, kernel_facets] : _exterior_facet_integrals)
590 std::vector<std::int32_t>& facets = kernel_facets.second;
593 auto f_to_c = topology.connectivity(tdim - 1, tdim);
595 auto c_to_f = topology.connectivity(tdim, tdim - 1);
597 for (std::int32_t f : boundary_facets)
600 std::array<std::int32_t, 2> pair
601 = get_cell_local_facet_pairs<1>(f, f_to_c->links(f), *c_to_f)[0];
602 facets.insert(facets.end(), pair.cbegin(), pair.cend());
609 for (
auto& [domain_id, kernel_facets] : _interior_facet_integrals)
613 std::vector<std::int32_t>& facets = kernel_facets.second;
616 mesh.topology_mutable().create_connectivity(tdim - 1, tdim);
617 auto f_to_c = topology.connectivity(tdim - 1, tdim);
619 mesh.topology_mutable().create_connectivity(tdim, tdim - 1);
620 auto c_to_f =
mesh.topology().connectivity(tdim, tdim - 1);
624 assert(topology.index_map(tdim - 1));
625 const int num_facets = topology.index_map(tdim - 1)->size_local();
626 facets.reserve(num_facets);
627 for (
int f = 0; f < num_facets; ++f)
629 if (f_to_c->num_links(f) == 2)
631 const std::array<std::array<std::int32_t, 2>, 2> pairs
632 = get_cell_local_facet_pairs<2>(f, f_to_c->links(f), *c_to_f);
633 facets.insert(facets.end(), pairs[0].cbegin(), pairs[0].cend());
634 facets.insert(facets.end(), pairs[1].cbegin(), pairs[1].cend());
642 std::vector<std::shared_ptr<const FunctionSpace>> _function_spaces;
645 std::vector<std::shared_ptr<const Function<T>>> _coefficients;
648 std::vector<std::shared_ptr<const Constant<T>>> _constants;
651 std::shared_ptr<const mesh::Mesh> _mesh;
654 std::map<int, std::pair<kern, std::vector<std::int32_t>>> _cell_integrals;
657 std::map<int, std::pair<kern, std::vector<std::int32_t>>>
658 _exterior_facet_integrals;
661 std::map<int, std::pair<kern, std::vector<std::int32_t>>>
662 _interior_facet_integrals;
665 bool _needs_facet_permutations;
Constant value which can be attached to a Form. Constants may be scalar (rank 0), vector (rank 1),...
Definition: Constant.h:20
This class represents a function in a finite element function space , given by.
Definition: Function.h:45
This class provides a static adjacency list data structure. It is commonly used to store directed gra...
Definition: AdjacencyList.h:26
std::span< T > links(int node)
Get the links (edges) for given node.
Definition: AdjacencyList.h:111
void cells(la::SparsityPattern &pattern, const mesh::Topology &topology, const std::array< const std::reference_wrapper< const DofMap >, 2 > &dofmaps)
Iterate over cells and insert entries into sparsity pattern.
Definition: sparsitybuild.cpp:18
Finite element method functionality.
Definition: assemble_matrix_impl.h:25
IntegralType
Type of integral.
Definition: Form.h:31
@ interior_facet
Interior facet.
@ exterior_facet
Exterior facet.
std::vector< std::int32_t > exterior_facet_indices(const Topology &topology)
Compute the indices of all exterior facets that are owned by the caller.
Definition: utils.cpp:570