Files
OpenFOAM-12/applications/utilities/postProcessing/postProcess/postProcess.C
Henry Weller 8cb90dab4f functionObjectList: Improved error message if a functionObject configuration file is not foundation
e.g. if streamlinesLines is specified rather than the correct streamlinesLine in
the controlDict::functions list:

functions
{
    #includeFunc streamlinesLines
    (
        funcName=streamlines,
        start=(-0.0205 0.001 0.00001),
        end=(-0.0205 0.0251 0.00001),
        nPoints=10,
        fields=(p k U)
    )
}

the following error message is generated providing a list of all the
functionObject configuration files available and the run stops rather than
continuing regardless of the error:

--> FOAM FATAL IO ERROR:
Cannot find functionObject configuration file streamlinesLines

Available configured functionObjects:
88
(
CourantNo
Lambda2
MachNo
PecletNo
Q
Qdot
XiReactionRate
add
age
boundaryProbes
cellMax
cellMaxMag
cellMin
cellMinMag
components
cutPlaneSurface
ddt
div
divide
dsmcFields
enstrophy
faceZoneAverage
faceZoneFlowRate
fieldAverage
flowType
forceCoeffsCompressible
forceCoeffsIncompressible
forcesCompressible
forcesIncompressible
grad
graphCell
graphCellFace
graphFace
graphLayerAverage
graphUniform
interfaceHeight
internalProbes
isoSurface
log
mag
magSqr
moments
multiply
particles
patchAverage
patchDifference
patchFlowRate
patchIntegrate
patchSurface
phaseForces
phaseMap
phaseScalarTransport
probes
randomise
residuals
scalarTransport
scale
shearStress
sizeDistribution
staticPressureIncompressible
stopAtClockTime
stopAtFile
streamFunction
streamlinesLine
streamlinesPatch
streamlinesPoints
streamlinesSphere
subtract
surfaceInterpolate
time
timeStep
totalEnthalpy
totalPressureCompressible
totalPressureIncompressible
triSurfaceDifference
triSurfaceVolumetricFlowRate
turbulenceFields
turbulenceIntensity
uniform
vorticity
wallHeatFlux
wallHeatTransferCoeff
wallShearStress
writeCellCentres
writeCellVolumes
writeObjects
writeVTK
yPlus
)

file: /home/dm2/henry/OpenFOAM/OpenFOAM-dev/tutorials/incompressible/simpleFoam/pitzDaily/system/controlDict/functions

    From function static bool Foam::functionObjectList::readFunctionObject(const Foam::string&, Foam::dictionary&, const Foam::Pair<Foam::string>&, const Foam::word&)
    in file db/functionObjects/functionObjectList/functionObjectList.C at line 250.

FOAM exiting
2022-01-26 12:02:07 +00:00

