tensor.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "tensor.H"
27 #include "cubicEqn.H"
28 #include "mathematicalConstants.H"
29 
30 using namespace Foam::constant::mathematical;
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 template<>
35 const char* const Foam::tensor::vsType::typeName = "tensor";
36 
37 template<>
38 const char* const Foam::tensor::vsType::componentNames[] =
39 {
40  "xx", "xy", "xz",
41  "yx", "yy", "yz",
42  "zx", "zy", "zz"
43 };
44 
45 template<>
47 
48 template<>
50 
51 template<>
53 
54 template<>
56 
57 template<>
59 
60 template<>
62 
63 template<>
65 
66 template<>
68 (
69  1, 0, 0,
70  0, 1, 0,
71  0, 0, 1
72 );
73 
74 
75 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76 
78 {
79  // Coefficients of the characteristic cubic polynomial (a = 1)
80  const scalar b =
81  - t.xx() - t.yy() - t.zz();
82  const scalar c =
83  t.xx()*t.yy() + t.xx()*t.zz() + t.yy()*t.zz()
84  - t.xy()*t.yx() - t.yz()*t.zy() - t.zx()*t.xz();
85  const scalar d =
86  - t.xx()*t.yy()*t.zz()
87  - t.xy()*t.yz()*t.zx() - t.xz()*t.zy()*t.yx()
88  + t.xx()*t.yz()*t.zy() + t.yy()*t.zx()*t.xz() + t.zz()*t.xy()*t.yx();
89 
90  // Solve
91  Roots<3> roots = cubicEqn(1, b, c, d).roots();
92 
93  // Check the root types
95  forAll(roots, i)
96  {
97  switch (roots.type(i))
98  {
99  case rootType::real:
100  lambda[i] = roots[i];
101  break;
102  case rootType::complex:
104  << "Complex eigenvalues detected for tensor: " << t
105  << endl;
106  lambda[i] = 0;
107  break;
108  case rootType::posInf:
109  lambda[i] = vGreat;
110  break;
111  case rootType::negInf:
112  lambda[i] = - vGreat;
113  break;
114  case rootType::nan:
116  << "Eigenvalue calculation failed for tensor: " << t
117  << exit(FatalError);
118  }
119  }
120 
121  // Sort the eigenvalues into ascending order
122  if (lambda.x() > lambda.y())
123  {
124  Swap(lambda.x(), lambda.y());
125  }
126  if (lambda.y() > lambda.z())
127  {
128  Swap(lambda.y(), lambda.z());
129  }
130  if (lambda.x() > lambda.y())
131  {
132  Swap(lambda.x(), lambda.y());
133  }
134 
135  return lambda;
136 }
137 
138 
140 {
141  // Coefficients of the characteristic cubic polynomial (a = 1)
142  const scalar b =
143  - t.xx() - t.yy() - t.zz();
144  const scalar c =
145  t.xx()*t.yy() + t.xx()*t.zz() + t.yy()*t.zz()
146  - t.xy()*t.yx() - t.yz()*t.zy() - t.zx()*t.xz();
147  const scalar d =
148  - t.xx()*t.yy()*t.zz()
149  - t.xy()*t.yz()*t.zx() - t.xz()*t.zy()*t.yx()
150  + t.xx()*t.yz()*t.zy() + t.yy()*t.zx()*t.xz() + t.zz()*t.xy()*t.yx();
151 
152  // Solve
153  Roots<3> roots = cubicEqn(1, b, c, d).roots(true);
154 
155  // Check the root types
157  forAll(roots, i)
158  {
159  switch (roots.type(i))
160  {
161  case rootType::real:
162  lambda[i] = roots[i];
163  break;
164  case rootType::complex:
166  << "Complex eigenvalues detected for symmTensor: " << t
167  << exit(FatalError);
168  break;
169  case rootType::posInf:
170  lambda[i] = vGreat;
171  break;
172  case rootType::negInf:
173  lambda[i] = - vGreat;
174  break;
175  case rootType::nan:
177  << "Eigenvalue calculation failed for symmTensor: " << t
178  << exit(FatalError);
179  }
180  }
181 
182  // Sort the eigenvalues into ascending order
183  if (lambda.x() > lambda.y())
184  {
185  Swap(lambda.x(), lambda.y());
186  }
187  if (lambda.y() > lambda.z())
188  {
189  Swap(lambda.y(), lambda.z());
190  }
191  if (lambda.x() > lambda.y())
192  {
193  Swap(lambda.x(), lambda.y());
194  }
195 
196  return lambda;
197 }
198 
199 
200 template<class TensorType>
202 (
203  const TensorType& T,
204  const scalar lambda,
205  const vector& direction1,
206  const vector& direction2
207 )
208 {
209  // Construct the linear system for this eigenvalue
210  TensorType A(T - lambda*I);
211 
212  // Determinants of the 2x2 sub-matrices used to find the eigenvectors
213  scalar sd0, sd1, sd2;
214  scalar magSd0, magSd1, magSd2;
215 
216  // Sub-determinants for a unique eigenvalue
217  sd0 = A.yy()*A.zz() - A.yz()*A.zy();
218  sd1 = A.zz()*A.xx() - A.zx()*A.xz();
219  sd2 = A.xx()*A.yy() - A.xy()*A.yx();
220  magSd0 = mag(sd0);
221  magSd1 = mag(sd1);
222  magSd2 = mag(sd2);
223 
224  // Evaluate the eigenvector using the largest sub-determinant
225  if (magSd0 >= magSd1 && magSd0 >= magSd2 && magSd0 > small)
226  {
227  vector ev
228  (
229  1,
230  (A.yz()*A.zx() - A.zz()*A.yx())/sd0,
231  (A.zy()*A.yx() - A.yy()*A.zx())/sd0
232  );
233 
234  return ev/mag(ev);
235  }
236  else if (magSd1 >= magSd2 && magSd1 > small)
237  {
238  vector ev
239  (
240  (A.xz()*A.zy() - A.zz()*A.xy())/sd1,
241  1,
242  (A.zx()*A.xy() - A.xx()*A.zy())/sd1
243  );
244 
245  return ev/mag(ev);
246  }
247  else if (magSd2 > small)
248  {
249  vector ev
250  (
251  (A.xy()*A.yz() - A.yy()*A.xz())/sd2,
252  (A.yx()*A.xz() - A.xx()*A.yz())/sd2,
253  1
254  );
255 
256  return ev/mag(ev);
257  }
258 
259  // Sub-determinants for a repeated eigenvalue
260  sd0 = A.yy()*direction1.z() - A.yz()*direction1.y();
261  sd1 = A.zz()*direction1.x() - A.zx()*direction1.z();
262  sd2 = A.xx()*direction1.y() - A.xy()*direction1.x();
263  magSd0 = mag(sd0);
264  magSd1 = mag(sd1);
265  magSd2 = mag(sd2);
266 
267  // Evaluate the eigenvector using the largest sub-determinant
268  if (magSd0 >= magSd1 && magSd0 >= magSd2 && magSd0 > small)
269  {
270  vector ev
271  (
272  1,
273  (A.yz()*direction1.x() - direction1.z()*A.yx())/sd0,
274  (direction1.y()*A.yx() - A.yy()*direction1.x())/sd0
275  );
276 
277  return ev/mag(ev);
278  }
279  else if (magSd1 >= magSd2 && magSd1 > small)
280  {
281  vector ev
282  (
283  (direction1.z()*A.zy() - A.zz()*direction1.y())/sd1,
284  1,
285  (A.zx()*direction1.y() - direction1.x()*A.zy())/sd1
286  );
287 
288  return ev/mag(ev);
289  }
290  else if (magSd2 > small)
291  {
292  vector ev
293  (
294  (A.xy()*direction1.z() - direction1.y()*A.xz())/sd2,
295  (direction1.x()*A.xz() - A.xx()*direction1.z())/sd2,
296  1
297  );
298 
299  return ev/mag(ev);
300  }
301 
302  // Triple eigenvalue
303  return direction1^direction2;
304 }
305 
306 
308 {
309  vector Ux(1, 0, 0), Uy(0, 1, 0), Uz(0, 0, 1);
310 
311  Ux = eigenVector(T, lambdas.x(), Uy, Uz);
312  Uy = eigenVector(T, lambdas.y(), Uz, Ux);
313  Uz = eigenVector(T, lambdas.z(), Ux, Uy);
314 
315  return tensor(Ux, Uy, Uz);
316 }
317 
318 
320 {
321  vector Ux(1, 0, 0), Uy(0, 1, 0), Uz(0, 0, 1);
322 
323  Ux = eigenVector(T, lambdas.x(), Uy, Uz);
324  Uy = eigenVector(T, lambdas.y(), Uz, Ux);
325  Uz = eigenVector(T, lambdas.z(), Ux, Uy);
326 
327  return tensor(Ux, Uy, Uz);
328 }
329 
330 
332 {
333  const vector lambdas(eigenValues(T));
334  return eigenVectors(T, lambdas);
335 }
336 
337 
339 {
340  const vector lambdas(eigenValues(T));
341  return eigenVectors(T, lambdas);
342 }
343 
344 
345 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Templated storage for the roots of polynomial equations, plus flags to indicate the nature of the roo...
Definition: Roots.H:67
void type(const direction i, const rootType t)
Set the type of the i-th root.
Definition: RootsI.H:120
const Cmpt & xx() const
Definition: SymmTensorI.H:87
const Cmpt & yx() const
Definition: SymmTensorI.H:105
const Cmpt & yz() const
Definition: SymmTensorI.H:117
const Cmpt & xz() const
Definition: SymmTensorI.H:99
const Cmpt & zz() const
Definition: SymmTensorI.H:135
const Cmpt & xy() const
Definition: SymmTensorI.H:93
const Cmpt & zx() const
Definition: SymmTensorI.H:123
const Cmpt & zy() const
Definition: SymmTensorI.H:129
const Cmpt & yy() const
Definition: SymmTensorI.H:111
const Cmpt & xx() const
Definition: TensorI.H:163
const Cmpt & yx() const
Definition: TensorI.H:184
static const Tensor I
Definition: Tensor.H:83
const Cmpt & yz() const
Definition: TensorI.H:198
const Cmpt & xz() const
Definition: TensorI.H:177
const Cmpt & zz() const
Definition: TensorI.H:219
const Cmpt & xy() const
Definition: TensorI.H:170
const Cmpt & zx() const
Definition: TensorI.H:205
const Cmpt & zy() const
Definition: TensorI.H:212
const Cmpt & yy() const
Definition: TensorI.H:191
static const Form zero
Definition: VectorSpace.H:118
static const char *const componentNames[]
Definition: VectorSpace.H:117
static const Form one
Definition: VectorSpace.H:119
static const Form nan
Definition: VectorSpace.H:124
static const Form rootMax
Definition: VectorSpace.H:122
static const Form max
Definition: VectorSpace.H:120
static const Form rootMin
Definition: VectorSpace.H:123
static Form uniform(const Cmpt &s)
Return a VectorSpace with all elements = s.
Definition: VectorSpaceI.H:168
static const Form min
Definition: VectorSpace.H:121
static const char *const typeName
Definition: VectorSpace.H:116
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
Cubic equation of the form a*x^3 + b*x^2 + c*x + d = 0.
Definition: cubicEqn.H:52
Roots< 3 > roots(const bool real=false) const
Get the roots.
Definition: cubicEqn.C:32
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
volScalarField & b
Definition: createFields.H:27
dimensionedScalar lambda(viscosity->lookup("lambda"))
#define WarningInFunction
Report a warning using Foam::Warning.
mathematical constants.
const dimensionedScalar c
Speed of light in a vacuum.
static const coefficient A("A", dimPressure, 611.21)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
vector eigenVector(const TensorType &T, const scalar lambda, const vector &direction1, const vector &direction2)
void eigenValues(LagrangianPatchField< vector > &f, const LagrangianPatchField< tensor > &f1)
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
static const Identity< scalar > I
Definition: Identity.H:93
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void eigenVectors(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
error FatalError
void Swap(T &a, T &b)
Definition: Swap.H:43