Go to the documentation of this file.
41 namespace regionModels
43 namespace thermalBaffleModels
104 const polyPatch& ppCoupled = rbm[patchi];
106 forAll(ppCoupled, localFacei)
158 thermalBaffle::thermalBaffle
160 const word& modelType,
166 nNonOrthCorr_(
solution().get<label>(
"nNonOrthCorr")),
199 dict.subDict(
"radiation"),
209 thermalBaffle::thermalBaffle
211 const word& modelType,
216 nNonOrthCorr_(
solution().get<label>(
"nNonOrthCorr")),
266 void thermalBaffle::init()
276 <<
"the boundary field of qs is "
277 << qsb <<
" and " <<
nl
278 <<
"the field 'thickness' is " <<
thickness_.size() <<
nl
340 const label patchi = coupledPatches[i];
348 *
thermo_->alpha().boundaryField()[patchi]
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.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Selector class for relaxation factors, solver type and solution.
A class for handling words, derived from Foam::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.
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
dimensionedScalar delta_
Baffle mesh thickness.
static constexpr const zero Zero
Global zero (0)
const dimensionSet dimEnergy
bool constantThickness_
Is thickness constant.
volScalarField qs_
Surface energy source / [J/m2/s].
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
volScalarField & h_
Enthalpy/internal energy.
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)
void solveEnergy()
Solve energy equation.
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)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
virtual const solidThermo & thermo() const
Return const reference to the solidThermo.
scalarField thickness_
Baffle physical thickness.
#define forAll(list, i)
Loop across all elements in list.
labelListList boundaryFaceCells_
Global cell IDs.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Fundamental solid thermodynamic properties.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
const dimensionSet dimArea(sqr(dimLength))
static autoPtr< radiationModel > New(const volScalarField &T)
Return a reference to the selected radiation model.
messageStream Info
Information stream (uses stdout - output is on the master only)
#define DebugInFunction
Report an information message using Foam::Info.
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).
const fvMesh & regionMesh() const
Return the region mesh database.
addToRunTimeSelectionTable(thermalBaffleModel, noThermo, mesh)
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
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.
const dictionary & solution() const
Return the solution dictionary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Macros for easy insertion into run-time selection tables.
Mesh data needed to do the Finite Volume discretisation.
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.
virtual void evolveRegion()
Evolve the thermal baffle.
Ostream & indent(Ostream &os)
Indent stream.
autoPtr< radiation::radiationModel > radiation_
Pointer to radiation model.
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.
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].
virtual bool read()
Read control parameters IO dictionary.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A special matrix type and solver, designed for finite volume solutions of scalar equations....
SolverPerformance< Type > solve(const dictionary &)
Solve returning the solution statistics.
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))
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
label nNonOrthCorr_
Number of non orthogonal correctors.