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 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
33using namespace Foam::constant::mathematical;
34
35// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36
38{
39 // Return diagonal if T is effectively diagonal tensor
40 if
41 (
42 (
43 sqr(T.xy()) + sqr(T.xz()) + sqr(T.yz())
44 + sqr(T.yx()) + sqr(T.zx()) + sqr(T.zy())
45 ) < ROOTSMALL
46 )
47 {
48 return Vector<complex>
49 (
50 complex(T.xx()), complex(T.yy()), complex(T.zz())
51 );
52 }
53
54 // Coefficients of the characteristic cubic polynomial (a = 1)
55 const scalar b = - T.xx() - T.yy() - T.zz();
56 const scalar c =
57 T.xx()*T.yy() + T.xx()*T.zz() + T.yy()*T.zz()
58 - T.xy()*T.yx() - T.yz()*T.zy() - T.zx()*T.xz();
59 const scalar d =
60 - T.xx()*T.yy()*T.zz()
61 - T.xy()*T.yz()*T.zx() - T.xz()*T.zy()*T.yx()
62 + T.xx()*T.yz()*T.zy() + T.yy()*T.zx()*T.xz() + T.zz()*T.xy()*T.yx();
63
64 // Determine the roots of the characteristic cubic polynomial
65 const Roots<3> roots(cubicEqn(1, b, c, d).roots());
66 // Check the root types
67 bool isComplex = false;
68 forAll(roots, i)
69 {
70 switch (roots.type(i))
71 {
72 case roots::complex:
73 isComplex = true;
74 break;
75 case roots::posInf:
76 case roots::negInf:
77 case roots::nan:
79 << "Eigenvalue computation fails for tensor: " << T
80 << "due to the not-a-number root = " << roots[i]
81 << endl;
82 case roots::real:
83 break;
84 }
85 }
86
87 if (isComplex)
88 {
89 return
91 (
92 complex(roots[0], 0),
93 complex(roots[1], roots[2]),
94 complex(roots[1], -roots[2])
95 );
96 }
97
98 return
100 (
101 complex(roots[0], 0),
102 complex(roots[1], 0),
103 complex(roots[2], 0)
104 );
105}
106
107
109(
110 const tensor& T,
111 const complex& eVal,
112 const Vector<complex>& standardBasis1,
113 const Vector<complex>& standardBasis2
114)
115{
116 // Construct the characteristic equation system for this eigenvalue
117 Tensor<complex> A(Zero);
118 forAll(A, i)
119 {
120 A[i] = complex(T[i], 0);
121 }
122 A.xx() -= eVal;
123 A.yy() -= eVal;
124 A.zz() -= eVal;
125
126 // Determinants of the 2x2 sub-matrices used to find the eigenvectors
127 // Sub-determinants for a unique eigenvenvalue
128 complex sd0 = A.yy()*A.zz() - A.yz()*A.zy();
129 complex sd1 = A.zz()*A.xx() - A.zx()*A.xz();
130 complex sd2 = A.xx()*A.yy() - A.xy()*A.yx();
131 scalar magSd0 = mag(sd0);
132 scalar magSd1 = mag(sd1);
133 scalar magSd2 = mag(sd2);
134
135 // Evaluate the eigenvector using the largest sub-determinant
136 if (magSd0 >= magSd1 && magSd0 >= magSd2 && magSd0 > SMALL)
137 {
138 const Vector<complex> eVec
139 (
140 complex(1, 0),
141 (A.yz()*A.zx() - A.zz()*A.yx())/sd0,
142 (A.zy()*A.yx() - A.yy()*A.zx())/sd0
143 );
144
145 #ifdef FULLDEBUG
146 if (mag(eVec) < SMALL)
147 {
149 << "Eigenvector magnitude should be non-zero:"
150 << "mag(eigenvector) = " << mag(eVec)
151 << abort(FatalError);
152 }
153 #endif
154
155 return eVec/mag(eVec);
156 }
157 else if (magSd1 >= magSd2 && magSd1 > SMALL)
158 {
159 const Vector<complex> eVec
160 (
161 (A.xz()*A.zy() - A.zz()*A.xy())/sd1,
162 complex(1, 0),
163 (A.zx()*A.xy() - A.xx()*A.zy())/sd1
164 );
165
166 #ifdef FULLDEBUG
167 if (mag(eVec) < SMALL)
168 {
170 << "Eigenvector magnitude should be non-zero:"
171 << "mag(eigenvector) = " << mag(eVec)
172 << abort(FatalError);
173 }
174 #endif
175
176 return eVec/mag(eVec);
177 }
178 else if (magSd2 > SMALL)
179 {
180 const Vector<complex> eVec
181 (
182 (A.xy()*A.yz() - A.yy()*A.xz())/sd2,
183 (A.yx()*A.xz() - A.xx()*A.yz())/sd2,
184 complex(1, 0)
185 );
186
187 #ifdef FULLDEBUG
188 if (mag(eVec) < SMALL)
189 {
191 << "Eigenvector magnitude should be non-zero:"
192 << "mag(eigenvector) = " << mag(eVec)
193 << abort(FatalError);
194 }
195 #endif
196
197 return eVec/mag(eVec);
198 }
199
200 // Sub-determinants for a repeated eigenvalue
201 sd0 = A.yy()*standardBasis1.z() - A.yz()*standardBasis1.y();
202 sd1 = A.zz()*standardBasis1.x() - A.zx()*standardBasis1.z();
203 sd2 = A.xx()*standardBasis1.y() - A.xy()*standardBasis1.x();
204 magSd0 = mag(sd0);
205 magSd1 = mag(sd1);
206 magSd2 = mag(sd2);
207
208 // Evaluate the eigenvector using the largest sub-determinant
209 if (magSd0 >= magSd1 && magSd0 >= magSd2 && magSd0 > SMALL)
210 {
211 const Vector<complex> eVec
212 (
213 complex(1, 0),
214 (A.yz()*standardBasis1.x() - standardBasis1.z()*A.yx())/sd0,
215 (standardBasis1.y()*A.yx() - A.yy()*standardBasis1.x())/sd0
216 );
217
218 #ifdef FULLDEBUG
219 if (mag(eVec) < SMALL)
220 {
222 << "Eigenvector magnitude should be non-zero:"
223 << "mag(eigenvector) = " << mag(eVec)
224 << abort(FatalError);
225 }
226 #endif
227
228 return eVec/mag(eVec);
229 }
230 else if (magSd1 >= magSd2 && magSd1 > SMALL)
231 {
232 const Vector<complex> eVec
233 (
234 (standardBasis1.z()*A.zy() - A.zz()*standardBasis1.y())/sd1,
235 complex(1, 0),
236 (A.zx()*standardBasis1.y() - standardBasis1.x()*A.zy())/sd1
237 );
238
239 #ifdef FULLDEBUG
240 if (mag(eVec) < SMALL)
241 {
243 << "Eigenvector magnitude should be non-zero:"
244 << "mag(eigenvector) = " << mag(eVec)
245 << abort(FatalError);
246 }
247 #endif
248
249 return eVec/mag(eVec);
250 }
251 else if (magSd2 > SMALL)
252 {
253 const Vector<complex> eVec
254 (
255 (A.xy()*standardBasis1.z() - standardBasis1.y()*A.xz())/sd2,
256 (standardBasis1.x()*A.xz() - A.xx()*standardBasis1.z())/sd2,
257 complex(1, 0)
258 );
259
260 #ifdef FULLDEBUG
261 if (mag(eVec) < SMALL)
262 {
264 << "Eigenvector magnitude should be non-zero:"
265 << "mag(eigenvector) = " << mag(eVec)
266 << abort(FatalError);
267 }
268 #endif
269
270 return eVec/mag(eVec);
271 }
272
273 // Triple eigenvalue
274 return standardBasis1^standardBasis2;
275}
276
277
279(
280 const tensor& T,
281 const Vector<complex>& eVals
282)
283{
284 Vector<complex> Ux(complex(1, 0), Zero, Zero);
285 Vector<complex> Uy(Zero, complex(1, 0), Zero);
286 Vector<complex> Uz(Zero, Zero, complex(1, 0));
287
288 Ux = eigenVector(T, eVals.x(), Uy, Uz);
289 Uy = eigenVector(T, eVals.y(), Uz, Ux);
290 Uz = eigenVector(T, eVals.z(), Ux, Uy);
291
292 return Tensor<complex>(Ux, Uy, Uz);
293}
294
295
297{
298 const Vector<complex> eVals(eigenValues(T));
299
300 return eigenVectors(T, eVals);
301}
302
303
304// ************************************************************************* //
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Templated storage for the roots of polynomial equations, plus flags to indicate the nature of the roo...
Definition: Roots.H:73
void type(const direction i, const roots::type t)
Set the type of the i-th root.
Definition: RootsI.H:122
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:65
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
A complex number, similar to the C++ complex type.
Definition: complex.H:83
Container to encapsulate various operations for cubic equation of the forms with real coefficients:
Definition: cubicEqn.H:115
const volScalarField & T
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
#define WarningInFunction
Report a warning using Foam::Warning.
Mathematical constants.
dimensionedTensor eigenVectors(const dimensionedSymmTensor &dt)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedVector eigenValues(const dimensionedSymmTensor &dt)
vector eigenVector(const symmTensor &T, const scalar eVal, const vector &standardBasis1, const vector &standardBasis2)
Definition: symmTensor.C:152
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
volScalarField & b
Definition: createFields.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333