functionObjects: layerAverage: Replacment for postChannel

This function generates plots of fields averaged over the layers in the
mesh. It is a generalised replacement for the postChannel utility, which
has been removed. An example of this function's usage is as follows:

    layerAverage1
    {
        type            layerAverage;
        libs            ("libfieldFunctionObjects.so");

        writeControl    writeTime;

        setFormat       raw;

        // Patches and/or zones from which layers extrude
        patches         (bottom);
        zones           (quarterPlane threeQuartersPlane);

        // Spatial component against which to plot
        component       y;

        // Is the geometry symmetric around the centre layer?
        symmetric       true;

        // Fields to average and plot
        fields          (pMean pPrime2Mean UMean UPrime2Mean k);
    }
This commit is contained in:
Will Bainbridge
2021-12-07 15:32:43 +00:00
parent 31fee136e1
commit 053eed714d
19 changed files with 829 additions and 935 deletions

View File

@ -81,4 +81,6 @@ cylindrical/cylindricalFunctionObject.C
uniform/uniform.C
layerAverage/layerAverage.C
LIB = $(FOAM_LIBBIN)/libfieldFunctionObjects

View File

@ -0,0 +1,420 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2021 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/>.
\*---------------------------------------------------------------------------*/
#include "layerAverage.H"
#include "regionSplit.H"
#include "syncTools.H"
#include "volFields.H"
#include "writeFile.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace functionObjects
{
defineTypeNameAndDebug(layerAverage, 0);
addToRunTimeSelectionTable(functionObject, layerAverage, dictionary);
}
template<>
const char*
Foam::NamedEnum<Foam::vector::components, 3>::names[] =
{"x", "y", "z"};
}
const Foam::NamedEnum<Foam::vector::components, 3>
Foam::functionObjects::layerAverage::axisNames_;
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void Foam::functionObjects::layerAverage::walkOppositeFaces
(
const labelList& startFaces,
const boolList& startFaceIntoOwners,
boolList& blockedFace,
boolList& cellIsLayer
)
{
// Initialise the front
DynamicList<label> frontFaces(startFaces);
DynamicList<bool> frontFaceIntoOwners(startFaceIntoOwners);
DynamicList<label> newFrontFaces(frontFaces.size());
DynamicList<bool> newFrontFaceIntoOwners(frontFaceIntoOwners.size());
// Set the start and end faces as blocked
UIndirectList<bool>(blockedFace, startFaces) = true;
// Iterate until the front is empty
while (returnReduce(frontFaces.size(), sumOp<label>()) > 0)
{
// Transfer front faces across couplings
boolList bndFaceIsFront(mesh_.nFaces() - mesh_.nInternalFaces(), false);
forAll(frontFaces, i)
{
const label facei = frontFaces[i];
const label intoOwner = frontFaceIntoOwners[i];
if (!mesh_.isInternalFace(facei) && !intoOwner)
{
bndFaceIsFront[facei - mesh_.nInternalFaces()] = true;
}
}
syncTools::swapBoundaryFaceList(mesh_, bndFaceIsFront);
// Set new front faces
forAll(bndFaceIsFront, i)
{
const label facei = mesh_.nInternalFaces() + i;
if (bndFaceIsFront[i] && !blockedFace[facei])
{
frontFaces.append(facei);
frontFaceIntoOwners.append(true);
}
}
// Transfer across cells
newFrontFaces.clear();
newFrontFaceIntoOwners.clear();
forAll(frontFaces, i)
{
const label facei = frontFaces[i];
const label intoOwner = frontFaceIntoOwners[i];
const label celli =
intoOwner
? mesh_.faceOwner()[facei]
: mesh_.isInternalFace(facei) ? mesh_.faceNeighbour()[facei] : -1;
if (celli != -1)
{
cellIsLayer[celli] = true;
const label oppositeFacei =
mesh_.cells()[celli]
.opposingFaceLabel(facei, mesh_.faces());
const bool oppositeIntoOwner =
mesh_.faceOwner()[oppositeFacei] != celli;
if (oppositeFacei == -1)
{
FatalErrorInFunction
<< "Cannot find face on cell " << mesh_.cells()[celli]
<< " opposing face " << mesh_.faces()[facei]
<< ". Mesh is not layered." << exit(FatalError);
}
else if (!blockedFace[oppositeFacei])
{
blockedFace[oppositeFacei] = true;
newFrontFaces.append(oppositeFacei);
newFrontFaceIntoOwners.append(oppositeIntoOwner);
}
}
}
frontFaces.transfer(newFrontFaces);
frontFaceIntoOwners.transfer(newFrontFaceIntoOwners);
}
}
void Foam::functionObjects::layerAverage::calcLayers()
{
DynamicList<label> startFaces;
DynamicList<bool> startFaceIntoOwners;
forAll(patchIDs_, i)
{
const polyPatch& pp = mesh_.boundaryMesh()[patchIDs_[i]];
startFaces.append(identity(pp.size()) + pp.start());
startFaceIntoOwners.append(boolList(pp.size(), true));
}
forAll(zoneIDs_, i)
{
const faceZone& fz = mesh_.faceZones()[zoneIDs_[i]];
startFaces.append(fz);
startFaceIntoOwners.append(fz.flipMap());
}
// Identify faces that separate the layers
boolList blockedFace(mesh_.nFaces(), false);
boolList cellIsLayer(mesh_.nCells(), false);
walkOppositeFaces
(
startFaces,
startFaceIntoOwners,
blockedFace,
cellIsLayer
);
// Do analysis for connected layers
regionSplit rs(mesh_, blockedFace);
nLayers_ = rs.nRegions();
cellLayer_.transfer(rs);
// Get rid of regions that are not layers
label layeri0 = labelMax, layeri1 = -labelMax;
forAll(cellLayer_, celli)
{
if (cellIsLayer[celli])
{
layeri0 = min(layeri0, cellLayer_[celli]);
layeri1 = max(layeri1, cellLayer_[celli]);
}
else
{
cellLayer_[celli] = -1;
}
}
reduce(layeri0, minOp<label>());
reduce(layeri1, maxOp<label>());
nLayers_ = layeri0 != labelMax ? 1 + layeri1 - layeri0 : 0;
forAll(cellLayer_, celli)
{
if (cellLayer_[celli] != -1)
{
cellLayer_[celli] -= layeri0;
}
}
// Report
if (nLayers_ != 0)
{
Info<< " Detected " << nLayers_ << " layers" << nl << endl;
}
else
{
WarningInFunction
<< "No layers detected" << endl;
}
// Sum number of entries per layer
layerCount_ = sum(scalarField(mesh_.nCells(), 1));
// Average the cell centres
const pointField layerCentres(sum(mesh_.cellCentres())/layerCount_);
// Sort by the direction component
x_ = layerCentres.component(axis_);
sortMap_ = identity(layerCentres.size());
stableSort(sortMap_, UList<scalar>::less(x_));
inplaceReorder(sortMap_, x_);
// If symmetric, keep only half of the coordinates
if (symmetric_)
{
x_.setSize(nLayers_/2);
}
}
template<>
Foam::vector
Foam::functionObjects::layerAverage::symmetricCoeff<Foam::vector>() const
{
vector result = vector::one;
result[axis_] = -1;
return result;
}
template<>
Foam::symmTensor
Foam::functionObjects::layerAverage::symmetricCoeff<Foam::symmTensor>() const
{
return sqr(symmetricCoeff<vector>());
}
template<>
Foam::tensor
Foam::functionObjects::layerAverage::symmetricCoeff<Foam::tensor>() const
{
return symmetricCoeff<symmTensor>();
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::functionObjects::layerAverage::layerAverage
(
const word& name,
const Time& runTime,
const dictionary& dict
)
:
fvMeshFunctionObject(name, runTime, dict)
{
read(dict);
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::functionObjects::layerAverage::~layerAverage()
{}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
bool Foam::functionObjects::layerAverage::read(const dictionary& dict)
{
Info<< type() << " " << name() << ":" << nl;
patchIDs_ =
mesh_.boundaryMesh().patchSet
(
dict.lookupOrDefault<wordReList>("patches", wordReList())
).toc();
zoneIDs_ =
findStrings
(
dict.lookupOrDefault<wordReList>("zones", wordReList()),
mesh_.faceZones().names()
);
if (patchIDs_.empty() && zoneIDs_.empty())
{
WarningInFunction
<< "No patches or zones specified" << endl;
}
symmetric_ = dict.lookupOrDefault<bool>("symmetric", false);
axis_ = axisNames_.read(dict.lookup("axis"));
fields_ = dict.lookup<wordList>("fields");
formatter_ = setWriter::New(dict.lookup("setFormat"), dict);
calcLayers();
return true;
}
Foam::wordList Foam::functionObjects::layerAverage::fields() const
{
return fields_;
}
bool Foam::functionObjects::layerAverage::execute()
{
return true;
}
bool Foam::functionObjects::layerAverage::write()
{
// Create list of available fields
wordList fieldNames;
forAll(fields_, fieldi)
{
if
(
false
#define FoundTypeField(Type, nullArg) \
|| foundObject<VolField<Type>>(fields_[fieldi])
FOR_ALL_FIELD_TYPES(FoundTypeField)
#undef FoundTypeField
)
{
fieldNames.append(fields_[fieldi]);
}
else
{
cannotFindObject(fields_[fieldi]);
}
}
// Collapse the fields
#define DeclareTypeValueSets(Type, nullArg) \
PtrList<Field<Type>> Type##ValueSets(fieldNames.size());
FOR_ALL_FIELD_TYPES(DeclareTypeValueSets);
#undef DeclareTypeValueSets
forAll(fieldNames, fieldi)
{
#define CollapseTypeFields(Type, nullArg) \
if (mesh_.foundObject<VolField<Type>>(fieldNames[fieldi])) \
{ \
const VolField<Type>& field = \
mesh_.lookupObject<VolField<Type>>(fieldNames[fieldi]); \
\
Type##ValueSets.set \
( \
fieldi, \
average(field.primitiveField()) \
); \
}
FOR_ALL_FIELD_TYPES(CollapseTypeFields);
#undef CollapseTypeFields
}
// Output directory
fileName outputPath =
mesh_.time().globalPath()/writeFile::outputPrefix/name();
if (mesh_.name() != fvMesh::defaultRegion)
{
outputPath = outputPath/mesh_.name();
}
// Write
formatter_->write
(
outputPath/mesh_.time().timeName(),
typeName,
coordSet(true, axisNames_[axis_], x_),
fieldNames
#define TypeValueSetsParameter(Type, nullArg) , Type##ValueSets
FOR_ALL_FIELD_TYPES(TypeValueSetsParameter)
#undef TypeValueSetsParameter
);
return true;
}
void Foam::functionObjects::layerAverage::updateMesh(const mapPolyMesh&)
{
Info<< type() << " " << name() << ":" << nl;
calcLayers();
}
void Foam::functionObjects::layerAverage::movePoints(const polyMesh&)
{
Info<< type() << " " << name() << ":" << nl;
calcLayers();
}
// ************************************************************************* //

View File

@ -0,0 +1,243 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2021 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/>.
Class
Foam::layerAverage
Description
Generates plots of fields averaged over the layers in the mesh
Example of function object specification:
\verbatim
layerAverage1
{
type layerAverage;
libs ("libfieldFunctionObjects.so");
writeControl writeTime;
setFormat raw;
patches (bottom);
zones (quarterPlane threeQuartersPlane);
component y;
symmetric true;
fields (pMean pPrime2Mean UMean UPrime2Mean k);
}
\endverbatim
Usage
\table
Property | Description | Required | Default value
type | Type name: layerAverage | yes |
setFormat | Output plotting format | yes |
patches | Patches that layers extrude from | no | ()
zones | Face zones that the layers extrude from | no | ()
component | Component of position to plot against | yes |
symmetric | Is the geometry symmetric around the centre layer? \
| no | false
field | Fields to average and plot | yes |
\endtable
SourceFiles
layerAverage.C
layerAverageTemplates.C
\*---------------------------------------------------------------------------*/
#ifndef layerAverage_H
#define layerAverage_H
#include "fvMeshFunctionObject.H"
#include "setWriter.H"
#include "boolList.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace functionObjects
{
/*---------------------------------------------------------------------------*\
Class layerAverage Declaration
\*---------------------------------------------------------------------------*/
class layerAverage
:
public fvMeshFunctionObject
{
// Private Static Data
//- Names of vector components
static const NamedEnum<vector::components, 3> axisNames_;
// Private Data
//- Patches which form the start of the layers
labelList patchIDs_;
//- Zones which form the start of the layers
labelList zoneIDs_;
//- Zones on which the layers are considered to end
labelList endZoneIDs_;
//- Is the case symmetric?
bool symmetric_;
//- The direction over which to plot the results
vector::components axis_;
//- Per cell the global layer
label nLayers_;
//- Per cell the global layer
labelList cellLayer_;
//- Per global layer the number of cells
scalarField layerCount_;
//- From sorted layer back to unsorted global layer
labelList sortMap_;
//- Sorted component of cell centres
scalarField x_;
//- Fields to sample
wordList fields_;
//- Set formatter
autoPtr<setWriter> formatter_;
// Private Member Functions
//- Walk through layers marking faces that separate layers and cells
// that are within layers
void walkOppositeFaces
(
const labelList& startFaces,
const boolList& startFaceIntoOwners,
boolList& blockedFace,
boolList& cellIsLayer
);
//- Create the layer information, the sort map, and the scalar axis
void calcLayers();
//- Return the coefficient to multiply onto symmetric values
template<class T>
T symmetricCoeff() const;
//- Sum field per layer
template<class T>
Field<T> sum(const Field<T>& cellField) const;
//- Average a field per layer
template<class T>
Field<T> average(const Field<T>& cellField) const;
public:
//- Runtime type information
TypeName("layerAverage");
// Constructors
//- Construct from Time and dictionary
layerAverage
(
const word& name,
const Time& runTime,
const dictionary& dict
);
//- Disallow default bitwise copy construction
layerAverage(const layerAverage&) = delete;
//- Destructor
virtual ~layerAverage();
// Member Functions
//- Read the field average data
virtual bool read(const dictionary&);
//- Return the list of fields required
virtual wordList fields() const;
//- Do nothing
virtual bool execute();
//- Calculate and write the graphs
virtual bool write();
//- Update for changes of mesh
virtual void updateMesh(const mapPolyMesh&);
//- Update for mesh point-motion
virtual void movePoints(const polyMesh&);
// Member Operators
//- Disallow default bitwise assignment
void operator=(const layerAverage&) = delete;
};
template<>
vector layerAverage::symmetricCoeff<vector>() const;
template<>
symmTensor layerAverage::symmetricCoeff<symmTensor>() const;
template<>
tensor layerAverage::symmetricCoeff<tensor>() const;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace functionObjects
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "layerAverageTemplates.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,89 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2021 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/>.
\*---------------------------------------------------------------------------*/
#include "layerAverage.H"
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
template<class T>
T Foam::functionObjects::layerAverage::symmetricCoeff() const
{
return pTraits<T>::one;
}
template<class T>
Foam::Field<T> Foam::functionObjects::layerAverage::sum
(
const Field<T>& cellField
) const
{
Field<T> layerField(nLayers_, Zero);
forAll(cellLayer_, celli)
{
if (cellLayer_[celli] != -1)
{
layerField[cellLayer_[celli]] += cellField[celli];
}
}
Pstream::listCombineGather(layerField, plusEqOp<T>());
Pstream::listCombineScatter(layerField);
return layerField;
}
template<class T>
Foam::Field<T> Foam::functionObjects::layerAverage::average
(
const Field<T>& cellField
) const
{
// Sum, average and sort
Field<T> layerField(sum(cellField)/layerCount_, sortMap_);
// Handle symmetry
if (symmetric_)
{
const T coeff = symmetricCoeff<T>();
for (label i=0; i<nLayers_/2; i++)
{
const label j = nLayers_ - i - 1;
layerField[i] =
(layerField[i] + cmptMultiply(coeff, layerField[j]))/2;
}
layerField.setSize(nLayers_/2);
}
return layerField;
}
// ************************************************************************* //