MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (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-2020 Manav Bhatia and MAST authors
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 
110 
116  bool if_hanging_node(const libMesh::Node* n) const;
117 
122  bool if_elem_has_boundary() const;
123 
129  bool if_intersection_through_elem() const;
130 
137  bool if_elem_on_positive_phi() const;
138 
145  bool if_elem_on_negative_phi() const;
146 
151  bool if_elem_has_positive_phi_region() const;
152 
153 
158  bool if_elem_has_negative_phi_region() const;
159 
160  const std::vector<const libMesh::Elem*>&
162 
163  const std::vector<const libMesh::Elem*>&
165 
166  void
167  get_corner_nodes_on_negative_phi(std::set<const libMesh::Node*>& nodes) const;
168 
174  unsigned int
175  get_side_on_interface(const libMesh::Elem& e) const;
176 
177  bool
178  has_side_on_interface(const libMesh::Elem& e) const;
179 
180  const libMesh::Point&
181  get_nondimensional_coordinate_for_node(const libMesh::Node& n) const;
182 
183 
188  bool if_node_is_new(const libMesh::Node& node) const;
189 
190 
196  bool if_interior_node(const libMesh::Node& node) const;
197 
198 
204  std::pair<const libMesh::Node*, const libMesh::Node*>
205  get_bounding_nodes_for_node(const libMesh::Node& node) const;
206 
207 
212  void
213  get_material_sides_without_intersection(std::set<unsigned int>& sides) const;
214 
215 
216  void get_nearest_intersection_point(const libMesh::Point& p,
217  libMesh::Point& pt);
218 
219  protected:
220 
227  std::unique_ptr<libMesh::Elem>
228  _first_order_elem(const libMesh::Elem& e);
229 
236  const libMesh::Elem& e,
237  const Real t);
238 
239 
241  ( const libMesh::Elem& e,
242  std::vector<std::pair<libMesh::Point, libMesh::Point> >& side_nondim_points,
243  std::map<const libMesh::Node*, libMesh::Point>& node_coord_map);
244 
245  void
247  const libMesh::Elem& e,
248  const Real t,
249  const std::map<const libMesh::Node*, std::pair<Real, bool> >&
250  node_phi_vals);
251 
256  Real
257  _find_intersection_on_straight_edge(const libMesh::Point& p0,
258  const libMesh::Point& p1,
259  const MAST::FieldFunction<Real>& phi,
260  const Real t);
261 
263 
264  unsigned int _max_iters;
265  unsigned int _max_mesh_elem_id;
266  unsigned int _max_mesh_node_id;
267  const unsigned int _max_elem_divs;
268  const libMesh::Elem* _elem;
270 
272 
278 
284 
286 
287  unsigned int _node_num_on_boundary;
288 
289  unsigned int _edge_num_on_boundary;
290 
291  std::vector<const libMesh::Elem*> _positive_phi_elems;
292 
293  std::vector<const libMesh::Elem*> _negative_phi_elems;
294 
295  std::map<const libMesh::Elem*, int> _elem_sides_on_interface;
296 
297  std::vector<libMesh::Node*> _new_nodes;
298  std::vector<libMesh::Elem*> _new_elems;
299  std::map<const libMesh::Node*, libMesh::Point> _node_local_coords;
300  std::map<const libMesh::Node*, std::pair<Real, bool> > _node_phi_vals;
301  std::set<const libMesh::Node*> _interior_nodes;
302  std::map<const libMesh::Node*, std::pair<const libMesh::Node*, const libMesh::Node*>> _bounding_nodes;
303  std::set<const libMesh::Node*> _hanging_node;
304  };
305 
306 }
307 
308 
309 
310 
311 #endif // __mast_level_set_intersection_h__
312 
unsigned int edge_on_boundary() const
MAST::LevelSet2DIntersectionMode _mode
const MAST::FieldFunction< Real > * _phi
const std::vector< const libMesh::Elem * > & get_sub_elems_positive_phi() const
unsigned int get_side_on_interface(const libMesh::Elem &e) const
Real _find_intersection_on_straight_edge(const libMesh::Point &p0, const libMesh::Point &p1, const MAST::FieldFunction< Real > &phi, const Real t)
const libMesh::Elem & elem() const
void get_material_sides_without_intersection(std::set< unsigned int > &sides) const
identifies the sides of the element that are completely on the material side without any intersection...
MAST::LevelSet2DIntersectionMode get_intersection_mode() const
bool _if_elem_on_negative_phi
true if element is completely on the negative side of level set with no intersection ...
Real get_node_phi_value(const libMesh::Node *n) 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_node_is_new(const libMesh::Node &node) const
identifies if the node from the subelements is a new node or an existing node from the parent element...
bool if_hanging_node(const libMesh::Node *n) const
The case of two adjacent edges results in a new node on an edge that is not coincident with the level...
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
std::map< const libMesh::Node *, std::pair< const libMesh::Node *, const libMesh::Node * > > _bounding_nodes
std::pair< const libMesh::Node *, const libMesh::Node * > get_bounding_nodes_for_node(const libMesh::Node &node) const
for new nodes required to create the subelements this method returns the nodes on an edge that bound ...
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)
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
std::vector< const libMesh::Elem * > _positive_phi_elems
std::map< const libMesh::Node *, libMesh::Point > _node_local_coords
std::vector< libMesh::Node * > _new_nodes
bool if_interior_node(const libMesh::Node &node) const
identifies if the new node is on an edge along the level-set method in the interior of the element (a...
bool has_side_on_interface(const libMesh::Elem &e) const
void get_corner_nodes_on_negative_phi(std::set< const libMesh::Node *> &nodes) const
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.
unsigned int node_on_boundary() const
const std::vector< const libMesh::Elem * > & get_sub_elems_negative_phi() const
std::vector< const libMesh::Elem * > _negative_phi_elems
const libMesh::Point & get_nondimensional_coordinate_for_node(const libMesh::Node &n) const
std::set< const libMesh::Node * > _interior_nodes
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 get_nearest_intersection_point(const libMesh::Point &p, libMesh::Point &pt)
void clear()
clears the data structures
std::set< const libMesh::Node * > _hanging_node