rhoReactingFoam.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) 2011-2019 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 Application
25  rhoReactingFoam
26 
27 Description
28  Solver for combustion with chemical reactions using density based
29  thermodynamics package.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "fvCFD.H"
34 #include "dynamicFvMesh.H"
35 #include "rhoReactionThermo.H"
36 #include "CombustionModel.H"
38 #include "multivariateScheme.H"
39 #include "pimpleControl.H"
40 #include "pressureControl.H"
41 #include "CorrectPhi.H"
42 #include "fvOptions.H"
43 #include "localEulerDdtScheme.H"
44 #include "fvcSmooth.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 int main(int argc, char *argv[])
49 {
50  #include "postProcess.H"
51 
52  #include "setRootCaseLists.H"
53  #include "createTime.H"
54  #include "createDynamicFvMesh.H"
55  #include "createDyMControls.H"
56  #include "initContinuityErrs.H"
57  #include "createFields.H"
58  #include "createFieldRefs.H"
59  #include "createRhoUfIfPresent.H"
60 
61  turbulence->validate();
62 
63  if (!LTS)
64  {
65  #include "compressibleCourantNo.H"
66  #include "setInitialDeltaT.H"
67  }
68 
69  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71  Info<< "\nStarting time loop\n" << endl;
72 
73  while (runTime.run())
74  {
75  #include "readDyMControls.H"
76 
77  // Store divrhoU from the previous mesh so that it can be mapped
78  // and used in correctPhi to ensure the corrected phi has the
79  // same divergence
80  autoPtr<volScalarField> divrhoU;
81  if (correctPhi)
82  {
83  divrhoU = new volScalarField
84  (
85  "divrhoU",
87  );
88  }
89 
90  if (LTS)
91  {
92  #include "setRDeltaT.H"
93  }
94  else
95  {
96  #include "compressibleCourantNo.H"
97  #include "setDeltaT.H"
98  }
99 
100  runTime++;
101 
102  Info<< "Time = " << runTime.timeName() << nl << endl;
103 
104  // --- Pressure-velocity PIMPLE corrector loop
105  while (pimple.loop())
106  {
107  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
108  {
109  // Store momentum to set rhoUf for introduced faces.
110  autoPtr<volVectorField> rhoU;
111  if (rhoUf.valid())
112  {
113  rhoU = new volVectorField("rhoU", rho*U);
114  }
115 
116  // Do any mesh changes
117  mesh.update();
118 
119  if (mesh.changing())
120  {
121  MRF.update();
122 
123  if (correctPhi)
124  {
125  // Calculate absolute flux
126  // from the mapped surface velocity
127  phi = mesh.Sf() & rhoUf();
128 
129  #include "../../../compressible/rhoPimpleFoam/correctPhi.H"
130 
131  // Make the fluxes relative to the mesh-motion
133  }
134 
135  if (checkMeshCourantNo)
136  {
137  #include "meshCourantNo.H"
138  }
139  }
140  }
141 
142  if (pimple.firstPimpleIter() && !pimple.simpleRho())
143  {
144  #include "rhoEqn.H"
145  }
146 
147  #include "UEqn.H"
148  #include "YEqn.H"
149  #include "EEqn.H"
150 
151  // --- Pressure corrector loop
152  while (pimple.correct())
153  {
154  if (pimple.consistent())
155  {
156  #include "../../../compressible/rhoPimpleFoam/pcEqn.H"
157  }
158  else
159  {
160  #include "../../../compressible/rhoPimpleFoam/pEqn.H"
161  }
162  }
163 
164  if (pimple.turbCorr())
165  {
166  turbulence->correct();
167  }
168  }
169 
170  rho = thermo.rho();
171 
172  runTime.write();
173 
174  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
175  << " ClockTime = " << runTime.elapsedClockTime() << " s"
176  << nl << endl;
177  }
178 
179  Info<< "End\n" << endl;
180 
181  return 0;
182 }
183 
184 
185 // ************************************************************************* //
pimpleNoLoopControl & pimple
surfaceScalarField & phi
IOMRFZoneList & MRF
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Info<< "Reading field U\"<< 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(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
correctPhi
checkMeshCourantNo
rhoReactionThermo & thermo
Definition: createFields.H:28
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Creates and initialises the velocity field rhoUf if required.
dynamicFvMesh & mesh
static const char nl
Definition: Ostream.H:265
bool LTS
Definition: createRDeltaT.H:1
moveMeshOuterCorrectors
Calculates and outputs the mean and maximum Courant Numbers.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
U
Definition: pEqn.H:72
messageStream Info
Execute application functionObjects to post-process existing results.
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:75
autoPtr< surfaceVectorField > rhoUf
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...