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-2023 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  const tensor& T,
142  const scalar lambda,
143  const vector& direction1,
144  const vector& direction2
145 )
146 {
147  // Construct the linear system for this eigenvalue
148  tensor A(T - lambda*I);
149 
150  // Determinants of the 2x2 sub-matrices used to find the eigenvectors
151  scalar sd0, sd1, sd2;
152  scalar magSd0, magSd1, magSd2;
153 
154  // Sub-determinants for a unique eigenvalue
155  sd0 = A.yy()*A.zz() - A.yz()*A.zy();
156  sd1 = A.zz()*A.xx() - A.zx()*A.xz();
157  sd2 = A.xx()*A.yy() - A.xy()*A.yx();
158  magSd0 = mag(sd0);
159  magSd1 = mag(sd1);
160  magSd2 = mag(sd2);
161 
162  // Evaluate the eigenvector using the largest sub-determinant
163  if (magSd0 >= magSd1 && magSd0 >= magSd2 && magSd0 > small)
164  {
165  vector ev
166  (
167  1,
168  (A.yz()*A.zx() - A.zz()*A.yx())/sd0,
169  (A.zy()*A.yx() - A.yy()*A.zx())/sd0
170  );
171 
172  return ev/mag(ev);
173  }
174  else if (magSd1 >= magSd2 && magSd1 > small)
175  {
176  vector ev
177  (
178  (A.xz()*A.zy() - A.zz()*A.xy())/sd1,
179  1,
180  (A.zx()*A.xy() - A.xx()*A.zy())/sd1
181  );
182 
183  return ev/mag(ev);
184  }
185  else if (magSd2 > small)
186  {
187  vector ev
188  (
189  (A.xy()*A.yz() - A.yy()*A.xz())/sd2,
190  (A.yx()*A.xz() - A.xx()*A.yz())/sd2,
191  1
192  );
193 
194  return ev/mag(ev);
195  }
196 
197  // Sub-determinants for a repeated eigenvalue
198  sd0 = A.yy()*direction1.z() - A.yz()*direction1.y();
199  sd1 = A.zz()*direction1.x() - A.zx()*direction1.z();
200  sd2 = A.xx()*direction1.y() - A.xy()*direction1.x();
201  magSd0 = mag(sd0);
202  magSd1 = mag(sd1);
203  magSd2 = mag(sd2);
204 
205  // Evaluate the eigenvector using the largest sub-determinant
206  if (magSd0 >= magSd1 && magSd0 >= magSd2 && magSd0 > small)
207  {
208  vector ev
209  (
210  1,
211  (A.yz()*direction1.x() - direction1.z()*A.yx())/sd0,
212  (direction1.y()*A.yx() - A.yy()*direction1.x())/sd0
213  );
214 
215  return ev/mag(ev);
216  }
217  else if (magSd1 >= magSd2 && magSd1 > small)
218  {
219  vector ev
220  (
221  (direction1.z()*A.zy() - A.zz()*direction1.y())/sd1,
222  1,
223  (A.zx()*direction1.y() - direction1.x()*A.zy())/sd1
224  );
225 
226  return ev/mag(ev);
227  }
228  else if (magSd2 > small)
229  {
230  vector ev
231  (
232  (A.xy()*direction1.z() - direction1.y()*A.xz())/sd2,
233  (direction1.x()*A.xz() - A.xx()*direction1.z())/sd2,
234  1
235  );
236 
237  return ev/mag(ev);
238  }
239 
240  // Triple eigenvalue
241  return direction1^direction2;
242 }
243 
244 
246 {
247  vector Ux(1, 0, 0), Uy(0, 1, 0), Uz(0, 0, 1);
248 
249  Ux = eigenVector(T, lambdas.x(), Uy, Uz);
250  Uy = eigenVector(T, lambdas.y(), Uz, Ux);
251  Uz = eigenVector(T, lambdas.z(), Ux, Uy);
252 
253  return tensor(Ux, Uy, Uz);
254 }
255 
256 
258 {
259  const vector lambdas(eigenValues(T));
260 
261  return eigenVectors(T, lambdas);
262 }
263 
264 
266 {
267  return eigenValues(tensor(T));
268 }
269 
270 
272 (
273  const symmTensor& T,
274  const scalar lambda,
275  const vector& direction1,
276  const vector& direction2
277 )
278 {
279  return eigenVector(tensor(T), lambda, direction1, direction2);
280 }
281 
282 
284 {
285  return eigenVectors(tensor(T), lambdas);
286 }
287 
288 
290 {
291  return eigenVectors(tensor(T));
292 }
293 
294 
295 // ************************************************************************* //
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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: 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
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:25
dimensionedScalar lambda(viscosity->lookup("lambda"))
#define WarningInFunction
Report a warning using Foam::Warning.
mathematical constants.
const dimensionedScalar c
Speed of light in a vacuum.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedTensor eigenVectors(const dimensionedTensor &dt)
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
static const Identity< scalar > I
Definition: Identity.H:93
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedVector eigenValues(const dimensionedTensor &dt)
error FatalError
void Swap(T &a, T &b)
Definition: Swap.H:43
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
vector eigenVector(const tensor &T, const scalar lambda, const vector &direction1, const vector &direction2)
Definition: tensor.C:140