temperatureCoupledBase.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) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2022 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 
29 #include "temperatureCoupledBase.H"
30 #include "volFields.H"
31 #include "fluidThermo.H"
32 #include "solidThermo.H"
34 #include "multiphaseInterSystem.H"
35 
36 // * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
37 
38 const Foam::Enum
39 <
41 >
43 {
44  { KMethodType::mtFluidThermo, "fluidThermo" },
45  { KMethodType::mtSolidThermo, "solidThermo" },
46  { KMethodType::mtDirectionalSolidThermo, "directionalSolidThermo" },
47  { KMethodType::mtLookup, "lookup" },
48  { KMethodType::mtFunction, "function" }
49 };
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const fvPatch& patch,
57  const KMethodType method
58 )
59 :
60  patch_(patch),
61  method_(method),
62  kappaName_(),
63  alphaName_(),
64  alphaAniName_(),
65  kappaFunction1_(nullptr),
66  alphaFunction1_(nullptr)
67 {
68  switch (method_)
69  {
71  case mtLookup:
72  case mtFunction:
73  {
75  << "Cannot construct kappaMethod: "
76  << KMethodTypeNames_[method_] << " without a dictionary"
77  << abort(FatalError);
78  break;
79  }
80  default:
81  {
82  break;
83  }
84  }
85 }
86 
87 
89 (
90  const fvPatch& patch,
91  const KMethodType method,
92  const word& kappaName,
93  const word& alphaName,
94  const word& alphaAniName
95 )
96 :
97  patch_(patch),
98  method_(method),
99  kappaName_(kappaName),
100  alphaName_(alphaName),
101  alphaAniName_(alphaAniName),
102  kappaFunction1_(nullptr),
103  alphaFunction1_(nullptr)
104 {
105  switch (method_)
106  {
107  case mtFunction:
108  {
110  << "Cannot construct kappaMethod: "
111  << KMethodTypeNames_[method_] << " without a dictionary"
112  << abort(FatalError);
113  break;
114  }
115  default:
116  {
117  break;
118  }
119  }
120 }
121 
122 
124 (
125  const fvPatch& patch,
126  const dictionary& dict
127 )
128 :
129  patch_(patch),
130  method_(KMethodTypeNames_.get("kappaMethod", dict)),
131  kappaName_(dict.getOrDefault<word>("kappa", word::null)),
132  alphaName_(dict.getOrDefault<word>("alpha", word::null)),
133  alphaAniName_(dict.getOrDefault<word>("alphaAni", word::null)),
134  kappaFunction1_(nullptr),
135  alphaFunction1_(nullptr)
136 {
137  switch (method_)
138  {
140  {
141  if (!dict.found("alphaAni"))
142  {
144  << "Did not find entry 'alphaAni'"
145  " required for 'kappaMethod' "
147  << exit(FatalIOError);
148  }
149 
150  break;
151  }
152 
153  case mtLookup:
154  {
155  if (!dict.found("kappa"))
156  {
158  << "Did not find entry 'kappa'"
159  " required for 'kappaMethod' "
161  << "Please set 'kappa' to the name of"
162  " a volScalar or volSymmTensor field" << nl
163  << exit(FatalIOError);
164  }
165 
166  break;
167  }
168 
169  case mtFunction:
170  {
171  kappaFunction1_ = PatchFunction1<scalar>::New
172  (
173  patch_.patch(),
174  "kappaValue",
175  dict
176  );
177  alphaFunction1_ = PatchFunction1<scalar>::New
178  (
179  patch_.patch(),
180  "alphaValue",
181  dict
182  );
183 
184  break;
185  }
186 
187  default:
188  {
189  break;
190  }
191  }
192 }
193 
194 
196 (
197  const temperatureCoupledBase& base
198 )
199 :
200  temperatureCoupledBase(base.patch_, base)
201 {}
202 
203 
205 (
206  const fvPatch& patch,
207  const temperatureCoupledBase& base
208 )
209 :
210  patch_(patch),
211  method_(base.method_),
212  kappaName_(base.kappaName_),
213  alphaName_(base.alphaName_),
214  alphaAniName_(base.alphaAniName_),
215  kappaFunction1_(base.kappaFunction1_.clone(patch_.patch())),
216  alphaFunction1_(base.alphaFunction1_.clone(patch_.patch()))
217 {}
218 
219 
220 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221 
223 (
224  const fvPatchFieldMapper& mapper
225 )
226 {
227  if (kappaFunction1_)
228  {
229  kappaFunction1_().autoMap(mapper);
230  }
231  if (alphaFunction1_)
232  {
233  alphaFunction1_().autoMap(mapper);
234  }
235 }
236 
237 
239 (
240  const fvPatchField<scalar>& ptf,
241  const labelList& addr
242 )
243 {
244  const auto* tcb = isA<temperatureCoupledBase>(ptf);
245 
246  if (tcb)
247  {
248  if (kappaFunction1_)
249  {
250  kappaFunction1_().rmap(tcb->kappaFunction1_(), addr);
251  }
252  if (alphaFunction1_)
253  {
254  alphaFunction1_().rmap(tcb->alphaFunction1_(), addr);
255  }
256  }
257 }
258 
259 
261 (
262  const scalarField& Tp
263 ) const
264 {
265  const fvMesh& mesh = patch_.boundaryMesh().mesh();
266  const label patchi = patch_.index();
267 
268  switch (method_)
269  {
270  case mtFluidThermo:
271  {
273 
274  {
275  const auto* ptr =
277  (
279  );
280 
281  if (ptr)
282  {
283  return ptr->kappaEff(patchi);
284  }
285  }
286 
287  {
288  const auto* ptr =
289  mesh.cfindObject<fluidThermo>(basicThermo::dictName);
290 
291  if (ptr)
292  {
293  return ptr->kappa(patchi);
294  }
295  }
296 
297  {
298  const auto* ptr =
299  mesh.cfindObject<basicThermo>(basicThermo::dictName);
300 
301  if (ptr)
302  {
303  return ptr->kappa(patchi);
304  }
305  }
306 
307  {
308  const auto* ptr =
309  mesh.cfindObject<multiphaseInterSystem>
310  (
312  );
313 
314  if (ptr)
315  {
316  return ptr->kappaEff(patchi);
317  }
318  }
319 
321  << "Using kappaMethod " << KMethodTypeNames_[method_]
322  << ", but thermo package not available\n"
323  << exit(FatalError);
324 
325  break;
326  }
327 
328  case mtSolidThermo:
329  {
330  const solidThermo& thermo =
332 
333  return thermo.kappa(patchi);
334 
335  break;
336  }
337 
338  case mtDirectionalSolidThermo:
339  {
340  const solidThermo& thermo =
342 
343  const symmTensorField& alphaAni =
344  patch_.lookupPatchField<volSymmTensorField, scalar>
345  (
346  alphaAniName_
347  );
348 
349  const scalarField& pp = thermo.p().boundaryField()[patchi];
350 
351  const symmTensorField kappa(alphaAni*thermo.Cp(pp, Tp, patchi));
352 
353  const vectorField n(patch_.nf());
354 
355  return n & kappa & n;
356  }
357 
358  case mtLookup:
359  {
360  {
361  const auto* ptr =
362  mesh.cfindObject<volScalarField>(kappaName_);
363 
364  if (ptr)
365  {
366  return patch_.patchField<volScalarField>(*ptr);
367  }
368  }
369 
370  {
371  const auto* ptr =
372  mesh.cfindObject<volSymmTensorField>(kappaName_);
373 
374  if (ptr)
375  {
376  const symmTensorField& wallValues =
377  patch_.patchField<volSymmTensorField>(*ptr);
378 
379  const vectorField n(patch_.nf());
380 
381  return n & wallValues & n;
382  }
383  }
384 
385 
387  << "Did not find field '" << kappaName_
388  << "' on mesh " << mesh.name()
389  << " patch " << patch_.name() << nl
390  << "Please set 'kappa' to the name of"
391  " a volScalar or volSymmTensor field"
392  ", or use another method" << nl
393  << " " << flatOutput(KMethodTypeNames_.sortedToc()) << nl
394  << exit(FatalError);
395 
396  break;
397  }
398 
399  case KMethodType::mtFunction:
400  {
401  const auto& tm = patch_.patch().boundaryMesh().mesh().time();
402  return kappaFunction1_->value(tm.timeOutputValue());
403  break;
404  }
405 
406  default:
407  {
409  << "Unimplemented method " << KMethodTypeNames_[method_] << nl
410  << "Please set 'kappaMethod' to one of "
411  << flatOutput(KMethodTypeNames_.sortedToc()) << nl
412  << "If kappaMethod=lookup, also set 'kappa' to the name of"
413  " a volScalar or volSymmTensor field" << nl
414  << exit(FatalError);
415 
416  break;
417  }
418  }
419 
420  return scalarField();
421 }
422 
423 
425 (
426  const scalarField& Tp
427 ) const
428 {
429  const fvMesh& mesh = patch_.boundaryMesh().mesh();
430  const label patchi = patch_.index();
431 
432  switch (method_)
433  {
434  case mtFluidThermo:
435  {
437 
438  {
439  const auto* ptr =
441  (
443  );
444 
445  if (ptr)
446  {
447  return ptr->alphaEff(patchi);
448  }
449  }
450 
451  {
452  const auto* ptr =
453  mesh.cfindObject<fluidThermo>(basicThermo::dictName);
454 
455  if (ptr)
456  {
457  return ptr->alpha(patchi);
458  }
459  }
460 
461  {
462  const auto* ptr =
463  mesh.cfindObject<basicThermo>(basicThermo::dictName);
464 
465  if (ptr)
466  {
467  return ptr->alpha(patchi);
468  }
469  }
470 
471  {
472  const auto* ptr =
473  mesh.cfindObject<basicThermo>("phaseProperties");
474 
475  if (ptr)
476  {
477  return ptr->alpha(patchi);
478  }
479  }
480 
482  << "Using kappaMethod " << KMethodTypeNames_[method_]
483  << ", but thermo package not available\n"
484  << exit(FatalError);
485 
486  break;
487  }
488 
489  case mtSolidThermo:
490  {
491  const solidThermo& thermo =
493 
494  return thermo.alpha(patchi);
495  break;
496  }
497 
498  case mtDirectionalSolidThermo:
499  {
500  const symmTensorField& alphaAni =
501  patch_.lookupPatchField<volSymmTensorField, scalar>
502  (
503  alphaAniName_
504  );
505 
506  const vectorField n(patch_.nf());
507 
508  return n & alphaAni & n;
509  }
510 
511  case mtLookup:
512  {
513  {
514  const auto* ptr =
515  mesh.cfindObject<volScalarField>(alphaName_);
516 
517  if (ptr)
518  {
519  return patch_.patchField<volScalarField>(*ptr);
520  }
521  }
522 
523  {
524  const auto* ptr =
525  mesh.cfindObject<volSymmTensorField>(alphaName_);
526 
527  if (ptr)
528  {
529  const symmTensorField& wallValues =
530  patch_.patchField<volSymmTensorField>(*ptr);
531 
532  const vectorField n(patch_.nf());
533 
534  return n & wallValues & n;
535  }
536  }
537 
539  << "Did not find field '" << alphaName_
540  << "' on mesh " << mesh.name()
541  << " patch " << patch_.name() << nl
542  << "Please set 'alpha' to the name of"
543  " a volScalar or volSymmTensor field"
544  ", or use another method" << nl
545  << " " << flatOutput(KMethodTypeNames_.sortedToc()) << nl
546  << exit(FatalError);
547 
548  break;
549  }
550 
551  case KMethodType::mtFunction:
552  {
553  const auto& tm = patch_.patch().boundaryMesh().mesh().time();
554  return alphaFunction1_->value(tm.timeOutputValue());
555  break;
556  }
557 
558  default:
559  {
561  << "Unimplemented method " << KMethodTypeNames_[method_] << nl
562  << "Please set 'kappaMethod' to one of "
563  << flatOutput(KMethodTypeNames_.sortedToc()) << nl
564  << "If kappaMethod=lookup, also set 'alpha' to the name of"
565  " a volScalar or volSymmTensor field" << nl
566  << exit(FatalError);
567 
568  break;
569  }
570  }
571 
572  return scalarField();
573 }
574 
575 
576 void Foam::temperatureCoupledBase::write(Ostream& os) const
577 {
578  os.writeEntry("kappaMethod", KMethodTypeNames_[method_]);
579  if (!kappaName_.empty())
580  {
581  os.writeEntry("kappa", kappaName_);
582  }
583  if (!alphaAniName_.empty())
584  {
585  os.writeEntry("alphaAni", alphaAniName_);
586  }
587  if (!alphaName_.empty())
588  {
589  os.writeEntry("alpha", alphaName_);
590  }
591  if (kappaFunction1_)
592  {
593  kappaFunction1_->writeData(os);
594  }
595  if (alphaFunction1_)
596  {
597  alphaFunction1_->writeData(os);
598  }
599 }
600 
601 
602 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:88
dictionary dict
const fvPatch & patch_
Underlying patch.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
virtual tmp< volScalarField > kappa() const =0
Thermal diffusivity for temperature of mixture [J/m/s/K].
Field< symmTensor > symmTensorField
Specialisation of Field<T> for symmTensor.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
static const word phasePropertiesName
Default name of the phase properties dictionary.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
autoPtr< PatchFunction1< scalar > > alphaFunction1_
Function1 for alpha.
KMethodType
Type of supplied Kappa.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:68
virtual tmp< scalarField > alpha(const scalarField &Tp) const
Given patch temperature calculate corresponding alphaEff field.
virtual const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Definition: basicThermo.C:627
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:312
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
void write(Ostream &os) const
Write.
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:603
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
psiReactionThermo & thermo
Definition: createFields.H:28
autoPtr< PatchFunction1< scalar > > kappaFunction1_
Function1 for kappa.
dynamicFvMesh & mesh
const polyMesh & mesh() const noexcept
Return the mesh reference.
static const word propertiesName
Default name of the turbulence properties dictionary.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
temperatureCoupledBase(const fvPatch &patch, const KMethodType method=KMethodType::mtFluidThermo)
Default construct from patch, using fluidThermo (default) or specified method.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
OBJstream os(runTime.globalPath()/outputName)
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
const KMethodType method_
How to get K.
Common functions used in temperature coupled boundaries.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
const word & name() const
Return reference to name.
Definition: fvMesh.H:388
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:607
const std::string patch
OpenFOAM patch number as a std::string.
static const Enum< KMethodType > KMethodTypeNames_
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const polyPatch & patch() const noexcept
Return the polyPatch.
Definition: fvPatch.H:200
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...