All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
setRDeltaT.H
Go to the documentation of this file.
1 {
2  volScalarField& rDeltaT = trDeltaT.ref();
3 
4  const dictionary& pimpleDict = pimple.dict();
5 
6  const scalar maxCo
7  (
8  pimpleDict.lookupOrDefault<scalar>("maxCo", 0.8)
9  );
10 
11  const scalar rDeltaTSmoothingCoeff
12  (
13  pimpleDict.lookupOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.02)
14  );
15 
16  const scalar rDeltaTDampingCoeff
17  (
18  pimpleDict.lookupOrDefault<scalar>("rDeltaTDampingCoeff", 1.0)
19  );
20 
21  const scalar maxDeltaT
22  (
23  pimpleDict.lookupOrDefault<scalar>("maxDeltaT", great)
24  );
25 
26  const scalar minDeltaT
27  (
28  pimpleDict.lookupOrDefault<scalar>("minDeltaT", small)
29  );
30 
31  const volScalarField rDeltaT0("rDeltaT0", rDeltaT);
32 
33  // Set the reciprocal time-step from the local Courant number
34  // and maximum and minimum time-steps
35  rDeltaT.ref() = min
36  (
38  max
39  (
41  fvc::surfaceSum(mag(phi))()()
42  /((2*maxCo)*mesh.V())
43  )
44  );
45 
46  // Update the boundary values of the reciprocal time-step
47  rDeltaT.correctBoundaryConditions();
48 
49  Info<< "Flow time scale min/max = "
50  << gMin(1/rDeltaT.primitiveField())
51  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
52 
53  if (rDeltaTSmoothingCoeff < 1.0)
54  {
56  }
57 
58  Info<< "Smoothed flow time scale min/max = "
59  << gMin(1/rDeltaT.primitiveField())
60  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
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  (
68  && runTime.timeIndex() > runTime.startTimeIndex() + 1
69  )
70  {
71  rDeltaT =
72  rDeltaT0
73  *max(rDeltaT/rDeltaT0, scalar(1) - rDeltaTDampingCoeff);
74 
75  Info<< "Damped flow time scale min/max = "
76  << gMin(1/rDeltaT.primitiveField())
77  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
78  }
79 }
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const scalar minDeltaT(pimpleDict.lookupOrDefault< scalar >("minDeltaT", small))
pimpleNoLoopControl & pimple
Type gMin(const FieldField< Field, Type > &f)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fvMesh & mesh
void smooth(volScalarField &field, const scalar coeff)
Definition: fvcSmooth.C:35
const volScalarField rDeltaT0("rDeltaT0", rDeltaT)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionSet dimTime
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
phi
Definition: correctPhi.H:3
const scalar maxCo(pimpleDict.lookupOrDefault< scalar >("maxCo", 0.9))
Type gMax(const FieldField< Field, Type > &f)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const scalar rDeltaTSmoothingCoeff(pimpleDict.lookupOrDefault< scalar >("rDeltaTSmoothingCoeff", 0.1))
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
const scalar rDeltaTDampingCoeff(pimpleDict.lookupOrDefault< scalar >("rDeltaTDampingCoeff", 1.0))
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< volScalarField > trDeltaT
Definition: createRDeltaT.H:3
const scalar maxDeltaT(pimpleDict.lookupOrDefault< scalar >("maxDeltaT", great))