kinematicSingleLayer.H
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-2019 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 Class
25  Foam::regionModels::surfaceFilmModels::kinematicSingleLayer
26 
27 Description
28  Kinematic form of single-cell layer surface film model
29 
30 SourceFiles
31  kinematicSingleLayer.C
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #ifndef kinematicSingleLayer_H
36 #define kinematicSingleLayer_H
37 
38 #include "surfaceFilmRegionModel.H"
39 #include "fvMesh.H"
40 #include "volFields.H"
41 #include "surfaceFields.H"
42 #include "fvMatrices.H"
43 
44 #include "injectionModelList.H"
45 #include "transferModelList.H"
46 #include "forceList.H"
47 #include "filmTurbulenceModel.H"
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 namespace Foam
52 {
53 namespace regionModels
54 {
55 namespace surfaceFilmModels
56 {
57 
58 // Forward class declarations
59 class filmThermoModel;
60 
61 /*---------------------------------------------------------------------------*\
62  Class kinematicSingleLayer Declaration
63 \*---------------------------------------------------------------------------*/
64 
66 :
68 {
69 protected:
70 
71  // Protected data
72 
73  // Solution parameters
74 
75  //- Momentum predictor
77 
78  //- Number of outer correctors
80 
81  //- Number of PISO-like correctors
82  label nCorr_;
83 
84  //- Number of non-orthogonal correctors
86 
87  //- Cumulative continuity error
88  scalar cumulativeContErr_;
89 
90  //- Small delta
92 
93  //- Film thickness above which Courant number calculation in valid
94  scalar deltaCoLimit_;
95 
96 
97  // Thermo properties
98 
99  // Fields
100 
101  //- Density [kg/m^3]
103 
104  //- Dynamic viscosity [Pa.s]
106 
107  //- Surface tension [m/s^2]
109 
110 
111  // Fields
112 
113  //- Film thickness [m]
115 
116  //- Film coverage indicator, 1 = covered, 0 = uncovered []
118 
119  //- Velocity - mean [m/s]
121 
122  //- Velocity - surface [m/s]
124 
125  //- Velocity - wall [m/s]
127 
128  //- Film thickness*density (helper field) [kg/m^2]
130 
131  //- Mass flux (includes film thickness) [kg m/s]
133 
134 
135  // Transfer fields
136 
137  //- Film mass available for transfer to the primary region
139 
140  //- Film mass available for transfer to cloud
142 
143  //- Parcel diameters originating from film to cloud
145 
146 
147  // Source term fields
148 
149  // Film region - registered to the film region mesh
150  // Note: need boundary value mapped from primary region, and then
151  // pushed into the patch internal field
152 
153  //- Momementum [kg/m/s^2]
155 
156  //- Pressure [Pa]
158 
159  //- Mass [kg/m^2/s]
161 
162 
163  // Primary region - registered to the primary region mesh
164  // Internal use only - not read-in
165 
166  //- Momementum [kg/m/s^2]
168 
169  //- Pressure [Pa]
171 
172  //- Mass [kg/m^2/s]
174 
175 
176  // Fields mapped from primary region - registered to the film region
177  // Note: need both boundary AND patch internal fields to be mapped
178 
179  //- Velocity [m/s]
181 
182  //- Pressure [Pa]
184 
185  //- Density [kg/m^3]
187 
188  //- Viscosity [Pa.s]
190 
191 
192  // Sub-models
193 
194  //- Film thermo model
196 
197  //- Available mass for transfer via sub-models
199 
200  //- Cloud injection
202 
203  //- Transfer with the continuous phase
205 
206  //- Turbulence model
208 
209  //- List of film forces
211 
212 
213  // Checks
214 
215  //- Cumulative mass added via sources [kg]
216  scalar addedMassTotal_;
217 
218 
219  // Protected member functions
220 
221  //- Read control parameters from dictionary
222  virtual bool read();
223 
224  //- Correct the thermo fields
225  virtual void correctThermoFields();
226 
227  //- Reset source term fields
228  virtual void resetPrimaryRegionSourceTerms();
229 
230  //- Transfer thermo fields from the primary region to the film region
231  virtual void transferPrimaryRegionThermoFields();
232 
233  //- Transfer source fields from the primary region to the film region
234  virtual void transferPrimaryRegionSourceFields();
235 
236  //- Explicit pressure source contribution
237  virtual tmp<volScalarField> pu();
238 
239  //- Implicit pressure source coefficient
240  virtual tmp<volScalarField> pp();
241 
242  //- Correct film coverage field
243  virtual void correctAlpha();
244 
245  //- Update the film sub-models
246  virtual void updateSubmodels();
247 
248  //- Continuity check
249  virtual void continuityCheck();
250 
251  //- Update film surface velocities
252  virtual void updateSurfaceVelocities();
253 
254  //- Constrain a film region master/slave boundaries of a field to a
255  // given value
256  template<class Type>
257  void constrainFilmField
258  (
259  Type& field,
260  const typename Type::cmptType& value
261  );
262 
263 
264  // Equations
265 
266  //- Solve continuity equation
267  virtual void solveContinuity();
268 
269  //- Solve for film velocity
271  (
272  const volScalarField& pu,
273  const volScalarField& pp
274  );
275 
276  //- Solve coupled velocity-thickness equations
277  virtual void solveThickness
278  (
279  const volScalarField& pu,
280  const volScalarField& pp,
281  const fvVectorMatrix& UEqn
282  );
283 
284 
285 public:
286 
287  //- Runtime type information
288  TypeName("kinematicSingleLayer");
289 
290 
291  // Constructors
292 
293  //- Construct from components
295  (
296  const word& modelType,
297  const fvMesh& mesh,
298  const dimensionedVector& g,
299  const word& regionType,
300  const bool readFields = true
301  );
302 
303  //- Disallow default bitwise copy construction
305 
306 
307  //- Destructor
308  virtual ~kinematicSingleLayer();
309 
310 
311  // Member Functions
312 
313  // Solution parameters
314 
315  //- Courant number evaluation
316  virtual scalar CourantNumber() const;
317 
318  //- Return the momentum predictor
319  inline const Switch& momentumPredictor() const;
320 
321  //- Return the number of outer correctors
322  inline label nOuterCorr() const;
323 
324  //- Return the number of PISO correctors
325  inline label nCorr() const;
326 
327  //- Return the number of non-orthogonal correctors
328  inline label nNonOrthCorr() const;
329 
330  //- Return small delta
331  inline const dimensionedScalar& deltaSmall() const;
332 
333 
334  // Thermo properties
335 
336  //- Return const access to the dynamic viscosity [Pa.s]
337  inline const volScalarField& mu() const;
338 
339  //- Return const access to the surface tension [kg/s^2]
340  inline const volScalarField& sigma() const;
341 
342 
343  // Fields
344 
345  //- Return const access to the film thickness [m]
346  inline const volScalarField& delta() const;
347 
348  //- Return the film coverage, 1 = covered, 0 = uncovered []
349  inline const volScalarField& alpha() const;
350 
351  //- Return the film velocity [m/s]
352  virtual const volVectorField& U() const;
353 
354  //- Return the film surface velocity [m/s]
355  virtual const volVectorField& Us() const;
356 
357  //- Return the film wall velocity [m/s]
358  virtual const volVectorField& Uw() const;
359 
360  //- Return the film thickness*density (helper field) [kg/m^3]
361  virtual const volScalarField& deltaRho() const;
362 
363  //- Return the film flux [kg m/s]
364  virtual const surfaceScalarField& phi() const;
365 
366  //- Return the film density [kg/m^3]
367  virtual const volScalarField& rho() const;
368 
369  //- Return the film mean temperature [K]
370  virtual const volScalarField& T() const;
371 
372  //- Return the film surface temperature [K]
373  virtual const volScalarField& Ts() const;
374 
375  //- Return the film wall temperature [K]
376  virtual const volScalarField& Tw() const;
377 
378  //- Return the film surface enthalpy [J/kg]
379  virtual const volScalarField& hs() const;
380 
381  //- Return the film specific heat capacity [J/kg/K]
382  virtual const volScalarField& Cp() const;
383 
384  //- Return the film thermal conductivity [W/m/K]
385  virtual const volScalarField& kappa() const;
386 
387 
388  // Transfer fields - to the primary region
389 
390  //- Return mass transfer source - Eulerian phase only
391  virtual tmp<volScalarField> primaryMassTrans() const;
392 
393  //- Return the film mass available for transfer to cloud
394  virtual const volScalarField& cloudMassTrans() const;
395 
396  //- Return the parcel diameters originating from film to cloud
397  virtual const volScalarField& cloudDiameterTrans() const;
398 
399 
400  // External helper functions
401 
402  //- External hook to add sources to the film
403  virtual void addSources
404  (
405  const label patchi, // patchi on primary region
406  const label facei, // facei of patchi
407  const scalar massSource, // [kg]
408  const vector& momentumSource, // [kg m/s] (tang'l momentum)
409  const scalar pressureSource, // [kg m/s] (normal momentum)
410  const scalar energySource = 0 // [J]
411  );
412 
413 
414  // Source fields (read/write access)
415 
416  // Primary region
417 
418  //- Momementum [kg/m/s^2]
419  inline volVectorField& USpPrimary();
420 
421  //- Pressure [Pa]
422  inline volScalarField& pSpPrimary();
423 
424  //- Mass [kg/m^2/s]
425  inline volScalarField& rhoSpPrimary();
426 
427 
428  // Film region
429 
430  //- Momentum [kg/m/s^2]
431  inline volVectorField& USp();
432 
433  //- Pressure [Pa]
434  inline volScalarField& pSp();
435 
436  //- Mass [kg/m^2/s]
437  inline volScalarField& rhoSp();
438 
439  //- Momentum [kg/m/s^2]
440  inline const volVectorField& USp() const;
441 
442  //- Pressure [Pa]
443  inline const volScalarField& pSp() const;
444 
445  //- Mass [kg/m^2/s]
446  inline const volScalarField& rhoSp() const;
447 
448 
449  // Fields mapped from primary region
450 
451  //- Velocity [m/s]
452  inline const volVectorField& UPrimary() const;
453 
454  //- Pressure [Pa]
455  inline const volScalarField& pPrimary() const;
456 
457  //- Density [kg/m^3]
458  inline const volScalarField& rhoPrimary() const;
459 
460  //- Viscosity [Pa.s]
461  inline const volScalarField& muPrimary() const;
462 
463 
464  // Sub-models
465 
466  //- Film thermo
467  inline const filmThermoModel& filmThermo() const;
468 
469  //- Injection
470  inline injectionModelList& injection();
471 
472  //- Transfer
473  inline transferModelList& transfer();
474 
475  //- Turbulence
476  inline const filmTurbulenceModel& turbulence() const;
477 
478 
479  // Helper functions
480 
481  //- Return the current film mass
482  inline tmp<volScalarField> mass() const;
483 
484  //- Return the change in film mass due to sources/sinks
485  inline tmp<volScalarField> deltaMass() const;
486 
487  //- Return the gravity normal-to-patch component contribution
488  inline tmp<volScalarField> gNorm() const;
489 
490  //- Return the gravity normal-to-patch component contribution
491  // Clipped so that only non-zero if g & nHat_ < 0
492  inline tmp<volScalarField> gNormClipped() const;
493 
494  //- Return the gravity tangential component contributions
495  inline tmp<volVectorField> gTan() const;
496 
497 
498  // Evolution
499 
500  //- Pre-evolve film hook
501  virtual void preEvolveRegion();
502 
503  //- Evolve the film equations
504  virtual void evolveRegion();
505 
506 
507  // Source fields
508 
509  // Mapped into primary region
510 
511  //- Return total mass source - Eulerian phase only
512  virtual tmp<volScalarField::Internal> Srho() const;
513 
514  //- Return mass source for specie i - Eulerian phase only
516  (
517  const label i
518  ) const;
519 
520  //- Return enthalpy source - Eulerian phase only
521  virtual tmp<volScalarField::Internal> Sh() const;
522 
523 
524  // I-O
525 
526  //- Provide some feedback
527  virtual void info();
528 
529 
530  // Member Operators
531 
532  //- Disallow default bitwise assignment
533  void operator=(const kinematicSingleLayer&) = delete;
534 };
535 
536 
537 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
538 
539 } // End namespace surfaceFilmModels
540 } // End namespace regionModels
541 } // End namespace Foam
542 
543 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
544 
545 #ifdef NoRepository
547 #endif
548 
549 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
550 
551 #include "kinematicSingleLayerI.H"
552 
553 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
554 
555 #endif
556 
557 // ************************************************************************* //
virtual bool read()
Read control parameters from dictionary.
Foam::surfaceFields.
autoPtr< filmThermoModel > filmThermo_
Film thermo model.
List container for film injection models.
virtual const volVectorField & U() const
Return the film velocity [m/s].
autoPtr< filmTurbulenceModel > turbulence_
Turbulence model.
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 const volScalarField & deltaRho() const
Return the film thickness*density (helper field) [kg/m^3].
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
volScalarField primaryMassTrans_
Film mass available for transfer to the primary region.
Kinematic form of single-cell layer surface film model.
label nOuterCorr() const
Return the number of outer correctors.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
label nCorr() const
Return the number of PISO correctors.
void operator=(const kinematicSingleLayer &)=delete
Disallow default bitwise assignment.
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
List container for film sources.
Definition: forceList.H:53
virtual const volVectorField & Us() const
Return the film surface velocity [m/s].
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
virtual tmp< volScalarField::Internal > Srho() const
Return total mass source - Eulerian phase only.
tmp< volScalarField > gNorm() const
Return the gravity normal-to-patch component contribution.
scalarField availableMass_
Available mass for transfer via sub-models.
volScalarField cloudDiameterTrans_
Parcel diameters originating from film to cloud.
scalar addedMassTotal_
Cumulative mass added via sources [kg].
const volScalarField & rhoPrimary() const
Density [kg/m^3].
virtual tmp< volScalarField > primaryMassTrans() const
Return mass transfer source - Eulerian phase only.
const Switch & momentumPredictor() const
Return the momentum predictor.
virtual void solveThickness(const volScalarField &pu, const volScalarField &pp, const fvVectorMatrix &UEqn)
Solve coupled velocity-thickness equations.
virtual void evolveRegion()
Evolve the film equations.
virtual void solveContinuity()
Solve continuity equation.
tmp< volScalarField > deltaMass() const
Return the change in film mass due to sources/sinks.
virtual tmp< volScalarField > pp()
Implicit pressure source coefficient.
volScalarField deltaRho_
Film thickness*density (helper field) [kg/m^2].
const volScalarField & muPrimary() const
Viscosity [Pa.s].
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
const filmThermoModel & filmThermo() const
Film thermo.
const volVectorField & UPrimary() const
Velocity [m/s].
dynamicFvMesh & mesh
const volScalarField & mu() const
Return const access to the dynamic viscosity [Pa.s].
virtual tmp< fvVectorMatrix > solveMomentum(const volScalarField &pu, const volScalarField &pp)
Solve for film velocity.
virtual const volScalarField & cloudMassTrans() const
Return the film mass available for transfer to cloud.
tmp< volVectorField > gTan() const
Return the gravity tangential component contributions.
tmp< volScalarField > mass() const
Return the current film mass.
A class for handling words, derived from string.
Definition: word.H:59
virtual const volScalarField & cloudDiameterTrans() const
Return the parcel diameters originating from film to cloud.
tmp< volScalarField > gNormClipped() const
Return the gravity normal-to-patch component contribution.
const dimensionedScalar & deltaSmall() const
Return small delta.
virtual const surfaceScalarField & phi() const
Return the film flux [kg m/s].
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
label nNonOrthCorr_
Number of non-orthogonal correctors.
virtual const volScalarField & kappa() const
Return the film thermal conductivity [W/m/K].
virtual tmp< volScalarField > pu()
Explicit pressure source contribution.
const volScalarField & pPrimary() const
Pressure [Pa].
List container for film transfer models.
const dimensionedVector & g() const
Return the acceleration due to gravity.
virtual void correctAlpha()
Correct film coverage field.
virtual const volScalarField & hs() const
Return the film surface enthalpy [J/kg].
virtual void resetPrimaryRegionSourceTerms()
Reset source term fields.
const volScalarField & delta() const
Return const access to the film thickness [m].
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.
virtual void updateSubmodels()
Update the film sub-models.
label patchi
virtual const volVectorField & Uw() const
Return the film wall velocity [m/s].
surfaceScalarField phi_
Mass flux (includes film thickness) [kg m/s].
const volScalarField & sigma() const
Return const access to the surface tension [kg/s^2].
label nNonOrthCorr() const
Return the number of non-orthogonal correctors.
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 [].
virtual void updateSurfaceVelocities()
Update film surface velocities.
virtual const volScalarField & Tw() const
Return the film wall temperature [K].
A special matrix type and solver, designed for finite volume solutions of scalar equations.
TypeName("kinematicSingleLayer")
Runtime type information.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
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.
virtual const volScalarField & rho() const
Return the film density [kg/m^3].
volScalarField cloudMassTrans_
Film mass available for transfer to cloud.
void constrainFilmField(Type &field, const typename Type::cmptType &value)
Constrain a film region master/slave boundaries of a field to a.
A class for managing temporary objects.
Definition: PtrList.H:53
virtual const volScalarField & Cp() const
Return the film specific heat capacity [J/kg/K].
const filmTurbulenceModel & turbulence() const
Turbulence.
virtual const volScalarField & T() const
Return the film mean temperature [K].
virtual tmp< volScalarField::Internal > Sh() const
Return enthalpy source - Eulerian phase only.
virtual scalar CourantNumber() const
Courant number evaluation.
scalar deltaCoLimit_
Film thickness above which Courant number calculation in valid.
const volScalarField & alpha() const
Return the film coverage, 1 = covered, 0 = uncovered [].
Namespace for OpenFOAM.
virtual void correctThermoFields()
Correct the thermo fields.