openfoam/applications/utilities/postProcessing/Field/magGrad/magGrad.C
2008-04-15 18:56:58 +01:00

129 lines
3.8 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 1991-2007 OpenCFD Ltd.
\\/ 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 2 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, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Application
magGrad
Description
Calculates and writes the scalar magnitude of a scalar or vector field
at each time
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// Main program:
int main(int argc, char *argv[])
{
argList::validArgs.append("field");
# include "addTimeOptions.H"
# include "setRootCase.H"
# include "createTime.H"
// Get times list
instantList Times = runTime.times();
// set startTime and endTime depending on -time and -latestTime options
# include "checkTimeOptions.H"
const word fieldName(args.additionalArgs()[0]);
runTime.setTime(Times[startTime], startTime);
# include "createMesh.H"
for (label i=startTime; i<endTime; i++)
{
runTime.setTime(Times[i], i);
Info<< "Time = " << runTime.timeName() << endl;
IOobject header
(
fieldName,
runTime.timeName(),
mesh,
IOobject::MUST_READ
);
// Check U exists
if (header.headerOk())
{
mesh.readUpdate();
if (header.headerClassName() == volVectorField::typeName)
{
Info<< " Reading " << fieldName << endl;
volVectorField U(header, mesh);
volScalarField magGrad
(
IOobject
(
"magGrad" + fieldName,
runTime.timeName(),
mesh,
IOobject::NO_READ
),
mag(fvc::grad(U))
);
Info<< " Calculating " << magGrad.name() << endl;
magGrad.write();
}
else if (header.headerClassName() == volScalarField::typeName)
{
Info<< " Reading " << fieldName << endl;
volScalarField U(header, mesh);
volScalarField magGrad
(
IOobject
(
"magGrad" + fieldName,
runTime.timeName(),
mesh,
IOobject::NO_READ
),
mag(fvc::grad(U))
);
Info<< " Calculating " << magGrad.name() << endl;
magGrad.write();
}
}
else
{
Info<< " No " << fieldName << endl;
}
}
return(0);
}
// ************************************************************************* //