26 #ifndef GRINS_RAYFIRE_MESH_H 
   27 #define GRINS_RAYFIRE_MESH_H 
   30 #include "libmesh/quadrature.h" 
   31 #include "libmesh/mesh.h" 
   76     RayfireMesh(libMesh::Point& origin, libMesh::Real theta);
 
   86     RayfireMesh(libMesh::Point& origin, libMesh::Real theta, libMesh::Real phi);
 
   92     RayfireMesh(
const GetPot & input, 
const std::string & qoi_string);
 
  107     void init(
const libMesh::MeshBase& mesh_base);
 
  131     void reinit(
const libMesh::MeshBase& mesh_base);
 
  148     libMesh::UniquePtr<libMesh::Mesh> 
_mesh;
 
  172     const libMesh::Elem* 
get_next_elem(
const libMesh::Elem* cur_elem, libMesh::Point& start_point, libMesh::Point& end_point, 
bool same_parent = 
false);
 
  175     bool check_valid_point(libMesh::Point& intersection_point, libMesh::Point& start_point, 
const libMesh::Elem& edge_elem, libMesh::Point& next_point);
 
  178     const libMesh::Elem* 
get_correct_neighbor(libMesh::Point & start_point, libMesh::Point & end_point, 
const libMesh::Elem * cur_elem, 
unsigned int side, 
bool same_parent);
 
  188     const libMesh::Elem* 
get_start_elem(
const libMesh::MeshBase& mesh_base);
 
  191     bool validate_edge(
const libMesh::Point & start_point, 
const libMesh::Point & end_point, 
const libMesh::Elem * side_elem, 
const libMesh::Elem * neighbor);
 
  194     bool rayfire_in_elem(
const libMesh::Point& end_point, 
const libMesh::Elem* elem);
 
  233     bool newton_solve_intersection(libMesh::Point& initial_point, 
const libMesh::Elem* edge_elem, libMesh::Point& intersection_point);
 
  236     void refine(
const libMesh::Elem* main_elem, libMesh::Elem* rayfire_elem);
 
  239     void coarsen(
const libMesh::Elem* rayfire_elem);
 
  244 #endif //GRINS_RAYFIRE_MESH_H 
void coarsen(const libMesh::Elem *rayfire_elem)
Coarsening of a rayfire element whose main mesh counterpart was coarsened. 
const libMesh::Elem * get_next_elem(const libMesh::Elem *cur_elem, libMesh::Point &start_point, libMesh::Point &end_point, bool same_parent=false)
Calculate the intersection point. 
void refine(const libMesh::Elem *main_elem, libMesh::Elem *rayfire_elem)
Refinement of a rayfire element whose main mesh counterpart was refined. 
bool validate_edge(const libMesh::Point &start_point, const libMesh::Point &end_point, const libMesh::Elem *side_elem, const libMesh::Elem *neighbor)
Ensures the rayfire doesn't wander into the middle of an elem. 
const libMesh::Elem * map_to_rayfire_elem(const libMesh::dof_id_type elem_id)
This function takes in an elem_id on the main mesh and returns an elem from the 1D rayfire mesh...
libMesh::UniquePtr< libMesh::Mesh > _mesh
Internal 1D mesh of EDGE2 elements. 
libMesh::Elem * get_rayfire_elem(const libMesh::dof_id_type elem_id)
Private function to get a rayfire elem from main_mesh elem ID. 
bool check_valid_point(libMesh::Point &intersection_point, libMesh::Point &start_point, const libMesh::Elem &edge_elem, libMesh::Point &next_point)
Ensure the calculated intersection point is on the edge_elem and is not the start_point. 
libMesh::Point _origin
Origin point. 
const libMesh::Elem * get_correct_neighbor(libMesh::Point &start_point, libMesh::Point &end_point, const libMesh::Elem *cur_elem, unsigned int side, bool same_parent)
Knowing the end_point, get the appropraite next elem along the path. 
const libMesh::Elem * get_start_elem(const libMesh::MeshBase &mesh_base)
Get the correct starting elem corresponding to _origin. 
void init(const libMesh::MeshBase &mesh_base)
Initialization. 
std::map< libMesh::dof_id_type, libMesh::dof_id_type > _elem_id_map
Map of main mesh elem_id to rayfire mesh elem_id. 
const unsigned int _dim
Dimension of the main mesh. 
void reinit(const libMesh::MeshBase &mesh_base)
Checks for refined and coarsened main mesh elements along the rayfire path. 
bool newton_solve_intersection(libMesh::Point &initial_point, const libMesh::Elem *edge_elem, libMesh::Point &intersection_point)
Iterative solver for calculating the intersection point of the rayfire on the side of an elem...
void elem_ids_in_rayfire(std::vector< libMesh::dof_id_type > &id_vector) const 
Returns a vector of elem IDs that are currently along the rayfire. 
libMesh::Real _phi
Rayfire Spherical azimuthal angle (in radians) 
bool rayfire_in_elem(const libMesh::Point &end_point, const libMesh::Elem *elem)
Walks a short distance along the rayfire and checks if elem contains that point. 
libMesh::Real _theta
Rayfire Spherical polar angle (in radians) 
RayfireMesh()
User should never call the default constructor. 
void check_origin_on_boundary(const libMesh::Elem *start_elem)
Ensure the supplied origin is on a boundary of the mesh.