Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
cell_weights.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2018 - 2023 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_distributed_cell_weights_h
17#define dealii_distributed_cell_weights_h
18
19#include <deal.II/base/config.h>
20
22
24
25
27
28namespace parallel
29{
99 template <int dim, int spacedim = dim>
101 {
102 public:
113 using WeightingFunction = std::function<
114 unsigned int(const typename DoFHandler<dim, spacedim>::cell_iterator &,
116
122 CellWeights() = default;
123
132 CellWeights(const DoFHandler<dim, spacedim> &dof_handler,
133 const WeightingFunction & weighting_function);
134
140 ~CellWeights();
141
148 void
149 reinit(const DoFHandler<dim, spacedim> &dof_handler,
150 const WeightingFunction & weighting_function);
151
160 static std::function<unsigned int(
161 const typename ::Triangulation<dim, spacedim>::cell_iterator &cell,
162 const typename ::Triangulation<dim, spacedim>::CellStatus status)>
164 const WeightingFunction &weighting_function);
165
174 static WeightingFunction
175 constant_weighting(const unsigned int factor = 1);
176
186 static WeightingFunction
187 ndofs_weighting(const std::pair<float, float> &coefficients);
188
198 static WeightingFunction
199 ndofs_weighting(const std::vector<std::pair<float, float>> &coefficients);
200
205 private:
210 boost::signals2::connection connection;
211
219 static unsigned int
221 const typename ::Triangulation<dim, spacedim>::cell_iterator &cell,
222 const typename ::Triangulation<dim, spacedim>::CellStatus status,
223 const DoFHandler<dim, spacedim> & dof_handler,
225 const WeightingFunction & weighting_function);
226 };
227} // namespace parallel
228
229
231
232#endif
static unsigned int weighting_callback(const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const typename ::Triangulation< dim, spacedim >::CellStatus status, const DoFHandler< dim, spacedim > &dof_handler, const parallel::TriangulationBase< dim, spacedim > &triangulation, const WeightingFunction &weighting_function)
static std::function< unsigned int(const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const typename ::Triangulation< dim, spacedim >::CellStatus status)> make_weighting_callback(const DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function)
static WeightingFunction constant_weighting(const unsigned int factor=1)
std::function< unsigned int(const typename DoFHandler< dim, spacedim >::cell_iterator &, const FiniteElement< dim, spacedim > &)> WeightingFunction
static WeightingFunction ndofs_weighting(const std::pair< float, float > &coefficients)
void reinit(const DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function)
boost::signals2::connection connection
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
typename ActiveSelector::cell_iterator cell_iterator
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation