Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/dolfinx/v0.9.0/cpp/doxygen/dd/d4f/DirichletBC_8h_source.html
DOLFINx 0.7.3
DOLFINx C++ interface
Loading...
Searching...
No Matches
DirichletBC.h
1// Copyright (C) 2007-2021 Michal Habera, Anders Logg, Garth N. Wells
2// and Jørgen S.Dokken
3//
4// This file is part of DOLFINx (https://www.fenicsproject.org)
5//
6// SPDX-License-Identifier: LGPL-3.0-or-later
7
8#pragma once
9
10#include "Constant.h"
11#include "DofMap.h"
12#include "Function.h"
13#include "FunctionSpace.h"
14#include <array>
15#include <concepts>
16#include <dolfinx/common/types.h>
17#include <functional>
18#include <memory>
19#include <span>
20#include <type_traits>
21#include <utility>
22#include <variant>
23#include <vector>
24
25namespace dolfinx::fem
26{
27
49std::vector<std::int32_t>
50locate_dofs_topological(mesh::Topology& topology, const DofMap& dofmap, int dim,
51 std::span<const std::int32_t> entities,
52 bool remote = true);
53
76std::array<std::vector<std::int32_t>, 2> locate_dofs_topological(
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);
80
92template <std::floating_point T, typename U>
93std::vector<std::int32_t> locate_dofs_geometrical(const FunctionSpace<T>& V,
94 U marker_fn)
95{
96 // FIXME: Calling V.tabulate_dof_coordinates() is very expensive,
97 // especially when we usually want the boundary dofs only. Add
98 // interface that computes dofs coordinates only for specified cell.
99
100 assert(V.element());
101 if (V.element()->is_mixed())
102 {
103 throw std::runtime_error(
104 "Cannot locate dofs geometrically for mixed space. Use subspaces.");
105 }
106
107 // Compute dof coordinates
108 const std::vector<T> dof_coordinates = V.tabulate_dof_coordinates(true);
109
110 using cmdspan3x_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
111 const T,
112 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
113 std::size_t, 3, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>;
114
115 // Compute marker for each dof coordinate
116 cmdspan3x_t x(dof_coordinates.data(), 3, dof_coordinates.size() / 3);
117 const std::vector<std::int8_t> marked_dofs = marker_fn(x);
118
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)
122 {
123 if (marked_dofs[i])
124 dofs.push_back(i);
125 }
126
127 return dofs;
128}
129
143template <std::floating_point T, typename U>
144std::array<std::vector<std::int32_t>, 2> locate_dofs_geometrical(
145 const std::array<std::reference_wrapper<const FunctionSpace<T>>, 2>& V,
146 U marker_fn)
147{
148 // FIXME: Calling V.tabulate_dof_coordinates() is very expensive,
149 // especially when we usually want the boundary dofs only. Add
150 // interface that computes dofs coordinates only for specified cell.
151
152 // Get function spaces
153 const FunctionSpace<T>& V0 = V.at(0).get();
154 const FunctionSpace<T>& V1 = V.at(1).get();
155
156 // Get mesh
157 auto mesh = V0.mesh();
158 assert(mesh);
159 assert(V1.mesh());
160 if (mesh != V1.mesh())
161 throw std::runtime_error("Meshes are not the same.");
162 const int tdim = mesh->topology()->dim();
163
164 assert(V0.element());
165 assert(V1.element());
166 if (*V0.element() != *V1.element())
167 throw std::runtime_error("Function spaces must have the same element.");
168
169 // Compute dof coordinates
170 const std::vector<T> dof_coordinates = V1.tabulate_dof_coordinates(true);
171
172 using cmdspan3x_t = MDSPAN_IMPL_STANDARD_NAMESPACE::mdspan<
173 const T,
174 MDSPAN_IMPL_STANDARD_NAMESPACE::extents<
175 std::size_t, 3, MDSPAN_IMPL_STANDARD_NAMESPACE::dynamic_extent>>;
176
177 // Evaluate marker for each dof coordinate
178 cmdspan3x_t x(dof_coordinates.data(), 3, dof_coordinates.size() / 3);
179 const std::vector<std::int8_t> marked_dofs = marker_fn(x);
180
181 // Get dofmaps
182 std::shared_ptr<const DofMap> dofmap0 = V0.dofmap();
183 assert(dofmap0);
184 const int bs0 = dofmap0->bs();
185 std::shared_ptr<const DofMap> dofmap1 = V1.dofmap();
186 assert(dofmap1);
187 const int bs1 = dofmap1->bs();
188
189 const int element_bs = dofmap0->element_dof_layout().block_size();
190 assert(element_bs == dofmap1->element_dof_layout().block_size());
191
192 // Iterate over cells
193 auto topology = mesh->topology();
194 assert(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)
197 {
198 // Get cell dofmaps
199 auto cell_dofs0 = dofmap0->cell_dofs(c);
200 auto cell_dofs1 = dofmap1->cell_dofs(c);
201
202 // Loop over cell dofs and add to bc_dofs if marked.
203 for (std::size_t i = 0; i < cell_dofs1.size(); ++i)
204 {
205 if (marked_dofs[cell_dofs1[i]])
206 {
207 // Unroll over blocks
208 for (int k = 0; k < element_bs; ++k)
209 {
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});
218 }
219 }
220 }
221 }
222
223 // Remove duplicates
224 std::sort(bc_dofs.begin(), bc_dofs.end());
225 bc_dofs.erase(std::unique(bc_dofs.begin(), bc_dofs.end()), bc_dofs.end());
226
227 // Copy to separate array
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]; });
234
235 return dofs;
236}
237
248template <dolfinx::scalar T,
249 std::floating_point U = dolfinx::scalar_value_type_t<T>>
251{
252private:
255 std::size_t num_owned(const DofMap& dofmap,
256 std::span<const std::int32_t> dofs)
257 {
258 int bs = dofmap.index_map_bs();
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);
263 }
264
266 static std::vector<std::int32_t>
267 unroll_dofs(std::span<const std::int32_t> dofs, int bs)
268 {
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;
274 }
275
276public:
294 template <typename S, typename X,
295 typename
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>>
300 DirichletBC(const S& g, X&& dofs, std::shared_ptr<const FunctionSpace<U>> V)
301 : DirichletBC(std::make_shared<Constant<T>>(g), dofs, V)
302 {
303 }
304
320 template <typename X>
321 requires std::is_convertible_v<std::remove_cvref_t<X>,
322 std::vector<std::int32_t>>
323 DirichletBC(std::shared_ptr<const Constant<T>> g, X&& dofs,
324 std::shared_ptr<const FunctionSpace<U>> V)
325 : _function_space(V), _g(g), _dofs0(std::forward<X>(dofs)),
326 _owned_indices0(num_owned(*V->dofmap(), _dofs0))
327 {
328 assert(g);
329 assert(V);
330 if (g->shape.size() != V->element()->value_shape().size())
331 {
332 throw std::runtime_error(
333 "Rank mis-match between Constant and function space in DirichletBC");
334 }
335
336 if (g->value.size() != _function_space->dofmap()->bs())
337 {
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.");
342 }
343
344 if (!V->element()->interpolation_ident())
345 {
346 throw std::runtime_error(
347 "Constant can be used only with point-evaluation elements");
348 }
349
350 // Unroll _dofs0 if dofmap block size > 1
351 if (const int bs = V->dofmap()->bs(); bs > 1)
352 {
353 _owned_indices0 *= bs;
354 _dofs0 = unroll_dofs(_dofs0, bs);
355 }
356 }
357
370 template <typename X>
371 requires std::is_convertible_v<std::remove_cvref_t<X>,
372 std::vector<std::int32_t>>
373 DirichletBC(std::shared_ptr<const Function<T, U>> g, X&& dofs)
374 : _function_space(g->function_space()), _g(g),
375 _dofs0(std::forward<X>(dofs)),
376 _owned_indices0(num_owned(*_function_space->dofmap(), _dofs0))
377 {
378 assert(_function_space);
379
380 // Unroll _dofs0 if dofmap block size > 1
381 if (const int bs = _function_space->dofmap()->bs(); bs > 1)
382 {
383 _owned_indices0 *= bs;
384 _dofs0 = unroll_dofs(_dofs0, bs);
385 }
386 }
387
409 template <typename X>
410 DirichletBC(std::shared_ptr<const Function<T, U>> g, X&& V_g_dofs,
411 std::shared_ptr<const FunctionSpace<U>> V)
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))
416 {
417 }
418
421 DirichletBC(const DirichletBC& bc) = default;
422
425 DirichletBC(DirichletBC&& bc) = default;
426
428 ~DirichletBC() = default;
429
432 DirichletBC& operator=(const DirichletBC& bc) = default;
433
436
439 std::shared_ptr<const FunctionSpace<U>> function_space() const
440 {
441 return _function_space;
442 }
443
446 std::variant<std::shared_ptr<const Function<T, U>>,
447 std::shared_ptr<const Constant<T>>>
448 value() const
449 {
450 return _g;
451 }
452
459 std::pair<std::span<const std::int32_t>, std::int32_t> dof_indices() const
460 {
461 return {_dofs0, _owned_indices0};
462 }
463
475 void set(std::span<T> x, T scale = 1) const
476 {
477 if (std::holds_alternative<std::shared_ptr<const Function<T, U>>>(_g))
478 {
479 auto g = std::get<std::shared_ptr<const Function<T, U>>>(_g);
480 assert(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)
485 {
486 if (_dofs0[i] < x_size)
487 {
488 assert(dofs1_g[i] < (std::int32_t)values.size());
489 x[_dofs0[i]] = scale * values[dofs1_g[i]];
490 }
491 }
492 }
493 else if (std::holds_alternative<std::shared_ptr<const Constant<T>>>(_g))
494 {
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)
501 {
502 if (dof < x_size)
503 x[dof] = scale * value[dof % bs];
504 });
505 }
506 }
507
512 void set(std::span<T> x, std::span<const T> x0, T scale = 1) const
513 {
514 if (std::holds_alternative<std::shared_ptr<const Function<T, U>>>(_g))
515 {
516 auto g = std::get<std::shared_ptr<const Function<T, U>>>(_g);
517 assert(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)
523 {
524 if (_dofs0[i] < x_size)
525 {
526 assert(dofs1_g[i] < (std::int32_t)values.size());
527 x[_dofs0[i]] = scale * (values[dofs1_g[i]] - x0[_dofs0[i]]);
528 }
529 }
530 }
531 else if (std::holds_alternative<std::shared_ptr<const Constant<T>>>(_g))
532 {
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)
538 {
539 if (dof < (std::int32_t)x.size())
540 x[dof] = scale * (value[dof % bs] - x0[dof]);
541 });
542 }
543 }
544
553 void dof_values(std::span<T> values) const
554 {
555 if (std::holds_alternative<std::shared_ptr<const Function<T, U>>>(_g))
556 {
557 auto g = std::get<std::shared_ptr<const Function<T, U>>>(_g);
558 assert(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]];
563 }
564 else if (std::holds_alternative<std::shared_ptr<const Constant<T>>>(_g))
565 {
566 auto g = std::get<std::shared_ptr<const Constant<T>>>(_g);
567 assert(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];
572 }
573 }
574
581 void mark_dofs(std::span<std::int8_t> markers) const
582 {
583 for (std::size_t i = 0; i < _dofs0.size(); ++i)
584 {
585 assert(_dofs0[i] < (std::int32_t)markers.size());
586 markers[_dofs0[i]] = true;
587 }
588 }
589
590private:
591 // The function space (possibly a sub function space)
592 std::shared_ptr<const FunctionSpace<U>> _function_space;
593
594 // The function
595 std::variant<std::shared_ptr<const Function<T, U>>,
596 std::shared_ptr<const Constant<T>>>
597 _g;
598
599 // Dof indices (_dofs0) in _function_space and (_dofs1_g) in the
600 // space of _g. _dofs1_g may be empty if _dofs0 can be re-used
601 std::vector<std::int32_t> _dofs0, _dofs1_g;
602
603 // The first _owned_indices in _dofs are owned by this process
604 std::int32_t _owned_indices0 = -1;
605};
606} // namespace dolfinx::fem
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