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