setRDeltaT.H
Go to the documentation of this file.
1 {
2  volScalarField& rDeltaT = trDeltaT.ref();
3 
4  const dictionary& pimpleDict = pimple.dict();
5 
6  scalar maxCo
7  (
8  pimpleDict.lookupOrDefault<scalar>("maxCo", 0.8)
9  );
10 
12  (
13  pimpleDict.lookupOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.02)
14  );
15 
16  scalar rDeltaTDampingCoeff
17  (
18  pimpleDict.lookupOrDefault<scalar>("rDeltaTDampingCoeff", 1.0)
19  );
20 
21  scalar maxDeltaT
22  (
23  pimpleDict.lookupOrDefault<scalar>("maxDeltaT", great)
24  );
25 
26  volScalarField rDeltaT0("rDeltaT0", rDeltaT);
27 
28  // Set the reciprocal time-step from the local Courant number
29  rDeltaT.ref() = max
30  (
32  fvc::surfaceSum(mag(phi))()()
33  /((2*maxCo)*mesh.V())
34  );
35 
36  // Update the boundary values of the reciprocal time-step
37  rDeltaT.correctBoundaryConditions();
38 
39  Info<< "Flow time scale min/max = "
40  << gMin(1/rDeltaT.primitiveField())
41  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
42 
43  if (rDeltaTSmoothingCoeff < 1.0)
44  {
46  }
47 
48  Info<< "Smoothed flow time scale min/max = "
49  << gMin(1/rDeltaT.primitiveField())
50  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
51 
52  // Limit rate of change of time scale
53  // - reduce as much as required
54  // - only increase at a fraction of old time scale
55  if
56  (
58  && runTime.timeIndex() > runTime.startTimeIndex() + 1
59  )
60  {
61  rDeltaT =
62  rDeltaT0
63  *max(rDeltaT/rDeltaT0, scalar(1) - rDeltaTDampingCoeff);
64 
65  Info<< "Damped flow time scale min/max = "
66  << gMin(1/rDeltaT.primitiveField())
67  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
68  }
69 }
pimpleNoLoopControl & pimple
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Type gMin(const FieldField< Field, Type > &f)
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
scalar rDeltaTDampingCoeff(pimpleDict.lookupOrDefault< scalar >("rDeltaTDampingCoeff", 1.0))
void smooth(volScalarField &field, const scalar coeff)
Definition: fvcSmooth.C:35
phi
Definition: pEqn.H:104
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
scalar maxCo(pimpleDict.lookupOrDefault< scalar >("maxCo", 0.9))
scalar rDeltaTSmoothingCoeff(pimpleDict.lookupOrDefault< scalar >("rDeltaTSmoothingCoeff", 0.1))
Type gMax(const FieldField< Field, Type > &f)
volScalarField rDeltaT0("rDeltaT0", rDeltaT)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
scalar maxDeltaT(pimpleDict.lookupOrDefault< scalar >("maxDeltaT", great))
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< volScalarField > trDeltaT
Definition: createRDeltaT.H:3