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-2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 "contactHeatFluxSource.H"
29 #include "faMatrices.H"
31 #include "volFields.H"
32 #include "famSup.H"
34 
35 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace fa
40 {
41  defineTypeNameAndDebug(contactHeatFluxSource, 0);
42  addToRunTimeSelectionTable(option, contactHeatFluxSource, dictionary);
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const word& sourceName,
52  const word& modelType,
53  const dictionary& dict,
54  const fvMesh& mesh
55 )
56 :
57  fa::faceSetOption(sourceName, modelType, dict, mesh),
58  TName_(dict.getOrDefault<word>("T", "T")),
59  TprimaryName_(dict.get<word>("Tprimary")),
60  Tprimary_(mesh_.lookupObject<volScalarField>(TprimaryName_)),
61  thicknessLayers_(),
62  kappaLayers_(),
63  contactRes_(0),
64  curTimeIndex_(-1),
65  coupling_()
66 {
67  fieldNames_.resize(1, TName_);
68 
70 
71  read(dict);
72 }
73 
74 
75 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
76 
78 Foam::fa::contactHeatFluxSource::htc() const
79 {
81  (
82  "htc_" + option::name(),
83  regionMesh(),
85  );
86  auto& htc = thtc.ref();
87 
88  // Set up values per coupled patch
89 
90  PtrList<scalarField> patchValues(coupling_.size());
91 
92  forAll(coupling_, patchi)
93  {
94  const auto* tempCoupled = coupling_.get(patchi);
95 
96  if (tempCoupled)
97  {
98  const fvPatch& p = tempCoupled->patch();
99 
100  patchValues.set
101  (
102  patchi,
103  (
104  tempCoupled->kappa(p.patchInternalField(Tprimary_))
105  * p.deltaCoeffs()
106  )
107  );
108  }
109  }
110 
111  vsm().mapToSurface<scalar>(patchValues, htc.field());
112 
113  if (contactRes_ != 0)
114  {
115  htc.field() += contactRes_;
116  }
117 
118  // Zero htc for non-mapped faces
120 
121  return thtc;
122 }
123 
124 
126 (
127  const areaScalarField& h,
128  const areaScalarField& rhoCph,
129  faMatrix<scalar>& eqn,
130  const label fieldi
131 )
132 {
133  if (isActive())
134  {
135  DebugInfo
136  << name() << ": applying source to "
137  << eqn.psi().name() << endl;
138 
139  if (curTimeIndex_ != mesh().time().timeIndex())
140  {
142 
143  // Wall temperature - mapped from primary field to finite-area
145  (
146  "Tw_" + option::name(),
147  regionMesh(),
149  );
150 
151  vsm().mapInternalToSurface<scalar>(Tprimary_, Twall.ref().field());
152 
153  eqn += -fam::Sp(htcw(), eqn.psi()) + htcw()*Twall;
155  curTimeIndex_ = mesh().time().timeIndex();
156  }
157  }
158 }
159 
160 
161 bool Foam::fa::contactHeatFluxSource::read(const dictionary& dict)
162 {
163  if (fa::option::read(dict))
164  {
165  coeffs_.readIfPresent("T", TName_);
166 
167  contactRes_ = 0;
168 
169  if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
170  {
171  dict.readEntry("kappaLayers", kappaLayers_);
172 
173  // Calculate effective thermal resistance by harmonic averaging
174 
175  forAll(thicknessLayers_, iLayer)
176  {
177  contactRes_ += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
178  }
179 
180  if (thicknessLayers_.size())
181  {
182  contactRes_ = scalar(1)/contactRes_;
183  }
184  }
185 
186 
187  // Set up coupling (per-patch) for referenced polyPatches (sorted order)
188  // - size is maxPolyPatch+1
189 
190  const labelList& patches = regionMesh().whichPolyPatches();
191 
192  coupling_.clear();
193  coupling_.resize(patches.empty() ? 0 : (patches.last()+1));
194 
195  for (const label patchi : patches)
196  {
197  const fvPatch& p = mesh_.boundary()[patchi];
198 
199  coupling_.set
200  (
201  patchi,
202  new temperatureCoupling(p, dict)
203  );
204  }
205 
206  return true;
207  }
208 
209  return false;
210 }
211 
212 
213 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
addToRunTimeSelectionTable(option, limitHeight, dictionary)
dictionary dict
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: faOption.H:171
const volSurfaceMapping vsm(aMesh)
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
zeroField Sp(const Foam::zero, const GeometricField< Type, faPatchField, areaMesh > &)
A no-op source.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
defineTypeNameAndDebug(limitHeight, 0)
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: faOptionIO.C:47
const word & name() const noexcept
Return const access to the source name.
Definition: faOptionI.H:23
#define DebugInfo
Report an information message using Foam::Info.
const dimensionSet dimPower
Calculate the finiteArea matrix for implicit and explicit sources.
virtual void addSup(const areaScalarField &h, const areaScalarField &rho, faMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to compressible momentum equation.
const dimensionedScalar h
Planck constant.
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: faOption.C:38
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
contactHeatFluxSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
const polyBoundaryMesh & patches
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...
A special matrix type and solver, designed for finite area solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: faMatricesFwd.H:37
void subsetFilter(List< Type > &field) const
Zero all non-selected locations within field.
const GeometricField< Type, faPatchField, areaMesh > & psi() const
Definition: faMatrix.H:352
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
virtual bool read(const dictionary &dict)
Read source dictionary.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
label timeIndex
Definition: getTimeIndex.H:24
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127