Reference documentation for deal.II version 9.3.3
\(\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 - 2021 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
25namespace parallel
26{
27 template <int dim, int spacedim>
29 const ::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 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 [factor](const typename DoFHandler<dim, spacedim>::cell_iterator &,
68 const FiniteElement<dim, spacedim> &) -> unsigned int {
69 return factor;
70 };
71 }
72
73
74
75 template <int dim, int spacedim>
78 const std::pair<float, float> &coefficients)
79 {
80 return [coefficients](
82 const FiniteElement<dim, spacedim> &future_fe) -> unsigned int {
83 const float result =
84 std::trunc(coefficients.first *
85 std::pow(future_fe.n_dofs_per_cell(), coefficients.second));
86
87 Assert(result >= 0. &&
88 result <=
89 static_cast<float>(std::numeric_limits<unsigned int>::max()),
91 "Cannot cast determined weight for this cell to unsigned int!"));
92
93 return static_cast<unsigned int>(result);
94 };
95 }
96
97
98
99 template <int dim, int spacedim>
102 const std::vector<std::pair<float, float>> &coefficients)
103 {
104 return [coefficients](
106 const FiniteElement<dim, spacedim> &future_fe) -> unsigned int {
107 float result = 0;
108 for (const auto &pair : coefficients)
109 result +=
110 pair.first * std::pow(future_fe.n_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()),
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 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,
143 "parallel::CellWeights requires a parallel::TriangulationBase object."));
144
145 return
146 [&dof_handler, tria, weighting_function](
147 const typename ::Triangulation<dim, spacedim>::cell_iterator
148 & cell,
149 const typename ::Triangulation<dim, spacedim>::CellStatus status)
150 -> unsigned int {
152 cell,
153 status,
154 std::cref(dof_handler),
155 std::cref(*tria),
156 weighting_function);
157 };
158 }
159
160
161
162 template <int dim, int spacedim>
163 unsigned int
165 const typename ::Triangulation<dim, spacedim>::cell_iterator &cell_,
166 const typename ::Triangulation<dim, spacedim>::CellStatus status,
167 const DoFHandler<dim, spacedim> & dof_handler,
170 &weighting_function)
171 {
172 // Check if we are still working with the correct combination of
173 // Triangulation and DoFHandler.
174 Assert(&triangulation == &(dof_handler.get_triangulation()),
176 "Triangulation associated with the DoFHandler has changed!"));
177 (void)triangulation;
178
179 // Skip if the DoFHandler has not been initialized yet.
180 if (dof_handler.get_fe_collection().size() == 0)
181 return 0;
182
183 // Convert cell type from Triangulation to DoFHandler to be able
184 // to access the information about the degrees of freedom.
185 const typename DoFHandler<dim, spacedim>::cell_iterator cell(*cell_,
186 &dof_handler);
187
188 // Determine which FiniteElement object will be present on this cell after
189 // refinement and will thus specify the number of degrees of freedom.
190 unsigned int fe_index = numbers::invalid_unsigned_int;
191 switch (status)
192 {
196 fe_index = cell->future_fe_index();
197 break;
198
200#ifdef DEBUG
201 for (const auto &child : cell->child_iterators())
202 Assert(child->is_active() && child->coarsen_flag_set(),
204 dim>::ExcInconsistentCoarseningFlags());
205#endif
206
207 fe_index = ::internal::hp::DoFHandlerImplementation::
208 dominated_future_fe_on_children<dim, spacedim>(cell);
209 break;
210
211 default:
212 Assert(false, ExcInternalError());
213 break;
214 }
215
216 // Return the cell weight determined by the function of choice.
217 return weighting_function(cell, dof_handler.get_fe(fe_index));
218 }
219
220
221
222 // ---------- deprecated functions ----------
223
224 template <int dim, int spacedim>
226 const ::DoFHandler<dim, spacedim> &dof_handler)
227 : dof_handler(&dof_handler)
229 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(
230 &(dof_handler.get_triangulation())))
231 {
232 Assert(
233 triangulation != nullptr,
235 "parallel::CellWeights requires a parallel::TriangulationBase object."));
236
238 }
239
240
241
242 template <int dim, int spacedim>
243 void
245 const unsigned int factor)
246 {
247 connection.disconnect();
248
249 connection = triangulation->signals.cell_weight.connect(
250 [this,
251 factor](const typename Triangulation<dim, spacedim>::cell_iterator &cell,
252 const typename Triangulation<dim, spacedim>::CellStatus status)
253 -> unsigned int {
254 return this->weighting_callback(cell,
255 status,
256 std::cref(*(this->dof_handler)),
257 std::cref(*(this->triangulation)),
258 this->constant_weighting(factor));
259 });
260 }
261
262
263
264 template <int dim, int spacedim>
265 void
267 const unsigned int factor)
268 {
269 connection.disconnect();
270
271 connection = triangulation->signals.cell_weight.connect(
272 [this,
273 factor](const typename Triangulation<dim, spacedim>::cell_iterator &cell,
274 const typename Triangulation<dim, spacedim>::CellStatus status)
275 -> unsigned int {
276 return this->weighting_callback(cell,
277 status,
278 std::cref(*(this->dof_handler)),
279 std::cref(*(this->triangulation)),
280 this->ndofs_weighting({factor, 1}));
281 });
282 }
283
284
285
286 template <int dim, int spacedim>
287 void
289 const unsigned int factor)
290 {
291 connection.disconnect();
292
293 connection = triangulation->signals.cell_weight.connect(
294 [this,
295 factor](const typename Triangulation<dim, spacedim>::cell_iterator &cell,
296 const typename Triangulation<dim, spacedim>::CellStatus status)
297 -> unsigned int {
298 return this->weighting_callback(cell,
299 status,
300 std::cref(*(this->dof_handler)),
301 std::cref(*(this->triangulation)),
302 this->ndofs_weighting({factor, 2}));
303 });
304 }
305
306
307
308 template <int dim, int spacedim>
309 void
311 const std::function<
312 unsigned int(const FiniteElement<dim, spacedim> &,
314 custom_function)
315 {
316 connection.disconnect();
317
318 const std::function<
319 unsigned int(const typename DoFHandler<dim, spacedim>::cell_iterator &,
321 converted_function =
322 [&custom_function](
323 const typename DoFHandler<dim, spacedim>::cell_iterator &cell,
324 const FiniteElement<dim, spacedim> &future_fe) -> unsigned int {
325 return custom_function(future_fe, cell);
326 };
327
328 connection = triangulation->signals.cell_weight.connect(
329 [this, converted_function](
331 const typename Triangulation<dim, spacedim>::CellStatus status)
332 -> unsigned int {
333 return this->weighting_callback(cell,
334 status,
335 std::cref(*(this->dof_handler)),
336 std::cref(*(this->triangulation)),
337 converted_function);
338 });
339 }
340} // namespace parallel
341
342
343// explicit instantiations
344#include "cell_weights.inst"
345
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
const Triangulation< dim, spacedim > & get_triangulation() const
Signals signals
Definition: tria.h:2295
unsigned int size() const
Definition: collection.h:109
void register_constant_weighting(const unsigned int factor=1000)
void register_custom_weighting(const std::function< unsigned int(const FiniteElement< dim, spacedim > &, const typename DoFHandler< dim, spacedim >::cell_iterator &)> custom_function)
void register_ndofs_weighting(const unsigned int factor=1000)
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)
void register_ndofs_squared_weighting(const unsigned int factor=1000)
SmartPointer< const parallel::TriangulationBase< dim, spacedim >, CellWeights > triangulation
Definition: cell_weights.h:274
static WeightingFunction constant_weighting(const unsigned int factor=1000)
Definition: cell_weights.cc:65
std::function< unsigned int(const typename DoFHandler< dim, spacedim >::cell_iterator &, const FiniteElement< dim, spacedim > &)> WeightingFunction
Definition: cell_weights.h:106
static WeightingFunction ndofs_weighting(const std::pair< float, float > &coefficients)
Definition: cell_weights.cc:77
CellWeights(const ::DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function)
Definition: cell_weights.cc:28
void reinit(const DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function)
Definition: cell_weights.cc:47
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:618
static const unsigned int invalid_unsigned_int
Definition: types.h:196
TriangulationBase< dim, spacedim > Triangulation
Definition: tria_base.h:396
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation