Reference documentation for deal.II version 9.4.1
|
#include <deal.II/non_matching/mesh_classifier.h>
Public Member Functions | |
template<class VectorType > | |
MeshClassifier (const DoFHandler< dim > &level_set_dof_handler, const VectorType &level_set) | |
MeshClassifier (const Triangulation< dim > &triangulation, const Function< dim > &level_set, const FiniteElement< dim > &element) | |
void | reclassify () |
LocationToLevelSet | location_to_level_set (const typename Triangulation< dim >::cell_iterator &cell) const |
LocationToLevelSet | location_to_level_set (const typename Triangulation< dim >::cell_iterator &cell, const unsigned int face_index) const |
Private Member Functions | |
void | initialize () |
LocationToLevelSet | determine_face_location_to_levelset (const typename Triangulation< dim >::active_cell_iterator &cell, const unsigned int face_index) |
Private Attributes | |
const SmartPointer< const Triangulation< dim > > | triangulation |
const std::unique_ptr< internal::MeshClassifierImplementation::LevelSetDescription< dim > > | level_set_description |
std::vector< LocationToLevelSet > | cell_locations |
std::vector< LocationToLevelSet > | face_locations |
std::vector< std::array< LAPACKFullMatrix< double >, GeometryInfo< dim >::faces_per_cell > > | lagrange_to_bernstein_face |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
void | check_no_subscribers () const noexcept |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
static std::mutex | mutex |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Class responsible for determining how the active cells and faces of a triangulation relate to the sign of a level set function. When calling the reclassify() function each of the active cells and faces are categorized as one of the values of LocationToLevelSet: inside, outside or intersected, depending on the sign of the level set function over the cell/face. This information is typically required in immersed/cut finite element methods, both when distributing degrees of freedom over the triangulation and when the system is assembled. The given class would then be used in the following way:
The level set function can either be described as a discrete function by a (DoFHandler, Vector)-pair or as a general Function. In the case of a discrete function, LocationToLevelSet for a given face is determined by looking at the local degrees of freedom on the face. Since the Lagrange basis functions are not positive definite, positive/negative definite dof values do not imply that the interpolated function is positive/negative definite. Thus, to classify a face this class internally transforms the local dofs to a basis spanning the same polynomial space over the face but where definite dof values imply a definite function. Currently, only the case of FE_Q-elements is implemented, where we internally change basis to FE_Bernstein. For cells, LocationToLevelSet is determined from the faces of the cell. That is, if all faces of the cell are inside/outside the LocationToLevelSet of the cell is set to inside/outside. LocationToLevelSet of the cell is set to intersected if at least one face is intersected or if its faces have different LocationToLevelSet. Note that, this procedure will incorrectly classify the cell as inside/outside, if the mesh refinement is so low that the whole zero-contour is contained in a single cell (so that none of its faces are intersected).
When the level set function is described as a Function, the level set function is locally interpolated to an FE_Q element and we proceed in the same way as for the discrete level set function.
Definition at line 109 of file mesh_classifier.h.
NonMatching::MeshClassifier< dim >::MeshClassifier | ( | const DoFHandler< dim > & | level_set_dof_handler, |
const VectorType & | level_set | ||
) |
Constructor. Takes a level set function described as a DoFHandler and a Vector. The triangulation attached to DoFHandler is the one that will be classified.
Definition at line 315 of file mesh_classifier.cc.
NonMatching::MeshClassifier< dim >::MeshClassifier | ( | const Triangulation< dim > & | triangulation, |
const Function< dim > & | level_set, | ||
const FiniteElement< dim > & | element | ||
) |
Constructor. Takes the triangulation that should be classified, a level set function described as a Function, and a scalar element that we interpolate the Function to in order to classify each cell/face.
Definition at line 345 of file mesh_classifier.cc.
void NonMatching::MeshClassifier< dim >::reclassify |
Perform the classification of the non artificial cells and faces in the triangulation.
Definition at line 363 of file mesh_classifier.cc.
LocationToLevelSet NonMatching::MeshClassifier< dim >::location_to_level_set | ( | const typename Triangulation< dim >::cell_iterator & | cell | ) | const |
Return how the incoming cell is located relative to the level set function.
Definition at line 433 of file mesh_classifier.cc.
LocationToLevelSet NonMatching::MeshClassifier< dim >::location_to_level_set | ( | const typename Triangulation< dim >::cell_iterator & | cell, |
const unsigned int | face_index | ||
) | const |
Return how a face of the incoming cell is located relative to the level set function.
Definition at line 448 of file mesh_classifier.cc.
|
private |
For each element in the hp::FECollection returned by level_set_description, sets up the local transformation matrices.
Definition at line 465 of file mesh_classifier.cc.
|
private |
Computes how the face with the given index on the incoming cell is located relative to the level set function.
Definition at line 400 of file mesh_classifier.cc.
|
private |
Pointer to the triangulation that should be classified.
Definition at line 176 of file mesh_classifier.h.
|
private |
Pointer to an object that describes what we need to know about the level set function. The underlying object will be of different type depending on whether the level set function is discrete (DoFHandler, Vector) or described by a Function.
Definition at line 186 of file mesh_classifier.h.
|
private |
A vector that stores how each active cell is located relative to the level set function, based on the cells active index.
Definition at line 192 of file mesh_classifier.h.
|
private |
A vector that stores how each active face is located relative to the level set function, based on the face's global index.
Definition at line 198 of file mesh_classifier.h.
|
private |
For each element in the hp::FECollection returned by the LevelSetDescription, and for each local face, this vector stores a transformation matrix to a basis where positive/negative definite face dofs implies that the underlying function is positive/negative definite over the face.
Definition at line 209 of file mesh_classifier.h.