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