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.getOrDefault<scalar>("maxCo", 0.2)
9  );
10 
11  scalar maxDeltaT
12  (
13  pimpleDict.getOrDefault<scalar>("maxDeltaT", GREAT)
14  );
15 
17  (
18  pimpleDict.getOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.02)
19  );
20 
21  surfaceScalarField maxPhi("maxPhi", phi);
22 
24  {
26  }
27 
28  // Set the reciprocal time-step from the local Courant number
29  rDeltaT.ref() = max
30  (
31  1/dimensionedScalar("maxDeltaT", dimTime, maxDeltaT),
33  /((2*maxCo)*mesh.V())
34  );
35 
36  // Update tho boundary values of the reciprocal time-step
37  rDeltaT.correctBoundaryConditions();
38 
40 
41  Info<< "Flow time scale min/max = "
42  << gMin(1/rDeltaT.primitiveField())
43  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
44 }
Foam::fvc::surfaceSum
tmp< GeometricField< Type, fvPatchField, volMesh > > surfaceSum(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcSurfaceIntegrate.C:135
rDeltaTSmoothingCoeff
scalar rDeltaTSmoothingCoeff(pimpleDict.getOrDefault< scalar >("rDeltaTSmoothingCoeff", 0.1))
maxPhi
surfaceScalarField maxPhi("maxPhi", phi)
phasei
label phasei
Definition: pEqn.H:27
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
maxCo
scalar maxCo(pimpleDict.get< scalar >("maxCo"))
pimpleDict
const dictionary & pimpleDict
Definition: setRDeltaT.H:32
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
maxDeltaT
scalar maxDeltaT(pimpleDict.getOrDefault< scalar >("maxDeltaT", GREAT))
forAll
forAll(phases, phasei)
Definition: setRDeltaT.H:23
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
trDeltaT
tmp< volScalarField > trDeltaT
Definition: createRDeltaT.H:3
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
phases
multiphaseSystem::phaseModelList & phases
Definition: createFields.H:12
Foam::fvc::smooth
void smooth(volScalarField &field, const scalar coeff)
Definition: fvcSmooth.C:44