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-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 "externalHeatFluxSource.H"
29 #include "fam.H"
30 #include "faScalarMatrix.H"
34 
36 
37 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace fa
42 {
43  defineTypeNameAndDebug(externalHeatFluxSource, 0);
44  addToRunTimeSelectionTable(option, externalHeatFluxSource, dictionary);
45 }
46 }
47 
48 
49 const Foam::Enum
50 <
52 >
54 ({
55  { operationMode::fixedPower, "power" },
56  { operationMode::fixedHeatFlux, "flux" },
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 )
70 :
71  fa::faceSetOption(sourceName, modelType, dict, m),
72  mode_(operationModeNames.get("mode", dict)),
73  TName_(dict.getOrDefault<word>("T", "T")),
74  Q_(nullptr),
75  q_(nullptr),
76  h_(nullptr),
77  Ta_(nullptr),
78  emissivity_(dict.getOrDefault<scalar>("emissivity", 0))
79 {
80  fieldNames_.resize(1, TName_);
81 
83 
85 }
86 
87 
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 
91 (
92  const areaScalarField& h,
93  const areaScalarField& rho,
94  faMatrix<scalar>& eqn,
95  const label fieldi
96 )
97 {
98  if (isActive())
99  {
100  DebugInfo
101  << name() << ": applying source to "
102  << eqn.psi().name() << endl;
103 
104  scalar qflux = 0;
105 
106  const scalar timeVal = mesh_.time().timeOutputValue();
107 
108  switch (mode_)
109  {
110  case fixedPower:
111  {
112  // From [W] to [W/m2]
113  qflux = Q_->value(timeVal)/(faceSetOption::A() + VSMALL);
114  break;
115  }
116 
117  case fixedHeatFlux:
118  {
119  qflux = q_->value(timeVal);
120  break;
121  }
122 
123  default:
124  {
125  break;
126  }
127  }
128 
129  switch (mode_)
130  {
131  case fixedPower:
132  case fixedHeatFlux:
133  {
135  (
136  "Q",
137  regionMesh(),
139  );
140  auto& Q = tQ.ref();
141 
143  {
144  UIndirectList<scalar>(Q.field(), faceSetOption::faces())
145  = qflux;
146  }
147  else
148  {
149  Q.field() = qflux;
150  }
151 
152  eqn += Q;
153 
154  break;
155  }
156 
157  case fixedHeatTransferCoeff:
158  {
159  const dimensionedScalar Ta
160  (
161  "Ta",
163  Ta_->value(timeVal)
164  );
165 
167  (
168  "h",
169  regionMesh(),
171  (
172  "h",
174  h_->value(timeVal)
175  )
176  );
177  auto& hp = thp.ref();
178 
179  DimensionedField<scalar, areaMesh> hpTa(hp*Ta);
180 
181  if (emissivity_ > 0)
182  {
183  hp -= emissivity_*sigma.value()*pow3(eqn.psi());
184  }
185 
186  // Zero htc for non-mapped faces
187  faceSetOption::subsetFilter(hp.field());
188  faceSetOption::subsetFilter(hpTa.field());
189 
190  eqn -= fam::SuSp(hp, eqn.psi()) - hpTa;
191  break;
192  }
193  }
194  }
195 }
196 
197 
198 bool Foam::fa::externalHeatFluxSource::read(const dictionary& dict)
199 {
200  if (fa::option::read(dict))
201  {
202  dict.readIfPresent("T", TName_);
203  dict.readIfPresent("emissivity", emissivity_);
204 
205  mode_ = operationModeNames.get("mode", dict);
206 
207  switch (mode_)
208  {
209  case fixedPower:
210  {
211  Q_ = Function1<scalar>::New("Q", dict, &mesh_);
212  break;
213  }
214  case fixedHeatFlux:
215  {
216  Q_ = Function1<scalar>::New("q", dict, &mesh_);
217  break;
218  }
219  case fixedHeatTransferCoeff:
220  {
221  h_ = Function1<scalar>::New("h", dict, &mesh_);
222  Ta_ = Function1<scalar>::New("Ta", dict, &mesh_);
223  break;
224  }
225  }
226 
227  return true;
228  }
229 
230  return false;
231 }
232 
233 
234 // ************************************************************************* //
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)
virtual bool read(const dictionary &dict)
Read source dictionary.
dictionary dict
Namespace of functions to calculate implicit derivatives returning a matrix. Time derivatives are cal...
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: faOption.H:171
externalHeatFluxSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
zeroField SuSp(const Foam::zero, const GeometricField< Type, faPatchField, areaMesh > &)
A no-op source.
virtual void addSup(const areaScalarField &h, const areaScalarField &rho, faMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to compressible momentum equation.
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Macros for easy insertion into run-time selection tables.
scalar A() const noexcept
Return const access to the total face area.
bool useSubMesh() const noexcept
True if sub-selection should be used.
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 labelList & faces() const noexcept
Return const access to the local finite-area face selection.
#define DebugInfo
Report an information message using Foam::Info.
const dimensionSet dimPower
operationMode
Options for the heat transfer condition mode.
static const Enum< operationMode > operationModeNames
Names for operationMode.
const dimensionedScalar h
Planck constant.
dimensionedScalar pow3(const dimensionedScalar &ds)
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: faOption.C:38
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
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
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127