Files
OpenFOAM-12/src/finiteVolume/cfdTools/general/solutionControl/convergenceControl/singleRegionConvergenceControl/singleRegionConvergenceControl.H
Will Bainbridge 4c8122783a solutionControl: Multi-region and PIMPLE time-loop control
The solution controls have been rewritten for use in multi-region
solvers, and PIMPLE fluid/solid solution controls have been implemented
within this framework.

PIMPLE also now has time-loop convergence control which can be used to
end the simulation once a certain initial residual is reached. This
allows a PIMPLE solver to run with equivalent convergence control to a
SIMPLE solver. Corrector loop convergence control is still available,
and can be used at the same time as the time-loop control.

The "residualControl" sub-dictionary of PIMPLE contains the residual
values required on the first solve of a time-step for the simulation to
end. This behaviour is the same as SIMPLE. The
"outerCorrectorResidualControl" sub-dictionary contains the tolerances
required for the corrector loop to exit. An example specification with
both types of control active is shown below.

PIMPLE
{
    // ...

    residualControl
    {
        p               1e-3;
        U               1e-4;
        "(k|epsilon|omega)" 1e-3;
    }

    outerCorrectorResidualControl
    {
        U
        {
            tolerance       1e-4;
            relTol          0.1;
        }
        "(k|epsilon|omega)"
        {
            tolerance       1e-3;
            relTol          0.1;
        }
    }
}

Note that existing PIMPLE "residualControl" entries will need to be
renamed "outerCorrectorResidualControl".

Application within a solver has also changed slightly. In order to have
convergence control for the time loop as a whole, the
solutionControl::loop(Time&) method (or the equivalent run method) must
be used; i.e.,

    while (simple.loop(runTime))
    {
        Info<< "Time = " << runTime.timeName() << nl << endl;

        // solve ...
    }

or,

    while (pimple.run(runTime))
    {
        // pre-time-increment operations ...

        runTime ++;
        Info<< "Time = " << runTime.timeName() << nl << endl;

        // solve ...
    }
2018-02-01 16:44:07 +00:00

116 lines
3.1 KiB
C++

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2018 OpenFOAM Foundation
\\/ 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::singleRegionConvergenceControl
Description
Single-region-specific derivation of the convergence control class
SourceFiles
singleRegionConvergenceControl.C
\*---------------------------------------------------------------------------*/
#ifndef singleRegionConvergenceControl_H
#define singleRegionConvergenceControl_H
#include "convergenceControl.H"
#include "singleRegionSolutionControl.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
/*---------------------------------------------------------------------------*\
Class singleRegionConvergenceControl Declaration
\*---------------------------------------------------------------------------*/
class singleRegionConvergenceControl
:
public convergenceControl
{
protected:
// Protected data
//- Reference to the mesh
const fvMesh& mesh_;
//- List of residual data per field
List<residualData> residualControl_;
public:
// Static Data Members
//- Run-time type information
TypeName("singleRegionConvergenceControl");
// Constructors
//- Construct from a solution control
singleRegionConvergenceControl
(
const singleRegionSolutionControl& control
);
//- Destructor
virtual ~singleRegionConvergenceControl();
// Member Functions
// IO
//- Read residual controls
bool readResidualControls();
//- Print the residual controls
void printResidualControls() const;
// Evolution
//- Return true if residual controls are present
virtual bool hasResidualControls() const;
//- Return true if all convergence checks are satisfied
virtual bool criteriaSatisfied() const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //