Flow123d  DF_patch_fe_data_tables-dc3adc1
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 "system/arena_resource.hh"
26 #include <armadillo>
27 #include <Eigen/Core>
28 #include <Eigen/Dense>
29 #include <unsupported/Eigen/MatrixFunctions>
30 
31 
32 /// Definitions of Eigen structures
33 typedef Eigen::Array<double,Eigen::Dynamic,1> ArrayDbl;
34 typedef Eigen::Array<uint,Eigen::Dynamic,1> ArrayInt;
35 typedef Eigen::Vector<ArrayDbl,Eigen::Dynamic> TableDbl;
36 typedef Eigen::Vector<ArrayInt,Eigen::Dynamic> TableInt;
37 
38 
39 namespace eigen_tools {
40 
41 /********************************************
42  * Matrix operations on ArrayDbl type *
43  ********************************************/
44 
45 /**
46  * @brief Calculates determinant of a rectangular matrix.
47  */
48 template<int m, int n>
49 ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,m,n> &A);
50 
51 
52 
53 inline Eigen::Matrix<ArrayDbl,1,1> normal_matrix(const Eigen::Matrix<ArrayDbl,1,2> &A) {
54  Eigen::Matrix<ArrayDbl,1,1> res;
55  res(0,0) = A(0,0)*A(0,0)+A(0,1)*A(0,1);
56  return res;
57 }
58 
59 inline Eigen::Matrix<ArrayDbl,1,1> normal_matrix(const Eigen::Matrix<ArrayDbl,2,1> &A) {
60  Eigen::Matrix<ArrayDbl,1,1> res;
61  res(0,0) = A(0,0)*A(0,0)+A(1,0)*A(1,0);
62  return res;
63 }
64 
65 inline Eigen::Matrix<ArrayDbl,1,1> normal_matrix(const Eigen::Matrix<ArrayDbl,1,3> &A) {
66  Eigen::Matrix<ArrayDbl,1,1> res;
67  res(0,0) = A(0,0)*A(0,0)+A(0,1)*A(0,1)+A(0,2)*A(0,2);
68  return res;
69 }
70 
71 inline Eigen::Matrix<ArrayDbl,1,1> normal_matrix(const Eigen::Matrix<ArrayDbl,3,1> &A) {
72  Eigen::Matrix<ArrayDbl,1,1> res;
73  res(0,0) = A(0,0)*A(0,0)+A(1,0)*A(1,0)+A(2,0)*A(2,0);
74  return res;
75 }
76 
77 inline Eigen::Matrix<ArrayDbl,2,2> normal_matrix(const Eigen::Matrix<ArrayDbl,2,3> &A) {
78  Eigen::Matrix<ArrayDbl,2,2> res;
79  res(0,0) = A(0,0)*A(0,0)+A(0,1)*A(0,1)+A(0,2)*A(0,2);
80  res(0,1) = A(0,0)*A(1,0)+A(0,1)*A(1,1)+A(0,2)*A(1,2);
81  res(1,0) = A(1,0)*A(0,0)+A(1,1)*A(0,1)+A(1,2)*A(0,2);
82  res(1,1) = A(1,0)*A(1,0)+A(1,1)*A(1,1)+A(1,2)*A(1,2);
83  return res;
84 }
85 
86 inline Eigen::Matrix<ArrayDbl,2,2> normal_matrix(const Eigen::Matrix<ArrayDbl,3,2> &A) {
87  Eigen::Matrix<ArrayDbl,2,2> res;
88  res(0,0) = A(0,0)*A(0,0)+A(1,0)*A(1,0)+A(2,0)*A(2,0);
89  res(0,1) = A(0,0)*A(0,1)+A(1,0)*A(1,1)+A(2,0)*A(2,1);
90  res(1,0) = A(0,1)*A(0,0)+A(1,1)*A(1,0)+A(2,1)*A(2,0);
91  res(1,1) = A(0,1)*A(0,1)+A(1,1)*A(1,1)+A(2,1)*A(2,1);
92  return res;
93 }
94 
95 
96 
97 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,1,1> &M)
98 {
99  return M(0,0);
100 }
101 
102 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,2,2> &M)
103 {
104  return M(0,0)*M(1,1) - M(1,0)*M(0,1);
105 }
106 
107 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,3,3> &M)
108 {
109  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) )
110  - ( 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) );
111 }
112 
113 template<> inline ArrayDbl determinant(FMT_UNUSED const Eigen::Matrix<ArrayDbl,0,3> &M)
114 {
115  return ArrayDbl();
116 }
117 
118 template<> inline ArrayDbl determinant(FMT_UNUSED const Eigen::Matrix<ArrayDbl,3,0> &M)
119 {
120  return ArrayDbl();
121 }
122 
123 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,1,2> &M)
124 {
125  return determinant(normal_matrix(M)).sqrt();
126 }
127 
128 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,2,1> &M)
129 {
130  return determinant(normal_matrix(M)).sqrt();
131 }
132 
133 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,1,3> &M)
134 {
135  return determinant(normal_matrix(M)).sqrt();
136 }
137 
138 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,3,1> &M)
139 {
140  return determinant(normal_matrix(M)).sqrt();
141 }
142 
143 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,2,3> &M)
144 {
145  return determinant(normal_matrix(M)).sqrt();
146 }
147 
148 template<> inline ArrayDbl determinant(const Eigen::Matrix<ArrayDbl,3,2> &M)
149 {
150  return determinant(normal_matrix(M)).sqrt();
151 }
152 
153 
154 /**
155  * @brief Calculates inverse of rectangular matrix or pseudoinverse of non-rectangular matrix.
156  */
157 template<int m, int n>
158 Eigen::Matrix<ArrayDbl,n,m> inverse(const Eigen::Matrix<ArrayDbl,m,n> &A) {
159  // only for cases m > n
160  return inverse(normal_matrix(A)) * A.transpose();
161 }
162 
163 
164 template<> inline Eigen::Matrix<ArrayDbl,1,1> inverse<1,1>(const Eigen::Matrix<ArrayDbl,1,1> &A)
165 {
166  Eigen::Matrix<ArrayDbl,1,1> B;
167  B(0,0) = A(0,0).inverse(); // 1/A(0,0)
168  return B;
169 }
170 
171 template<> inline Eigen::Matrix<ArrayDbl,2,2> inverse<2,2>(const Eigen::Matrix<ArrayDbl,2,2> &A)
172 {
173  Eigen::Matrix<ArrayDbl,2,2> B;
174  ArrayDbl det = determinant(A);
175 
176  B(0,0) = A(1,1) / det;
177  B(0,1) = A(0,1) / det * (-1);
178  B(1,0) = A(1,0) / det * (-1);
179  B(1,1) = A(0,0) / det;
180  return B;
181 }
182 
183 template<> inline Eigen::Matrix<ArrayDbl,3,3> inverse<3,3>(const Eigen::Matrix<ArrayDbl,3,3> &A)
184 {
185  Eigen::Matrix<ArrayDbl,3,3> B;
186 
187  B(0,0) = A(1,1)*A(2,2) - A(2,1)*A(1,2);
188  B(1,0) = A(2,0)*A(1,2) - A(1,0)*A(2,2);
189  B(2,0) = A(1,0)*A(2,1) - A(2,0)*A(1,1);
190 
191  ArrayDbl det = A(0,0)*B(0,0) + A(0,1)*B(1,0) + A(0,2)*B(2,0);
192  B(0,0) = B(0,0) / det;
193  B(1,0) = B(1,0) / det;
194  B(2,0) = B(2,0) / det;
195 
196  B(0,1) = (A(2,1)*A(0,2) - A(0,1)*A(2,2)) / det;
197  B(1,1) = (A(0,0)*A(2,2) - A(2,0)*A(0,2)) / det;
198  B(2,1) = (A(2,0)*A(0,1) - A(0,0)*A(2,1)) / det;
199 
200  B(0,2) = (A(0,1)*A(1,2) - A(1,1)*A(0,2)) / det;
201  B(1,2) = (A(1,0)*A(0,2) - A(0,0)*A(1,2)) / det;
202  B(2,2) = (A(0,0)*A(1,1) - A(1,0)*A(0,1)) / det;
203 
204  return B;
205 }
206 
207 template<> inline Eigen::Matrix<ArrayDbl,2,1> inverse<1,2>(const Eigen::Matrix<ArrayDbl,1,2> &A)
208 {
209  return A.transpose() * inverse(normal_matrix(A));
210 }
211 
212 template<> inline Eigen::Matrix<ArrayDbl,3,1> inverse<1,3>(const Eigen::Matrix<ArrayDbl,1,3> &A)
213 {
214  return A.transpose() * inverse(normal_matrix(A));
215 }
216 
217 template<> inline Eigen::Matrix<ArrayDbl,3,2> inverse<2,3>(const Eigen::Matrix<ArrayDbl,2,3> &A)
218 {
219  return A.transpose() * inverse(normal_matrix(A));
220 }
221 
222 
223 /// Resize vector of Eigen::Array to given size
224 template<class ET>
225 void resize_table(typename Eigen::Vector<ET,Eigen::Dynamic> &table, uint size) {
226  for (uint i=0; i<table.rows(); ++i) {
227  table(i).resize(size);
228  table(i).setZero(size,1);
229  }
230 }
231 
232 
233 template<unsigned int spacedim, unsigned int dim>
234 Eigen::Matrix<ArrayDbl, spacedim, dim> jacobian(const Eigen::Matrix<ArrayDbl, spacedim, dim+1> &coords)
235 {
236  Eigen::Matrix<ArrayDbl, spacedim, dim> jac;
237  for (unsigned int i=0; i<spacedim; i++)
238  for (unsigned int j=0; j<dim; j++)
239  jac(i,j) = coords(i,j+1) - coords(i,0);
240  return jac;
241 }
242 
243 }
244 
245 namespace eigen_arena_tools {
246 /********************************************
247  * Matrix operations on ArenaVec type *
248  ********************************************/
249 
250 /**
251  * @brief Calculates determinant of a rectangular matrix.
252  */
253 template<int m, int n>
254 ArenaVec<double> determinant(const Eigen::Matrix<ArenaVec<double>,m,n> &A);
255 
256 
257 
258 inline Eigen::Matrix<ArenaVec<double>,1,1> normal_matrix(const Eigen::Matrix<ArenaVec<double>,1,2> &A) {
259  Eigen::Matrix<ArenaVec<double>,1,1> res;
260  res(0,0) = A(0,0)*A(0,0)+A(0,1)*A(0,1);
261  return res;
262 }
263 
264 inline Eigen::Matrix<ArenaVec<double>,1,1> normal_matrix(const Eigen::Matrix<ArenaVec<double>,2,1> &A) {
265  Eigen::Matrix<ArenaVec<double>,1,1> res;
266  res(0,0) = A(0,0)*A(0,0)+A(1,0)*A(1,0);
267  return res;
268 }
269 
270 inline Eigen::Matrix<ArenaVec<double>,1,1> normal_matrix(const Eigen::Matrix<ArenaVec<double>,1,3> &A) {
271  Eigen::Matrix<ArenaVec<double>,1,1> res;
272  res(0,0) = A(0,0)*A(0,0)+A(0,1)*A(0,1)+A(0,2)*A(0,2);
273  return res;
274 }
275 
276 inline Eigen::Matrix<ArenaVec<double>,1,1> normal_matrix(const Eigen::Matrix<ArenaVec<double>,3,1> &A) {
277  Eigen::Matrix<ArenaVec<double>,1,1> res;
278  res(0,0) = A(0,0)*A(0,0)+A(1,0)*A(1,0)+A(2,0)*A(2,0);
279  return res;
280 }
281 
282 inline Eigen::Matrix<ArenaVec<double>,2,2> normal_matrix(const Eigen::Matrix<ArenaVec<double>,2,3> &A) {
283  Eigen::Matrix<ArenaVec<double>,2,2> res;
284  res(0,0) = A(0,0)*A(0,0)+A(0,1)*A(0,1)+A(0,2)*A(0,2);
285  res(0,1) = A(0,0)*A(1,0)+A(0,1)*A(1,1)+A(0,2)*A(1,2);
286  res(1,0) = A(1,0)*A(0,0)+A(1,1)*A(0,1)+A(1,2)*A(0,2);
287  res(1,1) = A(1,0)*A(1,0)+A(1,1)*A(1,1)+A(1,2)*A(1,2);
288  return res;
289 }
290 
291 inline Eigen::Matrix<ArenaVec<double>,2,2> normal_matrix(const Eigen::Matrix<ArenaVec<double>,3,2> &A) {
292  Eigen::Matrix<ArenaVec<double>,2,2> res;
293  res(0,0) = A(0,0)*A(0,0)+A(1,0)*A(1,0)+A(2,0)*A(2,0);
294  res(0,1) = A(0,0)*A(0,1)+A(1,0)*A(1,1)+A(2,0)*A(2,1);
295  res(1,0) = A(0,1)*A(0,0)+A(1,1)*A(1,0)+A(2,1)*A(2,0);
296  res(1,1) = A(0,1)*A(0,1)+A(1,1)*A(1,1)+A(2,1)*A(2,1);
297  return res;
298 }
299 
300 
301 
302 template<> inline ArenaVec<double> determinant(const Eigen::Matrix<ArenaVec<double>,1,1> &M)
303 {
304  return M(0,0);
305 }
306 
307 template<> inline ArenaVec<double> determinant(const Eigen::Matrix<ArenaVec<double>,2,2> &M)
308 {
309  return M(0,0)*M(1,1) - M(1,0)*M(0,1);
310 }
311 
312 template<> inline ArenaVec<double> determinant(const Eigen::Matrix<ArenaVec<double>,3,3> &M)
313 {
314  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) )
315  - ( 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) );
316 }
317 
318 template<> inline ArenaVec<double> determinant(FMT_UNUSED const Eigen::Matrix<ArenaVec<double>,0,3> &M)
319 {
320  return ArenaVec<double>();
321 }
322 
323 template<> inline ArenaVec<double> determinant(FMT_UNUSED const Eigen::Matrix<ArenaVec<double>,3,0> &M)
324 {
325  return ArenaVec<double>();
326 }
327 
328 template<> inline ArenaVec<double> determinant(const Eigen::Matrix<ArenaVec<double>,1,2> &M)
329 {
330  return determinant(normal_matrix(M)).sqrt();
331 }
332 
333 template<> inline ArenaVec<double> determinant(const Eigen::Matrix<ArenaVec<double>,2,1> &M)
334 {
335  return determinant(normal_matrix(M)).sqrt();
336 }
337 
338 template<> inline ArenaVec<double> determinant(const Eigen::Matrix<ArenaVec<double>,1,3> &M)
339 {
340  return determinant(normal_matrix(M)).sqrt();
341 }
342 
343 template<> inline ArenaVec<double> determinant(const Eigen::Matrix<ArenaVec<double>,3,1> &M)
344 {
345  return determinant(normal_matrix(M)).sqrt();
346 }
347 
348 template<> inline ArenaVec<double> determinant(const Eigen::Matrix<ArenaVec<double>,2,3> &M)
349 {
350  return determinant(normal_matrix(M)).sqrt();
351 }
352 
353 template<> inline ArenaVec<double> determinant(const Eigen::Matrix<ArenaVec<double>,3,2> &M)
354 {
355  return determinant(normal_matrix(M)).sqrt();
356 }
357 
358 
359 /**
360  * @brief Calculates inverse of rectangular matrix or pseudoinverse of non-rectangular matrix.
361  */
362 template<int m, int n>
363 Eigen::Matrix<ArenaVec<double>,n,m> inverse(const Eigen::Matrix<ArenaVec<double>,m,n> &A) {
364  // only for cases m > n
365  return inverse(normal_matrix(A)) * A.transpose();
366 }
367 
368 
369 template<> inline Eigen::Matrix<ArenaVec<double>,1,1> inverse<1,1>(const Eigen::Matrix<ArenaVec<double>,1,1> &A)
370 {
371  Eigen::Matrix<ArenaVec<double>,1,1> B;
372  B(0,0) = A(0,0).inverse(); // 1/A(0,0)
373  return B;
374 }
375 
376 template<> inline Eigen::Matrix<ArenaVec<double>,2,2> inverse<2,2>(const Eigen::Matrix<ArenaVec<double>,2,2> &A)
377 {
378  Eigen::Matrix<ArenaVec<double>,2,2> B;
379  ArenaVec<double> det = determinant(A);
380 
381  B(0,0) = A(1,1) / det;
382  B(0,1) = A(0,1) / det * (-1);
383  B(1,0) = A(1,0) / det * (-1);
384  B(1,1) = A(0,0) / det;
385  return B;
386 }
387 
388 template<> inline Eigen::Matrix<ArenaVec<double>,3,3> inverse<3,3>(const Eigen::Matrix<ArenaVec<double>,3,3> &A)
389 {
390  Eigen::Matrix<ArenaVec<double>,3,3> B;
391 
392  B(0,0) = A(1,1)*A(2,2) - A(2,1)*A(1,2);
393  B(1,0) = A(2,0)*A(1,2) - A(1,0)*A(2,2);
394  B(2,0) = A(1,0)*A(2,1) - A(2,0)*A(1,1);
395 
396  ArenaVec<double> det = A(0,0)*B(0,0) + A(0,1)*B(1,0) + A(0,2)*B(2,0);
397  B(0,0) = B(0,0) / det;
398  B(1,0) = B(1,0) / det;
399  B(2,0) = B(2,0) / det;
400 
401  B(0,1) = (A(2,1)*A(0,2) - A(0,1)*A(2,2)) / det;
402  B(1,1) = (A(0,0)*A(2,2) - A(2,0)*A(0,2)) / det;
403  B(2,1) = (A(2,0)*A(0,1) - A(0,0)*A(2,1)) / det;
404 
405  B(0,2) = (A(0,1)*A(1,2) - A(1,1)*A(0,2)) / det;
406  B(1,2) = (A(1,0)*A(0,2) - A(0,0)*A(1,2)) / det;
407  B(2,2) = (A(0,0)*A(1,1) - A(1,0)*A(0,1)) / det;
408 
409  return B;
410 }
411 
412 template<> inline Eigen::Matrix<ArenaVec<double>,2,1> inverse<1,2>(const Eigen::Matrix<ArenaVec<double>,1,2> &A)
413 {
414  return A.transpose() * inverse(normal_matrix(A));
415 }
416 
417 template<> inline Eigen::Matrix<ArenaVec<double>,3,1> inverse<1,3>(const Eigen::Matrix<ArenaVec<double>,1,3> &A)
418 {
419  return A.transpose() * inverse(normal_matrix(A));
420 }
421 
422 template<> inline Eigen::Matrix<ArenaVec<double>,3,2> inverse<2,3>(const Eigen::Matrix<ArenaVec<double>,2,3> &A)
423 {
424  return A.transpose() * inverse(normal_matrix(A));
425 }
426 
427 
428 } // closing namespace eigen_tools
429 
430 //template <unsigned int size>
431 //class VectorCol {
432 //public:
433 // typename Eigen::Matrix<double,size,1> data;
434 //
435 // VectorCol() {}
436 //
437 // VectorCol(int n) {}
438 //
439 // inline double & operator[](std::size_t item) {
440 // return data[item];
441 // }
442 //
443 // inline double & operator()(std::size_t item) {
444 // return data[item];
445 // }
446 //
447 // inline VectorCol<size> operator+(const VectorCol<size> &other) const {
448 // VectorCol<size> res;
449 // res.data = this->data + other.data;
450 // return res;
451 // }
452 //
453 // inline VectorCol<size> operator-(const VectorCol<size> &other) const {
454 // VectorCol<size> res;
455 // res.data = this->data - other.data;
456 // return res;
457 // }
458 //
459 // inline VectorCol<size> operator*(const double &coef) const {
460 // VectorCol<size> res;
461 // res.data = this->data * coef;
462 // return res;
463 // }
464 //
465 // inline VectorCol<size> operator/(const double &coef) const {
466 // VectorCol<size> res;
467 // res.data = this->data / coef;
468 // return res;
469 // }
470 //
471 // inline VectorCol<size> operator*(const VectorCol<size> &other) const {
472 // VectorCol<size> res;
473 // for (unsigned int i=0; i<size; ++i)
474 // res.data[i] = this->data[i] * other.data[i];
475 // return res;
476 // }
477 //
478 // inline VectorCol<size> operator/(const VectorCol<size> &other) const {
479 // VectorCol<size> res;
480 // for (unsigned int i=0; i<size; ++i)
481 // res.data[i] = this->data[i] / other.data[i];
482 // return res;
483 // }
484 //
485 //};
486 
487 
488 #endif /* PATCH_DATA_TABLE_HH_ */
Definitions of ASSERTS.
Eigen::Vector< ArrayDbl, Eigen::Dynamic > TableDbl
Definition: eigen_tools.hh:35
Eigen::Array< uint, Eigen::Dynamic, 1 > ArrayInt
Definition: eigen_tools.hh:34
Eigen::Array< double, Eigen::Dynamic, 1 > ArrayDbl
Definitions of Eigen structures.
Definition: eigen_tools.hh:33
Eigen::Vector< ArrayInt, Eigen::Dynamic > TableInt
Definition: eigen_tools.hh:36
unsigned int uint
Eigen::Matrix< ArenaVec< double >, n, m > inverse(const Eigen::Matrix< ArenaVec< double >, m, n > &A)
Calculates inverse of rectangular matrix or pseudoinverse of non-rectangular matrix.
Definition: eigen_tools.hh:363
Eigen::Matrix< ArenaVec< double >, 3, 2 > inverse< 2, 3 >(const Eigen::Matrix< ArenaVec< double >, 2, 3 > &A)
Definition: eigen_tools.hh:422
Eigen::Matrix< ArenaVec< double >, 1, 1 > normal_matrix(const Eigen::Matrix< ArenaVec< double >, 1, 2 > &A)
Definition: eigen_tools.hh:258
Eigen::Matrix< ArenaVec< double >, 1, 1 > inverse< 1, 1 >(const Eigen::Matrix< ArenaVec< double >, 1, 1 > &A)
Definition: eigen_tools.hh:369
Eigen::Matrix< ArenaVec< double >, 2, 2 > inverse< 2, 2 >(const Eigen::Matrix< ArenaVec< double >, 2, 2 > &A)
Definition: eigen_tools.hh:376
ArenaVec< double > determinant(const Eigen::Matrix< ArenaVec< double >, m, n > &A)
Calculates determinant of a rectangular matrix.
Eigen::Matrix< ArenaVec< double >, 3, 3 > inverse< 3, 3 >(const Eigen::Matrix< ArenaVec< double >, 3, 3 > &A)
Definition: eigen_tools.hh:388
Eigen::Matrix< ArenaVec< double >, 3, 1 > inverse< 1, 3 >(const Eigen::Matrix< ArenaVec< double >, 1, 3 > &A)
Definition: eigen_tools.hh:417
Eigen::Matrix< ArenaVec< double >, 2, 1 > inverse< 1, 2 >(const Eigen::Matrix< ArenaVec< double >, 1, 2 > &A)
Definition: eigen_tools.hh:412
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:158
Eigen::Matrix< ArrayDbl, 2, 1 > inverse< 1, 2 >(const Eigen::Matrix< ArrayDbl, 1, 2 > &A)
Definition: eigen_tools.hh:207
Eigen::Matrix< ArrayDbl, 3, 3 > inverse< 3, 3 >(const Eigen::Matrix< ArrayDbl, 3, 3 > &A)
Definition: eigen_tools.hh:183
Eigen::Matrix< ArrayDbl, 1, 1 > normal_matrix(const Eigen::Matrix< ArrayDbl, 1, 2 > &A)
Definition: eigen_tools.hh:53
Eigen::Matrix< ArrayDbl, spacedim, dim > jacobian(const Eigen::Matrix< ArrayDbl, spacedim, dim+1 > &coords)
Definition: eigen_tools.hh:234
Eigen::Matrix< ArrayDbl, 3, 1 > inverse< 1, 3 >(const Eigen::Matrix< ArrayDbl, 1, 3 > &A)
Definition: eigen_tools.hh:212
Eigen::Matrix< ArrayDbl, 1, 1 > inverse< 1, 1 >(const Eigen::Matrix< ArrayDbl, 1, 1 > &A)
Definition: eigen_tools.hh:164
void resize_table(typename Eigen::Vector< ET, Eigen::Dynamic > &table, uint size)
Resize vector of Eigen::Array to given size.
Definition: eigen_tools.hh:225
Eigen::Matrix< ArrayDbl, 3, 2 > inverse< 2, 3 >(const Eigen::Matrix< ArrayDbl, 2, 3 > &A)
Definition: eigen_tools.hh:217
ArrayDbl determinant(const Eigen::Matrix< ArrayDbl, m, n > &A)
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:171
#define FMT_UNUSED
Definition: posix.h:75