Loading...
Searching...
No Matches
thermalShell.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2019-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "thermalShell.H"
29#include "fam.H"
30#include "faMatrices.H"
31#include "fvPatchFields.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace regionModels
39{
40
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
44
46
47// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48
49bool thermalShell::init(const dictionary& dict)
50{
51 if (thickness_ > 0)
52 {
54 }
55
56 this->solution().readEntry("nNonOrthCorr", nNonOrthCorr_);
57
58 return true;
59}
60
61
62tmp<areaScalarField> thermalShell::qr()
63{
64 auto taqr =
66 (
68 (
69 "tqr",
71 regionMesh().thisDb()
72 ),
73 regionMesh(),
75 );
76
77 if (!qrName_.empty() && qrName_ != "none")
78 {
79 vsm().mapToSurface<scalar>
80 (
82 taqr.ref().primitiveFieldRef()
83 );
84 }
86 return taqr;
87}
88
89
90// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
91
93{
95
96 const areaScalarField rhoCph(Cp()*rho()*h_);
97
99 (
100 fam::ddt(rhoCph, T_)
102 ==
103 qs_
104 + qr()
105 + faOptions()(h_, rhoCph, T_)
106 );
107
108 TEqn.relax();
109
111
112 TEqn.solve();
115}
116
117
118// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
119
121(
122 const word& modelType,
123 const fvMesh& mesh,
124 const dictionary& dict
125)
126:
127 thermalShellModel(modelType, mesh, dict),
128 hName_(dict.getOrDefault<word>("h", suffixed("h"))),
129 qsName_(dict.getOrDefault<word>("qs", suffixed("qs"))),
130 qrName_(dict.getOrDefault<word>("qr", "none")),
131 nNonOrthCorr_(1),
132 // Only need/want thermal solid properties
133 thermo_(dict.subDict("thermo"), solidProperties::THERMAL),
134 h_
135 (
137 (
138 hName_,
139 regionMesh().time().timeName(),
140 regionMesh().thisDb(),
141 IOobject::MUST_READ,
142 IOobject::AUTO_WRITE
143 ),
144 regionMesh()
145 ),
146 qs_
147 (
149 (
150 qsName_,
151 regionMesh().time().timeName(),
152 regionMesh().thisDb(),
153 IOobject::READ_IF_PRESENT,
154 IOobject::NO_WRITE
155 ),
156 regionMesh(),
158 ),
159 thickness_(dict.getOrDefault<scalar>("thickness", 0))
161 init(dict);
162}
163
165// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
166
168{}
169
170
172{
173 nNonOrthCorr_ = solution().get<label>("nNonOrthCorr");
174
175 for (int nonOrth = 0; nonOrth <= nNonOrthCorr_; ++nonOrth)
176 {
178 }
179
180 Info<< T_.name() << " min/max = " << gMinMax(T_) << endl;
181}
182
183
185{
187 (
188 "Cps",
201 "rhos",
218 (
220 thermo_.kappa()
221 ),
223 );
224}
225
226
228{}
229
230
231// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232
233} // End namespace regionModels
234} // End namespace Foam
235
236// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< scalar, faPatchField, areaMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=faPatchField< scalar >::calculatedType())
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
@ NO_REGISTER
Do not request registration (bool: false).
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
void correct(GeometricField< Type, faPatchField, areaMesh > &field)
Apply correction to field.
void constrain(faMatrix< Type > &eqn)
Apply constraints to equation.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
const dictionary & solution() const
Return the solution dictionary.
word suffixed(const char *base) const
Return the concatenation of base and the suffix hint.
const Time & time() const noexcept
Return the reference to the time database.
const volSurfaceMapping & vsm() const
Return mapping between surface and volume fields.
const faMesh & regionMesh() const
Return the region mesh database.
const fvMesh & primaryMesh() const noexcept
Return the reference to the primary mesh database.
Intermediate class for thermal-shell finite-area models.
thermalShellModel(const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from type name and mesh and dict.
areaScalarField T_
Shell temperature.
Foam::fa::options & faOptions() noexcept
Return faOptions.
thermalShell(const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components and dict.
areaScalarField h_
Shell thickness field [m].
const tmp< areaScalarField > rho() const
Return density [kg/m3].
scalar thickness_
Uniform shell thickness [m].
const word qsName_
Name of surface energy source (default: "qs" + suffix).
void solveEnergy()
Solve energy equation.
areaScalarField qs_
External surface energy source [J/m2/s].
const word qrName_
Name of (volume) radiative flux field (default: none).
virtual void preEvolveRegion()
Pre-evolve thermal baffle.
solidProperties thermo_
Solid properties.
const tmp< areaScalarField > Cp() const
Return the shell specific heat capacity [J/kg/K].
const word hName_
Name of shell thickness [height] field (default: "hs" + suffix).
virtual void info()
Provide some feedback.
label nNonOrthCorr_
Number of non orthogonal correctors.
const tmp< areaScalarField > kappa() const
Return thermal conductivity [W/m/K].
virtual void evolveRegion()
Evolve the thermal baffle.
The thermophysical, mechanical properties of a solid.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
void mapToSurface(const GeometricBoundaryField< Type, fvPatchField, volMesh > &, Field< Type > &result) const
Map volume boundary fields as area field.
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
Namespace of functions to calculate implicit derivatives returning a matrix. Time derivatives are cal...
fvScalarMatrix TEqn(fvm::ddt(T)+fvm::div(phi, T) - fvm::laplacian(alphaEff, T)==radiation->ST(rhoCpRef, T)+fvOptions(T))
word timeName
Definition getTimeIndex.H:3
#define DebugInFunction
Report an information message using Foam::Info.
tmp< faMatrix< Type > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf)
tmp< faMatrix< Type > > ddt(const GeometricField< Type, faPatchField, areaMesh > &vf)
Definition famDdt.C:41
Namespace for OpenFOAM.
const dimensionSet dimPower
const dimensionSet dimEnergy
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
faMatrix< scalar > faScalarMatrix
Definition faMatrices.H:46
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
dictionary dict