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.cc
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 
18 
20 
21 
23 
24 
25 namespace parallel
26 {
27  template <int dim, int spacedim>
29  const hp::DoFHandler<dim, spacedim> &dof_handler,
30  const WeightingFunction & weighting_function)
31  {
32  reinit(dof_handler, weighting_function);
33  }
34 
35 
36 
37  template <int dim, int spacedim>
39  {
40  connection.disconnect();
41  }
42 
43 
44 
45  template <int dim, int spacedim>
46  void
48  const hp::DoFHandler<dim, spacedim> &dof_handler,
50  &weighting_function)
51  {
52  connection.disconnect();
53  connection = dof_handler.get_triangulation().signals.cell_weight.connect(
54  make_weighting_callback(dof_handler, weighting_function));
55  }
56 
57 
58 
59  // ---------- static member functions ----------
60 
61  // ---------- selection of weighting functions ----------
62 
63  template <int dim, int spacedim>
66  {
67  return
68  [factor](const typename hp::DoFHandler<dim, spacedim>::cell_iterator &,
69  const FiniteElement<dim, spacedim> &) -> unsigned int {
70  return factor;
71  };
72  }
73 
74 
75 
76  template <int dim, int spacedim>
79  const std::pair<float, float> &coefficients)
80  {
81  return [coefficients](
83  const FiniteElement<dim, spacedim> &future_fe) -> unsigned int {
84  const float result =
85  std::trunc(coefficients.first *
86  std::pow(future_fe.dofs_per_cell, coefficients.second));
87 
88  Assert(result >= 0. &&
89  result <=
90  static_cast<float>(std::numeric_limits<unsigned int>::max()),
91  ExcMessage(
92  "Cannot cast determined weight for this cell to unsigned int!"));
93 
94  return static_cast<unsigned int>(result);
95  };
96  }
97 
98 
99 
100  template <int dim, int spacedim>
103  const std::vector<std::pair<float, float>> &coefficients)
104  {
105  return [coefficients](
107  const FiniteElement<dim, spacedim> &future_fe) -> unsigned int {
108  float result = 0;
109  for (const auto &pair : coefficients)
110  result += pair.first * std::pow(future_fe.dofs_per_cell, pair.second);
111  result = std::trunc(result);
112 
113  Assert(result >= 0. &&
114  result <=
115  static_cast<float>(std::numeric_limits<unsigned int>::max()),
116  ExcMessage(
117  "Cannot cast determined weight for this cell to unsigned int!"));
118 
119  return static_cast<unsigned int>(result);
120  };
121  }
122 
123 
124 
125  // ---------- handling callback functions ----------
126 
127  template <int dim, int spacedim>
128  std::function<unsigned int(
129  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
130  const typename Triangulation<dim, spacedim>::CellStatus status)>
132  const hp::DoFHandler<dim, spacedim> &dof_handler,
134  &weighting_function)
135  {
137  dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
138  &(dof_handler.get_triangulation()));
139 
140  Assert(
141  tria != nullptr,
142  ExcMessage(
143  "parallel::CellWeights requires a parallel::TriangulationBase object."));
144 
145  return [&dof_handler, tria, weighting_function](
146  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
147  const typename Triangulation<dim, spacedim>::CellStatus status)
148  -> unsigned int {
150  status,
151  std::cref(
152  dof_handler),
153  std::cref(*tria),
154  weighting_function);
155  };
156  }
157 
158 
159 
160  template <int dim, int spacedim>
161  unsigned int
163  const typename Triangulation<dim, spacedim>::cell_iterator &cell_,
164  const typename Triangulation<dim, spacedim>::CellStatus status,
165  const hp::DoFHandler<dim, spacedim> & dof_handler,
168  &weighting_function)
169  {
170  // Check if we are still working with the correct combination of
171  // Triangulation and DoFHandler.
172  AssertThrow(&triangulation == &(dof_handler.get_triangulation()),
173  ExcMessage(
174  "Triangulation associated with the DoFHandler has changed!"));
175 
176  // Convert cell type from Triangulation to DoFHandler to be able
177  // to access the information about the degrees of freedom.
179  *cell_, &dof_handler);
180 
181  // Determine which FiniteElement object will be present on this cell after
182  // refinement and will thus specify the number of degrees of freedom.
183  unsigned int fe_index = numbers::invalid_unsigned_int;
184  switch (status)
185  {
189  fe_index = cell->future_fe_index();
190  break;
191 
193  {
194  std::set<unsigned int> fe_indices_children;
195  for (unsigned int child_index = 0; child_index < cell->n_children();
196  ++child_index)
197  {
198  const auto &child = cell->child(child_index);
199  Assert(child->is_active() && child->coarsen_flag_set(),
201  dim>::ExcInconsistentCoarseningFlags());
202 
203  fe_indices_children.insert(child->future_fe_index());
204  }
205  Assert(!fe_indices_children.empty(), ExcInternalError());
206 
207  fe_index =
208  dof_handler.get_fe_collection().find_dominated_fe_extended(
209  fe_indices_children, /*codim=*/0);
210 
212  typename ::hp::FECollection<
213  dim>::ExcNoDominatedFiniteElementAmongstChildren());
214  }
215  break;
216 
217  default:
218  Assert(false, ExcInternalError());
219  break;
220  }
221 
222  // Return the cell weight determined by the function of choice.
223  return weighting_function(cell, dof_handler.get_fe(fe_index));
224  }
225 
226 
227 
228  // ---------- deprecated functions ----------
229 
230  template <int dim, int spacedim>
232  const hp::DoFHandler<dim, spacedim> &dof_handler)
233  : dof_handler(&dof_handler)
234  , triangulation(
235  dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
236  &(dof_handler.get_triangulation())))
237  {
238  Assert(
239  triangulation != nullptr,
240  ExcMessage(
241  "parallel::CellWeights requires a parallel::TriangulationBase object."));
242 
244  }
245 
246 
247 
248  template <int dim, int spacedim>
249  void
251  const unsigned int factor)
252  {
253  connection.disconnect();
254 
255  connection = triangulation->signals.cell_weight.connect(
256  [this,
257  factor](const typename Triangulation<dim, spacedim>::cell_iterator &cell,
258  const typename Triangulation<dim, spacedim>::CellStatus status)
259  -> unsigned int {
260  return this->weighting_callback(cell,
261  status,
262  std::cref(*(this->dof_handler)),
263  std::cref(*(this->triangulation)),
264  this->constant_weighting(factor));
265  });
266  }
267 
268 
269 
270  template <int dim, int spacedim>
271  void
273  const unsigned int factor)
274  {
275  connection.disconnect();
276 
277  connection = triangulation->signals.cell_weight.connect(
278  [this,
279  factor](const typename Triangulation<dim, spacedim>::cell_iterator &cell,
280  const typename Triangulation<dim, spacedim>::CellStatus status)
281  -> unsigned int {
282  return this->weighting_callback(cell,
283  status,
284  std::cref(*(this->dof_handler)),
285  std::cref(*(this->triangulation)),
286  this->ndofs_weighting({factor, 1}));
287  });
288  }
289 
290 
291 
292  template <int dim, int spacedim>
293  void
295  const unsigned int factor)
296  {
297  connection.disconnect();
298 
299  connection = triangulation->signals.cell_weight.connect(
300  [this,
301  factor](const typename Triangulation<dim, spacedim>::cell_iterator &cell,
302  const typename Triangulation<dim, spacedim>::CellStatus status)
303  -> unsigned int {
304  return this->weighting_callback(cell,
305  status,
306  std::cref(*(this->dof_handler)),
307  std::cref(*(this->triangulation)),
308  this->ndofs_weighting({factor, 2}));
309  });
310  }
311 
312 
313 
314  template <int dim, int spacedim>
315  void
317  const std::function<unsigned int(
320  custom_function)
321  {
322  connection.disconnect();
323 
324  const std::function<unsigned int(
327  converted_function =
328  [&custom_function](
330  const FiniteElement<dim, spacedim> &future_fe) -> unsigned int {
331  return custom_function(future_fe, cell);
332  };
333 
334  connection = triangulation->signals.cell_weight.connect(
335  [this, converted_function](
336  const typename Triangulation<dim, spacedim>::cell_iterator &cell,
337  const typename Triangulation<dim, spacedim>::CellStatus status)
338  -> unsigned int {
339  return this->weighting_callback(cell,
340  status,
341  std::cref(*(this->dof_handler)),
342  std::cref(*(this->triangulation)),
343  converted_function);
344  });
345  }
346 } // namespace parallel
347 
348 
349 // explicit instantiations
350 #include "cell_weights.inst"
351 
parallel::CellWeights::constant_weighting
static WeightingFunction constant_weighting(const unsigned int factor=1000)
Definition: cell_weights.cc:65
parallel::CellWeights::register_ndofs_weighting
void register_ndofs_weighting(const unsigned int factor=1000)
Definition: cell_weights.cc:272
hp::DoFHandler::get_fe_collection
const hp::FECollection< dim, spacedim > & get_fe_collection() const
hp::DoFHandler::get_triangulation
const Triangulation< dim, spacedim > & get_triangulation() const
hp::DoFHandler::get_fe
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index) const
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
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
internal::reinit
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:621
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
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
cell_weights.h
parallel::Triangulation
TriangulationBase< dim, spacedim > Triangulation
Definition: tria_base.h:302
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
dof_accessor.h
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
numbers::invalid_unsigned_int
static const unsigned int invalid_unsigned_int
Definition: types.h:191
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
triangulation
const typename ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
Definition: p4est_wrappers.cc:69
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
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
parallel
Definition: distributed.h:416
Utilities::MPI::max
T max(const T &t, const MPI_Comm &mpi_communicator)
int