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