234 lines
6.7 KiB
C++

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2016-2022 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/>.
Application
postProcess
Description
Execute the set of functionObjects specified in the selected dictionary
(which defaults to system/controlDict) or on the command-line for the
selected set of times on the selected set of fields.
\*---------------------------------------------------------------------------*/
#include "argList.H"
#include "timeSelector.H"
#include "ReadFields.H"
#include "volFields.H"
#include "surfaceFields.H"
#include "pointFields.H"
#include "uniformDimensionedFields.H"
using namespace Foam;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#define ReadFields(GeoFieldType) \
readFields<GeoFieldType>(mesh, objects, requiredFields, storedObjects);
#define ReadPointFields(GeoFieldType) \
readFields<GeoFieldType>(pMesh, objects, requiredFields, storedObjects);
#define ReadUniformFields(FieldType) \
readUniformFields<FieldType> \
(constantObjects, requiredFields, storedObjects);
void executeFunctionObjects
(
const argList& args,
const Time& runTime,
fvMesh& mesh,
const HashSet<word>& requiredFields0,
functionObjectList& functions,
bool lastTime
)
{
Info<< nl << "Reading fields:" << endl;
// Maintain a stack of the stored objects to clear after executing
// the functionObjects
LIFOStack<regIOobject*> storedObjects;
// Read objects in time directory
IOobjectList objects(mesh, runTime.timeName());
HashSet<word> requiredFields(requiredFields0);
forAll(functions, i)
{
requiredFields.insert(functions[i].fields());
}
// Read volFields
ReadFields(volScalarField);
ReadFields(volVectorField);
ReadFields(volSphericalTensorField);
ReadFields(volSymmTensorField);
ReadFields(volTensorField);
// Read internal fields
ReadFields(volScalarField::Internal);
ReadFields(volVectorField::Internal);
ReadFields(volSphericalTensorField::Internal);
ReadFields(volSymmTensorField::Internal);
ReadFields(volTensorField::Internal);
// Read surface fields
ReadFields(surfaceScalarField);
ReadFields(surfaceVectorField);
ReadFields(surfaceSphericalTensorField);
ReadFields(surfaceSymmTensorField);
ReadFields(surfaceTensorField);
// Read point fields.
const pointMesh& pMesh = pointMesh::New(mesh);
ReadPointFields(pointScalarField)
ReadPointFields(pointVectorField);
ReadPointFields(pointSphericalTensorField);
ReadPointFields(pointSymmTensorField);
ReadPointFields(pointTensorField);
// Read uniform dimensioned fields
IOobjectList constantObjects(mesh, runTime.constant());
ReadUniformFields(uniformDimensionedScalarField);
ReadUniformFields(uniformDimensionedVectorField);
ReadUniformFields(uniformDimensionedSphericalTensorField);
ReadUniformFields(uniformDimensionedSymmTensorField);
ReadUniformFields(uniformDimensionedTensorField);
Info<< nl << "Executing functionObjects" << endl;
// Execute the functionObjects in post-processing mode
functions.execute();
// Execute the functionObject 'end()' function for the last time
if (lastTime)
{
functions.end();
}
while (!storedObjects.empty())
{
storedObjects.pop()->checkOut();
}
}
int main(int argc, char *argv[])
{
Foam::timeSelector::addOptions();
#include "addRegionOption.H"
#include "addFunctionObjectOptions.H"
// Set functionObject post-processing mode
functionObject::postProcess = true;
#include "setRootCase.H"
if (args.optionFound("list"))
{
Info<< nl
<< "Available configured functionObjects:"
<< functionObjectList::list()
<< endl;
return 0;
}
#include "createTime.H"
Foam::instantList timeDirs = Foam::timeSelector::select0(runTime, args);
#include "createNamedMesh.H"
// Initialise the set of selected fields from the command-line options
HashSet<word> requiredFields;
if (args.optionFound("fields"))
{
args.optionLookup("fields")() >> requiredFields;
}
if (args.optionFound("field"))
{
requiredFields.insert(args.optionLookup("field")());
}
// Externally stored dictionary for functionObjectList
// if not constructed from runTime
dictionary functionsControlDict("controlDict");
// Construct functionObjectList
autoPtr<functionObjectList> functionsPtr
(
functionObjectList::New
(
args,
runTime,
functionsControlDict
)
);
forAll(timeDirs, timei)
{
runTime.setTime(timeDirs[timei], timei);
Info<< "Time = " << runTime.userTimeName() << endl;
if (mesh.readUpdate() != polyMesh::UNCHANGED)
{
// Update functionObjectList if mesh changes
functionsPtr = functionObjectList::New
(
args,
runTime,
functionsControlDict
);
}
FatalIOError.throwExceptions();
try
{
executeFunctionObjects
(
args,
runTime,
mesh,
requiredFields,
functionsPtr(),
timei == timeDirs.size()-1
);
}
catch (IOerror& err)
{
Warning<< err << endl;
}
Info<< endl;
}
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //