EEqns.H
Go to the documentation of this file.
1 for (int Ecorr=0; Ecorr<nEnergyCorrectors; Ecorr++)
2 {
3  fluid.correctEnergyTransport();
4 
5  autoPtr<phaseSystem::heatTransferTable>
6  heatTransferPtr(fluid.heatTransfer());
7 
8  phaseSystem::heatTransferTable& heatTransfer = heatTransferPtr();
9 
11  {
12  phaseModel& phase = phases[phasei];
13 
14  const volScalarField& alpha = phase;
15  const volScalarField& rho = phase.rho();
16  const volVectorField& U = phase.U();
17 
18  tmp<fvScalarMatrix> EEqn(phase.heEqn());
19 
20  if (EEqn.valid())
21  {
22  EEqn =
23  (
24  EEqn
25  ==
26  *heatTransfer[phase.name()]
27  + alpha*rho*(U&g)
28  + fvOptions(alpha, rho, phase.thermo().he())
29  );
30 
31  EEqn->relax();
32  fvOptions.constrain(EEqn());
33  EEqn->solve();
34  }
35  }
36 }
37 
38 fluid.correctThermo();
39 
41 {
42  phaseModel& phase = phases[phasei];
43 
44  Info<< phase.name() << " min/max T "
45  << min(phase.thermo().T()).value()
46  << " - "
47  << max(phase.thermo().T()).value()
48  << endl;
49 }
fluid
multiphaseSystem & fluid
Definition: createFields.H:10
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
forAll
forAll(phases, phasei)
Definition: EEqns.H:40
g
const dimensionedVector & g
Definition: setRegionFluidFields.H:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
U
U
Definition: pEqn.H:46
Foam::Info
messageStream Info
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
rho
rho
Definition: pEqn.H:3
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
EEqn
fvScalarMatrix EEqn(fvm::div(phi, he)+(he.name()=="e" ? fvc::div(phi, volScalarField("Ekp", 0.5 *magSqr(U)+p/rho)) :fvc::div(phi, volScalarField("K", 0.5 *magSqr(U)))) - fvm::laplacian(turb.alphaEff(), he)==rho *(U &g)+rad.Sh(thermo)+fvOptions(rho, he))
phasei
label phasei
Definition: pEqn.H:37
phases
multiphaseSystem::phaseModelList & phases
Definition: createFields.H:11
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:16