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-2024 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 "isothermalFluid.H"
27 #include "fvcFlux.H"
28 #include "fvcSmooth.H"
29 #include "fvcSurfaceIntegrate.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 void Foam::solvers::isothermalFluid::setRDeltaT()
34 {
35  const volScalarField& psi = thermo.psi();
36 
37  volScalarField& rDeltaT = trDeltaT.ref();
38 
39  const dictionary& pimpleDict = pimple.dict();
40 
41  const scalar maxCo
42  (
43  pimpleDict.lookupOrDefault<scalar>("maxCo", 0.8)
44  );
45 
46  const scalar rDeltaTSmoothingCoeff
47  (
48  pimpleDict.lookupOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.02)
49  );
50 
51  const volScalarField rDeltaT0("rDeltaT0", rDeltaT);
52 
53  // Set the reciprocal time-step from the local Courant number
54  rDeltaT.internalFieldRef() =
55  fvc::surfaceSum(mag(phi))()()/((2*maxCo)*mesh.V()*rho());
56 
57  // Set the reciprocal time-step from the local acoustic Courant number
58  if (pimple.transonic())
59  {
61 
62  rDeltaT.internalFieldRef() =
63  max
64  (
65  rDeltaT(),
66  fvc::surfaceSum(mag(phid))()()/((2*maxCo)*mesh.V()*psi())
67  );
68  }
69 
70  // Clip to user-defined maximum and minimum time-steps
71  scalar minRDeltaT = gMin(rDeltaT.primitiveField());
72  if (pimpleDict.found("maxDeltaT") || minRDeltaT < rootVSmall)
73  {
74  const scalar clipRDeltaT = 1/pimpleDict.lookup<scalar>("maxDeltaT");
75  rDeltaT.max(clipRDeltaT);
76  minRDeltaT = max(minRDeltaT, clipRDeltaT);
77  }
78  if (pimpleDict.found("minDeltaT"))
79  {
80  const scalar clipRDeltaT = 1/pimpleDict.lookup<scalar>("minDeltaT");
81  rDeltaT.min(clipRDeltaT);
82  minRDeltaT = min(minRDeltaT, clipRDeltaT);
83  }
84 
85  Info<< "Flow time scale min/max = "
86  << gMin(1/rDeltaT.primitiveField()) << ", " << 1/minRDeltaT << endl;
87 
88  // Update the boundary values of the reciprocal time-step
89  rDeltaT.correctBoundaryConditions();
90 
91  if (rDeltaTSmoothingCoeff < 1.0)
92  {
93  fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
94  }
95 
96  Info<< "Smoothed flow time scale min/max = "
97  << gMin(1/rDeltaT.primitiveField())
98  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
99 
100  // Limit rate of change of time scale
101  // - reduce as much as required
102  // - only increase at a fraction of old time scale
103  if
104  (
105  pimpleDict.found("rDeltaTDampingCoeff")
107  )
108  {
109  // Damping coefficient (1-0)
110  const scalar rDeltaTDampingCoeff
111  (
112  pimpleDict.lookup<scalar>("rDeltaTDampingCoeff")
113  );
114 
115  rDeltaT =
116  rDeltaT0
117  *max(rDeltaT/rDeltaT0, scalar(1) - rDeltaTDampingCoeff);
118 
119  Info<< "Damped flow time scale min/max = "
120  << gMin(1/rDeltaT.primitiveField())
121  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
122  }
123 }
124 
125 
126 // ************************************************************************* //
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
virtual label startTimeIndex() const
Return start time index.
Definition: Time.C:785
bool transonic() const
Flag to indicate to solve using transonic algorithm.
virtual const volScalarField & psi() const =0
Compressibility [s^2/m^2].
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:107
const Time & runTime
Time.
Definition: solver.H:104
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
const surfaceScalarField & phi
Mass-flux field.
tmp< volScalarField > trDeltaT
Optional LTS reciprocal time-step field.
const volVectorField & U
Velocity field.
const volScalarField & rho
Reference to the continuity density field.
const fluidThermo & thermo
Reference to the fluid thermophysical properties.
scalar maxCo
Calculate the face-flux of the given field.
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.
const volScalarField & psi
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
void smooth(volScalarField &field, const scalar coeff)
Definition: fvcSmooth.C:35
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
SurfaceField< scalar > surfaceScalarField
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
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)