ENH: patchAverage,patchIntegrate: generalise

This commit is contained in:
mattijs
2012-01-25 15:32:30 +00:00
parent 242143be1b
commit 40f4146ce2
2 changed files with 206 additions and 70 deletions

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2012 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -32,17 +32,59 @@ Description
#include "fvCFD.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
template<class FieldType>
void printAverage
(
const fvMesh& mesh,
const IOobject& fieldHeader,
const scalar area,
const label patchI,
bool& done
)
{
if (!done && fieldHeader.headerClassName() == FieldType::typeName)
{
Info<< " Reading " << fieldHeader.headerClassName() << " "
<< fieldHeader.name() << endl;
FieldType field(fieldHeader, mesh);
typename FieldType::value_type sumField =
pTraits<typename FieldType::value_type>::zero;
if (area > 0)
{
sumField = gSum
(
mesh.magSf().boundaryField()[patchI]
* field.boundaryField()[patchI]
) / area;
}
Info<< " Average of " << fieldHeader.headerClassName()
<< " over patch "
<< mesh.boundary()[patchI].name()
<< '[' << patchI << ']' << " = "
<< sumField << endl;
done = true;
}
}
// Main program:
int main(int argc, char *argv[])
{
timeSelector::addOptions();
#include "addRegionOption.H"
argList::validArgs.append("fieldName");
argList::validArgs.append("patchName");
# include "setRootCase.H"
# include "createTime.H"
instantList timeDirs = timeSelector::select0(runTime, args);
# include "createMesh.H"
# include "createNamedMesh.H"
const word fieldName = args[1];
const word patchName = args[2];
@ -52,7 +94,7 @@ int main(int argc, char *argv[])
runTime.setTime(timeDirs[timeI], timeI);
Info<< "Time = " << runTime.timeName() << endl;
IOobject fieldHeader
IOobject io
(
fieldName,
runTime.timeName(),
@ -61,7 +103,7 @@ int main(int argc, char *argv[])
);
// Check field exists
if (fieldHeader.headerOk())
if (io.headerOk())
{
mesh.readUpdate();
@ -72,32 +114,21 @@ int main(int argc, char *argv[])
<< "Unable to find patch " << patchName << nl
<< exit(FatalError);
}
scalar area = gSum(mesh.magSf().boundaryField()[patchI]);
if (fieldHeader.headerClassName() == "volScalarField")
{
Info<< " Reading volScalarField " << fieldName << endl;
volScalarField field(fieldHeader, mesh);
bool done = false;
printAverage<volScalarField>(mesh, io, area, patchI, done);
printAverage<volVectorField>(mesh, io, area, patchI, done);
printAverage<volSphericalTensorField>(mesh, io, area, patchI, done);
printAverage<volSymmTensorField>(mesh, io, area, patchI, done);
printAverage<volTensorField>(mesh, io, area, patchI, done);
scalar area = gSum(mesh.magSf().boundaryField()[patchI]);
scalar sumField = 0;
if (area > 0)
{
sumField = gSum
(
mesh.magSf().boundaryField()[patchI]
* field.boundaryField()[patchI]
) / area;
}
Info<< " Average of " << fieldName << " over patch "
<< patchName << '[' << patchI << ']' << " = "
<< sumField << endl;
}
else
if (!done)
{
FatalError
<< "Only possible to average volScalarFields "
<< "Only possible to average volFields."
<< " Field " << fieldName << " is of type "
<< io.headerClassName()
<< nl << exit(FatalError);
}
}