mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-11-28 03:28:01 +00:00
pimpleControl: Added optional 'solveFlow' control
sprayFoam: Added support for the optional 'solveFlow' control to allow
simulation of the spray evolution with all sub-models in a 'frozen'
flow-field.
This commit is contained in:
@ -75,34 +75,37 @@ int main(int argc, char *argv[])
|
|||||||
|
|
||||||
parcels.evolve();
|
parcels.evolve();
|
||||||
|
|
||||||
#include "rhoEqn.H"
|
if (pimple.solveFlow())
|
||||||
|
|
||||||
// --- Pressure-velocity PIMPLE corrector loop
|
|
||||||
while (pimple.loop())
|
|
||||||
{
|
{
|
||||||
#include "UEqn.H"
|
#include "rhoEqn.H"
|
||||||
#include "YEqn.H"
|
|
||||||
#include "EEqn.H"
|
|
||||||
|
|
||||||
// --- Pressure corrector loop
|
// --- Pressure-velocity PIMPLE corrector loop
|
||||||
while (pimple.correct())
|
while (pimple.loop())
|
||||||
{
|
{
|
||||||
#include "pEqn.H"
|
#include "UEqn.H"
|
||||||
|
#include "YEqn.H"
|
||||||
|
#include "EEqn.H"
|
||||||
|
|
||||||
|
// --- Pressure corrector loop
|
||||||
|
while (pimple.correct())
|
||||||
|
{
|
||||||
|
#include "pEqn.H"
|
||||||
|
}
|
||||||
|
|
||||||
|
if (pimple.turbCorr())
|
||||||
|
{
|
||||||
|
turbulence->correct();
|
||||||
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
if (pimple.turbCorr())
|
rho = thermo.rho();
|
||||||
|
|
||||||
|
if (runTime.write())
|
||||||
{
|
{
|
||||||
turbulence->correct();
|
combustion->dQ()().write();
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
rho = thermo.rho();
|
|
||||||
|
|
||||||
if (runTime.write())
|
|
||||||
{
|
|
||||||
combustion->dQ()().write();
|
|
||||||
}
|
|
||||||
|
|
||||||
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
|
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
|
||||||
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
|
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
|
||||||
<< nl << endl;
|
<< nl << endl;
|
||||||
|
|||||||
@ -40,8 +40,9 @@ void Foam::pimpleControl::read()
|
|||||||
{
|
{
|
||||||
solutionControl::read(false);
|
solutionControl::read(false);
|
||||||
|
|
||||||
// Read solution controls
|
|
||||||
const dictionary& pimpleDict = dict();
|
const dictionary& pimpleDict = dict();
|
||||||
|
|
||||||
|
solveFlow_ = pimpleDict.lookupOrDefault<Switch>("solveFlow", true);
|
||||||
nCorrPIMPLE_ = pimpleDict.lookupOrDefault<label>("nOuterCorrectors", 1);
|
nCorrPIMPLE_ = pimpleDict.lookupOrDefault<label>("nOuterCorrectors", 1);
|
||||||
nCorrPISO_ = pimpleDict.lookupOrDefault<label>("nCorrectors", 1);
|
nCorrPISO_ = pimpleDict.lookupOrDefault<label>("nCorrectors", 1);
|
||||||
turbOnFinalIterOnly_ =
|
turbOnFinalIterOnly_ =
|
||||||
@ -123,6 +124,7 @@ bool Foam::pimpleControl::criteriaSatisfied()
|
|||||||
Foam::pimpleControl::pimpleControl(fvMesh& mesh, const word& dictName)
|
Foam::pimpleControl::pimpleControl(fvMesh& mesh, const word& dictName)
|
||||||
:
|
:
|
||||||
solutionControl(mesh, dictName),
|
solutionControl(mesh, dictName),
|
||||||
|
solveFlow_(true),
|
||||||
nCorrPIMPLE_(0),
|
nCorrPIMPLE_(0),
|
||||||
nCorrPISO_(0),
|
nCorrPISO_(0),
|
||||||
corrPISO_(0),
|
corrPISO_(0),
|
||||||
|
|||||||
@ -69,6 +69,9 @@ protected:
|
|||||||
|
|
||||||
// Solution controls
|
// Solution controls
|
||||||
|
|
||||||
|
//- Flag to indicate whether to solve for the flow
|
||||||
|
bool solveFlow_;
|
||||||
|
|
||||||
//- Maximum number of PIMPLE correctors
|
//- Maximum number of PIMPLE correctors
|
||||||
label nCorrPIMPLE_;
|
label nCorrPIMPLE_;
|
||||||
|
|
||||||
@ -131,22 +134,25 @@ public:
|
|||||||
//- PIMPLE loop
|
//- PIMPLE loop
|
||||||
virtual bool loop();
|
virtual bool loop();
|
||||||
|
|
||||||
//- Pressure corrector loop
|
//- Pressure corrector loop control
|
||||||
inline bool correct();
|
inline bool correct();
|
||||||
|
|
||||||
//- Helper function to identify when to store the intial residuals
|
//- Return true to store the intial residuals
|
||||||
inline bool storeInitialResiduals() const;
|
inline bool storeInitialResiduals() const;
|
||||||
|
|
||||||
//- Helper function to identify first PIMPLE (outer) iteration
|
//- Return true for first PIMPLE (outer) iteration
|
||||||
inline bool firstIter() const;
|
inline bool firstIter() const;
|
||||||
|
|
||||||
//- Helper function to identify final PIMPLE (outer) iteration
|
//- Return true fore final PIMPLE (outer) iteration
|
||||||
inline bool finalIter() const;
|
inline bool finalIter() const;
|
||||||
|
|
||||||
//- Helper function to identify final inner iteration
|
//- Return true for final inner iteration
|
||||||
inline bool finalInnerIter() const;
|
inline bool finalInnerIter() const;
|
||||||
|
|
||||||
//- Helper function to identify whether to solve for turbulence
|
//- Return true to solve for flow
|
||||||
|
inline bool solveFlow() const;
|
||||||
|
|
||||||
|
//- Return true to solve for turbulence
|
||||||
inline bool turbCorr() const;
|
inline bool turbCorr() const;
|
||||||
};
|
};
|
||||||
|
|
||||||
|
|||||||
@ -2,7 +2,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) 2011-2013 OpenFOAM Foundation
|
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
|
||||||
\\/ M anipulation |
|
\\/ M anipulation |
|
||||||
-------------------------------------------------------------------------------
|
-------------------------------------------------------------------------------
|
||||||
License
|
License
|
||||||
@ -92,6 +92,12 @@ inline bool Foam::pimpleControl::finalInnerIter() const
|
|||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
inline bool Foam::pimpleControl::solveFlow() const
|
||||||
|
{
|
||||||
|
return solveFlow_;
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
inline bool Foam::pimpleControl::turbCorr() const
|
inline bool Foam::pimpleControl::turbCorr() const
|
||||||
{
|
{
|
||||||
return !turbOnFinalIterOnly_ || finalIter();
|
return !turbOnFinalIterOnly_ || finalIter();
|
||||||
|
|||||||
Reference in New Issue
Block a user