solveSolid.H
Go to the documentation of this file.
1 {
2  while (pimple.correctNonOrthogonal())
3  {
4  fvScalarMatrix hEqn
5  (
6  fvm::ddt(betav*rho, h)
7  - (
8  thermo.isotropic()
9  ? fvm::laplacian(betav*thermo.alpha(), h, "laplacian(alpha,h)")
10  : fvm::laplacian(betav*taniAlpha(), h, "laplacian(alpha,h)")
11  )
12  ==
13  fvOptions(rho, h)
14  );
15 
16  hEqn.relax();
17 
18  fvOptions.constrain(hEqn);
19 
20  hEqn.solve(mesh.solver(h.select(pimples.finalIter())));
21 
22  fvOptions.correct(h);
23  }
24 }
25 
26 thermo.correct();
27 
28 Info<< "Min/max T:" << min(thermo.T()).value() << ' '
29  << max(thermo.T()).value() << endl;
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
fv::options & fvOptions
pimpleNoLoopControl & pimple
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
rhoReactionThermo & thermo
Definition: createFields.H:28
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< volSymmTensorField > taniAlpha
dynamicFvMesh & mesh
const volScalarField & betav
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
volScalarField & h
Planck constant.
messageStream Info