ENH: fvSolution: allow Function1 for all scalars

TUT: demonstrate some ramping

- compressible/rhoPimpleFoam/RAS/angledDuct
This commit is contained in:
mattijs
2021-09-30 13:37:08 +01:00
committed by Mark Olesen
parent ba95b89e5f
commit 53e8f2a807
4 changed files with 91 additions and 33 deletions

View File

@ -27,6 +27,8 @@ License
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
#include "solution.H" #include "solution.H"
#include "HashPtrTable.H"
#include "Function1.H"
#include "Time.H" #include "Time.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
@ -55,23 +57,27 @@ void Foam::solution::read(const dictionary& dict)
if (dict.found("relaxationFactors")) if (dict.found("relaxationFactors"))
{ {
const dictionary& relaxDict(dict.subDict("relaxationFactors")); const dictionary& relaxDict = dict.subDict("relaxationFactors");
if (relaxDict.found("fields") || relaxDict.found("equations")) if (relaxDict.found("fields") || relaxDict.found("equations"))
{ {
if (relaxDict.found("fields")) if (relaxDict.found("fields"))
{ {
fieldRelaxDict_ = relaxDict.subDict("fields"); fieldRelaxDict_ = relaxDict.subDict("fields");
fieldRelaxCache_.clear();
} }
if (relaxDict.found("equations")) if (relaxDict.found("equations"))
{ {
eqnRelaxDict_ = relaxDict.subDict("equations"); eqnRelaxDict_ = relaxDict.subDict("equations");
eqnRelaxCache_.clear();
} }
} }
else else
{ {
// backwards compatibility // backwards compatibility
fieldRelaxDict_.clear(); fieldRelaxDict_.clear();
fieldRelaxCache_.clear();
for (const word& e : relaxDict.toc()) for (const word& e : relaxDict.toc())
{ {
@ -85,17 +91,38 @@ void Foam::solution::read(const dictionary& dict)
{ {
fieldRelaxDict_.add(e, value); fieldRelaxDict_.add(e, value);
} }
} }
eqnRelaxDict_ = relaxDict; eqnRelaxDict_ = relaxDict;
eqnRelaxCache_.clear();
} }
fieldRelaxDefault_ =
fieldRelaxDict_.getOrDefault<scalar>("default", 0);
eqnRelaxDefault_ = fieldRelaxDefault_ = Function1<scalar>::NewIfPresent
eqnRelaxDict_.getOrDefault<scalar>("default", 0); (
"default",
fieldRelaxDict_
);
if (!fieldRelaxDefault_)
{
fieldRelaxDefault_.reset
(
new Function1Types::Constant<scalar>("default", 0)
);
}
eqnRelaxDefault_ = Function1<scalar>::NewIfPresent
(
"default",
eqnRelaxDict_
);
if (!eqnRelaxDefault_)
{
eqnRelaxDefault_.reset
(
new Function1Types::Constant<scalar>("default", 0)
);
}
DebugInfo DebugInfo
<< "Relaxation factors:" << nl << "Relaxation factors:" << nl
@ -141,8 +168,6 @@ Foam::solution::solution
caching_(false), caching_(false),
fieldRelaxDict_(), fieldRelaxDict_(),
eqnRelaxDict_(), eqnRelaxDict_(),
fieldRelaxDefault_(0),
eqnRelaxDefault_(0),
solvers_() solvers_()
{ {
if if
@ -157,6 +182,13 @@ Foam::solution::solution
} }
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
// A non-default destructor since we had incomplete types in the header
Foam::solution::~solution()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
Foam::label Foam::solution::upgradeSolverDict Foam::label Foam::solution::upgradeSolverDict
@ -257,11 +289,17 @@ Foam::scalar Foam::solution::fieldRelaxationFactor(const word& name) const
if (fieldRelaxDict_.found(name)) if (fieldRelaxDict_.found(name))
{ {
return fieldRelaxDict_.get<scalar>(name); return Function1<scalar>::New
(
fieldRelaxCache_, // cache
name,
fieldRelaxDict_,
keyType::REGEX
)().value(time().timeOutputValue());
} }
else if (fieldRelaxDefault_ > SMALL) else if (fieldRelaxDefault_)
{ {
return fieldRelaxDefault_; return fieldRelaxDefault_().value(time().timeOutputValue());
} }
FatalIOErrorInFunction(fieldRelaxDict_) FatalIOErrorInFunction(fieldRelaxDict_)
@ -279,11 +317,17 @@ Foam::scalar Foam::solution::equationRelaxationFactor(const word& name) const
if (eqnRelaxDict_.found(name)) if (eqnRelaxDict_.found(name))
{ {
return eqnRelaxDict_.get<scalar>(name); return Function1<scalar>::New
(
eqnRelaxCache_, // cache
name,
eqnRelaxDict_,
keyType::REGEX
)().value(time().timeOutputValue());
} }
else if (eqnRelaxDefault_ > SMALL) else if (eqnRelaxDefault_)
{ {
return eqnRelaxDefault_; return eqnRelaxDefault_().value(time().timeOutputValue());
} }
FatalIOErrorInFunction(eqnRelaxDict_) FatalIOErrorInFunction(eqnRelaxDict_)

View File

@ -39,12 +39,16 @@ SourceFiles
#define solution_H #define solution_H
#include "IOdictionary.H" #include "IOdictionary.H"
#include "HashPtrTable.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam namespace Foam
{ {
// Forward Declarations
template<class Type> class Function1;
/*---------------------------------------------------------------------------*\ /*---------------------------------------------------------------------------*\
Class solution Declaration Class solution Declaration
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
@ -64,14 +68,20 @@ class solution
//- Dictionary of relaxation factors for all the fields //- Dictionary of relaxation factors for all the fields
dictionary fieldRelaxDict_; dictionary fieldRelaxDict_;
//- Cache of Function1s in above dictionary
mutable HashPtrTable<Function1<scalar>> fieldRelaxCache_;
//- Dictionary of relaxation factors for all the equations //- Dictionary of relaxation factors for all the equations
dictionary eqnRelaxDict_; dictionary eqnRelaxDict_;
//- Cache of Function1s in above dictionary
mutable HashPtrTable<Function1<scalar>> eqnRelaxCache_;
//- Optional default relaxation factor for all the fields //- Optional default relaxation factor for all the fields
scalar fieldRelaxDefault_; autoPtr<Function1<scalar>> fieldRelaxDefault_;
//- Optional default relaxation factor for all the equations //- Optional default relaxation factor for all the equations
scalar eqnRelaxDefault_; autoPtr<Function1<scalar>> eqnRelaxDefault_;
//- Dictionary of solver parameters for all the fields //- Dictionary of solver parameters for all the fields
dictionary solvers_; dictionary solvers_;
@ -125,6 +135,10 @@ public:
{} {}
//- Destructor
virtual ~solution();
// Member Functions // Member Functions
// Access // Access
@ -132,15 +146,6 @@ public:
//- Return true if the given field should be cached //- Return true if the given field should be cached
bool cache(const word& name) const; bool cache(const word& name) const;
//- Helper for printing cache message
template<class FieldType>
static void cachePrintMessage
(
const char* message,
const word& name,
const FieldType& vf
);
//- Return true if the relaxation factor is given for the field //- Return true if the relaxation factor is given for the field
bool relaxField(const word& name) const; bool relaxField(const word& name) const;
@ -168,6 +173,18 @@ public:
//- Read the solution dictionary //- Read the solution dictionary
bool read(); bool read();
// Other
//- Helper for printing cache message
template<class FieldType>
static void cachePrintMessage
(
const char* message,
const word& name,
const FieldType& vf
);
}; };

View File

@ -27,8 +27,7 @@ License
#include "solution.H" #include "solution.H"
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
template<class FieldType> template<class FieldType>
void Foam::solution::cachePrintMessage void Foam::solution::cachePrintMessage
@ -48,7 +47,4 @@ void Foam::solution::cachePrintMessage
} }
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// ************************************************************************* // // ************************************************************************* //

View File

@ -79,9 +79,10 @@ relaxationFactors
} }
equations equations
{ {
"U.*" 0.7; // "(U|e|k|epsilon).*" 0.7;
"e.*" 0.7;
"(k|epsilon).*" 0.7; // Demonstrate some ramping
"(U|e|k|epsilon).*" table ((0 0.4) (0.5 0.7));
} }
} }