DOLFINx 0.11.0.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
FunctionSpace.h
1// Copyright (C) 2008-2023 Anders Logg and 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 "CoordinateElement.h"
10#include "DofMap.h"
11#include "FiniteElement.h"
12#include <boost/uuid/uuid.hpp>
13#include <boost/uuid/uuid_generators.hpp>
14#include <concepts>
15#include <cstddef>
16#include <cstdint>
17#include <dolfinx/common/IndexMap.h>
18#include <dolfinx/mesh/Geometry.h>
19#include <dolfinx/mesh/Mesh.h>
20#include <dolfinx/mesh/Topology.h>
21#include <map>
22#include <memory>
23#include <vector>
24
25namespace dolfinx::fem
26{
32template <std::floating_point T>
34{
35public:
37 using geometry_type = T;
38
45 std::shared_ptr<const FiniteElement<geometry_type>> element,
46 std::shared_ptr<const DofMap> dofmap)
47 : _mesh(mesh), _elements{element}, _dofmaps{std::move(dofmap)},
48 _id(boost::uuids::random_generator()()), _root_space_id(_id)
49 {
50 // Do nothing
51 }
52
62 std::shared_ptr<const mesh::Mesh<geometry_type>> mesh,
63 std::vector<std::shared_ptr<const FiniteElement<geometry_type>>> elements,
64 std::vector<std::shared_ptr<const DofMap>> dofmaps)
65 : _mesh(mesh), _elements(elements), _dofmaps(std::move(dofmaps)),
66 _id(boost::uuids::random_generator()()), _root_space_id(_id)
67 {
68 std::vector<mesh::CellType> cell_types = mesh->topology()->cell_types();
69 std::size_t num_cell_types = cell_types.size();
70 if (elements.size() != num_cell_types)
71 {
72 throw std::runtime_error(
73 "Number of elements must match number of cell types");
74 }
75 if (_dofmaps.size() != num_cell_types)
76 {
77 throw std::runtime_error(
78 "Number of dofmaps must match number of cell types");
79 }
80 for (std::size_t i = 0; i < num_cell_types; ++i)
81 {
82 if (elements.at(i)->cell_type() != cell_types.at(i))
83 throw std::runtime_error(
84 "Element cell types must match mesh cell types");
85 }
86 }
87
88 // Copy constructor (deleted)
89 FunctionSpace(const FunctionSpace& V) = delete;
90
93
95 virtual ~FunctionSpace() = default;
96
97 // Assignment operator (delete)
98 FunctionSpace& operator=(const FunctionSpace& V) = delete;
99
102
111 FunctionSpace sub(const std::vector<int>& component) const
112 {
113 assert(_mesh);
114 assert(_elements.size() > 0);
115 assert(_dofmaps.size() > 0);
116
117 // Check that component is valid
118 if (component.empty())
119 throw std::runtime_error("Component must be non-empty");
120
121 // Extract sub-element
122 std::vector<std::shared_ptr<const FiniteElement<geometry_type>>>
123 sub_elements;
124 for (auto e : _elements)
125 {
126 assert(e);
127 sub_elements.push_back(e->extract_sub_element(component));
128 }
129 // Extract sub dofmap
130 std::vector<std::shared_ptr<const DofMap>> sub_dofmaps;
131 for (auto d : _dofmaps)
132 {
133 assert(d);
134 sub_dofmaps.push_back(
135 std::make_shared<const DofMap>(d->extract_sub_dofmap(component)));
136 }
137
138 // Create new sub space
139 FunctionSpace sub_space(_mesh, sub_elements, sub_dofmaps);
140
141 // Set root space id and component w.r.t. root
142 sub_space._root_space_id = _root_space_id;
143 sub_space._component = _component;
144 sub_space._component.insert(sub_space._component.end(), component.begin(),
145 component.end());
146 return sub_space;
147 }
148
153 bool contains(const FunctionSpace& V) const
154 {
155 if (this == std::addressof(V)) // Spaces are the same (same memory address)
156 return true;
157 else if (_root_space_id != V._root_space_id) // Different root spaces
158 return false;
159 else if (_component.size()
160 > V._component.size()) // V is a superspace of *this
161 {
162 return false;
163 }
164 else if (!std::equal(_component.begin(), _component.end(),
165 V._component.begin())) // Components of 'this' are not
166 // the same as the leading
167 // components of V
168 {
169 return false;
170 }
171 else // Ok, V is really our subspace
172 return true;
173 }
174
178 std::pair<FunctionSpace, std::vector<std::vector<std::int32_t>>>
179 collapse() const
180 {
181 spdlog::debug("FunctionSpace::collapse");
182 if (_component.empty())
183 throw std::runtime_error("Function space is not a subspace");
184
185 // Create collapsed DofMap
186 std::vector<std::shared_ptr<const DofMap>> collapsed_dofmaps;
187 std::vector<std::vector<std::int32_t>> collapsed_dofs;
188 for (auto d : _dofmaps)
189 {
190 assert(d);
191 spdlog::debug("Call DofMap::collapse");
192 auto [_collapsed_dofmap, _collapsed_dofs]
193 = d->collapse(_mesh->comm(), *_mesh->topology());
194 spdlog::debug("Got collapsed dofmap");
195 collapsed_dofmaps.push_back(
196 std::make_shared<DofMap>(std::move(_collapsed_dofmap)));
197 collapsed_dofs.push_back(_collapsed_dofs);
198 }
199
200 return {FunctionSpace(_mesh, _elements, collapsed_dofmaps),
201 std::move(collapsed_dofs)};
202 }
203
207 std::vector<int> component() const { return _component; }
208
211 bool symmetric() const
212 {
213 if (_elements.front())
214 return _elements.front()->symmetric();
215 return false;
216 }
217
229 std::vector<geometry_type> tabulate_dof_coordinates(bool transpose) const
230 {
231 if (!_component.empty())
232 {
233 throw std::runtime_error("Cannot tabulate coordinates for a "
234 "FunctionSpace that is a subspace.");
235 }
236
237 assert(_elements.front());
238 if (_elements.front()->is_mixed())
239 {
240 throw std::runtime_error(
241 "Cannot tabulate coordinates for a mixed FunctionSpace.");
242 }
243
244 // Geometric dimension
245 assert(_mesh);
246 assert(_elements.front());
247 const std::size_t gdim = _mesh->geometry().dim();
248 const int tdim = _mesh->topology()->dim();
249
250 // Get dofmap local size
251 assert(_dofmaps.front());
252 std::shared_ptr<const common::IndexMap> index_map
253 = _dofmaps.front()->index_map;
254 assert(index_map);
255 const int index_map_bs = _dofmaps.front()->index_map_bs();
256 const int dofmap_bs = _dofmaps.front()->bs();
257
258 const int element_block_size = _elements.front()->block_size();
259 const std::size_t scalar_dofs
260 = _elements.front()->space_dimension() / element_block_size;
261 const std::int32_t num_dofs
262 = index_map_bs * (index_map->size_local() + index_map->num_ghosts())
263 / dofmap_bs;
264
265 std::size_t num_cell_types = _elements.size();
266
267 // Array to hold coordinates to return
268 const std::size_t shape_c0 = transpose ? 3 : num_dofs;
269 const std::size_t shape_c1 = transpose ? num_dofs : 3;
270 std::vector<geometry_type> coords(shape_c0 * shape_c1, 0);
271 using mdspan2_t = md::mdspan<geometry_type, md::dextents<std::size_t, 2>>;
272 using cmdspan4_t
273 = md::mdspan<const geometry_type, md::dextents<std::size_t, 4>>;
274 std::span<const geometry_type> x_g = _mesh->geometry().x();
275 std::vector<geometry_type> x_b(scalar_dofs * gdim);
276 mdspan2_t x(x_b.data(), scalar_dofs, gdim);
277
278 for (std::size_t i = 0; i < num_cell_types; ++i)
279 {
280 // Get the dof coordinates on the reference element
281 if (!_elements[i]->interpolation_ident())
282 {
283 throw std::runtime_error(
284 "Cannot evaluate dof coordinates - this element "
285 "does not have pointwise evaluation.");
286 }
287 const auto [X, Xshape] = _elements[i]->interpolation_points();
288
289 // Get coordinate map
291 = _mesh->geometry().cmaps()[i];
292
293 // Prepare cell geometry
294 auto x_dofmap = _mesh->geometry().dofmap(i);
295 const std::size_t num_dofs_g = cmap.dim();
296 std::vector<geometry_type> coordinate_dofs_b(num_dofs_g * gdim);
297 mdspan2_t coordinate_dofs(coordinate_dofs_b.data(), num_dofs_g, gdim);
298
299 auto map = _mesh->topology()->index_maps(tdim)[i];
300 assert(map);
301 const int num_cells = map->size_local() + map->num_ghosts();
302
303 std::span<const std::uint32_t> cell_info;
304 if (_elements[i]->needs_dof_transformations())
305 {
306 _mesh->topology_mutable()->create_entity_permutations();
307 cell_info = std::span(_mesh->topology()->get_cell_permutation_info());
308 }
309
310 const std::array<std::size_t, 4> phi_shape
311 = cmap.tabulate_shape(0, Xshape[0]);
312 std::vector<geometry_type> phi_b(std::reduce(
313 phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
314 cmdspan4_t phi_full(phi_b.data(), phi_shape);
315 cmap.tabulate(0, X, Xshape, phi_b);
316 auto phi
317 = md::submdspan(phi_full, 0, md::full_extent, md::full_extent, 0);
318
319 // TODO: Check transform
320 // Basis function reference-to-conforming transformation function
321 auto apply_dof_transformation
322 = _elements[i]->template dof_transformation_fn<geometry_type>(
324
325 for (int c = 0; c < num_cells; ++c)
326 {
327 // Extract cell geometry 'dofs'
328 auto x_dofs = md::submdspan(x_dofmap, c, md::full_extent);
329 for (std::size_t j = 0; j < x_dofs.size(); ++j)
330 for (std::size_t k = 0; k < gdim; ++k)
331 coordinate_dofs(j, k) = x_g[3 * x_dofs[j] + k];
332
333 // Tabulate dof coordinates on cell
334 cmap.push_forward(x, coordinate_dofs, phi);
335 apply_dof_transformation(
336 x_b, std::span(cell_info.data(), cell_info.size()), c, x.extent(1));
337
338 // Get cell dofmap
339 auto dofs = _dofmaps[i]->cell_dofs(c);
340
341 // Copy dof coordinates into vector
342 if (!transpose)
343 {
344 for (std::size_t j = 0; j < dofs.size(); ++j)
345 for (std::size_t k = 0; k < gdim; ++k)
346 coords[dofs[j] * 3 + k] = x(j, k);
347 }
348 else
349 {
350 for (std::size_t j = 0; j < dofs.size(); ++j)
351 for (std::size_t k = 0; k < gdim; ++k)
352 coords[k * num_dofs + dofs[j]] = x(j, k);
353 }
354 }
355 }
356
357 return coords;
358 }
359
361 std::shared_ptr<const mesh::Mesh<geometry_type>> mesh() const
362 {
363 return _mesh;
364 }
365
367 std::shared_ptr<const FiniteElement<geometry_type>> element() const
368 {
369 if (_elements.size() > 1)
370 {
371 throw std::runtime_error(
372 "FunctionSpace has multiple elements, call `elements` instead.");
373 }
374
375 return elements(0);
376 }
377
379 std::shared_ptr<const FiniteElement<geometry_type>>
380 elements(int cell_type_idx) const
381 {
382 return _elements.at(cell_type_idx);
383 }
384
386 std::shared_ptr<const DofMap> dofmap() const
387 {
388 if (_dofmaps.size() > 1)
389 {
390 throw std::runtime_error(
391 "FunctionSpace has multiple dofmaps, call `dofmaps` instead.");
392 }
393
394 return dofmaps(0);
395 }
396
398 std::shared_ptr<const DofMap> dofmaps(int cell_type_idx) const
399 {
400 return _dofmaps.at(cell_type_idx);
401 }
402
403private:
404 // The mesh
405 std::shared_ptr<const mesh::Mesh<geometry_type>> _mesh;
406
407 // The finite element
408 std::vector<std::shared_ptr<const FiniteElement<geometry_type>>> _elements;
409
410 // The dofmap
411 std::vector<std::shared_ptr<const DofMap>> _dofmaps;
412
413 // The component w.r.t. to root space
414 std::vector<int> _component;
415
416 // Unique identifier for the space and for its root space
417 boost::uuids::uuid _id;
418 boost::uuids::uuid _root_space_id;
419};
420
431template <dolfinx::scalar T>
432std::array<std::vector<std::shared_ptr<const FunctionSpace<T>>>, 2>
434 const std::vector<
435 std::vector<std::array<std::shared_ptr<const FunctionSpace<T>>, 2>>>& V)
436{
437 assert(!V.empty());
438 std::vector<std::shared_ptr<const FunctionSpace<T>>> spaces0(V.size(),
439 nullptr);
440 std::vector<std::shared_ptr<const FunctionSpace<T>>> spaces1(V.front().size(),
441 nullptr);
442
443 // Loop over rows
444 for (std::size_t i = 0; i < V.size(); ++i)
445 {
446 // Loop over columns
447 for (std::size_t j = 0; j < V[i].size(); ++j)
448 {
449 auto& V0 = V[i][j][0];
450 auto& V1 = V[i][j][1];
451 if (V0 and V1)
452 {
453 if (!spaces0[i])
454 spaces0[i] = V0;
455 else
456 {
457 if (spaces0[i] != V0)
458 throw std::runtime_error("Mismatched test space for row.");
459 }
460
461 if (!spaces1[j])
462 spaces1[j] = V1;
463 else
464 {
465 if (spaces1[j] != V1)
466 throw std::runtime_error("Mismatched trial space for column.");
467 }
468 }
469 }
470 }
471
472 // Check that there are no null entries
473 if (std::find(spaces0.begin(), spaces0.end(), nullptr) != spaces0.end())
474 throw std::runtime_error("Could not deduce all block test spaces.");
475 if (std::find(spaces1.begin(), spaces1.end(), nullptr) != spaces1.end())
476 throw std::runtime_error("Could not deduce all block trial spaces.");
477
478 return {spaces0, spaces1};
479}
480
482template <typename U, typename V, typename W>
483FunctionSpace(U mesh, V element, W dofmap)
484 -> FunctionSpace<typename std::remove_cvref<
485 typename U::element_type>::type::geometry_type::value_type>;
486
487} // namespace dolfinx::fem
Degree-of-freedom map representations and tools.
Definition CoordinateElement.h:38
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
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:205
static void push_forward(U &&x, const V &cell_geometry, const W &phi)
Compute physical coordinates x for points X in the reference configuration.
Definition CoordinateElement.h:194
Model of a finite element.
Definition FiniteElement.h:57
This class represents a finite element function space defined by a mesh, a finite element,...
Definition FunctionSpace.h:34
std::vector< int > component() const
Get the component with respect to the root superspace.
Definition FunctionSpace.h:207
bool symmetric() const
Indicate whether this function space represents a symmetric 2-tensor.
Definition FunctionSpace.h:211
FunctionSpace & operator=(FunctionSpace &&V)=default
Move assignment operator.
std::shared_ptr< const FiniteElement< geometry_type > > elements(int cell_type_idx) const
The finite elements.
Definition FunctionSpace.h:380
std::shared_ptr< const DofMap > dofmap() const
The dofmap.
Definition FunctionSpace.h:386
FunctionSpace(std::shared_ptr< const mesh::Mesh< geometry_type > > mesh, std::vector< std::shared_ptr< const FiniteElement< geometry_type > > > elements, std::vector< std::shared_ptr< const DofMap > > dofmaps)
Create function space for given mesh, elements and degree-of-freedom maps.
Definition FunctionSpace.h:61
FunctionSpace sub(const std::vector< int > &component) const
Create a subspace (view) for a specific component.
Definition FunctionSpace.h:111
bool contains(const FunctionSpace &V) const
Check whether V is subspace of this, or this itself.
Definition FunctionSpace.h:153
std::shared_ptr< const mesh::Mesh< geometry_type > > mesh() const
The mesh.
Definition FunctionSpace.h:361
virtual ~FunctionSpace()=default
Destructor.
std::shared_ptr< const FiniteElement< geometry_type > > element() const
The finite element.
Definition FunctionSpace.h:367
FunctionSpace(std::shared_ptr< const mesh::Mesh< geometry_type > > mesh, std::shared_ptr< const FiniteElement< geometry_type > > element, std::shared_ptr< const DofMap > dofmap)
Create function space for given mesh, element and degree-of-freedom map.
Definition FunctionSpace.h:44
T geometry_type
Geometry type of the Mesh that the FunctionSpace is defined on.
Definition FunctionSpace.h:37
std::vector< geometry_type > tabulate_dof_coordinates(bool transpose) const
Tabulate the physical coordinates of all dofs on this process.
Definition FunctionSpace.h:229
std::shared_ptr< const DofMap > dofmaps(int cell_type_idx) const
The dofmaps.
Definition FunctionSpace.h:398
FunctionSpace(FunctionSpace &&V)=default
Move constructor.
std::pair< FunctionSpace, std::vector< std::vector< std::int32_t > > > collapse() const
Definition FunctionSpace.h:179
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition Mesh.h:23
Finite element method functionality.
Definition assemble_expression_impl.h:23
@ transpose
Transpose.
Definition FiniteElement.h:28
@ standard
Standard.
Definition FiniteElement.h:27
std::array< std::vector< std::shared_ptr< const FunctionSpace< T > > >, 2 > common_function_spaces(const std::vector< std::vector< std::array< std::shared_ptr< const FunctionSpace< T > >, 2 > > > &V)
Extract FunctionSpaces for (0) rows blocks and (1) columns blocks from a rectangular array of (test,...
Definition FunctionSpace.h:433
FunctionSpace(U mesh, V element, W dofmap) -> FunctionSpace< typename std::remove_cvref< typename U::element_type >::type::geometry_type::value_type >
Type deduction.
Mesh data structures and algorithms on meshes.
Definition DofMap.h:32