43namespace thermalBaffleModels
104 const polyPatch& ppCoupled = rbm[patchi];
106 forAll(ppCoupled, localFacei)
160 const word& modelType,
166 nNonOrthCorr_(
solution().get<label>(
"nNonOrthCorr")),
199 dict.subDict(
"radiation"),
211 const word& modelType,
216 nNonOrthCorr_(
solution().get<label>(
"nNonOrthCorr")),
266void thermalBaffle::init()
276 <<
"the boundary field of qs is "
277 << qsb <<
" and " <<
nl
340 const label patchi = coupledPatches[i];
348 *
thermo_->alpha().boundaryField()[patchi]
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
void size(const label n)
Older name for setAddressableSize.
label size() const noexcept
The number of elements in the list.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
const Type & value() const
Return const reference to value.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
Mesh data needed to do the Finite Volume discretisation.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
A patch is a list of labels that address the faces in the global face list.
Top level model for radiation modelling.
Switch moveMesh_
Flag to allow mesh movement.
labelListList boundaryFaceCells_
Global cell IDs.
virtual bool read()
Read control parameters from dictionary.
const Time & time() const
Return the reference to the time database.
const dictionary & solution() const
Return the solution dictionary.
const fvMesh & regionMesh() const
Return the region mesh database.
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
bool constantThickness_
Is thickness constant.
scalarField thickness_
Baffle physical thickness.
bool oneD_
Is it one dimension.
dimensionedScalar delta_
Baffle mesh thickness.
virtual const solidThermo & thermo() const
Return const reference to the solidThermo.
virtual const volScalarField & rho() const
Return density [Kg/m3].
volScalarField qs_
Surface energy source / [J/m2/s].
volScalarField & h_
Enthalpy/internal energy.
virtual const volScalarField & kappa() const
Return thermal conductivity [W/m/K].
virtual const volScalarField & T() const
Return temperature [K].
void solveEnergy()
Solve energy equation.
virtual const tmp< volScalarField > Cp() const
Return the film specific heat capacity [J/kg/K].
volScalarField Q_
Volumetric energy source / [J/m3/s].
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
autoPtr< radiation::radiationModel > radiation_
Pointer to radiation model.
virtual const volScalarField & kappaRad() const
Return solid absortivity [1/m].
autoPtr< solidThermo > thermo_
Solid thermo.
virtual void info()
Provide some feedback.
virtual ~thermalBaffle()
Destructor.
label nNonOrthCorr_
Number of non orthogonal correctors.
virtual bool read()
Read control parameters IO dictionary.
virtual void evolveRegion()
Evolve the thermal baffle.
Fundamental solid thermodynamic properties.
Selector class for relaxation factors, solver type and solution.
A class for managing temporary objects.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Calculate the divergence of the given field.
#define DebugInFunction
Report an information message using Foam::Info.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Type gSum(const FieldField< Field, Type > &f)
const dimensionSet dimEnergy
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
const dimensionSet dimArea(sqr(dimLength))
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Ostream & indent(Ostream &os)
Indent stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
static constexpr const zero Zero
Global zero (0)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimVolume(pow3(dimLength))
errorManipArg< error, int > exit(error &err, const int errNo=1)
constexpr char nl
The newline '\n' character (0x0a)
#define forAll(list, i)
Loop across all elements in list.