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)
47 : _mesh(
mesh), _element(element), _dofmap(
dofmap),
48 _id(boost::uuids::random_generator()()), _root_space_id(_id)
84 throw std::runtime_error(
"Component must be non-empty");
87 auto element = this->_element->extract_sub_element(
component);
91 = std::make_shared<DofMap>(_dofmap->extract_sub_dofmap(
component));
97 sub_space._root_space_id = _root_space_id;
98 sub_space._component = _component;
99 sub_space._component.insert(sub_space._component.end(),
component.begin(),
110 if (
this == std::addressof(V))
112 else if (_root_space_id != V._root_space_id)
114 else if (_component.size()
115 > V._component.size())
119 else if (!std::equal(_component.begin(), _component.end(),
120 V._component.begin()))
133 std::pair<FunctionSpace, std::vector<std::int32_t>>
collapse()
const
135 if (_component.empty())
136 throw std::runtime_error(
"Function space is not a subspace");
139 auto [_collapsed_dofmap, collapsed_dofs]
140 = _dofmap->collapse(_mesh->comm(), *_mesh->topology());
141 auto collapsed_dofmap
142 = std::make_shared<DofMap>(std::move(_collapsed_dofmap));
145 std::move(collapsed_dofs)};
151 std::vector<int>
component()
const {
return _component; }
158 return _element->symmetric();
175 if (!_component.empty())
177 throw std::runtime_error(
"Cannot tabulate coordinates for a "
178 "FunctionSpace that is a subspace.");
182 if (_element->is_mixed())
184 throw std::runtime_error(
185 "Cannot tabulate coordinates for a mixed FunctionSpace.");
191 const std::size_t gdim = _mesh->geometry().dim();
192 const int tdim = _mesh->topology()->dim();
196 std::shared_ptr<const common::IndexMap> index_map = _dofmap->index_map;
198 const int index_map_bs = _dofmap->index_map_bs();
199 const int dofmap_bs = _dofmap->bs();
201 const int element_block_size = _element->block_size();
202 const std::size_t scalar_dofs
203 = _element->space_dimension() / element_block_size;
204 const std::int32_t num_dofs
205 = index_map_bs * (index_map->size_local() + index_map->num_ghosts())
209 if (!_element->interpolation_ident())
211 throw std::runtime_error(
"Cannot evaluate dof coordinates - this element "
212 "does not have pointwise evaluation.");
214 const auto [X, Xshape] = _element->interpolation_points();
220 auto x_dofmap = _mesh->geometry().dofmap();
221 const std::size_t num_dofs_g = cmap.
dim();
222 std::span<const geometry_type> x_g = _mesh->geometry().x();
225 const std::size_t shape_c0 =
transpose ? 3 : num_dofs;
226 const std::size_t shape_c1 =
transpose ? num_dofs : 3;
227 std::vector<geometry_type> coords(shape_c0 * shape_c1, 0);
229 using mdspan2_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
231 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 2>>;
232 using cmdspan4_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
234 MDSPAN_IMPL_STANDARD_NAMESPACE::dextents<std::size_t, 4>>;
237 std::vector<geometry_type> x_b(scalar_dofs * gdim);
238 mdspan2_t x(x_b.data(), scalar_dofs, gdim);
240 std::vector<geometry_type> coordinate_dofs_b(num_dofs_g * gdim);
241 mdspan2_t coordinate_dofs(coordinate_dofs_b.data(), num_dofs_g, gdim);
243 auto map = _mesh->topology()->index_map(tdim);
245 const int num_cells = map->size_local() + map->num_ghosts();
247 std::span<const std::uint32_t> cell_info;
248 if (_element->needs_dof_transformations())
250 _mesh->topology_mutable()->create_entity_permutations();
251 cell_info = std::span(_mesh->topology()->get_cell_permutation_info());
254 const std::array<std::size_t, 4> phi_shape
256 std::vector<geometry_type> phi_b(
257 std::reduce(phi_shape.begin(), phi_shape.end(), 1, std::multiplies{}));
258 cmdspan4_t phi_full(phi_b.data(), phi_shape);
260 auto phi = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
261 phi_full, 0, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent,
262 MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent, 0);
266 auto apply_dof_transformation
267 = _element->template dof_transformation_fn<geometry_type>(
270 for (
int c = 0; c < num_cells; ++c)
273 auto x_dofs = MDSPAN_IMPL_STANDARD_NAMESPACE::submdspan(
274 x_dofmap, c, MDSPAN_IMPL_STANDARD_NAMESPACE::full_extent);
275 for (std::size_t i = 0; i < x_dofs.size(); ++i)
276 for (std::size_t j = 0; j < gdim; ++j)
277 coordinate_dofs(i, j) = x_g[3 * x_dofs[i] + j];
281 apply_dof_transformation(
282 x_b, std::span(cell_info.data(), cell_info.size()), c, x.extent(1));
285 auto dofs = _dofmap->cell_dofs(c);
290 for (std::size_t i = 0; i < dofs.size(); ++i)
291 for (std::size_t j = 0; j < gdim; ++j)
292 coords[dofs[i] * 3 + j] = x(i, j);
296 for (std::size_t i = 0; i < dofs.size(); ++i)
297 for (std::size_t j = 0; j < gdim; ++j)
298 coords[j * num_dofs + dofs[i]] = x(i, j);
306 std::shared_ptr<const mesh::Mesh<geometry_type>>
mesh()
const
312 std::shared_ptr<const FiniteElement<geometry_type>>
element()
const
318 std::shared_ptr<const DofMap>
dofmap()
const {
return _dofmap; }
322 std::shared_ptr<const mesh::Mesh<geometry_type>> _mesh;
325 std::shared_ptr<const FiniteElement<geometry_type>> _element;
328 std::shared_ptr<const DofMap> _dofmap;
331 std::vector<int> _component;
334 boost::uuids::uuid _id;
335 boost::uuids::uuid _root_space_id;
348template <dolfinx::scalar T>
349std::array<std::vector<std::shared_ptr<const FunctionSpace<T>>>, 2>
355 std::vector<std::shared_ptr<const FunctionSpace<T>>> spaces0(V.size(),
357 std::vector<std::shared_ptr<const FunctionSpace<T>>> spaces1(V.front().size(),
361 for (std::size_t i = 0; i < V.size(); ++i)
364 for (std::size_t j = 0; j < V[i].size(); ++j)
366 auto& V0 = V[i][j][0];
367 auto& V1 = V[i][j][1];
374 if (spaces0[i] != V0)
375 throw std::runtime_error(
"Mismatched test space for row.");
382 if (spaces1[j] != V1)
383 throw std::runtime_error(
"Mismatched trial space for column.");
390 if (std::find(spaces0.begin(), spaces0.end(),
nullptr) != spaces0.end())
391 throw std::runtime_error(
"Could not deduce all block test spaces.");
392 if (std::find(spaces1.begin(), spaces1.end(),
nullptr) != spaces1.end())
393 throw std::runtime_error(
"Could not deduce all block trial spaces.");
395 return {spaces0, spaces1};
399template <
typename U,
typename V,
typename W>
402 typename U::element_type>::type::geometry_type::value_type>;
Degree-of-freedom map representations and tools.
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:206
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:188
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 vtk_utils.h:32
std::vector< int > component() const
Get the component with respect to the root superspace.
Definition FunctionSpace.h:151
bool symmetric() const
Indicate whether this function space represents a symmetric 2-tensor.
Definition FunctionSpace.h:155
FunctionSpace & operator=(FunctionSpace &&V)=default
Move assignment operator.
std::pair< FunctionSpace, std::vector< std::int32_t > > collapse() const
Definition FunctionSpace.h:133
std::shared_ptr< const DofMap > dofmap() const
The dofmap.
Definition FunctionSpace.h:318
FunctionSpace sub(const std::vector< int > &component) const
Create a subspace (view) for a specific component.
Definition FunctionSpace.h:76
bool contains(const FunctionSpace &V) const
Check whether V is subspace of this, or this itself.
Definition FunctionSpace.h:108
std::shared_ptr< const mesh::Mesh< geometry_type > > mesh() const
The mesh.
Definition FunctionSpace.h:306
virtual ~FunctionSpace()=default
Destructor.
std::shared_ptr< const FiniteElement< geometry_type > > element() const
The finite element.
Definition FunctionSpace.h:312
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:173
FunctionSpace(FunctionSpace &&V)=default
Move constructor.
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_matrix_impl.h:26
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:350
FunctionSpace(U mesh, V element, W dofmap) -> FunctionSpace< typename std::remove_cvref< typename U::element_type >::type::geometry_type::value_type >
Type deduction.