TEqns.H
Go to the documentation of this file.
1 {
2  volScalarField kByCp1("kByCp1", alpha1*(k1/Cp1/rho1 + sqr(Ct)*nut2/Prt));
3  volScalarField kByCp2("kByCp2", alpha2*(k2/Cp2/rho2 + nut2/Prt));
4 
6  (
7  fvm::ddt(alpha1, T1)
8  + fvm::div(alphaPhi1, T1)
9  - fvm::laplacian(kByCp1, T1)
10  ==
11  heatTransferCoeff*T2/Cp1/rho1
12  - fvm::Sp(heatTransferCoeff/Cp1/rho1, T1)
13  + alpha1*Dp1Dt/Cp1/rho1
14  );
15 
17  (
18  fvm::ddt(alpha2, T2)
19  + fvm::div(alphaPhi2, T2)
20  - fvm::laplacian(kByCp2, T2)
21  ==
22  heatTransferCoeff*T1/Cp2/rho2
23  - fvm::Sp(heatTransferCoeff/Cp2/rho2, T2)
24  + alpha2*Dp2Dt/Cp2/rho2
25  );
26 
27  T1Eqn.relax();
28  T1Eqn.solve();
29 
30  T2Eqn.relax();
31  T2Eqn.solve();
32 
33  // Update compressibilities
34  psi1 = 1.0/(R1*T1);
35  psi2 = 1.0/(R2*T2);
36 }
psi2
Definition: TEqns.H:35
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
const surfaceScalarField & alphaPhi2
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
const surfaceScalarField & alphaPhi1
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
volScalarField kByCp2("kByCp2", alpha2 *(k2/Cp2/rho2+nut2/Prt))
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const volScalarField & alpha1
alpha2
Definition: alphaEqn.H:115
dimensionedScalar Prt("Prt", dimless, laminarTransport)
psi1
Definition: TEqns.H:34
const dimensionedScalar & rho2
Definition: createFields.H:40
const dimensionedScalar & rho1
Definition: createFields.H:39
fvScalarMatrix T2Eqn(fvm::ddt(alpha2, T2)+fvm::div(alphaPhi2, T2) - fvm::laplacian(kByCp2, T2)==heatTransferCoeff *T1/Cp2/rho2 - fvm::Sp(heatTransferCoeff/Cp2/rho2, T2)+alpha2 *Dp2Dt/Cp2/rho2)
fvScalarMatrix T1Eqn(fvm::ddt(alpha1, T1)+fvm::div(alphaPhi1, T1) - fvm::laplacian(kByCp1, T1)==heatTransferCoeff *T2/Cp1/rho1 - fvm::Sp(heatTransferCoeff/Cp1/rho1, T1)+alpha1 *Dp1Dt/Cp1/rho1)
zeroField Sp
Definition: alphaSuSp.H:2