Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.1.0/v0.9.0/cpp
DOLFINx  0.1.0
DOLFINx C++ interface
Form.h
1 // Copyright (C) 2019-2020 Garth N. Wells and Chris Richardson
2 //
3 // This file is part of DOLFINx (https://www.fenicsproject.org)
4 //
5 // SPDX-License-Identifier: LGPL-3.0-or-later
6 
7 #pragma once
8 
9 #include <dolfinx/fem/FunctionSpace.h>
10 #include <dolfinx/mesh/Mesh.h>
11 #include <dolfinx/mesh/MeshTags.h>
12 #include <functional>
13 #include <memory>
14 #include <string>
15 #include <vector>
16 
17 namespace dolfinx::fem
18 {
19 
20 template <typename T>
21 class Constant;
22 template <typename T>
23 class Function;
24 
26 enum class IntegralType : std::int8_t
27 {
28  cell = 0,
29  exterior_facet = 1,
30  interior_facet = 2,
31  vertex = 3
32 };
33 
57 
58 template <typename T>
59 class Form
60 {
61 public:
75  Form(const std::vector<std::shared_ptr<const fem::FunctionSpace>>&
77  const std::map<
79  std::pair<
80  std::vector<std::pair<
81  int, std::function<void(
82  T*, const T*, const T*, const double*, const int*,
83  const std::uint8_t*, const std::uint32_t)>>>,
84  const mesh::MeshTags<int>*>>& integrals,
85  const std::vector<std::shared_ptr<const fem::Function<T>>>& coefficients,
86  const std::vector<std::shared_ptr<const fem::Constant<T>>>& constants,
88  const std::shared_ptr<const mesh::Mesh>& mesh = nullptr)
89  : _function_spaces(function_spaces), _coefficients(coefficients),
90  _constants(constants), _mesh(mesh),
91  _needs_permutation_data(needs_permutation_data)
92  {
93  // Extract _mesh from fem::FunctionSpace, and check they are the same
94  if (!_mesh and !function_spaces.empty())
95  _mesh = function_spaces[0]->mesh();
96  for (const auto& V : function_spaces)
97  if (_mesh != V->mesh())
98  throw std::runtime_error("Incompatible mesh");
99  if (!_mesh)
100  throw std::runtime_error("No mesh could be associated with the Form.");
101 
102  // Store kernels, looping over integrals by domain type (dimension)
103  for (auto& integral_type : integrals)
104  {
105  // Add key to map
106  const IntegralType type = integral_type.first;
107  auto it = _integrals.emplace(
108  type, std::map<int, std::pair<kern, std::vector<std::int32_t>>>());
109 
110  // Loop over integrals kernels
111  for (auto& integral : integral_type.second.first)
112  it.first->second.insert({integral.first, {integral.second, {}}});
113 
114  // FIXME: do this neatly via a static function
115  // Set domains for integral type
116  if (integral_type.second.second)
117  {
118  assert(_mesh == integral_type.second.second->mesh());
119  set_domains(type, *integral_type.second.second);
120  }
121  }
122 
123  // FIXME: do this neatly via a static function
124  // Set markers for default integrals
125  set_default_domains(*_mesh);
126  }
127 
129  Form(const Form& form) = delete;
130 
132  Form(Form&& form) = default;
133 
135  virtual ~Form() = default;
136 
140  int rank() const { return _function_spaces.size(); }
141 
144  std::shared_ptr<const mesh::Mesh> mesh() const { return _mesh; }
145 
148  const std::vector<std::shared_ptr<const fem::FunctionSpace>>&
150  {
151  return _function_spaces;
152  }
153 
159  const std::function<void(T*, const T*, const T*, const double*, const int*,
160  const std::uint8_t*, const std::uint32_t)>&
161  kernel(IntegralType type, int i) const
162  {
163  auto it0 = _integrals.find(type);
164  if (it0 == _integrals.end())
165  throw std::runtime_error("No kernels for requested type.");
166  auto it1 = it0->second.find(i);
167  if (it1 == it0->second.end())
168  throw std::runtime_error("No kernel for requested domain index.");
169 
170  return it1->second.first;
171  }
172 
175  std::set<IntegralType> integral_types() const
176  {
177  std::set<IntegralType> set;
178  for (auto& type : _integrals)
179  set.insert(type.first);
180  return set;
181  }
182 
186  int num_integrals(IntegralType type) const
187  {
188  if (auto it = _integrals.find(type); it == _integrals.end())
189  return 0;
190  else
191  return it->second.size();
192  }
193 
200  std::vector<int> integral_ids(IntegralType type) const
201  {
202  std::vector<int> ids;
203  if (auto it = _integrals.find(type); it != _integrals.end())
204  {
205  for (auto& kernel : it->second)
206  ids.push_back(kernel.first);
207  }
208  return ids;
209  }
210 
217  const std::vector<std::int32_t>& domains(IntegralType type, int i) const
218  {
219  auto it0 = _integrals.find(type);
220  if (it0 == _integrals.end())
221  throw std::runtime_error("No mesh entities for requested type.");
222  auto it1 = it0->second.find(i);
223  if (it1 == it0->second.end())
224  throw std::runtime_error("No mesh entities for requested domain index.");
225  return it1->second.second;
226  }
227 
229  const std::vector<std::shared_ptr<const fem::Function<T>>>
230  coefficients() const
231  {
232  return _coefficients;
233  }
234 
238  bool needs_permutation_data() const { return _needs_permutation_data; }
239 
243  std::vector<int> coefficient_offsets() const
244  {
245  std::vector<int> n{0};
246  for (const auto& c : _coefficients)
247  {
248  if (!c)
249  throw std::runtime_error("Not all form coefficients have been set.");
250  n.push_back(n.back() + c->function_space()->element()->space_dimension());
251  }
252  return n;
253  }
254 
256  const std::vector<std::shared_ptr<const fem::Constant<T>>>& constants() const
257  {
258  return _constants;
259  }
260 
262  using scalar_type = T;
263 
264 private:
270  void set_domains(IntegralType type, const mesh::MeshTags<int>& marker)
271  {
272  auto it0 = _integrals.find(type);
273  assert(it0 != _integrals.end());
274 
275  std::shared_ptr<const mesh::Mesh> mesh = marker.mesh();
276  const mesh::Topology& topology = mesh->topology();
277  const int tdim = topology.dim();
278  int dim = tdim;
279  if (type == IntegralType::exterior_facet
280  or type == IntegralType::interior_facet)
281  {
282  dim = tdim - 1;
283  mesh->topology_mutable().create_connectivity(dim, tdim);
284  }
285  else if (type == IntegralType::vertex)
286  dim = 0;
287 
288  if (dim != marker.dim())
289  {
290  throw std::runtime_error("Invalid MeshTags dimension:"
291  + std::to_string(marker.dim()));
292  }
293 
294  // Get all integrals for considered entity type
295  std::map<int, std::pair<kern, std::vector<std::int32_t>>>& integrals
296  = it0->second;
297 
298  // Get mesh tag data
299  const std::vector<int>& values = marker.values();
300  const std::vector<std::int32_t>& tagged_entities = marker.indices();
301  assert(topology.index_map(dim));
302  const auto entity_end
303  = std::lower_bound(tagged_entities.begin(), tagged_entities.end(),
304  topology.index_map(dim)->size_local());
305 
306  if (dim == tdim - 1)
307  {
308  auto f_to_c = topology.connectivity(tdim - 1, tdim);
309  assert(f_to_c);
310  if (type == IntegralType::exterior_facet)
311  {
312  // Only need to consider shared facets when there are no ghost
313  // cells
314  assert(topology.index_map(tdim));
315  std::set<std::int32_t> fwd_shared;
316  if (topology.index_map(tdim)->num_ghosts() == 0)
317  {
318  fwd_shared.insert(
319  topology.index_map(tdim - 1)->shared_indices().array().begin(),
320  topology.index_map(tdim - 1)->shared_indices().array().end());
321  }
322 
323  for (auto f = tagged_entities.begin(); f != entity_end; ++f)
324  {
325  // All "owned" facets connected to one cell, that are not
326  // shared, should be external
327  if (f_to_c->num_links(*f) == 1
328  and fwd_shared.find(*f) == fwd_shared.end())
329  {
330  const std::size_t i = std::distance(tagged_entities.cbegin(), f);
331  if (auto it = integrals.find(values[i]); it != integrals.end())
332  it->second.second.push_back(*f);
333  }
334  }
335  }
336  else if (type == IntegralType::interior_facet)
337  {
338  for (auto f = tagged_entities.begin(); f != entity_end; ++f)
339  {
340  if (f_to_c->num_links(*f) == 2)
341  {
342  const std::size_t i = std::distance(tagged_entities.cbegin(), f);
343  if (auto it = integrals.find(values[i]); it != integrals.end())
344  it->second.second.push_back(*f);
345  }
346  }
347  }
348  }
349  else
350  {
351  // For cell and vertex integrals use all markers (but not on ghost
352  // entities)
353  for (auto e = tagged_entities.begin(); e != entity_end; ++e)
354  {
355  const std::size_t i = std::distance(tagged_entities.cbegin(), e);
356  if (auto it = integrals.find(values[i]); it != integrals.end())
357  it->second.second.push_back(*e);
358  }
359  }
360  }
361 
367  void set_default_domains(const mesh::Mesh& mesh)
368  {
369  const mesh::Topology& topology = mesh.topology();
370  const int tdim = topology.dim();
371 
372  // Cells. If there is a default integral, define it on all owned cells
373  if (auto kernels = _integrals.find(IntegralType::cell);
374  kernels != _integrals.end())
375  {
376  if (auto it = kernels->second.find(-1); it != kernels->second.end())
377  {
378  std::vector<std::int32_t>& active_entities = it->second.second;
379  const int num_cells = topology.index_map(tdim)->size_local();
380  active_entities.resize(num_cells);
381  std::iota(active_entities.begin(), active_entities.end(), 0);
382  }
383  }
384 
385  // Exterior facets. If there is a default integral, define it only
386  // on owned surface facets.
387  if (auto kernels = _integrals.find(IntegralType::exterior_facet);
388  kernels != _integrals.end())
389  {
390  if (auto it = kernels->second.find(-1); it != kernels->second.end())
391  {
392  std::vector<std::int32_t>& active_entities = it->second.second;
393  active_entities.clear();
394 
395  // Get number of facets owned by this process
396  mesh.topology_mutable().create_connectivity(tdim - 1, tdim);
397  auto f_to_c = topology.connectivity(tdim - 1, tdim);
398  assert(topology.index_map(tdim - 1));
399  std::set<std::int32_t> fwd_shared_facets;
400 
401  // Only need to consider shared facets when there are no ghost cells
402  if (topology.index_map(tdim)->num_ghosts() == 0)
403  {
404  fwd_shared_facets.insert(
405  topology.index_map(tdim - 1)->shared_indices().array().begin(),
406  topology.index_map(tdim - 1)->shared_indices().array().end());
407  }
408 
409  const int num_facets = topology.index_map(tdim - 1)->size_local();
410  for (int f = 0; f < num_facets; ++f)
411  {
412  if (f_to_c->num_links(f) == 1
413  and fwd_shared_facets.find(f) == fwd_shared_facets.end())
414  {
415  active_entities.push_back(f);
416  }
417  }
418  }
419  }
420 
421  // Interior facets. If there is a default integral, define it only on
422  // owned interior facets.
423  if (auto kernels = _integrals.find(IntegralType::interior_facet);
424  kernels != _integrals.end())
425  {
426  if (auto it = kernels->second.find(-1); it != kernels->second.end())
427  {
428  std::vector<std::int32_t>& active_entities = it->second.second;
429 
430  // Get number of facets owned by this process
431  mesh.topology_mutable().create_connectivity(tdim - 1, tdim);
432  assert(topology.index_map(tdim - 1));
433  const int num_facets = topology.index_map(tdim - 1)->size_local();
434  auto f_to_c = topology.connectivity(tdim - 1, tdim);
435  active_entities.clear();
436  active_entities.reserve(num_facets);
437  for (int f = 0; f < num_facets; ++f)
438  {
439  if (f_to_c->num_links(f) == 2)
440  active_entities.push_back(f);
441  }
442  }
443  }
444  }
445 
446  // Function spaces (one for each argument)
447  std::vector<std::shared_ptr<const fem::FunctionSpace>> _function_spaces;
448 
449  // Form coefficients
450  std::vector<std::shared_ptr<const fem::Function<T>>> _coefficients;
451 
452  // Constants associated with the Form
453  std::vector<std::shared_ptr<const fem::Constant<T>>> _constants;
454 
455  // The mesh
456  std::shared_ptr<const mesh::Mesh> _mesh;
457 
458  using kern
459  = std::function<void(T*, const T*, const T*, const double*, const int*,
460  const std::uint8_t*, const std::uint32_t)>;
461  std::map<IntegralType,
462  std::map<int, std::pair<kern, std::vector<std::int32_t>>>>
463  _integrals;
464 
465  // True if permutation data needs to be passed into these integrals
466  bool _needs_permutation_data;
467 };
468 } // namespace dolfinx::fem
dolfinx::fem::Form::Form
Form(const std::vector< std::shared_ptr< const fem::FunctionSpace >> &function_spaces, const std::map< IntegralType, std::pair< std::vector< std::pair< int, std::function< void(T *, const T *, const T *, const double *, const int *, const std::uint8_t *, const std::uint32_t)>>>, const mesh::MeshTags< int > * >> &integrals, const std::vector< std::shared_ptr< const fem::Function< T >>> &coefficients, const std::vector< std::shared_ptr< const fem::Constant< T >>> &constants, bool needs_permutation_data, const std::shared_ptr< const mesh::Mesh > &mesh=nullptr)
Create form.
Definition: Form.h:75
dolfinx::fem::Form::coefficient_offsets
std::vector< int > coefficient_offsets() const
Offset for each coefficient expansion array on a cell. Used to pack data for multiple coefficients in...
Definition: Form.h:243
dolfinx::fem::IntegralType
IntegralType
Type of integral.
Definition: Form.h:26
dolfinx::fem::Function
This class represents a function in a finite element function space , given by.
Definition: Form.h:23
dolfinx::fem::Form::constants
const std::vector< std::shared_ptr< const fem::Constant< T > > > & constants() const
Access constants.
Definition: Form.h:256
dolfinx::fem::Form::integral_types
std::set< IntegralType > integral_types() const
Get types of integrals in the form.
Definition: Form.h:175
dolfinx::fem::Form::rank
int rank() const
Rank of the form (bilinear form = 2, linear form = 1, functional = 0, etc)
Definition: Form.h:140
dolfinx::fem::Form::num_integrals
int num_integrals(IntegralType type) const
Number of integrals of given type.
Definition: Form.h:186
dolfinx::mesh::Topology::dim
int dim() const
Return topological dimension.
Definition: Topology.cpp:162
dolfinx::mesh::MeshTags::dim
int dim() const
Return topological dimension of tagged entities.
Definition: MeshTags.h:105
dolfinx::fem::Form::~Form
virtual ~Form()=default
Destructor.
dolfinx::fem::Form
Class for variational forms.
Definition: assembler.h:22
dolfinx::mesh::MeshTags::mesh
std::shared_ptr< const Mesh > mesh() const
Return mesh.
Definition: MeshTags.h:108
dolfinx::mesh::MeshTags
A MeshTags are used to associate mesh entities with values. The entity index (local to process) ident...
Definition: XDMFFile.h:40
dolfinx::mesh::Topology::index_map
std::shared_ptr< const common::IndexMap > index_map(int dim) const
Get the IndexMap that described the parallel distribution of the mesh entities.
Definition: Topology.cpp:171
dolfinx::fem::Form::coefficients
const std::vector< std::shared_ptr< const fem::Function< T > > > coefficients() const
Access coefficients.
Definition: Form.h:230
dolfinx::mesh::Mesh
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition: Mesh.h:55
dolfinx::mesh::Topology::connectivity
std::shared_ptr< const graph::AdjacencyList< std::int32_t > > connectivity(int d0, int d1) const
Return connectivity from entities of dimension d0 to entities of dimension d1.
Definition: Topology.cpp:263
dolfinx::fem::Form::kernel
const std::function< void(T *, const T *, const T *, const double *, const int *, const std::uint8_t *, const std::uint32_t)> & kernel(IntegralType type, int i) const
Get the function for 'kernel' for integral i of given type.
Definition: Form.h:161
dolfinx::fem::Form::integral_ids
std::vector< int > integral_ids(IntegralType type) const
Get the IDs for integrals (kernels) for given integral type. The IDs correspond to the domain IDs whi...
Definition: Form.h:200
dolfinx::mesh::MeshTags::indices
const std::vector< std::int32_t > & indices() const
Indices of tagged mesh entities (local-to-process). The indices are sorted.
Definition: MeshTags.h:99
dolfinx::fem::Form::mesh
std::shared_ptr< const mesh::Mesh > mesh() const
Extract common mesh for the form.
Definition: Form.h:144
dolfinx::mesh::MeshTags::values
const std::vector< T > & values() const
Values attached to mesh entities.
Definition: MeshTags.h:102
dolfinx::fem::Form::scalar_type
T scalar_type
Scalar type (T).
Definition: Form.h:262
dolfinx::fem::Form::needs_permutation_data
bool needs_permutation_data() const
Get bool indicating whether permutation data needs to be passed into these integrals.
Definition: Form.h:238
dolfinx::fem::Form::domains
const std::vector< std::int32_t > & domains(IntegralType type, int i) const
Get the list of mesh entity indices for the ith integral (kernel) for the given domain type,...
Definition: Form.h:217
dolfinx::fem
Finite element method functionality.
Definition: assemble_matrix_impl.h:22
dolfinx::mesh::Topology
Topology stores the topology of a mesh, consisting of mesh entities and connectivity (incidence relat...
Definition: Topology.h:56
dolfinx::fem::Constant
A constant value which can be attached to a Form. Constants may be scalar (rank 0),...
Definition: Constant.h:18
dolfinx::fem::Form::function_spaces
const std::vector< std::shared_ptr< const fem::FunctionSpace > > & function_spaces() const
Return function spaces for all arguments.
Definition: Form.h:149