MAST
MAST::LevelSetIntersection Class Reference

#include <level_set_intersection.h>

Public Member Functions

 LevelSetIntersection ()
 
virtual ~LevelSetIntersection ()
 
void init (const MAST::FieldFunction< Real > &phi, const libMesh::Elem &e, const Real t, unsigned int max_elem_id, unsigned int max_node_id)
 
void clear ()
 clears the data structures More...
 
const libMesh::Elem & elem () const
 
MAST::LevelSet2DIntersectionMode get_intersection_mode () const
 
unsigned int node_on_boundary () const
 
unsigned int edge_on_boundary () const
 
Real get_positive_phi_volume_fraction () const
 
Real get_node_phi_value (const libMesh::Node *n) const
 
bool if_elem_has_boundary () const
 
bool if_intersection_through_elem () const
 
bool if_elem_on_positive_phi () const
 
bool if_elem_on_negative_phi () const
 
bool if_elem_has_positive_phi_region () const
 
const std::vector< const libMesh::Elem * > & get_sub_elems_positive_phi () const
 
const std::vector< const libMesh::Elem * > & get_sub_elems_negative_phi () const
 
void get_corner_nodes_on_negative_phi (std::set< const libMesh::Node * > &nodes) const
 
unsigned int get_side_on_interface (const libMesh::Elem &e) const
 
bool has_side_on_interface (const libMesh::Elem &e) const
 
const libMesh::Point & get_nondimensional_coordinate_for_node (const libMesh::Node &n) const
 

Protected Member Functions

std::unique_ptr< libMesh::Elem > _first_order_elem (const libMesh::Elem &e)
 creates a first order element from the given high-order element. More...
 
void _init_on_first_order_ref_elem (const MAST::FieldFunction< Real > &phi, const libMesh::Elem &e, const Real t)
 initializes on a reference element that is a first-order counterpart of the given high-order element. More...
 
void _add_node_local_coords (const libMesh::Elem &e, std::vector< std::pair< libMesh::Point, libMesh::Point > > &side_nondim_points, std::map< const libMesh::Node *, libMesh::Point > &node_coord_map)
 
void _find_quad4_intersections (const MAST::FieldFunction< Real > &phi, const libMesh::Elem &e, const Real t, const std::map< const libMesh::Node *, std::pair< Real, bool > > &node_phi_vals)
 
Real _find_intersection_on_straight_edge (const libMesh::Point &p0, const libMesh::Point &p1, const MAST::FieldFunction< Real > &phi, const Real t)
 

Protected Attributes

Real _tol
 
unsigned int _max_iters
 
unsigned int _max_mesh_elem_id
 
unsigned int _max_mesh_node_id
 
const unsigned int _max_elem_divs
 
const libMesh::Elem * _elem
 
bool _initialized
 
bool _if_elem_on_positive_phi
 true if element is completely on the positive side of level set with no intersection More...
 
bool _if_elem_on_negative_phi
 true if element is completely on the negative side of level set with no intersection More...
 
MAST::LevelSet2DIntersectionMode _mode
 
unsigned int _node_num_on_boundary
 
unsigned int _edge_num_on_boundary
 
std::vector< const libMesh::Elem * > _positive_phi_elems
 
std::vector< const libMesh::Elem * > _negative_phi_elems
 
std::map< const libMesh::Elem *, int > _elem_sides_on_interface
 
std::vector< libMesh::Node * > _new_nodes
 
std::vector< libMesh::Elem * > _new_elems
 
std::map< const libMesh::Node *, libMesh::Point > _node_local_coords
 
std::map< const libMesh::Node *, std::pair< Real, bool > > _node_phi_vals
 

Detailed Description

Definition at line 54 of file level_set_intersection.h.

Constructor & Destructor Documentation

MAST::LevelSetIntersection::LevelSetIntersection ( )

Definition at line 31 of file level_set_intersection.cpp.

MAST::LevelSetIntersection::~LevelSetIntersection ( )
virtual

Definition at line 48 of file level_set_intersection.cpp.

Here is the call graph for this function:

Member Function Documentation

void MAST::LevelSetIntersection::_add_node_local_coords ( const libMesh::Elem &  e,
std::vector< std::pair< libMesh::Point, libMesh::Point > > &  side_nondim_points,
std::map< const libMesh::Node *, libMesh::Point > &  node_coord_map 
)
protected

Definition at line 594 of file level_set_intersection.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

