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