Loading...
Searching...
No Matches
humidityTemperatureCoupledMixedFvPatchScalarField.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) 2015-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30#include "fvPatchFieldMapper.H"
31#include "volFields.H"
32#include "surfaceFields.H"
33#include "mappedPatchBase.H"
35
36// * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
37
38const Foam::Enum
39<
41>
42Foam::humidityTemperatureCoupledMixedFvPatchScalarField::massModeTypeNames_
43({
44 { massTransferMode::mtConstantMass, "constantMass" },
45 { massTransferMode::mtCondensation, "condensation" },
46 { massTransferMode::mtEvaporation, "evaporation" },
47 {
48 massTransferMode::mtCondensationAndEvaporation,
49 "condensationAndEvaporation"
50 },
51});
52
53
54// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
55
56Foam::scalar Foam::humidityTemperatureCoupledMixedFvPatchScalarField::Sh
57(
58 const scalar Re,
59 const scalar Sc
60) const
61{
62 if (Re < 5.0E+05)
63 {
64 return 0.664*sqrt(Re)*cbrt(Sc);
65 }
66 else
67 {
68 return 0.037*pow(Re, 0.8)*cbrt(Sc);
69 }
70}
71
72
73Foam::scalar
74Foam::humidityTemperatureCoupledMixedFvPatchScalarField::htcCondensation
75(
76 const scalar Tsat,
77 const scalar Re
78) const
79{
80 // (BLID:Eq. 10.52-10.53)
81 if (Tsat > 295.15 && Tsat < 373.15)
82 {
83 return 51104 + 2044*(Tsat - 273.15);
84 }
85 else
86 {
87 return 255510;
88 }
89}
90
91
93Foam::humidityTemperatureCoupledMixedFvPatchScalarField::thicknessField
94(
95 const word& fieldName,
96 const fvMesh& mesh
97)
98{
100
101 if (!ptr)
102 {
103 ptr = new volScalarField
104 (
106 (
107 fieldName,
108 mesh.time().timeName(),
109 mesh,
113 ),
114 mesh,
116 );
117
118 ptr->store();
119 }
120
121 return *ptr;
122}
123
124
125
126// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
127
130(
131 const fvPatch& p,
133)
134:
135 mixedFvPatchScalarField(p, iF),
136 temperatureCoupledBase(patch(), temperatureCoupledBase::mtFluidThermo),
137 mode_(mtConstantMass),
138 pName_("p"),
139 UName_("U"),
140 rhoName_("rho"),
141 muName_("thermo:mu"),
142 TnbrName_("T"),
143 qrNbrName_("none"),
144 qrName_("none"),
145 specieName_("none"),
146 liquid_(nullptr),
147 liquidDict_(nullptr),
148 mass_(patch().size(), Zero),
149 Tvap_(0.0),
150 myKDelta_(patch().size(), Zero),
151 dmHfg_(patch().size(), Zero),
152 mpCpTp_(patch().size(), Zero),
153 Mcomp_(0.0),
154 L_(0.0),
155 fluid_(false),
156 cp_(patch().size(), Zero),
157 thickness_(patch().size(), Zero),
158 rho_(patch().size(), Zero),
159 thicknessLayers_(0),
160 kappaLayers_(0)
162 this->refValue() = 0.0;
163 this->refGrad() = 0.0;
164 this->valueFraction() = 1.0;
165}
166
167
170(
172 const fvPatch& p,
174 const fvPatchFieldMapper& mapper
175)
176:
177 mixedFvPatchScalarField(psf, p, iF, mapper),
178 temperatureCoupledBase(patch(), psf),
179 mode_(psf.mode_),
180 pName_(psf.pName_),
181 UName_(psf.UName_),
182 rhoName_(psf.rhoName_),
183 muName_(psf.muName_),
184 TnbrName_(psf.TnbrName_),
185 qrNbrName_(psf.qrNbrName_),
186 qrName_(psf.qrName_),
187 specieName_(psf.specieName_),
188 liquid_(psf.liquid_.clone()),
189 liquidDict_(psf.liquidDict_),
190 mass_(psf.mass_, mapper),
191 Tvap_(psf.Tvap_),
192 myKDelta_(psf.myKDelta_, mapper),
193 dmHfg_(psf.dmHfg_, mapper),
194 mpCpTp_(psf.mpCpTp_, mapper),
195 Mcomp_(psf.Mcomp_),
196 L_(psf.L_),
197 fluid_(psf.fluid_),
198 cp_(psf.cp_, mapper),
199 thickness_(psf.thickness_, mapper),
200 rho_(psf.rho_, mapper),
201 thicknessLayers_(psf.thicknessLayers_),
202 kappaLayers_(psf.kappaLayers_)
203{}
204
205
208(
209 const fvPatch& p,
211 const dictionary& dict
212)
213:
214 mixedFvPatchScalarField(p, iF),
216 mode_(mtCondensationAndEvaporation),
217 pName_(dict.getOrDefault<word>("p", "p")),
218 UName_(dict.getOrDefault<word>("U", "U")),
219 rhoName_(dict.getOrDefault<word>("rho", "rho")),
220 muName_(dict.getOrDefault<word>("mu", "thermo:mu")),
221 TnbrName_(dict.getOrDefault<word>("Tnbr", "T")),
222 qrNbrName_(dict.getOrDefault<word>("qrNbr", "none")),
223 qrName_(dict.getOrDefault<word>("qr", "none")),
224 specieName_(dict.getOrDefault<word>("specie", "none")),
225 liquid_(nullptr),
226 liquidDict_(),
227 mass_(patch().size(), Zero),
228 Tvap_(0.0),
229 myKDelta_(patch().size(), Zero),
230 dmHfg_(patch().size(), Zero),
231 mpCpTp_(patch().size(), Zero),
232 Mcomp_(0.0),
233 L_(0.0),
234 fluid_(false),
235 cp_(patch().size(), Zero),
236 thickness_(patch().size(), Zero),
237 rho_(patch().size(), Zero),
238 thicknessLayers_(0),
239 kappaLayers_(0)
240{
241 if (!isA<mappedPatchBase>(this->patch().patch()))
242 {
244 << "\n patch type '" << p.type()
245 << "' not type '" << mappedPatchBase::typeName << "'"
246 << "\n for patch " << p.name()
247 << " of field " << internalField().name()
248 << " in file " << internalField().objectPath()
249 << exit(FatalIOError);
250 }
251
252 this->readValueEntry(dict, IOobjectOption::MUST_READ);
253
254 if (massModeTypeNames_.readIfPresent("mode", dict, mode_))
255 {
256 fluid_ = true;
257 }
258
259 if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
260 {
261 dict.readEntry("kappaLayers", kappaLayers_);
262 }
263
264 if (fluid_)
265 {
266 switch(mode_)
267 {
268 case mtConstantMass:
269 {
270 thickness_ = scalarField("thickness", dict, p.size());
271 cp_ = scalarField("cp", dict, p.size());
272 rho_ = scalarField("rho", dict, p.size());
273
274 break;
275 }
276 case mtCondensation:
277 case mtEvaporation:
279 {
280 dict.readEntry("carrierMolWeight", Mcomp_);
281 dict.readEntry("L", L_);
282 dict.readEntry("Tvap", Tvap_);
283 liquidDict_ = dict.subDict("liquid");
284 liquid_ =
285 liquidProperties::New(liquidDict_.subDict(specieName_));
286
287 if (dict.found("thickness"))
288 {
289 scalarField& Tp = *this;
290 const scalarField& magSf = patch().magSf();
291
292 // Assume initially standard pressure for rho calculation
293 scalar pf = 1e5;
294 thickness_ = scalarField("thickness", dict, p.size());
295 forAll(thickness_, i)
296 {
297 mass_[i] =
298 thickness_[i]*liquid_->rho(pf, Tp[i])*magSf[i];
299 }
300 }
301 fluid_ = true;
302
303 break;
304 }
305 default:
306 {
308 << "Did not find mode " << mode_
309 << " on patch " << patch().name()
310 << nl
311 << "Please set 'mode' to one of "
312 << massModeTypeNames_.sortedToc()
313 << exit(FatalIOError);
314 }
315 }
316 }
317
318
319 if (this->readMixedEntries(dict))
320 {
321 // Full restart
322 }
323 else
324 {
325 // Start from user entered data. Assume fixedValue.
326 refValue() = *this;
327 refGrad() = 0.0;
328 valueFraction() = 1.0;
329 }
330}
331
332
335(
337 const DimensionedField<scalar, volMesh>& iF
338)
339:
340 mixedFvPatchScalarField(psf, iF),
341 temperatureCoupledBase(patch(), psf),
342 mode_(psf.mode_),
343 pName_(psf.pName_),
344 UName_(psf.UName_),
345 rhoName_(psf.rhoName_),
346 muName_(psf.muName_),
347 TnbrName_(psf.TnbrName_),
348 qrNbrName_(psf.qrNbrName_),
349 qrName_(psf.qrName_),
350 specieName_(psf.specieName_),
351 liquid_(psf.liquid_.clone()),
352 liquidDict_(psf.liquidDict_),
353 mass_(psf.mass_),
354 Tvap_(psf.Tvap_),
355 myKDelta_(psf.myKDelta_),
356 dmHfg_(psf.dmHfg_),
357 mpCpTp_(psf.mpCpTp_),
358 Mcomp_(psf.Mcomp_),
359 L_(psf.L_),
360 fluid_(psf.fluid_),
361 cp_(psf.cp_),
362 thickness_(psf.thickness_),
363 rho_(psf.rho_),
364 thicknessLayers_(psf.thicknessLayers_),
365 kappaLayers_(psf.kappaLayers_)
366{}
367
368
369// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
370
372(
373 const fvPatchFieldMapper& m
374)
375{
376 mixedFvPatchScalarField::autoMap(m);
378
379 if (fluid_)
380 {
381 mass_.autoMap(m);
382 myKDelta_.autoMap(m);
383 dmHfg_.autoMap(m);
384 mpCpTp_.autoMap(m);
385 cp_.autoMap(m);
386 thickness_.autoMap(m);
387 rho_.autoMap(m);
388 }
389}
390
391
393(
394 const fvPatchScalarField& ptf,
395 const labelList& addr
396)
397{
398 mixedFvPatchScalarField::rmap(ptf, addr);
399
402 (
403 ptf
404 );
405
406 temperatureCoupledBase::rmap(tiptf, addr);
407 if (fluid_)
408 {
409 mass_.rmap(tiptf.mass_, addr);
410 myKDelta_.rmap(tiptf.myKDelta_, addr);
411 dmHfg_.rmap(tiptf.dmHfg_, addr);
412 mpCpTp_.rmap(tiptf.mpCpTp_, addr);
413 cp_.rmap(tiptf.cp_, addr);
414 thickness_.rmap(tiptf.thickness_, addr);
415 rho_.rmap(tiptf.rho_, addr);
416 }
417}
418
419
421{
422 if (updated())
423 {
424 return;
425 }
426
427 // Get the coupling information from the mappedPatchBase
428 const mappedPatchBase& mpp =
430
431 const scalarField& magSf = patch().magSf();
432
433 const label nbrPatchI = mpp.samplePolyPatch().index();
434 const polyMesh& mesh = patch().boundaryMesh().mesh();
435 const polyMesh& nbrMesh = mpp.sampleMesh();
436 const fvPatch& nbrPatch =
437 refCast<const fvMesh>(nbrMesh).boundary()[nbrPatchI];
438
439 const auto& nbrField = refCast
440 <
442 >
443 (
444 nbrPatch.lookupPatchField<volScalarField>(TnbrName_)
445 );
446
447 // Swap to obtain full local values of neighbour internal field
448 scalarField nbrIntFld(nbrField.patchInternalField());
449 mpp.distribute(nbrIntFld);
450
451
452 scalarField& Tp = *this;
453
454 const volScalarField& T =
455 static_cast<const volScalarField&>(internalField());
456
457 const scalarField TpOld(T.oldTime().boundaryField()[patch().index()]);
458
459 scalarField Tin(patchInternalField());
460
461 const scalarField K(this->kappa(*this));
462
463 // Neighbour kappa done separately because we need kappa solid for the
464 // htc correlation
465 scalarField nbrK(nbrField.kappa(*this));
466 mpp.distribute(nbrK);
467
468 scalarField nrbDeltaCoeffs(nbrPatch.deltaCoeffs());
469 mpp.distribute(nrbDeltaCoeffs);
470
471 scalarField KDeltaNbr(nbrField.kappa(*this)*nbrPatch.deltaCoeffs());
472 mpp.distribute(KDeltaNbr);
473
474 myKDelta_ = K*patch().deltaCoeffs();
475
476 if (thicknessLayers_.size() > 0)
477 {
478 myKDelta_ = 1.0/myKDelta_;
479 forAll(thicknessLayers_, iLayer)
480 {
481 myKDelta_ += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
482 }
483 myKDelta_ = 1.0/myKDelta_;
484 }
485
486 scalarField dm(patch().size(), Zero);
487
488 // Fluid Side
489 if (fluid_)
490 {
491 scalarField Yvp(patch().size(), Zero);
492 const scalar dt = mesh.time().deltaTValue();
493
494 const scalarField myDelta(patch().deltaCoeffs());
495
496 if (mode_ != mtConstantMass)
497 {
498 scalarField cp(patch().size(), Zero);
499 scalarField hfg(patch().size(), Zero);
500 scalarField htc(patch().size(), GREAT);
501 scalarField liquidRho(patch().size(), Zero);
502 scalarField Tdew(patch().size(), Zero);
503 scalarField RH(patch().size(), Zero);
504
507 (
508 refCast
509 <
511 >
512 (
513 patch().lookupPatchField<volScalarField>(specieName_)
514 )
515 );
516
517 const auto& pp =
518 patch().lookupPatchField<volScalarField>(pName_);
519
520 const auto& Up =
521 patch().lookupPatchField<volVectorField>(UName_);
522
523 const auto& rhop =
524 patch().lookupPatchField<volScalarField>(rhoName_);
525
526 const auto& mup =
527 patch().lookupPatchField<volScalarField>(muName_);
528
529 const vectorField Ui(Up.patchInternalField());
530 const scalarField Yi(Yp.patchInternalField());
531
532 forAll(Tp, faceI)
533 {
534 const scalar Tf = Tp[faceI];
535 const scalar Tint = Tin[faceI];
536 const vector Uf = Ui[faceI];
537 const scalar pf = pp[faceI];
538
539 const scalar muf = mup[faceI];
540 const scalar rhof = rhop[faceI];
541 const scalar nuf = muf/rhof;
542 const scalar pSat = liquid_->pv(pf, Tint);
543 const scalar Mv = liquid_->W();
544 const scalar TSat = liquid_->pvInvert(pSat);
545 const scalar Re = mag(Uf)*L_/nuf;
546
547 cp[faceI] = liquid_->Cp(pf, Tf);
548 hfg[faceI] = liquid_->hl(pf, Tf);
549
550 // Calculate relative humidity
551 const scalar invMwmean =
552 Yi[faceI]/Mv + (1.0 - Yi[faceI])/Mcomp_;
553 const scalar Xv = Yi[faceI]/invMwmean/Mv;
554 RH[faceI] = min(Xv*pf/pSat, 1.0);
555
556 scalar RHmin = 0.01;
557 Tdew[faceI] = 0.0;
558
559 if (RH[faceI] > RHmin)
560 {
561 scalar b = 243.5;
562 scalar c = 17.65;
563 scalar TintDeg = Tint - 273;
564 Tdew[faceI] =
565 b*(log(RH[faceI]) + (c*TintDeg)/(b + TintDeg))
566 /(c - log(RH[faceI]) - ((c*TintDeg)/(b + TintDeg)))
567 + 273;
568 }
569
570 if
571 (
572 Tf < Tdew[faceI]
573 && RH[faceI] > RHmin
574 && (
575 mode_ == mtCondensation
576 || mode_ == mtCondensationAndEvaporation
577 )
578 )
579 {
580 htc[faceI] = htcCondensation(TSat, Re)*nbrK[faceI]/L_;
581
582 scalar htcTotal =
583 1.0/((1.0/myKDelta_[faceI]) + (1.0/htc[faceI]));
584
585 // Heat flux [W] (>0 heat is converted into mass)
586 const scalar q = (Tint - Tf)*htcTotal*magSf[faceI];
587
588 // Mass flux rate [Kg/s/m2]
589 dm[faceI] = q/hfg[faceI]/magSf[faceI];
590
591 mass_[faceI] += q/hfg[faceI]*dt;
592
593 // -dYp/dn = q/Dab (fixedGradient)
594 const scalar Dab = liquid_->D(pf, Tf);
595 Yvp[faceI] =
596 -min(dm[faceI]/Dab/rhof, Yi[faceI]*myDelta[faceI]);
597 }
598 else if
599 (
600 Tf > Tvap_
601 && mass_[faceI] > 0.0
602 && (
603 mode_ == mtEvaporation
604 || mode_ == mtCondensationAndEvaporation
605 )
606 )
607 {
608 const scalar Dab = liquid_->D(pf, Tf);
609
610 const scalar Sc = nuf/Dab;
611 const scalar Sh = this->Sh(Re, Sc);
612
613 const scalar Ys = Mv*pSat/(Mv*pSat + Mcomp_*(pf - pSat));
614
615 // Mass transfer coefficient [m/s]
616 const scalar hm = Dab*Sh/L_;
617
618 const scalar Yinf = max(Yi[faceI], 0.0);
619
620 // Mass flux rate [Kg/s/m2]
621 dm[faceI] = -rhof*hm*max((Ys - Yinf), 0.0)/(1.0 - Ys);
622
623 // Set fixedGradient for carrier species.
624 Yvp[faceI] = -dm[faceI]/Dab/rhof;
625
626 // Total mass accumulated [Kg]
627 mass_[faceI] += dm[faceI]*magSf[faceI]*dt;
628
629 htc[faceI] = htcCondensation(TSat, Re)*nbrK[faceI]/L_;
630 }
631 else if (Tf > Tdew[faceI] && Tf < Tvap_ && mass_[faceI] > 0.0)
632 {
633 htc[faceI] = htcCondensation(TSat, Re)*nbrK[faceI]/L_;
634 }
635 else if (mass_[faceI] == 0.0)
636 {
637 // Do nothing
638 }
639
640 liquidRho[faceI] = liquid_->rho(pf, Tf);
641 }
642
643 mass_ = max(mass_, scalar(0));
644
645 Yp.gradient() = Yvp;
646
647 // Output film delta (e.g. H2OThickness) [m]
648 const word fieldName(specieName_ + "Thickness");
649
650 scalarField& pDelta =
651 thicknessField
652 (
653 fieldName,
655 ).boundaryFieldRef()[patch().index()];
656
657
658 pDelta = mass_/liquidRho/magSf;
659
660 // Weight myKDelta and htc
661 myKDelta_ = 1.0/((1.0/myKDelta_) + (1.0/htc));
662
663 mpCpTp_ = mass_*cp/dt/magSf;
664
665 // Heat flux due to change of phase [W/m2]
666 dmHfg_ = dm*hfg;
667
668 if (debug)
669 {
670 // Output RH and Tdew
671 scalarField& bRH =
672 thicknessField
673 (
674 "RH",
676 ).boundaryFieldRef()[patch().index()];
677 bRH = RH;
678
679 scalarField& bTdew =
680 thicknessField
681 (
682 "Tdew",
684 ).boundaryFieldRef()[patch().index()];
685 bTdew = Tdew;
686 }
687 }
688 else
689 {
690 // Inertia term [W/K/m2]
691 mpCpTp_ = thickness_*rho_*cp_/dt;
692 }
693 }
694
695 scalarField mpCpTpNbr(patch().size(), Zero);
696 scalarField dmHfgNbr(patch().size(), Zero);
697
698 if (!fluid_)
699 {
700 mpCpTpNbr = nbrField.mpCpTp();
701 mpp.distribute(mpCpTpNbr);
702
703 dmHfgNbr = nbrField.dmHfg();
704 mpp.distribute(dmHfgNbr);
705 }
706
707 // Obtain Rad heat (qr)
708 scalarField qr(Tp.size(), Zero);
709 if (qrName_ != "none")
710 {
711 qr = patch().lookupPatchField<volScalarField>(qrName_);
712 }
713
714 scalarField qrNbr(Tp.size(), Zero);
715 if (qrNbrName_ != "none")
716 {
717 qrNbr = nbrPatch.lookupPatchField<volScalarField>(qrNbrName_);
718 mpp.distribute(qrNbr);
719 }
720
721 const scalarField dmHfg(dmHfgNbr + dmHfg_);
722
723 const scalarField mpCpdt(mpCpTpNbr + mpCpTp_);
724
725 // qr > 0 (heat up the wall)
726 scalarField alpha(KDeltaNbr + mpCpdt - (qr + qrNbr)/Tp);
727
728 valueFraction() = alpha/(alpha + myKDelta_);
729
730 refValue() = (KDeltaNbr*nbrIntFld + mpCpdt*TpOld + dmHfg)/alpha;
731
732 mixedFvPatchScalarField::updateCoeffs();
733
734 if (debug && fluid_)
735 {
736 scalar Qdm = gWeightedSum(magSf, dm);
737 scalar QMass = gSum(mass_);
738 scalar Qt = gSum(myKDelta_*(Tp - Tin)*magSf);
739 scalar QtSolid = gSum(KDeltaNbr*(Tp - nbrIntFld)*magSf);
740
742
743 Info<< mesh.name() << ':'
744 << patch().name() << ':'
745 << internalField().name() << " <- "
746 << nbrMesh.name() << ':'
747 << nbrPatch.name() << ':'
748 << internalField().name() << " :" << nl
749 << " Total mass flux [Kg/s] : " << Qdm << nl
750 << " Total mass on the wall [Kg] : " << QMass << nl
751 << " Total heat (>0 leaving the wall to the fluid) [W] : "
752 << Qt << nl
753 << " Total heat (>0 leaving the wall to the solid) [W] : "
754 << QtSolid << nl
755 << " wall temperature "
756 << " min:" << limits.min()
757 << " max:" << limits.max()
758 << " avg:" << gAverage(*this)
759 << endl;
760 }
761}
762
763
765(
766 Ostream& os
767) const
768{
770 os.writeEntryIfDifferent<word>("p", "p", pName_);
771 os.writeEntryIfDifferent<word>("U", "U", UName_);
772 os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
773 os.writeEntryIfDifferent<word>("mu", "thermo:mu", muName_);
774 os.writeEntryIfDifferent<word>("Tnbr", "T", TnbrName_);
775 os.writeEntryIfDifferent<word>("qrNbr", "none", qrNbrName_);
776 os.writeEntryIfDifferent<word>("qr", "none", qrName_);
777
778 if (fluid_)
779 {
780 os.writeEntry("mode", massModeTypeNames_[mode_]);
781
782 os.writeEntryIfDifferent<word>("specie", "none", specieName_);
783
784 os.writeEntry("carrierMolWeight", Mcomp_);
785
786 os.writeEntry("L", L_);
787 os.writeEntry("Tvap", Tvap_);
788 os.writeEntry("fluid", fluid_);
789 mass_.writeEntry("mass", os);
790
791 if (mode_ == mtConstantMass)
792 {
793 cp_.writeEntry("cp", os);
794 rho_.writeEntry("rho", os);
795 }
796
797 thickness_.writeEntry("thickness", os);
798 word liq = "liquid";
799 os << token::TAB << token::TAB << liq;
800 liquidDict_.write(os);
801 }
802
803 if (thicknessLayers_.size())
804 {
805 thicknessLayers_.writeEntry("thicknessLayers", os);
806 kappaLayers_.writeEntry("kappaLayers", os);
807 }
808
811
812
813// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
814
815namespace Foam
816{
818 (
821 );
822}
823
824
825// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
Macros for easy insertion into run-time selection tables.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
bool readIfPresent(const word &key, const dictionary &dict, EnumType &val, const bool warnOnly=false) const
Find an entry if present, and assign to T val.
Definition EnumI.H:111
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A min/max value pair with additional methods. In addition to conveniently storing values,...
Definition MinMax.H:106
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
scalar deltaTValue() const noexcept
Return time step value.
Definition TimeStateI.H:49
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition Time.C:714
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Mixed boundary condition for temperature to be used at the coupling interface between fluid solid reg...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
This boundary condition supplies a fixed gradient condition, such that the patch values are calculate...
virtual Field< Type > & gradient()
Return gradient at boundary.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const word & name() const
Return reference to name.
Definition fvMesh.H:387
A FieldMapper for finite-volume patch fields.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual const word & name() const
Return name.
Definition fvPatch.H:210
const scalarField & deltaCoeffs() const
Return the face - cell distance coefficient except for coupled patches for which the cell-centre to c...
Definition fvPatch.C:183
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup the named field from the local registry and return the patch field corresponding to this patch...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
humidityTemperatureCoupledMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
static autoPtr< liquidProperties > New(const word &name)
Return a pointer to a new liquidProperties created from name.
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.
virtual void write(Ostream &) const
Write.
Type * getObjectPtr(const word &name, const bool recursive=false) const
Return non-const pointer to the object of the given Type, using a const-cast to have it behave like a...
label index() const noexcept
The index of this patch in the boundaryMesh.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Common functions used in temperature coupled boundaries.
virtual tmp< scalarField > alpha(const scalarField &Tp) const
Given patch temperature calculate corresponding alphaEff field.
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
void write(Ostream &os) const
Write.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
temperatureCoupledBase(const fvPatch &patch, const KMethodType method=KMethodType::mtFluidThermo)
Default construct from patch, using fluidThermo (default) or specified method.
@ TAB
Tab [isspace].
Definition token.H:142
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
auto limits
Definition setRDeltaT.H:186
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > Uf
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Namespace for OpenFOAM.
Type gAverage(const FieldField< Field, Type > &f, const label comm)
The global arithmetic average of a FieldField.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
Type gSum(const FieldField< Field, Type > &f)
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar log(const dimensionedScalar &ds)
Type gWeightedSum(const UList< scalar > &weights, const UList< Type > &fld, const label comm)
The global weighted sum (integral) of a field, using the mag() of the weights.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
MinMax< Type > gMinMax(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
scalarField Re(const UList< complex > &cmplx)
Extract real component.
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition POSIX.C:1065
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fvPatchField< scalar > fvPatchScalarField
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha
dictionary dict
const volScalarField & cp
volScalarField & b
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.