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-2024 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 "incompressibleFluid.H"
27 #include "fvcSmooth.H"
28 #include "fvcSurfaceIntegrate.H"
29 
30 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
31 
33 {
34  volScalarField& rDeltaT = trDeltaT.ref();
35 
36  const dictionary& pimpleDict = pimple.dict();
37 
38  const scalar maxCo
39  (
40  pimpleDict.lookupOrDefault<scalar>("maxCo", 0.8)
41  );
42 
43  const scalar rDeltaTSmoothingCoeff
44  (
45  pimpleDict.lookupOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.02)
46  );
47 
48  const scalar rDeltaTDampingCoeff
49  (
50  pimpleDict.lookupOrDefault<scalar>("rDeltaTDampingCoeff", 1.0)
51  );
52 
53  const volScalarField rDeltaT0("rDeltaT0", rDeltaT);
54 
55  // Set the reciprocal time-step from the local Courant number
56  rDeltaT.internalFieldRef() =
57  fvc::surfaceSum(mag(phi))()()/((2*maxCo)*mesh.V());
58 
59  // Clip to user-defined maximum and minimum time-steps
60  scalar minRDeltaT = gMin(rDeltaT.primitiveField());
61  scalar maxRDeltaT = gMax(rDeltaT.primitiveField());
62  if (pimpleDict.found("maxDeltaT") || minRDeltaT < rootVSmall)
63  {
64  const scalar clipRDeltaT = 1/pimpleDict.lookup<scalar>("maxDeltaT");
65  rDeltaT.max(clipRDeltaT);
66  minRDeltaT = max(minRDeltaT, clipRDeltaT);
67  maxRDeltaT = max(maxRDeltaT, clipRDeltaT);
68  }
69  if (pimpleDict.found("minDeltaT") || maxRDeltaT > rootVGreat)
70  {
71  const scalar clipRDeltaT = 1/pimpleDict.lookup<scalar>("minDeltaT");
72  rDeltaT.min(clipRDeltaT);
73  minRDeltaT = min(minRDeltaT, clipRDeltaT);
74  maxRDeltaT = min(maxRDeltaT, clipRDeltaT);
75  }
76 
77  Info<< "Flow time scale min/max = "
78  << 1/maxRDeltaT << ", " << 1/minRDeltaT << endl;
79 
80  // Update the boundary values of the reciprocal time-step
81  rDeltaT.correctBoundaryConditions();
82 
83  if (rDeltaTSmoothingCoeff < 1.0)
84  {
85  fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
86  }
87 
88  Info<< "Smoothed flow time scale min/max = "
89  << gMin(1/rDeltaT.primitiveField())
90  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
91 
92  // Limit rate of change of time scale
93  // - reduce as much as required
94  // - only increase at a fraction of old time scale
95  if
96  (
97  rDeltaTDampingCoeff < 1.0
99  )
100  {
101  rDeltaT =
102  rDeltaT0
103  *max(rDeltaT/rDeltaT0, scalar(1) - rDeltaTDampingCoeff);
104 
105  Info<< "Damped flow time scale min/max = "
106  << gMin(1/rDeltaT.primitiveField())
107  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
108  }
109 }
110 
111 
112 // ************************************************************************* //
Generic GeometricField class.
void max(const dimensioned< Type > &)
Internal & internalFieldRef()
Return a reference to the dimensioned internal field.
const Internal::FieldType & primitiveField() const
Return a const-reference to the primitive field.
void min(const dimensioned< Type > &)
void correctBoundaryConditions()
Correct boundary field.
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
virtual label startTimeIndex() const
Return start time index.
Definition: Time.C:785
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
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:107
const Time & runTime
Time.
Definition: solver.H:104
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
const surfaceScalarField & phi
Reference to the volumetric-flux field.
tmp< volScalarField > trDeltaT
Optional LTS reciprocal time-step field.
void setRDeltaT()
Set rDeltaT for LTS.
Definition: setRDeltaT.C:32
scalar maxCo
Provides functions smooth spread and sweep which use the FvFaceCellWave algorithm to smooth and redis...
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
void smooth(volScalarField &field, const scalar coeff)
Definition: fvcSmooth.C:35
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
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)