GRINS-0.8.0
rayfire_mesh.h
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2017 Paul T. Bauman, Roy H. Stogner
7 // Copyright (C) 2010-2013 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 
26 #ifndef GRINS_RAYFIRE_MESH_H
27 #define GRINS_RAYFIRE_MESH_H
28 
29 // libMesh
30 #include "libmesh/quadrature.h"
31 #include "libmesh/mesh.h"
32 
33 // GRINS
34 #include "grins/qoi_base.h"
36 
37 namespace GRINS
38 {
39 
41 
66  {
67  public:
68 
70 
76  RayfireMesh(libMesh::Point& origin, libMesh::Real theta);
77 
79 
86  RayfireMesh(libMesh::Point& origin, libMesh::Real theta, libMesh::Real phi);
87 
89 
92  RayfireMesh(const GetPot & input, const std::string & qoi_string);
93 
95 
98  RayfireMesh(const RayfireMesh & original);
99 
101 
107  void init(const libMesh::MeshBase& mesh_base);
108 
115  const libMesh::Elem* map_to_rayfire_elem(const libMesh::dof_id_type elem_id);
116 
120  void elem_ids_in_rayfire(std::vector<libMesh::dof_id_type>& id_vector) const;
121 
131  void reinit(const libMesh::MeshBase& mesh_base);
132 
133 
134  private:
136  const unsigned int _dim;
137 
139  libMesh::Point _origin;
140 
142  libMesh::Real _theta;
143 
145  libMesh::Real _phi;
146 
148  libMesh::UniquePtr<libMesh::Mesh> _mesh;
149 
151  std::map<libMesh::dof_id_type,libMesh::dof_id_type> _elem_id_map;
152 
154  RayfireMesh();
155 
157 
164  libMesh::Elem* get_rayfire_elem(const libMesh::dof_id_type elem_id);
165 
167 
172  const libMesh::Elem* get_next_elem(const libMesh::Elem* cur_elem, libMesh::Point& start_point, libMesh::Point& end_point, bool same_parent = false);
173 
175  bool check_valid_point(libMesh::Point& intersection_point, libMesh::Point& start_point, const libMesh::Elem& edge_elem, libMesh::Point& next_point);
176 
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);
179 
181  void check_origin_on_boundary(const libMesh::Elem* start_elem);
182 
184 
188  const libMesh::Elem* get_start_elem(const libMesh::MeshBase& mesh_base);
189 
191  bool validate_edge(const libMesh::Point & start_point, const libMesh::Point & end_point, const libMesh::Elem * side_elem, const libMesh::Elem * neighbor);
192 
194  bool rayfire_in_elem(const libMesh::Point& end_point, const libMesh::Elem* elem);
195 
197 
233  bool newton_solve_intersection(libMesh::Point& initial_point, const libMesh::Elem* edge_elem, libMesh::Point& intersection_point);
234 
236  void refine(const libMesh::Elem* main_elem, libMesh::Elem* rayfire_elem);
237 
239  void coarsen(const libMesh::Elem* rayfire_elem);
240 
241  };
242 
243 }
244 #endif //GRINS_RAYFIRE_MESH_H
void coarsen(const libMesh::Elem *rayfire_elem)
Coarsening of a rayfire element whose main mesh counterpart was coarsened.
Definition: rayfire_mesh.C:731
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.
Definition: rayfire_mesh.C:323
void refine(const libMesh::Elem *main_elem, libMesh::Elem *rayfire_elem)
Refinement of a rayfire element whose main mesh counterpart was refined.
Definition: rayfire_mesh.C:618
RayfireMesh.
Definition: rayfire_mesh.h:65
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.
Definition: rayfire_mesh.C:392
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...
Definition: rayfire_mesh.C:181
libMesh::UniquePtr< libMesh::Mesh > _mesh
Internal 1D mesh of EDGE2 elements.
Definition: rayfire_mesh.h:148
libMesh::Elem * get_rayfire_elem(const libMesh::dof_id_type elem_id)
Private function to get a rayfire elem from main_mesh elem ID.
Definition: rayfire_mesh.C:308
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.
Definition: rayfire_mesh.C:355
GRINS namespace.
libMesh::Point _origin
Origin point.
Definition: rayfire_mesh.h:139
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.
Definition: rayfire_mesh.C:432
const libMesh::Elem * get_start_elem(const libMesh::MeshBase &mesh_base)
Get the correct starting elem corresponding to _origin.
Definition: rayfire_mesh.C:275
void init(const libMesh::MeshBase &mesh_base)
Initialization.
Definition: rayfire_mesh.C:104
std::map< libMesh::dof_id_type, libMesh::dof_id_type > _elem_id_map
Map of main mesh elem_id to rayfire mesh elem_id.
Definition: rayfire_mesh.h:151
const unsigned int _dim
Dimension of the main mesh.
Definition: rayfire_mesh.h:136
void reinit(const libMesh::MeshBase &mesh_base)
Checks for refined and coarsened main mesh elements along the rayfire path.
Definition: rayfire_mesh.C:198
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...
Definition: rayfire_mesh.C:541
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.
Definition: rayfire_mesh.C:187
libMesh::Real _phi
Rayfire Spherical azimuthal angle (in radians)
Definition: rayfire_mesh.h:145
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.
Definition: rayfire_mesh.C:372
libMesh::Real _theta
Rayfire Spherical polar angle (in radians)
Definition: rayfire_mesh.h:142
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.
Definition: rayfire_mesh.C:247

Generated on Tue Dec 19 2017 12:47:28 for GRINS-0.8.0 by  doxygen 1.8.9.1