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-2023 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  // and maximum and minimum time-steps
49  rDeltaT.ref() =
50  fvc::surfaceSum(mag(phi))()()/((2*maxCo)*mesh.V()*rho());
51  if (pimpleDict.found("maxDeltaT"))
52  {
53  rDeltaT.max(1/pimpleDict.lookup<scalar>("maxDeltaT"));
54  }
55  if (pimpleDict.found("minDeltaT"))
56  {
57  rDeltaT.min(1/pimpleDict.lookup<scalar>("minDeltaT"));
58  }
59 
60  Info<< " Flow = "
61  << 1/gMax(rDeltaT.primitiveField()) << ", "
62  << 1/gMin(rDeltaT.primitiveField()) << endl;
63  }
64 
65  // Maximum change in cell temperature per iteration
66  // (relative to previous value)
67  const scalar alphaTemp(pimpleDict.lookupOrDefault("alphaTemp", 0.05));
68 
69  // Heat release rate time scale
70  if (alphaTemp < 1)
71  {
73  (
74  mag(reaction->Qdot())/(alphaTemp*rho*thermo.Cp()*thermo.T())
75  );
76 
77  Info<< " Temperature = "
78  << 1/(gMax(rDeltaTT.field()) + vSmall) << ", "
79  << 1/(gMin(rDeltaTT.field()) + vSmall) << endl;
80 
81  rDeltaT.ref() = max(rDeltaT(), rDeltaTT);
82  }
83 
84  // Reaction rate time scale
85  if (pimpleDict.found("alphaY"))
86  {
87  // Maximum change in cell concentration per iteration
88  // (relative to reference value)
89  const scalar alphaY(pimpleDict.lookup<scalar>("alphaY"));
90 
91  const dictionary Yref(pimpleDict.subDict("Yref"));
92 
94  (
95  IOobject
96  (
97  "rDeltaTY",
98  runTime.name(),
99  mesh
100  ),
101  mesh,
102  dimensionedScalar(rDeltaT.dimensions(), 0)
103  );
104 
105  bool foundY = false;
106 
107  forAll(Y, i)
108  {
109  if (composition.solve(i))
110  {
111  volScalarField& Yi = Y_[i];
112 
113  if (Yref.found(Yi.name()))
114  {
115  foundY = true;
116  scalar Yrefi = Yref.lookup<scalar>(Yi.name());
117 
118  rDeltaTY.field() = max
119  (
120  mag
121  (
122  reaction->R(Yi)().source()
123  /((Yrefi*alphaY)*(rho*mesh.V()))
124  ),
125  rDeltaTY
126  );
127  }
128  }
129  }
130 
131  if (foundY)
132  {
133  Info<< " Composition = "
134  << 1/(gMax(rDeltaTY.field()) + vSmall) << ", "
135  << 1/(gMin(rDeltaTY.field()) + vSmall) << endl;
136 
137  rDeltaT.ref() = max(rDeltaT(), rDeltaTY);
138  }
139  else
140  {
142  << "Cannot find any active species in Yref " << Yref
143  << endl;
144  }
145  }
146 
147  // Update the boundary values of the reciprocal time-step
148  rDeltaT.correctBoundaryConditions();
149 
150  // Smoothing parameter (0-1) when smoothing iterations > 0
151  const scalar rDeltaTSmoothingCoeff
152  (
153  pimpleDict.lookupOrDefault<scalar>("rDeltaTSmoothingCoeff", 0.1)
154  );
155 
156  // Spatially smooth the time scale field
157  if (rDeltaTSmoothingCoeff < 1)
158  {
159  fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
160  }
161 
162  // Limit rate of change of time scale
163  // - reduce as much as required
164  // - only increase at a fraction of old time scale
165  if
166  (
167  pimpleDict.found("rDeltaTDampingCoeff")
169  )
170  {
171  // Damping coefficient (1-0)
172  const scalar rDeltaTDampingCoeff
173  (
174  pimpleDict.lookup<scalar>("rDeltaTDampingCoeff")
175  );
176 
177  rDeltaT = max
178  (
179  rDeltaT,
180  (scalar(1) - rDeltaTDampingCoeff)*rDeltaT0
181  );
182  }
183 
184  // Update the boundary values of the reciprocal time-step
185  rDeltaT.correctBoundaryConditions();
186 
187  Info<< " Overall = "
188  << 1/gMax(rDeltaT.primitiveField())
189  << ", " << 1/gMin(rDeltaT.primitiveField()) << endl;
190 }
191 
192 
193 // ************************************************************************* //
#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:773
bool solve(label speciei) const
Return true if the specie should be solved for.
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.
virtual const dictionary & dict() const
Return the solution dictionary.
pimpleNoLoopControl pimple
PIMPLE inner-loop controls.
Definition: solver.H:100
const Time & runTime
Time.
Definition: solver.H:97
const fvMesh & mesh
Region mesh.
Definition: solver.H:94
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.
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
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
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)