ENH: Added new fieldExtents function object

Description
    Calculates the spatial minimum and maximum extents of a field

    The extents are derived from the bound box limits after identifying
    the locations where field values exceed the user-supplied threshold
    value.

Usage
    Example of function object specification:

    fieldExtents1
    {
        type        fieldExtents;
        libs        ("libfieldFunctionObjects.so");
        ...
        writeToFile yes;
        log         yes;
        fields      (alpha);
        threshold   0.5;
        patches     ();
    }

    Where the entries comprise:

        Property      | Description              | Required   | Default
        type          | type name: fieldExtents  | yes        |
        writeToFile   | write extents data to file | no       | yes
        log           | write extents data to standard output | no | yes
        internalField | Process the internal field | no       | yes
        threshold     | Field value to identify extents boundary | yes |
        referencePosition | Reference position   | no         | (0 0 0)
        fields        | list of fields to process | yes       |
        patches       | list of patches to process | no       | <all>

    Output data is written to the file \<timeDir\>/fieldExtents.dat

Note
    For non-scalar fields, the magnitude of the field is employed and
    compared to the threshold value.
This commit is contained in:
Andrew Heather
2018-10-16 13:01:15 +01:00
parent 5d7b2329dd
commit 1eba709319
5 changed files with 542 additions and 0 deletions

View File

@ -3,6 +3,7 @@ fieldAverage/fieldAverageItem/fieldAverageItem.C
fieldAverage/fieldAverageItem/fieldAverageItemIO.C fieldAverage/fieldAverageItem/fieldAverageItemIO.C
fieldCoordinateSystemTransform/fieldCoordinateSystemTransform.C fieldCoordinateSystemTransform/fieldCoordinateSystemTransform.C
fieldExtents/fieldExtents.C
fieldMinMax/fieldMinMax.C fieldMinMax/fieldMinMax.C
fieldValues/fieldValue/fieldValue.C fieldValues/fieldValue/fieldValue.C

View File

