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.6.0
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 "Function.h"
12#include "FunctionSpace.h"
13#include <array>
14#include <concepts>
15#include <functional>
16#include <memory>
17#include <span>
18#include <type_traits>
19#include <utility>
20#include <variant>
21#include <vector>
22
23namespace dolfinx::fem
24{
25
47std::vector<std::int32_t>
48locate_dofs_topological(const FunctionSpace& V, int dim,
49 std::span<const std::int32_t> entities,
50 bool remote = true);
51
75std::array<std::vector<std::int32_t>, 2> locate_dofs_topological(
76 const std::array<std::reference_wrapper<const FunctionSpace>, 2>& V,
77 int dim, std::span<const std::int32_t> entities, bool remote = true);
78
90std::vector<std::int32_t> locate_dofs_geometrical(
91 const FunctionSpace& V,
92 const std::function<std::vector<std::int8_t>(
93 std::experimental::mdspan<
94 const double,
95 std::experimental::extents<std::size_t, 3,
96 std::experimental::dynamic_extent>>)>&
97 marker_fn);
98
112std::array<std::vector<std::int32_t>, 2> locate_dofs_geometrical(
113 const std::array<std::reference_wrapper<const FunctionSpace>, 2>& V,
114 const std::function<std::vector<std::int8_t>(
115 std::experimental::mdspan<
116 const double,
117 std::experimental::extents<std::size_t, 3,
118 std::experimental::dynamic_extent>>)>&
119 marker_fn);
120
131template <typename T>
133{
134private:
137 std::size_t num_owned(const FunctionSpace& V,
138 std::span<const std::int32_t> dofs)
139 {
140 int bs = V.dofmap()->index_map_bs();
141 std::int32_t map_size = V.dofmap()->index_map->size_local();
142 std::int32_t owned_size = bs * map_size;
143 auto it = std::lower_bound(dofs.begin(), dofs.end(), owned_size);
144 return std::distance(dofs.begin(), it);
145 }
146
148 std::vector<std::int32_t> unroll_dofs(std::span<const std::int32_t> dofs,
149 int bs)
150 {
151 std::vector<std::int32_t> dofs_unrolled(bs * dofs.size());
152 for (std::size_t i = 0; i < dofs.size(); ++i)
153 for (int k = 0; k < bs; ++k)
154 dofs_unrolled[bs * i + k] = bs * dofs[i] + k;
155 return dofs_unrolled;
156 }
157
158public:
176 template <typename S, std::convertible_to<std::vector<std::int32_t>> U,
177 typename
178 = std::enable_if_t<std::is_convertible_v<S, T>
179 or std::is_convertible_v<S, std::span<const T>>>>
180 DirichletBC(const S& g, U&& dofs, std::shared_ptr<const FunctionSpace> V)
181 : DirichletBC(std::make_shared<Constant<T>>(g), dofs, V)
182 {
183 }
184
200 template <std::convertible_to<std::vector<std::int32_t>> U>
201 DirichletBC(std::shared_ptr<const Constant<T>> g, U&& dofs,
202 std::shared_ptr<const FunctionSpace> V)
203 : _function_space(V), _g(g), _dofs0(std::forward<U>(dofs)),
204 _owned_indices0(num_owned(*V, _dofs0))
205 {
206 assert(g);
207 assert(V);
208 if (g->shape.size() != V->element()->value_shape().size())
209 {
210 throw std::runtime_error(
211 "Rank mis-match between Constant and function space in DirichletBC");
212 }
213
214 if (g->value.size() != _function_space->dofmap()->bs())
215 {
216 throw std::runtime_error(
217 "Creating a DirichletBC using a Constant is not supported when the "
218 "Constant size is not equal to the block size of the constrained "
219 "(sub-)space. Use a fem::Function to create the fem::DirichletBC.");
220 }
221
222 if (!V->element()->interpolation_ident())
223 {
224 throw std::runtime_error(
225 "Constant can be used only with point-evaluation elements");
226 }
227
228 // Unroll _dofs0 if dofmap block size > 1
229 if (const int bs = V->dofmap()->bs(); bs > 1)
230 {
231 _owned_indices0 *= bs;
232 _dofs0 = unroll_dofs(_dofs0, bs);
233 }
234 }
235
248 template <std::convertible_to<std::vector<std::int32_t>> U>
249 DirichletBC(std::shared_ptr<const Function<T>> g, U&& dofs)
250 : _function_space(g->function_space()), _g(g),
251 _dofs0(std::forward<U>(dofs)),
252 _owned_indices0(num_owned(*_function_space, _dofs0))
253 {
254 assert(_function_space);
255
256 // Unroll _dofs0 if dofmap block size > 1
257 if (const int bs = _function_space->dofmap()->bs(); bs > 1)
258 {
259 _owned_indices0 *= bs;
260 _dofs0 = unroll_dofs(_dofs0, bs);
261 }
262 }
263
285 template <typename U>
286 DirichletBC(std::shared_ptr<const Function<T>> g, U&& V_g_dofs,
287 std::shared_ptr<const FunctionSpace> V)
288 : _function_space(V), _g(g),
289 _dofs0(std::forward<typename U::value_type>(V_g_dofs[0])),
290 _dofs1_g(std::forward<typename U::value_type>(V_g_dofs[1])),
291 _owned_indices0(num_owned(*_function_space, _dofs0))
292 {
293 }
294
297 DirichletBC(const DirichletBC& bc) = default;
298
301 DirichletBC(DirichletBC&& bc) = default;
302
304 ~DirichletBC() = default;
305
308 DirichletBC& operator=(const DirichletBC& bc) = default;
309
312
315 std::shared_ptr<const FunctionSpace> function_space() const
316 {
317 return _function_space;
318 }
319
322 std::variant<std::shared_ptr<const Function<T>>,
323 std::shared_ptr<const Constant<T>>>
324 value() const
325 {
326 return _g;
327 }
328
335 std::pair<std::span<const std::int32_t>, std::int32_t> dof_indices() const
336 {
337 return {_dofs0, _owned_indices0};
338 }
339
351 void set(std::span<T> x, double scale = 1.0) const
352 {
353 if (std::holds_alternative<std::shared_ptr<const Function<T>>>(_g))
354 {
355 auto g = std::get<std::shared_ptr<const Function<T>>>(_g);
356 assert(g);
357 std::span<const T> values = g->x()->array();
358 auto dofs1_g = _dofs1_g.empty() ? std::span(_dofs0) : std::span(_dofs1_g);
359 std::int32_t x_size = x.size();
360 for (std::size_t i = 0; i < _dofs0.size(); ++i)
361 {
362 if (_dofs0[i] < x_size)
363 {
364 assert(dofs1_g[i] < (std::int32_t)values.size());
365 x[_dofs0[i]] = scale * values[dofs1_g[i]];
366 }
367 }
368 }
369 else if (std::holds_alternative<std::shared_ptr<const Constant<T>>>(_g))
370 {
371 auto g = std::get<std::shared_ptr<const Constant<T>>>(_g);
372 std::vector<T> value = g->value;
373 int bs = _function_space->dofmap()->bs();
374 std::int32_t x_size = x.size();
375 std::for_each(_dofs0.cbegin(), _dofs0.cend(),
376 [x_size, bs, scale, &value, &x](auto dof)
377 {
378 if (dof < x_size)
379 x[dof] = scale * value[dof % bs];
380 });
381 }
382 }
383
388 void set(std::span<T> x, std::span<const T> x0, double scale = 1.0) const
389 {
390 if (std::holds_alternative<std::shared_ptr<const Function<T>>>(_g))
391 {
392 auto g = std::get<std::shared_ptr<const Function<T>>>(_g);
393 assert(g);
394 std::span<const T> values = g->x()->array();
395 assert(x.size() <= x0.size());
396 auto dofs1_g = _dofs1_g.empty() ? std::span(_dofs0) : std::span(_dofs1_g);
397 std::int32_t x_size = x.size();
398 for (std::size_t i = 0; i < _dofs0.size(); ++i)
399 {
400 if (_dofs0[i] < x_size)
401 {
402 assert(dofs1_g[i] < (std::int32_t)values.size());
403 x[_dofs0[i]] = scale * (values[dofs1_g[i]] - x0[_dofs0[i]]);
404 }
405 }
406 }
407 else if (std::holds_alternative<std::shared_ptr<const Constant<T>>>(_g))
408 {
409 auto g = std::get<std::shared_ptr<const Constant<T>>>(_g);
410 const std::vector<T>& value = g->value;
411 std::int32_t bs = _function_space->dofmap()->bs();
412 std::for_each(_dofs0.begin(), _dofs0.end(),
413 [&x, &x0, &value, scale, bs](auto dof)
414 {
415 if (dof < (std::int32_t)x.size())
416 x[dof] = scale * (value[dof % bs] - x0[dof]);
417 });
418 }
419 }
420
429 void dof_values(std::span<T> values) const
430 {
431 if (std::holds_alternative<std::shared_ptr<const Function<T>>>(_g))
432 {
433 auto g = std::get<std::shared_ptr<const Function<T>>>(_g);
434 assert(g);
435 std::span<const T> g_values = g->x()->array();
436 auto dofs1_g = _dofs1_g.empty() ? std::span(_dofs0) : std::span(_dofs1_g);
437 for (std::size_t i = 0; i < dofs1_g.size(); ++i)
438 values[_dofs0[i]] = g_values[dofs1_g[i]];
439 }
440 else if (std::holds_alternative<std::shared_ptr<const Constant<T>>>(_g))
441 {
442 auto g = std::get<std::shared_ptr<const Constant<T>>>(_g);
443 assert(g);
444 const std::vector<T>& g_value = g->value;
445 const std::int32_t bs = _function_space->dofmap()->bs();
446 for (std::size_t i = 0; i < _dofs0.size(); ++i)
447 values[_dofs0[i]] = g_value[_dofs0[i] % bs];
448 }
449 }
450
457 void mark_dofs(std::span<std::int8_t> markers) const
458 {
459 for (std::size_t i = 0; i < _dofs0.size(); ++i)
460 {
461 assert(_dofs0[i] < (std::int32_t)markers.size());
462 markers[_dofs0[i]] = true;
463 }
464 }
465
466private:
467 // The function space (possibly a sub function space)
468 std::shared_ptr<const FunctionSpace> _function_space;
469
470 // The function
471 std::variant<std::shared_ptr<const Function<T>>,
472 std::shared_ptr<const Constant<T>>>
473 _g;
474
475 // Dof indices (_dofs0) in _function_space and (_dofs1_g) in the
476 // space of _g. _dofs1_g may be empty if _dofs0 can be re-used
477 std::vector<std::int32_t> _dofs0, _dofs1_g;
478
479 // The first _owned_indices in _dofs are owned by this process
480 std::int32_t _owned_indices0 = -1;
481};
482} // namespace dolfinx::fem
Constant value which can be attached to a Form. Constants may be scalar (rank 0), vector (rank 1),...
Definition: Constant.h:20
Object for setting (strong) Dirichlet boundary conditions.
Definition: DirichletBC.h:133
void dof_values(std::span< T > values) const
Definition: DirichletBC.h:429
void set(std::span< T > x, double scale=1.0) const
Set bc entries in x to scale * x_bc
Definition: DirichletBC.h:351
DirichletBC & operator=(const DirichletBC &bc)=default
Assignment operator.
DirichletBC(const DirichletBC &bc)=default
Copy constructor.
DirichletBC(std::shared_ptr< const Constant< T > > g, U &&dofs, std::shared_ptr< const FunctionSpace > V)
Create a representation of a Dirichlet boundary condition constrained by a fem::Constant.
Definition: DirichletBC.h:201
DirichletBC(DirichletBC &&bc)=default
Move constructor.
DirichletBC(std::shared_ptr< const Function< T > > g, U &&dofs)
Create a representation of a Dirichlet boundary condition where the space being constrained is the sa...
Definition: DirichletBC.h:249
std::variant< std::shared_ptr< const Function< T > >, std::shared_ptr< const Constant< T > > > value() const
Return boundary value function g.
Definition: DirichletBC.h:324
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:335
std::shared_ptr< const FunctionSpace > function_space() const
The function space to which boundary conditions are applied.
Definition: DirichletBC.h:315
DirichletBC(std::shared_ptr< const Function< T > > g, U &&V_g_dofs, std::shared_ptr< const FunctionSpace > V)
Create a representation of a Dirichlet boundary condition where the space being constrained and the f...
Definition: DirichletBC.h:286
~DirichletBC()=default
Destructor.
void set(std::span< T > x, std::span< const T > x0, double scale=1.0) const
Set bc entries in x to scale * (x0 - x_bc)
Definition: DirichletBC.h:388
DirichletBC & operator=(DirichletBC &&bc)=default
Move assignment operator.
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:457
DirichletBC(const S &g, U &&dofs, std::shared_ptr< const FunctionSpace > V)
Create a representation of a Dirichlet boundary condition constrained by a scalar- or vector-valued c...
Definition: DirichletBC.h:180
This class represents a finite element function space defined by a mesh, a finite element,...
Definition: FunctionSpace.h:31
std::shared_ptr< const DofMap > dofmap() const
The dofmap.
Definition: FunctionSpace.cpp:253
This class represents a function in a finite element function space , given by.
Definition: Function.h:43
Finite element method functionality.
Definition: assemble_matrix_impl.h:25
std::vector< std::int32_t > locate_dofs_topological(const FunctionSpace &V, 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:185
std::vector< std::int32_t > locate_dofs_geometrical(const FunctionSpace &V, const std::function< std::vector< std::int8_t >(std::experimental::mdspan< const double, std::experimental::extents< std::size_t, 3, std::experimental::dynamic_extent > >)> &marker_fn)
Finds degrees of freedom whose geometric coordinate is true for the provided marking function.
Definition: DirichletBC.cpp:455