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 | Copyright (C) 2011-2016 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 "mathematicalConstants.H"
28 
29 using namespace Foam::constant::mathematical;
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 template<>
34 const char* const Foam::tensor::vsType::typeName = "tensor";
35 
36 template<>
37 const char* const Foam::tensor::vsType::componentNames[] =
38 {
39  "xx", "xy", "xz",
40  "yx", "yy", "yz",
41  "zx", "zy", "zz"
42 };
43 
44 template<>
46 
47 template<>
49 
50 template<>
52 
53 template<>
55 
56 template<>
58 
59 template<>
61 
62 template<>
64 (
65  1, 0, 0,
66  0, 1, 0,
67  0, 0, 1
68 );
69 
70 
71 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72 
74 {
75  // The eigenvalues
76  scalar i, ii, iii;
77 
78  // diagonal matrix
79  if
80  (
81  (
82  mag(t.xy()) + mag(t.xz()) + mag(t.yx())
83  + mag(t.yz()) + mag(t.zx()) + mag(t.zy())
84  )
85  < SMALL
86  )
87  {
88  i = t.xx();
89  ii = t.yy();
90  iii = t.zz();
91  }
92 
93  // non-diagonal matrix
94  else
95  {
96  // Coefficients of the characteristic polynmial
97  // x^3 + a*x^2 + b*x + c = 0
98  scalar a =
99  - t.xx() - t.yy() - t.zz();
100 
101  scalar b =
102  t.xx()*t.yy() + t.xx()*t.zz() + t.yy()*t.zz()
103  - t.xy()*t.yx() - t.yz()*t.zy() - t.zx()*t.xz();
104 
105  scalar c =
106  - t.xx()*t.yy()*t.zz()
107  - t.xy()*t.yz()*t.zx() - t.xz()*t.zy()*t.yx()
108  + t.xx()*t.yz()*t.zy() + t.yy()*t.zx()*t.xz() + t.zz()*t.xy()*t.yx();
109 
110  // Auxillary variables
111  scalar aBy3 = a/3;
112 
113  scalar P = (a*a - 3*b)/9; // == -p_wikipedia/3
114  scalar PPP = P*P*P;
115 
116  scalar Q = (2*a*a*a - 9*a*b + 27*c)/54; // == q_wikipedia/2
117  scalar QQ = Q*Q;
118 
119  // Three identical roots
120  if (mag(P) < SMALL && mag(Q) < SMALL)
121  {
122  return vector(- aBy3, - aBy3, - aBy3);
123  }
124 
125  // Two identical roots and one distinct root
126  else if (mag(PPP/QQ - 1) < SMALL)
127  {
128  scalar sqrtP = sqrt(P);
129  scalar signQ = sign(Q);
130 
131  i = ii = signQ*sqrtP - aBy3;
132  iii = - 2*signQ*sqrtP - aBy3;
133  }
134 
135  // Three distinct roots
136  else if (PPP > QQ)
137  {
138  scalar sqrtP = sqrt(P);
139  scalar value = cos(acos(Q/sqrt(PPP))/3);
140  scalar delta = sqrt(3 - 3*value*value);
141 
142  i = - 2*sqrtP*value - aBy3;
143  ii = sqrtP*(value + delta) - aBy3;
144  iii = sqrtP*(value - delta) - aBy3;
145  }
146 
147  // One real root, two imaginary roots
148  // based on the above logic, PPP must be less than QQ
149  else
150  {
152  << "complex eigenvalues detected for tensor: " << t
153  << endl;
154 
155  if (mag(P) < SMALL)
156  {
157  i = cbrt(QQ/2);
158  }
159  else
160  {
161  scalar w = cbrt(- Q - sqrt(QQ - PPP));
162  i = w + P/w - aBy3;
163  }
164 
165  return vector(-VGREAT, i, VGREAT);
166  }
167  }
168 
169  // Sort the eigenvalues into ascending order
170  if (i > ii)
171  {
172  Swap(i, ii);
173  }
174 
175  if (ii > iii)
176  {
177  Swap(ii, iii);
178  }
179 
180  if (i > ii)
181  {
182  Swap(i, ii);
183  }
184 
185  return vector(i, ii, iii);
186 }
187 
188 
190 (
191  const tensor& t,
192  const scalar lambda
193 )
194 {
195  // Constantly rotating direction ensures different eigenvectors are
196  // generated when called sequentially with a multiple eigenvalue
197  static vector direction(1,0,0);
198  vector oldDirection(direction);
199  scalar temp = direction[2];
200  direction[2] = direction[1];
201  direction[1] = direction[0];
202  direction[0] = temp;
203 
204  // Construct the linear system for this eigenvalue
205  tensor A(t - lambda*I);
206 
207  // Determinants of the 2x2 sub-matrices used to find the eigenvectors
208  scalar sd0, sd1, sd2;
209  scalar magSd0, magSd1, magSd2;
210 
211  // Sub-determinants for a unique eivenvalue
212  sd0 = A.yy()*A.zz() - A.yz()*A.zy();
213  sd1 = A.zz()*A.xx() - A.zx()*A.xz();
214  sd2 = A.xx()*A.yy() - A.xy()*A.yx();
215  magSd0 = mag(sd0);
216  magSd1 = mag(sd1);
217  magSd2 = mag(sd2);
218 
219  // Evaluate the eigenvector using the largest sub-determinant
220  if (magSd0 >= magSd1 && magSd0 >= magSd2 && magSd0 > SMALL)
221  {
222  vector ev
223  (
224  1,
225  (A.yz()*A.zx() - A.zz()*A.yx())/sd0,
226  (A.zy()*A.yx() - A.yy()*A.zx())/sd0
227  );
228 
229  return ev/mag(ev);
230  }
231  else if (magSd1 >= magSd2 && magSd1 > SMALL)
232  {
233  vector ev
234  (
235  (A.xz()*A.zy() - A.zz()*A.xy())/sd1,
236  1,
237  (A.zx()*A.xy() - A.xx()*A.zy())/sd1
238  );
239 
240  return ev/mag(ev);
241  }
242  else if (magSd2 > SMALL)
243  {
244  vector ev
245  (
246  (A.xy()*A.yz() - A.yy()*A.xz())/sd2,
247  (A.yx()*A.xz() - A.xx()*A.yz())/sd2,
248  1
249  );
250 
251  return ev/mag(ev);
252  }
253 
254  // Sub-determinants for a repeated eigenvalue
255  sd0 = A.yy()*direction.z() - A.yz()*direction.y();
256  sd1 = A.zz()*direction.x() - A.zx()*direction.z();
257  sd2 = A.xx()*direction.y() - A.xy()*direction.x();
258  magSd0 = mag(sd0);
259  magSd1 = mag(sd1);
260  magSd2 = mag(sd2);
261 
262  // Evaluate the eigenvector using the largest sub-determinant
263  if (magSd0 >= magSd1 && magSd0 >= magSd2 && magSd0 > SMALL)
264  {
265  vector ev
266  (
267  1,
268  (A.yz()*direction.x() - direction.z()*A.yx())/sd0,
269  (direction.y()*A.yx() - A.yy()*direction.x())/sd0
270  );
271 
272  return ev/mag(ev);
273  }
274  else if (magSd1 >= magSd2 && magSd1 > SMALL)
275  {
276  vector ev
277  (
278  (direction.z()*A.zy() - A.zz()*direction.y())/sd1,
279  1,
280  (A.zx()*direction.y() - direction.x()*A.zy())/sd1
281  );
282 
283  return ev/mag(ev);
284  }
285  else if (magSd2 > SMALL)
286  {
287  vector ev
288  (
289  (A.xy()*direction.z() - direction.y()*A.xz())/sd2,
290  (direction.x()*A.xz() - A.xx()*direction.z())/sd2,
291  1
292  );
293 
294  return ev/mag(ev);
295  }
296 
297  // Triple eigenvalue
298  return oldDirection;
299 }
300 
301 
303 {
304  vector evals(eigenValues(t));
305 
306  tensor evs
307  (
308  eigenVector(t, evals.x()),
309  eigenVector(t, evals.y()),
310  eigenVector(t, evals.z())
311  );
312 
313  return evs;
314 }
315 
316 
318 {
319  return eigenValues(tensor(t));
320 }
321 
322 
323 Foam::vector Foam::eigenVector(const symmTensor& t, const scalar lambda)
324 {
325  return eigenVector(tensor(t), lambda);
326 }
327 
328 
330 {
331  return eigenVectors(tensor(t));
332 }
333 
334 
335 // ************************************************************************* //
scalar delta
dimensionedScalar sign(const dimensionedScalar &ds)
static const Form max
Definition: VectorSpace.H:112
static const char *const typeName
Definition: VectorSpace.H:108
const Cmpt & z() const
Definition: VectorI.H:87
dimensionedScalar acos(const dimensionedScalar &ds)
const Cmpt & yy() const
Definition: TensorI.H:181
uint8_t direction
Definition: direction.H:46
const Cmpt & x() const
Definition: VectorI.H:75
static const char *const componentNames[]
Definition: VectorSpace.H:109
dimensionedVector eigenValues(const dimensionedTensor &dt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const Cmpt & zy() const
Definition: TensorI.H:202
static const Form rootMax
Definition: VectorSpace.H:114
const Cmpt & yz() const
Definition: TensorI.H:188
static const Form one
Definition: VectorSpace.H:111
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
static const Form min
Definition: VectorSpace.H:113
const Cmpt & xz() const
Definition: TensorI.H:167
static const Form rootMin
Definition: VectorSpace.H:115
dimensionedScalar cos(const dimensionedScalar &ds)
void Swap(T &a, T &b)
Definition: Swap.H:43
static const Identity< scalar > I
Definition: Identity.H:93
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
mathematical constants.
dimensionedScalar cbrt(const dimensionedScalar &ds)
const Cmpt & xx() const
Definition: TensorI.H:153
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & xy() const
Definition: TensorI.H:160
const Cmpt & yx() const
Definition: TensorI.H:174
static const Tensor I
Definition: Tensor.H:80
const Cmpt & zz() const
Definition: TensorI.H:209
dimensionedTensor eigenVectors(const dimensionedTensor &dt)
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
static Tensor< scalar > uniform(const scalar &s)
Return a VectorSpace with all elements = s.
dimensioned< scalar > mag(const dimensioned< Type > &)
static const Form zero
Definition: VectorSpace.H:110
const Cmpt & zx() const
Definition: TensorI.H:195
vector eigenVector(const tensor &, const scalar lambda)
Definition: tensor.C:190
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51