filmPyrolysisRadiativeCoupledMixedFvPatchScalarField.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) 2013-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "mappedPatchBase.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
41 filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::
42 filmModel() const
43 {
44  HashTable<const filmModelType*> models
45  = db().time().lookupClass<filmModelType>();
46 
47  forAllConstIters(models, iter)
48  {
49  if (iter()->regionMesh().name() == filmRegionName_)
50  {
51  return *iter();
52  }
53  }
54 
55  DynamicList<word> modelNames;
56  forAllConstIters(models, iter)
57  {
58  modelNames.append(iter()->regionMesh().name());
59  }
60 
62  << "Unable to locate film region " << filmRegionName_
63  << ". Available regions include: " << modelNames
64  << abort(FatalError);
65 
66  return **models.begin();
67 }
68 
69 
72 filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::
73 pyrModel() const
74 {
75  HashTable<const pyrolysisModelType*> models =
77 
78  forAllConstIters(models, iter)
79  {
80  if (iter()->regionMesh().name() == pyrolysisRegionName_)
81  {
82  return *iter();
83  }
84  }
85 
86  DynamicList<word> modelNames;
87  forAllConstIters(models, iter)
88  {
89  modelNames.append(iter()->regionMesh().name());
90  }
91 
93  << "Unable to locate pyrolysis region " << pyrolysisRegionName_
94  << ". Available regions include: " << modelNames
95  << abort(FatalError);
96 
97  return **models.begin();
98 }
99 
100 
101 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
102 
105 (
106  const fvPatch& p,
108 )
109 :
110  mixedFvPatchScalarField(p, iF),
112  (
113  patch(),
114  "undefined",
115  "undefined",
116  "undefined-K",
117  "undefined-alpha"
118  ),
119  filmRegionName_("surfaceFilmProperties"),
120  pyrolysisRegionName_("pyrolysisProperties"),
121  TnbrName_("undefined-Tnbr"),
122  qrName_("undefined-qr"),
123  convectiveScaling_(1.0),
124  filmDeltaDry_(0.0),
125  filmDeltaWet_(0.0)
126 {
127  this->refValue() = 0.0;
128  this->refGrad() = 0.0;
129  this->valueFraction() = 1.0;
130 }
131 
132 
135 (
137  const fvPatch& p,
139  const fvPatchFieldMapper& mapper
140 )
141 :
142  mixedFvPatchScalarField(psf, p, iF, mapper),
144  filmRegionName_(psf.filmRegionName_),
145  pyrolysisRegionName_(psf.pyrolysisRegionName_),
146  TnbrName_(psf.TnbrName_),
147  qrName_(psf.qrName_),
148  convectiveScaling_(psf.convectiveScaling_),
149  filmDeltaDry_(psf.filmDeltaDry_),
150  filmDeltaWet_(psf.filmDeltaWet_)
151 {}
152 
153 
156 (
157  const fvPatch& p,
159  const dictionary& dict
160 )
161 :
162  mixedFvPatchScalarField(p, iF),
164  filmRegionName_
165  (
166  dict.getOrDefault<word>("filmRegion", "surfaceFilmProperties")
167  ),
168  pyrolysisRegionName_
169  (
170  dict.getOrDefault<word>("pyrolysisRegion", "pyrolysisProperties")
171  ),
172  TnbrName_(dict.lookup("Tnbr")),
173  qrName_(dict.lookup("qr")),
174  convectiveScaling_(dict.getOrDefault<scalar>("convectiveScaling", 1)),
175  filmDeltaDry_(dict.get<scalar>("filmDeltaDry")),
176  filmDeltaWet_(dict.get<scalar>("filmDeltaWet"))
177 {
178  if (!isA<mappedPatchBase>(this->patch().patch()))
179  {
181  << "' not type '" << mappedPatchBase::typeName << "'"
182  << "\n for patch " << p.name()
183  << " of field " << internalField().name()
184  << " in file " << internalField().objectPath()
185  << exit(FatalError);
186  }
187 
188  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
189 
190  if (dict.found("refValue"))
191  {
192  // Full restart
193  refValue() = scalarField("refValue", dict, p.size());
194  refGrad() = scalarField("refGradient", dict, p.size());
195  valueFraction() = scalarField("valueFraction", dict, p.size());
196  }
197  else
198  {
199  // Start from user entered data. Assume fixedValue.
200  refValue() = *this;
201  refGrad() = 0.0;
202  valueFraction() = 1.0;
203  }
204 }
205 
206 
209 (
212 )
213 :
214  mixedFvPatchScalarField(psf, iF),
216  filmRegionName_(psf.filmRegionName_),
217  pyrolysisRegionName_(psf.pyrolysisRegionName_),
218  TnbrName_(psf.TnbrName_),
219  qrName_(psf.qrName_),
220  convectiveScaling_(psf.convectiveScaling_),
221  filmDeltaDry_(psf.filmDeltaDry_),
222  filmDeltaWet_(psf.filmDeltaWet_)
223 {}
224 
225 
226 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
227 
229 (
230  const fvPatchFieldMapper& mapper
231 )
232 {
233  mixedFvPatchScalarField::autoMap(mapper);
235 }
236 
237 
239 (
240  const fvPatchField<scalar>& ptf,
241  const labelList& addr
242 )
243 {
244  mixedFvPatchScalarField::rmap(ptf, addr);
245 
247  refCast
248  <
250  >(ptf);
251 
252  temperatureCoupledBase::rmap(tiptf, addr);
253 }
254 
255 
257 {
258  if (updated())
259  {
260  return;
261  }
262 
263  // Get the coupling information from the mappedPatchBase
264  const mappedPatchBase& mpp =
265  refCast<const mappedPatchBase>(patch().patch());
266 
267  const label patchi = patch().index();
268  const label nbrPatchi = mpp.samplePolyPatch().index();
269  const polyMesh& mesh = patch().boundaryMesh().mesh();
270  const polyMesh& nbrMesh = mpp.sampleMesh();
271  const fvPatch& nbrPatch =
272  refCast<const fvMesh>(nbrMesh).boundary()[nbrPatchi];
273 
274  scalarField intFld(patchInternalField());
275 
277  nbrField =
278  refCast
279  <
281  >
282  (
283  nbrPatch.lookupPatchField<volScalarField, scalar>(TnbrName_)
284  );
285 
286  // Swap to obtain full local values of neighbour internal field
287  scalarField nbrIntFld(nbrField.patchInternalField());
288  mpp.distribute(nbrIntFld);
289 
290  scalarField& Tp = *this;
291 
292  const scalarField K(this->kappa(*this));
293  const scalarField nbrK(nbrField.kappa(*this));
294 
295  // Swap to obtain full local values of neighbour K*delta
296  scalarField KDeltaNbr(nbrK*nbrPatch.deltaCoeffs());
297  mpp.distribute(KDeltaNbr);
298 
299  scalarField myKDelta(K*patch().deltaCoeffs());
300 
301  scalarList Tfilm(patch().size(), Zero);
302  scalarList htcwfilm(patch().size(), Zero);
303  scalarList filmDelta(patch().size(), Zero);
304 
305  const pyrolysisModelType& pyrolysis = pyrModel();
306  const filmModelType& film = filmModel();
307 
308  // Obtain Rad heat (qr)
309  scalarField qr(patch().size(), Zero);
310 
311  label coupledPatchi = -1;
312  if (pyrolysisRegionName_ == mesh.name())
313  {
314  coupledPatchi = patchi;
315  if (qrName_ != "none")
316  {
317  qr = nbrPatch.lookupPatchField<volScalarField, scalar>(qrName_);
318  mpp.distribute(qr);
319  }
320  }
321  else if (pyrolysis.primaryMesh().name() == mesh.name())
322  {
323  coupledPatchi = nbrPatch.index();
324  if (qrName_ != "none")
325  {
326  qr = patch().lookupPatchField<volScalarField, scalar>(qrName_);
327  }
328  }
329  else
330  {
332  << type() << " condition is intended to be applied to either the "
333  << "primary or pyrolysis regions only"
334  << exit(FatalError);
335  }
336 
337  const label filmPatchi = pyrolysis.nbrCoupledPatchID(film, coupledPatchi);
338 
339  const scalarField htcw(film.htcw().h()().boundaryField()[filmPatchi]);
340 
341  // Obtain htcw
342  htcwfilm =
343  pyrolysis.mapRegionPatchField
344  (
345  film,
346  coupledPatchi,
347  filmPatchi,
348  htcw,
349  true
350  );
351 
352 
353  // Obtain Tfilm at the boundary through Ts.
354  // NOTE: Tf is not good as at the boundary it will retrieve Tp
355  Tfilm = film.Ts().boundaryField()[filmPatchi];
356  film.toPrimary(filmPatchi, Tfilm);
357 
358  // Obtain delta
359  filmDelta =
360  pyrolysis.mapRegionPatchField<scalar>
361  (
362  film,
363  "deltaf",
364  coupledPatchi,
365  true
366  );
367 
368  // Estimate wetness of the film (1: wet , 0: dry)
369  scalarField ratio
370  (
371  min
372  (
373  max
374  (
375  (filmDelta - filmDeltaDry_)/(filmDeltaWet_ - filmDeltaDry_),
376  scalar(0)
377  ),
378  scalar(1)
379  )
380  );
381 
382  scalarField qConv(ratio*htcwfilm*(Tfilm - Tp)*convectiveScaling_);
383 
384  scalarField qRad((1.0 - ratio)*qr);
385 
386  scalarField alpha(KDeltaNbr - (qRad + qConv)/Tp);
387 
388  valueFraction() = alpha/(alpha + (1.0 - ratio)*myKDelta);
389 
390  refValue() = ratio*Tfilm + (1.0 - ratio)*(KDeltaNbr*nbrIntFld)/alpha;
391 
392  mixedFvPatchScalarField::updateCoeffs();
393 
394  if (debug)
395  {
396  scalar Qc = gSum(qConv*patch().magSf());
397  scalar qr = gSum(qRad*patch().magSf());
398  scalar Qt = gSum((qConv + qRad)*patch().magSf());
399 
400  Info<< mesh.name() << ':'
401  << patch().name() << ':'
402  << this->internalField().name() << " <- "
403  << nbrMesh.name() << ':'
404  << nbrPatch.name() << ':'
405  << this->internalField().name() << " :" << nl
406  << " convective heat[W] : " << Qc << nl
407  << " radiative heat [W] : " << qr << nl
408  << " total heat [W] : " << Qt << nl
409  << " wall temperature "
410  << " min:" << gMin(*this)
411  << " max:" << gMax(*this)
412  << " avg:" << gAverage(*this)
413  << endl;
414  }
415 }
416 
417 
419 (
420  Ostream& os
421 ) const
422 {
425  (
426  "filmRegion",
427  "surfaceFilmProperties",
428  filmRegionName_
429  );
431  (
432  "pyrolysisRegion",
433  "pyrolysisProperties",
434  pyrolysisRegionName_
435  );
436  os.writeEntry("Tnbr", TnbrName_);
437  os.writeEntry("qr", qrName_);
438  os.writeEntry("convectiveScaling", convectiveScaling_);
439  os.writeEntry("filmDeltaDry", filmDeltaDry_);
440  os.writeEntry("filmDeltaWet", filmDeltaWet_);
442 }
443 
444 
445 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
446 
448 (
451 );
452 
453 
454 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
455 
456 } // End namespace Foam
457 
458 
459 // ************************************************************************* //
Foam::regionModels::pyrolysisModels::pyrolysisModel
Base class for pyrolysis models.
Definition: pyrolysisModel.H:61
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< scalar >
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:248
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::filmPyrolysisRadiativeCoupledMixedFvPatchScalarField
Mixed boundary condition for temperature, to be used in the flow and pyrolysis regions when a film re...
Definition: filmPyrolysisRadiativeCoupledMixedFvPatchScalarField.H:99
Foam::filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: filmPyrolysisRadiativeCoupledMixedFvPatchScalarField.C:419
Foam::mappedPatchBase::samplePolyPatch
const polyPatch & samplePolyPatch() const
Get the patch on the region.
Definition: mappedPatchBase.C:1662
Foam::temperatureCoupledBase
Common functions used in temperature coupled boundaries.
Definition: temperatureCoupledBase.H:135
Foam::regionModels::regionModel::time
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:40
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::fvPatch::name
virtual const word & name() const
Return name.
Definition: fvPatch.H:167
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:112
Foam::temperatureCoupledBase::autoMap
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
Definition: temperatureCoupledBase.C:172
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: filmPyrolysisRadiativeCoupledMixedFvPatchScalarField.C:229
Foam::filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::pyrolysisModelType
Foam::regionModels::pyrolysisModels::pyrolysisModel pyrolysisModelType
Definition: filmPyrolysisRadiativeCoupledMixedFvPatchScalarField.H:110
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
Foam::regionModels::regionModel::primaryMesh
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionModelI.H:34
Foam::regionModels::regionModel::toPrimary
void toPrimary(const label regionPatchi, List< Type > &regionField) const
Convert a local region field to the primary region.
Definition: regionModelTemplates.C:147
Foam::objectRegistry::lookupClass
HashTable< const Type * > lookupClass(const bool strict=false) const
Return all objects with a class satisfying isA<Type>
Foam::mappedPatchBase::sampleMesh
const polyMesh & sampleMesh() const
Get the region mesh.
Definition: mappedPatchBase.C:1649
Foam::filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::rmap
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: filmPyrolysisRadiativeCoupledMixedFvPatchScalarField.C:239
Foam::filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::K
tmp< scalarField > K() const
Get corresponding K field.
Foam::regionModels::surfaceFilmModels::thermoSingleLayer
Thermodynamic form of single-cell layer surface film model.
Definition: thermoSingleLayer.H:67
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::mappedPatchBase::distribute
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
Definition: mappedPatchBaseTemplates.C:30
Foam::regionModels::surfaceFilmModels::heatTransferModel::h
virtual tmp< volScalarField > h() const =0
Return the heat transfer coefficient [W/m2/K].
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::temperatureCoupledBase::alpha
virtual tmp< scalarField > alpha(const scalarField &Tp) const
Given patch temperature calculate corresponding alphaEff field.
Definition: temperatureCoupledBase.C:363
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::regionModels::regionModel::nbrCoupledPatchID
label nbrCoupledPatchID(const regionModel &nbrRegion, const label regionPatchi) const
Return the coupled patch ID paired with coupled patch.
Definition: regionModel.C:272
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:386
Foam::fvPatch::lookupPatchField
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
Definition: fvPatchFvMeshTemplates.C:34
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::fvPatch::index
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:197
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::temperatureCoupledBase::rmap
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
Definition: temperatureCoupledBase.C:188
Foam::patchIdentifier::index
label index() const noexcept
The index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:147
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::refCast
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::filmPyrolysisRadiativeCoupledMixedFvPatchScalarField
filmPyrolysisRadiativeCoupledMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: filmPyrolysisRadiativeCoupledMixedFvPatchScalarField.C:105
filmPyrolysisRadiativeCoupledMixedFvPatchScalarField.H
Foam::IOobject::name
const word & name() const noexcept
Return name.
Definition: IOobjectI.H:65
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::htcw
const heatTransferModel & htcw() const
Return const access to the (wall) heat transfer model.
Definition: thermoSingleLayerI.H:142
Foam::regionModels::regionModel::mapRegionPatchField
tmp< Foam::Field< Type > > mapRegionPatchField(const regionModel &nbrRegion, const label regionPatchi, const label nbrPatchi, const Field< Type > &nbrField, const bool flip=false) const
Map patch field from another region model to local patch.
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::fvPatch::deltaCoeffs
const scalarField & deltaCoeffs() const
Definition: fvPatch.C:196
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::temperatureCoupledBase::kappa
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
Definition: temperatureCoupledBase.C:210
Foam::filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: filmPyrolysisRadiativeCoupledMixedFvPatchScalarField.C:256
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::filmPyrolysisRadiativeCoupledMixedFvPatchScalarField::filmModelType
Foam::regionModels::surfaceFilmModels::thermoSingleLayer filmModelType
Definition: filmPyrolysisRadiativeCoupledMixedFvPatchScalarField.H:107
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::temperatureCoupledBase::write
void write(Ostream &os) const
Write.
Definition: temperatureCoupledBase.C:512
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::Ts
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
Definition: thermoSingleLayer.C:668
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::fvPatchField< scalar >::operator=
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:404
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::fvMesh::name
const word & name() const
Return reference to name.
Definition: fvMesh.H:300
mappedPatchBase.H