Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30:02+00:00
\(\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.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2018 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
17
19
20#include <limits>
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.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 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 [&dof_handler, tria, weighting_function](
146 const typename ::Triangulation<dim, spacedim>::cell_iterator
147 &cell,
148 const CellStatus status) -> 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 CellStatus status,
165 const 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 Assert(&triangulation == &(dof_handler.get_triangulation()),
174 "Triangulation associated with the DoFHandler has changed!"));
175 (void)triangulation;
176
177 // Skip if the DoFHandler has not been initialized yet.
178 if (dof_handler.get_fe_collection().size() == 0)
179 return 0;
180
181 // Convert cell type from Triangulation to DoFHandler to be able
182 // to access the information about the degrees of freedom.
183 const typename DoFHandler<dim, spacedim>::cell_iterator cell(*cell_,
184 &dof_handler);
185
186 // Determine which FiniteElement object will be present on this cell after
187 // refinement and will thus specify the number of degrees of freedom.
188 unsigned int fe_index = numbers::invalid_unsigned_int;
189 switch (status)
190 {
194 fe_index = cell->future_fe_index();
195 break;
196
198#ifdef DEBUG
199 for (const auto &child : cell->child_iterators())
200 Assert(child->is_active() && child->coarsen_flag_set(),
201 typename ::Triangulation<
202 dim>::ExcInconsistentCoarseningFlags());
203#endif
204
205 fe_index = ::internal::hp::DoFHandlerImplementation::
206 dominated_future_fe_on_children<dim, spacedim>(cell);
207 break;
208
209 default:
211 break;
212 }
213
214 // Return the cell weight determined by the function of choice.
215 return weighting_function(cell, dof_handler.get_fe(fe_index));
216 }
217} // namespace parallel
218
219
220// explicit instantiations
221#include "cell_weights.inst"
222
CellStatus
Definition cell_status.h:31
@ cell_will_be_refined
@ children_will_be_coarsened
const hp::FECollection< dim, spacedim > & get_fe_collection() const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const Triangulation< dim, spacedim > & get_triangulation() const
Signals signals
Definition tria.h:2509
unsigned int size() const
Definition collection.h:264
static unsigned int weighting_callback(const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const CellStatus status, const DoFHandler< dim, spacedim > &dof_handler, const parallel::TriangulationBase< dim, spacedim > &triangulation, 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)
static std::function< unsigned int(const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const CellStatus status)> make_weighting_callback(const DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
#define Assert(cond, exc)
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::cell_iterator cell_iterator
#define DEAL_II_ASSERT_UNREACHABLE()
static const unsigned int invalid_unsigned_int
Definition types.h:220
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation