Go to the documentation of this file.
39 namespace regionModels
41 namespace thermalBaffleModels
71 Info<<
"thermalBaffle::solveEnergy()" <<
endl;
105 const polyPatch& ppCoupled = rbm[patchI];
107 forAll(ppCoupled, localFaceI)
161 const word& modelType,
210 dict.subDict(
"radiation"),
222 const word& modelType,
297 <<
"the boundary field of Qs is "
298 << Qsb <<
" and " <<
nl
299 <<
"the field 'thickness' is " <<
thickness_.size() <<
nl
361 const label patchI = coupledPatches[i];
369 *
thermo_->alpha().boundaryField()[patchI]
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
bool oneD_
Is it one dimension.
virtual bool read()
Read control parameters from dictionary.
virtual const volScalarField & T() const
Return temperature [K].
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Selector class for relaxation factors, solver type and solution.
A class for handling words, derived from string.
virtual const volScalarField & kappa() const
Return thermal conductivity [W/m/K].
static autoPtr< solidThermo > New(const fvMesh &, const word &phaseName=word::null)
Return a pointer to a new solidThermo created from.
dimensionedScalar delta_
Baffle mesh thickness.
#define forAll(list, i)
Loop across all elements in list.
A class for managing temporary objects.
const dimensionSet dimEnergy
bool constantThickness_
Is thickness constant.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & h_
Enthalpy/internal energy.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Calculate the divergence of the given field.
virtual ~thermalBaffle()
Destructor.
const Time & time() const
Return the reference to the time database.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
void solveEnergy()
Solve energy equation.
volScalarField Qs_
Surface energy source / [J/m2/s].
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type & value() const
Return const reference to value.
volScalarField Q_
Volumetric energy source / [J/m3/s].
Type gSum(const FieldField< Field, Type > &f)
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual const solidThermo & thermo() const
Return const reference to the solidThermo.
scalarField thickness_
Baffle physical thickness.
labelListList boundaryFaceCells_
Global cell IDs.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
This function object calculates and outputs the second invariant of the velocity gradient tensor [1/s...
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Fundamental solid thermodynamic properties.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
const dimensionSet dimArea(sqr(dimLength))
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
static autoPtr< radiationModel > New(const volScalarField &T)
Return a reference to the selected radiation model.
A patch is a list of labels that address the faces in the global face list.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
addToRunTimeSelectionTable(thermalBaffleModel, noThermo, mesh)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dictionary & solution() const
Return the solution dictionary.
thermalBaffle(const thermalBaffle &)
Disallow default bitwise copy construct.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Macros for easy insertion into run-time selection tables.
cellIbMask *cellIbMaskExt *faceIbMask *cellIbMask boundaryField().evaluateCoupled()
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
virtual void evolveRegion()
Evolve the thermal baffle.
Ostream & indent(Ostream &os)
Indent stream.
autoPtr< radiation::radiationModel > radiation_
Pointer to radiation model.
void init()
Initialize thermalBaffle.
errorManipArg< error, int > exit(error &err, const int errNo=1)
virtual const volScalarField & kappaRad() const
Return solid absortivity [1/m].
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
defineTypeNameAndDebug(noThermo, 0)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
virtual const volScalarField & rho() const
Return density [Kg/m3].
virtual void info()
Provide some feedback.
virtual const tmp< volScalarField > Cp() const
Return the film specific heat capacity [J/kg/K].
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Traits class for primitives.
virtual bool read()
Read control parameters IO dictionary.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
label readLabel(Istream &is)
Switch moveMesh_
Flag to allow mesh movement.
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
autoPtr< solidThermo > thermo_
Solid thermo.
const dimensionSet dimVolume(pow3(dimLength))
Generic GeometricField class.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
label nNonOrthCorr_
Number of non orthogonal correctors.
stressControl lookup("compactNormalStress") >> compactNormalStress
const fvMesh & regionMesh() const
Return the region mesh database.