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-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  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 "dynamicFvMesh.H"
39 #include "fluidReactionThermo.H"
40 #include "combustionModel.H"
43 #include "multivariateScheme.H"
44 #include "pimpleControl.H"
45 #include "pressureReference.H"
47 #include "CorrectPhi.H"
48 #include "fvModels.H"
49 #include "fvConstraints.H"
50 #include "localEulerDdtScheme.H"
51 #include "fvcSmooth.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 int main(int argc, char *argv[])
56 {
57  #include "postProcess.H"
58 
59  #include "setRootCaseLists.H"
60  #include "createTime.H"
61  #include "createDynamicFvMesh.H"
62  #include "createDyMControls.H"
63  #include "initContinuityErrs.H"
64  #include "createFields.H"
65  #include "createFieldRefs.H"
66  #include "createRhoUfIfPresent.H"
67 
68  turbulence->validate();
69 
70  if (!LTS)
71  {
72  #include "compressibleCourantNo.H"
73  #include "setInitialDeltaT.H"
74  }
75 
76  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77 
78  Info<< "\nStarting time loop\n" << endl;
79 
80  while (pimple.run(runTime))
81  {
82  #include "readDyMControls.H"
83 
84  // Store divrhoU from the previous mesh so that it can be mapped
85  // and used in correctPhi to ensure the corrected phi has the
86  // same divergence
87  autoPtr<volScalarField> divrhoU;
88  if (correctPhi)
89  {
90  divrhoU = new volScalarField
91  (
92  "divrhoU",
94  );
95  }
96 
97  if (LTS)
98  {
99  #include "setRDeltaT.H"
100  }
101  else
102  {
103  #include "compressibleCourantNo.H"
104  #include "setDeltaT.H"
105  }
106 
107  runTime++;
108 
109  Info<< "Time = " << runTime.timeName() << nl << endl;
110 
111  // --- Pressure-velocity PIMPLE corrector loop
112  while (pimple.loop())
113  {
114  if (!pimple.flow())
115  {
116  if (pimple.models())
117  {
118  fvModels.correct();
119  }
120 
121  if (pimple.thermophysics())
122  {
123  #include "YEqn.H"
124  #include "EEqn.H"
125  }
126  }
127  else
128  {
129  if (pimple.firstPimpleIter() || moveMeshOuterCorrectors)
130  {
131  // Store momentum to set rhoUf for introduced faces.
132  autoPtr<volVectorField> rhoU;
133  if (rhoUf.valid())
134  {
135  rhoU = new volVectorField("rhoU", rho*U);
136  }
137 
139 
140  // Do any mesh changes
141  mesh.update();
142 
143  if (mesh.changing())
144  {
145  gh = (g & mesh.C()) - ghRef;
146  ghf = (g & mesh.Cf()) - ghRef;
147 
148  MRF.update();
149 
150  if (correctPhi)
151  {
152  #include "correctPhi.H"
153  }
154 
155  if (checkMeshCourantNo)
156  {
157  #include "meshCourantNo.H"
158  }
159  }
160  }
161 
162  if (pimple.firstPimpleIter() && !pimple.simpleRho())
163  {
164  #include "rhoEqn.H"
165  }
166 
167  if (pimple.models())
168  {
169  fvModels.correct();
170  }
171 
172  #include "UEqn.H"
173 
174  if (pimple.thermophysics())
175  {
176  #include "YEqn.H"
177  #include "EEqn.H"
178  }
179 
180  // --- Pressure corrector loop
181  while (pimple.correct())
182  {
183  #include "../../../heatTransfer/buoyantPimpleFoam/pEqn.H"
184  }
185 
186  if (pimple.turbCorr())
187  {
188  turbulence->correct();
189  thermophysicalTransport->correct();
190  }
191  }
192  }
193 
194  rho = thermo.rho();
195 
196  runTime.write();
197 
198  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
199  << " ClockTime = " << runTime.elapsedClockTime() << " s"
200  << nl << endl;
201  }
202 
203  Info<< "End\n" << endl;
204 
205  return 0;
206 }
207 
208 
209 // ************************************************************************* //
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
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:188
U
Definition: pEqn.H:72
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 FaceCellWave algorithm to smooth and redistr...