All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
createFaceFluxes.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 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 Global
29  createFaceFluxes
30 
31 Description
32  Creates the face-flux fields.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 //Gradients and divergence
39 //---------Start---------
40 surfaceVectorField gradPf
41 (
42  "gradPf", fvsc::grad(p)
43 );
44 
45 surfaceTensorField gradUf
46 (
47  "gradUf",
48  fvsc::grad(U)
49 );
50 
51 surfaceScalarField divUf
52 (
53  "divUf",
54  Foam::tr(gradUf)
55 );
56 
57 surfaceVectorField gradef
58 (
59  "gradef",
60  fvsc::grad(e)
61 );
62 
63 surfaceVectorField gradRhof
64 (
65  "gradRhof",
67 );
68 
69 //---------End---------
70 
71 //Continuity equation fluxes
72 //---------Start---------
73 
74 surfaceVectorField rhoUf
75 (
76  "rhoUf",
77  linearInterpolate(rho*U)
78 );
79 
80 surfaceVectorField rhofUf
81 (
82  "rhofUf",
83  rhof*Uf
84 );
85 
86 surfaceTensorField UrhoUf
87 (
88  "UrhoUf",
89  linearInterpolate(U*rhoU)
90 );
91 
92 phi = rhoUf & mesh.Sf();
93 
94 surfaceVectorField rhoW
95 (
96  "rhoW",
98 );
99 
100 surfaceScalarField phiw
101 (
102  "phiwStar",
103  mesh.Sf() & rhoW
104 );
105 
106 surfaceVectorField jm
107 (
108  "jm",
109  rhoUf - rhoW
110 );
111 
112 surfaceScalarField phiJm
113 (
114  "phiJm",
115  mesh.Sf() & jm
116 );
117 
118 //---------End---------
119 
120 // Fluxes for momentum balance equation
121 //---------Start---------
122 surfaceVectorField phiJmU
123 (
124  "phiJmU",
125  mesh.Sf() & (jm * Uf)
126 );
127 
128 surfaceVectorField phiP
129 (
130  "phiP",
131  pf*mesh.Sf()
132 );
133 
134 surfaceTensorField Pif
135 (
136  "Pif",
137  //QGD diffusive fluxes
138  tauQGDf * ((UrhoUf & gradUf) + Uf * gradPf)
139  +
140  tauQGDf * I * ( (Uf & gradPf) + (gammaf * pf * divUf) )
141 );
142 
143 surfaceVectorField phiPi
144 (
145  "phiPi",
146  mesh.Sf() & Pif
147 );
148 
149 autoPtr<surfaceTensorField> tauMCPtr;
150 tauMCPtr.reset
151 (
152  new surfaceTensorField
153  (
154  "tauMC",
155  0.0*muf*
156  (
157  Foam::T(gradUf)//Don't forget to transpose
158  -
159  I*divUf
160  )
161  )
162 );
163 
164 surfaceVectorField phiTauMC
165 (
166  "phiTauMC",
167  tauMCPtr() & mesh.Sf()
168 );
169 
170 //---------End---------
171 
172 // Fluxes for energy balance equation
173 //---------Start---------
174 surfaceScalarField phiJmH
175 (
176  "phiJmH",
177  phiJm * Hf
178 );
179 
180 surfaceVectorField qf
181 (
182  "qf",
183  -
184  tauQGDf*
185  (
186  ((linearInterpolate(rho*(U*U))) & gradef)
187  -
188  (pf * rhof * Uf * (Uf & gradRhof) / rhof / rhof)
189  )
190 );
191 
192 surfaceScalarField phiQ
193 (
194  "phiQ",
195  qf & mesh.Sf()
196 );
197 
198 surfaceScalarField phiPiU
199 (
200  "phiPiU",
201  (Pif & Uf) & mesh.Sf()
202 );
203 
204 autoPtr<surfaceVectorField> sigmaDotUPtr;
205 sigmaDotUPtr.reset
206 (
207  new surfaceVectorField
208  (
209  "sigmaDotU",
210  tauMCPtr() & linearInterpolate(U) //to be updated later
211  )
212 );
213 
214 surfaceScalarField phiSigmaDotU
215 (
216  "phiSigmaDotU",
217  sigmaDotUPtr() & mesh.Sf()
218 );
219 
220 //---------End---------
221 
222 
223 // ************************************************************************* //
surfaceScalarField phiSigmaDotU("phiSigmaDotU", sigmaDotUPtr()&mesh.Sf())
rhoUf
Definition: updateFields.H:14
gradPf
Definition: updateFluxes.H:33
pf
Definition: updateFields.H:21
const surfaceScalarField & tauQGDf
Definition: createFields.H:55
phiP
Definition: updateFluxes.H:51
muf
Definition: updateFields.H:38
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh)
gradRhof
Definition: updateFluxes.H:11
phiJmH
——–End———
Definition: updateFluxes.H:95
gradUf
Definition: updateFields.H:34
phiPiU
Definition: updateFluxes.H:115
gradef
Definition: updateFluxes.H:9
gammaf
Definition: updateFields.H:27
autoPtr< surfaceVectorField > sigmaDotUPtr
phiPi
Definition: updateFluxes.H:85
phiJm
Definition: updateFluxes.H:39
qf
Definition: updateFluxes.H:97
rhof
Definition: updateFields.H:27
autoPtr< surfaceTensorField > tauMCPtr
phiJmU
——–End———
Definition: updateFluxes.H:50
rhoW
——–End———
Definition: updateFluxes.H:22
rhofUf
Definition: updateFields.H:15
UrhoUf
Definition: updateFields.H:18
volScalarField & e
Definition: createFields.H:50
Hf
Definition: updateFields.H:35
phiTauMC
Definition: updateFluxes.H:82
divUf
Definition: updateFluxes.H:7
jm
Definition: updateFluxes.H:37
Uf
Definition: updateFields.H:30
phiQ
Definition: updateFluxes.H:113
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)
volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U)
volScalarField & p
Definition: createFields.H:52
phiw
Definition: updateFluxes.H:31
Pif
Definition: updateFluxes.H:53