hEqn.H
Go to the documentation of this file.
1 {
2  fvScalarMatrix hEqn
3  (
5  - (
6  thermo.isotropic()
7  ? fvm::laplacian(betav*thermo.alpha(), h, "laplacian(alpha,h)")
8  : fvm::laplacian(betav*taniAlpha(), h, "laplacian(alpha,h)")
9  )
10  ==
11  fvOptions(rho, h)
12  );
13 
15  {
16  surfaceScalarField phihMesh
17  (
19  );
20 
21  hEqn -= fvc::div(phihMesh);
22  }
23 
24  hEqn.relax();
25 
26  fvOptions.constrain(hEqn);
27 
28  hEqn.solve(); //mesh.solver(h.select(finalIter)));
29 
30  fvOptions.correct(h);
31 
32  thermo.correct();
33 
34  Info<< "Min/max T:" << min(thermo.T()).value() << ' '
35  << max(thermo.T()).value() << endl;
36 
37  radiation->correct();
38 }
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:40
fv::options & fvOptions
psiReactionThermo & thermo
Definition: createFields.H:28
tmp< volSymmTensorField > taniAlpha
dynamicFvMesh & mesh
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:40
const volScalarField & betav
bool meshFluxCorr(false)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
const dimensionedScalar h
Planck constant.
messageStream Info
Information stream (stdout output on master, null elsewhere)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))