All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
QGDEEqn.H
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 
13 License
14  This file is part of QGDsolver, based on OpenFOAM library.
15 
16  OpenFOAM is free software: you can redistribute it and/or modify it
17  under the terms of the GNU General Public License as published by
18  the Free Software Foundation, either version 3 of the License, or
19  (at your option) any later version.
20 
21  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
22  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
23  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
24  for more details.
25 
26  You should have received a copy of the GNU General Public License
27  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
28 
29 Description
30  Equation of energy for QGD solver
31 
32 SourceFiles
33  QGDFoam.C
34 
35 \*----------------------------------------------------------------------------*/
36 
37  fvScalarMatrix EEqn
38  (
39  fvm::ddt(rhoE)
40  + fvc::div(phiJmH)
41  + fvc::div(phiQ)
42  - fvc::div(phiPiU)
44  );
45 
46  EEqn.solve();
47 
48  // Correct energy
49  e = rhoE/rho - 0.5*magSqr(U);
50  e.correctBoundaryConditions();
51 
52  // Solve diffusive QGD & NS part
54  {
55  solve
56  (
57  fvm::ddt(rho, e) - fvc::ddt(rho,e)
58  - fvm::laplacian(alphauf, e)
59  ==
60  rhoESu
61  );
62 
63  rhoE = rho*(e + 0.5*magSqr(U));
64  }
65  else
66  {
67  solve
68  (
69  fvm::ddt(rho,e) - fvc::ddt(rhoE)
70  ==
71  rhoESu
72  );
73  }
74 
75  rhoE.boundaryFieldRef() == rho.boundaryField()*
76  (e.boundaryField() + 0.5*magSqr(U.boundaryField()));
77 //
78 //END-OF-FILE
79 //
80 
81 
surfaceScalarField phiSigmaDotU("phiSigmaDotU", sigmaDotUPtr()&mesh.Sf())
fvScalarMatrix EEqn(fvm::ddt(rhoE)+fvc::div(phiJmH)+fvc::div(phiQ)-fvc::div(phiPiU)-fvc::div(phiSigmaDotU))
fvScalarMatrix rhoESu(e, rho.dimensions()*e.dimensions()*dimVolume/dimTime)
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh)
Switch implicitDiffusion(thermo.implicitDiffusion())
phiJmH
——–End———
Definition: updateFluxes.H:95
EEqn solve()
phiPiU
Definition: updateFluxes.H:115
alphauf
Definition: updateFields.H:42
volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *e+rho *0.5 *magSqr(U))
volScalarField & e
Definition: createFields.H:50
tmp< surfaceScalarField > div(const volVectorField &vF)
phiQ
Definition: updateFluxes.H:113
volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho())