9 #include "FunctionSpace.h"
12 #include <dolfinx/mesh/Mesh.h>
13 #include <dolfinx/mesh/MeshTags.h>
81 const std::vector<std::shared_ptr<const fem::FunctionSpace>>&
86 std::vector<std::pair<
87 int, std::function<
void(T*,
const T*,
const T*,
const double*,
88 const int*,
const std::uint8_t*)>>>,
93 const std::shared_ptr<const mesh::Mesh>&
mesh =
nullptr)
103 if (_mesh != V->mesh())
104 throw std::runtime_error(
"Incompatible mesh");
107 throw std::runtime_error(
"No mesh could be associated with the Form.");
110 for (
auto& integral_type : integrals)
117 for (
auto& integral : integral_type.second.first)
118 _cell_integrals.insert({integral.first, {integral.second, {}}});
121 for (
auto& integral : integral_type.second.first)
123 _exterior_facet_integrals.insert(
124 {integral.first, {integral.second, {}}});
128 for (
auto& integral : integral_type.second.first)
130 _interior_facet_integrals.insert(
131 {integral.first, {integral.second, {}}});
136 if (integral_type.second.second)
138 assert(_mesh == integral_type.second.second->mesh());
139 set_domains(type, *integral_type.second.second);
145 set_default_domains(*_mesh);
160 int rank()
const {
return _function_spaces.size(); }
164 std::shared_ptr<const mesh::Mesh>
mesh()
const {
return _mesh; }
168 const std::vector<std::shared_ptr<const fem::FunctionSpace>>&
171 return _function_spaces;
179 const std::function<void(T*,
const T*,
const T*,
const double*,
const int*,
180 const std::uint8_t*)>&
186 return get_kernel_from_integrals(_cell_integrals, i);
188 return get_kernel_from_integrals(_exterior_facet_integrals, i);
190 return get_kernel_from_integrals(_interior_facet_integrals, i);
192 throw std::runtime_error(
193 "Cannot access kernel. Integral type not supported.");
201 std::set<IntegralType> set;
202 if (!_cell_integrals.empty())
204 if (!_exterior_facet_integrals.empty())
206 if (!_interior_facet_integrals.empty())
220 return _cell_integrals.size();
222 return _exterior_facet_integrals.size();
224 return _interior_facet_integrals.size();
226 throw std::runtime_error(
"Integral type not supported.");
238 std::vector<int> ids;
242 std::transform(_cell_integrals.cbegin(), _cell_integrals.cend(),
243 std::back_inserter(ids),
244 [](
auto& integral) { return integral.first; });
247 std::transform(_exterior_facet_integrals.cbegin(),
248 _exterior_facet_integrals.cend(), std::back_inserter(ids),
249 [](
auto& integral) { return integral.first; });
252 std::transform(_interior_facet_integrals.cbegin(),
253 _interior_facet_integrals.cend(), std::back_inserter(ids),
254 [](
auto& integral) { return integral.first; });
257 throw std::runtime_error(
258 "Cannot return IDs. Integral type not supported.");
270 auto it = _cell_integrals.find(i);
271 if (it == _cell_integrals.end())
272 throw std::runtime_error(
"No mesh entities for requested domain index.");
273 return it->second.second;
280 const std::vector<std::pair<std::int32_t, int>>&
283 auto it = _exterior_facet_integrals.find(i);
284 if (it == _exterior_facet_integrals.end())
285 throw std::runtime_error(
"No mesh entities for requested domain index.");
286 return it->second.second;
295 const std::vector<std::tuple<std::int32_t, int, std::int32_t, int>>&
298 auto it = _interior_facet_integrals.find(i);
299 if (it == _interior_facet_integrals.end())
300 throw std::runtime_error(
"No mesh entities for requested domain index.");
301 return it->second.second;
305 const std::vector<std::shared_ptr<const fem::Function<T>>>&
308 return _coefficients;
321 std::vector<int> n = {0};
322 for (
const auto& c : _coefficients)
325 throw std::runtime_error(
"Not all form coefficients have been set.");
326 n.push_back(n.back() + c->function_space()->element()->space_dimension());
332 const std::vector<std::shared_ptr<const fem::Constant<T>>>&
constants()
const
341 using kern = std::function<void(T*,
const T*,
const T*,
const double*,
342 const int*,
const std::uint8_t*)>;
349 template <
typename U>
350 const std::function<void(T*,
const T*,
const T*,
const double*,
const int*,
351 const std::uint8_t*)>&
352 get_kernel_from_integrals(
const U& integrals,
int i)
const
354 auto it = integrals.find(i);
355 if (it == integrals.end())
356 throw std::runtime_error(
"No kernel for requested domain index.");
357 return it->second.first;
366 template <
int num_cells>
367 static std::array<std::pair<std::int32_t, int>, num_cells>
368 get_cell_local_facet_pairs(
369 std::int32_t f,
const xtl::span<const std::int32_t>& cells,
373 assert(cells.size() == num_cells);
374 std::array<std::pair<std::int32_t, int>, num_cells> cell_local_facet_pairs;
375 for (
int c = 0; c < num_cells; ++c)
378 std::int32_t
cell = cells[c];
379 auto cell_facets = c_to_f.
links(cell);
380 auto facet_it = std::find(cell_facets.begin(), cell_facets.end(), f);
381 assert(facet_it != cell_facets.end());
382 int local_f = std::distance(cell_facets.begin(), facet_it);
383 cell_local_facet_pairs[c] = {
cell, local_f};
386 return cell_local_facet_pairs;
390 template <
typename iterator>
391 void set_cell_domains(
392 std::map<
int, std::pair<kern, std::vector<std::int32_t>>>& integrals,
393 const iterator& tagged_cells_begin,
const iterator& tagged_cells_end,
394 const std::vector<int>& tags)
397 for (
auto c = tagged_cells_begin; c != tagged_cells_end; ++c)
399 const std::size_t pos = std::distance(tagged_cells_begin, c);
400 if (
auto it = integrals.find(tags[pos]); it != integrals.end())
401 it->second.second.push_back(*c);
406 template <
typename iterator>
407 void set_exterior_facet_domains(
408 const mesh::Topology& topology,
409 std::map<
int, std::pair<kern, std::vector<std::pair<std::int32_t, int>>>>&
411 const iterator& tagged_facets_begin,
const iterator& tagged_facets_end,
412 const std::vector<int>& tags)
420 int tdim = topology.dim();
421 assert(topology.index_map(tdim));
422 std::set<std::int32_t> fwd_shared_facets;
423 if (topology.index_map(tdim)->num_ghosts() == 0)
425 const std::vector<std::int32_t>& fwd_indices
426 = topology.index_map(tdim - 1)->scatter_fwd_indices().array();
427 fwd_shared_facets.insert(fwd_indices.begin(), fwd_indices.end());
430 auto f_to_c = topology.connectivity(tdim - 1, tdim);
432 auto c_to_f = topology.connectivity(tdim, tdim - 1);
434 for (
auto f = tagged_facets_begin; f != tagged_facets_end; ++f)
440 if (f_to_c->num_links(*f) == 1
441 and fwd_shared_facets.find(*f) == fwd_shared_facets.end())
443 const std::size_t pos = std::distance(tagged_facets_begin, f);
444 if (
auto it = integrals.find(tags[pos]); it != integrals.end())
447 std::pair<std::int32_t, int> pair = get_cell_local_facet_pairs<1>(
448 *f, f_to_c->links(*f), *c_to_f)[0];
449 it->second.second.push_back(pair);
456 template <
typename iterator>
457 static void set_interior_facet_domains(
458 const mesh::Topology& topology,
460 std::pair<kern, std::vector<std::tuple<std::int32_t,
int,
461 std::int32_t,
int>>>>&
463 const iterator& tagged_facets_begin,
const iterator& tagged_facets_end,
464 const std::vector<int>& tags)
466 int tdim = topology.dim();
467 auto f_to_c = topology.connectivity(tdim - 1, tdim);
469 auto c_to_f = topology.connectivity(tdim, tdim - 1);
471 for (
auto f = tagged_facets_begin; f != tagged_facets_end; ++f)
473 if (f_to_c->num_links(*f) == 2)
475 const std::size_t pos = std::distance(tagged_facets_begin, f);
476 if (
auto it = integrals.find(tags[pos]); it != integrals.end())
478 std::array<std::pair<std::int32_t, int>, 2> pairs
479 = get_cell_local_facet_pairs<2>(*f, f_to_c->links(*f), *c_to_f);
480 it->second.second.emplace_back(pairs[0].first, pairs[0].second,
481 pairs[1].first, pairs[1].second);
493 void set_domains(
IntegralType type,
const mesh::MeshTags<int>& marker)
495 std::shared_ptr<const mesh::Mesh>
mesh = marker.mesh();
496 const mesh::Topology& topology =
mesh->topology();
497 const int tdim = topology.dim();
499 if (dim != marker.dim())
501 throw std::runtime_error(
"Invalid MeshTags dimension: "
502 + std::to_string(marker.dim()));
506 const std::vector<int>& tags = marker.values();
507 const std::vector<std::int32_t>& tagged_entities = marker.indices();
508 assert(topology.index_map(dim));
509 const auto entity_end
510 = std::lower_bound(tagged_entities.begin(), tagged_entities.end(),
511 topology.index_map(dim)->size_local());
515 set_cell_domains(_cell_integrals, tagged_entities.cbegin(), entity_end,
519 mesh->topology_mutable().create_connectivity(dim, tdim);
520 mesh->topology_mutable().create_connectivity(tdim, dim);
524 set_exterior_facet_domains(topology, _exterior_facet_integrals,
525 tagged_entities.cbegin(), entity_end, tags);
528 set_interior_facet_domains(topology, _interior_facet_integrals,
529 tagged_entities.cbegin(), entity_end, tags);
532 throw std::runtime_error(
533 "Cannot set domains. Integral type not supported.");
543 void set_default_domains(
const mesh::Mesh&
mesh)
545 const mesh::Topology& topology =
mesh.topology();
546 const int tdim = topology.dim();
550 for (
auto& [domain_id, kernel_cells] : _cell_integrals)
554 std::vector<std::int32_t>&
cells = kernel_cells.second;
555 const int num_cells = topology.index_map(tdim)->size_local();
556 cells.resize(num_cells);
563 for (
auto& [domain_id, kernel_facets] : _exterior_facet_integrals)
567 std::vector<std::pair<std::int32_t, int>>& facets
568 = kernel_facets.second;
571 mesh.topology_mutable().create_connectivity(tdim - 1, tdim);
572 auto f_to_c = topology.connectivity(tdim - 1, tdim);
575 mesh.topology_mutable().create_connectivity(tdim, tdim - 1);
576 auto c_to_f =
mesh.topology().connectivity(tdim, tdim - 1);
579 std::vector<std::int8_t> boundary_facet_markers =
581 for (std::size_t f = 0; f < boundary_facet_markers.size(); ++f)
583 if (boundary_facet_markers[f])
587 std::pair<std::int32_t, int> pair = get_cell_local_facet_pairs<1>(
588 f, f_to_c->links(f), *c_to_f)[0];
589 facets.push_back(pair);
597 for (
auto& [domain_id, kernel_facets] : _interior_facet_integrals)
601 std::vector<std::tuple<std::int32_t, int, std::int32_t, int>>& facets
602 = kernel_facets.second;
605 mesh.topology_mutable().create_connectivity(tdim - 1, tdim);
606 auto f_to_c = topology.connectivity(tdim - 1, tdim);
608 mesh.topology_mutable().create_connectivity(tdim, tdim - 1);
609 auto c_to_f =
mesh.topology().connectivity(tdim, tdim - 1);
613 assert(topology.index_map(tdim - 1));
614 const int num_facets = topology.index_map(tdim - 1)->size_local();
615 facets.reserve(num_facets);
616 for (
int f = 0; f < num_facets; ++f)
618 if (f_to_c->num_links(f) == 2)
620 std::array<std::pair<std::int32_t, int>, 2> pairs
621 = get_cell_local_facet_pairs<2>(f, f_to_c->links(f), *c_to_f);
622 facets.emplace_back(pairs[0].first, pairs[0].second, pairs[1].first,
631 std::vector<std::shared_ptr<const fem::FunctionSpace>> _function_spaces;
634 std::vector<std::shared_ptr<const fem::Function<T>>> _coefficients;
637 std::vector<std::shared_ptr<const fem::Constant<T>>> _constants;
640 std::shared_ptr<const mesh::Mesh> _mesh;
643 std::map<int, std::pair<kern, std::vector<std::int32_t>>> _cell_integrals;
646 std::map<int, std::pair<kern, std::vector<std::pair<std::int32_t, int>>>>
647 _exterior_facet_integrals;
650 std::map<int, std::pair<kern, std::vector<std::tuple<std::int32_t, int,
651 std::int32_t,
int>>>>
652 _interior_facet_integrals;
655 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:21
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:46
xtl::span< T > links(int node)
Get the links (edges) for given node.
Definition: AdjacencyList.h:131
void cells(la::SparsityPattern &pattern, const mesh::Topology &topology, const std::array< const std::reference_wrapper< const fem::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:24
IntegralType
Type of integral.
Definition: Form.h:30
@ interior_facet
Interior facet.
@ exterior_facet
Exterior facet.
std::vector< std::int8_t > compute_boundary_facets(const Topology &topology)
Compute marker for owned facets that are on the exterior of the domain, i.e. are connected to only on...
Definition: Topology.cpp:702