tensor.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2020-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "tensor.H"
30 #include "cubicEqn.H"
32 #include "SVD.H"
33 
34 using namespace Foam::constant::mathematical;
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
39 {
40  // Return diagonal if T is effectively diagonal tensor
41  if
42  (
43  (
44  sqr(T.xy()) + sqr(T.xz()) + sqr(T.yz())
45  + sqr(T.yx()) + sqr(T.zx()) + sqr(T.zy())
46  ) < ROOTSMALL
47  )
48  {
49  return Vector<complex>
50  (
51  complex(T.xx()), complex(T.yy()), complex(T.zz())
52  );
53  }
54 
55  // Coefficients of the characteristic cubic polynomial (a = 1)
56  const scalar b = - T.xx() - T.yy() - T.zz();
57  const scalar c =
58  T.xx()*T.yy() + T.xx()*T.zz() + T.yy()*T.zz()
59  - T.xy()*T.yx() - T.yz()*T.zy() - T.zx()*T.xz();
60  const scalar d =
61  - T.xx()*T.yy()*T.zz()
62  - T.xy()*T.yz()*T.zx() - T.xz()*T.zy()*T.yx()
63  + T.xx()*T.yz()*T.zy() + T.yy()*T.zx()*T.xz() + T.zz()*T.xy()*T.yx();
64 
65  // Determine the roots of the characteristic cubic polynomial
66  const Roots<3> roots(cubicEqn(1, b, c, d).roots());
67  // Check the root types
68  bool isComplex = false;
69  forAll(roots, i)
70  {
71  switch (roots.type(i))
72  {
73  case roots::complex:
74  isComplex = true;
75  break;
76  case roots::posInf:
77  case roots::negInf:
78  case roots::nan:
80  << "Eigenvalue computation fails for tensor: " << T
81  << "due to the not-a-number root = " << roots[i]
82  << endl;
83  case roots::real:
84  break;
85  }
86  }
87 
88  if (isComplex)
89  {
90  return
92  (
93  complex(roots[0], 0),
94  complex(roots[1], roots[2]),
95  complex(roots[1], -roots[2])
96  );
97  }
98 
99  return
101  (
102  complex(roots[0], 0),
103  complex(roots[1], 0),
104  complex(roots[2], 0)
105  );
106 }
107 
108 
110 (
111  const tensor& T,
112  const complex& eVal,
113  const Vector<complex>& standardBasis1,
114  const Vector<complex>& standardBasis2
115 )
116 {
117  // Construct the characteristic equation system for this eigenvalue
118  Tensor<complex> A(Zero);
119  forAll(A, i)
120  {
121  A[i] = complex(T[i], 0);
122  }
123  A.xx() -= eVal;
124  A.yy() -= eVal;
125  A.zz() -= eVal;
126 
127  // Determinants of the 2x2 sub-matrices used to find the eigenvectors
128  // Sub-determinants for a unique eigenvenvalue
129  complex sd0 = A.yy()*A.zz() - A.yz()*A.zy();
130  complex sd1 = A.zz()*A.xx() - A.zx()*A.xz();
131  complex sd2 = A.xx()*A.yy() - A.xy()*A.yx();
132  scalar magSd0 = mag(sd0);
133  scalar magSd1 = mag(sd1);
134  scalar magSd2 = mag(sd2);
135 
136  // Evaluate the eigenvector using the largest sub-determinant
137  if (magSd0 >= magSd1 && magSd0 >= magSd2 && magSd0 > SMALL)
138  {
139  const Vector<complex> eVec
140  (
141  complex(1, 0),
142  (A.yz()*A.zx() - A.zz()*A.yx())/sd0,
143  (A.zy()*A.yx() - A.yy()*A.zx())/sd0
144  );
145 
146  #ifdef FULLDEBUG
147  if (mag(eVec) < SMALL)
148  {
150  << "Eigenvector magnitude should be non-zero:"
151  << "mag(eigenvector) = " << mag(eVec)
152  << abort(FatalError);
153  }
154  #endif
155 
156  return eVec/mag(eVec);
157  }
158  else if (magSd1 >= magSd2 && magSd1 > SMALL)
159  {
160  const Vector<complex> eVec
161  (
162  (A.xz()*A.zy() - A.zz()*A.xy())/sd1,
163  complex(1, 0),
164  (A.zx()*A.xy() - A.xx()*A.zy())/sd1
165  );
166 
167  #ifdef FULLDEBUG
168  if (mag(eVec) < SMALL)
169  {
171  << "Eigenvector magnitude should be non-zero:"
172  << "mag(eigenvector) = " << mag(eVec)
173  << abort(FatalError);
174  }
175  #endif
176 
177  return eVec/mag(eVec);
178  }
179  else if (magSd2 > SMALL)
180  {
181  const Vector<complex> eVec
182  (
183  (A.xy()*A.yz() - A.yy()*A.xz())/sd2,
184  (A.yx()*A.xz() - A.xx()*A.yz())/sd2,
185  complex(1, 0)
186  );
187 
188  #ifdef FULLDEBUG
189  if (mag(eVec) < SMALL)
190  {
192  << "Eigenvector magnitude should be non-zero:"
193  << "mag(eigenvector) = " << mag(eVec)
194  << abort(FatalError);
195  }
196  #endif
197 
198  return eVec/mag(eVec);
199  }
200 
201  // Sub-determinants for a repeated eigenvalue
202  sd0 = A.yy()*standardBasis1.z() - A.yz()*standardBasis1.y();
203  sd1 = A.zz()*standardBasis1.x() - A.zx()*standardBasis1.z();
204  sd2 = A.xx()*standardBasis1.y() - A.xy()*standardBasis1.x();
205  magSd0 = mag(sd0);
206  magSd1 = mag(sd1);
207  magSd2 = mag(sd2);
208 
209  // Evaluate the eigenvector using the largest sub-determinant
210  if (magSd0 >= magSd1 && magSd0 >= magSd2 && magSd0 > SMALL)
211  {
212  const Vector<complex> eVec
213  (
214  complex(1, 0),
215  (A.yz()*standardBasis1.x() - standardBasis1.z()*A.yx())/sd0,
216  (standardBasis1.y()*A.yx() - A.yy()*standardBasis1.x())/sd0
217  );
218 
219  #ifdef FULLDEBUG
220  if (mag(eVec) < SMALL)
221  {
223  << "Eigenvector magnitude should be non-zero:"
224  << "mag(eigenvector) = " << mag(eVec)
225  << abort(FatalError);
226  }
227  #endif
228 
229  return eVec/mag(eVec);
230  }
231  else if (magSd1 >= magSd2 && magSd1 > SMALL)
232  {
233  const Vector<complex> eVec
234  (
235  (standardBasis1.z()*A.zy() - A.zz()*standardBasis1.y())/sd1,
236  complex(1, 0),
237  (A.zx()*standardBasis1.y() - standardBasis1.x()*A.zy())/sd1
238  );
239 
240  #ifdef FULLDEBUG
241  if (mag(eVec) < SMALL)
242  {
244  << "Eigenvector magnitude should be non-zero:"
245  << "mag(eigenvector) = " << mag(eVec)
246  << abort(FatalError);
247  }
248  #endif
249 
250  return eVec/mag(eVec);
251  }
252  else if (magSd2 > SMALL)
253  {
254  const Vector<complex> eVec
255  (
256  (A.xy()*standardBasis1.z() - standardBasis1.y()*A.xz())/sd2,
257  (standardBasis1.x()*A.xz() - A.xx()*standardBasis1.z())/sd2,
258  complex(1, 0)
259  );
260 
261  #ifdef FULLDEBUG
262  if (mag(eVec) < SMALL)
263  {
265  << "Eigenvector magnitude should be non-zero:"
266  << "mag(eigenvector) = " << mag(eVec)
267  << abort(FatalError);
268  }
269  #endif
270 
271  return eVec/mag(eVec);
272  }
274  // Triple eigenvalue
275  return standardBasis1^standardBasis2;
276 }
277 
278 
280 (
281  const tensor& T,
282  const Vector<complex>& eVals
283 )
284 {
285  Vector<complex> Ux(complex(1, 0), Zero, Zero);
286  Vector<complex> Uy(Zero, complex(1, 0), Zero);
287  Vector<complex> Uz(Zero, Zero, complex(1, 0));
288 
289  Ux = eigenVector(T, eVals.x(), Uy, Uz);
290  Uy = eigenVector(T, eVals.y(), Uz, Ux);
291  Uz = eigenVector(T, eVals.z(), Ux, Uy);
292 
293  return Tensor<complex>(Ux, Uy, Uz);
294 }
295 
296 
298 {
299  const Vector<complex> eVals(eigenValues(T));
300 
301  return eigenVectors(T, eVals);
302 }
303 
304 
306 {
307  const scalar detval = t.det();
308 
309  if (detval < ROOTVSMALL)
310  {
311  // Fall back to pseudo inverse
312  scalarRectangularMatrix mat(3, 3);
313 
314  mat(0,0) = t.xx(); mat(0,1) = t.xy(); mat(0,2) = t.xz();
315  mat(1,0) = t.yx(); mat(1,1) = t.yy(); mat(1,2) = t.yz();
316  mat(2,0) = t.zx(); mat(2,1) = t.zy(); mat(2,2) = t.zz();
317 
318  mat = SVDinv(mat);
319 
320  return tensor
321  (
322  mat(0,0), mat(0,1), mat(0,2),
323  mat(1,0), mat(1,1), mat(1,2),
324  mat(2,0), mat(2,1), mat(2,2)
325  );
326  }
327 
328  return Foam::inv(t, detval);
329 }
330 
331 
332 // ************************************************************************* //
vector eigenVector(const symmTensor &T, const scalar eVal, const vector &standardBasis1, const vector &standardBasis2)
Return a real eigenvector corresponding to a given real eigenvalue of a given symmTensor.
Definition: symmTensor.C:146
scalarRectangularMatrix SVDinv(const scalarRectangularMatrix &A, scalar minCondition=0)
Return the inverse of matrix A using SVD.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Tensor< scalar > tensor
Definition: symmTensor.H:57
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Templated storage for the roots of polynomial equations, plus flags to indicate the nature of the roo...
Definition: Roots.H:68
Cmpt det() const
The determinate.
Definition: SymmTensorI.H:235
dimensionedTensor eigenVectors(const dimensionedSymmTensor &dt)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Container to encapsulate various operations for cubic equation of the forms with real coefficients: ...
Definition: cubicEqn.H:108
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Mathematical constants.
symmTensor pinv(const symmTensor &st)
Return inverse of a given symmTensor, and fall back to pseudo-inverse if the symmTensor is singular...
Definition: symmTensor.C:334
errorManip< error > abort(error &err)
Definition: errorManip.H:139
Templated 3D Vector derived from VectorSpace adding construction from 3 components, element access using x(), y() and z() member functions and the inner-product (dot-product) and cross-product operators.
Definition: Vector.H:58
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
#define WarningInFunction
Report a warning using Foam::Warning.
void type(const direction i, const roots::type t)
Set the type of the i-th root.
Definition: RootsI.H:115
const dimensionedScalar c
Speed of light in a vacuum.
dimensionedVector eigenValues(const dimensionedSymmTensor &dt)
A templated (3 x 3) tensor of objects of <T> derived from MatrixSpace.
Definition: complexI.H:266
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Tensor of scalars, i.e. Tensor<scalar>.
RectangularMatrix< scalar > scalarRectangularMatrix
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127