ENH: residuals function object - extended to write residual fields

Residual fields can be written using the new 'writeFields' entry, e.g.

    functions
    {
        residual
        {
            type            residuals;
            libs            ("libutilityFunctionObjects.so");
            fields          (".*");
            writeControl    writeTime;
            writeFields     true;
        }
    }

Fields currently correspond to the initial residual for the last solver
iteration.
This commit is contained in:
Andrew Heather
2018-01-16 12:13:11 +00:00
parent 08193a50fa
commit 6312f80918
14 changed files with 265 additions and 44 deletions

View File

@ -79,6 +79,69 @@ void Foam::functionObjects::residuals::writeFileHeader(Ostream& os)
}
void Foam::functionObjects::residuals::createField(const word& fieldName)
{
if (!writeFields_)
{
return;
}
const word residualName("initialResidual:" + fieldName);
if (!mesh_.foundObject<IOField<scalar>>(residualName))
{
IOField<scalar>* fieldPtr =
new IOField<scalar>
(
IOobject
(
residualName,
mesh_.time().timeName(),
mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE
),
Field<scalar>(mesh_.nCells(), scalar(0))
);
fieldPtr->store();
}
}
void Foam::functionObjects::residuals::writeField(const word& fieldName) const
{
const word residualName("initialResidual:" + fieldName);
const IOField<scalar>* residualPtr =
mesh_.lookupObjectPtr<IOField<scalar>>(residualName);
if (residualPtr)
{
volScalarField residual
(
IOobject
(
residualName,
mesh_.time().timeName(),
mesh_,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
mesh_,
dimensionedScalar("0", dimless, 0),
zeroGradientFvPatchField<scalar>::typeName
);
residual.primitiveFieldRef() = *residualPtr;
residual.correctBoundaryConditions();
residual.write();
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::functionObjects::residuals::residuals
@ -90,7 +153,9 @@ Foam::functionObjects::residuals::residuals
:
fvMeshFunctionObject(name, runTime, dict),
writeFile(obr_, name, typeName, dict),
fieldSet_(mesh_)
fieldSet_(mesh_),
writeFields_(false),
initialised_(false)
{
read(dict);
}
@ -109,6 +174,9 @@ bool Foam::functionObjects::residuals::read(const dictionary& dict)
if (fvMeshFunctionObject::read(dict))
{
fieldSet_.read(dict);
writeFields_ = dict.lookupOrDefault("writeFields", false);
return true;
}
@ -118,32 +186,48 @@ bool Foam::functionObjects::residuals::read(const dictionary& dict)
bool Foam::functionObjects::residuals::execute()
{
return true;
}
bool Foam::functionObjects::residuals::write()
{
if (Pstream::master())
// Note: delaying initialisation until after first iteration so that
// we can find wildcard fields
if (!initialised_)
{
writeFileHeader(file());
writeTime(file());
for (const word& fieldName : fieldSet_.selection())
if (writeFields_)
{
writeResidual<scalar>(fieldName);
writeResidual<vector>(fieldName);
writeResidual<sphericalTensor>(fieldName);
writeResidual<symmTensor>(fieldName);
writeResidual<tensor>(fieldName);
for (const word& fieldName : fieldSet_.selection())
{
initialiseField<scalar>(fieldName);
initialiseField<vector>(fieldName);
initialiseField<sphericalTensor>(fieldName);
initialiseField<symmTensor>(fieldName);
initialiseField<tensor>(fieldName);
}
}
file() << endl;
initialised_ = true;
}
return true;
}
bool Foam::functionObjects::residuals::write()
{
writeTime(file());
for (const word& fieldName : fieldSet_.selection())
{
writeResidual<scalar>(fieldName);
writeResidual<vector>(fieldName);
writeResidual<sphericalTensor>(fieldName);
writeResidual<symmTensor>(fieldName);
writeResidual<tensor>(fieldName);
}
file() << endl;
return true;
}
// ************************************************************************* //

View File

@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2015-2016 OpenFOAM Foundation
\\/ M anipulation | Copyright (C) 2015-2017 OpenCFD Ltd.
\\/ M anipulation | Copyright (C) 2015-2018 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -87,17 +87,33 @@ protected:
//- Fields to write residuals
solverFieldSelection fieldSet_;
//- Flag to write the residual as a vol field
bool writeFields_;
//- Initialisation flag
bool initialised_;
// Protected Member Functions
//- Output file header information
void writeFileHeader(Ostream& os);
//- Create and store a residual field on the mesh database
void createField(const word& fieldName);
//- Write a residual field
void writeField(const word& fieldName) const;
//- Output file header information per primitive type value
template<class Type>
void writeFileHeader(Ostream& os, const word& fileName) const;
//- Calculate the field min/max
//- Initialise a residual field
template<class Type>
void initialiseField(const word& fieldName);
//- Calculate the residual
template<class Type>
void writeResidual(const word& fieldName);

View File

@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2015-2016 OpenFOAM Foundation
\\/ M anipulation | Copyright (C) 2015-2016 OpenCFD Ltd.
\\/ M anipulation | Copyright (C) 2015-2018 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -26,6 +26,7 @@ License
#include "residuals.H"
#include "volFields.H"
#include "ListOps.H"
#include "zeroGradientFvPatchField.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -60,12 +61,43 @@ void Foam::functionObjects::residuals::writeFileHeader
}
template<class Type>
void Foam::functionObjects::residuals::initialiseField(const word& fieldName)
{
typedef GeometricField<Type, fvPatchField, volMesh> volFieldType;
if (foundObject<volFieldType>(fieldName))
{
const Foam::dictionary& solverDict = mesh_.solverPerformanceDict();
if (solverDict.found(fieldName))
{
typename pTraits<Type>::labelType validComponents
(
mesh_.validComponents<Type>()
);
for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; ++cmpt)
{
if (component(validComponents, cmpt) != -1)
{
const word resultName =
fieldName + word(pTraits<Type>::componentNames[cmpt]);
createField(resultName);
}
}
}
}
}
template<class Type>
void Foam::functionObjects::residuals::writeResidual(const word& fieldName)
{
typedef GeometricField<Type, fvPatchField, volMesh> fieldType;
typedef GeometricField<Type, fvPatchField, volMesh> volFieldType;
if (foundObject<fieldType>(fieldName))
if (foundObject<volFieldType>(fieldName))
{
const Foam::dictionary& solverDict = mesh_.solverPerformanceDict();
@ -83,7 +115,7 @@ void Foam::functionObjects::residuals::writeResidual(const word& fieldName)
mesh_.validComponents<Type>()
);
for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; ++cmpt)
{
if (component(validComponents, cmpt) != -1)
{
@ -94,6 +126,8 @@ void Foam::functionObjects::residuals::writeResidual(const word& fieldName)
const word resultName =
fieldName + word(pTraits<Type>::componentNames[cmpt]);
setResult(resultName, r);
writeField(resultName);
}
}
}