This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/utilities/postProcessing/velocityField/Q/Q.C

114 lines
3.2 KiB
C++
Raw Normal View History

/*---------------------------------------------------------------------------*\
========= |
2013-12-11 16:09:41 +00:00
\\ / F ield | foam-extend: Open Source CFD
\\ / O peration | Version: 4.1
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
2013-12-11 16:09:41 +00:00
This file is part of foam-extend.
2013-12-11 16:09:41 +00:00
foam-extend is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
2013-12-11 16:09:41 +00:00
Free Software Foundation, either version 3 of the License, or (at your
option) any later version.
2013-12-11 16:09:41 +00:00
foam-extend 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
2013-12-11 16:09:41 +00:00
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
Application
Q
Description
Calculates and writes the second invariant of the velocity gradient tensor.
2010-08-26 14:22:03 +00:00
The -noWrite option just outputs the max/min values without writing
the field.
\*---------------------------------------------------------------------------*/
#include "calc.H"
#include "fvc.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
{
2010-08-26 14:22:03 +00:00
bool writeResults = !args.optionFound("noWrite");
IOobject Uheader
(
"U",
runTime.timeName(),
mesh,
IOobject::MUST_READ
);
if (Uheader.headerOk())
{
Info<< " Reading U" << endl;
volVectorField U(Uheader, mesh);
volTensorField gradU = fvc::grad(U);
volScalarField Q
(
IOobject
(
"Q",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
0.5*(sqr(tr(gradU)) - tr(((gradU)&(gradU))))
);
/*
// This is a second way of calculating Q, that delivers results
// very close, but not identical to the first approach.
volSymmTensorField S = symm(gradU); // symmetric part of tensor
volTensorField W = skew(gradU); // anti-symmetric part
volScalarField SS = S&&S;
volScalarField WW = W&&W;
volScalarField Q
(
IOobject
(
"Q",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
0.5*(WW - SS)
);
*/
Info<< "mag(Q) max/min : "
<< max(Q).value() << " "
<< min(Q).value() << endl;
if (writeResults)
{
Q.write();
}
}
else
{
Info<< " No U" << endl;
}
2010-08-26 14:22:03 +00:00
Info<< "\nEnd\n" << endl;
}
2010-08-26 14:22:03 +00:00
// ************************************************************************* //