Real MAST::LevelSetIntersection::_find_intersection_on_straight_edge ( const libMesh::Point &  p0,
const libMesh::Point &  p1,
const MAST::FieldFunction< Real > &  phi,
const Real  t 
)
protected
Returns
the value along a straight edge where the level-set function zero.

Definition at line 1556 of file level_set_intersection.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::LevelSetIntersection::_find_quad4_intersections ( const MAST::FieldFunction< Real > &  phi,
const libMesh::Elem &  e,
const Real  t,
const std::map< const libMesh::Node *, std::pair< Real, bool > > &  node_phi_vals 
)
protected

Definition at line 640 of file level_set_intersection.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

std::unique_ptr< libMesh::Elem > MAST::LevelSetIntersection::_first_order_elem ( const libMesh::Elem &  e)
protected

creates a first order element from the given high-order element.

For a QUAD9 a QUAD4 is obtained by only using the corner nodes. Note that this does not create any new nodes. The element can be deleted after use.

Definition at line 260 of file level_set_intersection.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::LevelSetIntersection::_init_on_first_order_ref_elem ( const MAST::FieldFunction< Real > &  phi,
const libMesh::Elem &  e,
const Real  t 
)
protected

initializes on a reference element that is a first-order counterpart of the given high-order element.

For two-dimensional elements this is a QUAD4.

Definition at line 290 of file level_set_intersection.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

void MAST::LevelSetIntersection::clear ( )

clears the data structures

Definition at line 184 of file level_set_intersection.cpp.

Here is the caller graph for this function:

unsigned int MAST::LevelSetIntersection::edge_on_boundary ( ) const
Returns
the edge number on the element if the mode is COLINEAR_EDGE, otherwise throws an error.

Definition at line 85 of file level_set_intersection.cpp.

const libMesh::Elem & MAST::LevelSetIntersection::elem ( ) const
Returns
a reference to the element on which the intersection is defined

Definition at line 55 of file level_set_intersection.cpp.

Here is the caller graph for this function:

void MAST::LevelSetIntersection::get_corner_nodes_on_negative_phi ( std::set< const libMesh::Node * > &  nodes) const

Definition at line 525 of file level_set_intersection.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

MAST::LevelSet2DIntersectionMode MAST::LevelSetIntersection::get_intersection_mode ( ) const
Returns
mode of intersection

Definition at line 64 of file level_set_intersection.cpp.

Real MAST::LevelSetIntersection::get_node_phi_value ( const libMesh::Node *  n) const
Returns
value of phi on element node. This can only be the node in the

Definition at line 170 of file level_set_intersection.cpp.

Here is the caller graph for this function:

const libMesh::Point & MAST::LevelSetIntersection::get_nondimensional_coordinate_for_node ( const libMesh::Node &  n) const

Definition at line 1619 of file level_set_intersection.cpp.

Here is the caller graph for this function:

Real MAST::LevelSetIntersection::get_positive_phi_volume_fraction ( ) const
Returns
the area/volume fraction of element on positive phi of the element

Definition at line 149 of file level_set_intersection.cpp.

Here is the caller graph for this function:

unsigned int MAST::LevelSetIntersection::get_side_on_interface ( const libMesh::Elem &  e) const
Returns
the id of side that is on the interface. In case the element does not have a side on the interface, then a negative value is returned.

Definition at line 580 of file level_set_intersection.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

const std::vector< const libMesh::Elem * > & MAST::LevelSetIntersection::get_sub_elems_negative_phi ( ) const

Definition at line 517 of file level_set_intersection.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

const std::vector< const libMesh::Elem * > & MAST::LevelSetIntersection::get_sub_elems_positive_phi ( ) const

Definition at line 509 of file level_set_intersection.cpp.

Here is the caller graph for this function:

bool MAST::LevelSetIntersection::has_side_on_interface ( const libMesh::Elem &  e) const

Definition at line 566 of file level_set_intersection.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

bool MAST::LevelSetIntersection::if_elem_has_boundary ( ) const
Returns
true if the intersection is through the element, or has colinear edge.

Definition at line 126 of file level_set_intersection.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

bool MAST::LevelSetIntersection::if_elem_has_positive_phi_region ( ) const
Returns
true if there is any portion of the element (interior or edge) that is on the positive side of the level set function.

Definition at line 116 of file level_set_intersection.cpp.

Here is the caller graph for this function:

