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