MAST
level_set_intersection.h
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2019 Manav Bhatia
4  *
5  * This library is free software; you can redistribute it and/or
6  * modify it under the terms of the GNU Lesser General Public
7  * License as published by the Free Software Foundation; either
8  * version 2.1 of the License, or (at your option) any later version.
9  *
10  * This library is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13  * Lesser General Public License for more details.
14  *
15  * You should have received a copy of the GNU Lesser General Public
16  * License along with this library; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 
20 #ifndef __mast_level_set_intersection_h__
21 #define __mast_level_set_intersection_h__
22 
23 // C++ includes
24 #include <vector>
25 
26 // MAST includes
27 #include "base/mast_data_types.h"
28 
29 // libMesh includes
30 #include "libmesh/elem.h"
31 
32 
33 namespace MAST {
34 
36 
37  THROUGH_NODE, // level set passes through node
38  COLINEAR_EDGE, // level set coliniear with edge of element
39  ADJACENT_EDGES, // level set passes through two adjacent edges
40  OPPOSITE_EDGES, // level set passes through two opposite edges
41  OPPOSITE_NODES, // level set passes through diagonally opposite nodes
42  NODE_AND_EDGE, // level set passes through a node and edge
43  TWO_ADJACENT_EDGES, // level set passes through four edges
44  NODE_AND_TWO_EDGES, // level set passes through a node and two edges
46  };
47 
48 
49  // Forward declerations
50  template <typename ValType> class FieldFunction;
51 
52 
53 
55 
56  public:
57 
59 
60  virtual ~LevelSetIntersection();
61 
62  void init(const MAST::FieldFunction<Real>& phi,
63  const libMesh::Elem& e,
64  const Real t,
65  unsigned int max_elem_id,
66  unsigned int max_node_id);
67 
71  void clear();
72 
77  const libMesh::Elem& elem() const;
78 
83  get_intersection_mode() const;
84 
89  unsigned int node_on_boundary() const;
90 
95  unsigned int edge_on_boundary() const;
96 
97 
103 
108  Real get_node_phi_value(const libMesh::Node* n) const;
109 
114  bool if_elem_has_boundary() const;
115 
121  bool if_intersection_through_elem() const;
122 
129  bool if_elem_on_positive_phi() const;
130 
137  bool if_elem_on_negative_phi() const;
138 
143  bool if_elem_has_positive_phi_region() const;
144 
145  const std::vector<const libMesh::Elem*>&
147 
148  const std::vector<const libMesh::Elem*>&
150 
151  void
152  get_corner_nodes_on_negative_phi(std::set<const libMesh::Node*>& nodes) const;
153 
159  unsigned int
160  get_side_on_interface(const libMesh::Elem& e) const;
161 
162  bool
163  has_side_on_interface(const libMesh::Elem& e) const;
164 
165  const libMesh::Point&
166  get_nondimensional_coordinate_for_node(const libMesh::Node& n) const;
167 
168  protected:
169 
176  std::unique_ptr<libMesh::Elem>
177  _first_order_elem(const libMesh::Elem& e);
178 
185  const libMesh::Elem& e,
186  const Real t);
187 
188 
190  ( const libMesh::Elem& e,
191  std::vector<std::pair<libMesh::Point, libMesh::Point> >& side_nondim_points,
192  std::map<const libMesh::Node*, libMesh::Point>& node_coord_map);
193 
194  void
196  const libMesh::Elem& e,
197  const Real t,
198  const std::map<const libMesh::Node*, std::pair<Real, bool> >&
199  node_phi_vals);
200 
205  Real
206  _find_intersection_on_straight_edge(const libMesh::Point& p0,
207  const libMesh::Point& p1,
208  const MAST::FieldFunction<Real>& phi,
209  const Real t);
210 
212 
213  unsigned int _max_iters;
214  unsigned int _max_mesh_elem_id;
215  unsigned int _max_mesh_node_id;
216  const unsigned int _max_elem_divs;
217  const libMesh::Elem* _elem;
219 
225 
231 
233 
234  unsigned int _node_num_on_boundary;
235 
236  unsigned int _edge_num_on_boundary;
237 
238  std::vector<const libMesh::Elem*> _positive_phi_elems;
239 
240  std::vector<const libMesh::Elem*> _negative_phi_elems;
241 
242  std::map<const libMesh::Elem*, int> _elem_sides_on_interface;
243 
244  std::vector<libMesh::Node*> _new_nodes;
245  std::vector<libMesh::Elem*> _new_elems;
246  std::map<const libMesh::Node*, libMesh::Point> _node_local_coords;
247  std::map<const libMesh::Node*, std::pair<Real, bool> > _node_phi_vals;
248  };
249 
250 }
251 
252 
253 
254 
255 #endif // __mast_level_set_intersection_h__
256 
MAST::LevelSet2DIntersectionMode _mode
Real _find_intersection_on_straight_edge(const libMesh::Point &p0, const libMesh::Point &p1, const MAST::FieldFunction< Real > &phi, const Real t)
bool _if_elem_on_negative_phi
true if element is completely on the negative side of level set with no intersection ...
unsigned int edge_on_boundary() const
unsigned int get_side_on_interface(const libMesh::Elem &e) const
std::unique_ptr< libMesh::Elem > _first_order_elem(const libMesh::Elem &e)
creates a first order element from the given high-order element.
bool _if_elem_on_positive_phi
true if element is completely on the positive side of level set with no intersection ...
std::map< const libMesh::Elem *, int > _elem_sides_on_interface
unsigned int node_on_boundary() const
const libMesh::Point & get_nondimensional_coordinate_for_node(const libMesh::Node &n) const
const std::vector< const libMesh::Elem * > & get_sub_elems_negative_phi() const
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)
MAST::LevelSet2DIntersectionMode get_intersection_mode() const
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)
libMesh::Real Real
void init(const MAST::FieldFunction< Real > &phi, const libMesh::Elem &e, const Real t, unsigned int max_elem_id, unsigned int max_node_id)
std::vector< libMesh::Elem * > _new_elems
const libMesh::Elem & elem() const
std::vector< const libMesh::Elem * > _positive_phi_elems
Real get_node_phi_value(const libMesh::Node *n) const
std::map< const libMesh::Node *, libMesh::Point > _node_local_coords
std::vector< libMesh::Node * > _new_nodes
This creates the base class for functions that have a saptial and temporal dependence, and provide sensitivity operations with respect to the functions and parameters.
bool has_side_on_interface(const libMesh::Elem &e) const
void get_corner_nodes_on_negative_phi(std::set< const libMesh::Node * > &nodes) const
const std::vector< const libMesh::Elem * > & get_sub_elems_positive_phi() const
std::vector< const libMesh::Elem * > _negative_phi_elems
std::map< const libMesh::Node *, std::pair< Real, bool > > _node_phi_vals
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...
void clear()
clears the data structures