bool MAST::LevelSetIntersection::if_elem_on_negative_phi ( ) const
Returns
true if the element is entirely on the negative side of the level set without any intersection. This will return true only if the mode is NO_INTERSECTION and all nodal phi values are negative.

Definition at line 106 of file level_set_intersection.cpp.

Here is the caller graph for this function:

bool MAST::LevelSetIntersection::if_elem_on_positive_phi ( ) const
Returns
true if the element is entirely on the positive side of the level set without any intersection. This will return true only if the mode is NO_INTERSECTION and all nodal phi values are positive.

Definition at line 97 of file level_set_intersection.cpp.

bool MAST::LevelSetIntersection::if_intersection_through_elem ( ) const
Returns
true if the intersection passes through the interior of the element. Note that COLINEAR_EDGE or THROUGH_NODE do not qualify for this.

Definition at line 136 of file level_set_intersection.cpp.

Here is the caller graph for this function:

void MAST::LevelSetIntersection::init ( const MAST::FieldFunction< Real > &  phi,
const libMesh::Elem &  e,
const Real  t,
unsigned int  max_elem_id,
unsigned int  max_node_id 
)

Definition at line 223 of file level_set_intersection.cpp.

Here is the call graph for this function:

Here is the caller graph for this function:

unsigned int MAST::LevelSetIntersection::node_on_boundary ( ) const
Returns
the node number on the element if the mode is THROUDH_NODE, otherwise throws an error.

Definition at line 74 of file level_set_intersection.cpp.

Member Data Documentation

unsigned int MAST::LevelSetIntersection::_edge_num_on_boundary
protected

Definition at line 236 of file level_set_intersection.h.

const libMesh::Elem* MAST::LevelSetIntersection::_elem
protected

Definition at line 217 of file level_set_intersection.h.

std::map<const libMesh::Elem*, int> MAST::LevelSetIntersection::_elem_sides_on_interface
protected

Definition at line 242 of file level_set_intersection.h.

bool MAST::LevelSetIntersection::_if_elem_on_negative_phi
protected

true if element is completely on the negative side of level set with no intersection

Definition at line 230 of file level_set_intersection.h.

bool MAST::LevelSetIntersection::_if_elem_on_positive_phi
protected

true if element is completely on the positive side of level set with no intersection

Definition at line 224 of file level_set_intersection.h.

bool MAST::LevelSetIntersection::_initialized
protected

Definition at line 218 of file level_set_intersection.h.

const unsigned int MAST::LevelSetIntersection::_max_elem_divs
protected

Definition at line 216 of file level_set_intersection.h.

unsigned int MAST::LevelSetIntersection::_max_iters
protected

Definition at line 213 of file level_set_intersection.h.

unsigned int MAST::LevelSetIntersection::_max_mesh_elem_id
protected

Definition at line 214 of file level_set_intersection.h.

unsigned int MAST::LevelSetIntersection::_max_mesh_node_id
protected

Definition at line 215 of file level_set_intersection.h.

MAST::LevelSet2DIntersectionMode MAST::LevelSetIntersection::_mode
protected

Definition at line 232 of file level_set_intersection.h.

std::vector<const libMesh::Elem*> MAST::LevelSetIntersection::_negative_phi_elems
protected

Definition at line 240 of file level_set_intersection.h.

std::vector<libMesh::Elem*> MAST::LevelSetIntersection::_new_elems
protected

Definition at line 245 of file level_set_intersection.h.

std::vector<libMesh::Node*> MAST::LevelSetIntersection::_new_nodes
protected

Definition at line 244 of file level_set_intersection.h.

std::map<const libMesh::Node*, libMesh::Point> MAST::LevelSetIntersection::_node_local_coords
protected

Definition at line 246 of file level_set_intersection.h.

unsigned int MAST::LevelSetIntersection::_node_num_on_boundary
protected

Definition at line 234 of file level_set_intersection.h.

std::map<const libMesh::Node*, std::pair<Real, bool> > MAST::LevelSetIntersection::_node_phi_vals
protected

Definition at line 247 of file level_set_intersection.h.

std::vector<const libMesh::Elem*> MAST::LevelSetIntersection::_positive_phi_elems
protected

Definition at line 238 of file level_set_intersection.h.

Real MAST::LevelSetIntersection::_tol
protected

Definition at line 211 of file level_set_intersection.h.


The documentation for this class was generated from the following files: