Note: this is documentation for an old release. View the latest documentation at docs.fenicsproject.org/v0.3.0/v0.9.0/cpp
DOLFINx  0.3.0
DOLFINx C++ interface
DirichletBC.h
1 // Copyright (C) 2007-2020 Michal Habera, Anders Logg and Garth N. Wells
2 //
3 // This file is part of DOLFINx (https://www.fenicsproject.org)
4 //
5 // SPDX-License-Identifier: LGPL-3.0-or-later
6 
7 #pragma once
8 
9 #include <array>
10 #include <dolfinx/fem/Function.h>
11 #include <dolfinx/fem/FunctionSpace.h>
12 #include <dolfinx/la/utils.h>
13 #include <functional>
14 #include <memory>
15 #include <vector>
16 #include <xtensor/xtensor.hpp>
17 #include <xtl/xspan.hpp>
18 
19 namespace dolfinx::mesh
20 {
21 class Mesh;
22 }
23 
24 namespace dolfinx::fem
25 {
26 
50 std::array<std::vector<std::int32_t>, 2> locate_dofs_topological(
51  const std::array<std::reference_wrapper<const fem::FunctionSpace>, 2>& V,
52  const int dim, const xtl::span<const std::int32_t>& entities,
53  bool remote = true);
54 
76 std::vector<std::int32_t>
77 locate_dofs_topological(const fem::FunctionSpace& V, const int dim,
78  const xtl::span<const std::int32_t>& entities,
79  bool remote = true);
80 
94 std::array<std::vector<std::int32_t>, 2> locate_dofs_geometrical(
95  const std::array<std::reference_wrapper<const fem::FunctionSpace>, 2>& V,
96  const std::function<xt::xtensor<bool, 1>(const xt::xtensor<double, 2>&)>&
97  marker_fn);
98 
110 std::vector<std::int32_t> locate_dofs_geometrical(
111  const fem::FunctionSpace& V,
112  const std::function<xt::xtensor<bool, 1>(const xt::xtensor<double, 2>&)>&
113  marker_fn);
114 
125 
126 template <typename T>
128 {
129 
130 public:
145  template <typename U>
146  DirichletBC(const std::shared_ptr<const fem::Function<T>>& g, U&& dofs)
147  : _function_space(g->function_space()), _g(g),
148  _dofs0(std::forward<U>(dofs))
149  {
150  const int owned_size0 = _function_space->dofmap()->index_map->size_local();
151  auto it = std::lower_bound(_dofs0.begin(), _dofs0.end(), owned_size0);
152  const int map0_bs = _function_space->dofmap()->index_map_bs();
153  _owned_indices0 = map0_bs * std::distance(_dofs0.begin(), it);
154 
155  const int bs = _function_space->dofmap()->bs();
156  if (bs > 1)
157  {
158  // Unroll for the block size
159  const std::vector<std::int32_t> dof_tmp = _dofs0;
160  _dofs0.resize(bs * dof_tmp.size());
161  for (std::size_t i = 0; i < dof_tmp.size(); ++i)
162  {
163  for (int k = 0; k < bs; ++k)
164  _dofs0[bs * i + k] = bs * dof_tmp[i] + k;
165  }
166  }
167 
168  // TODO: allows single dofs array (let one point to the other)
169  _dofs1_g = _dofs0;
170  }
171 
189  DirichletBC(const std::shared_ptr<const fem::Function<T>>& g,
190  const std::array<std::vector<std::int32_t>, 2>& V_g_dofs,
191  std::shared_ptr<const fem::FunctionSpace> V)
192  : _function_space(V), _g(g), _dofs0(V_g_dofs[0]), _dofs1_g(V_g_dofs[1])
193  {
194  assert(_dofs0.size() == _dofs1_g.size());
195  assert(_function_space);
196  assert(_g);
197 
198  const int map0_bs = _function_space->dofmap()->index_map_bs();
199  const int map0_size = _function_space->dofmap()->index_map->size_local();
200  const int owned_size0 = map0_bs * map0_size;
201  auto it0 = std::lower_bound(_dofs0.begin(), _dofs0.end(), owned_size0);
202  _owned_indices0 = std::distance(_dofs0.begin(), it0);
203  }
204 
207  DirichletBC(const DirichletBC& bc) = default;
208 
211  DirichletBC(DirichletBC&& bc) = default;
212 
214  ~DirichletBC() = default;
215 
218  DirichletBC& operator=(const DirichletBC& bc) = default;
219 
222 
225  std::shared_ptr<const fem::FunctionSpace> function_space() const
226  {
227  return _function_space;
228  }
229 
232  std::shared_ptr<const fem::Function<T>> value() const { return _g; }
233 
240  std::pair<xtl::span<const std::int32_t>, std::int32_t> dof_indices() const
241  {
242  return {_dofs0, _owned_indices0};
243  }
244 
256  void set(xtl::span<T> x, double scale = 1.0) const
257  {
258  assert(_g);
259  const std::vector<T>& g = _g->x()->array();
260  for (std::size_t i = 0; i < _dofs0.size(); ++i)
261  {
262  if (_dofs0[i] < (std::int32_t)x.size())
263  {
264  assert(_dofs1_g[i] < (std::int32_t)g.size());
265  x[_dofs0[i]] = scale * g[_dofs1_g[i]];
266  }
267  }
268  }
269 
274  void set(xtl::span<T> x, const xtl::span<const T>& x0,
275  double scale = 1.0) const
276  {
277  assert(_g);
278  const std::vector<T>& g = _g->x()->array();
279  assert(x.size() <= x0.size());
280  for (std::size_t i = 0; i < _dofs0.size(); ++i)
281  {
282  if (_dofs0[i] < (std::int32_t)x.size())
283  {
284  assert(_dofs1_g[i] < (std::int32_t)g.size());
285  x[_dofs0[i]] = scale * (g[_dofs1_g[i]] - x0[_dofs0[i]]);
286  }
287  }
288  }
289 
298  void dof_values(xtl::span<T> values) const
299  {
300  assert(_g);
301  const std::vector<T>& g = _g->x()->array();
302  for (std::size_t i = 0; i < _dofs1_g.size(); ++i)
303  values[_dofs0[i]] = g[_dofs1_g[i]];
304  }
305 
312  void mark_dofs(std::vector<bool>& markers) const
313  {
314  for (std::size_t i = 0; i < _dofs0.size(); ++i)
315  {
316  assert(_dofs0[i] < (std::int32_t)markers.size());
317  markers[_dofs0[i]] = true;
318  }
319  }
320 
321 private:
322  // The function space (possibly a sub function space)
323  std::shared_ptr<const fem::FunctionSpace> _function_space;
324 
325  // The function
326  std::shared_ptr<const fem::Function<T>> _g;
327 
328  // Dof indices (_dofs0) in _function_space and (_dofs1_g) in the
329  // space of _g
330  std::vector<std::int32_t> _dofs0, _dofs1_g;
331 
332  // The first _owned_indices in _dofs are owned by this process
333  int _owned_indices0 = -1;
334  int _owned_indices1 = -1;
335 };
336 } // namespace dolfinx::fem
Interface for setting (strong) Dirichlet boundary conditions.
Definition: DirichletBC.h:128
DirichletBC & operator=(const DirichletBC &bc)=default
Assignment operator.
DirichletBC & operator=(DirichletBC &&bc)=default
Move assignment operator.
DirichletBC(DirichletBC &&bc)=default
Move constructor.
~DirichletBC()=default
Destructor.
void dof_values(xtl::span< T > values) const
Definition: DirichletBC.h:298
DirichletBC(const std::shared_ptr< const fem::Function< T >> &g, const std::array< std::vector< std::int32_t >, 2 > &V_g_dofs, std::shared_ptr< const fem::FunctionSpace > V)
Create a representation of a Dirichlet boundary condition where the space being constrained and the f...
Definition: DirichletBC.h:189
void set(xtl::span< T > x, const xtl::span< const T > &x0, double scale=1.0) const
Set bc entries in x to scale * (x0 - x_bc)
Definition: DirichletBC.h:274
DirichletBC(const DirichletBC &bc)=default
Copy constructor.
void set(xtl::span< T > x, double scale=1.0) const
Set bc entries in x to scale * x_bc
Definition: DirichletBC.h:256
std::shared_ptr< const fem::FunctionSpace > function_space() const
The function space to which boundary conditions are applied.
Definition: DirichletBC.h:225
std::pair< xtl::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:240
DirichletBC(const std::shared_ptr< const fem::Function< T >> &g, U &&dofs)
Create a representation of a Dirichlet boundary condition where the space being constrained is the sa...
Definition: DirichletBC.h:146
std::shared_ptr< const fem::Function< T > > value() const
Return boundary value function g.
Definition: DirichletBC.h:232
void mark_dofs(std::vector< bool > &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:312
This class represents a finite element function space defined by a mesh, a finite element,...
Definition: FunctionSpace.h:31
This class represents a function in a finite element function space , given by.
Definition: Function.h:47
A Mesh consists of a set of connected and numbered mesh topological entities, and geometry data.
Definition: Mesh.h:53
Finite element method functionality.
Definition: assemble_matrix_impl.h:23
std::array< std::vector< std::int32_t >, 2 > locate_dofs_topological(const std::array< std::reference_wrapper< const fem::FunctionSpace >, 2 > &V, const int dim, const xtl::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:279
std::array< std::vector< std::int32_t >, 2 > locate_dofs_geometrical(const std::array< std::reference_wrapper< const fem::FunctionSpace >, 2 > &V, const std::function< xt::xtensor< bool, 1 >(const xt::xtensor< double, 2 > &)> &marker_fn)
Finds degrees of freedom whose geometric coordinate is true for the provided marking function.
Definition: DirichletBC.cpp:485
Mesh data structures and algorithms on meshes.
Definition: DirichletBC.h:20