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 "multicomponentFluid.H"
27 #include "fvcSmooth.H"
28 
29 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
30 
31 void Foam::solvers::multicomponentFluid::setRDeltaT()
32 {
33  volScalarField& rDeltaT = trDeltaT.ref();
34 
35  const dictionary& pimpleDict = pimple.dict();
36 
37  Info<< "Time scales min/max:" << endl;
38 
39  // Cache old reciprocal time scale field
40  const volScalarField rDeltaT0("rDeltaT0", rDeltaT);
41 
42  // Flow time scale
43  {
44  // Maximum flow Courant number
45  const scalar maxCo(pimpleDict.lookup<scalar>("maxCo"));
46 
47  // Set the reciprocal time-step from the local Courant number
48  rDeltaT.internalFieldRef() =
49  fvc::surfaceSum(mag(phi))()()/((2*maxCo)*mesh.V()*rho());
50 
51  // Clip to user-defined maximum and minimum time-steps
52  scalar minRDeltaT = gMin(rDeltaT.primitiveField());
53  if (pimpleDict.found("maxDeltaT") || minRDeltaT < rootVSmall)
54  {
55  const scalar clipRDeltaT = 1/pimpleDict.lookup<scalar>("maxDeltaT");
56  rDeltaT.max(clipRDeltaT);
57  minRDeltaT = max(minRDeltaT, clipRDeltaT);
58  }
59  if (pimpleDict.found("minDeltaT"))
60  {
61  const scalar clipRDeltaT = 1/pimpleDict.lookup<scalar>("minDeltaT");
62  rDeltaT.min(clipRDeltaT);
63  minRDeltaT = min(minRDeltaT, clipRDeltaT);
64  }
65 
66  Info<< "Flow time scale min/max = "
67  << gMin(1/rDeltaT.primitiveField()) << ", " << 1/minRDeltaT << endl;
68  }
69 
70  // Maximum change in cell temperature per iteration
71  // (relative to previous value)
72  const scalar alphaTemp(pimpleDict.lookupOrDefault("alphaTemp", 0.05));
73 
74  // Heat release rate time scale
75  if (alphaTemp < 1)
76  {
78  (
79  mag(reaction->Qdot())/(alphaTemp*rho*thermo.Cp()*thermo.T())
80  );
81 
82  Info<< " Temperature = "
83  << 1/(gMax(rDeltaTT.primitiveField()) + vSmall) << ", "
84  << 1/(gMin(rDeltaTT.primitiveField()) + vSmall) << endl;
85 
86  rDeltaT.internalFieldRef() = max(rDeltaT(), rDeltaTT);
87  }
88 
89  // Reaction rate time scale
90  if (pimpleDict.found("alphaY"))
91  {
92  // Maximum change in cell concentration per iteration
93  // (relative to reference value)
94  const scalar alphaY(pimpleDict.lookup<scalar>("alphaY"));
95 
96  const dictionary Yref(pimpleDict.subDict("Yref"));
97 
99  (
100  IOobject
101  (
102  "rDeltaTY",
103  runTime.name(),
104  mesh
105  ),
106  mesh,
107  dimensionedScalar(rDeltaT.dimensions(), 0)
108  );
109 
110  bool foundY = false;
111 
112  forAll(Y, i)
113  {
114  if (thermo_.solveSpecie(i))
115  {
116  volScalarField& Yi = Y_[i];
117 
118  if (Yref.found(Yi.name()))
119  {
120  foundY = true;
121  scalar Yrefi = Yref.lookup<scalar>(Yi.name());
122 
123  rDeltaTY.primitiveFieldRef() = max
124  (
125  mag
126  (
127  reaction->R(Yi)().source()
128  /((Yrefi*alphaY)*(rho*mesh.V()))
129  ),
130  rDeltaTY
131  );
132  }
133  }
134  }
135 
136  if (foundY)
137  {
138  Info<< " Composition = "
139  << 1/(gMax(rDeltaTY.primitiveField()) + vSmall) << ", "
140  << 1/(gMin(rDeltaTY.primitiveField()) + vSmall) << endl;
141 
142  rDeltaT.internalFieldRef() = max(rDeltaT(), rDeltaTY);
143  }
144  else
145  {
147  << "Cannot find any active species in Yref " << Yref
148  << endl;
149  }
150  }
151 
152  // Update the boundary values of the reciprocal time-step
153  rDeltaT.correctBoundaryConditions();
154 
155  // Smoothing parameter (0-1) when smoothing iterations > 0
156  const scalar rDeltaTSmoothingCoeff
157  (
158  pimpleDict.lookupOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.1)
159  );
160 
161  // Spatially smooth the time scale field
162  if (rDeltaTSmoothingCoeff < 1)
163  {
164  fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
165  }
166 
167  // Limit rate of change of time scale
168  // - reduce as much as required
169  // - only increase at a fraction of old time scale
170  if
171  (
172  pimpleDict.found("rDeltaTDampingCoeff")
174  )
175  {
176  // Damping coefficient (1-0)
177  const scalar rDeltaTDampingCoeff
178  (
179  pimpleDict.lookup<scalar>("rDeltaTDampingCoeff")
180  );
181 
182  rDeltaT = max
183  (
184  rDeltaT,
185  (scalar(1) - rDeltaTDampingCoeff)*rDeltaT0
186  );
187  }
188 
189  // Update the boundary values of the reciprocal time-step
190  rDeltaT.correctBoundaryConditions();
191 
192  Info<< " Overall = "
193  << 1/gMax(rDeltaT.primitiveField())
194  << ", " << 1/gMin(rDeltaT.primitiveField()) << endl;
195 }
196 
197 
198 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
DimensionedField< Type, GeoMesh > Internal
Type of the internal field from which this GeometricField is derived.
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:167
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
virtual label startTimeIndex() const
Return start time index.
Definition: Time.C:785
virtual const volScalarField & T() const =0
Temperature [K].
virtual const volScalarField & Cp() const =0
Heat capacity at constant pressure [J/kg/K].
const word & name() const
Return const reference to name.
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
bool solveSpecie(const label speciei) const
Should the given specie be solved for? I.e., is it active and.
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 volScalarField & rho
Reference to the continuity density field.
const PtrList< volScalarField > & Y
Reference to the composition.
fluidMulticomponentThermo & thermo_
PtrList< volScalarField > & Y_
autoPtr< combustionModel > reaction
const fluidMulticomponentThermo & thermo
Reference to the fluid thermophysical properties.
scalar maxCo
Provides functions smooth spread and sweep which use the FvFaceCellWave algorithm to smooth and redis...
#define WarningInFunction
Report a warning using Foam::Warning.
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
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)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.