All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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) 2023-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 "shockFluid.H"
27 #include "fvcSmooth.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 void Foam::solvers::shockFluid::setRDeltaT(const surfaceScalarField& amaxSf)
32 {
33  volScalarField& rDeltaT = trDeltaT.ref();
34 
35  const dictionary& pimpleDict = pimple.dict();
36 
37  const scalar maxCo
38  (
39  pimpleDict.lookupOrDefault<scalar>("maxCo", 0.8)
40  );
41 
42  const scalar rDeltaTSmoothingCoeff
43  (
44  pimpleDict.lookupOrDefault<scalar>
45  (
46  "rDeltaTSmoothingCoeff",
47  0.02
48  )
49  );
50 
51  // Set the reciprocal time-step from the local Courant number
52  rDeltaT.internalFieldRef() =
53  fvc::surfaceSum(amaxSf)()()/((2*maxCo)*mesh.V());
54 
55  // Clip to user-defined maximum and minimum time-steps
56  scalar minRDeltaT = gMin(rDeltaT.primitiveField());
57  if (pimpleDict.found("maxDeltaT") || minRDeltaT < rootVSmall)
58  {
59  const scalar clipRDeltaT = 1/pimpleDict.lookup<scalar>("maxDeltaT");
60  rDeltaT.max(clipRDeltaT);
61  minRDeltaT = max(minRDeltaT, clipRDeltaT);
62  }
63  if (pimpleDict.found("minDeltaT"))
64  {
65  const scalar clipRDeltaT = 1/pimpleDict.lookup<scalar>("minDeltaT");
66  rDeltaT.min(clipRDeltaT);
67  minRDeltaT = min(minRDeltaT, clipRDeltaT);
68  }
69 
70  Info<< "Flow time scale min/max = "
71  << gMin(1/rDeltaT.primitiveField()) << ", " << 1/minRDeltaT << endl;
72 
73  // Update the boundary values of the reciprocal time-step
74  rDeltaT.correctBoundaryConditions();
75 
76  fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
77 
78  Info<< "Smoothed flow time scale min/max = "
79  << gMin(1/rDeltaT.primitiveField())
80  << ", " << gMax(1/rDeltaT.primitiveField()) << endl;
81 }
82 
83 
84 // ************************************************************************* //
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 fvMesh & mesh
Region mesh.
Definition: solver.H:101
tmp< volScalarField > trDeltaT
Optional LTS reciprocal time-step field.
Definition: shockFluid.H:152
scalar maxCo
Provides functions smooth spread and sweep which use the FvFaceCellWave algorithm to smooth and redis...
tmp< VolField< Type > > surfaceSum(const SurfaceField< Type > &ssf)
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
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Type gMin(const FieldField< Field, Type > &f)
Type gMax(const FieldField< Field, Type > &f)