XiDyMFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2015 OpenFOAM Foundation
9  Copyright (C) 2015-2016 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  XiFoam
29 
30 Group
31  grpCombustionSolvers grpMovingMeshSolvers
32 
33 Description
34  Solver for compressible premixed/partially-premixed combustion with
35  turbulence modelling.
36 
37  Combusting RANS code using the b-Xi two-equation model.
38  Xi may be obtained by either the solution of the Xi transport
39  equation or from an algebraic expression. Both approaches are
40  based on Gulder's flame speed correlation which has been shown
41  to be appropriate by comparison with the results from the
42  spectral model.
43 
44  Strain effects are encorporated directly into the Xi equation
45  but not in the algebraic approximation. Further work need to be
46  done on this issue, particularly regarding the enhanced removal rate
47  caused by flame compression. Analysis using results of the spectral
48  model will be required.
49 
50  For cases involving very lean Propane flames or other flames which are
51  very strain-sensitive, a transport equation for the laminar flame
52  speed is present. This equation is derived using heuristic arguments
53  involving the strain time scale and the strain-rate at extinction.
54  the transport velocity is the same as that for the Xi equation.
55 
56 \*---------------------------------------------------------------------------*/
57 
58 #include "fvCFD.H"
59 #include "dynamicFvMesh.H"
60 #include "psiuReactionThermo.H"
62 #include "laminarFlameSpeed.H"
63 #include "ignition.H"
64 #include "Switch.H"
65 #include "pimpleControl.H"
66 #include "CorrectPhi.H"
67 #include "fvOptions.H"
68 
69 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71 int main(int argc, char *argv[])
72 {
73  argList::addNote
74  (
75  "Solver for compressible premixed/partially-premixed combustion with"
76  " turbulence modelling."
77  );
78 
79  #include "postProcess.H"
80 
81  #include "setRootCaseLists.H"
82  #include "createTime.H"
83  #include "createDynamicFvMesh.H"
84  #include "createControl.H"
85  #include "readCombustionProperties.H"
86  #include "readGravitationalAcceleration.H"
87  #include "createFields.H"
88  #include "createFieldRefs.H"
89  #include "initContinuityErrs.H"
90  #include "createRhoUf.H"
91  #include "createControls.H"
92  #include "initContinuityErrs.H"
93  #include "compressibleCourantNo.H"
94  #include "setInitialDeltaT.H"
95 
96  turbulence->validate();
97 
98  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
99 
100  Info<< "\nStarting time loop\n" << endl;
101 
102  while (runTime.run())
103  {
104  #include "createTimeControls.H"
105 
106  {
107  // Store divrhoU from the previous mesh so that it can be mapped
108  // and used in correctPhi to ensure the corrected phi has the
109  // same divergence
110  volScalarField divrhoU
111  (
112  "divrhoU",
114  );
115 
116  #include "compressibleCourantNo.H"
117  #include "setDeltaT.H"
118 
119  ++runTime;
120 
121  Info<< "Time = " << runTime.timeName() << nl << endl;
122 
123  // Store momentum to set rhoUf for introduced faces.
124  volVectorField rhoU("rhoU", rho*U);
125 
126  // Do any mesh changes
127  mesh.update();
128 
129  if (mesh.changing() && correctPhi)
130  {
131  // Calculate absolute flux from the mapped surface velocity
132  phi = mesh.Sf() & rhoUf;
133 
134  #include "correctPhi.H"
135 
136  // Make the fluxes relative to the mesh-motion
138  }
139  }
140 
141  if (mesh.changing() && checkMeshCourantNo)
142  {
143  #include "meshCourantNo.H"
144  }
145 
146  #include "rhoEqn.H"
147  Info<< "rho min/max : " << min(rho).value() << " " << max(rho).value()
148  << endl;
149 
150  // --- Pressure-velocity PIMPLE corrector loop
151  while (pimple.loop())
152  {
153  #include "UEqn.H"
154 
155  #include "ftEqn.H"
156  #include "bEqn.H"
157  #include "EauEqn.H"
158  #include "EaEqn.H"
159 
160  if (!ign.ignited())
161  {
162  thermo.heu() == thermo.he();
163  }
164 
165  // --- Pressure corrector loop
166  while (pimple.correct())
167  {
168  #include "pEqn.H"
169  }
170 
171  if (pimple.turbCorr())
172  {
173  turbulence->correct();
174  }
175  }
176 
177  rho = thermo.rho();
178 
179  runTime.write();
180 
181  runTime.printExecutionTime(Info);
182  }
183 
184  Info<< "End\n" << endl;
185 
186  return 0;
187 }
188 
189 
190 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
psiuReactionThermo.H
fvOptions.H
turbulence
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
laminarFlameSpeed.H
correctPhi
correctPhi
Definition: readDyMControls.H:3
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
rho
rho
Definition: readInitialConditions.H:88
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
pimpleControl.H
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:77
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
setRootCaseLists.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
ignition.H
Switch.H
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
postProcess.H
Execute application functionObjects to post-process existing results.
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
U
U
Definition: pEqn.H:72
rhoUf
autoPtr< surfaceVectorField > rhoUf
Definition: createRhoUfIfPresent.H:33
Foam::nl
constexpr char nl
Definition: Ostream.H:404
createTimeControls.H
Read the control parameters used by setDeltaT.
meshCourantNo.H
Calculates and outputs the mean and maximum Courant Numbers.
createRhoUf.H
Creates and initialises the velocity velocity field rhoUf.
CorrectPhi.H
createTime.H
dynamicFvMesh.H
fvCFD.H
checkMeshCourantNo
checkMeshCourantNo
Definition: readDyMControls.H:9
createDynamicFvMesh.H
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
turbulentFluidThermoModel.H