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-2024 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 "IOField.H"
36 #include "mappedPatchFieldBase.H"
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace compressible
43 {
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
48 turbulentTemperatureRadCoupledMixedFvPatchScalarField::getOrCreateField
49 (
50  const word& fieldName
51 ) const
52 {
53  const fvMesh& mesh = patch().boundaryMesh().mesh();
54 
55  auto* ptr = mesh.getObjectPtr<volScalarField>(fieldName);
56 
57  if (!ptr)
58  {
59  ptr = new volScalarField
60  (
61  IOobject
62  (
63  fieldName,
64  mesh.time().timeName(),
65  mesh.thisDb(),
69  ),
70  mesh,
72  );
73  regIOobject::store(ptr);
74  }
75 
76  return *ptr;
77 }
78 
79 
80 void turbulentTemperatureRadCoupledMixedFvPatchScalarField::storeHTCFields
81 (
82  const word& prefix,
83  const scalarField& shtc,
84  const scalarField& shtcPatch
85 )
86 const
87 {
88  volScalarField& htc =
89  getOrCreateField(IOobject::scopedName(prefix, "htc"));
90  htc.boundaryFieldRef()[patch().index()] = shtc;
91 
92  volScalarField& htcPatch =
93  getOrCreateField(IOobject::scopedName(prefix, "htcPatch"));
94  htcPatch.boundaryFieldRef()[patch().index()] = shtcPatch;
95 }
96 
97 
98 void turbulentTemperatureRadCoupledMixedFvPatchScalarField::writeFileHeader
99 (
100  Ostream& os
101 )
102 {
103  writeCommented(os, "Time");
104  writeTabbed(os, "Q_[W]");
105  writeTabbed(os, "q_[W/m^2]");
106  writeTabbed(os, "HTCavg_[W/m^2/K]");
107  writeTabbed(os, "patchHTCavg_[W/m^2/K]");
108  writeTabbed(os, "TpMin_[K]");
109  writeTabbed(os, "TpMax_[K]");
110  writeTabbed(os, "TpAvg_[K]");
111  writeTabbed(os, "TpNbrMin_[K]");
112  writeTabbed(os, "TpNbrMax_[K]");
113  writeTabbed(os, "TpNbrAvg_[K]");
114 
115  os << endl;
116 
117  writtenHeader_ = true;
118  updateHeader_ = false;
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
123 
126 (
127  const fvPatch& p,
129 )
130 :
131  mixedFvPatchScalarField(p, iF),
132  temperatureCoupledBase(patch()), // default method (fluidThermo)
133  mappedPatchFieldBase<scalar>
134  (
135  mappedPatchFieldBase<scalar>::mapper(p, iF),
136  *this
137  ),
138  functionObjects::writeFile
139  (
140  db(),
141  "turbulentTemperatureRadCoupledMixed",
142  "undefined",
143  false
144  ),
145  TnbrName_("undefined-Tnbr"),
146  qrNbrName_("undefined-qrNbr"),
147  qrName_("undefined-qr"),
148  logInterval_(-1),
149  executionIndex_(0),
150  thermalInertia_(false),
151  verbose_(false),
152  prefix_()
153 {
154  this->refValue() = Zero;
155  this->refGrad() = Zero;
156  this->valueFraction() = 1.0;
157  this->source() = 0.0;
158 }
159 
160 
163 (
165  const fvPatch& p,
167  const fvPatchFieldMapper& mapper
168 )
169 :
170  mixedFvPatchScalarField(psf, p, iF, mapper),
172  mappedPatchFieldBase<scalar>
173  (
174  mappedPatchFieldBase<scalar>::mapper(p, iF),
175  *this,
176  psf
177  ),
178  functionObjects::writeFile(psf),
179  TnbrName_(psf.TnbrName_),
180  qrNbrName_(psf.qrNbrName_),
181  qrName_(psf.qrName_),
182  thicknessLayers_(psf.thicknessLayers_),
183  thicknessLayer_(psf.thicknessLayer_.clone(p.patch())),
184  kappaLayers_(psf.kappaLayers_),
185  kappaLayer_(psf.kappaLayer_.clone(p.patch())),
186  logInterval_(psf.logInterval_),
187  executionIndex_(psf.executionIndex_),
188  thermalInertia_(psf.thermalInertia_),
189  verbose_(psf.verbose_),
190  prefix_(psf.prefix_)
191 {}
192 
193 
196 (
197  const fvPatch& p,
199  const dictionary& dict
200 )
201 :
202  mixedFvPatchScalarField(p, iF),
204  mappedPatchFieldBase<scalar>
205  (
206  mappedPatchFieldBase<scalar>::mapper(p, iF),
207  *this,
208  dict
209  ),
210  functionObjects::writeFile
211  (
212  db(),
213  "turbulentTemperatureRadCoupledMixed",
214  patch().name(),
215  false
216  ),
217  TnbrName_(dict.getOrDefault<word>("Tnbr", "T")),
218  qrNbrName_(dict.getOrDefault<word>("qrNbr", "none")),
219  qrName_(dict.getOrDefault<word>("qr", "none")),
220  logInterval_(dict.getOrDefault<scalar>("logInterval", -1)),
221  executionIndex_(0),
222  thermalInertia_(dict.getOrDefault<Switch>("thermalInertia", false)),
223  verbose_(dict.getOrDefault<bool>("verbose", false)),
224  prefix_(dict.getOrDefault<word>("prefix", "multiWorld"))
225 {
226  if (!isA<mappedPatchBase>(this->patch().patch()))
227  {
229  << "' not type '" << mappedPatchBase::typeName << "'"
230  << "\n for patch " << p.name()
231  << " of field " << internalField().name()
232  << " in file " << internalField().objectPath()
233  << exit(FatalError);
234  }
235 
236  //const auto* eptr = dict.findEntry("thicknessLayers", keyType::LITERAL);
237  //if (eptr)
238  //{
239  // // Detect either a list (parsed as a scalarList) or
240  // // a single entry (parsed as a PatchFunction1) or
241  //
242  // if
243  // (
244  // eptr->isStream()
245  // && eptr->stream().peek().isPunctuation(token::BEGIN_LIST)
246  // )
247  // {
248  // // Backwards compatibility
249  // thicknessLayers_ = dict.get<scalarList>("thicknessLayers");
250  // kappaLayers_ = dict.get<scalarList>("kappaLayers");
251  //
252  // if (thicknessLayers_.size() != kappaLayers_.size())
253  // {
254  // FatalIOErrorInFunction(dict) << "Inconstent sizes :"
255  // << "thicknessLayers:" << thicknessLayers_
256  // << "kappaLayers:" << kappaLayers_
257  // << exit(FatalIOError);
258  // }
259  // }
260  // else
261  // {
262  // thicknessLayer_ = PatchFunction1<scalar>::New
263  // (
264  // p.patch(),
265  // "thicknessLayers",
266  // dict
267  // );
268  // kappaLayer_ = PatchFunction1<scalar>::New
269  // (
270  // p.patch(),
271  // "kappaLayers",
272  // dict
273  // );
274  // }
275  //}
276 
277  // Read list of layers
278  if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
279  {
280  dict.readEntry("kappaLayers", kappaLayers_);
281  }
282  // Read single additional layer as PatchFunction1
283  thicknessLayer_ = PatchFunction1<scalar>::NewIfPresent
284  (
285  p.patch(),
286  "thicknessLayer",
287  dict
288  );
289  if (thicknessLayer_)
290  {
291  kappaLayer_ = PatchFunction1<scalar>::New
292  (
293  p.patch(),
294  "kappaLayer",
295  dict
296  );
297  }
298 
299  this->readValueEntry(dict, IOobjectOption::MUST_READ);
300 
301  if (this->readMixedEntries(dict))
302  {
303  // Full restart
304  }
305  else
306  {
307  // Start from user entered data. Assume fixedValue.
308  refValue() = *this;
309  refGrad() = Zero;
310  valueFraction() = 1.0;
311  }
312 
313  bool boolVal(false);
314  if (dict.readIfPresent("useImplicit", boolVal))
315  {
316  this->useImplicit(boolVal);
317  }
318  if (dict.found("source"))
319  {
320  source() = scalarField("source", dict, p.size());
321  }
322  else
323  {
324  source() = 0.0;
325  }
328 }
329 
330 
333 (
336 )
337 :
338  mixedFvPatchScalarField(psf, iF),
340  mappedPatchFieldBase<scalar>
341  (
342  mappedPatchFieldBase<scalar>::mapper(patch(), iF),
343  *this,
344  psf
345  ),
346  functionObjects::writeFile(psf),
347  TnbrName_(psf.TnbrName_),
348  qrNbrName_(psf.qrNbrName_),
349  qrName_(psf.qrName_),
350  thicknessLayers_(psf.thicknessLayers_),
351  thicknessLayer_(psf.thicknessLayer_.clone(patch().patch())),
352  kappaLayers_(psf.kappaLayers_),
353  kappaLayer_(psf.kappaLayer_.clone(patch().patch())),
354  logInterval_(psf.logInterval_),
355  executionIndex_(psf.executionIndex_),
356  thermalInertia_(psf.thermalInertia_),
357  verbose_(psf.verbose_),
358  prefix_(psf.prefix_)
359 {}
360 
361 
364 (
366 )
367 :
368  mixedFvPatchScalarField(psf),
370  mappedPatchFieldBase<scalar>
371  (
372  mappedPatchFieldBase<scalar>::mapper(patch(), psf.internalField()),
373  *this,
374  psf
375  ),
376  functionObjects::writeFile(psf),
377  TnbrName_(psf.TnbrName_),
378  qrNbrName_(psf.qrNbrName_),
379  qrName_(psf.qrName_),
380  thicknessLayers_(psf.thicknessLayers_),
381  thicknessLayer_(psf.thicknessLayer_.clone(patch().patch())),
382  kappaLayers_(psf.kappaLayers_),
383  kappaLayer_(psf.kappaLayer_.clone(patch().patch())),
384  logInterval_(psf.logInterval_),
385  executionIndex_(psf.executionIndex_),
386  thermalInertia_(psf.thermalInertia_),
387  verbose_(psf.verbose_),
388  prefix_(psf.prefix_)
389 {}
390 
391 
392 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
393 
395 (
396  const fvPatchFieldMapper& mapper
397 )
398 {
399  mixedFvPatchScalarField::autoMap(mapper);
401  //mappedPatchFieldBase<scalar>::autoMap(mapper);
402  if (thicknessLayer_)
403  {
404  thicknessLayer_().autoMap(mapper);
405  kappaLayer_().autoMap(mapper);
406  }
407 }
408 
409 
411 (
412  const fvPatchField<scalar>& ptf,
413  const labelList& addr
414 )
415 {
416  mixedFvPatchScalarField::rmap(ptf, addr);
417 
419  refCast
420  <
422  >(ptf);
423 
424  temperatureCoupledBase::rmap(tiptf, addr);
425  //mappedPatchFieldBase<scalar>::rmap(ptf, addr);
426  if (thicknessLayer_)
427  {
428  thicknessLayer_().rmap(tiptf.thicknessLayer_(), addr);
429  kappaLayer_().rmap(tiptf.kappaLayer_(), addr);
430  }
431 }
432 
433 
436 (
437  const scalarField& Tp
438 ) const
439 {
440  // Get kappa from relevant thermo
442 
443  return tk;
444 }
445 
446 
448 {
449  if (updated())
450  {
451  return;
452  }
453 
454  const polyMesh& mesh = patch().boundaryMesh().mesh();
455 
456  // Since we're inside initEvaluate/evaluate there might be processor
457  // comms underway. Change the tag we use.
458  const int oldTag = UPstream::incrMsgType();
459 
460  // Get the coupling information from the mappedPatchBase
461  const label patchi = patch().index();
462  const mappedPatchBase& mpp =
464  (
465  patch(),
466  this->internalField()
467  );
468 
469  const scalarField Tc(patchInternalField());
470  const scalarField& Tp = *this;
471 
472  const scalarField kappaTp(kappa(Tp));
473  const scalarField KDelta(kappaTp*patch().deltaCoeffs());
474 
475 
476  scalarField TcNbr;
477  scalarField TpNbr;
478  scalarField KDeltaNbr;
479 
480  if (mpp.sameWorld())
481  {
482  const polyMesh& nbrMesh = mpp.sampleMesh();
483  const label samplePatchi = mpp.samplePolyPatch().index();
484  const fvPatch& nbrPatch =
485  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
486 
487  const auto& nbrField = refCast
489  (
490  nbrPatch.lookupPatchField<volScalarField>(TnbrName_)
491  );
492 
493  // Swap to obtain full local values of neighbour K*delta
494  TcNbr = nbrField.patchInternalField();
495  TpNbr = nbrField;
496  KDeltaNbr = nbrField.kappa(nbrField)*nbrPatch.deltaCoeffs();
497  }
498  else
499  {
500  // Different world so use my region,patch. Distribution below will
501  // do the reordering.
502  TcNbr = patchInternalField();
503  TpNbr = Tp;
504  KDeltaNbr = KDelta;
505  }
506  distribute(this->internalField().name() + "_value", TcNbr);
507  distribute(this->internalField().name() + "_patchValue", TpNbr);
508  distribute(this->internalField().name() + "_weights", KDeltaNbr);
509 
510  scalarField KDeltaC(this->size(), GREAT);
511  if (thicknessLayer_ || thicknessLayers_.size())
512  {
513  // Harmonic averaging
514  {
515  KDeltaC = 0.0;
516 
517  if (thicknessLayer_)
518  {
519  const scalar t = db().time().timeOutputValue();
520  KDeltaC +=
521  thicknessLayer_().value(t)
522  /kappaLayer_().value(t);
523 
524  }
525  if (thicknessLayers_.size())
526  {
527  forAll(thicknessLayers_, iLayer)
528  {
529  KDeltaC += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
530  }
531  }
532  KDeltaC = 1.0/(KDeltaC + SMALL);
533  }
534  }
535 
536  scalarField alpha(kappaTp*(1 + KDeltaNbr/KDeltaC)*patch().deltaCoeffs());
537 
538  scalarField qr(Tp.size(), Zero);
539  if (qrName_ != "none")
540  {
541  qr = patch().lookupPatchField<volScalarField>(qrName_);
542  }
543 
544  scalarField qrNbr(Tp.size(), Zero);
545  if (qrNbrName_ != "none")
546  {
547  if (mpp.sameWorld())
548  {
549  const polyMesh& nbrMesh = mpp.sampleMesh();
550  const label samplePatchi = mpp.samplePolyPatch().index();
551  const fvPatch& nbrPatch =
552  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
553  qrNbr = nbrPatch.lookupPatchField<volScalarField>(qrNbrName_);
554  }
555  else
556  {
557  qrNbr = patch().lookupPatchField<volScalarField>(qrNbrName_);
558  }
559  distribute(qrNbrName_, qrNbr);
560  }
561 
562  // inertia therm
563  if (thermalInertia_ && !mpp.sameWorld())
564  {
566  << "thermalInertia not supported in combination with multi-world"
567  << exit(FatalError);
568  }
569  if (thermalInertia_)
570  {
571  const scalar dt = mesh.time().deltaTValue();
572  scalarField mCpDtNbr;
573 
574  {
575  const polyMesh& nbrMesh = mpp.sampleMesh();
576 
577  const basicThermo* thermo =
578  nbrMesh.findObject<basicThermo>(basicThermo::dictName);
579 
580  if (thermo)
581  {
582  const label samplePatchi = mpp.samplePolyPatch().index();
583  const fvPatch& nbrPatch =
584  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
585  const scalarField& ppn =
586  thermo->p().boundaryField()[samplePatchi];
587  const scalarField& Tpn =
588  thermo->T().boundaryField()[samplePatchi];
589 
590  mCpDtNbr =
591  (
592  thermo->Cp(ppn, Tpn, samplePatchi)
593  * thermo->rho(samplePatchi)
594  / nbrPatch.deltaCoeffs()/dt
595  );
596 
597  mpp.distribute(mCpDtNbr);
598  }
599  else
600  {
601  mCpDtNbr.setSize(Tp.size(), Zero);
602  }
603  }
604 
605  scalarField mCpDt;
606 
607  // Local inertia therm
608  {
609  const basicThermo* thermo =
610  mesh.findObject<basicThermo>(basicThermo::dictName);
611 
612  if (thermo)
613  {
614  const scalarField& pp = thermo->p().boundaryField()[patchi];
615 
616  mCpDt =
617  (
618  thermo->Cp(pp, Tp, patchi)
619  * thermo->rho(patchi)
620  / patch().deltaCoeffs()/dt
621  );
622  }
623  else
624  {
625  // Issue warning?
626  mCpDt.setSize(Tp.size(), Zero);
627  }
628  }
629 
630  const volScalarField& T =
631  this->db().lookupObject<volScalarField>
632  (
633  this->internalField().name()
634  );
635 
636  const fvPatchField<scalar>& TpOld = T.oldTime().boundaryField()[patchi];
637 
638  scalarField alpha(KDeltaNbr + mCpDt + mCpDtNbr);
639 
640  valueFraction() = alpha/(alpha + KDelta);
641  scalarField c(KDeltaNbr*TcNbr + (mCpDt + mCpDtNbr)*TpOld);
642  refValue() = c/alpha;
643  refGrad() = (qr + qrNbr)/kappaTp;
644  }
645  else
646  {
647 
648  valueFraction() = KDeltaNbr/(KDeltaNbr + alpha);
649  refValue() = TcNbr;
650  refGrad() = (qr + qrNbr)/kappaTp;
651  }
652 
653  source() = Zero;
654 
655  // If useImplicit is true we need the source term associated with this BC
656  if (this->useImplicit())
657  {
658  source() =
659  alphaSfDelta()*
660  (
661  valueFraction()*deltaH()
662  + (qr + qrNbr)/beta()
663  );
664  }
665 
666  mixedFvPatchScalarField::updateCoeffs();
667 
668 
669  if (verbose_)
670  {
671  // Calculate heat-transfer rate and heat flux
672  const scalar Q = gSum(kappaTp*patch().magSf()*snGrad());
673  const scalar magSf = gSum(patch().magSf());
674  const scalar q = Q/max(magSf, SMALL);
675 
676 
677  // Calculate heat-transfer coeff based on the first definition
678  // [W/m^2] = [W/m/K K * 1/m]
679  const scalarField qField
680  (
681  kappaTp*snGrad()
682  );
683  const scalarField deltaT(TcNbr - Tc);
684  scalarField htc(deltaT.size(), Zero);
685 
686  forAll(deltaT, i)
687  {
688  if (mag(deltaT[i]) > SMALL)
689  {
690  htc[i] = qField[i]/deltaT[i];
691  }
692  }
693  const scalar aveHtc = gSum(htc*patch().magSf())/max(magSf, SMALL);
694 
695 
696  // Calculate heat-transfer coeff based on the second definition
697  const scalarField deltaTPatch(TpNbr - Tp);
698  scalarField htcPatch(deltaTPatch.size(), Zero);
699 
700  forAll(deltaTPatch, i)
701  {
702  if (mag(deltaTPatch[i]) > SMALL)
703  {
704  htcPatch[i] = qField[i]/deltaTPatch[i];
705  }
706  }
707  const scalar aveHtcPatch =
708  gSum(htcPatch*patch().magSf())/max(magSf, SMALL);
709 
710  // Calculate various averages of temperature
711  const scalarMinMax TpMinMax = gMinMax(Tp);
712  const scalar TpAvg = gAverage(Tp);
713  const scalarMinMax TpNbrMinMax = gMinMax(TpNbr);
714  const scalar TpNbrAvg = gAverage(TpNbr);
715 
716 
717  Info<< nl
718  << patch().boundaryMesh().mesh().name() << ':'
719  << patch().name() << ':'
720  << this->internalField().name() << " <- "
721  << mpp.sampleRegion() << ':'
722  << mpp.samplePatch() << ':'
723  << this->internalField().name() << " :" << nl
724  << " Heat transfer rate [W]:" << Q << nl
725  << " Area [m^2]:" << magSf << nl
726  << " Heat flux [W/m^2]:" << q << nl
727  << " Area-averaged heat-transfer coefficient [W/m^2/K]:"
728  << aveHtc << nl
729  << " Area-averaged patch heat-transfer coefficient [W/m^2/K]:"
730  << aveHtcPatch << nl
731  << " Wall temperature [K]"
732  << " min:" << TpMinMax.min()
733  << " max:" << TpMinMax.max()
734  << " avg:" << TpAvg << nl
735  << " Neighbour wall temperature [K]"
736  << " min:" << TpNbrMinMax.min()
737  << " max:" << TpNbrMinMax.max()
738  << " avg:" << TpNbrAvg
739  << nl << endl;
740 
741 
742  // Handle data for file output
743  if (canResetFile())
744  {
745  resetFile(patch().name());
746  }
747 
748  if (canWriteHeader())
749  {
750  writeFileHeader(file());
751  }
752 
753  if (canWriteToFile() && writeFile())
754  {
755  file()
756  << db().time().timeOutputValue() << token::TAB
757  << Q << token::TAB
758  << q << token::TAB
759  << aveHtc << token::TAB
760  << aveHtcPatch << token::TAB
761  << TpMinMax.min() << token::TAB
762  << TpMinMax.max() << token::TAB
763  << TpAvg << token::TAB
764  << TpNbrMinMax.min() << token::TAB
765  << TpNbrMinMax.max() << token::TAB
766  << TpNbrAvg << token::TAB
767  << endl;
768  }
769 
770 
771  // Store htc fields as patch fields of a volScalarField
772  storeHTCFields(prefix_, htc, htcPatch);
773  }
774 
775  UPstream::msgType(oldTag); // Restore tag
776 }
777 
778 
780 (
781  fvMatrix<scalar>& m,
782  const label iMatrix,
783  const direction cmpt
784 )
785 {
787  << "This T BC does not support energy coupling "
788  << "It is implemented on he field "
789  << abort(FatalError);
790 }
791 
792 
793 tmp<Field<scalar>> turbulentTemperatureRadCoupledMixedFvPatchScalarField::coeffs
794 (
795  fvMatrix<scalar>& matrix,
796  const Field<scalar>& coeffs,
797  const label mat
798 ) const
799 {
801  << "This BC does not support energy coupling "
802  << "Use compressible::turbulentTemperatureRadCoupledMixed "
803  << "which has more functionalities and it can handle "
804  << "the assemble coupled option for energy. "
805  << abort(FatalError);
806 
807  return nullptr;
808 }
809 
810 
811 tmp<scalarField>
812 turbulentTemperatureRadCoupledMixedFvPatchScalarField::alphaSfDelta() const
813 {
814  return (alpha(*this)*patch().deltaCoeffs()*patch().magSf());
815 }
816 
817 
818 tmp<scalarField> turbulentTemperatureRadCoupledMixedFvPatchScalarField::
819 beta() const
820 {
821  const mappedPatchBase& mpp =
822  refCast<const mappedPatchBase>(patch().patch());
823 
824  if (!mpp.sameWorld())
825  {
827  << "coupled energy not supported in combination with multi-world"
828  << exit(FatalError);
829  }
830 
831  const label samplePatchi = mpp.samplePolyPatch().index();
832  const polyMesh& nbrMesh = mpp.sampleMesh();
833 
834  const fvPatch& nbrPatch =
835  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
836 
838  nbrField = refCast
840  (
841  nbrPatch.lookupPatchField<volScalarField>(TnbrName_)
842  );
843 
844  // Swap to obtain full local values of neighbour internal field
845  scalarField TcNbr(nbrField.patchInternalField());
846  mpp.distribute(TcNbr);
847 
848  scalarField alphaDeltaNbr
849  (
850  nbrField.alpha(TcNbr)*nbrPatch.deltaCoeffs()
851  );
852  mpp.distribute(alphaDeltaNbr);
853 
854  scalarField alphaDelta
855  (
856  alpha(*this)*patch().deltaCoeffs()
857  );
858 
859  return (alphaDeltaNbr + alphaDelta);
860 }
861 
862 
863 tmp<scalarField> turbulentTemperatureRadCoupledMixedFvPatchScalarField::
864 deltaH() const
865 {
866  const mappedPatchBase& mpp =
867  refCast<const mappedPatchBase>(patch().patch());
868 
869  if (!mpp.sameWorld())
870  {
872  << "coupled energy not supported in combination with multi-world"
873  << exit(FatalError);
874  }
875 
876  const polyMesh& nbrMesh = mpp.sampleMesh();
877 
878  const basicThermo* nbrThermo =
879  nbrMesh.cfindObject<basicThermo>(basicThermo::dictName);
880 
881  const polyMesh& mesh = patch().boundaryMesh().mesh();
882 
883  const basicThermo* localThermo =
884  mesh.cfindObject<basicThermo>(basicThermo::dictName);
885 
886 
887  if (nbrThermo && localThermo)
888  {
889  const label patchi = patch().index();
890  const scalarField& pp = localThermo->p().boundaryField()[patchi];
891  const scalarField& Tp = *this;
892 
893  const mappedPatchBase& mpp =
894  refCast<const mappedPatchBase>(patch().patch());
895 
896  const label patchiNrb = mpp.samplePolyPatch().index();
897 
898  const scalarField& ppNbr = nbrThermo->p().boundaryField()[patchiNrb];
899  //const scalarField& TpNbr = nbrThermo->T().boundaryField()[patchiNrb];
900 
901  // Use this Tp to evaluate he jump as this is updated while doing
902  // updateCoeffs on boundary fields which set T values on patches
903  // then non consistent Tp and Tpnbr could be used from different
904  // updated values (specially when T changes drastically between time
905  // steps/
906  return
907  (
908  - localThermo->he(pp, Tp, patchi)
909  + nbrThermo->he(ppNbr, Tp, patchiNrb)
910  );
911  }
912  else
913  {
915  << "Can't find thermos on mapped patch "
916  << " method, but thermo package not available"
917  << exit(FatalError);
918  }
919 
920  return tmp<scalarField>::New(patch().size(), Zero);
921 }
922 
923 
924 bool turbulentTemperatureRadCoupledMixedFvPatchScalarField::writeFile()
925 {
926  if (!verbose_ || (logInterval_ <= 0))
927  {
928  return false;
929  }
930 
931  const auto& time = patch().boundaryMesh().mesh().time();
932 
933  const scalar t = time.timeOutputValue();
934  const scalar ts = time.startTime().value();
935  const scalar deltaT = time.deltaTValue();
936 
937  const label executionIndex = label
938  (
939  (
940  (t - ts)
941  + 0.5*deltaT
942  )
943  /logInterval_
944  );
945 
946  bool write = false;
947  if (executionIndex > executionIndex_)
948  {
949  executionIndex_ = executionIndex;
950  write = true;
951  }
952 
953  return write;
954 }
955 
956 
958 (
959  Ostream& os
960 ) const
961 {
963 
964  os.writeEntryIfDifferent<word>("Tnbr", "T", TnbrName_);
965  os.writeEntryIfDifferent<word>("qrNbr", "none", qrNbrName_);
966  os.writeEntryIfDifferent<word>("qr", "none", qrName_);
967  os.writeEntry<scalar>("logInterval", logInterval_);
968 
969  if (thermalInertia_)
970  {
971  os.writeEntry("thermalInertia", thermalInertia_);
972  }
973  os.writeEntryIfDifferent<bool>("verbose", false, verbose_);
974  os.writeEntryIfDifferent<word>("prefix", "multiWorld", prefix_);
975 
976  if (thicknessLayer_)
977  {
978  thicknessLayer_().writeData(os);
979  kappaLayer_().writeData(os);
980  }
981  if (thicknessLayers_.size())
982  {
983  thicknessLayers_.writeEntry("thicknessLayers", os);
984  kappaLayers_.writeEntry("kappaLayers", os);
985  }
986 
987  // Write writeFile entries
988  os.writeEntry<label>("writePrecision", writePrecision_);
989  os.writeEntry<bool>("updateHeader", updateHeader_);
990  os.writeEntry<bool>("writeToFile", writeToFile_);
991  os.writeEntry<bool>("useUserTime", useUserTime_);
992 
995 }
996 
997 
998 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
999 
1001 (
1003  turbulentTemperatureRadCoupledMixedFvPatchScalarField
1004 );
1005 
1006 
1007 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1008 
1009 } // End namespace compressible
1010 } // End namespace Foam
1011 
1012 
1013 // ************************************************************************* //
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.
virtual OFstream & file()
Return access to the file (if only 1)
Definition: writeFile.C:264
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
uint8_t direction
Definition: direction.H:46
virtual bool canWriteToFile() const
Flag to allow writing to the file.
Definition: writeFile.C:287
MinMax< scalar > scalarMinMax
A scalar min/max range.
Definition: MinMax.H:97
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition: UPstream.H:1274
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
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:129
static const mappedPatchBase & mapper(const fvPatch &p, const DimensionedField< scalar, volMesh > &iF)
Check that patch is of correct type.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
Tab [isspace].
Definition: token.H:129
const polyMesh & sampleMesh() const
Get the region mesh.
virtual void write(Ostream &os) const
Write.
virtual void resetFile(const word &name)
Reset internal file pointer to new file with new name.
Definition: writeFile.C:165
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition: typeInfo.H:172
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:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
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, suitable for multiregion cases...
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
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:1252
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.
bool writtenHeader_
Flag to identify whether the header has been written.
Definition: writeFile.H:157
const dimensionSet dimless
Dimensionless.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
psiReactionThermo & thermo
Definition: createFields.H:28
virtual bool canResetFile() const
Flag to allow resetting the file.
Definition: writeFile.C:293
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition: writeFile.C:313
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
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 expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual bool canWriteHeader() const
Flag to allow writing the header.
Definition: writeFile.C:299
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:206
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
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
label writePrecision_
Write precision.
Definition: writeFile.H:141
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:56
bool compressible
Definition: pEqn.H:2
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
OBJstream os(runTime.globalPath()/outputName)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Functionality for sampling fields using mappedPatchBase. Every call to mappedField() returns a sample...
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
bool writeToFile_
Flag to enable/disable writing to file.
Definition: writeFile.H:146
Common functions used in temperature coupled boundaries.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Type gAverage(const FieldField< Field, Type > &f)
const dimensionedScalar c
Speed of light in a vacuum.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Nothing to be read.
Automatically write from objectRegistry::writeObject()
bool useUserTime_
Flag to use the specified user time, e.g. CA deg instead of seconds. Default = true.
Definition: writeFile.H:163
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.
bool updateHeader_
Flag to update the header, e.g. on mesh changes. Default is true.
Definition: writeFile.H:152
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:30
virtual void write(Ostream &) const
Write.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
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.
Request registration (bool: true)
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
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition: writeFile.C:329