solveSolid.H
Go to the documentation of this file.
1 if (finalIter)
2 {
3  mesh.data::add("finalIteration", true);
4 }
5 
6 {
7  for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
8  {
9  fvScalarMatrix hEqn
10  (
12  - (
13  thermo.isotropic()
14  ? fvm::laplacian(betav*thermo.alpha(), h, "laplacian(alpha,h)")
15  : fvm::laplacian(betav*taniAlpha(), h, "laplacian(alpha,h)")
16  )
17  ==
18  fvOptions(rho, h)
19  );
20 
21  hEqn.relax();
22 
23  fvOptions.constrain(hEqn);
24 
25  hEqn.solve(mesh.solver(h.select(finalIter)));
26 
27  fvOptions.correct(h);
28  }
29 }
30 
31 thermo.correct();
32 
33 Info<< "Min/max T:" << min(thermo.T()).value() << ' '
34  << max(thermo.T()).value() << endl;
35 
36 if (finalIter)
37 {
38  mesh.data::remove("finalIteration");
39 }
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
const int nNonOrthCorr
fv::options & fvOptions
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
psiReactionThermo & thermo
Definition: createFields.H:31
tmp< volSymmTensorField > taniAlpha
dynamicFvMesh & mesh
const volScalarField & betav
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const dimensionedScalar h
Planck constant.
Definition: createFields.H:6
messageStream Info