Reference documentation for deal.II version 9.2.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\}}\)
cell_weights.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2018 - 2020 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 
23 #include <deal.II/hp/dof_handler.h>
24 
25 
27 
28 namespace parallel
29 {
86  template <int dim, int spacedim = dim>
88  {
89  public:
100  using WeightingFunction = std::function<unsigned int(
103 
113  const WeightingFunction & weighting_function);
114 
120  ~CellWeights();
121 
128  void
130  const WeightingFunction & weighting_function);
131 
140  static std::function<unsigned int(
141  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
142  const typename Triangulation<dim, spacedim>::CellStatus status)>
144  const WeightingFunction &weighting_function);
145 
154  static WeightingFunction
155  constant_weighting(const unsigned int factor = 1000);
156 
166  static WeightingFunction
167  ndofs_weighting(const std::pair<float, float> &coefficients);
168 
178  static WeightingFunction
179  ndofs_weighting(const std::vector<std::pair<float, float>> &coefficients);
180 
198 
204  DEAL_II_DEPRECATED void
205  register_constant_weighting(const unsigned int factor = 1000);
206 
213  DEAL_II_DEPRECATED void
214  register_ndofs_weighting(const unsigned int factor = 1000);
215 
222  DEAL_II_DEPRECATED void
223  register_ndofs_squared_weighting(const unsigned int factor = 1000);
224 
237  DEAL_II_DEPRECATED void
239  const std::function<unsigned int(
242  custom_function);
243 
248  private:
260 
271 
280  boost::signals2::connection connection;
281 
288  static unsigned int
290  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
291  const typename Triangulation<dim, spacedim>::CellStatus status,
294  const WeightingFunction &weighting_function);
295  };
296 } // namespace parallel
297 
298 
300 
301 #endif
parallel::CellWeights::constant_weighting
static WeightingFunction constant_weighting(const unsigned int factor=1000)
Definition: cell_weights.cc:65
parallel::CellWeights
Definition: cell_weights.h:87
parallel::CellWeights::register_ndofs_weighting
void register_ndofs_weighting(const unsigned int factor=1000)
Definition: cell_weights.cc:272
parallel::CellWeights::connection
boost::signals2::connection connection
Definition: cell_weights.h:280
hp::DoFHandler
Definition: dof_handler.h:203
Triangulation< dim, dim >::CellStatus
CellStatus
Definition: tria.h:1969
parallel::CellWeights::make_weighting_callback
static std::function< unsigned int(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status)> make_weighting_callback(const hp::DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function)
Definition: cell_weights.cc:131
parallel::CellWeights::reinit
void reinit(const hp::DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function)
Definition: cell_weights.cc:47
parallel::CellWeights::~CellWeights
~CellWeights()
Definition: cell_weights.cc:38
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
parallel::CellWeights::weighting_callback
static unsigned int weighting_callback(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status, const hp::DoFHandler< dim, spacedim > &dof_handler, const parallel::TriangulationBase< dim, spacedim > &triangulation, const WeightingFunction &weighting_function)
Definition: cell_weights.cc:162
DEAL_II_DEPRECATED
#define DEAL_II_DEPRECATED
Definition: config.h:98
parallel::CellWeights::CellWeights
CellWeights(const hp::DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function)
Definition: cell_weights.cc:28
FiniteElement
Definition: fe.h:648
hp::DoFHandler::cell_iterator
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:325
SmartPointer
Definition: smartpointer.h:68
parallel::CellWeights::register_constant_weighting
void register_constant_weighting(const unsigned int factor=1000)
Definition: cell_weights.cc:250
parallel::CellWeights::ndofs_weighting
static WeightingFunction ndofs_weighting(const std::pair< float, float > &coefficients)
Definition: cell_weights.cc:78
parallel::CellWeights::WeightingFunction
std::function< unsigned int(const typename hp::DoFHandler< dim, spacedim >::cell_iterator &, const FiniteElement< dim, spacedim > &)> WeightingFunction
Definition: cell_weights.h:102
parallel::TriangulationBase
Definition: tria_base.h:78
parallel::CellWeights::register_ndofs_squared_weighting
void register_ndofs_squared_weighting(const unsigned int factor=1000)
Definition: cell_weights.cc:294
config.h
parallel::CellWeights::dof_handler
SmartPointer< const ::hp::DoFHandler< dim, spacedim >, CellWeights > dof_handler
Definition: cell_weights.h:259
parallel::CellWeights::triangulation
SmartPointer< const parallel::TriangulationBase< dim, spacedim >, CellWeights > triangulation
Definition: cell_weights.h:270
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
TriaIterator
Definition: tria_iterator.h:578
parallel::CellWeights::register_custom_weighting
void register_custom_weighting(const std::function< unsigned int(const FiniteElement< dim, spacedim > &, const typename hp::DoFHandler< dim, spacedim >::cell_iterator &)> custom_function)
Definition: cell_weights.cc:316
parallel
Definition: distributed.h:416
tria_base.h
int
dof_handler.h