Flow123d  release_2.2.0-914-gf1a3a4f
intersection_local.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_local.hh
15  * @brief Classes representing final intersection objects.
16  * @author Viktor Fris, Pavel Exner
17  *
18  */
19 
20 #ifndef INTERSECTION_LOCAL_H_
21 #define INTERSECTION_LOCAL_H_
22 
23 #include <armadillo>
24 #include <mesh/mesh_types.hh>
25 #include "system/system.hh"
26 
27 
28 
29 //forwward declare
30 template<unsigned int, unsigned int> class IntersectionPointAux;
31 template<unsigned int, unsigned int> class IntersectionAux;
32 
33 template<unsigned int, unsigned int> class IntersectionPoint;
34 template<unsigned int, unsigned int> class IntersectionLocal;
35 template<unsigned int dimA, unsigned int dimB> std::ostream& operator<<(std::ostream& os, const IntersectionLocal<dimA,dimB>& il);
36 template<unsigned int dimA, unsigned int dimB> std::ostream& operator<<(std::ostream& os, const IntersectionPoint<dimA,dimB>& ip);
37 
38 
39 
40 /** @brief Common base for intersection object.
41  *
42  * This base class provides unification of all intersection objects.
43  * The common part is a pair of component element index and bulk element index.
44  *
45  * The derived class @p IntersectionLocal<dimA,dimB> differs in coordinates length
46  * according to @p dimA and @p dimB.
47  *
48  */
50 {
51 protected:
52  /// Index of intersecting element in the component.
53  unsigned int component_element_idx_;
54  /// Index of intersecting element in the bulk.
55  unsigned int bulk_element_idx_;
56 
57 public:
59  /// Constructor taking in element indices.
60  IntersectionLocalBase(unsigned int component_element_idx,
61  unsigned int bulk_element_idx);
63 
64  unsigned int component_ele_idx() const; ///< Returns index of component element.
65  unsigned int bulk_ele_idx() const; ///< Returns index of bulk element.
66 
67  virtual double compute_measure() const =0;
68 };
69 
70 /// First = element index, Second = pointer to intersection object.
71 typedef std::pair<unsigned int, IntersectionLocalBase*> ILpair;
72 
73 
74 
75 inline unsigned int IntersectionLocalBase::component_ele_idx() const
76 { return component_element_idx_; }
77 
78 inline unsigned int IntersectionLocalBase::bulk_ele_idx() const
79 { return bulk_element_idx_; }
80 
81 
82 
83 /** @brief Class represents intersection of two elements.
84  *
85  * It contains indices of intersecting elements (inherited from base class)
86  * and vector of intersection points which provides barycentric coordinates
87  * on both elements.
88  */
89 template<unsigned int dimA, unsigned int dimB>
91 {
92  /// Vector of intersection points.
94 
95 public:
96 
97  /// Default constructor.
99  /// Constructor taking in element indices.
100  IntersectionLocal(unsigned int component_element_idx, unsigned int bulk_element_idx);
101  /// Copy constructor.
103  /// Destructor.
105 
106  ///@name Getters.
107  //@{
108  /// Returns intersection points by a reference.
110 
111  /// Returns intersection points by a constant reference.
112  const std::vector<IntersectionPoint<dimA,dimB>> &points() const;
113 
114  /// Returns intersection point of given @p index.
115  const IntersectionPoint<dimA,dimB> &operator[](unsigned int index) const;
116 
117  unsigned int size() const; ///< Returns number of intersection points.
118  //@}
119 
120  /// Computes the relative measure of intersection object.
121  double compute_measure() const override;
122 
123 
124  /// Friend output operator.
125  friend std::ostream& operator<< <>(std::ostream& os, const IntersectionLocal<dimA,dimB>& intersection);
126 };
127 
128 /********************************************* IMPLEMENTATION ***********************************************/
129 
130 template<unsigned int dimA, unsigned int dimB>
132 { return i_points_; }
133 
134 template<unsigned int dimA, unsigned int dimB>
136 { return i_points_; }
137 
138 template<unsigned int dimA, unsigned int dimB>
140 { ASSERT_DBG(index < i_points_.size());
141  return i_points_[index]; }
142 
143 template<unsigned int dimA, unsigned int dimB>
144 inline unsigned int IntersectionLocal<dimA,dimB>::size() const
145 { return i_points_.size(); }
146 
147 
148 
149 
150 /** @brief Class represents an intersection point of simplex<N> and simplex<M>.
151  * It contains barycentric coordinates of the point on both simplices.
152  */
153 template<unsigned int dimA, unsigned int dimB>
154 class IntersectionPoint {
155 
156  arma::vec::fixed<dimA> comp_coords_; ///< Local coordinates of an IP on simplex<dimA>.
157  arma::vec::fixed<dimB> bulk_coords_; ///< Local coordinates of an IP on simplex<dimB>.
158 
159 public:
160 
161  IntersectionPoint(); ///< Default constructor.
162  ~IntersectionPoint(); ///< Destructor.
163 
164  /// Constructs IP from the auxiliary one coming from computation.
166  /**
167  * Constructor taking local coordinates on simplices as input parameters.
168  * @param comp_coords local coordinates of IP in Simplex<dimA>
169  * @param bulk_coords local coordinates of IP in Simplex<dimB>
170  */
171  IntersectionPoint(const arma::vec::fixed<dimA> &comp_coords, const arma::vec::fixed<dimB> &bulk_coords);
172 
173  ///@name Getters.
174  //@{
175  /// Returns local coordinates in the Simplex<N>.
176  const arma::vec::fixed<dimA> &comp_coords() const;
177 
178  /// Returns local coordinates in the Simplex<M>.
179  const arma::vec::fixed<dimB> &bulk_coords() const;
180  //@}
181 
182  /**
183  * Computes the real coordinates.
184  * comp_ele is component element
185  */
186 
187  arma::vec3 coords(ElementFullIter comp_ele) const;
188 
189  /// Friend output operator.
190  friend std::ostream& operator<< <>(std::ostream& os, const IntersectionPoint<dimA,dimB>& IP);
191 };
192 
193 
194 /********************************************* IMPLEMENTATION ***********************************************/
195 
196 template<unsigned int dimA, unsigned int dimB>
197 const arma::vec::fixed< dimA >& IntersectionPoint<dimA,dimB>::comp_coords() const
198 { return comp_coords_; }
199 
200 template<unsigned int dimA, unsigned int dimB>
201 const arma::vec::fixed< dimB >& IntersectionPoint<dimA,dimB>::bulk_coords() const
202 { return bulk_coords_; }
203 
204 
205 
206 #endif /* INTERSECTION_LOCAL_H_ */
Common base for intersection object.
arma::vec::fixed< dimB > bulk_coords_
Local coordinates of an IP on simplex<dimB>.
Class represents intersection of two elements.
Class represents an intersection point of simplex<N> and simplex<M>. It contains barycentric coordina...
unsigned int size() const
Returns number of intersection points.
unsigned int bulk_ele_idx() const
Returns index of bulk element.
unsigned int component_ele_idx() const
Returns index of component element.
const arma::vec::fixed< dimA > & comp_coords() const
Returns local coordinates in the Simplex<N>.
const arma::vec::fixed< dimB > & bulk_coords() const
Returns local coordinates in the Simplex<M>.
std::pair< unsigned int, IntersectionLocalBase * > ILpair
First = element index, Second = pointer to intersection object.
arma::vec::fixed< dimA > comp_coords_
Local coordinates of an IP on simplex<dimA>.
unsigned int component_element_idx_
Index of intersecting element in the component.
std::vector< IntersectionPoint< dimA, dimB > > i_points_
Vector of intersection points.
#define ASSERT_DBG(expr)
Definition: asserts.hh:349
Internal auxiliary class represents an intersection point of simplex<N> and simplex<M>.
std::vector< IntersectionPoint< dimA, dimB > > & points()
Returns intersection points by a reference.
Internal auxiliary class representing intersection object of simplex<dimA> and simplex<dimB>.
virtual double compute_measure() const =0
const IntersectionPoint< dimA, dimB > & operator[](unsigned int index) const
Returns intersection point of given index.
unsigned int bulk_element_idx_
Index of intersecting element in the bulk.