Flow123d  JS_before_hm-978-gf0793cd
intersection_point_aux.hh
Go to the documentation of this file.
1 /*!
2  *
3 * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file intersection_point_aux.hh
15  * @brief Internal class representing intersection point.
16  * @author Viktor Fris, Pavel Exner
17  *
18  */
19 
20 #ifndef INTERSECTIONPOINT_H_
21 #define INTERSECTIONPOINT_H_
22 
23 #include <armadillo>
24 template <int spacedim> class ElementAccessor;
25 
26 
27 
28 static const double geometry_epsilon = 1e-9;
29 
30 //http://en.cppreference.com/w/cpp/types/numeric_limits/epsilon
31 // static const double rounding_epsilon = 8*std::numeric_limits<double>::epsilon();
32 
33 //forward declare
34 template<unsigned int N, unsigned int M> class IntersectionPointAux;
35 template<unsigned int N, unsigned int M> std::ostream& operator<<(std::ostream& os, const IntersectionPointAux<N,M>& IP);
36 
37 
38 enum class IntersectionResult {
39  positive = 1,
40  negative = 0,
41  degenerate = 2,
42  none = 3
43 };
44 
46  return int(a) < int(b);
47 }
48 
49 /** @brief Internal auxiliary class represents an intersection point of simplex<N> and simplex<M>.
50  *
51  * It contains barycentric coordinates of the point on both simplices.
52  * Further, it contains topology information of the intersection point.
53  * Namely its orientation (according to Plucker coordinates product),
54  * pathologic flag, local indices and dimension of element objects (node/line/triangle).
55  * Refered as 'IP' in further documentation.
56  */
57 template<unsigned int N, unsigned int M>
59 
60  arma::vec::fixed<N+1> local_bcoords_A_; ///< Barycentric coordinates of an IP on simplex<N>.
61  arma::vec::fixed<M+1> local_bcoords_B_; ///< Barycentric coordinates of an IP on simplex<M>.
62 
63  /** @brief Local indices of element objects that intersects.
64  *
65  * Can be indices of point (vertex), lines (edge), triangle (side) in intersecting element,
66  * depending on the dimensions \p dim_A_, resp. \p dim_B_.
67  *
68  * Examples:
69  * dim_A=0, N=2 -> idx_A is local node index of triangle
70  * dim_B=2, N=3 -> idx_B is local side index of tetrahedron
71  * dim_B=1, N=3 -> idx_B is local line index of tetrahedron
72  */
73  //@{
74  unsigned int idx_A_;
75  unsigned int idx_B_;
76  //@}
77 
78  /** @brief Intersection result according to Plucker products.
79  *
80  * All Plucker products > 0 then @b positive.
81  * All Plucker products < 0 then @b negative.
82  * @b degenerate case.
83  * @b none - no intersection.
84  */
86 
87  unsigned int dim_A_; ///< Dimension of the object A of intersection. Equal \p N, by default.
88  unsigned int dim_B_; ///< Dimension of the object B of intersection. Equal \p M, by default.
89 
90 public:
91 
92  IntersectionPointAux(); ///< Default constructor.
93  ~IntersectionPointAux(){}; ///< Destructor.
94 
95  /** @brief Constructor taking barycentric coordinates on simplices as input parameters.
96  * @param lcA barycentric coordinates of IP in Simplex<N>
97  * @param lcB barycentric coordinates of IP in Simplex<M>
98  * @param dim_A dimension of object A
99  * @param dim_B dimension of object B
100  */
101  IntersectionPointAux(const arma::vec::fixed<N+1> &lcA, const arma::vec::fixed<M+1> &lcB,
102  unsigned int dim_A = N, unsigned int dim_B = M);
103 
104  /** @brief Constructor interpolates the second bary coords
105  * of IntersectionPointAux<N,M-1> to IntersectionPointAux<N,M>
106  *
107  * Allowed only from dimension \p M 1 to 2 and from 2 to 3.
108  * @param IP intersection point of lower dimension of object B
109  * @param idx_B is the index of object B of IntersectionPointAux<N,M-1> in object B of IntersectionPointAux<N,M>
110  */
111  IntersectionPointAux(const IntersectionPointAux<N,M-1> &IP, unsigned int idx_B);
112 
113  /** @brief Constructor interpolates the second bary coords
114  * of IntersectionPointAux<N,M-2> to IntersectionPointAux<N,M>
115  *
116  * Allowed only from dimension \p M 1 to 3.
117  * @param IP intersection point of lower dimension of object B
118  * @param idx_B is the index of object B of IntersectionPointAux<N,M-2> in object B of IntersectionPointAux<N,M>
119  */
120  IntersectionPointAux(const IntersectionPointAux<N,M-2> &IP, unsigned int idx_B);
121 
122  /// Resets the object to default values.
123  void clear();
124 
125  /// Switches the object A and B.
126  IntersectionPointAux<M, N> switch_objects() const;
127 
128  ///@name Setters.
129  //@{
130  /** @brief Setter for coordinates.
131  * @param lcA braycentric coordinates of A
132  * @param lcB braycentric coordinates of B
133  */
134  void set_coordinates(const arma::vec::fixed<N+1> &lcA, const arma::vec::fixed<M+1> &lcB);
135 
136  /// Setter for topology data. See description of @p idx_A_ and @p dim_A_.
137  void set_topology(unsigned int idx_A, unsigned int dim_A,
138  unsigned int idx_B, unsigned int dim_B);
139 
140  void set_topology_A(unsigned int idx, unsigned int dim_A); ///< Sets the topology of object A in Simplex<N>.
141  void set_topology_B(unsigned int idx, unsigned int dim_B); ///< Sets the topology of object B in Simplex<M>.
142  void set_result(IntersectionResult result); ///< Setter orientation flag.
143 
144  //@}
145 
146 
147  ///@name Getters.
148  //@{
149  /// Returns barycentric coordinates in the Simplex<N>.
150  const arma::vec::fixed<N+1> &local_bcoords_A() const;
151 
152  /// Returns barycentric coordinates in the Simplex<M>.
153  const arma::vec::fixed<M+1> &local_bcoords_B() const;
154 
155  unsigned int dim_A() const; ///< Returns dimension of object A.
156  unsigned int dim_B() const; ///< Returns dimension of object B.
157  unsigned int idx_A() const; ///< Returns the index of Simplex<N>.
158  unsigned int idx_B() const; ///< Returns the index of Simplex<M>.
159  /// Result: 0 - negative sign, 1 - positive sign, 2 - degenerate (zero for all sides), 3 - none
160  IntersectionResult result() const; ///< Returns the orientation.
161  //@}
162 
163  /// Computes real coordinates of IP, given the element @p ele in which IP lies.
164  arma::vec::fixed<3> coords(ElementAccessor<3> ele) const;
165 
166  /// Returns true, if @p other intersection point has the same topology.
167  bool topology_equal(const IntersectionPointAux<N,M> &other) const;
168 
169  /** @brief Comparison operator for sorting the IPs in convex hull tracing algorithm.
170  * Compares the points by x-coordinate (in case of a tie, compares by y-coordinate).
171  */
172  //bool operator<(const IntersectionPointAux<N,M> &ip) const;
173 
174  /// Friend output operator.
175  friend std::ostream& operator<< <>(std::ostream& os, const IntersectionPointAux<N,M>& IP);
176 };
177 
178 
179 /********************************************* IMPLEMENTATION ***********************************************/
180 
181 template<unsigned int N, unsigned int M>
182 void IntersectionPointAux<N,M>::set_coordinates(const arma::vec::fixed< N + 1 >& lcA, const arma::vec::fixed< M + 1 >& lcB)
183 { local_bcoords_A_ = lcA;
184  local_bcoords_B_ = lcB; }
185 
186 template<unsigned int N, unsigned int M>
187 void IntersectionPointAux<N,M>::set_topology(unsigned int idx_A,unsigned int dim_A, unsigned int idx_B, unsigned int dim_B)
188 { idx_A_ = idx_A;
189  idx_B_ = idx_B;
190  dim_A_ = dim_A;
191  dim_B_ = dim_B;
192 }
193 
194 template<unsigned int N, unsigned int M>
196 { result_ = res; }
197 
198 template<unsigned int N, unsigned int M>
200 { return dim_A_; }
201 
202 template<unsigned int N, unsigned int M>
204 { return dim_B_; }
205 
206 template<unsigned int N, unsigned int M>
207 const arma::vec::fixed< N + 1 >& IntersectionPointAux<N,M>::local_bcoords_A() const
208 { return local_bcoords_A_; }
209 
210 template<unsigned int N, unsigned int M>
211 const arma::vec::fixed< M + 1 >& IntersectionPointAux<N,M>::local_bcoords_B() const
212 { return local_bcoords_B_; }
213 
214 template<unsigned int N, unsigned int M>
215 void IntersectionPointAux<N,M>::set_topology_A(unsigned int idx_A, unsigned int dim_A)
216 { idx_A_ = idx_A;
217  dim_A_ = dim_A; }
218 
219 template<unsigned int N, unsigned int M>
220 void IntersectionPointAux<N,M>::set_topology_B(unsigned int idx_B, unsigned int dim_B)
221 { idx_B_ = idx_B;
222  dim_B_ = dim_B; }
223 
224 template<unsigned int N, unsigned int M>
226 { return idx_A_; }
227 
228 template<unsigned int N, unsigned int M>
230 { return idx_B_; }
231 
232 template<unsigned int N, unsigned int M>
234 { return result_; }
235 
236 #endif /* INTERSECTIONPOINT_H_ */
unsigned int dim_A() const
Returns dimension of object A.
void set_topology_B(unsigned int idx, unsigned int dim_B)
Sets the topology of object B in Simplex<M>.
IntersectionResult result() const
Result: 0 - negative sign, 1 - positive sign, 2 - degenerate (zero for all sides), 3 - none.
arma::vec::fixed< N+1 > local_bcoords_A_
Barycentric coordinates of an IP on simplex<N>.
IntersectionResult result_
Intersection result according to Plucker products.
void set_topology_A(unsigned int idx, unsigned int dim_A)
Sets the topology of object A in Simplex<N>.
arma::vec::fixed< M+1 > local_bcoords_B_
Barycentric coordinates of an IP on simplex<M>.
static const double geometry_epsilon
const arma::vec::fixed< N+1 > & local_bcoords_A() const
Returns barycentric coordinates in the Simplex<N>.
const arma::vec::fixed< M+1 > & local_bcoords_B() const
Returns barycentric coordinates in the Simplex<M>.
unsigned int idx_A() const
Returns the index of Simplex<N>.
unsigned int dim_B() const
Returns dimension of object B.
void set_coordinates(const arma::vec::fixed< N+1 > &lcA, const arma::vec::fixed< M+1 > &lcB)
Setter for coordinates.
unsigned int idx_A_
Local indices of element objects that intersects.
unsigned int dim_A_
Dimension of the object A of intersection. Equal N, by default.
unsigned int dim_B_
Dimension of the object B of intersection. Equal M, by default.
Internal auxiliary class represents an intersection point of simplex<N> and simplex<M>.
unsigned int idx_B() const
void set_result(IntersectionResult result)
Setter orientation flag.
bool operator<(IntersectionResult a, IntersectionResult b)
void set_topology(unsigned int idx_A, unsigned int dim_A, unsigned int idx_B, unsigned int dim_B)
Setter for topology data. See description of idx_A_ and dim_A_.