QGDsolver
The open source CFD toolbox
Main Page
Modules
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Friends
Macros
Groups
QHDFoam
QHDFoam.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
OpenFOAM is free software: you can redistribute it and/or modify it
15
under the terms of the GNU General Public License as published by
16
the Free Software Foundation, either version 3 of the License, or
17
(at your option) any later version.
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
You should have received a copy of the GNU General Public License
23
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24
Application
25
QHDFoam
26
Description
27
Solver for unsteady 3D turbulent flow of incompressible viscous fluid
28
governed by quasi-hydrodynamic dynamic (QHD) equations.
29
30
QHD system of equations has been developed by scientific group from
31
Keldysh Institute of Applied Mathematics,
32
see http://elizarova.imamod.ru/selection-of-papers.html
33
A comprehensive description of QGD equations and their applications
34
can be found here:
35
\verbatim
36
Elizarova, T.G.
37
"Quasi-Gas Dynamic equations"
38
Springer, 2009
39
\endverbatim
40
A brief of theory on QGD and QHD system of equations:
41
\verbatim
42
Elizarova, T.G. and Sheretov, Y.V.
43
"Theoretical and numerical analysis of quasi-gasdynamic and quasi-hydrodynamic
44
equations"
45
J. Computational Mathematics and Mathematical Physics, vol. 41, no. 2, pp 219-234,
46
2001
47
\endverbatim
48
49
Developed by UniCFD group (www.unicfd.ru) of ISP RAS (www.ispras.ru).
50
\*---------------------------------------------------------------------------*/
51
52
#include "fvCFD.H"
53
#include "
QHD.H
"
54
#include "turbulentFluidThermoModel.H"
55
#include "turbulentTransportModel.H"
56
57
58
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59
60
int
main
(
int
argc,
char
*argv[])
61
{
62
#define NO_CONTROL
63
#include "postProcess.H"
64
65
#include "setRootCase.H"
66
#include "createTime.H"
67
#include "createMesh.H"
68
#include "
createFields.H
"
69
#include "
createFaceFields.H
"
70
#include "
createFaceFluxes.H
"
71
#include "createTimeControls.H"
72
73
turbulence
->validate();
74
75
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76
77
// Courant numbers used to adjust the time-step
78
scalar CoNum = 0.0;
79
scalar meanCoNum = 0.0;
80
81
Info<<
"\nStarting time loop\n"
<< endl;
82
83
while
(runTime.run())
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 "
QHDCourantNo.H
"
106
#include "
setDeltaT-QGDQHD.H
"
107
108
runTime++;
109
110
Info<<
"Time = "
<< runTime.timeName() << nl << endl;
111
112
// --- Store old time values
113
U
.oldTime();
114
T
.oldTime();
115
turbulence
->correct();
116
117
#include "
QHDpEqn.H
"
118
119
#include "
QHDUEqn.H
"
120
121
#include "
QHDTEqn.H
"
122
123
if
(
p
.needReference())
124
{
125
p
+= dimensionedScalar
126
(
127
"p"
,
128
p
.dimensions(),
129
pRefValue - getRefCellValue(
p
, pRefCell)
130
);
131
}
132
133
runTime.write();
134
135
Info<<
"ExecutionTime = "
<< runTime.elapsedCpuTime() <<
" s"
136
<<
" ClockTime = "
<< runTime.elapsedClockTime() <<
" s"
137
<< nl << endl;
138
139
}
140
141
142
143
Info<<
"End\n"
<< endl;
144
145
return
0;
146
}
147
148
// ************************************************************************* //
createFaceFields.H
QHDCourantNo.H
Calculates the mean and maximum wave speed based Courant Numbers.
main
int main(int argc, char *argv[])
Definition:
particlesQGDFoam.C:45
U
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh)
createFaceFluxes.H
setDeltaT-QGDQHD.H
Reset the timestep to maintain a constant maximum courant Number. Reduction of time-step is immediate...
createFields.H
QHDTEqn.H
Solver for unsteady 3D turbulent flow of incompressible fluid governed by quasi-hydrodynamic dynamic ...
updateFluxes.H
turbulence
Info<< "Thermo corrected"<< endl;autoPtr< compressible::turbulenceModel > turbulence
Definition:
createFields.H:59
T
volScalarField & T
Definition:
createFields.H:53
QHD.H
Includation for QHD solver. Equations of state for QHD based on density. Class for fvsc namespace...
QHDUEqn.H
Solution of momentum equation for QHD solver.
QHDpEqn.H
Solution of continuity equation for QGD solver.
updateFields.H
p
volScalarField & p
Definition:
createFields.H:52
Generated by
1.8.5