thermoSingleLayer.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-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "thermoSingleLayer.H"
30 #include "fvcDdt.H"
31 #include "fvcDiv.H"
32 #include "fvcLaplacian.H"
33 #include "fvcFlux.H"
34 #include "fvm.H"
36 #include "mixedFvPatchFields.H"
38 #include "mapDistribute.H"
39 #include "constants.H"
41 
42 // Sub-models
43 #include "filmThermoModel.H"
44 #include "filmViscosityModel.H"
45 #include "heatTransferModel.H"
46 #include "phaseChangeModel.H"
47 #include "filmRadiationModel.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 namespace regionModels
54 {
55 namespace surfaceFilmModels
56 {
57 
58 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
59 
60 defineTypeNameAndDebug(thermoSingleLayer, 0);
61 
62 addToRunTimeSelectionTable(surfaceFilmRegionModel, thermoSingleLayer, mesh);
63 
64 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
65 
66 wordList thermoSingleLayer::hsBoundaryTypes()
67 {
68  wordList bTypes(T_.boundaryField().types());
69  forAll(bTypes, patchi)
70  {
71  if
72  (
73  T_.boundaryField()[patchi].fixesValue()
74  || isA<mixedFvPatchScalarField>(T_.boundaryField()[patchi])
75  || isA<mappedFieldFvPatchField<scalar>>(T_.boundaryField()[patchi])
76  )
77  {
78  bTypes[patchi] = fixedValueFvPatchField<scalar>::typeName;
79  }
80  }
81 
82  return bTypes;
83 }
84 
85 
86 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
87 
89 {
90  // No additional properties to read
92 }
93 
94 
96 {
98 
100 
101  hsSpPrimary_ == dimensionedScalar(hsSp_.dimensions(), Zero);
102 }
103 
104 
106 {
107  rho_ == filmThermo_->rho();
108  sigma_ == filmThermo_->sigma();
109  Cp_ == filmThermo_->Cp();
110  kappa_ == filmThermo_->kappa();
111 }
112 
113 
115 {
117 
118  volScalarField::Boundary& hsBf = hs_.boundaryFieldRef();
119 
120  forAll(hsBf, patchi)
121  {
122  const fvPatchField<scalar>& Tp = T_.boundaryField()[patchi];
124  {
125  hsBf[patchi] == hs(Tp, patchi);
126  }
127  }
128 }
129 
130 
132 {
134 
135  // Push boundary film temperature into wall temperature internal field
136  for (label i=0; i<intCoupledPatchIDs_.size(); i++)
137  {
138  label patchi = intCoupledPatchIDs_[i];
139  const polyPatch& pp = regionMesh().boundaryMesh()[patchi];
140  UIndirectList<scalar>(Tw_, pp.faceCells()) =
141  T_.boundaryField()[patchi];
142  }
144 
145  // Update film surface temperature
146  Ts_ = T_;
148 }
149 
150 
152 {
154 
156 
157  // Update primary region fields on local region via direct mapped (coupled)
158  // boundary conditions
160  forAll(YPrimary_, i)
161  {
162  YPrimary_[i].correctBoundaryConditions();
163  }
164 }
165 
166 
168 {
170 
172 
173  volScalarField::Boundary& hsSpPrimaryBf =
175 
176  // Convert accumulated source terms into per unit area per unit time
177  const scalar deltaT = time_.deltaTValue();
178  forAll(hsSpPrimaryBf, patchi)
179  {
180  scalarField rpriMagSfdeltaT
181  (
182  (1.0/deltaT)/primaryMesh().magSf().boundaryField()[patchi]
183  );
184 
185  hsSpPrimaryBf[patchi] *= rpriMagSfdeltaT;
186  }
187 
188  // Retrieve the source fields from the primary region via direct mapped
189  // (coupled) boundary conditions
190  // - fields require transfer of values for both patch AND to push the
191  // values into the first layer of internal cells
193 }
194 
195 
197 {
198  if (hydrophilic_)
199  {
200  const scalar hydrophilicDry = hydrophilicDryScale_*deltaWet_;
201  const scalar hydrophilicWet = hydrophilicWetScale_*deltaWet_;
202 
203  forAll(alpha_, i)
204  {
205  if ((alpha_[i] < 0.5) && (delta_[i] > hydrophilicWet))
206  {
207  alpha_[i] = 1.0;
208  }
209  else if ((alpha_[i] > 0.5) && (delta_[i] < hydrophilicDry))
210  {
211  alpha_[i] = 0.0;
212  }
213  }
214 
216  }
217  else
218  {
219  alpha_ ==
221  }
222 }
223 
224 
226 {
228 
229  // Update heat transfer coefficient sub-models
230  htcs_->correct();
231  htcw_->correct();
232 
233  // Update radiation
234  radiation_->correct();
235 
236  // Update injection model - mass returned is mass available for injection
238 
239  phaseChange_->correct
240  (
241  time_.deltaTValue(),
245  );
246 
247  const volScalarField rMagSfDt((1/time().deltaT())/magSf());
248 
249  // Vapour recoil pressure
250  pSp_ -= sqr(rMagSfDt*primaryMassTrans_)/(2*rhoPrimary_);
251 
252  // Update transfer model - mass returned is mass available for transfer
254 
255  // Update source fields
256  rhoSp_ += rMagSfDt*(cloudMassTrans_ + primaryMassTrans_);
258 
259  turbulence_->correct();
260 }
261 
262 
264 {
265  return
266  (
267  // Heat-transfer to the primary region
268  - fvm::Sp(htcs_->h()/Cp_, hs)
269  + htcs_->h()*(hs/Cp_ + alpha_*(TPrimary_ - T_))
270 
271  // Heat-transfer to the wall
272  - fvm::Sp(htcw_->h()/Cp_, hs)
273  + htcw_->h()*(hs/Cp_ + alpha_*(Tw_- T_))
274  );
275 }
276 
277 
279 {
281 
282  dimensionedScalar residualDeltaRho
283  (
284  "residualDeltaRho",
285  deltaRho_.dimensions(),
286  1e-10
287  );
288 
289  solve
290  (
292  + fvm::div(phi_, hs_)
293  ==
294  - hsSp_
295  + q(hs_)
296  + radiation_->Shs()
297  );
298 
300 
301  // Evaluate viscosity from user-model
302  viscosity_->correct(pPrimary_, T_);
303 
304  // Update film wall and surface temperatures
306 }
307 
308 
309 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
310 
311 thermoSingleLayer::thermoSingleLayer
312 (
313  const word& modelType,
314  const fvMesh& mesh,
315  const dimensionedVector& g,
316  const word& regionType,
317  const bool readFields
318 )
319 :
320  kinematicSingleLayer(modelType, mesh, g, regionType, false),
321  thermo_(mesh.lookupObject<SLGThermo>("SLGThermo")),
322  Cp_
323  (
324  IOobject
325  (
326  "Cp",
327  time().timeName(),
328  regionMesh(),
331  ),
332  regionMesh(),
334  zeroGradientFvPatchScalarField::typeName
335  ),
336  kappa_
337  (
338  IOobject
339  (
340  "kappa",
341  time().timeName(),
342  regionMesh(),
345  ),
346  regionMesh(),
348  zeroGradientFvPatchScalarField::typeName
349  ),
350 
351  T_
352  (
353  IOobject
354  (
355  "Tf",
356  time().timeName(),
357  regionMesh(),
360  ),
361  regionMesh()
362  ),
363  Ts_
364  (
365  IOobject
366  (
367  "Tsf",
368  time().timeName(),
369  regionMesh(),
372  ),
373  T_,
374  zeroGradientFvPatchScalarField::typeName
375  ),
376  Tw_
377  (
378  IOobject
379  (
380  "Twf",
381  time().timeName(),
382  regionMesh(),
385  ),
386  T_,
387  zeroGradientFvPatchScalarField::typeName
388  ),
389  hs_
390  (
391  IOobject
392  (
393  "hf",
394  time().timeName(),
395  regionMesh(),
398  ),
399  regionMesh(),
401  hsBoundaryTypes()
402  ),
403 
404  primaryEnergyTrans_
405  (
406  IOobject
407  (
408  "primaryEnergyTrans",
409  time().timeName(),
410  regionMesh(),
413  ),
414  regionMesh(),
416  zeroGradientFvPatchScalarField::typeName
417  ),
418 
419  deltaWet_(coeffs_.get<scalar>("deltaWet")),
420  hydrophilic_(coeffs_.get<bool>("hydrophilic")),
421  hydrophilicDryScale_(0.0),
422  hydrophilicWetScale_(0.0),
423 
424  hsSp_
425  (
426  IOobject
427  (
428  "hsSp",
429  time().timeName(),
430  regionMesh(),
433  ),
434  regionMesh(),
436  this->mappedPushedFieldPatchTypes<scalar>()
437  ),
438 
439  hsSpPrimary_
440  (
441  IOobject
442  (
443  hsSp_.name(), // Must have same name as hSp_ to enable mapping
444  time().timeName(),
445  primaryMesh(),
448  ),
449  primaryMesh(),
450  dimensionedScalar(hsSp_.dimensions(), Zero)
451  ),
452 
453  TPrimary_
454  (
455  IOobject
456  (
457  "T", // Same name as T on primary region to enable mapping
458  time().timeName(),
459  regionMesh(),
462  ),
463  regionMesh(),
465  this->mappedFieldAndInternalPatchTypes<scalar>()
466  ),
467 
468  YPrimary_(),
469 
470  viscosity_(filmViscosityModel::New(*this, coeffs(), mu_)),
471  htcs_
472  (
473  heatTransferModel::New(*this, coeffs().subDict("upperSurfaceModels"))
474  ),
475  htcw_
476  (
477  heatTransferModel::New(*this, coeffs().subDict("lowerSurfaceModels"))
478  ),
479  phaseChange_(phaseChangeModel::New(*this, coeffs())),
480  radiation_(filmRadiationModel::New(*this, coeffs())),
481  Tmin_(-VGREAT),
482  Tmax_(VGREAT)
483 {
484  if (coeffs().readIfPresent("Tmin", Tmin_))
485  {
486  Info<< " limiting minimum temperature to " << Tmin_ << endl;
487  }
488 
489  if (coeffs().readIfPresent("Tmax", Tmax_))
490  {
491  Info<< " limiting maximum temperature to " << Tmax_ << endl;
492  }
493 
494  if (thermo_.hasMultiComponentCarrier())
495  {
496  YPrimary_.setSize(thermo_.carrier().species().size());
497 
498  forAll(thermo_.carrier().species(), i)
499  {
500  YPrimary_.set
501  (
502  i,
503  new volScalarField
504  (
505  IOobject
506  (
507  thermo_.carrier().species()[i],
508  time().timeName(),
509  regionMesh(),
512  ),
513  regionMesh(),
515  pSp_.boundaryField().types()
516  )
517  );
518  }
519  }
520 
521  if (hydrophilic_)
522  {
523  coeffs_.readEntry("hydrophilicDryScale", hydrophilicDryScale_);
524  coeffs_.readEntry("hydrophilicWetScale", hydrophilicWetScale_);
525  }
526 
527  if (readFields)
528  {
529  transferPrimaryRegionThermoFields();
530 
531  correctAlpha();
532 
533  correctThermoFields();
534 
535  // Update derived fields
536  hs_ == hs(T_);
537 
538  deltaRho_ == delta_*rho_;
539 
541  (
542  IOobject
543  (
544  "phi",
545  time().timeName(),
546  regionMesh(),
549  false
550  ),
551  fvc::flux(deltaRho_*U_)
552  );
553 
554  phi_ == phi0;
555 
556  // Evaluate viscosity from user-model
557  viscosity_->correct(pPrimary_, T_);
558  }
559 }
560 
561 
562 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
563 
565 {}
566 
567 
568 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
569 
571 (
572  const label patchi,
573  const label facei,
574  const scalar massSource,
575  const vector& momentumSource,
576  const scalar pressureSource,
577  const scalar energySource
578 )
579 {
581  (
582  patchi,
583  facei,
584  massSource,
585  momentumSource,
586  pressureSource,
587  energySource
588  );
589 
590  DebugInfo
591  << " energy = " << energySource << nl << nl;
592 
593  hsSpPrimary_.boundaryFieldRef()[patchi][facei] -= energySource;
594 }
595 
596 
598 {
600 
603 }
604 
605 
607 {
609 
610  // Solve continuity for deltaRho_
611  solveContinuity();
612 
613  // Update sub-models to provide updated source contributions
614  updateSubmodels();
615 
616  // Solve continuity for deltaRho_
617  solveContinuity();
618 
619  for (int oCorr=1; oCorr<=nOuterCorr_; oCorr++)
620  {
621  // Explicit pressure source contribution
622  tmp<volScalarField> tpu(this->pu());
623 
624  // Implicit pressure source coefficient
625  tmp<volScalarField> tpp(this->pp());
626 
627  // Solve for momentum for U_
628  tmp<fvVectorMatrix> tUEqn = solveMomentum(tpu(), tpp());
629  fvVectorMatrix& UEqn = tUEqn.ref();
630 
631  // Solve energy for hs_ - also updates thermo
632  solveEnergy();
633 
634  // Film thickness correction loop
635  for (int corr=1; corr<=nCorr_; corr++)
636  {
637  // Solve thickness for delta_
638  solveThickness(tpu(), tpp(), UEqn);
639  }
640  }
641 
642  // Update deltaRho_ with new delta_
643  deltaRho_ == delta_*rho_;
644 
645  // Update temperature using latest hs_
646  T_ == T(hs_);
647 }
648 
649 
651 {
652  return Cp_;
653 }
654 
655 
657 {
658  return kappa_;
659 }
660 
661 
663 {
664  return T_;
665 }
666 
667 
669 {
670  return Ts_;
671 }
672 
673 
675 {
676  return Tw_;
677 }
678 
679 
681 {
682  return hs_;
683 }
684 
685 
687 {
689 
690  const scalarField& Tinternal = T_;
691 
692  Info<< indent << "min/mean/max(T) = "
693  << gMin(Tinternal) << ", "
694  << gAverage(Tinternal) << ", "
695  << gMax(Tinternal) << nl;
696 
697  phaseChange_->info(Info);
698 }
699 
700 
702 {
704  (
706  (
707  IOobject
708  (
709  typeName + ":Srho",
710  time().timeName(),
711  primaryMesh(),
714  false
715  ),
716  primaryMesh(),
718  )
719  );
720 
721  scalarField& Srho = tSrho.ref();
722  const scalarField& V = primaryMesh().V();
723  const scalar dt = time_.deltaTValue();
724 
726  {
727  const label filmPatchi = intCoupledPatchIDs()[i];
728 
729  scalarField patchMass =
730  primaryMassTrans_.boundaryField()[filmPatchi];
731 
732  toPrimary(filmPatchi, patchMass);
733 
734  const label primaryPatchi = primaryPatchIDs()[i];
735  const labelUList& cells =
736  primaryMesh().boundaryMesh()[primaryPatchi].faceCells();
737 
738  forAll(patchMass, j)
739  {
740  Srho[cells[j]] += patchMass[j]/(V[cells[j]]*dt);
741  }
742  }
743 
744  return tSrho;
745 }
746 
747 
749 (
750  const label i
751 ) const
752 {
753  const label vapId = thermo_.carrierId(filmThermo_->name());
754 
756  (
758  (
759  IOobject
760  (
761  typeName + ":Srho(" + Foam::name(i) + ")",
762  time_.timeName(),
763  primaryMesh(),
766  false
767  ),
768  primaryMesh(),
770  )
771  );
772 
773  if (vapId == i)
774  {
775  scalarField& Srho = tSrho.ref();
776  const scalarField& V = primaryMesh().V();
777  const scalar dt = time().deltaTValue();
778 
779  forAll(intCoupledPatchIDs_, i)
780  {
781  const label filmPatchi = intCoupledPatchIDs_[i];
782 
783  scalarField patchMass =
784  primaryMassTrans_.boundaryField()[filmPatchi];
785 
786  toPrimary(filmPatchi, patchMass);
787 
788  const label primaryPatchi = primaryPatchIDs()[i];
789  const labelUList& cells =
790  primaryMesh().boundaryMesh()[primaryPatchi].faceCells();
791 
792  forAll(patchMass, j)
793  {
794  Srho[cells[j]] += patchMass[j]/(V[cells[j]]*dt);
795  }
796  }
797  }
798 
799  return tSrho;
800 }
801 
802 
804 {
806  (
808  (
809  IOobject
810  (
811  typeName + ":Sh",
812  time().timeName(),
813  primaryMesh(),
816  false
817  ),
818  primaryMesh(),
820  )
821  );
822 
823  scalarField& Sh = tSh.ref();
824  const scalarField& V = primaryMesh().V();
825  const scalar dt = time_.deltaTValue();
826 
828  {
829  const label filmPatchi = intCoupledPatchIDs_[i];
830 
831  scalarField patchEnergy =
832  primaryEnergyTrans_.boundaryField()[filmPatchi];
833 
834  toPrimary(filmPatchi, patchEnergy);
835 
836  const label primaryPatchi = primaryPatchIDs()[i];
837  const labelUList& cells =
838  primaryMesh().boundaryMesh()[primaryPatchi].faceCells();
839 
840  forAll(patchEnergy, j)
841  {
842  Sh[cells[j]] += patchEnergy[j]/(V[cells[j]]*dt);
843  }
844  }
845 
846  return tSh;
847 }
848 
849 
850 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
851 
852 } // End namespace Foam
853 } // End namespace regionModels
854 } // End namespace surfaceFilmModels
855 
856 // ************************************************************************* //
Foam::fvPatchField< scalar >
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::solveEnergy
virtual void solveEnergy()
Solve energy equation.
Definition: thermoSingleLayer.C:278
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::Tw_
volScalarField Tw_
Temperature - wall [K].
Definition: thermoSingleLayer.H:108
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::read
virtual bool read()
Read control parameters from dictionary.
Definition: kinematicSingleLayer.C:59
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::read
virtual bool read()
Read control parameters from dictionary.
Definition: thermoSingleLayer.C:88
filmThermoModel.H
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::transferPrimaryRegionThermoFields
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
Definition: kinematicSingleLayer.C:94
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::resetPrimaryRegionSourceTerms
virtual void resetPrimaryRegionSourceTerms()
Reset source term fields.
Definition: thermoSingleLayer.C:95
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::info
virtual void info()
Provide some feedback.
Definition: thermoSingleLayer.C:686
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hs
virtual const volScalarField & hs() const
Return the film sensible enthalpy [J/kg].
Definition: thermoSingleLayer.C:680
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::correctHsForMappedT
virtual void correctHsForMappedT()
Correct sensible enthalpy for mapped temperature fields.
Definition: thermoSingleLayer.C:114
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::~thermoSingleLayer
virtual ~thermoSingleLayer()
Destructor.
Definition: thermoSingleLayer.C:564
Foam::regionModels::surfaceFilmModels::transferModelList::correct
virtual void correct(scalarField &availableMass, volScalarField &massToTransfer)
Correct kinematic transfers.
Definition: transferModelList.C:101
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::injection_
injectionModelList injection_
Cloud injection.
Definition: kinematicSingleLayer.H:214
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve film hook.
Definition: kinematicSingleLayer.C:853
Foam::SLGThermo
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:64
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::primaryMassTrans_
volScalarField primaryMassTrans_
Film mass available for transfer to the primary region.
Definition: kinematicSingleLayer.H:151
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::viscosity_
autoPtr< filmViscosityModel > viscosity_
Viscosity model.
Definition: thermoSingleLayer.H:168
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::htcs_
autoPtr< heatTransferModel > htcs_
Heat transfer coefficient between film surface and primary.
Definition: thermoSingleLayer.H:172
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::gAverage
Type gAverage(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:604
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::turbulence_
autoPtr< filmTurbulenceModel > turbulence_
Turbulence model.
Definition: kinematicSingleLayer.H:220
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::filmThermo_
autoPtr< filmThermoModel > filmThermo_
Film thermo model.
Definition: kinematicSingleLayer.H:208
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::solveContinuity
virtual void solveContinuity()
Solve continuity equation.
Definition: kinematicSingleLayer.C:261
Foam::regionModels::surfaceFilmModels::heatTransferModel::New
static autoPtr< heatTransferModel > New(surfaceFilmRegionModel &film, const dictionary &dict)
Return a reference to the selected phase change model.
Definition: heatTransferModelNew.C:43
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::primaryEnergyTrans_
volScalarField primaryEnergyTrans_
Film energy transfer.
Definition: thermoSingleLayer.H:117
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::phaseChange_
autoPtr< phaseChangeModel > phaseChange_
Phase change.
Definition: thermoSingleLayer.H:178
fvcDiv.H
Calculate the divergence of the given field.
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::Cp
virtual const volScalarField & Cp() const
Return the film specific heat capacity [J/kg/K].
Definition: thermoSingleLayer.C:650
Foam::regionModels::regionModel::time
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:40
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer
Kinematic form of single-cell layer surface film model.
Definition: kinematicSingleLayer.H:67
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rho_
volScalarField rho_
Density [kg/m3].
Definition: kinematicSingleLayer.H:115
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::sigma_
volScalarField sigma_
Surface tension [m/s2].
Definition: kinematicSingleLayer.H:121
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::kappa
virtual const volScalarField & kappa() const
Return the film thermal conductivity [W/m/K].
Definition: thermoSingleLayer.C:656
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::T_
volScalarField T_
Temperature - mean [K].
Definition: thermoSingleLayer.H:102
Foam::regionModels::surfaceFilmModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::Cp_
volScalarField Cp_
Specific heat capacity [J/kg/K].
Definition: thermoSingleLayer.H:96
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hydrophilicDryScale_
scalar hydrophilicDryScale_
Length scale applied to deltaWet_ to determine when a wet.
Definition: thermoSingleLayer.H:131
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::availableMass_
scalarField availableMass_
Available mass for transfer via sub-models.
Definition: kinematicSingleLayer.H:211
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pPrimary_
volScalarField pPrimary_
Pressure [Pa].
Definition: kinematicSingleLayer.H:196
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
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::correctThermoFields
virtual void correctThermoFields()
Correct the thermo fields.
Definition: thermoSingleLayer.C:105
thermoSingleLayer.H
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::evolveRegion
virtual void evolveRegion()
Evolve the film equations.
Definition: thermoSingleLayer.C:606
Foam::regionModels::regionModel::primaryMesh
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionModelI.H:34
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:62
Foam::regionModels::regionModel::toPrimary
void toPrimary(const label regionPatchi, List< Type > &regionField) const
Convert a local region field to the primary region.
Definition: regionModelTemplates.C:147
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::deltaWet_
scalar deltaWet_
Threshold film thickness beyond which the film is considered 'wet'.
Definition: thermoSingleLayer.H:121
mappedFieldFvPatchField.H
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::cloudDiameterTrans_
volScalarField cloudDiameterTrans_
Parcel diameters originating from film to cloud.
Definition: kinematicSingleLayer.H:157
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hsSp_
volScalarField hsSp_
Energy [J/m2/s].
Definition: thermoSingleLayer.H:145
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::deltaRho_
volScalarField deltaRho_
Film thickness*density (helper field) [kg/m2].
Definition: kinematicSingleLayer.H:142
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::transferPrimaryRegionThermoFields
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
Definition: thermoSingleLayer.C:151
tUEqn
tmp< fvVectorMatrix > tUEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::transferPrimaryRegionSourceFields
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
Definition: kinematicSingleLayer.C:107
Foam::Field< scalar >
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::delta_
volScalarField delta_
Film thickness [m].
Definition: kinematicSingleLayer.H:127
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::polyBoundaryMesh::faceCells
UPtrList< const labelUList > faceCells() const
Return a list of faceCells for each patch.
Definition: polyBoundaryMesh.C:311
Foam::regionModels::singleLayerRegion::magSf
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m2].
Definition: singleLayerRegion.C:223
Foam::regionModels::regionModel::regionMesh
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:64
Foam::IOobject::READ_IF_PRESENT
Definition: IOobject.H:187
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pp
virtual tmp< volScalarField > pp()
Implicit pressure source coefficient.
Definition: kinematicSingleLayer.C:184
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::info
virtual void info()
Provide some feedback.
Definition: kinematicSingleLayer.C:1060
Foam::constant::electromagnetic::phi0
const dimensionedScalar phi0
Magnetic flux quantum: default SI units: [Wb].
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::resetPrimaryRegionSourceTerms
virtual void resetPrimaryRegionSourceTerms()
Reset source term fields.
Definition: kinematicSingleLayer.C:84
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::addSources
virtual void addSources(const label patchi, const label facei, const scalar massSource, const vector &momentumSource, const scalar pressureSource, const scalar energySource)
External hook to add sources to the film.
Definition: thermoSingleLayer.C:571
Foam::regionModels::surfaceFilmModels::phaseChangeModel::New
static autoPtr< phaseChangeModel > New(surfaceFilmRegionModel &film, const dictionary &dict)
Return a reference to the selected phase change model.
Definition: phaseChangeModelNew.C:43
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::Ts_
volScalarField Ts_
Temperature - surface [K].
Definition: thermoSingleLayer.H:105
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::solveMomentum
virtual tmp< fvVectorMatrix > solveMomentum(const volScalarField &pu, const volScalarField &pp)
Solve for film velocity.
Definition: kinematicSingleLayer.C:293
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::updateSubmodels
virtual void updateSubmodels()
Update the film sub-models.
Definition: thermoSingleLayer.C:225
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::T
virtual const volScalarField & T() const
Return the film mean temperature [K].
Definition: thermoSingleLayer.C:662
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::YPrimary_
PtrList< volScalarField > YPrimary_
List of specie mass fractions [0-1].
Definition: thermoSingleLayer.H:162
fvm.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::TPrimary_
volScalarField TPrimary_
Temperature [K].
Definition: thermoSingleLayer.H:159
Foam::dimensioned< scalar >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:26
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::nCorr_
label nCorr_
Number of PISO-like correctors.
Definition: kinematicSingleLayer.H:95
Foam::indent
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:339
mixedFvPatchFields.H
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::q
virtual tmp< fvScalarMatrix > q(volScalarField &hs) const
Return the wall/surface heat transfer term for the enthalpy equation.
Definition: thermoSingleLayer.C:263
constants.H
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve film hook.
Definition: thermoSingleLayer.C:597
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pu
virtual tmp< volScalarField > pu()
Explicit pressure source contribution.
Definition: kinematicSingleLayer.C:162
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pSp_
volScalarField pSp_
Pressure [Pa].
Definition: kinematicSingleLayer.H:170
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::phi_
surfaceScalarField phi_
Mass flux (includes film thickness) [kg.m/s].
Definition: kinematicSingleLayer.H:145
Foam::regionModels::regionModel::primaryPatchIDs
const labelList & primaryPatchIDs() const
Return the list of patch IDs on the primary region coupled.
Definition: regionModelI.H:168
filmRadiationModel.H
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hsSpPrimary_
volScalarField hsSpPrimary_
Energy [J/m2/s].
Definition: thermoSingleLayer.H:152
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hs_
volScalarField hs_
Sensible enthalpy [J/kg].
Definition: thermoSingleLayer.H:111
Foam::regionModels::surfaceFilmModels::defineTypeNameAndDebug
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Foam::isA
const TargetType * isA(const Type &t)
Check if dynamic_cast to TargetType is possible.
Definition: typeInfo.H:197
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::addSources
virtual void addSources(const label patchi, const label facei, const scalar massSource, const vector &momentumSource, const scalar pressureSource, const scalar energySource=0)
External hook to add sources to the film.
Definition: kinematicSingleLayer.C:830
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::regionModels::regionModel::time_
const Time & time_
Reference to the time database.
Definition: regionModel.H:90
fvcLaplacian.H
Calculate the laplacian of the given field.
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:382
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::Srho
virtual tmp< volScalarField::Internal > Srho() const
Return total mass source - Eulerian phase only.
Definition: thermoSingleLayer.C:701
Foam::mappedFieldFvPatchField
This boundary condition provides a self-contained version of the mapped condition....
Definition: mappedFieldFvPatchField.H:114
mapDistribute.H
Foam::regionModels::surfaceFilmModels::filmRadiationModel::New
static autoPtr< filmRadiationModel > New(surfaceFilmRegionModel &film, const dictionary &dict)
Return a reference to the selected phase change model.
Definition: filmRadiationModelNew.C:43
Foam::regionModels::surfaceFilmModels::injectionModelList::correct
virtual void correct(scalarField &availableMass, volScalarField &massToInject, volScalarField &diameterToInject)
Correct.
Definition: injectionModelList.C:98
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::alpha_
volScalarField alpha_
Film coverage indicator, 1 = covered, 0 = uncovered [].
Definition: kinematicSingleLayer.H:130
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::Vector< scalar >
fvcFlux.H
Calculate the face-flux of the given field.
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoPrimary_
volScalarField rhoPrimary_
Density [kg/m3].
Definition: kinematicSingleLayer.H:199
Foam::UList< label >
filmViscosityModel.H
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::kappa_
volScalarField kappa_
Thermal conductivity [W/m/K].
Definition: thermoSingleLayer.H:99
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hydrophilic_
bool hydrophilic_
Activation flag.
Definition: thermoSingleLayer.H:127
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::radiation_
autoPtr< filmRadiationModel > radiation_
Radiation.
Definition: thermoSingleLayer.H:181
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::transfer_
transferModelList transfer_
Transfer with the continuous phase.
Definition: kinematicSingleLayer.H:217
Foam::regionModels::regionModel::intCoupledPatchIDs_
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
Definition: regionModel.H:114
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::correctAlpha
virtual void correctAlpha()
Correct film coverage field.
Definition: thermoSingleLayer.C:196
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hydrophilicWetScale_
scalar hydrophilicWetScale_
Length scale applied to deltaWet_ to determine when a dry.
Definition: thermoSingleLayer.H:135
fvcDdt.H
Calculate the first temporal derivative.
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::Sh
virtual tmp< volScalarField::Internal > Sh() const
Return enthalpy source - Eulerian phase only.
Definition: thermoSingleLayer.C:803
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::Tw
virtual const volScalarField & Tw() const
Return the film wall temperature [K].
Definition: thermoSingleLayer.C:674
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::htcw_
autoPtr< heatTransferModel > htcw_
Heat transfer coefficient between wall and film [W/m2/K].
Definition: thermoSingleLayer.H:175
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::Ts
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
Definition: thermoSingleLayer.C:668
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
zeroGradientFvPatchFields.H
Foam::regionModels::regionModel::intCoupledPatchIDs
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
Definition: regionModelI.H:175
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoSp_
volScalarField rhoSp_
Mass [kg/m2/s].
Definition: kinematicSingleLayer.H:173
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::solveThickness
virtual void solveThickness(const volScalarField &pu, const volScalarField &pp, fvVectorMatrix &UEqn)
Solve coupled velocity-thickness equations.
Definition: kinematicSingleLayer.C:348
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::cloudMassTrans_
volScalarField cloudMassTrans_
Film mass available for transfer to cloud.
Definition: kinematicSingleLayer.H:154
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::regionModels::surfaceFilmModels::filmViscosityModel::New
static autoPtr< filmViscosityModel > New(surfaceFilmRegionModel &film, const dictionary &dict, volScalarField &mu)
Return a reference to the selected phase change model.
Definition: filmViscosityModelNew.C:43
phaseChangeModel.H
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::nOuterCorr_
label nOuterCorr_
Number of outer correctors.
Definition: kinematicSingleLayer.H:92
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::updateSurfaceTemperatures
virtual void updateSurfaceTemperatures()
Correct the film surface and wall temperatures.
Definition: thermoSingleLayer.C:131
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::transferPrimaryRegionSourceFields
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
Definition: thermoSingleLayer.C:167
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::V
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Definition: fvMeshGeometry.C:179
Foam::IOobject::MUST_READ
Definition: IOobject.H:185