turbulentTemperatureRadCoupledMixedFvPatchScalarField.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-2017 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 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "mappedPatchBase.H"
34 #include "basicThermo.H"
35 #include "mappedPatchFieldBase.H"
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace compressible
42 {
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
48 (
49  const fvPatch& p,
51 )
52 :
53  mixedFvPatchScalarField(p, iF),
54  temperatureCoupledBase(patch()), // default method (fluidThermo)
56  (
58  *this
59  ),
60  TnbrName_("undefined-Tnbr"),
61  qrNbrName_("undefined-qrNbr"),
62  qrName_("undefined-qr"),
63  thermalInertia_(false)
64 {
65  this->refValue() = 0.0;
66  this->refGrad() = 0.0;
67  this->valueFraction() = 1.0;
68  this->source() = 0.0;
69 }
70 
71 
74 (
76  const fvPatch& p,
78  const fvPatchFieldMapper& mapper
79 )
80 :
81  mixedFvPatchScalarField(psf, p, iF, mapper),
83  mappedPatchFieldBase<scalar>
84  (
85  mappedPatchFieldBase<scalar>::mapper(p, iF),
86  *this,
87  psf
88  ),
89  TnbrName_(psf.TnbrName_),
90  qrNbrName_(psf.qrNbrName_),
91  qrName_(psf.qrName_),
92  thicknessLayers_(psf.thicknessLayers_),
93  thicknessLayer_(psf.thicknessLayer_.clone(p.patch())),
94  kappaLayers_(psf.kappaLayers_),
95  kappaLayer_(psf.kappaLayer_.clone(p.patch())),
96  thermalInertia_(psf.thermalInertia_)
97 {}
98 
99 
102 (
103  const fvPatch& p,
105  const dictionary& dict
106 )
107 :
108  mixedFvPatchScalarField(p, iF),
110  mappedPatchFieldBase<scalar>
111  (
112  mappedPatchFieldBase<scalar>::mapper(p, iF),
113  *this,
114  dict
115  ),
116  TnbrName_(dict.getOrDefault<word>("Tnbr", "T")),
117  qrNbrName_(dict.getOrDefault<word>("qrNbr", "none")),
118  qrName_(dict.getOrDefault<word>("qr", "none")),
119  thermalInertia_(dict.getOrDefault<Switch>("thermalInertia", false))
120 {
121  if (!isA<mappedPatchBase>(this->patch().patch()))
122  {
124  << "' not type '" << mappedPatchBase::typeName << "'"
125  << "\n for patch " << p.name()
126  << " of field " << internalField().name()
127  << " in file " << internalField().objectPath()
128  << exit(FatalError);
129  }
130 
131  //const auto* eptr = dict.findEntry("thicknessLayers");
132  //if (eptr)
133  //{
134  // // Detect either a list (parsed as a scalarList) or
135  // // a single entry (parsed as a PatchFunction1) or
136  //
137  // if
138  // (
139  // eptr->isStream()
140  // && eptr->stream().peek().isPunctuation(token::BEGIN_LIST)
141  // )
142  // {
143  // // Backwards compatibility
144  // thicknessLayers_ = dict.get<scalarList>("thicknessLayers");
145  // kappaLayers_ = dict.get<scalarList>("kappaLayers");
146  //
147  // if (thicknessLayers_.size() != kappaLayers_.size())
148  // {
149  // FatalIOErrorInFunction(dict) << "Inconstent sizes :"
150  // << "thicknessLayers:" << thicknessLayers_
151  // << "kappaLayers:" << kappaLayers_
152  // << exit(FatalIOError);
153  // }
154  // }
155  // else
156  // {
157  // thicknessLayer_ = PatchFunction1<scalar>::New
158  // (
159  // p.patch(),
160  // "thicknessLayers",
161  // dict
162  // );
163  // kappaLayer_ = PatchFunction1<scalar>::New
164  // (
165  // p.patch(),
166  // "kappaLayers",
167  // dict
168  // );
169  // }
170  //}
171 
172  // Read list of layers
173  if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
174  {
175  dict.readEntry("kappaLayers", kappaLayers_);
176  }
177  // Read single additional layer as PatchFunction1
178  thicknessLayer_ = PatchFunction1<scalar>::NewIfPresent
179  (
180  p.patch(),
181  "thicknessLayer",
182  dict
183  );
184  if (thicknessLayer_)
185  {
186  kappaLayer_ = PatchFunction1<scalar>::New
187  (
188  p.patch(),
189  "kappaLayer",
190  dict
191  );
192  }
193 
194  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
195 
196  if (dict.found("refValue"))
197  {
198  // Full restart
199  refValue() = scalarField("refValue", dict, p.size());
200  refGrad() = scalarField("refGradient", dict, p.size());
201  valueFraction() = scalarField("valueFraction", dict, p.size());
202  }
203  else
204  {
205  // Start from user entered data. Assume fixedValue.
206  refValue() = *this;
207  refGrad() = 0.0;
208  valueFraction() = 1.0;
209  }
210 
211  bool boolVal(false);
212  if (dict.readIfPresent("useImplicit", boolVal))
213  {
214  this->useImplicit(boolVal);
215  }
216  if (dict.found("source"))
217  {
218  source() = scalarField("source", dict, p.size());
219  }
220  else
221  {
222  source() = 0.0;
223  }
224 }
225 
226 
229 (
230  const turbulentTemperatureRadCoupledMixedFvPatchScalarField& psf,
232 )
233 :
234  mixedFvPatchScalarField(psf, iF),
236  mappedPatchFieldBase<scalar>
237  (
238  mappedPatchFieldBase<scalar>::mapper(patch(), iF),
239  *this,
240  psf
241  ),
242  TnbrName_(psf.TnbrName_),
243  qrNbrName_(psf.qrNbrName_),
244  qrName_(psf.qrName_),
245  thicknessLayers_(psf.thicknessLayers_),
246  thicknessLayer_(psf.thicknessLayer_.clone(patch().patch())),
247  kappaLayers_(psf.kappaLayers_),
248  kappaLayer_(psf.kappaLayer_.clone(patch().patch())),
249  thermalInertia_(psf.thermalInertia_)
250 {}
251 
252 
255 (
257 )
258 :
259  mixedFvPatchScalarField(psf),
261  mappedPatchFieldBase<scalar>
262  (
263  mappedPatchFieldBase<scalar>::mapper(patch(), psf.internalField()),
264  *this,
265  psf
266  ),
267  TnbrName_(psf.TnbrName_),
268  qrNbrName_(psf.qrNbrName_),
269  qrName_(psf.qrName_),
270  thicknessLayers_(psf.thicknessLayers_),
271  thicknessLayer_(psf.thicknessLayer_.clone(patch().patch())),
272  kappaLayers_(psf.kappaLayers_),
273  kappaLayer_(psf.kappaLayer_.clone(patch().patch())),
274  thermalInertia_(psf.thermalInertia_)
275 {}
276 
277 
278 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
279 
281 (
282  const fvPatchFieldMapper& mapper
283 )
284 {
285  mixedFvPatchScalarField::autoMap(mapper);
287  //mappedPatchFieldBase<scalar>::autoMap(mapper);
288  if (thicknessLayer_)
289  {
290  thicknessLayer_().autoMap(mapper);
291  kappaLayer_().autoMap(mapper);
292  }
293 }
294 
295 
297 (
298  const fvPatchField<scalar>& ptf,
299  const labelList& addr
300 )
301 {
302  mixedFvPatchScalarField::rmap(ptf, addr);
303 
305  refCast
306  <
308  >(ptf);
309 
310  temperatureCoupledBase::rmap(tiptf, addr);
311  //mappedPatchFieldBase<scalar>::rmap(ptf, addr);
312  if (thicknessLayer_)
313  {
314  thicknessLayer_().rmap(tiptf.thicknessLayer_(), addr);
315  kappaLayer_().rmap(tiptf.kappaLayer_(), addr);
316  }
317 }
318 
319 
322 (
323  const scalarField& Tp
324 ) const
325 {
326  // Get kappa from relevant thermo
328 
329  return tk;
330 }
331 
332 
334 {
335  if (updated())
336  {
337  return;
338  }
339 
340  const polyMesh& mesh = patch().boundaryMesh().mesh();
341 
342  // Since we're inside initEvaluate/evaluate there might be processor
343  // comms underway. Change the tag we use.
344  int oldTag = UPstream::msgType();
345  UPstream::msgType() = oldTag+1;
346 
347  // Get the coupling information from the mappedPatchBase
348  const label patchi = patch().index();
349  const mappedPatchBase& mpp =
351  (
352  patch(),
353  this->internalField()
354  );
355 
356  const scalarField Tc(patchInternalField());
357  const scalarField& Tp = *this;
358 
359  const scalarField kappaTp(kappa(Tp));
360  const scalarField KDelta(kappaTp*patch().deltaCoeffs());
361 
362 
363  scalarField TcNbr;
364  scalarField KDeltaNbr;
365 
366  if (mpp.sameWorld())
367  {
368  const polyMesh& nbrMesh = mpp.sampleMesh();
369  const label samplePatchi = mpp.samplePolyPatch().index();
370  const fvPatch& nbrPatch =
371  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
372 
373  const auto& nbrField = refCast
375  (
376  nbrPatch.lookupPatchField<volScalarField, scalar>(TnbrName_)
377  );
378 
379  // Swap to obtain full local values of neighbour K*delta
380  TcNbr = nbrField.patchInternalField();
381  KDeltaNbr = nbrField.kappa(nbrField)*nbrPatch.deltaCoeffs();
382  }
383  else
384  {
385  // Different world so use my region,patch. Distribution below will
386  // do the reordering.
387  TcNbr = patchInternalField();
388  KDeltaNbr = KDelta;
389  }
390  distribute(this->internalField().name() + "_value", TcNbr);
391  distribute(this->internalField().name() + "_weights", KDeltaNbr);
392 
393  scalarField KDeltaC(this->size(), GREAT);
394  if (thicknessLayer_ || thicknessLayers_.size())
395  {
396  // Harmonic averaging
397  {
398  KDeltaC = 0.0;
399 
400  if (thicknessLayer_)
401  {
402  const scalar t = db().time().timeOutputValue();
403  KDeltaC +=
404  thicknessLayer_().value(t)
405  /kappaLayer_().value(t);
406 
407  }
408  if (thicknessLayers_.size())
409  {
410  forAll(thicknessLayers_, iLayer)
411  {
412  KDeltaC += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
413  }
414  }
415  KDeltaC = 1.0/(KDeltaC + SMALL);
416  }
417  }
418 
419  scalarField alpha(kappaTp*(1 + KDeltaNbr/KDeltaC)*patch().deltaCoeffs());
420 
421  scalarField qr(Tp.size(), Zero);
422  if (qrName_ != "none")
423  {
424  qr = patch().lookupPatchField<volScalarField, scalar>(qrName_);
425  }
426 
427  scalarField qrNbr(Tp.size(), Zero);
428  if (qrNbrName_ != "none")
429  {
430  if (mpp.sameWorld())
431  {
432  const polyMesh& nbrMesh = mpp.sampleMesh();
433  const label samplePatchi = mpp.samplePolyPatch().index();
434  const fvPatch& nbrPatch =
435  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
436  qrNbr =
437  nbrPatch.lookupPatchField<volScalarField, scalar>(qrNbrName_);
438  }
439  else
440  {
441  qrNbr =
442  patch().lookupPatchField<volScalarField, scalar>(qrNbrName_);
443  }
444  distribute(qrNbrName_, qrNbr);
445  }
446 
447  // inertia therm
448  if (thermalInertia_ && !mpp.sameWorld())
449  {
451  << "thermalInertia not supported in combination with multi-world"
452  << exit(FatalError);
453  }
454  if (thermalInertia_)
455  {
456  const scalar dt = mesh.time().deltaTValue();
457  scalarField mCpDtNbr;
458 
459  {
460  const polyMesh& nbrMesh = mpp.sampleMesh();
461 
462  const basicThermo* thermo =
463  nbrMesh.findObject<basicThermo>(basicThermo::dictName);
464 
465  if (thermo)
466  {
467  const label samplePatchi = mpp.samplePolyPatch().index();
468  const fvPatch& nbrPatch =
469  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
470  const scalarField& ppn =
471  thermo->p().boundaryField()[samplePatchi];
472  const scalarField& Tpn =
473  thermo->T().boundaryField()[samplePatchi];
474 
475  mCpDtNbr =
476  (
477  thermo->Cp(ppn, Tpn, samplePatchi)
478  * thermo->rho(samplePatchi)
479  / nbrPatch.deltaCoeffs()/dt
480  );
481 
482  mpp.distribute(mCpDtNbr);
483  }
484  else
485  {
486  mCpDtNbr.setSize(Tp.size(), Zero);
487  }
488  }
489 
490  scalarField mCpDt;
491 
492  // Local inertia therm
493  {
494  const basicThermo* thermo =
495  mesh.findObject<basicThermo>(basicThermo::dictName);
496 
497  if (thermo)
498  {
499  const scalarField& pp = thermo->p().boundaryField()[patchi];
500 
501  mCpDt =
502  (
503  thermo->Cp(pp, Tp, patchi)
504  * thermo->rho(patchi)
505  / patch().deltaCoeffs()/dt
506  );
507  }
508  else
509  {
510  // Issue warning?
511  mCpDt.setSize(Tp.size(), Zero);
512  }
513  }
514 
515  const volScalarField& T =
516  this->db().lookupObject<volScalarField>
517  (
518  this->internalField().name()
519  );
520 
521  const fvPatchField<scalar>& TpOld = T.oldTime().boundaryField()[patchi];
522 
523  scalarField alpha(KDeltaNbr + mCpDt + mCpDtNbr);
524 
525  valueFraction() = alpha/(alpha + KDelta);
526  scalarField c(KDeltaNbr*TcNbr + (mCpDt + mCpDtNbr)*TpOld);
527  refValue() = c/alpha;
528  refGrad() = (qr + qrNbr)/kappaTp;
529  }
530  else
531  {
532 
533  valueFraction() = KDeltaNbr/(KDeltaNbr + alpha);
534  refValue() = TcNbr;
535  refGrad() = (qr + qrNbr)/kappaTp;
536  }
537 
538  source() = Zero;
539 
540  // If useImplicit is true we need the source term associated with this BC
541  if (this->useImplicit())
542  {
543  source() =
544  alphaSfDelta()*
545  (
546  valueFraction()*deltaH()
547  + (qr + qrNbr)/beta()
548  );
549  }
550 
551  mixedFvPatchScalarField::updateCoeffs();
552 
553  if (debug)
554  {
555  scalar Q = gSum(kappaTp*patch().magSf()*snGrad());
556 
557  Info<< patch().boundaryMesh().mesh().name() << ':'
558  << patch().name() << ':'
559  << this->internalField().name() << " <- "
560  << mpp.sampleRegion() << ':'
561  << mpp.samplePatch() << ':'
562  << this->internalField().name() << " :"
563  << " heat transfer rate:" << Q
564  << " walltemperature "
565  << " min:" << gMin(Tp)
566  << " max:" << gMax(Tp)
567  << " avg:" << gAverage(Tp)
568  << endl;
569  }
571  // Restore tag
572  UPstream::msgType() = oldTag;
573 }
574 
575 
577 (
578  fvMatrix<scalar>& m,
579  const label iMatrix,
580  const direction cmpt
581 )
582 {
584  << "This T BC does not support energy coupling "
585  << "It is implemented on he field "
586  << abort(FatalError);
587 }
588 
589 
590 tmp<Field<scalar>> turbulentTemperatureRadCoupledMixedFvPatchScalarField::coeffs
591 (
592  fvMatrix<scalar>& matrix,
593  const Field<scalar>& coeffs,
594  const label mat
595 ) const
596 {
598  << "This BC does not support energy coupling "
599  << "Use compressible::turbulentTemperatureRadCoupledMixed "
600  << "which has more functionalities and it can handle "
601  << "the assemble coupled option for energy. "
602  << abort(FatalError);
603 
604  return tmp<Field<scalar>>(new Field<scalar>());
605 }
606 
607 
608 tmp<scalarField>
609 turbulentTemperatureRadCoupledMixedFvPatchScalarField::alphaSfDelta() const
610 {
611  return (alpha(*this)*patch().deltaCoeffs()*patch().magSf());
612 }
613 
614 
615 tmp<scalarField> turbulentTemperatureRadCoupledMixedFvPatchScalarField::
616 beta() const
617 {
618  const mappedPatchBase& mpp =
619  refCast<const mappedPatchBase>(patch().patch());
620 
621  if (!mpp.sameWorld())
622  {
624  << "coupled energy not supported in combination with multi-world"
625  << exit(FatalError);
626  }
627 
628  const label samplePatchi = mpp.samplePolyPatch().index();
629  const polyMesh& nbrMesh = mpp.sampleMesh();
630 
631  const fvPatch& nbrPatch =
632  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
633 
635  nbrField = refCast
637  (
638  nbrPatch.lookupPatchField<volScalarField, scalar>(TnbrName_)
639  );
640 
641  // Swap to obtain full local values of neighbour internal field
642  scalarField TcNbr(nbrField.patchInternalField());
643  mpp.distribute(TcNbr);
644 
645  scalarField alphaDeltaNbr
646  (
647  nbrField.alpha(TcNbr)*nbrPatch.deltaCoeffs()
648  );
649  mpp.distribute(alphaDeltaNbr);
650 
651  scalarField alphaDelta
652  (
653  alpha(*this)*patch().deltaCoeffs()
654  );
655 
656  return (alphaDeltaNbr + alphaDelta);
657 }
658 
659 
660 tmp<scalarField> turbulentTemperatureRadCoupledMixedFvPatchScalarField::
661 deltaH() const
662 {
663  const mappedPatchBase& mpp =
664  refCast<const mappedPatchBase>(patch().patch());
665 
666  if (!mpp.sameWorld())
667  {
669  << "coupled energy not supported in combination with multi-world"
670  << exit(FatalError);
671  }
672 
673  const polyMesh& nbrMesh = mpp.sampleMesh();
674 
675  const basicThermo* nbrThermo =
676  nbrMesh.cfindObject<basicThermo>(basicThermo::dictName);
677 
678  const polyMesh& mesh = patch().boundaryMesh().mesh();
679 
680  const basicThermo* localThermo =
681  mesh.cfindObject<basicThermo>(basicThermo::dictName);
682 
683 
684  if (nbrThermo && localThermo)
685  {
686  const label patchi = patch().index();
687  const scalarField& pp = localThermo->p().boundaryField()[patchi];
688  const scalarField& Tp = *this;
689 
690  const mappedPatchBase& mpp =
691  refCast<const mappedPatchBase>(patch().patch());
692 
693  const label patchiNrb = mpp.samplePolyPatch().index();
694 
695  const scalarField& ppNbr = nbrThermo->p().boundaryField()[patchiNrb];
696  //const scalarField& TpNbr = nbrThermo->T().boundaryField()[patchiNrb];
697 
698  // Use this Tp to evaluate he jump as this is updated while doing
699  // updateCoeffs on boundary fields which set T values on patches
700  // then non consistent Tp and Tpnbr could be used from different
701  // updated values (specially when T changes drastically between time
702  // steps/
703  return
704  (
705  - localThermo->he(pp, Tp, patchi)
706  + nbrThermo->he(ppNbr, Tp, patchiNrb)
707  );
708  }
709  else
710  {
712  << "Can't find thermos on mapped patch "
713  << " method, but thermo package not available"
714  << exit(FatalError);
715  }
716 
717  return tmp<scalarField>::New(patch().size(), Zero);
718 }
719 
720 
722 (
723  Ostream& os
724 ) const
725 {
727 
728  //os.writeEntry("Tnbr", TnbrName_);
729  os.writeEntryIfDifferent<word>("Tnbr", "T", TnbrName_);
730 
731  //os.writeEntry("qrNbr", qrNbrName_);
732  os.writeEntryIfDifferent<word>("qrNbr", "none", qrNbrName_);
733  //os.writeEntry("qr", qrName_);
734  os.writeEntryIfDifferent<word>("qr", "none", qrName_);
735  if (thermalInertia_)
736  {
737  os.writeEntry("thermalInertia", thermalInertia_);
738  }
739 
740  if (thicknessLayer_)
741  {
742  thicknessLayer_().writeData(os);
743  kappaLayer_().writeData(os);
744  }
745  if (thicknessLayers_.size())
746  {
747  thicknessLayers_.writeEntry("thicknessLayers", os);
748  kappaLayers_.writeEntry("kappaLayers", os);
749  }
750 
753 }
754 
755 
756 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
757 
759 (
761  turbulentTemperatureRadCoupledMixedFvPatchScalarField
762 );
763 
764 
765 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
766 
767 } // End namespace compressible
768 } // End namespace Foam
769 
770 
771 // ************************************************************************* //
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
faceListList boundary
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
uint8_t direction
Definition: direction.H:46
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
static const mappedPatchBase & mapper(const fvPatch &p, const DimensionedField< Type, volMesh > &iF)
Check that patch is of correct type.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
Type gMin(const FieldField< Field, Type > &f)
const polyMesh & sampleMesh() const
Get the region mesh.
virtual void write(Ostream &os) const
Write.
Type & refCast(U &obj)
A dynamic_cast (for references) that generates FatalError on failed casts, uses the virtual type() me...
Definition: typeInfo.H:151
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
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.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, any/none. Also accepts 0/1 as a string and shortcuts t/f, y/n.
Definition: Switch.H:77
Mixed boundary condition for temperature and radiation heat transfer to be used for in multiregion ca...
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:806
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const AnyType *=nullptr) const
Lookup the named field from the local registry and return the patch field corresponding to this patch...
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
void write(Ostream &os) const
Write.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
psiReactionThermo & thermo
Definition: createFields.H:28
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
bool sameWorld() const
Is sample world the local world?
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
fvPatchField< scalar > fvPatchScalarField
virtual void manipulateMatrix(fvMatrix< scalar > &m, const label iMatrix, const direction cmpt)
Manipulate matrix.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
A FieldMapper for finite-volume patch fields.
turbulentTemperatureRadCoupledMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:203
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
errorManip< error > abort(error &err)
Definition: errorManip.H:139
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field. Override temperatureCoupledBase::kappa to in...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:55
bool compressible
Definition: pEqn.H:2
int debug
Static debugging option.
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
Common functions used in temperature coupled boundaries.
Type gAverage(const FieldField< Field, Type > &f)
const dimensionedScalar c
Speed of light in a vacuum.
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:352
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
static autoPtr< PatchFunction1< Type > > NewIfPresent(const polyPatch &pp, const word &entryName, const dictionary &dict, const bool faceValues=true)
An optional selector.
messageStream Info
Information stream (stdout output on master, null elsewhere)
const word & samplePatch() const
Patch (only if NEARESTPATCHFACE)
void writeEntry(Ostream &os) const
Write the UList with its compound type.
Definition: UListIO.C:31
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
label index() const noexcept
The index of this patch in the boundaryMesh.
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const word & sampleRegion() const
Region to sample.
const polyPatch & samplePolyPatch() const
Get the patch on the region.
void distribute(const word &fieldName, Field< T > &newValues) const
Wrapper for mapDistribute::distribute that knows about dabase mapping.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:40
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157