All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
reactingLagrangianQGDFoam.C
Go to the documentation of this file.
1 
2 /*---------------------------------------------------------------------------*\
3  ========= |
4  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5  \\ / O peration |
6  \\ / A nd | www.openfoam.com
7  \\/ M anipulation |
8 -------------------------------------------------------------------------------
9  Copyright (C) 2011-2016 OpenFOAM Foundation
10  Copyright (C) 2019 OpenCFD Ltd.
11  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
12 -------------------------------------------------------------------------------
13 License
14  This file is part of QGDsolver library, based on OpenFOAM+.
15 
16  OpenFOAM is free software: you can redistribute it and/or modify it
17  under the terms of the GNU General Public License as published by
18  the Free Software Foundation, either version 3 of the License, or
19  (at your option) any later version.
20 
21  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24  for more details.
25 
26  You should have received a copy of the GNU General Public License
27  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28 
29 Application
30  reactinLagrangianQGDFoam
31 
32 Description
33  Solver for unsteady 3D turbulent flow of perfect gas governed by
34  quasi-gas dynamic (QGD) equations at all Mach numbers (from 0 to
35  infinity) coupled with chemistry reactions.
36 
37  QGD system of equations has been developed by scientific group from
38  Keldysh Institute of Applied Mathematics,
39  see http://elizarova.imamod.ru/selection-of-papers.html
40 
41  Developed by UniCFD group (www.unicfd.ru) of ISP RAS (www.ispras.ru).
42 
43 
44 \*---------------------------------------------------------------------------*/
45 
46 #include "fvCFD.H"
47 #include "mcQGD.H"
48 #include "fvOptions.H"
49 #include "turbulentFluidThermoModel.H"
50 #include "basicReactingCloud.H"
51 #include "CombustionModel.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 int main(int argc, char *argv[])
56 {
57  #define NO_CONTROL
58  #include "postProcess.H"
59 
60  #include "setRootCase.H"
61  #include "createTime.H"
62  #include "createMesh.H"
63  #include "createFields.H"
64  Info << "Basic fields are created" << endl;
65  #include "createFaceFields.H"
66  #include "createFaceFluxes.H"
67  #include "createTimeControls.H"
68  #include "createFvOptions.H"
69 
70  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72  // Courant numbers used to adjust the time-step
73  scalar CoNum = 0.0;
74  scalar meanCoNum = 0.0;
75 
76  Info<< "\nStarting time loop\n" << endl;
77 
78  while (runTime.run())
79  {
80  /*
81  *
82  * Update QGD viscosity
83  *
84  */
85  turbulence->correct();
86 
87  /*
88  *
89  * Update fields
90  *
91  */
92  #include "updateFields.H"
93 
94  /*
95  *
96  * Update fluxes
97  *
98  */
99  #include "updateFluxes.H"
100 
101  /*
102  *
103  * Update time step
104  *
105  */
106  #include "readTimeControls.H"
107  #include "QGDCourantNo.H"
108  #include "setDeltaT-QGDQHD.H"
109 
110  runTime++;
111 
112  Info<< "Time = " << runTime.timeName() << nl << endl;
113 
114  parcels.evolve();
115 
116  // --- Store old time values
117  rho.oldTime();
118  rhoU.oldTime();
119  U.oldTime();
120  rhoE.oldTime();
121  e.oldTime();
122  forAll(Y,i)
123  {
124  Y[i].oldTime();
125  }
126 
127  // --- Solve density
128  #include "QGDRhoEqn.H"
129 
130  // --- Solve for mass fractions
131  #include "QGDYEqn.H"
132 
133  // --- Solve momentum
134  rhoUSu = parcels.SU(U);
135  #include "QGDUEqn.H"
136 
137  //--- Solve energy
138  rhoESu = parcels.Sh(e) + Qdot;
139  #include "addEnergyFluxes.H"
140  #include "QGDEEqn.H"
141 
142  if ( (min(T).value() <= 0.0) || (min(rho).value() <= 0.0) )
143  {
144  U.write();
145  T.write();
146  rho.write();
147  forAll(Y, i)
148  Y[i].write();
149  p.write();
150  }
151 
152  thermo.correct();
153 
154  // Correct pressure
155  p.ref() =
156  rho()
157  /psi();
158  p.correctBoundaryConditions();
159  rho.boundaryFieldRef() = psi.boundaryField()*p.boundaryField();
160 
161  runTime.write();
162  if (runTime.outputTime())
163  {
164  turbulence->alphaEff()().write();
165  }
166 
167  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
168  << " ClockTime = " << runTime.elapsedClockTime() << " s"
169  << nl << endl;
170  }
171 
172  Info<< "End\n" << endl;
173 
174  return 0;
175 }
176 
177 // ************************************************************************* //
fvScalarMatrix rhoESu(e, rho.dimensions()*e.dimensions()*dimVolume/dimTime)
int main(int argc, char *argv[])
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh)
Equation of energy for QGD solver.
Qdot
Definition: QGDYEqn.H:33
volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *e+rho *0.5 *magSqr(U))
Reset the timestep to maintain a constant maximum courant Number. Reduction of time-step is immediate...
forAll(Y, i)
Definition: QGDYEqn.H:36
Calculates the mean and maximum wave speed based Courant Numbers.
volScalarField & e
Definition: createFields.H:50
PtrList< volScalarField > & Y
Definition: createFields.H:29
rhoQGDThermo & thermo
Definition: createFields.H:47
fvVectorMatrix rhoUSu(U, rho.dimensions()*U.dimensions()*dimVolume/dimTime)
Solution of momentum equation for QGD solver.
const volScalarField & psi
Definition: createFields.H:20
volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho())
Info<< "Thermo corrected"<< endl;autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:59
volScalarField & T
Definition: createFields.H:53
Solution of momentum equation for QGD solver.
volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U)
volScalarField & p
Definition: createFields.H:52