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
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 (
105 IOobject
106 (
107 fieldName,
108 mesh.time().timeName(),
109 mesh,
112 ),
113 mesh,
115 );
116
117 ptr->store();
118 }
119
120 return *ptr;
121}
122
123
124
125// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
126
127Foam::humidityTemperatureCoupledMixedFvPatchScalarField::
128humidityTemperatureCoupledMixedFvPatchScalarField
129(
130 const fvPatch& p,
132)
133:
134 mixedFvPatchScalarField(p, iF),
136 (
137 patch(),
138 "fluidThermo",
139 "undefined",
140 "undefined-K",
141 "undefined-alpha"
142 ),
143 mode_(mtConstantMass),
144 pName_("p"),
145 UName_("U"),
146 rhoName_("rho"),
147 muName_("thermo:mu"),
148 TnbrName_("T"),
149 qrNbrName_("none"),
150 qrName_("none"),
151 specieName_("none"),
152 liquid_(nullptr),
153 liquidDict_(nullptr),
154 mass_(patch().size(), Zero),
155 Tvap_(0.0),
156 myKDelta_(patch().size(), Zero),
157 dmHfg_(patch().size(), Zero),
158 mpCpTp_(patch().size(), Zero),
159 Mcomp_(0.0),
160 L_(0.0),
161 fluid_(false),
162 cp_(patch().size(), Zero),
163 thickness_(patch().size(), Zero),
164 rho_(patch().size(), Zero),
165 thicknessLayers_(0),
166 kappaLayers_(0)
167{
168 this->refValue() = 0.0;
169 this->refGrad() = 0.0;
170 this->valueFraction() = 1.0;
171}
172
173
174Foam::humidityTemperatureCoupledMixedFvPatchScalarField::
175humidityTemperatureCoupledMixedFvPatchScalarField
176(
178 const fvPatch& p,
180 const fvPatchFieldMapper& mapper
181)
182:
183 mixedFvPatchScalarField(psf, p, iF, mapper),
184 temperatureCoupledBase(patch(), psf),
185 mode_(psf.mode_),
186 pName_(psf.pName_),
187 UName_(psf.UName_),
188 rhoName_(psf.rhoName_),
189 muName_(psf.muName_),
190 TnbrName_(psf.TnbrName_),
191 qrNbrName_(psf.qrNbrName_),
192 qrName_(psf.qrName_),
193 specieName_(psf.specieName_),
194 liquid_(psf.liquid_.clone()),
195 liquidDict_(psf.liquidDict_),
196 mass_(psf.mass_, mapper),
197 Tvap_(psf.Tvap_),
198 myKDelta_(psf.myKDelta_, mapper),
199 dmHfg_(psf.dmHfg_, mapper),
200 mpCpTp_(psf.mpCpTp_, mapper),
201 Mcomp_(psf.Mcomp_),
202 L_(psf.L_),
203 fluid_(psf.fluid_),
204 cp_(psf.cp_, mapper),
205 thickness_(psf.thickness_, mapper),
206 rho_(psf.rho_, mapper),
207 thicknessLayers_(psf.thicknessLayers_),
208 kappaLayers_(psf.kappaLayers_)
209{}
210
211
212Foam::humidityTemperatureCoupledMixedFvPatchScalarField::
213humidityTemperatureCoupledMixedFvPatchScalarField
214(
215 const fvPatch& p,
217 const dictionary& dict
218)
219:
220 mixedFvPatchScalarField(p, iF),
222 mode_(mtCondensationAndEvaporation),
223 pName_(dict.getOrDefault<word>("p", "p")),
224 UName_(dict.getOrDefault<word>("U", "U")),
225 rhoName_(dict.getOrDefault<word>("rho", "rho")),
226 muName_(dict.getOrDefault<word>("mu", "thermo:mu")),
227 TnbrName_(dict.getOrDefault<word>("Tnbr", "T")),
228 qrNbrName_(dict.getOrDefault<word>("qrNbr", "none")),
229 qrName_(dict.getOrDefault<word>("qr", "none")),
230 specieName_(dict.getOrDefault<word>("specie", "none")),
231 liquid_(nullptr),
232 liquidDict_(),
233 mass_(patch().size(), Zero),
234 Tvap_(0.0),
235 myKDelta_(patch().size(), Zero),
236 dmHfg_(patch().size(), Zero),
237 mpCpTp_(patch().size(), Zero),
238 Mcomp_(0.0),
239 L_(0.0),
240 fluid_(false),
241 cp_(patch().size(), Zero),
242 thickness_(patch().size(), Zero),
243 rho_(patch().size(), Zero),
244 thicknessLayers_(0),
245 kappaLayers_(0)
246{
247 if (!isA<mappedPatchBase>(this->patch().patch()))
248 {
250 << "\n patch type '" << p.type()
251 << "' not type '" << mappedPatchBase::typeName << "'"
252 << "\n for patch " << p.name()
253 << " of field " << internalField().name()
254 << " in file " << internalField().objectPath()
255 << exit(FatalIOError);
256 }
257
259
260 if (massModeTypeNames_.readIfPresent("mode", dict, mode_))
261 {
262 fluid_ = true;
263 }
264
265 if (dict.readIfPresent("thicknessLayers", thicknessLayers_))
266 {
267 dict.readEntry("kappaLayers", kappaLayers_);
268 }
269
270 if (fluid_)
271 {
272 switch(mode_)
273 {
274 case mtConstantMass:
275 {
276 thickness_ = scalarField("thickness", dict, p.size());
277 cp_ = scalarField("cp", dict, p.size());
278 rho_ = scalarField("rho", dict, p.size());
279
280 break;
281 }
282 case mtCondensation:
283 case mtEvaporation:
285 {
286 dict.readEntry("carrierMolWeight", Mcomp_);
287 dict.readEntry("L", L_);
288 dict.readEntry("Tvap", Tvap_);
289 liquidDict_ = dict.subDict("liquid");
290 liquid_ =
291 liquidProperties::New(liquidDict_.subDict(specieName_));
292
293 if (dict.found("thickness"))
294 {
295 scalarField& Tp = *this;
296 const scalarField& magSf = patch().magSf();
297
298 // Assume initially standard pressure for rho calculation
299 scalar pf = 1e5;
300 thickness_ = scalarField("thickness", dict, p.size());
301 forAll(thickness_, i)
302 {
303 mass_[i] =
304 thickness_[i]*liquid_->rho(pf, Tp[i])*magSf[i];
305 }
306 }
307 fluid_ = true;
308
309 break;
310 }
311 default:
312 {
314 << "Did not find mode " << mode_
315 << " on patch " << patch().name()
316 << nl
317 << "Please set 'mode' to one of "
318 << massModeTypeNames_.sortedToc()
319 << exit(FatalIOError);
320 }
321 }
322 }
323
324
325
326 if (dict.found("refValue"))
327 {
328 // Full restart
329 refValue() = scalarField("refValue", dict, p.size());
330 refGrad() = scalarField("refGradient", dict, p.size());
331 valueFraction() = scalarField("valueFraction", dict, p.size());
332 }
333 else
334 {
335 // Start from user entered data. Assume fixedValue.
336 refValue() = *this;
337 refGrad() = 0.0;
338 valueFraction() = 1.0;
339 }
340}
341
342
343Foam::humidityTemperatureCoupledMixedFvPatchScalarField::
344humidityTemperatureCoupledMixedFvPatchScalarField
345(
348)
349:
350 mixedFvPatchScalarField(psf, iF),
351 temperatureCoupledBase(patch(), psf),
352 mode_(psf.mode_),
353 pName_(psf.pName_),
354 UName_(psf.UName_),
355 rhoName_(psf.rhoName_),
356 muName_(psf.muName_),
357 TnbrName_(psf.TnbrName_),
358 qrNbrName_(psf.qrNbrName_),
359 qrName_(psf.qrName_),
360 specieName_(psf.specieName_),
361 liquid_(psf.liquid_.clone()),
362 liquidDict_(psf.liquidDict_),
363 mass_(psf.mass_),
364 Tvap_(psf.Tvap_),
365 myKDelta_(psf.myKDelta_),
366 dmHfg_(psf.dmHfg_),
367 mpCpTp_(psf.mpCpTp_),
368 Mcomp_(psf.Mcomp_),
369 L_(psf.L_),
370 fluid_(psf.fluid_),
371 cp_(psf.cp_),
372 thickness_(psf.thickness_),
373 rho_(psf.rho_),
374 thicknessLayers_(psf.thicknessLayers_),
375 kappaLayers_(psf.kappaLayers_)
376{}
377
378
379// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
380
382(
383 const fvPatchFieldMapper& m
384)
385{
386 mixedFvPatchScalarField::autoMap(m);
388
389 if (fluid_)
390 {
391 mass_.autoMap(m);
392 myKDelta_.autoMap(m);
393 dmHfg_.autoMap(m);
394 mpCpTp_.autoMap(m);
395 cp_.autoMap(m);
396 thickness_.autoMap(m);
397 rho_.autoMap(m);
398 }
399}
400
401
403(
404 const fvPatchScalarField& ptf,
405 const labelList& addr
406)
407{
408 mixedFvPatchScalarField::rmap(ptf, addr);
409
411 refCast<const humidityTemperatureCoupledMixedFvPatchScalarField>
412 (
413 ptf
414 );
415
416 temperatureCoupledBase::rmap(tiptf, addr);
417 if (fluid_)
418 {
419 mass_.rmap(tiptf.mass_, addr);
420 myKDelta_.rmap(tiptf.myKDelta_, addr);
421 dmHfg_.rmap(tiptf.dmHfg_, addr);
422 mpCpTp_.rmap(tiptf.mpCpTp_, addr);
423 cp_.rmap(tiptf.cp_, addr);
424 thickness_.rmap(tiptf.thickness_, addr);
425 rho_.rmap(tiptf.rho_, addr);
426 }
427}
428
429
431{
432 if (updated())
433 {
434 return;
435 }
436
437 // Get the coupling information from the mappedPatchBase
438 const mappedPatchBase& mpp =
439 refCast<const mappedPatchBase>(patch().patch());
440
441 const scalarField& magSf = patch().magSf();
442
443 const label nbrPatchI = mpp.samplePolyPatch().index();
444 const polyMesh& mesh = patch().boundaryMesh().mesh();
445 const polyMesh& nbrMesh = mpp.sampleMesh();
446 const fvPatch& nbrPatch =
447 refCast<const fvMesh>(nbrMesh).boundary()[nbrPatchI];
448
450 nbrField =
451 refCast
452 <
454 >
455 (
456 nbrPatch.lookupPatchField<volScalarField, scalar>(TnbrName_)
457 );
458
459 // Swap to obtain full local values of neighbour internal field
460 scalarField nbrIntFld(nbrField.patchInternalField());
461 mpp.distribute(nbrIntFld);
462
463
464 scalarField& Tp = *this;
465
466 const volScalarField& T =
467 static_cast<const volScalarField&>(internalField());
468
469 const scalarField TpOld(T.oldTime().boundaryField()[patch().index()]);
470
471 scalarField Tin(patchInternalField());
472
473 const scalarField K(this->kappa(*this));
474
475 // Neighbour kappa done separately because we need kappa solid for the
476 // htc correlation
477 scalarField nbrK(nbrField.kappa(*this));
478 mpp.distribute(nbrK);
479
480 scalarField nrbDeltaCoeffs(nbrPatch.deltaCoeffs());
481 mpp.distribute(nrbDeltaCoeffs);
482
483 scalarField KDeltaNbr(nbrField.kappa(*this)*nbrPatch.deltaCoeffs());
484 mpp.distribute(KDeltaNbr);
485
486 myKDelta_ = K*patch().deltaCoeffs();
487
488 if (thicknessLayers_.size() > 0)
489 {
490 myKDelta_ = 1.0/myKDelta_;
491 forAll(thicknessLayers_, iLayer)
492 {
493 myKDelta_ += thicknessLayers_[iLayer]/kappaLayers_[iLayer];
494 }
495 myKDelta_ = 1.0/myKDelta_;
496 }
497
498 scalarField dm(patch().size(), Zero);
499
500 // Fluid Side
501 if (fluid_)
502 {
503 scalarField Yvp(patch().size(), Zero);
504 const scalar dt = mesh.time().deltaTValue();
505
506 const scalarField myDelta(patch().deltaCoeffs());
507
508 if (mode_ != mtConstantMass)
509 {
510 scalarField cp(patch().size(), Zero);
511 scalarField hfg(patch().size(), Zero);
512 scalarField htc(patch().size(), GREAT);
513 scalarField liquidRho(patch().size(), Zero);
514 scalarField Tdew(patch().size(), Zero);
515 scalarField RH(patch().size(), Zero);
516
519 (
520 refCast
521 <
523 >
524 (
525 patch().lookupPatchField<volScalarField, scalar>
526 (
527 specieName_
528 )
529 )
530 );
531
532 const fvPatchScalarField& pp =
533 patch().lookupPatchField<volScalarField, scalar>(pName_);
534
535 const fvPatchVectorField& Up =
536 patch().lookupPatchField<volVectorField, vector>(UName_);
537
538 const fvPatchScalarField& rhop =
539 patch().lookupPatchField<volScalarField, scalar>(rhoName_);
540
541 const fvPatchScalarField& mup =
542 patch().lookupPatchField<volScalarField, scalar>(muName_);
543
544 const vectorField Ui(Up.patchInternalField());
545 const scalarField Yi(Yp.patchInternalField());
546
547 forAll(Tp, faceI)
548 {
549 const scalar Tf = Tp[faceI];
550 const scalar Tint = Tin[faceI];
551 const vector Uf = Ui[faceI];
552 const scalar pf = pp[faceI];
553
554 const scalar muf = mup[faceI];
555 const scalar rhof = rhop[faceI];
556 const scalar nuf = muf/rhof;
557 const scalar pSat = liquid_->pv(pf, Tint);
558 const scalar Mv = liquid_->W();
559 const scalar TSat = liquid_->pvInvert(pSat);
560 const scalar Re = mag(Uf)*L_/nuf;
561
562 cp[faceI] = liquid_->Cp(pf, Tf);
563 hfg[faceI] = liquid_->hl(pf, Tf);
564
565 // Calculate relative humidity
566 const scalar invMwmean =
567 Yi[faceI]/Mv + (1.0 - Yi[faceI])/Mcomp_;
568 const scalar Xv = Yi[faceI]/invMwmean/Mv;
569 RH[faceI] = min(Xv*pf/pSat, 1.0);
570
571 scalar RHmin = 0.01;
572 Tdew[faceI] = 0.0;
573
574 if (RH[faceI] > RHmin)
575 {
576 scalar b = 243.5;
577 scalar c = 17.65;
578 scalar TintDeg = Tint - 273;
579 Tdew[faceI] =
580 b*(log(RH[faceI]) + (c*TintDeg)/(b + TintDeg))
581 /(c - log(RH[faceI]) - ((c*TintDeg)/(b + TintDeg)))
582 + 273;
583 }
584
585 if
586 (
587 Tf < Tdew[faceI]
588 && RH[faceI] > RHmin
589 && (
590 mode_ == mtCondensation
591 || mode_ == mtCondensationAndEvaporation
592 )
593 )
594 {
595 htc[faceI] = htcCondensation(TSat, Re)*nbrK[faceI]/L_;
596
597 scalar htcTotal =
598 1.0/((1.0/myKDelta_[faceI]) + (1.0/htc[faceI]));
599
600 // Heat flux [W] (>0 heat is converted into mass)
601 const scalar q = (Tint - Tf)*htcTotal*magSf[faceI];
602
603 // Mass flux rate [Kg/s/m2]
604 dm[faceI] = q/hfg[faceI]/magSf[faceI];
605
606 mass_[faceI] += q/hfg[faceI]*dt;
607
608 // -dYp/dn = q/Dab (fixedGradient)
609 const scalar Dab = liquid_->D(pf, Tf);
610 Yvp[faceI] =
611 -min(dm[faceI]/Dab/rhof, Yi[faceI]*myDelta[faceI]);
612 }
613 else if
614 (
615 Tf > Tvap_
616 && mass_[faceI] > 0.0
617 && (
618 mode_ == mtEvaporation
619 || mode_ == mtCondensationAndEvaporation
620 )
621 )
622 {
623 const scalar Dab = liquid_->D(pf, Tf);
624
625 const scalar Sc = nuf/Dab;
626 const scalar Sh = this->Sh(Re, Sc);
627
628 const scalar Ys = Mv*pSat/(Mv*pSat + Mcomp_*(pf - pSat));
629
630 // Mass transfer coefficient [m/s]
631 const scalar hm = Dab*Sh/L_;
632
633 const scalar Yinf = max(Yi[faceI], 0.0);
634
635 // Mass flux rate [Kg/s/m2]
636 dm[faceI] = -rhof*hm*max((Ys - Yinf), 0.0)/(1.0 - Ys);
637
638 // Set fixedGradient for carrier species.
639 Yvp[faceI] = -dm[faceI]/Dab/rhof;
640
641 // Total mass accumulated [Kg]
642 mass_[faceI] += dm[faceI]*magSf[faceI]*dt;
643
644 htc[faceI] = htcCondensation(TSat, Re)*nbrK[faceI]/L_;
645 }
646 else if (Tf > Tdew[faceI] && Tf < Tvap_ && mass_[faceI] > 0.0)
647 {
648 htc[faceI] = htcCondensation(TSat, Re)*nbrK[faceI]/L_;
649 }
650 else if (mass_[faceI] == 0.0)
651 {
652 // Do nothing
653 }
654
655 liquidRho[faceI] = liquid_->rho(pf, Tf);
656 }
657
658 mass_ = max(mass_, scalar(0));
659
660 Yp.gradient() = Yvp;
661
662 // Output film delta (e.g. H2OThickness) [m]
663 const word fieldName(specieName_ + "Thickness");
664
665 scalarField& pDelta =
666 thicknessField
667 (
668 fieldName,
669 refCast<const fvMesh>(mesh)
670 ).boundaryFieldRef()[patch().index()];
671
672
673 pDelta = mass_/liquidRho/magSf;
674
675 // Weight myKDelta and htc
676 myKDelta_ = 1.0/((1.0/myKDelta_) + (1.0/htc));
677
678 mpCpTp_ = mass_*cp/dt/magSf;
679
680 // Heat flux due to change of phase [W/m2]
681 dmHfg_ = dm*hfg;
682
683 if (debug)
684 {
685 // Output RH and Tdew
686 scalarField& bRH =
687 thicknessField
688 (
689 "RH",
690 refCast<const fvMesh>(mesh)
691 ).boundaryFieldRef()[patch().index()];
692 bRH = RH;
693
694 scalarField& bTdew =
695 thicknessField
696 (
697 "Tdew",
698 refCast<const fvMesh>(mesh)
699 ).boundaryFieldRef()[patch().index()];
700 bTdew = Tdew;
701 }
702 }
703 else
704 {
705 // Inertia term [W/K/m2]
706 mpCpTp_ = thickness_*rho_*cp_/dt;
707 }
708 }
709
710 scalarField mpCpTpNbr(patch().size(), Zero);
711 scalarField dmHfgNbr(patch().size(), Zero);
712
713 if (!fluid_)
714 {
715 mpCpTpNbr = nbrField.mpCpTp();
716 mpp.distribute(mpCpTpNbr);
717
718 dmHfgNbr = nbrField.dmHfg();
719 mpp.distribute(dmHfgNbr);
720 }
721
722 // Obtain Rad heat (qr)
723 scalarField qr(Tp.size(), Zero);
724 if (qrName_ != "none")
725 {
726 qr = patch().lookupPatchField<volScalarField, scalar>(qrName_);
727 }
728
729 scalarField qrNbr(Tp.size(), Zero);
730 if (qrNbrName_ != "none")
731 {
732 qrNbr = nbrPatch.lookupPatchField<volScalarField, scalar>(qrNbrName_);
733 mpp.distribute(qrNbr);
734 }
735
736 const scalarField dmHfg(dmHfgNbr + dmHfg_);
737
738 const scalarField mpCpdt(mpCpTpNbr + mpCpTp_);
739
740 // qr > 0 (heat up the wall)
741 scalarField alpha(KDeltaNbr + mpCpdt - (qr + qrNbr)/Tp);
742
743 valueFraction() = alpha/(alpha + myKDelta_);
744
745 refValue() = (KDeltaNbr*nbrIntFld + mpCpdt*TpOld + dmHfg)/alpha;
746
747 mixedFvPatchScalarField::updateCoeffs();
748
749 if (debug && fluid_)
750 {
751 scalar Qdm = gSum(dm*magSf);
752 scalar QMass = gSum(mass_);
753 scalar Qt = gSum(myKDelta_*(Tp - Tin)*magSf);
754 scalar QtSolid = gSum(KDeltaNbr*(Tp - nbrIntFld)*magSf);
755
756 Info<< mesh.name() << ':'
757 << patch().name() << ':'
758 << internalField().name() << " <- "
759 << nbrMesh.name() << ':'
760 << nbrPatch.name() << ':'
761 << internalField().name() << " :" << nl
762 << " Total mass flux [Kg/s] : " << Qdm << nl
763 << " Total mass on the wall [Kg] : " << QMass << nl
764 << " Total heat (>0 leaving the wall to the fluid) [W] : "
765 << Qt << nl
766 << " Total heat (>0 leaving the wall to the solid) [W] : "
767 << QtSolid << nl
768 << " wall temperature "
769 << " min:" << gMin(*this)
770 << " max:" << gMax(*this)
771 << " avg:" << gAverage(*this)
772 << endl;
773 }
774}
775
776
778(
779 Ostream& os
780) const
781{
782 mixedFvPatchScalarField::write(os);
783 os.writeEntryIfDifferent<word>("p", "p", pName_);
784 os.writeEntryIfDifferent<word>("U", "U", UName_);
785 os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
786 os.writeEntryIfDifferent<word>("mu", "thermo:mu", muName_);
787 os.writeEntryIfDifferent<word>("Tnbr", "T", TnbrName_);
788 os.writeEntryIfDifferent<word>("qrNbr", "none", qrNbrName_);
789 os.writeEntryIfDifferent<word>("qr", "none", qrName_);
790
791 if (fluid_)
792 {
793 os.writeEntry("mode", massModeTypeNames_[mode_]);
794
795 os.writeEntryIfDifferent<word>("specie", "none", specieName_);
796
797 os.writeEntry("carrierMolWeight", Mcomp_);
798
799 os.writeEntry("L", L_);
800 os.writeEntry("Tvap", Tvap_);
801 os.writeEntry("fluid", fluid_);
802 mass_.writeEntry("mass", os);
803
804 if (mode_ == mtConstantMass)
805 {
806 cp_.writeEntry("cp", os);
807 rho_.writeEntry("rho", os);
808 }
809
810 thickness_.writeEntry("thickness", os);
811 word liq = "liquid";
812 os << token::TAB << token::TAB << liq;
813 liquidDict_.write(os);
814 }
815
816 if (thicknessLayers_.size())
817 {
818 thicknessLayers_.writeEntry("thicknessLayers", os);
819 kappaLayers_.writeEntry("kappaLayers", os);
820 }
821
823}
824
825
826// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
827
828namespace Foam
829{
831 (
834 );
835}
836
837
838// ************************************************************************* //
CGAL::Exact_predicates_exact_constructions_kernel K
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...
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
bool readIfPresent(const word &key, const dictionary &dict, EnumType &val) const
Find an entry if present, and assign to T val.
Definition: EnumI.H:132
List< word > sortedToc() const
The sorted list of enum names.
Definition: EnumI.H:70
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
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
scalar Sh() const
Sherwood number.
The TAB Method for Numerical Calculation of Spray Droplet Breakup.
Definition: TAB.H:69
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
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
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
This boundary condition supplies a fixed gradient condition, such that the patch values are calculate...
virtual Field< Type > & gradient()
Return gradient at boundary.
virtual bool write()
Write the output fields.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
const word & name() const
Return reference to name.
Definition: fvMesh.H:310
A FieldMapper for finite-volume patch fields.
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:237
virtual void operator=(const UList< scalar > &)
Definition: fvPatchField.C:408
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
virtual const word & name() const
Return name.
Definition: fvPatch.H:173
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
const scalarField & deltaCoeffs() const
Definition: fvPatch.C:196
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
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.
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.
Type * getObjectPtr(const word &name, const bool recursive=false) const
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.
virtual tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
const volScalarField & T
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > Uf
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Type gSum(const FieldField< Field, Type > &f)
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
const dimensionSet dimless
Dimensionless.
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
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.
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
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:33
IOerror FatalIOError
Type gAverage(const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dimensionedScalar cbrt(const dimensionedScalar &ds)
Type gMin(const FieldField< Field, Type > &f)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
volScalarField & alpha
dictionary dict
const volScalarField & cp
volScalarField & b
Definition: createFields.H:27
#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.
Foam::surfaceFields.