All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
SRFQHDFoam.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  SRFQHDFoam
30 
31 Description
32  Solver for unsteady 3D turbulent flow of incompressible fluid governed by
33  quasi-hydrodynamic dynamic (QHD) equations with support for single rotating
34  frame (SRF) models.
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  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 "QHD.H"
63 #include "turbulentFluidThermoModel.H"
64 #include "turbulentTransportModel.H"
65 #include "SRFModel.H"
66 
67 
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70 int main(int argc, char *argv[])
71 {
72  #define NO_CONTROL
73  #include "postProcess.H"
74 
75  #include "setRootCase.H"
76  #include "createTime.H"
77  #include "createMesh.H"
78  #include "createFields.H"
79  #include "createFaceFields.H"
80  #include "createFaceFluxes.H"
81  #include "createTimeControls.H"
82 
83  turbulence->validate();
84 
85  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 
87  // Courant numbers used to adjust the time-step
88  scalar CoNum = 0.0;
89  scalar meanCoNum = 0.0;
90 
91  Info<< "\nStarting time loop\n" << endl;
92 
93  while (runTime.run())
94  {
95  /*
96  *
97  * Update fields
98  *
99  */
100  #include "updateFields.H"
101 
102  /*
103  *
104  * Update fluxes
105  *
106  */
107  #include "updateFluxes.H"
108 
109  /*
110  *
111  * Update time step
112  *
113  */
114  #include "readTimeControls.H"
115  #include "QHDCourantNo.H"
116  #include "setDeltaT-QGDQHD.H"
117 
118  runTime++;
119 
120  Info<< "Time = " << runTime.timeName() << nl << endl;
121 
122  // --- Store old time values
123  U.oldTime();
124  T.oldTime();
125  turbulence->correct();
126 
127  #include "QHDpEqn.H"
128 
129  #include "QHDUEqn.H"
130 
131  if (p.needReference())
132  {
133  p += dimensionedScalar
134  (
135  "p",
136  p.dimensions(),
137  pRefValue - getRefCellValue(p, pRefCell)
138  );
139  }
140 
141  runTime.write();
142 
143  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
144  << " ClockTime = " << runTime.elapsedClockTime() << " s"
145  << nl << endl;
146 
147  }
148 
149 
150 
151  Info<< "End\n" << endl;
152 
153  return 0;
154 }
155 
156 
157 // ************************************************************************* //
Calculates the mean and maximum wave speed based Courant Numbers.
Updates fluxes for continuity equation.
Updates fields.
int main(int argc, char *argv[])
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh)
Reset the timestep to maintain a constant maximum courant Number. Reduction of time-step is immediate...
Info<< "Thermo corrected"<< endl;autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:59
volScalarField & T
Definition: createFields.H:53
Includation for QHD solver. Equations of state for QHD based on density. Class for fvsc namespace...
Solution of momentum equation for QHD solver.
Solution of continuity equation for QGD solver.
Creates the face fields: linear interpolation of fields from volumes to face centers.
volScalarField & p
Definition: createFields.H:52
Creates the face-flux fields.