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