All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
reactingFoam.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-2021 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  reactingFoam
26 
27 Description
28  Transient solver for turbulent flow of compressible reacting fluids with
29  optional mesh motion and mesh topology changes.
30 
31  Uses the flexible PIMPLE (PISO-SIMPLE) solution for time-resolved and
32  pseudo-transient simulations.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "fvCFD.H"
37 #include "dynamicFvMesh.H"
38 #include "fluidReactionThermo.H"
39 #include "combustionModel.H"
42 #include "multivariateScheme.H"
43 #include "pimpleControl.H"
44 #include "pressureReference.H"
45 #include "CorrectPhi.H"
46 #include "fvModels.H"
47 #include "fvConstraints.H"
48 #include "localEulerDdtScheme.H"
49 #include "fvcSmooth.H"
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 int main(int argc, char *argv[])
54 {
55  #include "postProcess.H"
56 
57  #include "setRootCaseLists.H"
58  #include "createTime.H"
59  #include "createDynamicFvMesh.H"
60  #include "createDyMControls.H"
61  #include "initContinuityErrs.H"
62  #include "createFields.H"
63  #include "createFieldRefs.H"
64  #include "createRhoUfIfPresent.H"
65 
66  turbulence->validate();
67 
68  if (!LTS)
69  {
70  #include "compressibleCourantNo.H"
71  #include "setInitialDeltaT.H"
72  }
73 
74  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75 
76  Info<< "\nStarting time loop\n" << endl;
77 
78  while (pimple.run(runTime))
79  {
80  #include "readDyMControls.H"
81 
82  // Store divrhoU from the previous mesh so that it can be mapped
83  // and used in correctPhi to ensure the corrected phi has the
84  // same divergence
85  autoPtr<volScalarField> divrhoU;
86  if (correctPhi)
87  {
88  divrhoU = new volScalarField
89  (
90  "divrhoU",
92  );
93  }
94 
95  if (LTS)
96  {
97  #include "setRDeltaT.H"
98  }
99  else
100  {
101  #include "compressibleCourantNo.H"
102  #include "setDeltaT.H"
103  }
104 
105  runTime++;
106 
107  Info<< "Time = " << runTime.timeName() << nl << endl;
108 
109  // --- Pressure-velocity PIMPLE corrector loop
110  while (pimple.loop())
111  {
112  if (!pimple.flow())
113  {
114  if (pimple.models())
115  {
116  fvModels.correct();
117  }
118 
119  if (pimple.thermophysics())
120  {
121  #include "YEqn.H"
122  #include "EEqn.H"
123  }
124  }
125  else
126  {
127  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
128  {
129  // Store momentum to set rhoUf for introduced faces.
130  autoPtr<volVectorField> rhoU;
131  if (rhoUf.valid())
132  {
133  rhoU = new volVectorField("rhoU", rho*U);
134  }
135 
137 
138  // Do any mesh changes
139  mesh.update();
140 
141  if (mesh.changing())
142  {
143  MRF.update();
144 
145  if (correctPhi)
146  {
147  #include "correctPhi.H"
148  }
149 
150  if (checkMeshCourantNo)
151  {
152  #include "meshCourantNo.H"
153  }
154  }
155  }
156 
157  if (pimple.firstPimpleIter() && !pimple.simpleRho())
158  {
159  #include "rhoEqn.H"
160  }
161 
162  if (pimple.models())
163  {
164  fvModels.correct();
165  }
166 
167  #include "UEqn.H"
168 
169  if (pimple.thermophysics())
170  {
171  #include "YEqn.H"
172  #include "EEqn.H"
173  }
174 
175  // --- Pressure corrector loop
176  while (pimple.correct())
177  {
178  #include "../../compressible/rhoPimpleFoam/pEqn.H"
179  }
180 
181  if (pimple.turbCorr())
182  {
183  turbulence->correct();
184  thermophysicalTransport->correct();
185  }
186  }
187  }
188 
189  rho = thermo.rho();
190 
191  runTime.write();
192 
193  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
194  << " ClockTime = " << runTime.elapsedClockTime() << " s"
195  << nl << endl;
196  }
197 
198  Info<< "End\n" << endl;
199 
200  return 0;
201 }
202 
203 
204 // ************************************************************************* //
fluidReactionThermo & thermo
Definition: createFields.H:28
pimpleNoLoopControl & pimple
IOMRFZoneList & MRF
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:316
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:240
correctPhi
checkMeshCourantNo
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Creates and initialises the velocity field rhoUf if required.
dynamicFvMesh & mesh
phi
Definition: correctPhi.H:3
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::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
static const char nl
Definition: Ostream.H:260
bool LTS
Definition: createRDeltaT.H:1
fluidReactionThermophysicalTransportModel & thermophysicalTransport
moveMeshOuterCorrectors
Foam::fvModels & fvModels
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.
autoPtr< surfaceVectorField > rhoUf
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...