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