thermoSingleLayer.H
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 -------------------------------------------------------------------------------
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 Class
27  Foam::regionModels::surfaceFilmModels::thermoSingleLayer
28 
29 Description
30  Thermodynamic form of single-cell layer surface film model
31 
32  Note: defining enthalpy as Cp(T - Tstd) - when using liquids from the
33  thermophysical library, their enthalpies are calculated similarly, where
34  Tstd = 298.15K. This is clearly non-conservative unless the heat-capacity
35  is constant and should be rewritten to use the standard thermodynamics
36  packages.
37 
38 SourceFiles
39  thermoSingleLayer.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef thermoSingleLayer_H
44 #define thermoSingleLayer_H
45 
46 #include "kinematicSingleLayer.H"
47 #include "SLGThermo.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 namespace regionModels
54 {
55 namespace surfaceFilmModels
56 {
57 
58 // Forward declaration of classes
59 class filmViscosityModel;
60 class heatTransferModel;
61 class phaseChangeModel;
62 class filmRadiationModel;
63 
64 /*---------------------------------------------------------------------------*\
65  Class thermoSingleLayer Declaration
66 \*---------------------------------------------------------------------------*/
67 
69 :
71 {
72  // Private member functions
73 
74  //- No copy construct
75  thermoSingleLayer(const thermoSingleLayer&) = delete;
76 
77  //- No copy assignment
78  void operator=(const thermoSingleLayer&) = delete;
79 
80  //- Return boundary types for sensible enthalpy field
81  wordList hsBoundaryTypes();
82 
83 
84 protected:
85 
86  // Protected data
87 
88  // Thermo properties
89 
90  //- Reference to the SLGThermo
91  const SLGThermo& thermo_;
92 
93 
94  // Fields
95 
96  //- Specific heat capacity [J/kg/K]
98 
99  //- Thermal conductivity [W/m/K]
101 
102  //- Temperature - mean [K]
104 
105  //- Temperature - surface [K]
107 
108  //- Temperature - wall [K]
110 
111  //- Sensible enthalpy [J/kg]
113 
114 
115  // Transfer fields - to the primary region
116 
117  //- Film energy transfer
119 
120 
121  //- Threshold film thickness beyond which the film is considered 'wet'
122  scalar deltaWet_;
123 
124 
125  // Hydrophilic/phobic properties
126 
127  //- Activation flag
128  bool hydrophilic_;
129 
130  //- Length scale applied to deltaWet_ to determine when a wet
131  // surface becomes dry, typically 0.5
132  scalar hydrophilicDryScale_;
133 
134  //- Length scale applied to deltaWet_ to determine when a dry
135  // surface becomes wet, typically 0.001
136  scalar hydrophilicWetScale_;
137 
138 
139  // Source term fields
140 
141  // Film region - registered to the film region mesh
142  // Note: need boundary value mapped from primary region, and then
143  // pushed into the patch internal field
144 
145  //- Energy [J/m2/s]
147 
148 
149  // Primary region - registered to the primary region mesh
150  // Internal use only - not read-in
151 
152  //- Energy [J/m2/s]
154 
155 
156  // Fields mapped from primary region - registered to the film region
157  // Note: need both boundary AND patch internal fields to be mapped
158 
159  //- Temperature [K]
161 
162  //- List of specie mass fractions [0-1]
164 
165 
166  // Sub-models
167 
168  //- Viscosity model
170 
171  //- Heat transfer coefficient between film surface and primary
172  // region [W/m2/K]
174 
175  //- Heat transfer coefficient between wall and film [W/m2/K]
177 
178  //- Phase change
180 
181  //- Radiation
183 
184 
185  // Limits
186 
187  //- Minimum temperature limit (optional)
188  scalar Tmin_;
189 
190  //- Maximum temperature limit (optional)
191  scalar Tmax_;
192 
193 
194  // Protected member functions
195 
196  //- Read control parameters from dictionary
197  virtual bool read();
198 
199  //- Correct the thermo fields
200  virtual void correctThermoFields();
201 
202  //- Correct sensible enthalpy for mapped temperature fields
203  virtual void correctHsForMappedT();
204 
205  //- Correct the film surface and wall temperatures
206  virtual void updateSurfaceTemperatures();
207 
208  //- Reset source term fields
209  virtual void resetPrimaryRegionSourceTerms();
210 
211  //- Transfer thermo fields from the primary region to the film region
212  virtual void transferPrimaryRegionThermoFields();
213 
214  //- Transfer source fields from the primary region to the film region
215  virtual void transferPrimaryRegionSourceFields();
216 
217  //- Correct film coverage field
218  virtual void correctAlpha();
219 
220  //- Update the film sub-models
221  virtual void updateSubmodels();
222 
223  //- Return the wall/surface heat transfer term for the enthalpy equation
224  virtual tmp<fvScalarMatrix> q(volScalarField& hs) const;
225 
226 
227  // Equations
228 
229  //- Solve energy equation
230  virtual void solveEnergy();
231 
232 
233 public:
234 
235  //- Runtime type information
236  TypeName("thermoSingleLayer");
237 
238 
239  // Constructors
240 
241  //- Construct from components
243  (
244  const word& modelType,
245  const fvMesh& mesh,
246  const dimensionedVector& g,
247  const word& regionType,
248  const bool readFields = true
249  );
250 
251 
252  //- Destructor
253  virtual ~thermoSingleLayer();
254 
255 
256  // Member Functions
257 
258  // Thermo properties
259 
260  //- Return const reference to the SLGThermo object
261  inline const SLGThermo& thermo() const;
262 
263 
264  // Fields
265 
266  //- Return the film specific heat capacity [J/kg/K]
267  virtual const volScalarField& Cp() const;
268 
269  //- Return the film thermal conductivity [W/m/K]
270  virtual const volScalarField& kappa() const;
271 
272  //- Return the film mean temperature [K]
273  virtual const volScalarField& T() const;
274 
275  //- Return the film surface temperature [K]
276  virtual const volScalarField& Ts() const;
277 
278  //- Return the film wall temperature [K]
279  virtual const volScalarField& Tw() const;
280 
281  //- Return the film sensible enthalpy [J/kg]
282  virtual const volScalarField& hs() const;
283 
284 
285  // Helper functions
286 
287  //- Return sensible enthalpy as a function of temperature
288  // for a patch
289  inline tmp<scalarField> hs
290  (
291  const scalarField& T,
292  const label patchi
293  ) const;
294 
295  //- Return sensible enthalpy as a function of temperature
296  inline tmp<volScalarField> hs
297  (
298  const volScalarField& T
299  ) const;
300 
301  //- Return temperature as a function of sensible enthalpy
302  inline tmp<volScalarField> T
303  (
304  const volScalarField& hs
305  ) const;
306 
307 
308  // Source fields (read/write access)
309 
310  //- External hook to add sources to the film
311  virtual void addSources
312  (
313  const label patchi, // patchi on primary region
314  const label facei, // facei of patchi
315  const scalar massSource, // [kg]
316  const vector& momentumSource, // [kg.m/s] (tangential momentum)
317  const scalar pressureSource, // [kg.m/s] (normal momentum)
318  const scalar energySource // [J]
319  );
320 
321 
322  // Source term fields
323 
324  // Film region
325 
326  //- Energy [J/m2/s]
327  inline const volScalarField& hsSp() const;
328 
329 
330  // Primary region
331 
332  //- Energy [J/m2/s]
333  inline const volScalarField& hsSpPrimary() const;
334 
335 
336  // Fields mapped from the primary region
337 
338  //- Temperature [K]
339  inline const volScalarField& TPrimary() const;
340 
341  //- Specie mass fractions [0-1]
342  inline const PtrList<volScalarField>& YPrimary() const;
343 
344 
345 
346  // Sub-models
347 
348  //- Return const access to the (surface) heat transfer model
349  inline const heatTransferModel& htcs() const;
350 
351  //- Return const access to the (wall) heat transfer model
352  inline const heatTransferModel& htcw() const;
353 
354  //- Return const access to the phase change model
355  inline const phaseChangeModel& phaseChange() const;
356 
357  //- Return const access to the radiation model
358  inline const filmRadiationModel& radiation() const;
359 
360 
361  // Derived fields (calculated on-the-fly)
362 
363  //- Return the convective heat energy from film to wall
364  inline tmp<scalarField> qconvw(const label patchi) const;
365 
366  //- Return the convective heat energy from primary region to film
367  inline tmp<scalarField> qconvp(const label patchi) const;
368 
369 
370  // Evolution
371 
372  //- Pre-evolve film hook
373  virtual void preEvolveRegion();
374 
375  //- Evolve the film equations
376  virtual void evolveRegion();
377 
378 
379  // Source fields
380 
381  // Mapped into primary region
382 
383  //- Return total mass source - Eulerian phase only
384  virtual tmp<volScalarField::Internal> Srho() const;
385 
386  //- Return mass source for specie i - Eulerian phase only
388  (
389  const label i
390  ) const;
391 
392  //- Return enthalpy source - Eulerian phase only
393  virtual tmp<volScalarField::Internal> Sh() const;
394 
395 
396  // I-O
397 
398  //- Provide some feedback
399  virtual void info();
400 };
401 
402 
403 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
404 
405 } // End namespace surfaceFilmModels
406 } // End namespace regionModels
407 } // End namespace Foam
408 
409 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
410 
411 #include "thermoSingleLayerI.H"
412 
413 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
414 
415 #endif
416 
417 // ************************************************************************* //
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::TypeName
TypeName("thermoSingleLayer")
Runtime type information.
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::thermoSingleLayer::read
virtual bool read()
Read control parameters from dictionary.
Definition: thermoSingleLayer.C:88
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::heatTransferModel
Base class for film heat transfer models.
Definition: heatTransferModel.H:56
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::YPrimary
const PtrList< volScalarField > & YPrimary() const
Specie mass fractions [0-1].
Definition: thermoSingleLayerI.H:130
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::~thermoSingleLayer
virtual ~thermoSingleLayer()
Destructor.
Definition: thermoSingleLayer.C:564
Foam::SLGThermo
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:64
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::regionModels::surfaceFilmModels::thermoSingleLayer::htcs
const heatTransferModel & htcs() const
Return const access to the (surface) heat transfer model.
Definition: thermoSingleLayerI.H:136
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::primaryEnergyTrans_
volScalarField primaryEnergyTrans_
Film energy transfer.
Definition: thermoSingleLayer.H:117
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::thermo
const SLGThermo & thermo() const
Return const reference to the SLGThermo object.
Definition: thermoSingleLayerI.H:44
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::phaseChange_
autoPtr< phaseChangeModel > phaseChange_
Phase change.
Definition: thermoSingleLayer.H:178
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::surfaceFilmModels::kinematicSingleLayer
Kinematic form of single-cell layer surface film model.
Definition: kinematicSingleLayer.H:67
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::kappa
virtual const volScalarField & kappa() const
Return the film thermal conductivity [W/m/K].
Definition: thermoSingleLayer.C:656
SLGThermo.H
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::T_
volScalarField T_
Temperature - mean [K].
Definition: thermoSingleLayer.H:102
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::thermoSingleLayer::qconvp
tmp< scalarField > qconvp(const label patchi) const
Return the convective heat energy from primary region to film.
Definition: thermoSingleLayerI.H:170
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::radiation
const filmRadiationModel & radiation() const
Return const access to the radiation model.
Definition: thermoSingleLayerI.H:154
Foam::regionModels::surfaceFilmModels::surfaceFilmRegionModel::g
const dimensionedVector & g() const
Return the acceleration due to gravity.
Definition: surfaceFilmRegionModelI.H:41
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::correctThermoFields
virtual void correctThermoFields()
Correct the thermo fields.
Definition: thermoSingleLayer.C:105
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::evolveRegion
virtual void evolveRegion()
Evolve the film equations.
Definition: thermoSingleLayer.C:606
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::deltaWet_
scalar deltaWet_
Threshold film thickness beyond which the film is considered 'wet'.
Definition: thermoSingleLayer.H:121
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hsSp_
volScalarField hsSp_
Energy [J/m2/s].
Definition: thermoSingleLayer.H:145
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::transferPrimaryRegionThermoFields
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
Definition: thermoSingleLayer.C:151
Foam::regionModels::surfaceFilmModels::thermoSingleLayer
Thermodynamic form of single-cell layer surface film model.
Definition: thermoSingleLayer.H:67
Foam::Field< scalar >
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::thermo_
const SLGThermo & thermo_
Reference to the SLGThermo.
Definition: thermoSingleLayer.H:90
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
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::thermoSingleLayer::Ts_
volScalarField Ts_
Temperature - surface [K].
Definition: thermoSingleLayer.H:105
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
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::YPrimary_
PtrList< volScalarField > YPrimary_
List of specie mass fractions [0-1].
Definition: thermoSingleLayer.H:162
kinematicSingleLayer.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::TPrimary_
volScalarField TPrimary_
Temperature [K].
Definition: thermoSingleLayer.H:159
Foam::dimensioned< vector >
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::Tmin_
scalar Tmin_
Minimum temperature limit (optional)
Definition: thermoSingleLayer.H:187
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
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve film hook.
Definition: thermoSingleLayer.C:597
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hsSp
const volScalarField & hsSp() const
Energy [J/m2/s].
Definition: thermoSingleLayerI.H:112
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::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::thermoSingleLayer::hsSpPrimary
const volScalarField & hsSpPrimary() const
Energy [J/m2/s].
Definition: thermoSingleLayerI.H:118
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::phaseChange
const phaseChangeModel & phaseChange() const
Return const access to the phase change model.
Definition: thermoSingleLayerI.H:148
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::Tmax_
scalar Tmax_
Maximum temperature limit (optional)
Definition: thermoSingleLayer.H:190
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::Srho
virtual tmp< volScalarField::Internal > Srho() const
Return total mass source - Eulerian phase only.
Definition: thermoSingleLayer.C:701
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::htcw
const heatTransferModel & htcw() const
Return const access to the (wall) heat transfer model.
Definition: thermoSingleLayerI.H:142
Foam::Vector< scalar >
Foam::List< word >
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::qconvw
tmp< scalarField > qconvw(const label patchi) const
Return the convective heat energy from film to wall.
Definition: thermoSingleLayerI.H:160
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::regionModels::surfaceFilmModels::thermoSingleLayer::correctAlpha
virtual void correctAlpha()
Correct film coverage field.
Definition: thermoSingleLayer.C:196
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::hydrophilicWetScale_
scalar hydrophilicWetScale_
Length scale applied to deltaWet_ to determine when a dry.
Definition: thermoSingleLayer.H:135
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::TPrimary
const volScalarField & TPrimary() const
Temperature [K].
Definition: thermoSingleLayerI.H:124
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
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::htcw_
autoPtr< heatTransferModel > htcw_
Heat transfer coefficient between wall and film [W/m2/K].
Definition: thermoSingleLayer.H:175
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::Ts
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
Definition: thermoSingleLayer.C:668
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::regionModels::surfaceFilmModels::filmRadiationModel
Base class for film radiation models.
Definition: filmRadiationModel.H:56
thermoSingleLayerI.H
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::updateSurfaceTemperatures
virtual void updateSurfaceTemperatures()
Correct the film surface and wall temperatures.
Definition: thermoSingleLayer.C:131
Foam::regionModels::surfaceFilmModels::thermoSingleLayer::transferPrimaryRegionSourceFields
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
Definition: thermoSingleLayer.C:167
Foam::regionModels::surfaceFilmModels::phaseChangeModel
Base class for surface film phase change models.
Definition: phaseChangeModel.H:57