mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
Adjust some utilities to use postCalc wrapper or timeSelector directly.
When required, also adjusted to use XXXApp.C for the source name. Adjusted some names in preparation for merge with master.
This commit is contained in:
@ -37,31 +37,20 @@ Description
|
||||
|
||||
int main(int argc, char *argv[])
|
||||
{
|
||||
timeSelector::addOptions();
|
||||
argList::validArgs.append("fieldName");
|
||||
argList::validArgs.append("patchName");
|
||||
|
||||
# include "addTimeOptions.H"
|
||||
# include "setRootCase.H"
|
||||
# include "createTime.H"
|
||||
instantList timeDirs = timeSelector::select0(runTime, args);
|
||||
# include "createMesh.H"
|
||||
|
||||
word fieldName(args.additionalArgs()[0]);
|
||||
word patchName(args.additionalArgs()[1]);
|
||||
|
||||
# include "createTime.H"
|
||||
|
||||
// Get times list
|
||||
instantList Times = runTime.times();
|
||||
|
||||
// set startTime and endTime depending on -time and -latestTime options
|
||||
# include "checkTimeOptions.H"
|
||||
|
||||
runTime.setTime(Times[startTime], startTime);
|
||||
|
||||
# include "createMesh.H"
|
||||
|
||||
for (label i=startTime; i<endTime; i++)
|
||||
forAll(timeDirs, timeI)
|
||||
{
|
||||
runTime.setTime(Times[i], i);
|
||||
|
||||
runTime.setTime(timeDirs[timeI], timeI);
|
||||
Info<< "Time = " << runTime.timeName() << endl;
|
||||
|
||||
IOobject fieldHeader
|
||||
@ -77,21 +66,40 @@ int main(int argc, char *argv[])
|
||||
{
|
||||
mesh.readUpdate();
|
||||
|
||||
Info<< " Reading field " << fieldName << endl;
|
||||
volScalarField field(fieldHeader, mesh);
|
||||
|
||||
label patchi = mesh.boundaryMesh().findPatchID(patchName);
|
||||
|
||||
if (patchi >= 0)
|
||||
if (patchi < 0)
|
||||
{
|
||||
FatalError
|
||||
<< "Unable to find patch " << patchName << nl
|
||||
<< exit(FatalError);
|
||||
}
|
||||
|
||||
if (fieldHeader.headerClassName() == "volScalarField")
|
||||
{
|
||||
Info<< " Reading volScalarField " << fieldName << endl;
|
||||
volScalarField field(fieldHeader, mesh);
|
||||
|
||||
scalar area = sum(mesh.magSf().boundaryField()[patchi]);
|
||||
scalar sumField = 0;
|
||||
|
||||
if (area > 0)
|
||||
{
|
||||
sumField = sum
|
||||
(
|
||||
mesh.magSf().boundaryField()[patchi]
|
||||
* field.boundaryField()[patchi]
|
||||
) / area;
|
||||
}
|
||||
|
||||
Info<< " Average of " << fieldName << " over patch "
|
||||
<< patchName << '[' << patchi << ']' << " = "
|
||||
<< sum
|
||||
(
|
||||
mesh.magSf().boundaryField()[patchi]
|
||||
*field.boundaryField()[patchi]
|
||||
)/sum(mesh.magSf().boundaryField()[patchi])
|
||||
<< endl;
|
||||
<< sumField << endl;
|
||||
}
|
||||
else
|
||||
{
|
||||
FatalError
|
||||
<< "Only possible to average volScalarFields "
|
||||
<< nl << exit(FatalError);
|
||||
}
|
||||
}
|
||||
else
|
||||
@ -104,8 +112,7 @@ int main(int argc, char *argv[])
|
||||
|
||||
Info<< "End\n" << endl;
|
||||
|
||||
return(0);
|
||||
return 0;
|
||||
}
|
||||
|
||||
|
||||
// ************************************************************************* //
|
||||
|
||||
Reference in New Issue
Block a user