tensor2D.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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "tensor2D.H"
29 #include "quadraticEqn.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 template<>
34 const char* const Foam::tensor2D::vsType::typeName = "tensor2D";
35 
36 template<>
37 const char* const Foam::tensor2D::vsType::componentNames[] =
38 {
39  "xx", "xy",
40  "yx", "yy"
41 };
42 
43 template<>
45 (
46  tensor2D::uniform(0)
47 );
48 
49 template<>
51 (
52  tensor2D::uniform(1)
53 );
54 
55 template<>
57 (
58  tensor2D::uniform(VGREAT)
59 );
60 
61 template<>
63 (
64  tensor2D::uniform(-VGREAT)
65 );
66 
67 template<>
69 (
70  tensor2D::uniform(ROOTVGREAT)
71 );
72 
73 template<>
75 (
76  tensor2D::uniform(-ROOTVGREAT)
77 );
78 
79 template<>
81 (
82  1, 0,
83  0, 1
84 );
85 
86 
87 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
88 
90 {
91  // Coefficients of the characteristic quadratic polynomial (a = 1)
92  const scalar b = - t.xx() - t.yy();
93  const scalar c = t.xx()*t.yy() - t.xy()*t.yx();
94 
95  // Solve
96  Roots<2> roots = quadraticEqn(1, b, c).roots();
97 
98  // Check the root types
99  vector2D lambda = vector2D::zero;
100  forAll(roots, i)
101  {
102  switch (roots.type(i))
103  {
104  case roots::real:
105  lambda[i] = roots[i];
106  break;
107  case roots::complex:
109  << "Complex eigenvalues detected for tensor: " << t
110  << endl;
111  lambda[i] = 0;
112  break;
113  case roots::posInf:
114  lambda[i] = VGREAT;
115  break;
116  case roots::negInf:
117  lambda[i] = - VGREAT;
118  break;
119  case roots::nan:
121  << "Eigenvalue calculation failed for tensor: " << t
122  << exit(FatalError);
123  }
124  }
125 
126  // Sort the eigenvalues into ascending order
127  if (lambda.x() > lambda.y())
128  {
129  Swap(lambda.x(), lambda.y());
130  }
131 
132  return lambda;
133 }
134 
135 
137 (
138  const tensor2D& T,
139  const scalar lambda,
140  const vector2D& direction1
141 )
142 {
143  // Construct the linear system for this eigenvalue
145 
146  // Evaluate the eigenvector using the largest divisor
147  if (mag(A.yy()) > mag(A.xx()) && mag(A.yy()) > SMALL)
148  {
149  vector2D ev(1, - A.yx()/A.yy());
150 
151  return ev/mag(ev);
152  }
153  else if (mag(A.xx()) > SMALL)
154  {
155  vector2D ev(- A.xy()/A.xx(), 1);
156 
157  return ev/mag(ev);
158  }
159 
160  // Repeated eigenvalue
161  return vector2D(- direction1.y(), direction1.x());
162 }
163 
164 
165 Foam::tensor2D Foam::eigenVectors(const tensor2D& T, const vector2D& lambdas)
166 {
167  vector2D Ux(1, 0), Uy(0, 1);
168 
169  Ux = eigenVector(T, lambdas.x(), Uy);
170  Uy = eigenVector(T, lambdas.y(), Ux);
171 
172  return tensor2D(Ux, Uy);
173 }
174 
175 
177 {
178  const vector2D lambdas(eigenValues(T));
179 
180  return eigenVectors(T, lambdas);
181 }
182 
183 
184 // ************************************************************************* //
Foam::VectorSpace::rootMax
static const Form rootMax
Definition: VectorSpace.H:119
Foam::roots::complex
Definition: Roots.H:57
Foam::roots::posInf
Definition: Roots.H:58
Foam::VectorSpace::one
static const Form one
Definition: VectorSpace.H:116
Foam::roots::nan
Definition: Roots.H:60
Foam::VectorSpace::componentNames
static const char *const componentNames[]
Definition: VectorSpace.H:114
Foam::eigenVectors
dimensionedTensor eigenVectors(const dimensionedTensor &dt)
Definition: dimensionedTensor.C:159
tensor2D.H
A
static const Foam::dimensionedScalar A("", Foam::dimPressure, 611.21)
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::Vector2D
Templated 2D Vector derived from VectorSpace adding construction from 2 components,...
Definition: Vector2D.H:55
Foam::Tensor2D::I
static const Tensor2D I
Definition: Tensor2D.H:78
Foam::Swap
void Swap(DynamicList< T, SizeMin1 > &a, DynamicList< T, SizeMin2 > &b)
Definition: DynamicListI.H:909
Foam::eigenVector
vector eigenVector(const tensor &T, const scalar lambda, const vector &direction1, const vector &direction2)
Definition: tensor.C:139
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::Tensor2D
Templated 2D tensor derived from VectorSpace adding construction from 4 components,...
Definition: Tensor2D.H:59
Foam::VectorSpace::min
static const Form min
Definition: VectorSpace.H:118
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::eigenValues
dimensionedVector eigenValues(const dimensionedTensor &dt)
Definition: dimensionedTensor.C:148
Foam::roots::negInf
Definition: Roots.H:59
Foam::roots::real
Definition: Roots.H:56
Foam::FatalError
error FatalError
lambda
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::VectorSpace::rootMin
static const Form rootMin
Definition: VectorSpace.H:120
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::tensor2D
Tensor2D< scalar > tensor2D
Tensor2D of scalars.
Definition: tensor2D.H:51
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::VectorSpace::max
static const Form max
Definition: VectorSpace.H:117
Foam::vector2D
Vector2D< scalar > vector2D
A 2D vector of scalars obtained from the generic Vector2D.
Definition: vector2D.H:51
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::VectorSpace::zero
static const Form zero
Definition: VectorSpace.H:115
quadraticEqn.H
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::I
static const Identity< scalar > I
Definition: Identity.H:95
Foam::VectorSpace::typeName
static const char *const typeName
Definition: VectorSpace.H:113