Loading...
Searching...
No Matches
externalHeatFluxSource.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
29#include "fam.H"
30#include "faScalarMatrix.H"
34
37// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
38
39namespace Foam
40{
41namespace fa
42{
45}
46}
47
48
49const Foam::Enum
50<
52>
54({
55 { operationMode::fixedPower, "power" },
57 { operationMode::fixedHeatTransferCoeff, "coefficient" },
58});
59
60
61// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
62
64(
65 const word& sourceName,
66 const word& modelType,
67 const dictionary& dict,
68 const fvMesh& m,
69 const word& defaultAreaName
70)
71:
72 fa::faceSetOption(sourceName, modelType, dict, m, defaultAreaName),
73 mode_(operationModeNames.get("mode", dict)),
74 TName_
75 (
76 dict.getOrDefaultCompat<word>
77 (
78 "Ts", {{"T", -2506}},
79 suffixed("Ts"), keyType::LITERAL
80 )
81 ),
82 emissivity_(dict.getOrDefault<scalar>("emissivity", 0))
83{
84 fieldNames_.resize(1, TName_);
85
88 read(dict);
89}
90
91
92// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93
95(
96 const areaScalarField& h,
97 const areaScalarField& rho,
99 const label fieldi
100)
101{
102 if (isActive())
103 {
105 << name() << ": applying source to "
106 << eqn.psi().name() << endl;
107
108 scalar qflux = 0;
109
110 const scalar timeVal = mesh_.time().timeOutputValue();
111
112 switch (mode_)
113 {
114 case fixedPower:
115 {
116 // From [W] to [W/m2]
117 qflux = Q_->value(timeVal)/(faceSetOption::A() + VSMALL);
118 break;
119 }
120
121 case fixedHeatFlux:
122 {
123 qflux = q_->value(timeVal);
124 break;
125 }
126
127 default:
128 {
129 break;
130 }
131 }
132
133 switch (mode_)
134 {
135 case fixedPower:
136 case fixedHeatFlux:
137 {
139 (
140 "Q",
141 regionMesh(),
143 );
144 auto& Q = tQ.ref();
145
147 {
148 UIndirectList<scalar>(Q.field(), faceSetOption::faces())
149 = qflux;
150 }
151 else
152 {
153 Q.field() = qflux;
154 }
155
156 eqn += Q;
157
158 break;
159 }
160
161 case fixedHeatTransferCoeff:
162 {
163 const dimensionedScalar Ta
164 (
165 "Ta",
167 Ta_->value(timeVal)
168 );
169
171 (
172 "h",
173 regionMesh(),
175 (
176 "h",
178 h_->value(timeVal)
179 )
180 );
181 auto& hp = thp.ref();
182
183 DimensionedField<scalar, areaMesh> hpTa(hp*Ta);
184
185 if (emissivity_ > 0)
186 {
187 hp -= emissivity_*sigma.value()*pow3(eqn.psi());
188 }
189
190 // Zero htc for non-mapped faces
191 faceSetOption::subsetFilter(hp.field());
192 faceSetOption::subsetFilter(hpTa.field());
193
194 eqn -= fam::SuSp(hp, eqn.psi()) - hpTa;
195 break;
196 }
197 }
198 }
199}
200
201
203{
205 {
206 dict.readIfPresent("T", TName_);
207 dict.readIfPresent("emissivity", emissivity_);
208
209 mode_ = operationModeNames.get("mode", dict);
210
211 Q_.reset(nullptr);
212 q_.reset(nullptr);
213 h_.reset(nullptr);
214 Ta_.reset(nullptr);
215
216 switch (mode_)
217 {
218 case fixedPower:
219 {
220 Q_ = Function1<scalar>::New("Q", dict, &mesh_);
221 break;
222 }
223 case fixedHeatFlux:
224 {
225 Q_ = Function1<scalar>::New("q", dict, &mesh_);
226 break;
227 }
228 case fixedHeatTransferCoeff:
229 {
230 h_ = Function1<scalar>::New("h", dict, &mesh_);
231 Ta_ = Function1<scalar>::New("Ta", dict, &mesh_);
232 break;
233 }
234 }
235
236 return true;
237 }
238
239 return false;
240}
241
242
243// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const DynamicField< Type > & field() const noexcept
Return const-reference to the primitive field values.
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const Field< Type > &iField)
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions, copy of internal field....
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A List with indirect addressing. Like IndirectList but does not store addressing.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A special matrix type and solver, designed for finite area solutions of scalar equations....
Definition faMatrix.H:108
const GeometricField< Type, faPatchField, areaMesh > & psi() const
Definition faMatrix.H:348
Applies a heat flux condition for a specified faMesh region to temperature on an external wall for co...
static const Enum< operationMode > operationModeNames
Names for operationMode.
virtual void addSup(const areaScalarField &h, const areaScalarField &rho, faMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to compressible momentum equation.
virtual bool read(const dictionary &dict)
Read source dictionary.
operationMode
Options for the heat transfer condition mode.
@ fixedHeatTransferCoeff
Fixed heat transfer coefficient.
externalHeatFluxSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh, const word &defaultAreaName=word())
Construct from explicit source name and mesh.
void subsetFilter(List< Type > &field) const
Zero all non-selected locations within field.
bool useSubMesh() const noexcept
True if sub-selection should be used.
virtual bool isActive()
Is the source active?
faceSetOption(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh, const word &defaultAreaName=word())
Construct from components.
const labelList & faces() const noexcept
Return const access to the local finite-area face selection.
scalar A() const noexcept
Return const access to the total face area.
Base abstract class for handling finite area options (i.e. faOption).
Definition faOption.H:149
virtual bool read(const dictionary &dict)
Read source dictionary.
const fvMesh & mesh_
Reference to the mesh database.
Definition faOption.H:175
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
const faMesh & regionMesh() const
Return the region mesh database (demand-driven).
Definition faOptionI.H:37
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
@ LITERAL
String literal.
Definition keyType.H:82
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
Namespace of functions to calculate implicit derivatives returning a matrix. Time derivatives are cal...
auto & name
#define DebugInfo
Report an information message using Foam::Info.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Namespace for finite-area.
Definition limitHeight.C:30
zeroField SuSp(const Foam::zero, const GeometricField< Type, faPatchField, areaMesh > &)
A no-op source.
Namespace for OpenFOAM.
const dimensionSet dimPower
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
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
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
volScalarField & h
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)