All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
mulesQHDFoam.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  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 Application
26  mulesQHDFoam
27 
28 Description
29  Solver for unsteady 3D turbulent flow of incompressible viscous fluid
30  governed by quasi-hydrodynamic dynamic (QHD) equations. The
31  temperature scalar transport equation is solved by MULES technology.
32 
33  QHD system of equations has been developed by scientific group from
34  Keldysh Institute of Applied Mathematics,
35  see http://elizarova.imamod.ru/selection-of-papers.html
36  A comprehensive description of QGD equations and their applications
37  can be found here:
38  \verbatim
39  Elizarova, T.G.
40  "Quasi-Gas Dynamic equations"
41  Springer, 2009
42  \endverbatim
43  A brief of theory on QGD and QHD system of equations:
44  \verbatim
45  Elizarova, T.G. and Sheretov, Y.V.
46  "Theoretical and numerical analysis of quasi-gasdynamic and quasi-hydrodynamic
47  equations"
48  J. Computational Mathematics and Mathematical Physics, vol. 41, no. 2, pp 219-234,
49  2001
50  \endverbatim
51  Developed by UniCFD group (www.unicfd.ru) of ISP RAS (www.ispras.ru).
52 
53  MULES: Multidimensional universal limiter for explicit solution.
54 
55  Solve a convective-only transport equation using an explicit universal
56  multi-dimensional limiter.
57  See https://www.openfoam.com/documentation/guides/latest/api/MULES_8H_source.html
58 
59 \*---------------------------------------------------------------------------*/
60 
61 #include "fvCFD.H"
62 #include "upwind.H"
63 #include "CMULES.H"
64 #include "EulerDdtScheme.H"
65 #include "localEulerDdtScheme.H"
66 #include "CrankNicolsonDdtScheme.H"
67 #include "subCycle.H"
68 #include "pimpleControl.H"
69 #include "fvOptions.H"
70 #include "CorrectPhi.H"
71 #include "localEulerDdtScheme.H"
72 #include "fvcSmooth.H"
73 #include "QHD.H"
74 #include "turbulentTransportModel.H"
75 #include "turbulentFluidThermoModel.H"
76 
77 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
78 
79 int main(int argc, char* argv[])
80 {
81 #define NO_CONTROL
82 #include "postProcess.H"
83 #include "setRootCase.H"
84 #include "createTime.H"
85 #include "createMesh.H"
86 #include "createFields.H"
87 #include "createFaceFields.H"
88 #include "createFaceFluxes.H"
89 #include "createTimeControls.H"
90  pimpleControl pimple(mesh);
91  turbulence->validate();
92  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
93  // Courant numbers used to adjust the time-step
94  scalar CoNum = 0.0;
95  scalar meanCoNum = 0.0;
96  Info << "\nStarting time loop\n" << endl;
97 
98  while (runTime.run())
99  {
100  /*
101  *
102  * Update fields
103  *
104  */
105 #include "updateFields.H"
106  /*
107  *
108  * Update fluxes
109  *
110  */
111 #include "updateFluxes.H"
112  /*
113  *
114  * Update time step
115  *
116  */
117 #include "readTimeControls.H"
118 #include "QHDCourantNo.H"
119 #include "setDeltaT-QGDQHD.H"
120  runTime++;
121  Info << "Time = " << runTime.timeName() << nl << endl;
122  // --- Store old time values
123  U.oldTime();
124  T.oldTime();
125  turbulence->correct();
126 #include "QHDpEqn.H"
127 #include "QHDUEqn.H"
128  phiTf = qgdFlux(phi, T, Tf);
129  // --- Solve T
130  {
131  volScalarField& alpha1 = T;
132 #include "alphaControls.H"
133 #include "MULESTEqn.H"
134  }
135 
136  if (p.needReference())
137  {
138  p += dimensionedScalar
139  (
140  "p",
141  p.dimensions(),
142  pRefValue - getRefCellValue(p, pRefCell)
143  );
144  }
145 
146  runTime.write();
147  Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
148  << " ClockTime = " << runTime.elapsedClockTime() << " s"
149  << nl << endl;
150  }
151 
152  Info << "End\n" << endl;
153  return 0;
154 }
155 // ************************************************************************* //
tmp< GeometricField< T, Foam::fvsPatchField, Foam::surfaceMesh > > qgdFlux(const GeometricField< scalar, Foam::fvsPatchField, Foam::surfaceMesh > &flux, const GeometricField< T, Foam::fvPatchField, Foam::volMesh > &psi, const GeometricField< T, Foam::fvsPatchField, Foam::surfaceMesh > &psif, const word fluxName)
Calculates the mean and maximum wave speed based Courant Numbers.
int main(int argc, char *argv[])
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh)
surfaceScalarField phiTf("phiTf", phi *Tf)
Reset the timestep to maintain a constant maximum courant Number. Reduction of time-step is immediate...
Tf
Definition: updateFields.H:33
Solves temperature scalar transport equation using MULES (Multidimensional universal limiter for expl...
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.
volScalarField & p
Definition: createFields.H:52