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 maxDeltaT
17  (
18  pimpleDict.lookupOrDefault<scalar>("maxDeltaT", great)
19  );
20 
21  const scalar minDeltaT
22  (
23  pimpleDict.lookupOrDefault<scalar>("minDeltaT", small)
24  );
25 
26  const volScalarField rDeltaT0("rDeltaT0", rDeltaT);
27 
28  // Set the reciprocal time-step from the local Courant number
29  // and maximum and minimum time-steps
30  rDeltaT.ref() = min
31  (
33  max
34  (
36  fvc::surfaceSum(mag(phi))()()
37  /((2*maxCo)*mesh.V()*rho())
38  )
39  );
40 
41  if (pimple.transonic())
42  {
44  (
45  "phid",
47  );
48 
49  rDeltaT.ref() = max
50  (
51  rDeltaT(),
52  fvc::surfaceSum(mag(phid))()()
53  /((2*maxCo)*mesh.V()*psi())
54  );
55  }
56 
57  // Update the boundary values of the reciprocal time-step
58  rDeltaT.correctBoundaryConditions();
59 
60  Info<< "Flow time scale min/max = "
61  << gMin(1/rDeltaT.primitiveField())
62  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
63 
64  if (rDeltaTSmoothingCoeff < 1.0)
65  {
67  }
68 
69  Info<< "Smoothed flow time scale min/max = "
70  << gMin(1/rDeltaT.primitiveField())
71  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
72 
73  // Limit rate of change of time scale
74  // - reduce as much as required
75  // - only increase at a fraction of old time scale
76  if
77  (
78  pimpleDict.found("rDeltaTDampingCoeff")
79  && runTime.timeIndex() > runTime.startTimeIndex() + 1
80  )
81  {
82  // Damping coefficient (1-0)
83  const scalar rDeltaTDampingCoeff
84  (
85  pimpleDict.lookup<scalar>("rDeltaTDampingCoeff")
86  );
87 
88  rDeltaT =
89  rDeltaT0
90  *max(rDeltaT/rDeltaT0, scalar(1) - rDeltaTDampingCoeff);
91 
92  Info<< "Damped flow time scale min/max = "
93  << gMin(1/rDeltaT.primitiveField())
94  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
95  }
96 }
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const scalar minDeltaT(pimpleDict.lookupOrDefault< scalar >("minDeltaT", small))
pimpleNoLoopControl & pimple
U
Definition: pEqn.H:72
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
const volScalarField & psi
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)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
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< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< volScalarField > trDeltaT
Definition: createRDeltaT.H:3
const scalar maxDeltaT(pimpleDict.lookupOrDefault< scalar >("maxDeltaT", great))