13#include "FunctionSpace.h"
16#include <dolfinx/common/types.h>
49std::vector<std::int32_t>
51 std::span<const std::int32_t> entities,
77 mesh::Topology& topology,
78 std::array<std::reference_wrapper<const DofMap>, 2> dofmaps,
int dim,
79 std::span<const std::int32_t> entities,
bool remote =
true);
92template <std::
floating_po
int T,
typename U>
103 throw std::runtime_error(
104 "Cannot locate dofs geometrically for mixed space. Use subspaces.");
110 using cmdspan3x_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
112 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
113 std::size_t, 3, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>;
116 cmdspan3x_t x(dof_coordinates.data(), 3, dof_coordinates.size() / 3);
117 const std::vector<std::int8_t> marked_dofs = marker_fn(x);
119 std::vector<std::int32_t> dofs;
120 dofs.reserve(std::count(marked_dofs.begin(), marked_dofs.end(),
true));
121 for (std::size_t i = 0; i < marked_dofs.size(); ++i)
143template <std::
floating_po
int T,
typename U>
157 auto mesh = V0.
mesh();
160 if (mesh != V1.
mesh())
161 throw std::runtime_error(
"Meshes are not the same.");
162 const int tdim = mesh->topology()->dim();
167 throw std::runtime_error(
"Function spaces must have the same element.");
172 using cmdspan3x_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
174 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
175 std::size_t, 3, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>;
178 cmdspan3x_t x(dof_coordinates.data(), 3, dof_coordinates.size() / 3);
179 const std::vector<std::int8_t> marked_dofs = marker_fn(x);
182 std::shared_ptr<const DofMap> dofmap0 = V0.
dofmap();
184 const int bs0 = dofmap0->bs();
185 std::shared_ptr<const DofMap> dofmap1 = V1.
dofmap();
187 const int bs1 = dofmap1->bs();
189 const int element_bs = dofmap0->element_dof_layout().block_size();
190 assert(element_bs == dofmap1->element_dof_layout().block_size());
193 auto topology = mesh->topology();
195 std::vector<std::array<std::int32_t, 2>> bc_dofs;
196 for (
int c = 0; c < topology->connectivity(tdim, 0)->num_nodes(); ++c)
199 auto cell_dofs0 = dofmap0->cell_dofs(c);
200 auto cell_dofs1 = dofmap1->cell_dofs(c);
203 for (std::size_t i = 0; i < cell_dofs1.size(); ++i)
205 if (marked_dofs[cell_dofs1[i]])
208 for (
int k = 0; k < element_bs; ++k)
210 const int local_pos = element_bs * i + k;
211 const std::div_t pos0 = std::div(local_pos, bs0);
212 const std::div_t pos1 = std::div(local_pos, bs1);
213 const std::int32_t dof_index0
214 = bs0 * cell_dofs0[pos0.quot] + pos0.rem;
215 const std::int32_t dof_index1
216 = bs1 * cell_dofs1[pos1.quot] + pos1.rem;
217 bc_dofs.push_back({dof_index0, dof_index1});
224 std::sort(bc_dofs.begin(), bc_dofs.end());
225 bc_dofs.erase(std::unique(bc_dofs.begin(), bc_dofs.end()), bc_dofs.end());
228 std::array dofs = {std::vector<std::int32_t>(bc_dofs.size()),
229 std::vector<std::int32_t>(bc_dofs.size())};
230 std::transform(bc_dofs.cbegin(), bc_dofs.cend(), dofs[0].begin(),
231 [](
auto dof) { return dof[0]; });
232 std::transform(bc_dofs.cbegin(), bc_dofs.cend(), dofs[1].begin(),
233 [](
auto dof) { return dof[1]; });
249 std::floating_point U = dolfinx::scalar_value_type_t<T>>
255 std::size_t num_owned(
const DofMap& dofmap,
256 std::span<const std::int32_t> dofs)
259 std::int32_t map_size = dofmap.
index_map->size_local();
260 std::int32_t owned_size = bs * map_size;
261 auto it = std::lower_bound(dofs.begin(), dofs.end(), owned_size);
262 return std::distance(dofs.begin(), it);
266 static std::vector<std::int32_t>
267 unroll_dofs(std::span<const std::int32_t> dofs,
int bs)
269 std::vector<std::int32_t> dofs_unrolled(bs * dofs.size());
270 for (std::size_t i = 0; i < dofs.size(); ++i)
271 for (
int k = 0; k < bs; ++k)
272 dofs_unrolled[bs * i + k] = bs * dofs[i] + k;
273 return dofs_unrolled;
294 template <
typename S,
typename X,
296 = std::enable_if_t<std::is_convertible_v<S, T>
297 or std::is_convertible_v<S, std::span<const T>>>>
298 requires std::is_convertible_v<std::remove_cvref_t<X>,
299 std::vector<std::int32_t>>
320 template <
typename X>
321 requires std::is_convertible_v<std::remove_cvref_t<X>,
322 std::vector<std::int32_t>>
325 : _function_space(V), _g(g), _dofs0(std::forward<X>(dofs)),
326 _owned_indices0(num_owned(*V->dofmap(), _dofs0))
330 if (g->shape.size() != V->element()->value_shape().size())
332 throw std::runtime_error(
333 "Rank mis-match between Constant and function space in DirichletBC");
336 if (g->value.size() != _function_space->dofmap()->bs())
338 throw std::runtime_error(
339 "Creating a DirichletBC using a Constant is not supported when the "
340 "Constant size is not equal to the block size of the constrained "
341 "(sub-)space. Use a fem::Function to create the fem::DirichletBC.");
344 if (!V->element()->interpolation_ident())
346 throw std::runtime_error(
347 "Constant can be used only with point-evaluation elements");
351 if (
const int bs = V->dofmap()->bs(); bs > 1)
353 _owned_indices0 *= bs;
354 _dofs0 = unroll_dofs(_dofs0, bs);
370 template <
typename X>
371 requires std::is_convertible_v<std::remove_cvref_t<X>,
372 std::vector<std::int32_t>>
375 _dofs0(std::forward<X>(dofs)),
376 _owned_indices0(num_owned(*_function_space->dofmap(), _dofs0))
378 assert(_function_space);
381 if (
const int bs = _function_space->dofmap()->bs(); bs > 1)
383 _owned_indices0 *= bs;
384 _dofs0 = unroll_dofs(_dofs0, bs);
409 template <
typename X>
412 : _function_space(V), _g(g),
413 _dofs0(std::forward<typename X::value_type>(V_g_dofs[0])),
414 _dofs1_g(std::forward<typename X::value_type>(V_g_dofs[1])),
415 _owned_indices0(num_owned(*_function_space->dofmap(), _dofs0))
441 return _function_space;
446 std::variant<std::shared_ptr<const Function<T, U>>,
447 std::shared_ptr<const Constant<T>>>
459 std::pair<std::span<const std::int32_t>, std::int32_t>
dof_indices()
const
461 return {_dofs0, _owned_indices0};
475 void set(std::span<T> x, T scale = 1)
const
477 if (std::holds_alternative<std::shared_ptr<
const Function<T, U>>>(_g))
479 auto g = std::get<std::shared_ptr<const Function<T, U>>>(_g);
481 std::span<const T> values = g->x()->array();
482 auto dofs1_g = _dofs1_g.empty() ? std::span(_dofs0) : std::span(_dofs1_g);
483 std::int32_t x_size = x.size();
484 for (std::size_t i = 0; i < _dofs0.size(); ++i)
486 if (_dofs0[i] < x_size)
488 assert(dofs1_g[i] < (std::int32_t)values.size());
489 x[_dofs0[i]] = scale * values[dofs1_g[i]];
493 else if (std::holds_alternative<std::shared_ptr<
const Constant<T>>>(_g))
495 auto g = std::get<std::shared_ptr<const Constant<T>>>(_g);
496 std::vector<T>
value = g->value;
497 int bs = _function_space->dofmap()->bs();
498 std::int32_t x_size = x.size();
499 std::for_each(_dofs0.cbegin(), _dofs0.cend(),
500 [x_size, bs, scale, &
value, &x](
auto dof)
503 x[dof] = scale * value[dof % bs];
512 void set(std::span<T> x, std::span<const T> x0, T scale = 1)
const
514 if (std::holds_alternative<std::shared_ptr<
const Function<T, U>>>(_g))
516 auto g = std::get<std::shared_ptr<const Function<T, U>>>(_g);
518 std::span<const T> values = g->x()->array();
519 assert(x.size() <= x0.size());
520 auto dofs1_g = _dofs1_g.empty() ? std::span(_dofs0) : std::span(_dofs1_g);
521 std::int32_t x_size = x.size();
522 for (std::size_t i = 0; i < _dofs0.size(); ++i)
524 if (_dofs0[i] < x_size)
526 assert(dofs1_g[i] < (std::int32_t)values.size());
527 x[_dofs0[i]] = scale * (values[dofs1_g[i]] - x0[_dofs0[i]]);
531 else if (std::holds_alternative<std::shared_ptr<
const Constant<T>>>(_g))
533 auto g = std::get<std::shared_ptr<const Constant<T>>>(_g);
534 const std::vector<T>&
value = g->value;
535 std::int32_t bs = _function_space->dofmap()->bs();
536 std::for_each(_dofs0.begin(), _dofs0.end(),
537 [&x, &x0, &
value, scale, bs](
auto dof)
539 if (dof < (std::int32_t)x.size())
540 x[dof] = scale * (value[dof % bs] - x0[dof]);
555 if (std::holds_alternative<std::shared_ptr<
const Function<T, U>>>(_g))
557 auto g = std::get<std::shared_ptr<const Function<T, U>>>(_g);
559 std::span<const T> g_values = g->x()->array();
560 auto dofs1_g = _dofs1_g.empty() ? std::span(_dofs0) : std::span(_dofs1_g);
561 for (std::size_t i = 0; i < dofs1_g.size(); ++i)
562 values[_dofs0[i]] = g_values[dofs1_g[i]];
564 else if (std::holds_alternative<std::shared_ptr<
const Constant<T>>>(_g))
566 auto g = std::get<std::shared_ptr<const Constant<T>>>(_g);
568 const std::vector<T>& g_value = g->value;
569 const std::int32_t bs = _function_space->dofmap()->bs();
570 for (std::size_t i = 0; i < _dofs0.size(); ++i)
571 values[_dofs0[i]] = g_value[_dofs0[i] % bs];
583 for (std::size_t i = 0; i < _dofs0.size(); ++i)
585 assert(_dofs0[i] < (std::int32_t)markers.size());
586 markers[_dofs0[i]] =
true;
592 std::shared_ptr<const FunctionSpace<U>> _function_space;
595 std::variant<std::shared_ptr<const Function<T, U>>,
596 std::shared_ptr<const Constant<T>>>
601 std::vector<std::int32_t> _dofs0, _dofs1_g;
604 std::int32_t _owned_indices0 = -1;
Degree-of-freedom map representations and tools.
Constant value which can be attached to a Form.
Definition Constant.h:23
Object for setting (strong) Dirichlet boundary conditions.
Definition DirichletBC.h:251
void dof_values(std::span< T > values) const
Definition DirichletBC.h:553
void set(std::span< T > x, T scale=1) const
Set bc entries in x to scale * x_bc
Definition DirichletBC.h:475
DirichletBC & operator=(const DirichletBC &bc)=default
Assignment operator.
DirichletBC(std::shared_ptr< const Constant< T > > g, X &&dofs, std::shared_ptr< const FunctionSpace< U > > V)
Create a representation of a Dirichlet boundary condition constrained by a fem::Constant.
Definition DirichletBC.h:323
std::variant< std::shared_ptr< const Function< T, U > >, std::shared_ptr< const Constant< T > > > value() const
Return boundary value function g.
Definition DirichletBC.h:448
DirichletBC(const DirichletBC &bc)=default
Copy constructor.
DirichletBC(DirichletBC &&bc)=default
Move constructor.
std::pair< std::span< const std::int32_t >, std::int32_t > dof_indices() const
Access dof indices (local indices, unrolled), including ghosts, to which a Dirichlet condition is app...
Definition DirichletBC.h:459
std::shared_ptr< const FunctionSpace< U > > function_space() const
The function space to which boundary conditions are applied.
Definition DirichletBC.h:439
void set(std::span< T > x, std::span< const T > x0, T scale=1) const
Set bc entries in x to scale * (x0 - x_bc)
Definition DirichletBC.h:512
~DirichletBC()=default
Destructor.
DirichletBC(std::shared_ptr< const Function< T, U > > g, X &&dofs)
Create a representation of a Dirichlet boundary condition where the space being constrained is the sa...
Definition DirichletBC.h:373
DirichletBC & operator=(DirichletBC &&bc)=default
Move assignment operator.
DirichletBC(const S &g, X &&dofs, std::shared_ptr< const FunctionSpace< U > > V)
Create a representation of a Dirichlet boundary condition constrained by a scalar- or vector-valued c...
Definition DirichletBC.h:300
void mark_dofs(std::span< std::int8_t > markers) const
Set markers[i] = true if dof i has a boundary condition applied. Value of markers[i] is not changed o...
Definition DirichletBC.h:581
DirichletBC(std::shared_ptr< const Function< T, U > > g, X &&V_g_dofs, std::shared_ptr< const FunctionSpace< U > > V)
Create a representation of a Dirichlet boundary condition where the space being constrained and the f...
Definition DirichletBC.h:410
Degree-of-freedom map.
Definition DofMap.h:76
std::shared_ptr< const common::IndexMap > index_map
Index map that describes the parallel distribution of the dofmap.
Definition DofMap.h:171
int index_map_bs() const
Block size associated with the index_map.
Definition DofMap.cpp:298
This class represents a finite element function space defined by a mesh, a finite element,...
Definition FunctionSpace.h:35
std::shared_ptr< const DofMap > dofmap() const
The dofmap.
Definition FunctionSpace.h:309
std::shared_ptr< const FiniteElement< T > > element() const
The finite element.
Definition FunctionSpace.h:306
std::vector< T > tabulate_dof_coordinates(bool transpose) const
Tabulate the physical coordinates of all dofs on this process.
Definition FunctionSpace.h:172
std::shared_ptr< const mesh::Mesh< T > > mesh() const
The mesh.
Definition FunctionSpace.h:303
This class represents a function in a finite element function space , given by.
Definition Function.h:45
This concept is used to constrain the a template type to floating point real or complex types....
Definition types.h:20
Finite element method functionality.
Definition assemble_matrix_impl.h:25
std::vector< std::int32_t > locate_dofs_topological(mesh::Topology &topology, const DofMap &dofmap, int dim, std::span< const std::int32_t > entities, bool remote=true)
Find degrees-of-freedom which belong to the provided mesh entities (topological). Note that degrees-o...
Definition DirichletBC.cpp:180
std::vector< std::int32_t > locate_dofs_geometrical(const FunctionSpace< T > &V, U marker_fn)
Find degrees of freedom whose geometric coordinate is true for the provided marking function.
Definition DirichletBC.h:93