All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
createFaceFluxes.H
Go to the documentation of this file.
1 
2 /*Field for debugging parallel execution*/
3 //surfaceVectorField cellNoGradf
4 //(
5 // "cellNoGradf", faceStencil.faceScalarGrad(cellNo)
6 //);
7 //
8 //cellNoGradf.write();
9 
10 //Gradients and divergence
11 //---------Start---------
12 surfaceVectorField gradPf
13 (
14  "gradPf", fvsc::grad(p)
15 );
16 
17 surfaceTensorField gradUf
18 (
19  "gradUf",
20  fvsc::grad(U)
21 );
22 
23 surfaceScalarField divUf
24 (
25  "divUf",
26  Foam::tr(gradUf)
27 );
28 
29 surfaceVectorField gradef
30 (
31  "gradef",
32  fvsc::grad(e)
33 );
34 
35 surfaceVectorField gradRhof
36 (
37  "gradRhof",
39 );
40 
41 //---------End---------
42 
43 //Continuity equation fluxes
44 //---------Start---------
45 
46 surfaceVectorField rhoUf
47 (
48  "rhoUf",
49  linearInterpolate(rho*U)
50 );
51 
52 surfaceVectorField rhofUf
53 (
54  "rhofUf",
55  rhof*Uf
56 );
57 
58 surfaceTensorField UrhoUf
59 (
60  "UrhoUf",
61  linearInterpolate(U*rhoU)
62 );
63 
64 phi = rhoUf & mesh.Sf();
65 
66 surfaceVectorField rhoW
67 (
68  "rhoW",
70 );
71 
72 surfaceScalarField phiw
73 (
74  "phiwStar",
75  mesh.Sf() & rhoW
76 );
77 
78 surfaceVectorField jm
79 (
80  "jm",
81  rhoUf - rhoW
82 );
83 
84 surfaceScalarField phiJm
85 (
86  "phiJm",
87  mesh.Sf() & jm
88 );
89 
90 //---------End---------
91 
92 // Fluxes for momentum balance equation
93 //---------Start---------
94 surfaceVectorField phiJmU
95 (
96  "phiJmU",
97  mesh.Sf() & (jm * Uf)
98 );
99 
100 surfaceVectorField phiP
101 (
102  "phiP",
103  pf*mesh.Sf()
104 );
105 
106 surfaceTensorField Pif
107 (
108  "Pif",
109  //QGD diffusive fluxes
110  tauQGDf * ((UrhoUf & gradUf) + Uf * gradPf)
111  +
112  tauQGDf * I * ( (Uf & gradPf) + (gammaf * pf * divUf) )
113 );
114 
115 surfaceVectorField phiPi
116 (
117  "phiPi",
118  mesh.Sf() & Pif
119 );
120 
121 autoPtr<surfaceTensorField> tauMCPtr;
122 tauMCPtr.reset
123 (
124  new surfaceTensorField
125  (
126  "tauMC",
127  0.0*muf*
128  (
129  Foam::T(gradUf)//Don't forget to transpose
130  -
132  )
133  )
134 );
135 
136 surfaceVectorField phiTauMC
137 (
138  "phiTauMC",
139  tauMCPtr() & mesh.Sf()
140 );
141 
142 //---------End---------
143 
144 // Fluxes for energy balance equation
145 //---------Start---------
146 surfaceScalarField phiJmH
147 (
148  "phiJmH",
149  phiJm * Hf
150 );
151 
152 surfaceVectorField qf
153 (
154  "qf",
155  -
156  tauQGDf*
157  (
158  ((linearInterpolate(rho*(U*U))) & gradef)
159  -
160  (pf * rhof * Uf * (Uf & gradRhof) / rhof / rhof)
161  )
162 );
163 
164 surfaceScalarField phiQ
165 (
166  "phiQ",
167  qf & mesh.Sf()
168 );
169 
170 surfaceScalarField phiPiU
171 (
172  "phiPiU",
173  (Pif & Uf) & mesh.Sf()
174 );
175 
176 autoPtr<surfaceVectorField> sigmaDotUPtr;
177 sigmaDotUPtr.reset
178 (
179  new surfaceVectorField
180  (
181  "sigmaDotU",
182  tauMCPtr() & linearInterpolate(U) //to be updated later
183  )
184 );
185 
186 surfaceScalarField phiSigmaDotU
187 (
188  "phiSigmaDotU",
189  sigmaDotUPtr() & mesh.Sf()
190 );
191 
192 //---------End---------
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