@ -0,0 +1,201 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2018 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 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 "fieldExtents.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace functionObjects
{
defineTypeNameAndDebug(fieldExtents, 0);
addToRunTimeSelectionTable(functionObject, fieldExtents, dictionary);
}
}
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
void Foam::functionObjects::fieldExtents::writeFileHeader(Ostream& os)
{
if (!fieldSet_.updateSelection())
{
return;
}
if (writtenHeader_)
{
writeBreak(os);
}
else
{
writeHeader(os, "Field extents");
writeHeaderValue(os, "Reference position", C0_);
writeHeaderValue(os, "Threshold", threshold_);
}
writeCommented(os, "Time");
forAllConstIters(fieldSet_.selection(), iter)
{
const word& fieldName = iter();
if (internalField_)
{
writeTabbed(os, fieldName + "_internal");
}
for (const label patchi : patchIDs_)
{
const word& patchName = mesh_.boundaryMesh()[patchi].name();
writeTabbed(os, fieldName + "_" + patchName);
}
}
os << endl;
writtenHeader_ = true;
}
template<>
Foam::tmp<Foam::volScalarField> Foam::functionObjects::fieldExtents::calcMask
(
const GeometricField<scalar, fvPatchField, volMesh>& field
) const
{
return
pos
(
field
- dimensionedScalar("t", field.dimensions(), threshold_)
);
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::functionObjects::fieldExtents::fieldExtents
(
const word& name,
const Time& runTime,
const dictionary& dict
)
:
fvMeshFunctionObject(name, runTime, dict),
writeFile(mesh_, name, typeName, dict),
internalField_(true),
threshold_(0),
C0_(Zero),
fieldSet_(mesh_),
patchIDs_()
{
read(dict);
// Note: delay creating the output file header to handle field names
// specified using regular expressions
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
bool Foam::functionObjects::fieldExtents::read(const dictionary& dict)
{
if (fvMeshFunctionObject::read(dict) && writeFile::read(dict))
{
threshold_ = dict.get<scalar>("threshold");
dict.readIfPresent<bool>("internalField", internalField_);
dict.readIfPresent<vector>("referencePosition", C0_);
patchIDs_.clear();
const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
wordReList patchNames;
if (dict.readIfPresent("patches", patchNames))
{
for (const wordRe& name : patchNames)
{
patchIDs_.insert(pbm.findIndices(name));
}
}
else
{
// Add all non-processor and non-empty patches
forAll(pbm, patchi)
{
const polyPatch& pp = pbm[patchi];
if (!isA<processorPolyPatch>(pp) && !isA<emptyPolyPatch>(pp))
{
patchIDs_.insert(patchi);
}
}
}
if (!internalField_ && patchIDs_.empty())
{
IOWarningInFunction(dict)
<< "No internal field or patches selected - no field extent "
<< "information will be generated" << endl;
}
fieldSet_.read(dict);
return true;
}
return false;
}
bool Foam::functionObjects::fieldExtents::execute()
{
return true;
}
bool Foam::functionObjects::fieldExtents::write()
{
writeFileHeader(file());
Log << type() << " " << name() << " write:" << nl;
for (const word& fieldName : fieldSet_.selection())
{
calcFieldExtents<scalar>(fieldName, true);
calcFieldExtents<vector>(fieldName);
calcFieldExtents<sphericalTensor>(fieldName);
calcFieldExtents<symmTensor>(fieldName);
calcFieldExtents<tensor>(fieldName);
}
Log << endl;
return true;
}
// ************************************************************************* //

View File

@ -0,0 +1,207 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2018 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 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::functionObjects::fieldExtents
Group
grpFieldFunctionObjects
Description
Calculates the spatial minimum and maximum extents of a field
The extents are derived from the bound box limits after identifying the
locations where field values exceed the user-supplied threshold value.
Usage
Example of function object specification:
\verbatim
fieldExtents1
{
type fieldExtents;
libs ("libfieldFunctionObjects.so");
...
writeToFile yes;
log yes;
fields (alpha);
threshold 0.5;
patches ();
}
\endverbatim
Where the entries comprise:
\table
Property | Description | Required | Default value
type | type name: fieldExtents | yes |
writeToFile | write extents data to file | no | yes
log | write extents data to standard output | no | yes
internalField | Process the internal field | no | yes
threshold | Field value to identify extents boundary | yes |
referencePosition | Reference position | no | (0 0 0)
fields | list of fields to process | yes |
patches | list of patches to process | no | \<all patches\>
\endtable
Output data is written to the file \<timeDir\>/fieldExtents.dat
Note
For non-scalar fields, the magnitude of the field is employed and compared
to the threshold value.
See also
Foam::functionObjects::fvMeshFunctionObject
Foam::functionObjects::writeFile
SourceFiles
fieldExtents.C
\*---------------------------------------------------------------------------*/
#ifndef functionObjects_fieldExtents_H
#define functionObjects_fieldExtents_H
#include "fvMeshFunctionObject.H"
#include "writeFile.H"
#include "volFieldSelection.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace functionObjects
{
/*---------------------------------------------------------------------------*\
Class fieldExtents Declaration
\*---------------------------------------------------------------------------*/
class fieldExtents
:
public fvMeshFunctionObject,
public writeFile
{
protected:
// Protected data
//- Flag to write the internal field extents
bool internalField_;
//- Threshold value
scalar threshold_;
//- Reference position; default = (0 0 0)
point C0_;
//- Fields to assess
volFieldSelection fieldSet_;
//- Patches to assess
labelHashSet patchIDs_;
// Protected Member Functions
//- Output file header information
virtual void writeFileHeader(Ostream& os);
//- Return the field mask
template<class Type>
tmp<volScalarField> calcMask
(
const GeometricField<Type, fvPatchField, volMesh>& field
) const;
//- Main calculation
template<class Type>
void calcFieldExtents
(
const word& fieldName,
const bool calcMag = false
);
//- No copy construct
fieldExtents(const fieldExtents&) = delete;
//- No copy assignment
void operator=(const fieldExtents&) = delete;
public:
//- Runtime type information
TypeName("fieldExtents");
// Constructors
//- Construct from Time and dictionary
fieldExtents
(
const word& name,
const Time& runTime,
const dictionary& dict
);
//- Destructor
virtual ~fieldExtents() = default;
// Member Functions
//- Read the field extents data
virtual bool read(const dictionary&);
//- Execute, currently does nothing
virtual bool execute();
//- Write the fieldExtents
virtual bool write();
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
template<>
tmp<volScalarField> fieldExtents::calcMask
(
const GeometricField<scalar, fvPatchField, volMesh>& field
) const;
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace functionObjects
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "fieldExtentsTemplates.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,120 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2018 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 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 "fieldExtents.H"
#include "volFields.H"
#include "boundBox.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
template<class Type>
Foam::tmp<Foam::volScalarField> Foam::functionObjects::fieldExtents::calcMask
(
const GeometricField<Type, fvPatchField, volMesh>& field
) const
{
return
pos
(
mag(field)
- dimensionedScalar("t", field.dimensions(), threshold_)
);
}
template<class Type>
void Foam::functionObjects::fieldExtents::calcFieldExtents
(
const word& fieldName,
const bool calcMag
)
{
typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
const VolFieldType* fieldPtr =
obr_.lookupObjectPtr<VolFieldType>(fieldName);
if (!fieldPtr)
{
return;
}
auto extents = [this](const scalarField& mask, const vectorField& C)
{
boundBox extents(boundBox::invertedBox);
forAll(mask, i)
{
if (mask[i] > 0.5)
{
extents.add(C[i] - C0_);
}
};
extents.reduce();
if (extents.empty())
{
extents.add(point::zero);
}
return extents;
};
Log << "field: " << fieldName << nl;
writeTime(file());
tmp<volScalarField> tmask = calcMask<Type>(*fieldPtr);
const volScalarField& mask = tmask();
// Internal field
if (internalField_)
{
boundBox bb(extents(mask, mesh_.C()));
Log << " internal field: " << bb << nl;
file() << bb;
this->setResult(fieldName + "_internal_min" , bb.min());
this->setResult(fieldName + "_internal_max", bb.max());
}
// Patches
for (const label patchi : patchIDs_)
{
const fvPatchScalarField& maskp = mask.boundaryField()[patchi];
boundBox bb(extents(maskp, maskp.patch().Cf()));
const word& patchName = maskp.patch().name();
Log << " patch " << patchName << ": " << bb << nl;
file() << bb;
this->setResult(fieldName + "_" + patchName + "_min", bb.min());
this->setResult(fieldName + "_" + patchName + "_max", bb.max());
}
Log << endl;
file() << endl;
}
// ************************************************************************* //

View File

@ -51,5 +51,18 @@ maxCo 0.3;
maxDeltaT 1e-03; maxDeltaT 1e-03;
functions
{
fieldExtents
{
type fieldExtents;
libs ("libfieldFunctionObjects.so");
region wallFilmRegion;
threshold 0.5;
fields (alpha);
internalField yes;
patches ();
}
}
// ************************************************************************* // // ************************************************************************* //