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
assemble_scalar_impl.h
1 // Copyright (C) 2019-2020 Garth N. Wells
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 "Form.h"
10 #include "utils.h"
11 #include <dolfinx/common/IndexMap.h>
12 #include <dolfinx/common/types.h>
13 #include <dolfinx/fem/Constant.h>
14 #include <dolfinx/fem/FunctionSpace.h>
15 #include <dolfinx/mesh/Geometry.h>
16 #include <dolfinx/mesh/Mesh.h>
17 #include <dolfinx/mesh/Topology.h>
18 #include <memory>
19 #include <vector>
20 
21 namespace dolfinx::fem::impl
22 {
23 
25 template <typename T>
26 T assemble_scalar(const fem::Form<T>& M);
27 
29 template <typename T>
30 T assemble_cells(
31  const mesh::Geometry& geometry,
32  const std::vector<std::int32_t>& active_cells,
33  const std::function<void(T*, const T*, const T*, const double*, const int*,
34  const std::uint8_t*, const std::uint32_t)>& fn,
35  const array2d<T>& coeffs, const std::vector<T>& constant_values,
36  const std::vector<std::uint32_t>& cell_info);
37 
39 template <typename T>
40 T assemble_exterior_facets(
41  const mesh::Mesh& mesh, const std::vector<std::int32_t>& active_cells,
42  const std::function<void(T*, const T*, const T*, const double*, const int*,
43  const std::uint8_t*, const std::uint32_t)>& fn,
44  const array2d<T>& coeffs, const std::vector<T>& constant_values,
45  const std::vector<std::uint32_t>& cell_info,
46  const std::vector<std::uint8_t>& perms);
47 
49 template <typename T>
50 T assemble_interior_facets(
51  const mesh::Mesh& mesh, const std::vector<std::int32_t>& active_cells,
52  const std::function<void(T*, const T*, const T*, const double*, const int*,
53  const std::uint8_t*, const std::uint32_t)>& fn,
54  const array2d<T>& coeffs, const std::vector<int>& offsets,
55  const std::vector<T>& constant_values,
56  const std::vector<std::uint32_t>& cell_info,
57  const std::vector<std::uint8_t>& perms);
58 
59 //-----------------------------------------------------------------------------
60 template <typename T>
61 T assemble_scalar(const fem::Form<T>& M)
62 {
63  std::shared_ptr<const mesh::Mesh> mesh = M.mesh();
64  assert(mesh);
65  const int tdim = mesh->topology().dim();
66  const std::int32_t num_cells
67  = mesh->topology().connectivity(tdim, 0)->num_nodes();
68 
69  // Prepare constants
70  const std::vector<std::shared_ptr<const fem::Constant<T>>>& constants
71  = M.constants();
72 
73  std::vector<T> constant_values;
74  for (auto const& constant : constants)
75  {
76  // Get underlying data array of this Constant
77  const std::vector<T>& array = constant->value;
78  constant_values.insert(constant_values.end(), array.begin(), array.end());
79  }
80 
81  // Prepare coefficients
82  const array2d<T> coeffs = pack_coefficients(M);
83 
84  const bool needs_permutation_data = M.needs_permutation_data();
85  if (needs_permutation_data)
86  mesh->topology_mutable().create_entity_permutations();
87  const std::vector<std::uint32_t>& cell_info
88  = needs_permutation_data ? mesh->topology().get_cell_permutation_info()
89  : std::vector<std::uint32_t>(num_cells);
90 
91  T value(0);
92  for (int i : M.integral_ids(IntegralType::cell))
93  {
94  const auto& fn = M.kernel(IntegralType::cell, i);
95  const std::vector<std::int32_t>& active_cells
96  = M.domains(IntegralType::cell, i);
97  value += fem::impl::assemble_cells(mesh->geometry(), active_cells, fn,
98  coeffs, constant_values, cell_info);
99  }
100 
101  if (M.num_integrals(IntegralType::exterior_facet) > 0
102  or M.num_integrals(IntegralType::interior_facet) > 0)
103  {
104  // FIXME: cleanup these calls? Some of these happen internally again.
105  mesh->topology_mutable().create_entities(tdim - 1);
106  mesh->topology_mutable().create_connectivity(tdim - 1, tdim);
107  mesh->topology_mutable().create_entity_permutations();
108 
109  const std::vector<std::uint8_t>& perms
110  = mesh->topology().get_facet_permutations();
111 
112  for (int i : M.integral_ids(IntegralType::exterior_facet))
113  {
114  const auto& fn = M.kernel(IntegralType::exterior_facet, i);
115  const std::vector<std::int32_t>& active_facets
116  = M.domains(IntegralType::exterior_facet, i);
117  value += fem::impl::assemble_exterior_facets(
118  *mesh, active_facets, fn, coeffs, constant_values, cell_info, perms);
119  }
120 
121  const std::vector<int> c_offsets = M.coefficient_offsets();
122  for (int i : M.integral_ids(IntegralType::interior_facet))
123  {
124  const auto& fn = M.kernel(IntegralType::interior_facet, i);
125  const std::vector<std::int32_t>& active_facets
126  = M.domains(IntegralType::interior_facet, i);
127  value += fem::impl::assemble_interior_facets(
128  *mesh, active_facets, fn, coeffs, c_offsets, constant_values,
129  cell_info, perms);
130  }
131  }
132 
133  return value;
134 }
135 //-----------------------------------------------------------------------------
136 template <typename T>
137 T assemble_cells(
138  const mesh::Geometry& geometry,
139  const std::vector<std::int32_t>& active_cells,
140  const std::function<void(T*, const T*, const T*, const double*, const int*,
141  const std::uint8_t*, const std::uint32_t)>& fn,
142  const array2d<T>& coeffs, const std::vector<T>& constant_values,
143  const std::vector<std::uint32_t>& cell_info)
144 {
145  const std::size_t gdim = geometry.dim();
146 
147  // Prepare cell geometry
148  const graph::AdjacencyList<std::int32_t>& x_dofmap = geometry.dofmap();
149 
150  // FIXME: Add proper interface for num coordinate dofs
151  const std::size_t num_dofs_g = x_dofmap.num_links(0);
152  const xt::xtensor<double, 2>& x_g = geometry.x();
153 
154  // Create data structures used in assembly
155  std::vector<double> coordinate_dofs(num_dofs_g * gdim);
156 
157  // Iterate over all cells
158  T value(0);
159  for (std::int32_t c : active_cells)
160  {
161  // Get cell coordinates/geometry
162  auto x_dofs = x_dofmap.links(c);
163  for (std::size_t i = 0; i < x_dofs.size(); ++i)
164  {
165  std::copy_n(xt::row(x_g, x_dofs[i]).begin(), gdim,
166  std::next(coordinate_dofs.begin(), i * gdim));
167  }
168 
169  auto coeff_cell = coeffs.row(c);
170  fn(&value, coeff_cell.data(), constant_values.data(),
171  coordinate_dofs.data(), nullptr, nullptr, cell_info[c]);
172  }
173 
174  return value;
175 }
176 //-----------------------------------------------------------------------------
177 template <typename T>
178 T assemble_exterior_facets(
179  const mesh::Mesh& mesh, const std::vector<std::int32_t>& active_facets,
180  const std::function<void(T*, const T*, const T*, const double*, const int*,
181  const std::uint8_t*, const std::uint32_t)>& fn,
182  const array2d<T>& coeffs, const std::vector<T>& constant_values,
183  const std::vector<std::uint32_t>& cell_info,
184  const std::vector<std::uint8_t>& perms)
185 {
186  const std::size_t gdim = mesh.geometry().dim();
187  const int tdim = mesh.topology().dim();
188 
189  // Prepare cell geometry
190  const graph::AdjacencyList<std::int32_t>& x_dofmap = mesh.geometry().dofmap();
191 
192  // FIXME: Add proper interface for num coordinate dofs
193  const std::size_t num_dofs_g = x_dofmap.num_links(0);
194  const xt::xtensor<double, 2>& x_g = mesh.geometry().x();
195 
196  // Creat data structures used in assembly
197  std::vector<double> coordinate_dofs(num_dofs_g * gdim);
198 
199  auto f_to_c = mesh.topology().connectivity(tdim - 1, tdim);
200  assert(f_to_c);
201  auto c_to_f = mesh.topology().connectivity(tdim, tdim - 1);
202  assert(c_to_f);
203 
204  // Iterate over all facets
205  T value(0);
206  for (std::int32_t facet : active_facets)
207  {
208  // Create attached cell
209  assert(f_to_c->num_links(facet) == 1);
210  const int cell = f_to_c->links(facet)[0];
211 
212  // Get local index of facet with respect to the cell
213  auto facets = c_to_f->links(cell);
214  auto it = std::find(facets.begin(), facets.end(), facet);
215  assert(it != facets.end());
216  const int local_facet = std::distance(facets.data(), it);
217 
218  // Get cell coordinates/geometry
219  auto x_dofs = x_dofmap.links(cell);
220  for (std::size_t i = 0; i < x_dofs.size(); ++i)
221  {
222  std::copy_n(xt::row(x_g, x_dofs[i]).begin(), gdim,
223  std::next(coordinate_dofs.begin(), i * gdim));
224  }
225 
226  auto coeff_cell = coeffs.row(cell);
227  fn(&value, coeff_cell.data(), constant_values.data(),
228  coordinate_dofs.data(), &local_facet,
229  &perms[cell * facets.size() + local_facet], cell_info[cell]);
230  }
231 
232  return value;
233 }
234 //-----------------------------------------------------------------------------
235 template <typename T>
236 T assemble_interior_facets(
237  const mesh::Mesh& mesh, const std::vector<std::int32_t>& active_facets,
238  const std::function<void(T*, const T*, const T*, const double*, const int*,
239  const std::uint8_t*, const std::uint32_t)>& fn,
240  const array2d<T>& coeffs, const std::vector<int>& offsets,
241  const std::vector<T>& constant_values,
242  const std::vector<std::uint32_t>& cell_info,
243  const std::vector<std::uint8_t>& perms)
244 {
245  const std::size_t gdim = mesh.geometry().dim();
246  const int tdim = mesh.topology().dim();
247 
248  // Prepare cell geometry
249  const graph::AdjacencyList<std::int32_t>& x_dofmap = mesh.geometry().dofmap();
250 
251  // FIXME: Add proper interface for num coordinate dofs
252  const std::size_t num_dofs_g = x_dofmap.num_links(0);
253  const xt::xtensor<double, 2>& x_g = mesh.geometry().x();
254 
255  // Creat data structures used in assembly
256  xt::xtensor<double, 3> coordinate_dofs({2, num_dofs_g, gdim});
257  std::vector<T> coeff_array(2 * offsets.back());
258  assert(offsets.back() == coeffs.shape[1]);
259 
260  auto f_to_c = mesh.topology().connectivity(tdim - 1, tdim);
261  assert(f_to_c);
262  auto c_to_f = mesh.topology().connectivity(tdim, tdim - 1);
263  assert(c_to_f);
264 
265  // Iterate over all facets
266  T value = 0;
267  for (std::int32_t f : active_facets)
268  {
269  // Create attached cell
270  auto cells = f_to_c->links(f);
271  assert(cells.size() == 2);
272 
273  const int facets_per_cell = c_to_f->num_links(cells[0]);
274 
275  // Get local index of facet with respect to the cell
276  std::array<int, 2> local_facet;
277  for (int i = 0; i < 2; ++i)
278  {
279  auto facets = c_to_f->links(cells[i]);
280  auto it = std::find(facets.begin(), facets.end(), f);
281  assert(it != facets.end());
282  local_facet[i] = std::distance(facets.begin(), it);
283  }
284 
285  // Get cell geometry
286  auto x_dofs0 = x_dofmap.links(cells[0]);
287  for (std::size_t i = 0; i < x_dofs0.size(); ++i)
288  {
289  std::copy_n(xt::view(x_g, x_dofs0[i]).begin(), gdim,
290  xt::view(coordinate_dofs, 0, i, xt::all()).begin());
291  }
292  auto x_dofs1 = x_dofmap.links(cells[1]);
293  for (std::size_t i = 0; i < x_dofs1.size(); ++i)
294  {
295  std::copy_n(xt::view(x_g, x_dofs1[i]).begin(), gdim,
296  xt::view(coordinate_dofs, 1, i, xt::all()).begin());
297  }
298 
299  // Layout for the restricted coefficients is flattened
300  // w[coefficient][restriction][dof]
301  auto coeff_cell0 = coeffs.row(cells[0]);
302  auto coeff_cell1 = coeffs.row(cells[1]);
303 
304  // Loop over coefficients
305  for (std::size_t i = 0; i < offsets.size() - 1; ++i)
306  {
307  // Loop over entries for coefficient i
308  const int num_entries = offsets[i + 1] - offsets[i];
309  std::copy_n(coeff_cell0.data() + offsets[i], num_entries,
310  std::next(coeff_array.begin(), 2 * offsets[i]));
311  std::copy_n(coeff_cell1.data() + offsets[i], num_entries,
312  std::next(coeff_array.begin(), offsets[i + 1] + offsets[i]));
313  }
314 
315  const std::array perm{perms[cells[0] * facets_per_cell + local_facet[0]],
316  perms[cells[1] * facets_per_cell + local_facet[1]]};
317  fn(&value, coeff_array.data(), constant_values.data(),
318  coordinate_dofs.data(), local_facet.data(), perm.data(),
319  cell_info[cells[0]]);
320  }
321 
322  return value;
323 }
324 //-----------------------------------------------------------------------------
325 
326 } // namespace dolfinx::fem::impl
dolfinx::fem::sparsitybuild::cells
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
dolfinx::fem::pack_coefficients
array2d< typename U::scalar_type > pack_coefficients(const U &u)
Pack coefficients of u of generic type U ready for assembly.
Definition: utils.h:376
dolfinx::fem::assemble_scalar
T assemble_scalar(const Form< T > &M)
Assemble functional into scalar. Caller is responsible for accumulation across processes.
Definition: assembler.h:33