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

@ -26,6 +26,9 @@ License
#include "lduMatrix.H" #include "lduMatrix.H"
#include "IOstreams.H" #include "IOstreams.H"
#include "Switch.H" #include "Switch.H"
#include "objectRegistry.H"
#include "IOField.H"
#include "Time.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
@ -312,6 +315,38 @@ const Foam::scalarField& Foam::lduMatrix::upper() const
} }
void Foam::lduMatrix::setResidualField
(
const Field<scalar>& residual,
const word& fieldName,
const bool initial
) const
{
if (!lduMesh_.hasDb())
{
return;
}
word lookupName;
if (initial)
{
lookupName = word("initialResidual:" + fieldName);
}
else
{
lookupName = word("residual:" + fieldName);
}
IOField<scalar>* residualPtr =
lduMesh_.thisDb().lookupObjectRefPtr<IOField<scalar>>(lookupName);
if (residualPtr)
{
*residualPtr = residual;
}
}
// * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * Friend Operators * * * * * * * * * * * * * //
Foam::Ostream& Foam::operator<<(Ostream& os, const lduMatrix& ldum) Foam::Ostream& Foam::operator<<(Ostream& os, const lduMatrix& ldum)

View File

@ -124,6 +124,7 @@ public:
profilingTrigger profiling_; profilingTrigger profiling_;
// Protected Member Functions // Protected Member Functions
//- Read the control parameters from the controlDict_ //- Read the control parameters from the controlDict_
@ -258,7 +259,7 @@ public:
) const = 0; ) const = 0;
//- Return the matrix norm used to normalise the residual for the //- Return the matrix norm used to normalise the residual for the
// stopping criterion //- stopping criterion
scalar normFactor scalar normFactor
( (
const scalarField& psi, const scalarField& psi,
@ -485,7 +486,7 @@ public:
// Member functions // Member functions
//- Read and reset the preconditioner parameters //- Read and reset the preconditioner parameters
// from the given stream //- from the given stream
virtual void read(const dictionary&) virtual void read(const dictionary&)
{} {}
@ -498,7 +499,7 @@ public:
) const = 0; ) const = 0;
//- Return wT the transpose-matrix preconditioned form of //- Return wT the transpose-matrix preconditioned form of
// residual rT. //- residual rT.
// This is only required for preconditioning asymmetric matrices. // This is only required for preconditioning asymmetric matrices.
virtual void preconditionT virtual void preconditionT
( (
@ -531,7 +532,7 @@ public:
lduMatrix(lduMatrix&, bool reuse); lduMatrix(lduMatrix&, bool reuse);
//- Construct given an LDU addressed mesh and an Istream //- Construct given an LDU addressed mesh and an Istream
// from which the coefficients are read //- from which the coefficients are read
lduMatrix(const lduMesh&, Istream&); lduMatrix(const lduMesh&, Istream&);
@ -669,7 +670,7 @@ public:
//- Initialise the update of interfaced interfaces //- Initialise the update of interfaced interfaces
// for matrix operations //- for matrix operations
void initMatrixInterfaces void initMatrixInterfaces
( (
const bool add, const bool add,
@ -691,6 +692,14 @@ public:
const direction cmpt const direction cmpt
) const; ) const;
//- Set the residual field using an IOField on the object registry
//- if it exists
void setResidualField
(
const Field<scalar>& residual,
const word& fieldName,
const bool initial
) const;
template<class Type> template<class Type>
tmp<Field<Type>> H(const Field<Type>&) const; tmp<Field<Type>> H(const Field<Type>&) const;

View File

@ -59,6 +59,8 @@ Foam::solverPerformance Foam::GAMGSolver::solve
// Calculate initial finest-grid residual field // Calculate initial finest-grid residual field
scalarField finestResidual(source - Apsi); scalarField finestResidual(source - Apsi);
matrix().setResidualField(finestResidual, fieldName_, true);
// Calculate normalised residual for convergence test // Calculate normalised residual for convergence test
solverPerf.initialResidual() = gSumMag solverPerf.initialResidual() = gSumMag
( (
@ -143,6 +145,8 @@ Foam::solverPerformance Foam::GAMGSolver::solve
); );
} }
matrix().setResidualField(finestResidual, fieldName_, false);
return solverPerf; return solverPerf;
} }

View File

@ -93,6 +93,8 @@ Foam::solverPerformance Foam::PBiCG::solve
scalarField rA(source - wA); scalarField rA(source - wA);
scalar* __restrict__ rAPtr = rA.begin(); scalar* __restrict__ rAPtr = rA.begin();
matrix().setResidualField(rA, fieldName_, true);
// --- Calculate normalisation factor // --- Calculate normalisation factor
const scalar normFactor = this->normFactor(psi, source, wA, pA); const scalar normFactor = this->normFactor(psi, source, wA, pA);
@ -218,6 +220,8 @@ Foam::solverPerformance Foam::PBiCG::solve
<< exit(FatalError); << exit(FatalError);
} }
matrix().setResidualField(rA, fieldName_, false);
return solverPerf; return solverPerf;
} }

