9#include "CoordinateElement.h"
11#include "FiniteElement.h"
12#include <boost/uuid/uuid.hpp>
13#include <boost/uuid/uuid_generators.hpp>
17#include <dolfinx/common/IndexMap.h>
18#include <dolfinx/mesh/Geometry.h>
19#include <dolfinx/mesh/Mesh.h>
20#include <dolfinx/mesh/Topology.h>
32template <std::
floating_po
int T>
46 std::shared_ptr<const DofMap>
dofmap)
48 _id(boost::uuids::random_generator()()), _root_space_id(_id)
64 std::vector<std::shared_ptr<const DofMap>>
dofmaps)
66 _id(boost::uuids::random_generator()()), _root_space_id(_id)
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)
72 throw std::runtime_error(
73 "Number of elements must match number of cell types");
75 if (_dofmaps.size() != num_cell_types)
77 throw std::runtime_error(
78 "Number of dofmaps must match number of cell types");
80 for (std::size_t i = 0; i < num_cell_types; ++i)
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");
114 assert(_elements.size() > 0);
115 assert(_dofmaps.size() > 0);
119 throw std::runtime_error(
"Component must be non-empty");
122 std::vector<std::shared_ptr<const FiniteElement<geometry_type>>>
124 for (
auto e : _elements)
127 sub_elements.push_back(e->extract_sub_element(
component));
130 std::vector<std::shared_ptr<const DofMap>> sub_dofmaps;
131 for (
auto d : _dofmaps)
134 sub_dofmaps.push_back(
135 std::make_shared<const DofMap>(d->extract_sub_dofmap(
component)));
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(),
155 if (
this == std::addressof(V))
157 else if (_root_space_id != V._root_space_id)
159 else if (_component.size()
160 > V._component.size())
164 else if (!std::equal(_component.begin(), _component.end(),
165 V._component.begin()))
178 std::pair<FunctionSpace, std::vector<std::vector<std::int32_t>>>
181 spdlog::debug(
"FunctionSpace::collapse");
182 if (_component.empty())
183 throw std::runtime_error(
"Function space is not a subspace");
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)
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);
201 std::move(collapsed_dofs)};
207 std::vector<int>
component()
const {
return _component; }
213 if (_elements.front())
214 return _elements.front()->symmetric();
231 if (!_component.empty())
233 throw std::runtime_error(
"Cannot tabulate coordinates for a "
234 "FunctionSpace that is a subspace.");
237 assert(_elements.front());
238 if (_elements.front()->is_mixed())
240 throw std::runtime_error(
241 "Cannot tabulate coordinates for a mixed FunctionSpace.");
246 assert(_elements.front());
247 const std::size_t gdim = _mesh->geometry().dim();
248 const int tdim = _mesh->topology()->dim();
251 assert(_dofmaps.front());
252 std::shared_ptr<const common::IndexMap> index_map
253 = _dofmaps.front()->index_map;
255 const int index_map_bs = _dofmaps.front()->index_map_bs();
256 const int dofmap_bs = _dofmaps.front()->bs();
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())
265 std::size_t num_cell_types = _elements.size();
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>>;
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);
278 for (std::size_t i = 0; i < num_cell_types; ++i)
281 if (!_elements[i]->interpolation_ident())
283 throw std::runtime_error(
284 "Cannot evaluate dof coordinates - this element "
285 "does not have pointwise evaluation.");
287 const auto [X, Xshape] = _elements[i]->interpolation_points();
291 = _mesh->geometry().cmaps()[i];
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);
299 auto map = _mesh->topology()->index_maps(tdim)[i];
301 const int num_cells = map->size_local() + map->num_ghosts();
303 std::span<const std::uint32_t> cell_info;
304 if (_elements[i]->needs_dof_transformations())
306 _mesh->topology_mutable()->create_entity_permutations();
307 cell_info = std::span(_mesh->topology()->get_cell_permutation_info());
310 const std::array<std::size_t, 4> phi_shape
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);
317 = md::submdspan(phi_full, 0, md::full_extent, md::full_extent, 0);
321 auto apply_dof_transformation
322 = _elements[i]->template dof_transformation_fn<geometry_type>(
325 for (
int c = 0; c < num_cells; ++c)
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];
335 apply_dof_transformation(
336 x_b, std::span(cell_info.data(), cell_info.size()), c, x.extent(1));
339 auto dofs = _dofmaps[i]->cell_dofs(c);
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);
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);
361 std::shared_ptr<const mesh::Mesh<geometry_type>>
mesh()
const
367 std::shared_ptr<const FiniteElement<geometry_type>>
element()
const
369 if (_elements.size() > 1)
371 throw std::runtime_error(
372 "FunctionSpace has multiple elements, call `elements` instead.");
379 std::shared_ptr<const FiniteElement<geometry_type>>
382 return _elements.at(cell_type_idx);
386 std::shared_ptr<const DofMap>
dofmap()
const
388 if (_dofmaps.size() > 1)
390 throw std::runtime_error(
391 "FunctionSpace has multiple dofmaps, call `dofmaps` instead.");
398 std::shared_ptr<const DofMap>
dofmaps(
int cell_type_idx)
const
400 return _dofmaps.at(cell_type_idx);
405 std::shared_ptr<const mesh::Mesh<geometry_type>> _mesh;
408 std::vector<std::shared_ptr<const FiniteElement<geometry_type>>> _elements;
411 std::vector<std::shared_ptr<const DofMap>> _dofmaps;
414 std::vector<int> _component;
417 boost::uuids::uuid _id;
418 boost::uuids::uuid _root_space_id;
431template <dolfinx::scalar T>
432std::array<std::vector<std::shared_ptr<const FunctionSpace<T>>>, 2>
438 std::vector<std::shared_ptr<const FunctionSpace<T>>> spaces0(V.size(),
440 std::vector<std::shared_ptr<const FunctionSpace<T>>> spaces1(V.front().size(),
444 for (std::size_t i = 0; i < V.size(); ++i)
447 for (std::size_t j = 0; j < V[i].size(); ++j)
449 auto& V0 = V[i][j][0];
450 auto& V1 = V[i][j][1];
457 if (spaces0[i] != V0)
458 throw std::runtime_error(
"Mismatched test space for row.");
465 if (spaces1[j] != V1)
466 throw std::runtime_error(
"Mismatched trial space for column.");
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.");
478 return {spaces0, spaces1};
482template <
typename U,
typename V,
typename W>
485 typename U::element_type>::type::geometry_type::value_type>;
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