27 template <
int dim,
int spacedim>
29 const ::DoFHandler<dim, spacedim> &dof_handler,
32 reinit(dof_handler, weighting_function);
37 template <
int dim,
int spacedim>
40 connection.disconnect();
45 template <
int dim,
int spacedim>
52 connection.disconnect();
54 make_weighting_callback(dof_handler, weighting_function));
63 template <
int dim,
int spacedim>
75 template <
int dim,
int spacedim>
78 const std::pair<float, float> &coefficients)
80 return [coefficients](
84 std::trunc(coefficients.first *
85 std::pow(future_fe.n_dofs_per_cell(), coefficients.second));
91 "Cannot cast determined weight for this cell to unsigned int!"));
93 return static_cast<unsigned int>(result);
99 template <
int dim,
int spacedim>
102 const std::vector<std::pair<float, float>> &coefficients)
104 return [coefficients](
108 for (
const auto &pair : coefficients)
110 pair.first *
std::pow(future_fe.n_dofs_per_cell(), pair.second);
111 result = std::trunc(result);
117 "Cannot cast determined weight for this cell to unsigned int!"));
119 return static_cast<unsigned int>(result);
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)>
143 "parallel::CellWeights requires a parallel::TriangulationBase object."));
146 [&dof_handler, tria, weighting_function](
147 const typename ::Triangulation<dim, spacedim>::cell_iterator
149 const typename ::Triangulation<dim, spacedim>::CellStatus status)
154 std::cref(dof_handler),
162 template <
int dim,
int spacedim>
165 const typename ::Triangulation<dim, spacedim>::cell_iterator &cell_,
166 const typename ::Triangulation<dim, spacedim>::CellStatus status,
176 "Triangulation associated with the DoFHandler has changed!"));
196 fe_index = cell->future_fe_index();
201 for (
const auto &child : cell->child_iterators())
202 Assert(child->is_active() && child->coarsen_flag_set(),
204 dim>::ExcInconsistentCoarseningFlags());
207 fe_index = ::internal::hp::DoFHandlerImplementation::
208 dominated_future_fe_on_children<dim, spacedim>(cell);
217 return weighting_function(cell, dof_handler.
get_fe(fe_index));
224 template <
int dim,
int spacedim>
226 const ::DoFHandler<dim, spacedim> &dof_handler)
227 : dof_handler(&dof_handler)
230 &(dof_handler.get_triangulation())))
235 "parallel::CellWeights requires a parallel::TriangulationBase object."));
242 template <
int dim,
int spacedim>
245 const unsigned int factor)
247 connection.disconnect();
254 return this->weighting_callback(cell,
256 std::cref(*(this->dof_handler)),
258 this->constant_weighting(factor));
264 template <
int dim,
int spacedim>
267 const unsigned int factor)
269 connection.disconnect();
276 return this->weighting_callback(cell,
278 std::cref(*(this->dof_handler)),
280 this->ndofs_weighting({factor, 1}));
286 template <
int dim,
int spacedim>
289 const unsigned int factor)
291 connection.disconnect();
298 return this->weighting_callback(cell,
300 std::cref(*(this->dof_handler)),
302 this->ndofs_weighting({factor, 2}));
308 template <
int dim,
int spacedim>
316 connection.disconnect();
325 return custom_function(future_fe, cell);
329 [
this, converted_function](
333 return this->weighting_callback(cell,
335 std::cref(*(this->dof_handler)),
344#include "cell_weights.inst"
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
unsigned int size() const
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
static WeightingFunction constant_weighting(const unsigned int factor=1000)
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)
CellWeights(const ::DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function)
void reinit(const DoFHandler< dim, spacedim > &dof_handler, const WeightingFunction &weighting_function)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
static const unsigned int invalid_unsigned_int
TriangulationBase< dim, spacedim > Triangulation
::VectorizedArray< Number, width > pow(const ::VectorizedArray< Number, width > &, const Number p)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation