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