47 namespace MeshClassifierImplementation
51 "The Triangulation has not been classified. You need to call the "
52 "reclassify()-function before using this function.");
56 "The incoming cell does not belong to the triangulation passed to "
64 template <
class VectorType>
68 const auto min_max_element =
69 std::minmax_element(local_levelset_values.begin(),
70 local_levelset_values.end());
72 if (*min_max_element.second < 0)
74 if (0 < *min_max_element.first)
86 template <
int dim,
class VectorType>
108 &cell)
const override;
117 const unsigned int face_index,
135 template <
int dim,
class VectorType>
138 const VectorType & level_set)
139 : dof_handler(&dof_handler)
140 , level_set(&level_set)
145 template <
int dim,
class VectorType>
149 return dof_handler->get_fe_collection();
154 template <
int dim,
class VectorType>
158 const unsigned int face_index,
162 &dof_handler->get_triangulation(),
167 const unsigned int n_dofs_per_face =
168 dof_handler->get_fe().n_dofs_per_face();
169 std::vector<types::global_dof_index> dof_indices(n_dofs_per_face);
170 cell_with_dofs->face(face_index)->get_dof_indices(dof_indices);
172 local_levelset_values.
reinit(dof_indices.size());
174 for (
unsigned int i = 0; i < dof_indices.size(); i++)
175 local_levelset_values[i] =
182 template <
int dim,
class VectorType>
188 &dof_handler->get_triangulation(),
193 return cell_with_dofs->active_fe_index();
225 &cell)
const override;
234 const unsigned int face_index,
262 : level_set(&level_set)
263 , fe_collection(element)
264 , fe_face_values(element,
266 element.get_unit_face_support_points()),
276 const unsigned int face_index,
280 fe_face_values.n_quadrature_points);
282 fe_face_values.reinit(cell, face_index);
283 const std::vector<Point<dim>> &points =
284 fe_face_values.get_quadrature_points();
286 for (
unsigned int i = 0; i < points.size(); i++)
287 local_levelset_values[i] = level_set->value(points[i]);
296 return fe_collection;
314 template <
class VectorType>
316 const VectorType & level_set)
318 , level_set_description(
319 std::make_unique<
internal::MeshClassifierImplementation::
320 DiscreteLevelSetDescription<dim, VectorType>>(
324#ifdef DEAL_II_WITH_LAPACK
327 for (
unsigned int i = 0; i < fe_collection.
size(); i++)
332 Assert(fe_collection[i].has_face_support_points(),
334 "The elements in the FECollection of the incoming DoFHandler "
335 "must have face support points."));
349 , level_set_description(
350 std::make_unique<
internal::MeshClassifierImplementation::
351 AnalyticLevelSetDescription<dim>>(level_set,
373 for (
const auto &cell :
triangulation->active_cell_iterators())
374 if (!cell->is_artificial())
377 determine_face_location_to_levelset(cell, 0);
379 face_locations[cell->face(0)->index()] = face0_location;
382 for (
unsigned int f = 1; f < GeometryInfo<dim>::faces_per_cell; ++f)
385 determine_face_location_to_levelset(cell, f);
387 face_locations[cell->face(f)->index()] = face_location;
389 if (face_location != face0_location)
392 cell_locations[cell->active_cell_index()] = cell_location;
402 const unsigned int face_index)
407 face_locations.at(cell->face(face_index)->index());
413 const unsigned int fe_index = level_set_description->active_fe_index(cell);
414 const unsigned int n_local_dofs =
415 lagrange_to_bernstein_face[fe_index][face_index].m();
418 level_set_description->get_local_level_set_values(cell,
420 local_levelset_values);
422 lagrange_to_bernstein_face[fe_index][face_index].solve(
423 local_levelset_values);
426 local_levelset_values);
441 return cell_locations.at(cell->active_cell_index());
450 const unsigned int face_index)
const
458 return face_locations.at(cell->face(face_index)->index());
468 level_set_description->get_fe_collection();
473 lagrange_to_bernstein_face.resize(fe_collection.
size());
475 for (
unsigned int i = 0; i < fe_collection.
size(); i++)
484 for (
unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; f++)
490 *fe_q, face_interpolation_matrix, f);
491 lagrange_to_bernstein_face[i][f].reinit(dofs_per_face);
492 lagrange_to_bernstein_face[i][f] = face_interpolation_matrix;
493 lagrange_to_bernstein_face[i][f].compute_lu_factorization();
500#include "mesh_classifier.inst"
const hp::FECollection< dim, spacedim > & get_fe_collection() const
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const override
unsigned int get_degree() const
const unsigned int dofs_per_face
unsigned int n_components() const
const unsigned int n_components
LocationToLevelSet location_to_level_set(const typename Triangulation< dim >::cell_iterator &cell) const
LocationToLevelSet determine_face_location_to_levelset(const typename Triangulation< dim >::active_cell_iterator &cell, const unsigned int face_index)
MeshClassifier(const DoFHandler< dim > &level_set_dof_handler, const VectorType &level_set)
AnalyticLevelSetDescription(const Function< dim > &level_set, const FiniteElement< dim > &element)
unsigned int active_fe_index(const typename Triangulation< dim >::active_cell_iterator &cell) const override
FEFaceValues< dim > fe_face_values
const hp::FECollection< dim > & get_fe_collection() const override
const hp::FECollection< dim > fe_collection
const SmartPointer< const Function< dim > > level_set
void get_local_level_set_values(const typename Triangulation< dim >::active_cell_iterator &cell, const unsigned int face_index, Vector< double > &local_levelset_values) override
void get_local_level_set_values(const typename Triangulation< dim >::active_cell_iterator &cell, const unsigned int face_index, Vector< double > &local_levelset_values) override
const SmartPointer< const DoFHandler< dim > > dof_handler
unsigned int active_fe_index(const typename Triangulation< dim >::active_cell_iterator &cell) const override
const hp::FECollection< dim > & get_fe_collection() const override
const SmartPointer< const VectorType > level_set
DiscreteLevelSetDescription(const DoFHandler< dim > &dof_handler, const VectorType &level_set)
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)
unsigned int size() const
unsigned int n_components() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcNeedsLAPACK()
static ::ExceptionBase & ExcReclassifyNotCalled()
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcTriangulationMismatch()
#define AssertThrow(cond, exc)
typename ActiveSelector::active_cell_iterator active_cell_iterator
@ update_quadrature_points
Transformed quadrature points.
LocationToLevelSet location_from_dof_signs(const VectorType &local_levelset_values)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation