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