DOLFINx 0.9.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
discreteoperators.h
1// Copyright (C) 2015-2022 Garth N. Wells, Jørgen S. Dokken
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 "DofMap.h"
10#include "FiniteElement.h"
11#include "FunctionSpace.h"
12#include <algorithm>
13#include <array>
14#include <concepts>
15#include <dolfinx/common/IndexMap.h>
16#include <dolfinx/common/math.h>
17#include <dolfinx/mesh/Mesh.h>
18#include <memory>
19#include <span>
20#include <vector>
21
22namespace dolfinx::fem
23{
24
51template <dolfinx::scalar T,
52 std::floating_point U = dolfinx::scalar_value_type_t<T>>
54 std::pair<std::reference_wrapper<const FiniteElement<U>>,
55 std::reference_wrapper<const DofMap>>
56 V0,
57 std::pair<std::reference_wrapper<const FiniteElement<U>>,
58 std::reference_wrapper<const DofMap>>
59 V1,
60 auto&& mat_set)
61{
62 auto& e0 = V0.first.get();
63 const DofMap& dofmap0 = V0.second.get();
64 auto& e1 = V1.first.get();
65 const DofMap& dofmap1 = V1.second.get();
66
67 using cmdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
68 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
69 using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
70 U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
71 using cmdspan4_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
72 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 4>>;
73
74 // Check elements
75 if (e0.map_type() != basix::maps::type::identity)
76 throw std::runtime_error("Wrong finite element space for V0.");
77 if (e0.block_size() != 1)
78 throw std::runtime_error("Block size is greater than 1 for V0.");
79 if (e0.reference_value_size() != 1)
80 throw std::runtime_error("Wrong value size for V0.");
81
82 if (e1.map_type() != basix::maps::type::covariantPiola)
83 throw std::runtime_error("Wrong finite element space for V1.");
84 if (e1.block_size() != 1)
85 throw std::runtime_error("Block size is greater than 1 for V1.");
86
87 // Get V0 (H(curl)) space interpolation points
88 const auto [X, Xshape] = e1.interpolation_points();
89
90 // Tabulate first order derivatives of Lagrange space at H(curl)
91 // interpolation points
92 const int ndofs0 = e0.space_dimension();
93 const int tdim = topology.dim();
94 std::vector<U> phi0_b((tdim + 1) * Xshape[0] * ndofs0 * 1);
95 cmdspan4_t phi0(phi0_b.data(), tdim + 1, Xshape[0], ndofs0, 1);
96 e0.tabulate(phi0_b, X, Xshape, 1);
97
98 // Reshape lagrange basis derivatives as a matrix of shape (tdim *
99 // num_points, num_dofs_per_cell)
100 cmdspan2_t dphi_reshaped(
101 phi0_b.data() + phi0.extent(3) * phi0.extent(2) * phi0.extent(1),
102 tdim * phi0.extent(1), phi0.extent(2));
103
104 // Get inverse DOF transform function
105 auto apply_inverse_dof_transform = e1.template dof_transformation_fn<T>(
107
108 // Generate cell permutations
110 const std::vector<std::uint32_t>& cell_info
111 = topology.get_cell_permutation_info();
112
113 // Create element kernel function
114
115 // Build the element interpolation matrix
116 std::vector<T> Ab(e1.space_dimension() * ndofs0);
117 {
118 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
119 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>
120 A(Ab.data(), e1.space_dimension(), ndofs0);
121 const auto [Pi, shape] = e1.interpolation_operator();
122 cmdspan2_t _Pi(Pi.data(), shape);
123 math::dot(_Pi, dphi_reshaped, A);
124 }
125
126 // Insert local interpolation matrix for each cell
127 auto cell_map = topology.index_map(tdim);
128 assert(cell_map);
129 std::int32_t num_cells = cell_map->size_local();
130 std::vector<T> Ae(Ab.size());
131 for (std::int32_t c = 0; c < num_cells; ++c)
132 {
133 std::ranges::copy(Ab, Ae.begin());
134 apply_inverse_dof_transform(Ae, cell_info, c, ndofs0);
135 mat_set(dofmap1.cell_dofs(c), dofmap0.cell_dofs(c), Ae);
136 }
137}
138
154template <dolfinx::scalar T, std::floating_point U>
156 const FunctionSpace<U>& V1, auto&& mat_set)
157{
158 // Get mesh
159 auto mesh = V0.mesh();
160 assert(mesh);
161
162 // Mesh dims
163 const int tdim = mesh->topology()->dim();
164 const int gdim = mesh->geometry().dim();
165
166 // Get elements
167 std::shared_ptr<const FiniteElement<U>> e0 = V0.element();
168 assert(e0);
169 std::shared_ptr<const FiniteElement<U>> e1 = V1.element();
170 assert(e1);
171
172 std::span<const std::uint32_t> cell_info;
173 if (e1->needs_dof_transformations() or e0->needs_dof_transformations())
174 {
175 mesh->topology_mutable()->create_entity_permutations();
176 cell_info = std::span(mesh->topology()->get_cell_permutation_info());
177 }
178
179 // Get dofmaps
180 auto dofmap0 = V0.dofmap();
181 assert(dofmap0);
182 auto dofmap1 = V1.dofmap();
183 assert(dofmap1);
184
185 // Get block sizes and dof transformation operators
186 const int bs0 = e0->block_size();
187 const int bs1 = e1->block_size();
188 auto apply_dof_transformation0
189 = e0->template dof_transformation_fn<U>(doftransform::standard, false);
190 auto apply_inverse_dof_transform1 = e1->template dof_transformation_fn<T>(
192
193 // Get sizes of elements
194 const std::size_t space_dim0 = e0->space_dimension();
195 const std::size_t space_dim1 = e1->space_dimension();
196 const std::size_t dim0 = space_dim0 / bs0;
197 const std::size_t value_size_ref0 = e0->reference_value_size() / bs0;
198 const std::size_t value_size0 = V0.value_size() / bs0;
199 const std::size_t value_size1 = V1.value_size() / bs1;
200
201 // Get geometry data
202 const CoordinateElement<U>& cmap = mesh->geometry().cmap();
203 auto x_dofmap = mesh->geometry().dofmap();
204 const std::size_t num_dofs_g = cmap.dim();
205 std::span<const U> x_g = mesh->geometry().x();
206
207 using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
208 U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
209 using cmdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
210 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
211 using cmdspan3_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
212 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>;
213 using cmdspan4_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
214 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 4>>;
215 using mdspan3_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
216 U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>;
217
218 // Evaluate coordinate map basis at reference interpolation points
219 const auto [X, Xshape] = e1->interpolation_points();
220 std::array<std::size_t, 4> phi_shape = cmap.tabulate_shape(1, Xshape[0]);
221 std::vector<U> phi_b(
222 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
223 cmdspan4_t phi(phi_b.data(), phi_shape);
224 cmap.tabulate(1, X, Xshape, phi_b);
225
226 // Evaluate V0 basis functions at reference interpolation points for V1
227 std::vector<U> basis_derivatives_reference0_b(Xshape[0] * dim0
228 * value_size_ref0);
229 cmdspan4_t basis_derivatives_reference0(basis_derivatives_reference0_b.data(),
230 1, Xshape[0], dim0, value_size_ref0);
231 e0->tabulate(basis_derivatives_reference0_b, X, Xshape, 0);
232
233 // Clamp values
234 std::ranges::transform(
235 basis_derivatives_reference0_b, basis_derivatives_reference0_b.begin(),
236 [atol = 1e-14](auto x) { return std::abs(x) < atol ? 0.0 : x; });
237
238 // Create working arrays
239 std::vector<U> basis_reference0_b(Xshape[0] * dim0 * value_size_ref0);
240 mdspan3_t basis_reference0(basis_reference0_b.data(), Xshape[0], dim0,
241 value_size_ref0);
242 std::vector<U> J_b(Xshape[0] * gdim * tdim);
243 mdspan3_t J(J_b.data(), Xshape[0], gdim, tdim);
244 std::vector<U> K_b(Xshape[0] * tdim * gdim);
245 mdspan3_t K(K_b.data(), Xshape[0], tdim, gdim);
246 std::vector<U> detJ(Xshape[0]);
247 std::vector<U> det_scratch(2 * tdim * gdim);
248
249 // Get the interpolation operator (matrix) `Pi` that maps a function
250 // evaluated at the interpolation points to the element degrees of
251 // freedom, i.e. dofs = Pi f_x
252 const auto [_Pi_1, pi_shape] = e1->interpolation_operator();
253 cmdspan2_t Pi_1(_Pi_1.data(), pi_shape);
254
255 bool interpolation_ident = e1->interpolation_ident();
256
257 using u_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
258 U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
259 using U_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
260 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
261 using J_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
262 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
263 using K_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
264 const U, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
265 auto push_forward_fn0
266 = e0->basix_element().template map_fn<u_t, U_t, J_t, K_t>();
267
268 // Basis values of Lagrange space unrolled for block size
269 // (num_quadrature_points, Lagrange dof, value_size)
270 std::vector<U> basis_values_b(Xshape[0] * bs0 * dim0 * V1.value_size());
271 mdspan3_t basis_values(basis_values_b.data(), Xshape[0], bs0 * dim0,
272 V1.value_size());
273 std::vector<U> mapped_values_b(Xshape[0] * bs0 * dim0 * V1.value_size());
274 mdspan3_t mapped_values(mapped_values_b.data(), Xshape[0], bs0 * dim0,
275 V1.value_size());
276
277 auto pull_back_fn1
278 = e1->basix_element().template map_fn<u_t, U_t, K_t, J_t>();
279
280 std::vector<U> coord_dofs_b(num_dofs_g * gdim);
281 mdspan2_t coord_dofs(coord_dofs_b.data(), num_dofs_g, gdim);
282 std::vector<U> basis0_b(Xshape[0] * dim0 * value_size0);
283 mdspan3_t basis0(basis0_b.data(), Xshape[0], dim0, value_size0);
284
285 // Buffers
286 std::vector<T> Ab(space_dim0 * space_dim1);
287 std::vector<T> local1(space_dim1);
288
289 // Iterate over mesh and interpolate on each cell
290 auto cell_map = mesh->topology()->index_map(tdim);
291 assert(cell_map);
292 std::int32_t num_cells = cell_map->size_local();
293 for (std::int32_t c = 0; c < num_cells; ++c)
294 {
295 // Get cell geometry (coordinate dofs)
296 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
297 x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
298 for (std::size_t i = 0; i < x_dofs.size(); ++i)
299 {
300 for (std::size_t j = 0; j < gdim; ++j)
301 coord_dofs(i, j) = x_g[3 * x_dofs[i] + j];
302 }
303
304 // Compute Jacobians and reference points for current cell
305 std::ranges::fill(J_b, 0);
306 for (std::size_t p = 0; p < Xshape[0]; ++p)
307 {
308 auto dphi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
309 phi, std::pair(1, tdim + 1), p,
310 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
311 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
312 J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
313 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
314 cmap.compute_jacobian(dphi, coord_dofs, _J);
315 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
316 K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
317 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
318 cmap.compute_jacobian_inverse(_J, _K);
319 detJ[p] = cmap.compute_jacobian_determinant(_J, det_scratch);
320 }
321
322 // Copy evaluated basis on reference, apply DOF transformations, and
323 // push forward to physical element
324 for (std::size_t k0 = 0; k0 < basis_reference0.extent(0); ++k0)
325 for (std::size_t k1 = 0; k1 < basis_reference0.extent(1); ++k1)
326 for (std::size_t k2 = 0; k2 < basis_reference0.extent(2); ++k2)
327 basis_reference0(k0, k1, k2)
328 = basis_derivatives_reference0(0, k0, k1, k2);
329 for (std::size_t p = 0; p < Xshape[0]; ++p)
330 {
331 apply_dof_transformation0(
332 std::span(basis_reference0.data_handle() + p * dim0 * value_size_ref0,
333 dim0 * value_size_ref0),
334 cell_info, c, value_size_ref0);
335 }
336
337 for (std::size_t p = 0; p < basis0.extent(0); ++p)
338 {
339 auto _u = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
340 basis0, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
341 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
342 auto _U = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
343 basis_reference0, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
344 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
345 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
346 K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
347 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
348 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
349 J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
350 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
351 push_forward_fn0(_u, _U, _J, detJ[p], _K);
352 }
353
354 // Unroll basis function for input space for block size
355 for (std::size_t p = 0; p < Xshape[0]; ++p)
356 for (std::size_t i = 0; i < dim0; ++i)
357 for (std::size_t j = 0; j < value_size0; ++j)
358 for (int k = 0; k < bs0; ++k)
359 basis_values(p, i * bs0 + k, j * bs0 + k) = basis0(p, i, j);
360
361 // Pull back the physical values to the reference of output space
362 for (std::size_t p = 0; p < basis_values.extent(0); ++p)
363 {
364 auto _u = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
365 basis_values, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
366 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
367 auto _U = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
368 mapped_values, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
369 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
370 auto _K = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
371 K, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
372 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
373 auto _J = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
374 J, p, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
375 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
376 pull_back_fn1(_U, _u, _K, 1.0 / detJ[p], _J);
377 }
378
379 // Apply interpolation matrix to basis values of V0 at the
380 // interpolation points of V1
381 if (interpolation_ident)
382 {
383 MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
384 T, MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 3>>
385 A(Ab.data(), Xshape[0], V1.value_size(), space_dim0);
386 for (std::size_t i = 0; i < mapped_values.extent(0); ++i)
387 for (std::size_t j = 0; j < mapped_values.extent(1); ++j)
388 for (std::size_t k = 0; k < mapped_values.extent(2); ++k)
389 A(i, k, j) = mapped_values(i, j, k);
390 }
391 else
392 {
393 for (std::size_t i = 0; i < mapped_values.extent(1); ++i)
394 {
395 auto values = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
396 mapped_values, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, i,
397 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
398 impl::interpolation_apply(Pi_1, values, std::span(local1), bs1);
399 for (std::size_t j = 0; j < local1.size(); j++)
400 Ab[space_dim0 * j + i] = local1[j];
401 }
402 }
403
404 apply_inverse_dof_transform1(Ab, cell_info, c, space_dim0);
405 mat_set(dofmap1->cell_dofs(c), dofmap0->cell_dofs(c), Ab);
406 }
407}
408
409} // namespace dolfinx::fem
Degree-of-freedom map representations and tools.
Definition XDMFFile.h:27
static void compute_jacobian(const U &dphi, const V &cell_geometry, W &&J)
Definition CoordinateElement.h:126
void tabulate(int nd, std::span< const T > X, std::array< std::size_t, 2 > shape, std::span< T > basis) const
Evaluate basis values and derivatives at set of points.
Definition CoordinateElement.cpp:55
static void compute_jacobian_inverse(const U &J, V &&K)
Compute the inverse of the Jacobian.
Definition CoordinateElement.h:135
std::array< std::size_t, 4 > tabulate_shape(std::size_t nd, std::size_t num_points) const
Shape of array to fill when calling tabulate.
Definition CoordinateElement.cpp:48
int dim() const
The dimension of the coordinate element space.
Definition CoordinateElement.cpp:206
static double compute_jacobian_determinant(const U &J, std::span< typename U::value_type > w)
Compute the determinant of the Jacobian.
Definition CoordinateElement.h:152
Degree-of-freedom map.
Definition DofMap.h:76
std::span< const std::int32_t > cell_dofs(std::int32_t c) const
Local-to-global mapping of dofs on a cell.
Definition DofMap.h:130
Model of a finite element.
Definition FiniteElement.h:38
This class represents a finite element function space defined by a mesh, a finite element,...
Definition vtk_utils.h:32
std::shared_ptr< const DofMap > dofmap() const
The dofmap.
Definition FunctionSpace.h:335
int value_size() const
Definition FunctionSpace.h:348
std::shared_ptr< const mesh::Mesh< geometry_type > > mesh() const
The mesh.
Definition FunctionSpace.h:323
std::shared_ptr< const FiniteElement< geometry_type > > element() const
The finite element.
Definition FunctionSpace.h:329
Topology stores the topology of a mesh, consisting of mesh entities and connectivity (incidence relat...
Definition Topology.h:44
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:815
void create_entity_permutations()
Compute entity permutations and reflections.
Definition Topology.cpp:910
const std::vector< std::uint32_t > & get_cell_permutation_info() const
Returns the permutation information.
Definition Topology.cpp:981
int dim() const noexcept
Return the topological dimension of the mesh.
Definition Topology.cpp:794
Definition types.h:20
Finite element method functionality.
Definition assemble_matrix_impl.h:26
void discrete_gradient(mesh::Topology &topology, std::pair< std::reference_wrapper< const FiniteElement< U > >, std::reference_wrapper< const DofMap > > V0, std::pair< std::reference_wrapper< const FiniteElement< U > >, std::reference_wrapper< const DofMap > > V1, auto &&mat_set)
Assemble a discrete gradient operator.
Definition discreteoperators.h:53
@ inverse_transpose
Transpose inverse.
void interpolation_matrix(const FunctionSpace< U > &V0, const FunctionSpace< U > &V1, auto &&mat_set)
Assemble an interpolation operator matrix.
Definition discreteoperators.h:155