setRDeltaT.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2022-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "VoFSolver.H"
27 #include "fvcSurfaceIntegrate.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 void Foam::solvers::VoFSolver::setRDeltaT()
32 {
33  volScalarField& rDeltaT = trDeltaT.ref();
34 
35  const dictionary& pimpleDict = pimple.dict();
36 
37  const scalar maxCo
38  (
39  pimpleDict.lookupOrDefault<scalar>("maxCo", 0.9)
40  );
41 
42  const volScalarField rDeltaT0("rDeltaT0", rDeltaT);
43 
44  // Set the reciprocal time-step from the local Courant number
45  // and maximum and minimum time-steps
46  rDeltaT.ref() = fvc::surfaceSum(mag(phi))()()/((2*maxCo)*mesh.V());
47  if (pimpleDict.found("maxDeltaT"))
48  {
49  rDeltaT.max(1/pimpleDict.lookup<scalar>("maxDeltaT"));
50  }
51  if (pimpleDict.found("minDeltaT"))
52  {
53  rDeltaT.min(1/pimpleDict.lookup<scalar>("minDeltaT"));
54  }
55 
56  Info<< "Flow time scale min/max = "
57  << gMin(1/rDeltaT.primitiveField())
58  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
59 
60  setInterfaceRDeltaT(rDeltaT);
61 
62  // Limit rate of change of time scale
63  // - reduce as much as required
64  // - only increase at a fraction of old time scale
65  if
66  (
67  pimpleDict.found("rDeltaTDampingCoeff")
69  )
70  {
71  // Damping coefficient (1-0)
72  const scalar rDeltaTDampingCoeff
73  (
74  pimpleDict.lookup<scalar>("rDeltaTDampingCoeff")
75  );
76 
77  rDeltaT = max
78  (
79  rDeltaT,
80  (scalar(1) - rDeltaTDampingCoeff)*rDeltaT0
81  );
82 
83  Info<< "Damped flow time scale min/max = "
84  << gMin(1/rDeltaT.primitiveField())
85  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
86  }
87 }
88 
89 
90 // ************************************************************************* //
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
virtual label startTimeIndex() const
Return start time index.
Definition: Time.C:773
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
virtual const dictionary & dict() const
Return the solution dictionary.
pimpleNoLoopControl pimple
PIMPLE inner-loop controls.
Definition: solver.H:100
const Time & runTime
Time.
Definition: solver.H:97
const fvMesh & mesh
Region mesh.
Definition: solver.H:94
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: VoFSolver.H:212
tmp< volScalarField > trDeltaT
Optional LTS reciprocal time-step field.
Definition: VoFSolver.H:141
virtual void setInterfaceRDeltaT(volScalarField &rDeltaT)=0
Adjust the rDeltaT in the vicinity of the interface.
scalar maxCo
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)