All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
updateFluxes.H
Go to the documentation of this file.
1 //Gradients and divergence
2 //---------Start---------
4 
5 divUf = tr(gradUf);
6 
8 
10 
12 
14 
15 //---------End---------
16 
17 //Continuity equation fluxes
18 //---------Start---------
19 
20 surfaceVectorField wHatf =
21  (tauQGDf / rhof)
22  *
23  (
24  ((rhof*Uf) & gradUf) + gradPf
25  );
26 surfaceVectorField wf =
27  wHatf + (tauQGDf/rhof)*(Uf * divRhoUf);
28 
29 jm = rhoLnf*(Uf - wf);
30 phiJm = mesh.Sf() & jm;
31 phi = mesh.Sf() & linearInterpolate(rhoU);
32 //---------End---------
33 
34 // Fluxes for momentum balance equation
35 //---------Start---------
36 
37 phiJmU = mesh.Sf() & (jm * Uf);
38 phiP = mesh.Sf() * pf;
39 
40 Info << "gammaf=" << gammaf << endl;
41 
42 //QGD diffusive fluxes
43 Pif =
44  rhof*Uf*wHatf
45  +
46  tauQGDf *
47  (
48  I*((Uf & gradPf) + (gammaf*p1f*divUf))
49  );
50 
52 {
53  Pif +=
54  muf*
55  (
56  gradUf
57 // +
58 // Foam::T(gradUf)
59 // -
60 // (2.0/3.0)*I*divUf
61  );
62 }
63 else
64 {
65  //tauMCPtr() = muf*linearInterpolate(Foam::T(fvc::grad(U)) - (2.0/3.0)*I*fvc::div(U));
66  //tauMCPtr() = muf*(Foam::T(gradUf) - (2.0/3.0)*I*divUf);
67  tauMCPtr() = linearInterpolate(turbulence->muEff() * dev2(Foam::T(fvc::grad(U))));
68  phiTauMC = mesh.Sf() & tauMCPtr();
69 }
70 
71 phiPi = mesh.Sf() & Pif;
72 
73 //---------End---------
74 
75 // Fluxes for energy balance equation
76 //---------Start---------
77 phih2 = (0.25*hQGDf*hQGDf*(fvc::snGrad(U) * fvc::snGrad(p)) & mesh.Sf());
78 
79 phiJmH = (E1f*(Uf - wf)) & mesh.Sf();
80 phih2.setOriented(phiJmH.oriented()());
81 phiJmH -= phih2;
82 
83 qf =
84  -tauQGDf*rhof*Uf*
85  (
86  (
87  (Uf&gradef)
88  -
89  (
90  Uf&
91  (
93  )
94  )
95  )
96  );
97 
99 {
100  qf -=
101  alphauf*gradef;
102 }
103 
104 phiQ = mesh.Sf() & qf;
105 
106 phiPiU = mesh.Sf() & (Pif & Uf);
107 
108 //---------End---------
gradPf
Definition: updateFluxes.H:33
pf
Definition: updateFields.H:21
const surfaceScalarField & tauQGDf
Definition: createFields.H:55
phiP
Definition: updateFluxes.H:51
divRhoUf
Definition: updateFluxes.H:15
const surfaceScalarField & hQGDf
Definition: createFields.H:54
muf
Definition: updateFields.H:38
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh)
Switch implicitDiffusion(thermo.implicitDiffusion())
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
phiPi
Definition: updateFluxes.H:85
phiJm
Definition: updateFluxes.H:39
qf
Definition: updateFluxes.H:97
surfaceVectorField wf
Definition: updateFluxes.H:32
alphauf
Definition: updateFields.H:42
surfaceVectorField wHatf
——–End———
Definition: updateFluxes.H:26
rhof
Definition: updateFields.H:27
p1f
Definition: updateFields.H:35
autoPtr< surfaceTensorField > tauMCPtr
phiJmU
——–End———
Definition: updateFluxes.H:50
volScalarField & e
Definition: createFields.H:50
tmp< surfaceScalarField > div(const volVectorField &vF)
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())
Info<< "Thermo corrected"<< endl;autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:59
volScalarField & T
Definition: createFields.H:53
tmp< surfaceVectorField > grad(const volScalarField &vF)
E1f
Definition: updateFields.H:34
volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U)
surfaceScalarField phih2("phih2",(0.25 *hQGDf *hQGDf *(fvc::snGrad(U)*fvc::snGrad(p))&mesh.Sf()))
volScalarField & p
Definition: createFields.H:52
surfaceScalarField rhoLnf
Pif
Definition: updateFluxes.H:53