Loading...
Searching...
No Matches
contactHeatFluxSource.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 "faMatrices.H"
31#include "volFields.H"
32#include "famSup.H"
34// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace fa
39{
42}
43}
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
49(
50 const word& sourceName,
51 const word& modelType,
52 const dictionary& dict,
53 const fvMesh& mesh,
54 const word& defaultAreaName
55)
56:
57 fa::faceSetOption(sourceName, modelType, dict, mesh, defaultAreaName),
58 TName_
59 (
60 dict.getOrDefaultCompat<word>
61 (
62 "Ts", {{"T", -2506}}, suffixed("Ts"), keyType::LITERAL
63 )
64 ),
65 TprimaryName_(dict.getOrDefault<word>("Tprimary", "T")),
66 Tprimary_(mesh_.lookupObject<volScalarField>(TprimaryName_)),
67 contactRes_(0),
68 curTimeIndex_(-1)
69{
70 fieldNames_.resize(1, TName_);
71
73
74 read(dict);
75}
76
77
78// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79
80Foam::tmp<Foam::DimensionedField<Foam::scalar, Foam::areaMesh>>
81Foam::fa::contactHeatFluxSource::htc() const
82{
84 (
85 "htc_" + option::name(),
86 regionMesh(),
88 );
89 auto& htc = thtc.ref();
90
91 // Set up values per coupled patch
92
93 PtrList<scalarField> patchValues(coupling_.size());
94
95 forAll(coupling_, patchi)
96 {
97 const auto* tempCoupled = coupling_.get(patchi);
98
99 if (tempCoupled)
100 {
101 const fvPatch& p = tempCoupled->patch();
102
103 patchValues.set
104 (
105 patchi,
106 (
107 tempCoupled->kappa(p.patchInternalField(Tprimary_))
108 * p.deltaCoeffs()
109 )
110 );
111 }
112 }
113
114 vsm().mapToSurface<scalar>(patchValues, htc.field());
115
116 if (contactRes_ != 0)
117 {
118 htc.field() += contactRes_;
119 }
120
121 // Zero htc for non-mapped faces
122 faceSetOption::subsetFilter(htc.field());
123
124 return thtc;
125}
126
127
129(
130 const areaScalarField& h,
131 const areaScalarField& rhoCph,
132 faMatrix<scalar>& eqn,
133 const label fieldi
134)
135{
136 if (!isActive())
137 {
138 return;
139 }
140
141 DebugInfo<< name() << ": applying source to " << eqn.psi().name() << endl;
142
143 if (curTimeIndex_ != mesh().time().timeIndex())
144 {
146
147 // Wall temperature - mapped from primary field to finite-area
149 (
150 "Tw_" + option::name(),
151 regionMesh(),
153 );
154
155 vsm().mapInternalToSurface<scalar>(Tprimary_, Twall.ref().field());
156
157 eqn += -fam::Sp(htcw(), eqn.psi()) + htcw()*Twall;
158
159 curTimeIndex_ = mesh().time().timeIndex();
160 }
161}
162
163
165{
167 {
168 coeffs_.readIfPresent("T", TName_);
169
170 contactRes_ = 0;
171
172 if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
173 {
174 dict.readEntry("kappaLayers", kappaLayers_);
175
176 // Calculate effective thermal resistance by harmonic averaging
177
178 forAll(thicknessLayers_, iLayer)
179 {
180 contactRes_ += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
181 }
182
183 if (thicknessLayers_.size())
184 {
185 contactRes_ = scalar(1)/contactRes_;
186 }
187 }
188
189
190 // Set up coupling (per-patch) for referenced polyPatches (sorted order)
191 // - size is maxPolyPatch+1
192
193 const labelList& patches = regionMesh().whichPolyPatches();
194
195 if (patches.empty())
196 {
197 coupling_.clear();
198 }
199 else
200 {
201 coupling_.resize_null(patches.back()+1);
202 }
203
204 for (const label patchi : patches)
205 {
206 const fvPatch& p = mesh_.boundary()[patchi];
207
208 coupling_.set
209 (
210 patchi,
211 new temperatureCoupling(p, dict)
212 );
213 }
214
215 return true;
216 }
217
218 return false;
219}
220
221
222// ************************************************************************* //
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< 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....
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
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 contact heat flux between specified faMesh and fvMesh within a specified region for compressi...
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.
contactHeatFluxSource(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.
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.
Base abstract class for handling finite area options (i.e. faOption).
Definition faOption.H:149
const fvMesh & mesh() const noexcept
Return const access to the volume mesh.
Definition faOption.H:385
virtual bool read(const dictionary &dict)
Read source dictionary.
const fvMesh & mesh_
Reference to the mesh database.
Definition faOption.H:175
dictionary coeffs_
Dictionary containing source coefficients.
Definition faOption.H:185
const word & name() const noexcept
The source name.
Definition faOption.H:380
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
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
@ LITERAL
String literal.
Definition keyType.H:82
A class for managing temporary objects.
Definition tmp.H:75
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
volScalarField & p
rDeltaTY field()
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
const volSurfaceMapping vsm(aMesh)
Calculate the finiteArea matrix for implicit and explicit sources.
auto & name
#define DebugInfo
Report an information message using Foam::Info.
Namespace for finite-area.
Definition limitHeight.C:30
zeroField Sp(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
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
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.
label timeIndex
dictionary dict
volScalarField & h
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299