Lambda2.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) 2013-2020 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 "Lambda2.H"
27 #include "fvcGrad.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace functionObjects
35 {
36  defineTypeNameAndDebug(Lambda2, 0);
37 
39  (
40  functionObject,
41  Lambda2,
42  dictionary
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 bool Foam::functionObjects::Lambda2::calc()
51 {
52  if (foundObject<volVectorField>(fieldName_))
53  {
54  const volVectorField& U = lookupObject<volVectorField>(fieldName_);
55  const tmp<volTensorField> tgradU(fvc::grad(U));
56  const volTensorField& gradU = tgradU();
57 
58  const volTensorField SSplusWW
59  (
60  (symm(gradU) & symm(gradU))
61  + (skew(gradU) & skew(gradU))
62  );
63 
64  return store
65  (
66  resultName_,
67  -eigenValues(SSplusWW)().component(vector::Y)
68  );
69  }
70  else
71  {
72  cannotFindObject<volVectorField>(fieldName_);
73 
74  return false;
75  }
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80 
82 (
83  const word& name,
84  const Time& runTime,
85  const dictionary& dict
86 )
87 :
88  fieldExpression(name, runTime, dict, typeName, "U")
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
93 
95 {}
96 
97 
98 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensionedTensor skew(const dimensionedTensor &dt)
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:63
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
dimensionedVector eigenValues(const dimensionedTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
Calculate the gradient of the given field.
virtual ~Lambda2()
Destructor.
Definition: Lambda2.C:94
A class for handling words, derived from string.
Definition: word.H:59
Lambda2(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: Lambda2.C:82
defineTypeNameAndDebug(Qdot, 0)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
Namespace for OpenFOAM.