All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
QHDUEqn.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 License
13  This file is part of 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 Group
29  grpQGDcommon
30 
31 Description
32  Solution of momentum equation for QHD solver
33 
34 \*---------------------------------------------------------------------------*/
35 
36  gradPf = fvsc::grad(p);
37  Wf = tauQGDf*((Uf & gradUf) + gradPf/rhof - BdFrcf);
38 
39  surfaceVectorField phiUfWf = mesh.Sf() & (Uf * Wf);
40  phiUfWf.setOriented(true);
41  phiUf = qgdFlux(phi,U,Uf);
42  phiUf.setOriented(true);
43  phiUf -= phiUfWf;
44 
45  // --- Solve U
47  {
48  solve
49  (
50  fvm::ddt(U)
51  +
53  -
54  fvm::laplacian(muf/rhof,U)
55  -
57  ==
58  -
59  fvc::grad(p)/rho
60  +
62  +
63  USu
64  );
65  }
66  else
67  {
68  solve
69  (
70  fvm::ddt(U)
71  +
73  -
74  fvc::laplacian(muf/rhof,U)
75  -
77  ==
78  -
79  fvc::grad(p)/rho
80  +
81  BdFrc
82  +
83  USu
84  );
85  }
86 
87 // ************************************************************************* //
gradPf
Definition: updateFluxes.H:33
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)
BdFrcf
Definition: updateFields.H:36
const surfaceScalarField & tauQGDf
Definition: createFields.H:55
muf
Definition: updateFields.H:38
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh)
Switch implicitDiffusion(thermo.implicitDiffusion())
EEqn solve()
gradUf
Definition: updateFields.H:34
tmp< GeometricField< T, Foam::fvsPatchField, Foam::surfaceMesh > > qgdInterpolate(const GeometricField< T, Foam::fvPatchField, Foam::volMesh > &psi)
rhof
Definition: updateFields.H:27
BdFrc
Definition: updateFields.H:35
tmp< surfaceScalarField > div(const volVectorField &vF)
fvVectorMatrix USu(U, U.dimensions()*dimVolume/dimTime)
Wf
Definition: QHDUEqn.H:31
Uf
Definition: updateFields.H:30
phiUf
Definition: QHDUEqn.H:35
volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho())
volScalarField & T
Definition: createFields.H:53
tmp< surfaceVectorField > grad(const volScalarField &vF)
surfaceVectorField phiUfWf
Definition: QHDUEqn.H:33
volScalarField & p
Definition: createFields.H:52