All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
zQGDFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2016-2018 OpenCFD Ltd.
7 -------------------------------------------------------------------------------
8  QGDsolver | Copyright (C) 2016-2018 ISP RAS (www.unicfd.ru)
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  QGDFoam
28 
29 Description
30  Solver for unsteady 3D turbulent flow of perfect gas governed by
31  quasi-gas dynamic (QGD) equations at high Mach numbers (from 2 to
32  infinity).
33 
34  QGD system of equations has been developed by scientific group from
35  Keldysh Institute of Applied Mathematics,
36  see http://elizarova.imamod.ru/selection-of-papers.html
37 
38  A comprehensive description of QGD equations and their applications
39  can be found here:
40  \verbatim
41  Elizarova, T.G.
42  "Quasi-Gas Dynamic equations"
43  Springer, 2009
44  \endverbatim
45 
46  A brief of theory on QGD and QHD system of equations:
47  \verbatim
48  Elizarova, T.G. and Sheretov, Y.V.
49  "Theoretical and numerical analysis of quasi-gasdynamic and quasi-hydrodynamic
50  equations"
51  J. Computational Mathematics and Mathematical Physics, vol. 41, no. 2, pp 219-234,
52  2001
53  \endverbatim
54 
55  Developed by UniCFD group (www.unicfd.ru) of ISP RAS (www.ispras.ru).
56 
57 
58 \*---------------------------------------------------------------------------*/
59 
60 #include "fvCFD.H"
61 #include "QGD.H"
62 #include "fvOptions.H"
63 #include "turbulentFluidThermoModel.H"
64 #include "directionInterpolate.H"
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 tmp<surfaceScalarField> logMean(const surfaceScalarField& a, const surfaceScalarField& b)
69 {
70  tmp<surfaceScalarField> tres (1.0 / (6.0*a) + 4.0 / (3.0*(a+b)) + 1.0 / (6.0*b));
71 
72  return tres;
73 }
74 
75 int main(int argc, char *argv[])
76 {
77  #define NO_CONTROL
78  #include "postProcess.H"
79 
80  #include "setRootCase.H"
81  #include "createTime.H"
82  #include "createMesh.H"
83  #include "createFields.H"
84  #include "createFaceFields.H"
85  #include "createFaceFluxes.H"
86  #include "createTimeControls.H"
87  #include "createFvOptions.H"
88 
89  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
90 
91  // Courant numbers used to adjust the time-step
92  scalar CoNum = 0.0;
93  scalar meanCoNum = 0.0;
94 
95  Info<< "\nStarting time loop\n" << endl;thermo.correct();
96 
97  while (runTime.run())
98  {
99  /*
100  *
101  * Update QGD viscosity
102  *
103  */
104  turbulence->correct();
105 
106  /*
107  *
108  * Update fields
109  *
110  */
111  #include "updateFields.H"
112 
113  /*
114  *
115  * Update fluxes
116  *
117  */
118  #include "updateFluxes.H"
119 
120  /*
121  *
122  * Update time step
123  *
124  */
125  #include "readTimeControls.H"
126  #include "QGDCourantNo.H"
127  #include "setDeltaT-QGDQHD.H"
128 
129  runTime++;
130 
131  Info<< "Time = " << runTime.timeName() << nl << endl;
132 
133  // --- Store old time values
134  rho.oldTime();
135  rhoU.oldTime();
136  U.oldTime();
137  rhoE.oldTime();
138  e.oldTime();
139 
140  // --- Solve density
141  #include "QGDRhoEqn.H"
142 
143  // --- Solve momentum
144  #include "QGDUEqn.H"
145 
146  //--- Solve energy
147  #include "QGDEEqn.H"
148 
149  if ( (min(e).value() <= 0.0) || (min(rho).value() <= 0.0) )
150  {
151  U.write();
152  e.write();
153  rho.write();
154  thermo.tauQGD().write();
155  }
156 
157  thermo.correct();
158 
159  // Correct pressure
160  p.ref() =
161  rho()
162  /psi();
163  p.correctBoundaryConditions();
164  rho.boundaryFieldRef() = psi.boundaryField()*p.boundaryField();
165 
166  if(runTime.write())
167  {
168  thermo.tauQGD().write();
169  }
170 
171  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
172  << " ClockTime = " << runTime.elapsedClockTime() << " s"
173  << nl << endl;
174  }
175 
176  Info<< "End\n" << endl;
177 
178  return 0;
179 }
180 
181 // ************************************************************************* //
tmp< surfaceScalarField > logMean(const surfaceScalarField &a, const surfaceScalarField &b)
Definition: zQGDFoam.C:59
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.
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...
Calculates the mean and maximum wave speed based Courant Numbers.
volScalarField & e
Definition: createFields.H:50
rhoQGDThermo & thermo
Definition: createFields.H:47
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
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