Go to the documentation of this file. 1 Info<<
"Reading thermal properties\n" <<
endl;
3 IOdictionary thermalProperties
10 IOobject::MUST_READ_IF_MODIFIED,
15 Switch thermalStress(thermalProperties.lookup(
"thermalStress"));
50 autoPtr<volScalarField>
CPtr;
61 const dictionary&
CDict(thermalProperties.subDict(
"C"));
76 dimensionSet(0, 2, -2 , -1, 0),
79 zeroGradientFvPatchField<scalar>::typeName
84 else if (
CType ==
"field")
86 CIO.readOpt() = IOobject::MUST_READ;
100 <<
"Valid type entries are uniform or field for C"
117 const dictionary&
kDict(thermalProperties.subDict(
"k"));
132 dimensionSet(1, 1, -3 , -1, 0),
135 zeroGradientFvPatchField<scalar>::typeName
140 else if (
kType ==
"field")
142 rhoKIO.readOpt() = IOobject::MUST_READ;
156 <<
"Valid type entries are uniform or field for K"
174 const dictionary&
alphaDict(thermalProperties.subDict(
"alpha"));
192 zeroGradientFvPatchField<scalar>::typeName
198 alphaIO.readOpt() = IOobject::MUST_READ;
212 <<
"Valid type entries are uniform or field for alpha"
218 Info<<
"Normalising k : k/rho\n" <<
endl;
221 Info<<
"Calculating thermal coefficients\n" <<
endl;
223 threeKalpha = threeK*
alpha;
word kType(kDict.lookup("type"))
autoPtr< volScalarField > alphaPtr
word alphaType(alphaDict.lookup("type"))
Ostream & endl(Ostream &os)
Add newline and flush stream.
word CType(CDict.lookup("type"))
IOobject alphaIO("alpha", runTime.timeName(0), mesh, IOobject::NO_READ, IOobject::NO_WRITE)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
autoPtr< volScalarField > rhoKPtr
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dictionary & alphaDict(thermalProperties.subDict("alpha"))
IOobject rhoKIO("k", runTime.timeName(0), mesh, IOobject::NO_READ, IOobject::NO_WRITE)
errorManip< error > abort(error &err)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
label k
Boltzmann constant.
const dictionary & CDict(thermalProperties.subDict("C"))
volScalarField & alpha
Fine-structure constant: default SI units: [].
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
IOobject CIO("C", runTime.timeName(0), mesh, IOobject::NO_READ, IOobject::NO_WRITE)
Info<< "Reading thermal properties\n"<< endl;IOdictionary thermalProperties(IOobject("thermalProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));Switch thermalStress(thermalProperties.lookup("thermalStress"));volScalarField threeKalpha(IOobject("threeKalpha", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedScalar("0", dimensionSet(0, 2, -2, -1, 0), 0.0));volScalarField DT(IOobject("DT", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedScalar("0", dimensionSet(0, 2, -1, 0, 0), 0.0));if(thermalStress){ autoPtr< volScalarField > CPtr
const dictionary & kDict(thermalProperties.subDict("k"))