turbulentTemperatureRadCoupledMixedFvPatchScalarField.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2017-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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"
36
37// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38
39namespace Foam
40{
41namespace compressible
42{
43
44// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45
46turbulentTemperatureRadCoupledMixedFvPatchScalarField::
47turbulentTemperatureRadCoupledMixedFvPatchScalarField
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 (
64 mappedPatchFieldBase<scalar>::mapper(p, iF),
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
79turbulentTemperatureRadCoupledMixedFvPatchScalarField::
80turbulentTemperatureRadCoupledMixedFvPatchScalarField
81(
83 const fvPatch& p,
85 const fvPatchFieldMapper& mapper
86)
87:
88 mixedFvPatchScalarField(psf, p, iF, mapper),
89 temperatureCoupledBase(patch(), psf),
91 (
92 mappedPatchFieldBase<scalar>::mapper(p, iF),
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
107turbulentTemperatureRadCoupledMixedFvPatchScalarField::
108turbulentTemperatureRadCoupledMixedFvPatchScalarField
109(
110 const fvPatch& p,
112 const dictionary& dict
113)
114:
115 mixedFvPatchScalarField(p, iF),
118 (
119 mappedPatchFieldBase<scalar>::mapper(p, iF),
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 layer as PatchFunction1
186 (
187 p.patch(),
188 "thicknessLayer",
189 dict
190 );
191 if (thicknessLayer_)
192 {
193 kappaLayer_ = PatchFunction1<scalar>::New
194 (
195 p.patch(),
196 "kappaLayer",
197 dict
198 );
199 }
200
202
203 if (dict.found("refValue"))
204 {
205 // Full restart
206 refValue() = scalarField("refValue", dict, p.size());
207 refGrad() = scalarField("refGradient", dict, p.size());
208 valueFraction() = scalarField("valueFraction", dict, p.size());
209 }
210 else
211 {
212 // Start from user entered data. Assume fixedValue.
213 refValue() = *this;
214 refGrad() = 0.0;
215 valueFraction() = 1.0;
216 }
217
218 bool boolVal(false);
219 if (dict.readIfPresent("useImplicit", boolVal))
220 {
221 this->useImplicit(boolVal);
222 }
223 if (dict.found("source"))
224 {
225 source() = scalarField("source", dict, p.size());
226 }
227 else
228 {
229 source() = 0.0;
230 }
231}
232
233
234turbulentTemperatureRadCoupledMixedFvPatchScalarField::
235turbulentTemperatureRadCoupledMixedFvPatchScalarField
236(
239)
240:
241 mixedFvPatchScalarField(psf, iF),
242 temperatureCoupledBase(patch(), psf),
244 (
245 mappedPatchFieldBase<scalar>::mapper(patch(), iF),
246 *this,
247 psf
248 ),
249 TnbrName_(psf.TnbrName_),
250 qrNbrName_(psf.qrNbrName_),
251 qrName_(psf.qrName_),
252 thicknessLayers_(psf.thicknessLayers_),
253 thicknessLayer_(psf.thicknessLayer_.clone(patch().patch())),
254 kappaLayers_(psf.kappaLayers_),
255 kappaLayer_(psf.kappaLayer_.clone(patch().patch())),
256 thermalInertia_(psf.thermalInertia_)
257{}
258
259
260turbulentTemperatureRadCoupledMixedFvPatchScalarField::
261turbulentTemperatureRadCoupledMixedFvPatchScalarField
262(
264)
265:
266 mixedFvPatchScalarField(psf),
267 temperatureCoupledBase(patch(), psf),
269 (
270 mappedPatchFieldBase<scalar>::mapper(patch(), psf.internalField()),
271 *this,
272 psf
273 ),
274 TnbrName_(psf.TnbrName_),
275 qrNbrName_(psf.qrNbrName_),
276 qrName_(psf.qrName_),
277 thicknessLayers_(psf.thicknessLayers_),
278 thicknessLayer_(psf.thicknessLayer_.clone(patch().patch())),
279 kappaLayers_(psf.kappaLayers_),
280 kappaLayer_(psf.kappaLayer_.clone(patch().patch())),
281 thermalInertia_(psf.thermalInertia_)
282{}
283
284
285// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
286
288(
289 const fvPatchFieldMapper& mapper
290)
291{
292 mixedFvPatchScalarField::autoMap(mapper);
294 //mappedPatchFieldBase<scalar>::autoMap(mapper);
295 if (thicknessLayer_)
296 {
297 thicknessLayer_().autoMap(mapper);
298 kappaLayer_().autoMap(mapper);
299 }
300}
301
302
304(
305 const fvPatchField<scalar>& ptf,
306 const labelList& addr
307)
308{
309 mixedFvPatchScalarField::rmap(ptf, addr);
310
312 refCast
313 <
315 >(ptf);
316
317 temperatureCoupledBase::rmap(tiptf, addr);
318 //mappedPatchFieldBase<scalar>::rmap(ptf, addr);
319 if (thicknessLayer_)
320 {
321 thicknessLayer_().rmap(tiptf.thicknessLayer_(), addr);
322 kappaLayer_().rmap(tiptf.kappaLayer_(), addr);
323 }
324}
325
326
329(
330 const scalarField& Tp
331) const
332{
333 // Get kappa from relevant thermo
335
336 return tk;
337}
338
339
341{
342 if (updated())
343 {
344 return;
345 }
346
347 const polyMesh& mesh = patch().boundaryMesh().mesh();
348
349 // Since we're inside initEvaluate/evaluate there might be processor
350 // comms underway. Change the tag we use.
351 int oldTag = UPstream::msgType();
352 UPstream::msgType() = oldTag+1;
353
354 // Get the coupling information from the mappedPatchBase
355 const label patchi = patch().index();
356 const mappedPatchBase& mpp =
358 (
359 patch(),
360 this->internalField()
361 );
362
363 const scalarField Tc(patchInternalField());
364 const scalarField& Tp = *this;
365
366 const scalarField kappaTp(kappa(Tp));
367 const scalarField KDelta(kappaTp*patch().deltaCoeffs());
368
369
370 scalarField TcNbr;
371 scalarField KDeltaNbr;
372
373 if (mpp.sameWorld())
374 {
375 const polyMesh& nbrMesh = mpp.sampleMesh();
376 const label samplePatchi = mpp.samplePolyPatch().index();
377 const fvPatch& nbrPatch =
378 refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
379
380 const auto& nbrField = refCast
382 (
383 nbrPatch.lookupPatchField<volScalarField, scalar>(TnbrName_)
384 );
385
386 // Swap to obtain full local values of neighbour K*delta
387 TcNbr = nbrField.patchInternalField();
388 KDeltaNbr = nbrField.kappa(nbrField)*nbrPatch.deltaCoeffs();
389 }
390 else
391 {
392 // Different world so use my region,patch. Distribution below will
393 // do the reordering.
394 TcNbr = patchInternalField();
395 KDeltaNbr = KDelta;
396 }
397 distribute(this->internalField().name() + "_value", TcNbr);
398 distribute(this->internalField().name() + "_weights", KDeltaNbr);
399
400 scalarField KDeltaC(this->size(), GREAT);
401 if (thicknessLayer_ || thicknessLayers_.size())
402 {
403 // Harmonic averaging
404 {
405 KDeltaC = 0.0;
406
407 if (thicknessLayer_)
408 {
409 const scalar t = db().time().timeOutputValue();
410 KDeltaC +=
411 thicknessLayer_().value(t)
412 /kappaLayer_().value(t);
413
414 }
415 if (thicknessLayers_.size())
416 {
417 forAll(thicknessLayers_, iLayer)
418 {
419 KDeltaC += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
420 }
421 }
422 KDeltaC = 1.0/(KDeltaC + SMALL);
423 }
424 }
425
426 scalarField alpha(kappaTp*(1 + KDeltaNbr/KDeltaC)*patch().deltaCoeffs());
427
428 scalarField qr(Tp.size(), Zero);
429 if (qrName_ != "none")
430 {
431 qr = patch().lookupPatchField<volScalarField, scalar>(qrName_);
432 }
433
434 scalarField qrNbr(Tp.size(), Zero);
435 if (qrNbrName_ != "none")
436 {
437 if (mpp.sameWorld())
438 {
439 const polyMesh& nbrMesh = mpp.sampleMesh();
440 const label samplePatchi = mpp.samplePolyPatch().index();
441 const fvPatch& nbrPatch =
442 refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
443 qrNbr =
444 nbrPatch.lookupPatchField<volScalarField, scalar>(qrNbrName_);
445 }
446 else
447 {
448 qrNbr =
449 patch().lookupPatchField<volScalarField, scalar>(qrNbrName_);
450 }
451 distribute(qrNbrName_, qrNbr);
452 }
453
454 // inertia therm
455 if (thermalInertia_ && !mpp.sameWorld())
456 {
458 << "thermalInertia not supported in combination with multi-world"
459 << exit(FatalError);
460 }
461 if (thermalInertia_)
462 {
463 const scalar dt = mesh.time().deltaTValue();
464 scalarField mCpDtNbr;
465
466 {
467 const polyMesh& nbrMesh = mpp.sampleMesh();
468
469 const basicThermo* thermo =
471
472 if (thermo)
473 {
474 const label samplePatchi = mpp.samplePolyPatch().index();
475 const fvPatch& nbrPatch =
476 refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
477 const scalarField& ppn =
478 thermo->p().boundaryField()[samplePatchi];
479 const scalarField& Tpn =
480 thermo->T().boundaryField()[samplePatchi];
481
482 mCpDtNbr =
483 (
484 thermo->Cp(ppn, Tpn, samplePatchi)
485 * thermo->rho(samplePatchi)
486 / nbrPatch.deltaCoeffs()/dt
487 );
488
489 mpp.distribute(mCpDtNbr);
490 }
491 else
492 {
493 mCpDtNbr.setSize(Tp.size(), Zero);
494 }
495 }
496
497 scalarField mCpDt;
498
499 // Local inertia therm
500 {
501 const basicThermo* thermo =
503
504 if (thermo)
505 {
506 const scalarField& pp = thermo->p().boundaryField()[patchi];
507
508 mCpDt =
509 (
510 thermo->Cp(pp, Tp, patchi)
511 * thermo->rho(patchi)
512 / patch().deltaCoeffs()/dt
513 );
514 }
515 else
516 {
517 // Issue warning?
518 mCpDt.setSize(Tp.size(), Zero);
519 }
520 }
521
522 const volScalarField& T =
523 this->db().lookupObject<volScalarField>
524 (
525 this->internalField().name()
526 );
527
528 const fvPatchField<scalar>& TpOld = T.oldTime().boundaryField()[patchi];
529
530 scalarField alpha(KDeltaNbr + mCpDt + mCpDtNbr);
531
532 valueFraction() = alpha/(alpha + KDelta);
533 scalarField c(KDeltaNbr*TcNbr + (mCpDt + mCpDtNbr)*TpOld);
534 refValue() = c/alpha;
535 refGrad() = (qr + qrNbr)/kappaTp;
536 }
537 else
538 {
539
540 valueFraction() = KDeltaNbr/(KDeltaNbr + alpha);
541 refValue() = TcNbr;
542 refGrad() = (qr + qrNbr)/kappaTp;
543 }
544
545 source() = Zero;
546
547 // If useImplicit is true we need the source term associated with this BC
548 if (this->useImplicit())
549 {
550 source() =
551 alphaSfDelta()*
552 (
553 valueFraction()*deltaH()
554 + (qr + qrNbr)/beta()
555 );
556 }
557
558 mixedFvPatchScalarField::updateCoeffs();
559
560 if (debug)
561 {
562 scalar Q = gSum(kappaTp*patch().magSf()*snGrad());
563
564 Info<< patch().boundaryMesh().mesh().name() << ':'
565 << patch().name() << ':'
566 << this->internalField().name() << " <- "
567 << mpp.sampleRegion() << ':'
568 << mpp.samplePatch() << ':'
569 << this->internalField().name() << " :"
570 << " heat transfer rate:" << Q
571 << " walltemperature "
572 << " min:" << gMin(Tp)
573 << " max:" << gMax(Tp)
574 << " avg:" << gAverage(Tp)
575 << endl;
576 }
577
578 // Restore tag
579 UPstream::msgType() = oldTag;
580}
581
582
584(
586 const label iMatrix,
587 const direction cmpt
588)
589{
591 << "This T BC does not support energy coupling "
592 << "It is implemented on he field "
593 << abort(FatalError);
594}
595
596
598(
599 fvMatrix<scalar>& matrix,
600 const Field<scalar>& coeffs,
601 const label mat
602) const
603{
605 << "This BC does not support energy coupling "
606 << "Use compressible::turbulentTemperatureRadCoupledMixed "
607 << "which has more functionalities and it can handle "
608 << "the assemble coupled option for energy. "
609 << abort(FatalError);
610
611 return tmp<Field<scalar>>(new Field<scalar>());
612}
613
614
616turbulentTemperatureRadCoupledMixedFvPatchScalarField::alphaSfDelta() const
617{
618 return (alpha(*this)*patch().deltaCoeffs()*patch().magSf());
619}
620
621
622tmp<scalarField> turbulentTemperatureRadCoupledMixedFvPatchScalarField::
623beta() const
624{
625 const mappedPatchBase& mpp =
626 refCast<const mappedPatchBase>(patch().patch());
627
628 if (!mpp.sameWorld())
629 {
631 << "coupled energy not supported in combination with multi-world"
632 << exit(FatalError);
633 }
634
635 const label samplePatchi = mpp.samplePolyPatch().index();
636 const polyMesh& nbrMesh = mpp.sampleMesh();
637
638 const fvPatch& nbrPatch =
639 refCast<const fvMesh>(nbrMesh).boundary()[samplePatchi];
640
642 nbrField = refCast
644 (
645 nbrPatch.lookupPatchField<volScalarField, scalar>(TnbrName_)
646 );
647
648 // Swap to obtain full local values of neighbour internal field
649 scalarField TcNbr(nbrField.patchInternalField());
650 mpp.distribute(TcNbr);
651
652 scalarField alphaDeltaNbr
653 (
654 nbrField.alpha(TcNbr)*nbrPatch.deltaCoeffs()
655 );
656 mpp.distribute(alphaDeltaNbr);
657
658 scalarField alphaDelta
659 (
660 alpha(*this)*patch().deltaCoeffs()
661 );
662
663 return (alphaDeltaNbr + alphaDelta);
664}
665
666
667tmp<scalarField> turbulentTemperatureRadCoupledMixedFvPatchScalarField::
668deltaH() const
669{
670 const mappedPatchBase& mpp =
671 refCast<const mappedPatchBase>(patch().patch());
672
673 if (!mpp.sameWorld())
674 {
676 << "coupled energy not supported in combination with multi-world"
677 << exit(FatalError);
678 }
679
680 const polyMesh& nbrMesh = mpp.sampleMesh();
681
682 const basicThermo* nbrThermo =
683 nbrMesh.cfindObject<basicThermo>(basicThermo::dictName);
684
685 const polyMesh& mesh = patch().boundaryMesh().mesh();
686
687 const basicThermo* localThermo =
689
690
691 if (nbrThermo && localThermo)
692 {
693 const label patchi = patch().index();
694 const scalarField& pp = localThermo->p().boundaryField()[patchi];
695 const scalarField& Tp = *this;
696
697 const mappedPatchBase& mpp =
698 refCast<const mappedPatchBase>(patch().patch());
699
700 const label patchiNrb = mpp.samplePolyPatch().index();
701
702 const scalarField& ppNbr = nbrThermo->p().boundaryField()[patchiNrb];
703 //const scalarField& TpNbr = nbrThermo->T().boundaryField()[patchiNrb];
704
705 // Use this Tp to evaluate he jump as this is updated while doing
706 // updateCoeffs on boundary fields which set T values on patches
707 // then non consistent Tp and Tpnbr could be used from different
708 // updated values (specially when T changes drastically between time
709 // steps/
710 return
711 (
712 - localThermo->he(pp, Tp, patchi)
713 + nbrThermo->he(ppNbr, Tp, patchiNrb)
714 );
715 }
716 else
717 {
719 << "Can't find thermos on mapped patch "
720 << " method, but thermo package not available"
721 << exit(FatalError);
722 }
723
724 return tmp<scalarField>::New(patch().size(), Zero);
725}
726
727
729(
730 Ostream& os
731) const
732{
733 mixedFvPatchScalarField::write(os);
734
735 //os.writeEntry("Tnbr", TnbrName_);
736 os.writeEntryIfDifferent<word>("Tnbr", "T", TnbrName_);
737
738 //os.writeEntry("qrNbr", qrNbrName_);
739 os.writeEntryIfDifferent<word>("qrNbr", "none", qrNbrName_);
740 //os.writeEntry("qr", qrName_);
741 os.writeEntryIfDifferent<word>("qr", "none", qrName_);
742 if (thermalInertia_)
743 {
744 os.writeEntry("thermalInertia", thermalInertia_);
745 }
746
747 if (thicknessLayer_)
748 {
749 thicknessLayer_().writeData(os);
750 kappaLayer_().writeData(os);
751 }
752 if (thicknessLayers_.size())
753 {
754 thicknessLayers_.writeEntry("thicknessLayers", os);
755 kappaLayers_.writeEntry("kappaLayers", os);
756 }
757
760}
761
762
763// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
764
766(
769);
770
771
772// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
773
774} // End namespace compressible
775} // End namespace Foam
776
777
778// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
void setSize(const label n)
Alias for resize()
Definition: List.H:218
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:251
static autoPtr< PatchFunction1< Type > > NewIfPresent(const polyPatch &pp, const word &entryName, const dictionary &dict, const bool faceValues=true)
An optional selector.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:78
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:43
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void writeEntry(Ostream &os) const
Write the UList with its compound type.
Definition: UListIO.C:38
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:556
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:66
static const word dictName
Definition: basicThermo.H:256
const dictionary & coeffs() const
Return const dictionary of the model.
Mixed boundary condition for temperature and radiation heat transfer to be used for in multiregion ca...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void manipulateMatrix(fvMatrix< scalar > &m, const label iMatrix, const direction cmpt)
Manipulate matrix.
turbulentTemperatureRadCoupledMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
virtual bool write()
Write the output fields.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:408
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
const scalarField & deltaCoeffs() const
Definition: fvPatch.C:196
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
const polyMesh & sampleMesh() const
Get the region mesh.
const polyPatch & samplePolyPatch() const
Get the patch on the region.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
const word & samplePatch() const
Patch (only if NEARESTPATCHFACE)
const word & sampleRegion() const
Region to sample.
bool sameWorld() const
Is sample world the local world?
Functionality for sampling fields using mappedPatchBase. Every call to mappedField() returns a sample...
static const mappedPatchBase & mapper(const fvPatch &p, const DimensionedField< scalar, volMesh > &iF)
Check that patch is of correct type.
void distribute(const word &fieldName, Field< T > &newValues) const
const Type * cfindObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
const Type * findObject(const word &name, const bool recursive=false) const
Return const pointer to the object of the given Type.
label index() const noexcept
The index of this patch in the boundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
Common functions used in temperature coupled boundaries.
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
const volScalarField & T
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
bool compressible
Definition: pEqn.H:2
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
errorManip< error > abort(error &err)
Definition: errorManip.H:144
uint8_t direction
Definition: direction.H:56
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Type gMax(const FieldField< Field, Type > &f)
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
static const char *const typeName
The type name used in ensight case files.