Flow123d  DF_patch_fe_data_tables-47aef41
eigen_tools.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 eigen_tools.hh
15  * @brief Store finite element data on the actual patch
16  * such as shape function values, gradients, Jacobian
17  * of the mapping from the reference cell etc.
18  * @author David Flanderka
19  */
20 
21 #ifndef PATCH_DATA_TABLE_HH_
22 #define PATCH_DATA_TABLE_HH_
23 
24 #include "system/asserts.hh"
25 #include <armadillo>
26 #include <Eigen/Core>
27 #include <Eigen/Dense>
28 
29 
30 /// Definitions of Eigen structures
31 typedef Eigen::Array<double,Eigen::Dynamic,1> ArrayDbl;
32 typedef Eigen::Array<uint,Eigen::Dynamic,1> ArrayInt;
33 typedef Eigen::Vector<ArrayDbl,Eigen::Dynamic> TableDbl;
34 typedef Eigen::Vector<ArrayInt,Eigen::Dynamic> TableInt;
35 
36 
37 namespace eigen_tools {
38 
39 /// Resize vector of Eigen::Array to given size
40 template<class ET>
41 void resize_table(typename Eigen::Vector<ET,Eigen::Dynamic> &table, uint size) {
42  for (uint i=0; i<table.rows(); ++i) {
43  table(i).resize(size);
44  table(i).setZero(size,1);
45  }
46 }
47 
48 
49 /**
50  * @brief Calculates determinant of a rectangular matrix.
51  */
52 template<class T>
53 ArrayDbl determinant(const T &M);
54 
55 
56 
57 inline Eigen::Matrix<ArrayDbl,1,1> normal_matrix(const Eigen::Matrix<ArrayDbl,1,2> &A) {
58  Eigen::Matrix<ArrayDbl,1,1> res;
59  res(0,0) = A(0,0)*A(0,0)+A(0,1)*A(0,1);
60  return res;
61 }
62 
63 inline Eigen::Matrix<ArrayDbl,1,1> normal_matrix(const Eigen::Matrix<ArrayDbl,2,1> &A) {
64  Eigen::Matrix<ArrayDbl,1,1> res;
65  res(0,0) = A(0,0)*A(0,0)+A(1,0)*A(1,0);
66  return res;
67 }
68 
69 inline Eigen::Matrix<ArrayDbl,1,1> normal_matrix(const Eigen::Matrix<ArrayDbl,1,3> &A) {
70  Eigen::Matrix<ArrayDbl,1,1> res;
71  res(0,0) = A(0,0)*A(0,0)+A(0,1)*A(0,1)+A(0,2)*A(0,2);
72  return res;
73 }
74 
75 inline Eigen::Matrix<ArrayDbl,1,1> normal_matrix(const Eigen::Matrix<ArrayDbl,3,1> &A) {
76  Eigen::Matrix<ArrayDbl,1,1> res;
77  res(0,0) = A(0,0)*A(0,0)+A(1,0)*A(1,0)+A(2,0)*A(2,0);
78  return res;
79 }
80 
81 inline Eigen::Matrix<ArrayDbl,2,2> normal_matrix(const Eigen::Matrix<ArrayDbl,2,3> &A) {
82  Eigen::Matrix<ArrayDbl,2,2> res;
83  res(0,0) = A(0,0)*A(0,0)+A(0,1)*A(0,1)+A(0,2)*A(0,2);
84  res(0,1) = A(0,0)*A(1,0)+A(0,1)*A(1,1)+A(0,2)*A(1,2);
85  res(1,0) = A(1,0)*A(0,0)+A(1,1)*A(0,1)+A(1,2)*A(0,2);
86  res(1,1) = A(1,0)*A(1,0)+A(1,1)*A(1,1)+A(1,2)*A(1,2);
87  return res;
88 }
89 
90 inline Eigen::Matrix<ArrayDbl,2,2> normal_matrix(const Eigen::Matrix<ArrayDbl,3,2> &A) {
91  Eigen::Matrix<ArrayDbl,2,2> res;
92  res(0,0) = A(0,0)*A(0,0)+A(1,0)*A(1,0)+A(2,0)*A(2,0);
93  res(0,1) = A(0,0)*A(0,1)+A(1,0)*A(1,1)+A(2,0)*A(2,1);
94  res(1,0) = A(0,1)*A(0,0)+A(1,1)*A(1,0)+A(2,1)*A(2,0);
95  res(1,1) = A(0,1)*A(0,1)+A(1,1)*A(1,1)+A(2,1)*A(2,1);
96  return res;
97 }
98 
99 
100 
101 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,1,1> &M)
102 {
103  return M(0,0);
104 }
105 
106 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,2,2> &M)
107 {
108  return M(0,0)*M(1,1) - M(1,0)*M(0,1);
109 }
110 
111 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,3,3> &M)
112 {
113  return ( M(0,0)*M(1,1)*M(2,2) + M(0,1)*M(1,2)*M(2,0) + M(0,2)*M(1,0)*M(2,1) )
114  - ( M(2,0)*M(1,1)*M(0,2) + M(2,1)*M(1,2)*M(0,0) + M(2,2)*M(1,0)*M(0,1) );
115 }
116 
117 template<> inline ArrayDbl determinant(FMT_UNUSED const Eigen::Matrix<ArrayDbl,0,3> &M)
118 {
119  return ArrayDbl();
120 }
121 
122 template<> inline ArrayDbl determinant(FMT_UNUSED const Eigen::Matrix<ArrayDbl,3,0> &M)
123 {
124  return ArrayDbl();
125 }
126 
127 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,1,2> &M)
128 {
129  return determinant(normal_matrix(M)).sqrt();
130 }
131 
132 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,2,1> &M)
133 {
134  return determinant(normal_matrix(M)).sqrt();
135 }
136 
137 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,1,3> &M)
138 {
139  return determinant(normal_matrix(M)).sqrt();
140 }
141 
142 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,3,1> &M)
143 {
144  return determinant(normal_matrix(M)).sqrt();
145 }
146 
147 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,2,3> &M)
148 {
149  return determinant(normal_matrix(M)).sqrt();
150 }
151 
152 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,3,2> &M)
153 {
154  return determinant(normal_matrix(M)).sqrt();
155 }
156 
157 
158 /**
159  * @brief Calculates inverse of rectangular matrix or pseudoinverse of non-rectangular matrix.
160  */
161 template<int m, int n>
162 Eigen::Matrix<ArrayDbl,n,m> inverse(const Eigen::Matrix<ArrayDbl,m,n> &A) {
163  // only for cases m > n
164  return inverse(normal_matrix(A)) * A.transpose();
165 }
166 
167 
168 template<> inline Eigen::Matrix<ArrayDbl,1,1> inverse<1,1>(const Eigen::Matrix<ArrayDbl,1,1> &A)
169 {
170  Eigen::Matrix<ArrayDbl,1,1> B;
171  B(0,0) = A(0,0).inverse(); // 1/A(0,0)
172  return B;
173 }
174 
175 template<> inline Eigen::Matrix<ArrayDbl,2,2> inverse<2,2>(const Eigen::Matrix<ArrayDbl,2,2> &A)
176 {
177  Eigen::Matrix<ArrayDbl,2,2> B;
178  ArrayDbl det = determinant(A);
179 
180  B(0,0) = A(1,1) / det;
181  B(0,1) = -A(0,1) / det;
182  B(1,0) = -A(1,0) / det;
183  B(1,1) = A(0,0) / det;
184  return B;
185 }
186 
187 template<> inline Eigen::Matrix<ArrayDbl,3,3> inverse<3,3>(const Eigen::Matrix<ArrayDbl,3,3> &A)
188 {
189  Eigen::Matrix<ArrayDbl,3,3> B;
190 
191  B(0,0) = A(1,1)*A(2,2) - A(2,1)*A(1,2);
192  B(1,0) = A(2,0)*A(1,2) - A(1,0)*A(2,2);
193  B(2,0) = A(1,0)*A(2,1) - A(2,0)*A(1,1);
194 
195  ArrayDbl det = A(0,0)*B(0,0) + A(0,1)*B(1,0) + A(0,2)*B(2,0);
196  B(0,0) = B(0,0) / det;
197  B(1,0) = B(1,0) / det;
198  B(2,0) = B(2,0) / det;
199 
200  B(0,1) = (A(2,1)*A(0,2) - A(0,1)*A(2,2)) / det;
201  B(1,1) = (A(0,0)*A(2,2) - A(2,0)*A(0,2)) / det;
202  B(2,1) = (A(2,0)*A(0,1) - A(0,0)*A(2,1)) / det;
203 
204  B(0,2) = (A(0,1)*A(1,2) - A(1,1)*A(0,2)) / det;
205  B(1,2) = (A(1,0)*A(0,2) - A(0,0)*A(1,2)) / det;
206  B(2,2) = (A(0,0)*A(1,1) - A(1,0)*A(0,1)) / det;
207 
208  return B;
209 }
210 
211 template<> inline Eigen::Matrix<ArrayDbl,2,1> inverse<1,2>(const Eigen::Matrix<ArrayDbl,1,2> &A)
212 {
213  return A.transpose() * inverse(normal_matrix(A));
214 }
215 
216 template<> inline Eigen::Matrix<ArrayDbl,3,1> inverse<1,3>(const Eigen::Matrix<ArrayDbl,1,3> &A)
217 {
218  return A.transpose() * inverse(normal_matrix(A));
219 }
220 
221 template<> inline Eigen::Matrix<ArrayDbl,3,2> inverse<2,3>(const Eigen::Matrix<ArrayDbl,2,3> &A)
222 {
223  return A.transpose() * inverse(normal_matrix(A));
224 }
225 
226 
227 template<unsigned int spacedim, unsigned int dim>
228 Eigen::Matrix<ArrayDbl, spacedim, dim> jacobian(const Eigen::Matrix<ArrayDbl, spacedim, dim+1> &coords)
229 {
230  Eigen::Matrix<ArrayDbl, spacedim, dim> jac;
231  for (unsigned int i=0; i<spacedim; i++)
232  for (unsigned int j=0; j<dim; j++)
233  jac(i,j) = coords(i,j+1) - coords(i,0);
234  return jac;
235 }
236 
237 } // closing namespace eigen_tools
238 
239 //template <unsigned int size>
240 //class VectorCol {
241 //public:
242 // typename Eigen::Matrix<double,size,1> data;
243 //
244 // VectorCol() {}
245 //
246 // VectorCol(int n) {}
247 //
248 // inline double & operator[](std::size_t item) {
249 // return data[item];
250 // }
251 //
252 // inline double & operator()(std::size_t item) {
253 // return data[item];
254 // }
255 //
256 // inline VectorCol<size> operator+(const VectorCol<size> &other) const {
257 // VectorCol<size> res;
258 // res.data = this->data + other.data;
259 // return res;
260 // }
261 //
262 // inline VectorCol<size> operator-(const VectorCol<size> &other) const {
263 // VectorCol<size> res;
264 // res.data = this->data - other.data;
265 // return res;
266 // }
267 //
268 // inline VectorCol<size> operator*(const double &coef) const {
269 // VectorCol<size> res;
270 // res.data = this->data * coef;
271 // return res;
272 // }
273 //
274 // inline VectorCol<size> operator/(const double &coef) const {
275 // VectorCol<size> res;
276 // res.data = this->data / coef;
277 // return res;
278 // }
279 //
280 // inline VectorCol<size> operator*(const VectorCol<size> &other) const {
281 // VectorCol<size> res;
282 // for (unsigned int i=0; i<size; ++i)
283 // res.data[i] = this->data[i] * other.data[i];
284 // return res;
285 // }
286 //
287 // inline VectorCol<size> operator/(const VectorCol<size> &other) const {
288 // VectorCol<size> res;
289 // for (unsigned int i=0; i<size; ++i)
290 // res.data[i] = this->data[i] / other.data[i];
291 // return res;
292 // }
293 //
294 //};
295 
296 
297 #endif /* PATCH_DATA_TABLE_HH_ */
Definitions of ASSERTS.
Eigen::Vector< ArrayDbl, Eigen::Dynamic > TableDbl
Definition: eigen_tools.hh:33
Eigen::Array< uint, Eigen::Dynamic, 1 > ArrayInt
Definition: eigen_tools.hh:32
Eigen::Array< double, Eigen::Dynamic, 1 > ArrayDbl
Definitions of Eigen structures.
Definition: eigen_tools.hh:31
Eigen::Vector< ArrayInt, Eigen::Dynamic > TableInt
Definition: eigen_tools.hh:34
unsigned int uint
Eigen::Matrix< ArrayDbl, n, m > inverse(const Eigen::Matrix< ArrayDbl, m, n > &A)
Calculates inverse of rectangular matrix or pseudoinverse of non-rectangular matrix.
Definition: eigen_tools.hh:162
Eigen::Matrix< ArrayDbl, 2, 1 > inverse< 1, 2 >(const Eigen::Matrix< ArrayDbl, 1, 2 > &A)
Definition: eigen_tools.hh:211
Eigen::Matrix< ArrayDbl, 3, 3 > inverse< 3, 3 >(const Eigen::Matrix< ArrayDbl, 3, 3 > &A)
Definition: eigen_tools.hh:187
Eigen::Matrix< ArrayDbl, 1, 1 > normal_matrix(const Eigen::Matrix< ArrayDbl, 1, 2 > &A)
Definition: eigen_tools.hh:57
Eigen::Matrix< ArrayDbl, spacedim, dim > jacobian(const Eigen::Matrix< ArrayDbl, spacedim, dim+1 > &coords)
Definition: eigen_tools.hh:228
Eigen::Matrix< ArrayDbl, 3, 1 > inverse< 1, 3 >(const Eigen::Matrix< ArrayDbl, 1, 3 > &A)
Definition: eigen_tools.hh:216
Eigen::Matrix< ArrayDbl, 1, 1 > inverse< 1, 1 >(const Eigen::Matrix< ArrayDbl, 1, 1 > &A)
Definition: eigen_tools.hh:168
void resize_table(typename Eigen::Vector< ET, Eigen::Dynamic > &table, uint size)
Resize vector of Eigen::Array to given size.
Definition: eigen_tools.hh:41
Eigen::Matrix< ArrayDbl, 3, 2 > inverse< 2, 3 >(const Eigen::Matrix< ArrayDbl, 2, 3 > &A)
Definition: eigen_tools.hh:221
ArrayDbl determinant(const T &M)
Calculates determinant of a rectangular matrix.
Eigen::Matrix< ArrayDbl, 2, 2 > inverse< 2, 2 >(const Eigen::Matrix< ArrayDbl, 2, 2 > &A)
Definition: eigen_tools.hh:175
#define FMT_UNUSED
Definition: posix.h:75