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-2024 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 {
37 
39  (
41  Lambda2,
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  store
65  (
67  -eigenValues(SSplusWW)().component(vector::Y)
68  );
69 
70  return true;
71  }
72  else
73  {
74  cannotFindObject<volVectorField>(fieldName_);
75 
76  return false;
77  }
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
82 
84 (
85  const word& name,
86  const Time& runTime,
87  const dictionary& dict
88 )
89 :
90  fieldExpression(name, runTime, dict, typeName, "U")
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
95 
97 {}
98 
99 
100 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base-class for Time/database functionObjects.
Calculates and outputs the second largest eigenvalue of the sum of the square of the symmetric and an...
Definition: Lambda2.H:59
virtual ~Lambda2()
Destructor.
Definition: Lambda2.C:96
Lambda2(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: Lambda2.C:84
word resultName_
Name of result field.
const word fieldName_
Name of field to process.
ObjectType & store(const tmp< ObjectType > &tfield)
Store the given field in the objectRegistry.
A class for handling words, derived from string.
Definition: word.H:62
Calculate the gradient of the given field.
U
Definition: pEqn.H:72
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
Namespace for OpenFOAM.
void skew(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
VolField< vector > volVectorField
Definition: volFieldsFwd.H:63
void eigenValues(LagrangianPatchField< vector > &f, const LagrangianPatchField< tensor > &f1)
VolField< tensor > volTensorField
Definition: volFieldsFwd.H:66
void component(LagrangianPatchField< typename LagrangianPatchField< Type >::cmptType > &sf, const LagrangianPatchField< Type > &f, const direction d)
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
dictionary dict