thermoSingleLayer.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "thermoSingleLayer.H"
27 #include "fvcDiv.H"
28 #include "fvcLaplacian.H"
29 #include "fvcFlux.H"
30 #include "fvm.H"
32 #include "mixedFvPatchFields.H"
34 #include "mapDistribute.H"
35 #include "constants.H"
37 
38 // Sub-models
39 #include "filmThermoModel.H"
40 #include "filmViscosityModel.H"
41 #include "heatTransferModel.H"
42 #include "phaseChangeModel.H"
43 #include "filmRadiationModel.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 namespace regionModels
50 {
51 namespace surfaceFilmModels
52 {
53 
54 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
55 
56 defineTypeNameAndDebug(thermoSingleLayer, 0);
57 
58 addToRunTimeSelectionTable(surfaceFilmRegionModel, thermoSingleLayer, mesh);
59 
60 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
61 
62 wordList thermoSingleLayer::hsBoundaryTypes()
63 {
64  wordList bTypes(T_.boundaryField().types());
65  forAll(bTypes, patchi)
66  {
67  if
68  (
69  T_.boundaryField()[patchi].fixesValue()
70  || isA<mixedFvPatchScalarField>(T_.boundaryField()[patchi])
71  || isA<mappedFieldFvPatchField<scalar>>(T_.boundaryField()[patchi])
72  )
73  {
75  }
76  }
77 
78  return bTypes;
79 }
80 
81 
82 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
83 
85 {
86  // No additional properties to read
88 }
89 
90 
92 {
93  if (debug)
94  {
96  }
97 
99 
101 }
102 
103 
105 {
106  rho_ == filmThermo_->rho();
107  sigma_ == filmThermo_->sigma();
108  Cp_ == filmThermo_->Cp();
109  kappa_ == filmThermo_->kappa();
110 }
111 
112 
114 {
116 
117  volScalarField::Boundary& hsBf = hs_.boundaryFieldRef();
118 
119  forAll(hsBf, patchi)
120  {
123  {
124  hsBf[patchi] == hs(Tp, patchi);
125  }
126  }
127 }
128 
129 
131 {
133 
134  // Push boundary film temperature into wall temperature internal field
135  for (label i=0; i<intCoupledPatchIDs_.size(); i++)
136  {
138  const polyPatch& pp = regionMesh().boundaryMesh()[patchi];
141  }
143 
144  // Update film surface temperature
145  Ts_ = T_;
147 }
148 
149 
151 {
152  if (debug)
153  {
154  InfoInFunction << endl;
155  }
156 
158 
159  // Update primary region fields on local region via direct mapped (coupled)
160  // boundary conditions
162  forAll(YPrimary_, i)
163  {
164  YPrimary_[i].correctBoundaryConditions();
165  }
166 }
167 
168 
170 {
171  if (debug)
172  {
173  InfoInFunction << endl;
174  }
175 
177 
178  volScalarField::Boundary& hsSpPrimaryBf =
180 
181  // Convert accummulated source terms into per unit area per unit time
182  const scalar deltaT = time_.deltaTValue();
183  forAll(hsSpPrimaryBf, patchi)
184  {
185  scalarField rpriMagSfdeltaT
186  (
187  (1.0/deltaT)/primaryMesh().magSf().boundaryField()[patchi]
188  );
189 
190  hsSpPrimaryBf[patchi] *= rpriMagSfdeltaT;
191  }
192 
193  // Retrieve the source fields from the primary region via direct mapped
194  // (coupled) boundary conditions
195  // - fields require transfer of values for both patch AND to push the
196  // values into the first layer of internal cells
198 }
199 
200 
202 {
203  if (hydrophilic_)
204  {
205  const scalar hydrophilicDry = hydrophilicDryScale_*deltaWet_;
206  const scalar hydrophilicWet = hydrophilicWetScale_*deltaWet_;
207 
208  forAll(alpha_, i)
209  {
210  if ((alpha_[i] < 0.5) && (delta_[i] > hydrophilicWet))
211  {
212  alpha_[i] = 1.0;
213  }
214  else if ((alpha_[i] > 0.5) && (delta_[i] < hydrophilicDry))
215  {
216  alpha_[i] = 0.0;
217  }
218  }
219 
221  }
222  else
223  {
224  alpha_ ==
226  }
227 }
228 
229 
231 {
232  if (debug)
233  {
234  InfoInFunction << endl;
235  }
236 
237  // Update heat transfer coefficient sub-models
238  htcs_->correct();
239  htcw_->correct();
240 
241  // Update radiation
242  radiation_->correct();
243 
244  // Update injection model - mass returned is mass available for injection
246 
247  phaseChange_->correct
248  (
249  time_.deltaTValue(),
253  );
254 
255  const volScalarField rMagSfDt((1/time().deltaT())/magSf());
256 
257  // Vapour recoil pressure
258  pSp_ -= sqr(rMagSfDt*primaryMassTrans_)/(2*rhoPrimary_);
259 
260  // Update transfer model - mass returned is mass available for transfer
262 
263  // Update source fields
264  rhoSp_ += rMagSfDt*(cloudMassTrans_ + primaryMassTrans_);
266 
267  turbulence_->correct();
268 }
269 
270 
272 {
274 
275  return
276  (
277  // Heat-transfer to the primary region
278  - fvm::Sp(htcs_->h()/Cp_, hs)
279  + htcs_->h()*(hs/Cp_ + alpha*(TPrimary_ - T_))
280 
281  // Heat-transfer to the wall
282  - fvm::Sp(htcw_->h()/Cp_, hs)
283  + htcw_->h()*(hs/Cp_ + alpha*(Tw_- T_))
284  );
285 }
286 
287 
289 {
290  if (debug)
291  {
292  InfoInFunction << endl;
293  }
294 
296 
297  solve
298  (
300  + fvm::div(phi_, hs_)
301  ==
302  - hsSp_
303  + q(hs_)
304  + radiation_->Shs()
305  );
306 
308 
309  // Evaluate viscosity from user-model
310  viscosity_->correct(pPrimary_, T_);
311 }
312 
313 
314 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
315 
317 (
318  const word& modelType,
319  const fvMesh& mesh,
320  const dimensionedVector& g,
321  const word& regionType,
322  const bool readFields
323 )
324 :
325  kinematicSingleLayer(modelType, mesh, g, regionType, false),
326  thermo_(mesh.lookupObject<SLGThermo>("SLGThermo")),
327  Cp_
328  (
329  IOobject
330  (
331  "Cp",
332  time().timeName(),
333  regionMesh(),
336  ),
337  regionMesh(),
339  zeroGradientFvPatchScalarField::typeName
340  ),
341  kappa_
342  (
343  IOobject
344  (
345  "kappa",
346  time().timeName(),
347  regionMesh(),
350  ),
351  regionMesh(),
353  zeroGradientFvPatchScalarField::typeName
354  ),
355 
356  T_
357  (
358  IOobject
359  (
360  "Tf",
361  time().timeName(),
362  regionMesh(),
365  ),
366  regionMesh()
367  ),
368  Ts_
369  (
370  IOobject
371  (
372  "Tsf",
373  time().timeName(),
374  regionMesh(),
377  ),
378  T_,
379  zeroGradientFvPatchScalarField::typeName
380  ),
381  Tw_
382  (
383  IOobject
384  (
385  "Twf",
386  time().timeName(),
387  regionMesh(),
390  ),
391  T_,
392  zeroGradientFvPatchScalarField::typeName
393  ),
394  hs_
395  (
396  IOobject
397  (
398  "hf",
399  time().timeName(),
400  regionMesh(),
403  ),
404  regionMesh(),
406  hsBoundaryTypes()
407  ),
408 
410  (
411  IOobject
412  (
413  "primaryEnergyTrans",
414  time().timeName(),
415  regionMesh(),
418  ),
419  regionMesh(),
421  zeroGradientFvPatchScalarField::typeName
422  ),
423 
424  deltaWet_(readScalar(coeffs_.lookup("deltaWet"))),
425  hydrophilic_(readBool(coeffs_.lookup("hydrophilic"))),
428 
429  hsSp_
430  (
431  IOobject
432  (
433  "hsSp",
434  time().timeName(),
435  regionMesh(),
438  ),
439  regionMesh(),
441  this->mappedPushedFieldPatchTypes<scalar>()
442  ),
443 
445  (
446  IOobject
447  (
448  hsSp_.name(), // Must have same name as hSp_ to enable mapping
449  time().timeName(),
450  primaryMesh(),
453  ),
454  primaryMesh(),
456  ),
457 
458  TPrimary_
459  (
460  IOobject
461  (
462  "T", // Same name as T on primary region to enable mapping
463  time().timeName(),
464  regionMesh(),
467  ),
468  regionMesh(),
470  this->mappedFieldAndInternalPatchTypes<scalar>()
471  ),
472 
473  YPrimary_(),
474 
476  htcs_
477  (
478  heatTransferModel::New(*this, coeffs().subDict("upperSurfaceModels"))
479  ),
480  htcw_
481  (
482  heatTransferModel::New(*this, coeffs().subDict("lowerSurfaceModels"))
483  ),
486  Tmin_(-vGreat),
487  Tmax_(vGreat)
488 {
489  if (coeffs().readIfPresent("Tmin", Tmin_))
490  {
491  Info<< " limiting minimum temperature to " << Tmin_ << endl;
492  }
493 
494  if (coeffs().readIfPresent("Tmax", Tmax_))
495  {
496  Info<< " limiting maximum temperature to " << Tmax_ << endl;
497  }
498 
500  {
501  YPrimary_.setSize(thermo_.carrier().species().size());
502 
503  forAll(thermo_.carrier().species(), i)
504  {
505  YPrimary_.set
506  (
507  i,
508  new volScalarField
509  (
510  IOobject
511  (
512  thermo_.carrier().species()[i],
513  time().timeName(),
514  regionMesh(),
517  ),
518  regionMesh(),
520  pSp_.boundaryField().types()
521  )
522  );
523  }
524  }
525 
526  if (hydrophilic_)
527  {
528  coeffs_.lookup("hydrophilicDryScale") >> hydrophilicDryScale_;
529  coeffs_.lookup("hydrophilicWetScale") >> hydrophilicWetScale_;
530  }
531 
532  if (readFields)
533  {
535 
536  correctAlpha();
537 
539 
540  // Update derived fields
541  hs_ == hs(T_);
542 
543  deltaRho_ == delta_*rho_;
544 
546  (
547  IOobject
548  (
549  "phi",
550  time().timeName(),
551  regionMesh(),
554  false
555  ),
557  );
558 
559  phi_ == phi0;
560 
561  // Evaluate viscosity from user-model
562  viscosity_->correct(pPrimary_, T_);
563  }
564 }
565 
566 
567 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
568 
570 {}
571 
572 
573 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
574 
576 (
577  const label patchi,
578  const label facei,
579  const scalar massSource,
580  const vector& momentumSource,
581  const scalar pressureSource,
582  const scalar energySource
583 )
584 {
586  (
587  patchi,
588  facei,
589  massSource,
590  momentumSource,
591  pressureSource,
592  energySource
593  );
594 
595  if (debug)
596  {
597  Info<< " energy = " << energySource << nl << endl;
598  }
599 
600  hsSpPrimary_.boundaryFieldRef()[patchi][facei] -= energySource;
601 }
602 
603 
605 {
606  if (debug)
607  {
608  InfoInFunction << endl;
609  }
610 
613 }
614 
615 
617 {
618  if (debug)
619  {
620  InfoInFunction << endl;
621  }
622 
623  // Update film coverage indicator
624  correctAlpha();
625 
626  // Update film wall and surface velocities
628 
629  // Update film wall and surface temperatures
631 
632  // Solve continuity for deltaRho_
633  solveContinuity();
634 
635  // Update sub-models to provide updated source contributions
636  updateSubmodels();
637 
638  // Solve continuity for deltaRho_
639  solveContinuity();
640 
641  for (int oCorr=1; oCorr<=nOuterCorr_; oCorr++)
642  {
643  // Explicit pressure source contribution
644  tmp<volScalarField> tpu(this->pu());
645 
646  // Implicit pressure source coefficient
647  tmp<volScalarField> tpp(this->pp());
648 
649  // Solve for momentum for U_
650  tmp<fvVectorMatrix> UEqn = solveMomentum(tpu(), tpp());
651 
652  // Solve energy for hs_ - also updates thermo
653  solveEnergy();
654 
655  // Film thickness correction loop
656  for (int corr=1; corr<=nCorr_; corr++)
657  {
658  // Solve thickness for delta_
659  solveThickness(tpu(), tpp(), UEqn());
660  }
661  }
662 
663  // Update deltaRho_ with new delta_
664  deltaRho_ == delta_*rho_;
665 
666  // Update temperature using latest hs_
667  T_ == T(hs_);
668 
669  // Reset source terms for next time integration
671 }
672 
673 
675 {
676  return Cp_;
677 }
678 
679 
681 {
682  return kappa_;
683 }
684 
685 
687 {
688  return T_;
689 }
690 
691 
693 {
694  return Ts_;
695 }
696 
697 
699 {
700  return Tw_;
701 }
702 
703 
705 {
706  return hs_;
707 }
708 
709 
711 {
713 
714  const scalarField& Tinternal = T_;
715 
716  Info<< indent << "min/mean/max(T) = "
717  << gMin(Tinternal) << ", "
718  << gAverage(Tinternal) << ", "
719  << gMax(Tinternal) << nl;
720 
721  phaseChange_->info(Info);
722 }
723 
724 
726 {
728  (
730  (
731  "thermoSingleLayer::Srho",
732  primaryMesh(),
734  )
735  );
736 
737  scalarField& Srho = tSrho.ref();
738  const scalarField& V = primaryMesh().V();
739  const scalar dt = time_.deltaTValue();
740 
742  {
743  const label filmPatchi = intCoupledPatchIDs()[i];
744 
745  scalarField patchMass =
746  primaryMassTrans_.boundaryField()[filmPatchi];
747 
748  toPrimary(filmPatchi, patchMass);
749 
750  const label primaryPatchi = primaryPatchIDs()[i];
751  const unallocLabelList& cells =
752  primaryMesh().boundaryMesh()[primaryPatchi].faceCells();
753 
754  forAll(patchMass, j)
755  {
756  Srho[cells[j]] += patchMass[j]/(V[cells[j]]*dt);
757  }
758  }
759 
760  return tSrho;
761 }
762 
763 
765 (
766  const label i
767 ) const
768 {
769  const label vapId = thermo_.carrierId(filmThermo_->name());
770 
772  (
774  (
775  "thermoSingleLayer::Srho(" + Foam::name(i) + ")",
776  primaryMesh(),
778  )
779  );
780 
781  if (vapId == i)
782  {
783  scalarField& Srho = tSrho.ref();
784  const scalarField& V = primaryMesh().V();
785  const scalar dt = time().deltaTValue();
786 
788  {
789  const label filmPatchi = intCoupledPatchIDs_[i];
790 
791  scalarField patchMass =
792  primaryMassTrans_.boundaryField()[filmPatchi];
793 
794  toPrimary(filmPatchi, patchMass);
795 
796  const label primaryPatchi = primaryPatchIDs()[i];
797  const unallocLabelList& cells =
798  primaryMesh().boundaryMesh()[primaryPatchi].faceCells();
799 
800  forAll(patchMass, j)
801  {
802  Srho[cells[j]] += patchMass[j]/(V[cells[j]]*dt);
803  }
804  }
805  }
806 
807  return tSrho;
808 }
809 
810 
812 {
814  (
816  (
817  "thermoSingleLayer::Sh",
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 unallocLabelList& cells =
837  primaryMesh().boundaryMesh()[primaryPatchIDs()[i]].faceCells();
838 
839  forAll(patchEnergy, j)
840  {
841  Sh[cells[j]] += patchEnergy[j]/(V[cells[j]]*dt);
842  }
843  }
844 
845  return tSh;
846 }
847 
848 
849 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
850 
851 } // end namespace Foam
852 } // end namespace regionModels
853 } // end namespace surfaceFilmModels
854 
855 // ************************************************************************* //
volScalarField Ts_
Temperature - surface [K].
virtual bool read()
Read control parameters from dictionary.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
autoPtr< filmThermoModel > filmThermo_
Film thermo model.
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
autoPtr< heatTransferModel > htcw_
Heat transfer coefficient between wall and film [W/m^2/K].
virtual void correctAlpha()
Correct film coverage field.
autoPtr< phaseChangeModel > phaseChange_
Phase change.
virtual const volScalarField & Cp() const
Return the film specific heat capacity [J/kg/K].
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
autoPtr< filmTurbulenceModel > turbulence_
Turbulence model.
thermoSingleLayer(const word &modelType, const fvMesh &mesh, const dimensionedVector &g, const word &regionType, const bool readFields=true)
Construct from components.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual void resetPrimaryRegionSourceTerms()
Reset source term fields.
const word & name() const
Return name.
Definition: IOobject.H:295
virtual void correct(scalarField &availableMass, volScalarField &massToInject, volScalarField &diameterToInject)
Correct.
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:226
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
volScalarField primaryMassTrans_
Film mass available for transfer to the primary region.
static autoPtr< filmRadiationModel > New(surfaceFilmRegionModel &film, const dictionary &dict)
Return a reference to the selected phase change model.
scalar Tmax_
Maximum temperature limit (optional)
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
Type gMin(const FieldField< Field, Type > &f)
static const char *const typeName
Definition: Field.H:105
const Boundary & boundaryField() const
Return const-reference to the boundary field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
virtual const volScalarField & kappa() const
Return the film thermal conductivity [W/m/K].
static autoPtr< phaseChangeModel > New(surfaceFilmRegionModel &film, const dictionary &dict)
Return a reference to the selected phase change model.
const speciesTable & species() const
Return the table of species.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
const dimensionedScalar phi0
Magnetic flux quantum: default SI units: [Wb].
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
bool hasMultiComponentCarrier() const
Thermo database has multi-component carrier flag.
Definition: SLGThermo.C:223
scalarField availableMass_
Available mass for transfer via sub-models.
volScalarField cloudDiameterTrans_
Parcel diameters originating from film to cloud.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
bool readBool(Istream &)
Definition: boolIO.C:60
autoPtr< heatTransferModel > htcs_
Heat transfer coefficient between film surface and primary.
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
static autoPtr< filmViscosityModel > New(surfaceFilmRegionModel &film, const dictionary &dict, volScalarField &mu)
Return a reference to the selected phase change model.
virtual void updateSurfaceTemperatures()
Correct the film surface and wall temperatures.
virtual tmp< volScalarField::Internal > Srho() const
Return total mass source - Eulerian phase only.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
virtual void solveThickness(const volScalarField &pu, const volScalarField &pp, const fvVectorMatrix &UEqn)
Solve coupled velocity-thickness equations.
Macros for easy insertion into run-time selection tables.
virtual void solveContinuity()
Solve continuity equation.
virtual tmp< volScalarField > pp()
Implicit pressure source coefficient.
volScalarField deltaRho_
Film thickness*density (helper field) [kg/m^2].
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:699
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
virtual tmp< fvScalarMatrix > q(volScalarField &hs) const
Return the wall/surface heat transfer term for the enthalpy equation.
virtual void correct(scalarField &availableMass, volScalarField &massToTransfer)
Correct kinematic transfers.
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const dimensionSet & dimensions() const
Return dimensions.
void toPrimary(const label regionPatchi, List< Type > &regionField) const
Convert a local region field to the primary region.
const fvMesh & primaryMesh() const
Return the reference to the primary mesh database.
Definition: regionModelI.H:31
virtual void updateSubmodels()
Update the film sub-models.
const cellShapeList & cells
virtual tmp< fvVectorMatrix > solveMomentum(const volScalarField &pu, const volScalarField &pp)
Solve for film velocity.
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:315
scalar hydrophilicDryScale_
Length scale applied to deltaWet_ to determine when a wet.
A class for handling words, derived from string.
Definition: word.H:59
Calculate the face-flux of the given field.
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
Calculate the laplacian of the given field.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
volScalarField kappa_
Thermal conductivity [W/m/K].
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:61
word timeName
Definition: getTimeIndex.H:3
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
Definition: SLGThermo.H:62
scalar hydrophilicWetScale_
Length scale applied to deltaWet_ to determine when a dry.
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
scalar Tmin_
Minimum temperature limit (optional)
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
virtual tmp< volScalarField > pu()
Explicit pressure source contribution.
const SLGThermo & thermo_
Reference to the SLGThermo.
Calculate the divergence of the given field.
const labelList & primaryPatchIDs() const
Return the list of patch IDs on the primary region coupled.
Definition: regionModelI.H:172
dimensionedScalar pos0(const dimensionedScalar &ds)
virtual void evolveRegion()
Evolve the film equations.
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.
static const char nl
Definition: Ostream.H:265
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Type gMax(const FieldField< Field, Type > &f)
virtual void resetPrimaryRegionSourceTerms()
Reset source term fields.
volScalarField primaryEnergyTrans_
Film energy transfer.
virtual void correctThermoFields()
Correct the thermo fields.
autoPtr< filmRadiationModel > radiation_
Radiation.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const dimensionSet dimEnergy
transferModelList transfer_
Transfer with the continuous phase.
kinematicSingleLayer(const word &modelType, const fvMesh &mesh, const dimensionedVector &g, const word &regionType, const bool readFields=true)
Construct from components.
dictionary coeffs_
Model coefficients dictionary.
Definition: regionModel.H:99
const basicSpecieMixture & carrier() const
Return reference to the gaseous components.
Definition: SLGThermo.C:108
virtual bool read()
Read control parameters from dictionary.
List< word > wordList
A List of words.
Definition: fileName.H:54
label carrierId(const word &cmptName, bool allowNotFound=false) const
Index of carrier component.
Definition: SLGThermo.C:148
virtual const volScalarField & magSf() const
Return the face area magnitudes / [m^2].
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
scalar deltaWet_
Threshold film thickness beyond which the film is considered &#39;wet&#39;.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
volScalarField Cp_
Specific heat capacity [J/kg/K].
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Type gAverage(const FieldField< Field, Type > &f)
const Time & time_
Reference to the time database.
Definition: regionModel.H:84
surfaceScalarField phi_
Mass flux (includes film thickness) [kg m/s].
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
volScalarField alpha_
Film coverage indicator, 1 = covered, 0 = uncovered [].
A List with indirect addressing.
Definition: fvMatrix.H:106
virtual void preEvolveRegion()
Pre-evolve film hook.
virtual void updateSurfaceVelocities()
Update film surface velocities.
PtrList< volScalarField > YPrimary_
List of specie mass fractions [0-1].
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
void correctBoundaryConditions()
Correct boundary field.
virtual const volScalarField & hs() const
Return the film sensible enthalpy [J/kg].
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
messageStream Info
autoPtr< filmViscosityModel > viscosity_
Viscosity model.
static autoPtr< heatTransferModel > New(surfaceFilmRegionModel &film, const dictionary &dict)
Return a reference to the selected phase change model.
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.
volScalarField cloudMassTrans_
Film mass available for transfer to cloud.
This boundary condition provides a self-contained version of the mapped condition. It does not use information on the patch; instead it holds thr data locally.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
labelList intCoupledPatchIDs_
List of patch IDs internally coupled with the primary region.
Definition: regionModel.H:111
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
const labelList & intCoupledPatchIDs() const
Return the list of patch IDs internally coupled with the.
Definition: regionModelI.H:179
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
const dictionary & coeffs() const
Return the model coefficients dictionary.
Definition: regionModelI.H:96
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
defineTypeNameAndDebug(kinematicSingleLayer, 0)
virtual void correctHsForMappedT()
Correct sensible enthalpy for mapped temperature fields.
volScalarField hs_
Sensible enthalpy [J/kg].
virtual tmp< volScalarField::Internal > Sh() const
Return enthalpy source - Eulerian phase only.
const volScalarField & alpha() const
Return the film coverage, 1 = covered, 0 = uncovered [].
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583
virtual const volScalarField & T() const
Return the film mean temperature [K].
#define InfoInFunction
Report an information message using Foam::Info.
virtual const volScalarField & Tw() const
Return the film wall temperature [K].