Newer
Older
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2013 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "Lambda2.H"
#include "volFields.H"
#include "dictionary.H"
#include "zeroGradientFvPatchFields.H"
#include "fvcGrad.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
defineTypeNameAndDebug(Lambda2, 0);
}
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::Lambda2::Lambda2
(
const word& name,
const objectRegistry& obr,
const dictionary& dict,
const bool loadFromFiles
)
:
name_(name),
obr_(obr),
active_(true),
UName_("U")
{
// Check if the available mesh is an fvMesh, otherwise deactivate
if (!isA<fvMesh>(obr_))
{
active_ = false;
WarningIn
(
"Lambda2::Lambda2"
"("
"const word&, "
"const objectRegistry&, "
"const dictionary&, "
"const bool"
")"
) << "No fvMesh available, deactivating." << nl
<< endl;
}
read(dict);
if (active_)
{
const fvMesh& mesh = refCast<const fvMesh>(obr_);
volScalarField* Lambda2Ptr
(
new volScalarField
(
IOobject
(
type(),
mesh.time().timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
mesh,
dimensionedScalar("0", dimless/sqr(dimTime), 0.0)
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
)
);
mesh.objectRegistry::store(Lambda2Ptr);
}
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::Lambda2::~Lambda2()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
void Foam::Lambda2::read(const dictionary& dict)
{
if (active_)
{
UName_ = dict.lookupOrDefault<word>("UName", "U");
}
}
void Foam::Lambda2::execute()
{
// Do nothing - only valid on write
}
void Foam::Lambda2::end()
{
// Do nothing - only valid on write
}
void Foam::Lambda2::timeSet()
{
// Do nothing - only valid on write
}
void Foam::Lambda2::write()
{
if (active_)
{
const fvMesh& mesh = refCast<const fvMesh>(obr_);
const volVectorField& U =
mesh.lookupObject<volVectorField>(UName_);
const volTensorField gradU(fvc::grad(U));
const volTensorField SSplusWW
(
(symm(gradU) & symm(gradU))
+ (skew(gradU) & skew(gradU))
);
volScalarField& Lambda2 =
const_cast<volScalarField&>
(
mesh.lookupObject<volScalarField>(type())
);
Lambda2 = -eigenValues(SSplusWW)().component(vector::Y);
Info<< type() << " " << name_ << " output:" << nl
<< " writing field " << Lambda2.name() << nl
<< endl;