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
QHDDyMFoam
QHDDyMFoam.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
QHDDyMFoam
30
31
Description
32
Solver for unsteady 3D turbulent flow of incompressible fluid governed by
33
quasi-hydrodynamic dynamic (QHD) equations with support for deforming
34
meshes.
35
36
QHD 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
#include "fvCFD.H"
45
#include "dynamicFvMesh.H"
46
#include "
QHD.H
"
47
#include "turbulentFluidThermoModel.H"
48
#include "turbulentTransportModel.H"
49
50
51
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52
53
int
main
(
int
argc,
char
*argv[])
54
{
55
#define NO_CONTROL
56
#include "postProcess.H"
57
58
#include "setRootCase.H"
59
#include "createTime.H"
60
#include "createDynamicFvMesh.H"
61
#include "
createFields.H
"
62
bool
checkMeshCourantNo
63
(
64
thermo
.subDict(
"QGD"
).lookupOrDefault(
"checkMeshCourantNo"
,
false
)
65
);
66
#include "
createFaceFields.H
"
67
#include "
createFaceFluxes.H
"
68
#include "createTimeControls.H"
69
70
turbulence
->validate();
71
72
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73
74
// Courant numbers used to adjust the time-step
75
scalar CoNum = 0.0;
76
scalar meanCoNum = 0.0;
77
78
Info<<
"\nStarting time loop\n"
<< endl;
79
80
while
(runTime.run())
81
{
82
checkMeshCourantNo =
thermo
.subDict(
"QGD"
).lookupOrDefault
83
(
84
"checkMeshCourantNo"
,
85
checkMeshCourantNo
86
);
87
88
#include "readTimeControls.H"
89
#include "
QHDCourantNo.H
"
90
#include "
setDeltaT-QGDQHD.H
"
91
92
/*
93
*
94
* Update time step
95
*
96
*/
97
runTime++;
98
99
Info<<
"Time = "
<< runTime.timeName() << nl << endl;
100
101
// --- Store old time values
102
U
.oldTime();
103
T
.oldTime();
104
105
/*
106
*
107
* Update the mesh
108
*
109
*/
110
mesh.update();
111
112
/*
113
*
114
* Update fields
115
*
116
*/
117
#include "
updateFields.H
"
118
119
/*
120
*
121
* Update fluxes
122
*
123
*/
124
#include "
updateFluxes.H
"
125
126
#include "
QHDpEqn.H
"
127
128
if
(mesh.changing())
129
{
130
fvc::makeRelative(
phi
,
U
);
131
132
if
(checkMeshCourantNo)
133
{
134
#include "meshCourantNo.H"
135
}
136
}
137
138
turbulence
->correct();
139
140
#include "
QHDUEqn.H
"
141
142
#include "
QHDTEqn.H
"
143
144
if
(
p
.needReference())
145
{
146
p
+= dimensionedScalar
147
(
148
"p"
,
149
p
.dimensions(),
150
pRefValue - getRefCellValue(
p
, pRefCell)
151
);
152
}
153
154
runTime.write();
155
156
Info<<
"ExecutionTime = "
<< runTime.elapsedCpuTime() <<
" s"
157
<<
" ClockTime = "
<< runTime.elapsedClockTime() <<
" s"
158
<< nl << endl;
159
160
}
161
162
163
164
Info<<
"End\n"
<< endl;
165
166
return
0;
167
}
168
169
170
// ************************************************************************* //
QHDCourantNo.H
Calculates the mean and maximum wave speed based Courant Numbers.
updateFluxes.H
Updates fluxes for continuity equation.
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)
createFields.H
setDeltaT-QGDQHD.H
Reset the timestep to maintain a constant maximum courant Number. Reduction of time-step is immediate...
turbulence
Info<< "Thermo corrected"<< endl;autoPtr< compressible::turbulenceModel > turbulence
Definition:
createFields.H:24
p
volScalarField & p
Definition:
createFields.H:14
QHDTEqn.H
Solver for unsteady 3D turbulent flow of incompressible fluid governed by quasi-hydrodynamic dynamic ...
updateFields.H
Updates fields.
phi
phi
Definition:
createFaceFluxes.H:70
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.
thermo
psiQGDThermo & thermo
Definition:
createFields.H:7
QHDpEqn.H
Solution of continuity equation for QGD solver.
T
const volScalarField & T
Definition:
createFields.H:15
createFaceFields.H
Creates the face fields: linear interpolation of fields from volumes to face centers.
createFaceFluxes.H
Creates the face-flux fields.
Generated by
1.8.5