Files
OpenFOAM-12/applications/modules/multiphaseEuler/setRDeltaT.C
Will Bainbridge 597121a4a7 multiphaseEuler: Library reorganisation
This change makes multiphaseEuler more consistent with other modules and
makes its sub-libraries less inter-dependent. Some left-over references
to multiphaseEulerFoam have also been removed.
2023-09-15 14:45:26 +01:00

93 lines
3.1 KiB
C++

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2022-2023 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/>.
\*---------------------------------------------------------------------------*/
#include "multiphaseEuler.H"
#include "fvcSmooth.H"
#include "fvcSurfaceIntegrate.H"
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void Foam::solvers::multiphaseEuler::setRDeltaT()
{
volScalarField& rDeltaT = trDeltaT.ref();
const dictionary& pimpleDict = pimple.dict();
const scalar maxCo
(
pimpleDict.lookupOrDefault<scalar>("maxCo", 0.2)
);
const scalar rDeltaTSmoothingCoeff
(
pimpleDict.lookupOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.02)
);
surfaceScalarField maxPhi("maxPhi", phi);
forAll(movingPhases, movingPhasei)
{
maxPhi = max(maxPhi, mag(movingPhases[movingPhasei].phi()));
}
// Set the reciprocal time-step from the local Courant number
rDeltaT.ref() = fvc::surfaceSum(maxPhi)()()/((2*maxCo)*mesh.V());
// Clip to user-defined maximum and minimum time-steps
scalar minRDeltaT = gMin(rDeltaT.primitiveField());
if (pimpleDict.found("maxDeltaT") || minRDeltaT < rootVSmall)
{
const scalar clipRDeltaT = 1/pimpleDict.lookup<scalar>("maxDeltaT");
rDeltaT.max(clipRDeltaT);
minRDeltaT = max(minRDeltaT, clipRDeltaT);
}
if (pimpleDict.found("minDeltaT"))
{
const scalar clipRDeltaT = 1/pimpleDict.lookup<scalar>("minDeltaT");
rDeltaT.min(clipRDeltaT);
minRDeltaT = min(minRDeltaT, clipRDeltaT);
}
Info<< "Flow time scale min/max = "
<< gMin(1/rDeltaT.primitiveField()) << ", " << 1/minRDeltaT << endl;
// Update the boundary values of the reciprocal time-step
rDeltaT.correctBoundaryConditions();
fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
Info<< "Smoothed flow time scale min/max = "
<< gMin(1/rDeltaT.primitiveField())
<< ", " << gMax(1/rDeltaT.primitiveField()) << endl;
if (faceMomentum)
{
trDeltaTf.ref() = fvc::interpolate(rDeltaT);
}
}
// ************************************************************************* //