View File

@ -96,6 +96,8 @@ Foam::solverPerformance Foam::PBiCGStab::solve
scalarField rA(source - yA); scalarField rA(source - yA);
scalar* __restrict__ rAPtr = rA.begin(); scalar* __restrict__ rAPtr = rA.begin();
matrix().setResidualField(rA, fieldName_, true);
// --- Calculate normalisation factor // --- Calculate normalisation factor
const scalar normFactor = this->normFactor(psi, source, yA, pA); const scalar normFactor = this->normFactor(psi, source, yA, pA);
@ -248,6 +250,8 @@ Foam::solverPerformance Foam::PBiCGStab::solve
); );
} }
matrix().setResidualField(rA, fieldName_, false);
return solverPerf; return solverPerf;
} }

View File

@ -96,6 +96,8 @@ Foam::solverPerformance Foam::PCG::solve
scalarField rA(source - wA); scalarField rA(source - wA);
scalar* __restrict__ rAPtr = rA.begin(); scalar* __restrict__ rAPtr = rA.begin();
matrix().setResidualField(rA, fieldName_, true);
// --- Calculate normalisation factor // --- Calculate normalisation factor
scalar normFactor = this->normFactor(psi, source, wA, pA); scalar normFactor = this->normFactor(psi, source, wA, pA);
@ -189,6 +191,8 @@ Foam::solverPerformance Foam::PCG::solve
); );
} }
matrix().setResidualField(rA, fieldName_, false);
return solverPerf; return solverPerf;
} }

View File

@ -113,6 +113,7 @@ Foam::solverPerformance Foam::smoothSolver::solve
else else
{ {
scalar normFactor = 0; scalar normFactor = 0;
scalarField residual;
{ {
scalarField Apsi(psi.size()); scalarField Apsi(psi.size());
@ -124,12 +125,13 @@ Foam::solverPerformance Foam::smoothSolver::solve
// Calculate normalisation factor // Calculate normalisation factor
normFactor = this->normFactor(psi, source, Apsi, temp); normFactor = this->normFactor(psi, source, Apsi, temp);
residual = source - Apsi;
matrix().setResidualField(residual, fieldName_, true);
// Calculate residual magnitude // Calculate residual magnitude
solverPerf.initialResidual() = gSumMag solverPerf.initialResidual() =
( gSumMag(residual, matrix().mesh().comm())/normFactor;
(source - Apsi)(),
matrix().mesh().comm()
)/normFactor;
solverPerf.finalResidual() = solverPerf.initialResidual(); solverPerf.finalResidual() = solverPerf.initialResidual();
} }
@ -170,9 +172,7 @@ Foam::solverPerformance Foam::smoothSolver::solve
nSweeps_ nSweeps_
); );
// Calculate the residual to check convergence residual =
solverPerf.finalResidual() = gSumMag
(
matrix_.residual matrix_.residual
( (
psi, psi,
@ -180,9 +180,11 @@ Foam::solverPerformance Foam::smoothSolver::solve
interfaceBouCoeffs_, interfaceBouCoeffs_,
interfaces_, interfaces_,
cmpt cmpt
)(), );
matrix().mesh().comm()
)/normFactor; // Calculate the residual to check convergence
solverPerf.finalResidual() =
gSumMag(residual, matrix().mesh().comm())/normFactor;
} while } while
( (
( (
@ -192,6 +194,8 @@ Foam::solverPerformance Foam::smoothSolver::solve
|| solverPerf.nIterations() < minIter_ || solverPerf.nIterations() < minIter_
); );
} }
matrix().setResidualField(residual, fieldName_, false);
} }
return solverPerf; return solverPerf;

View File

@ -74,6 +74,12 @@ public:
// Member Functions // Member Functions
//- Return true if thisDb() is a valid DB - here = false
bool hasDb() const
{
return true;
}
//- Return the object registry //- Return the object registry
const objectRegistry& thisDb() const const objectRegistry& thisDb() const
{ {

View File

@ -77,6 +77,9 @@ public:
// Access // Access
//- Return true if thisDb() is a valid DB
virtual bool hasDb() const = 0;
//- Return the object registry //- Return the object registry
virtual const objectRegistry& thisDb() const; virtual const objectRegistry& thisDb() const;

View File

@ -73,6 +73,7 @@ class lduPrimitiveMesh
//- Communicator to use for any parallel communication //- Communicator to use for any parallel communication
const label comm_; const label comm_;
// Private Member Functions // Private Member Functions
//- Get size of all meshes //- Get size of all meshes
@ -92,6 +93,7 @@ public:
// Declare name of the class and its debug switch // Declare name of the class and its debug switch
ClassName("lduPrimitiveMesh"); ClassName("lduPrimitiveMesh");
// Constructors // Constructors
//- Construct from components but without interfaces. Add interfaces //- Construct from components but without interfaces. Add interfaces
@ -170,6 +172,12 @@ public:
// Access // Access
//- Return true if thisDb() is a valid DB
virtual bool hasDb() const
{
return false;
}
//- Return ldu addressing //- Return ldu addressing
virtual const lduAddressing& lduAddr() const virtual const lduAddressing& lduAddr() const
{ {

View File

@ -239,6 +239,12 @@ public:
return polyMesh::time(); return polyMesh::time();
} }
//- Return true if thisDb() is a valid DB
virtual bool hasDb() const
{
return true;
}
//- Return the object registry - resolve conflict polyMesh/lduMesh //- Return the object registry - resolve conflict polyMesh/lduMesh
virtual const objectRegistry& thisDb() const virtual const objectRegistry& thisDb() const
{ {

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 * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::functionObjects::residuals::residuals Foam::functionObjects::residuals::residuals
@ -90,7 +153,9 @@ Foam::functionObjects::residuals::residuals
: :
fvMeshFunctionObject(name, runTime, dict), fvMeshFunctionObject(name, runTime, dict),
writeFile(obr_, name, typeName, dict), writeFile(obr_, name, typeName, dict),
fieldSet_(mesh_) fieldSet_(mesh_),
writeFields_(false),
initialised_(false)
{ {
read(dict); read(dict);
} }
@ -109,6 +174,9 @@ bool Foam::functionObjects::residuals::read(const dictionary& dict)
if (fvMeshFunctionObject::read(dict)) if (fvMeshFunctionObject::read(dict))
{ {
fieldSet_.read(dict); fieldSet_.read(dict);
writeFields_ = dict.lookupOrDefault("writeFields", false);
return true; return true;
} }
@ -118,16 +186,33 @@ bool Foam::functionObjects::residuals::read(const dictionary& dict)
bool Foam::functionObjects::residuals::execute() bool Foam::functionObjects::residuals::execute()
{ {
// Note: delaying initialisation until after first iteration so that
// we can find wildcard fields
if (!initialised_)
{
writeFileHeader(file());
if (writeFields_)
{
for (const word& fieldName : fieldSet_.selection())
{
initialiseField<scalar>(fieldName);
initialiseField<vector>(fieldName);
initialiseField<sphericalTensor>(fieldName);
initialiseField<symmTensor>(fieldName);
initialiseField<tensor>(fieldName);
}
}
initialised_ = true;
}
return true; return true;
} }
bool Foam::functionObjects::residuals::write() bool Foam::functionObjects::residuals::write()
{ {
if (Pstream::master())
{
writeFileHeader(file());
writeTime(file()); writeTime(file());
for (const word& fieldName : fieldSet_.selection()) for (const word& fieldName : fieldSet_.selection())
@ -140,7 +225,6 @@ bool Foam::functionObjects::residuals::write()
} }
file() << endl; file() << endl;
}
return true; return true;
} }

View File

@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | \\ / O peration |
\\ / A nd | Copyright (C) 2015-2016 OpenFOAM Foundation \\ / 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 License
This file is part of OpenFOAM. This file is part of OpenFOAM.
@ -87,17 +87,33 @@ protected:
//- Fields to write residuals //- Fields to write residuals
solverFieldSelection fieldSet_; solverFieldSelection fieldSet_;
//- Flag to write the residual as a vol field
bool writeFields_;
//- Initialisation flag
bool initialised_;
// Protected Member Functions // Protected Member Functions
//- Output file header information //- Output file header information
void writeFileHeader(Ostream& os); 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 //- Output file header information per primitive type value
template<class Type> template<class Type>
void writeFileHeader(Ostream& os, const word& fileName) const; 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> template<class Type>
void writeResidual(const word& fieldName); void writeResidual(const word& fieldName);

View File

@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | \\ / O peration |
\\ / A nd | Copyright (C) 2015-2016 OpenFOAM Foundation \\ / 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 License
This file is part of OpenFOAM. This file is part of OpenFOAM.
@ -26,6 +26,7 @@ License
#include "residuals.H" #include "residuals.H"
#include "volFields.H" #include "volFields.H"
#include "ListOps.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> template<class Type>
void Foam::functionObjects::residuals::writeResidual(const word& fieldName) 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(); const Foam::dictionary& solverDict = mesh_.solverPerformanceDict();
@ -83,7 +115,7 @@ void Foam::functionObjects::residuals::writeResidual(const word& fieldName)
mesh_.validComponents<Type>() 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) if (component(validComponents, cmpt) != -1)
{ {
@ -94,6 +126,8 @@ void Foam::functionObjects::residuals::writeResidual(const word& fieldName)
const word resultName = const word resultName =
fieldName + word(pTraits<Type>::componentNames[cmpt]); fieldName + word(pTraits<Type>::componentNames[cmpt]);
setResult(resultName, r); setResult(resultName, r);
writeField(resultName);
} }
} }
} }