All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
particlesQGDFoam.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  particlesQGDFoam
30 
31 Description
32  Solver for unsteady 3D turbulent flow of perfect gas governed by
33  quasi-gas dynamic (QGD) equations at all Mach numbers (from 0 to
34  infinity) coupled with particles motion.
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  Developed by UniCFD group (www.unicfd.ru) of ISP RAS (www.ispras.ru).
41 
42 
43 \*---------------------------------------------------------------------------*/
44 
45 #include "fvCFD.H"
46 #include "QGD.H"
47 #include "fvOptions.H"
48 #include "turbulentFluidThermoModel.H"
49 
50 #include "basicThermoCloud.H"
51 
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 
54 int main(int argc, char *argv[])
55 {
56  #define NO_CONTROL
57  #include "postProcess.H"
58 
59  #include "setRootCase.H"
60  #include "createTime.H"
61  #include "createMesh.H"
62  #include "createFields.H"
63  #include "createFaceFields.H"
64  #include "createFaceFluxes.H"
65  #include "createTimeControls.H"
66  #include "createFvOptions.H"
67 
68  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70  // Courant numbers used to adjust the time-step
71  scalar CoNum = 0.0;
72  scalar meanCoNum = 0.0;
73 
74  Info<< "\nStarting time loop\n" << endl;
75 
76  while (runTime.run())
77  {
78  /*
79  *
80  * Update QGD viscosity
81  *
82  */
83  turbulence->correct();
84 
85  /*
86  *
87  * Update fields
88  *
89  */
90  #include "updateFields.H"
91 
92  /*
93  *
94  * Update fluxes
95  *
96  */
97  #include "updateFluxes.H"
98 
99  /*
100  *
101  * Update time step
102  *
103  */
104  #include "readTimeControls.H"
105  #include "QGDCourantNo.H"
106  #include "setDeltaT-QGDQHD.H"
107 
108  runTime++;
109 
110  Info<< "Time = " << runTime.timeName() << nl << endl;
111 
112  parcels.evolve();
113 
114  // --- Store old time values
115  rho.oldTime();
116  rhoU.oldTime();
117  U.oldTime();
118  rhoE.oldTime();
119  e.oldTime();
120 
121  // --- Solve density
122  #include "QGDRhoEqn.H"
123 
124  // --- Solve momentum
125  rhoUSu = parcels.SU(U);
126  #include "QGDUEqn.H"
127 
128  //--- Solve energy
129  rhoESu = parcels.Sh(e);
130  #include "QGDEEqn.H"
131 
132  if ( (min(e).value() <= 0.0) || (min(rho).value() <= 0.0) )
133  {
134  U.write();
135  e.write();
136  rho.write();
137  }
138 
139  thermo.correct();
140 
141  // Correct pressure
142  p.ref() =
143  rho()
144  /psi();
145  p.correctBoundaryConditions();
146  rho.boundaryFieldRef() = psi.boundaryField()*p.boundaryField();
147 
148  runTime.write();
149 
150  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
151  << " ClockTime = " << runTime.elapsedClockTime() << " s"
152  << nl << endl;
153  }
154 
155  Info<< "End\n" << endl;
156 
157  return 0;
158 }
159 
160 // ************************************************************************* //
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.
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
fvVectorMatrix rhoUSu(U, rho.dimensions()*U.dimensions()*dimVolume/dimTime)
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