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 "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  // and maximum and minimum time-steps
57  rDeltaT.ref() = fvc::surfaceSum(mag(phi))()()/((2*maxCo)*mesh.V());
58  if (pimpleDict.found("maxDeltaT"))
59  {
60  rDeltaT.max(1/pimpleDict.lookup<scalar>("maxDeltaT"));
61  }
62  if (pimpleDict.found("minDeltaT"))
63  {
64  rDeltaT.min(1/pimpleDict.lookup<scalar>("minDeltaT"));
65  }
66 
67  // Update the boundary values of the reciprocal time-step
68  rDeltaT.correctBoundaryConditions();
69 
70  Info<< "Flow time scale min/max = "
71  << gMin(1/rDeltaT.primitiveField())
72  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
73 
74  if (rDeltaTSmoothingCoeff < 1.0)
75  {
76  fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
77  }
78 
79  Info<< "Smoothed flow time scale min/max = "
80  << gMin(1/rDeltaT.primitiveField())
81  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
82 
83  // Limit rate of change of time scale
84  // - reduce as much as required
85  // - only increase at a fraction of old time scale
86  if
87  (
88  rDeltaTDampingCoeff < 1.0
90  )
91  {
92  rDeltaT =
93  rDeltaT0
94  *max(rDeltaT/rDeltaT0, scalar(1) - rDeltaTDampingCoeff);
95 
96  Info<< "Damped flow time scale min/max = "
97  << gMin(1/rDeltaT.primitiveField())
98  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
99  }
100 }
101 
102 
103 // ************************************************************************* //
Generic GeometricField class.
void max(const dimensioned< Type > &)
Internal & ref()
Return a reference to the dimensioned internal field.
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal 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:773
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:860
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:659
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 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:251
messageStream Info
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)