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 | Copyright (C) 2011-2015 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::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 "surfaceFilmModel.H"
39 #include "fvMesh.H"
40 #include "volFields.H"
41 #include "surfaceFields.H"
42 #include "fvMatrices.H"
43 
44 #include "injectionModelList.H"
45 #include "forceList.H"
46 #include "filmTurbulenceModel.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 namespace regionModels
53 {
54 namespace surfaceFilmModels
55 {
56 
57 // Forward class declarations
58 class filmThermoModel;
59 
60 /*---------------------------------------------------------------------------*\
61  Class kinematicSingleLayer Declaration
62 \*---------------------------------------------------------------------------*/
63 
65 :
66  public surfaceFilmModel
67 {
68 private:
69 
70  // Private member functions
71 
72  //- Disallow default bitwise copy construct
74 
75  //- Disallow default bitwise assignment
76  void operator=(const kinematicSingleLayer&);
77 
78 
79 protected:
80 
81  // Protected data
82 
83  // Solution parameters
84 
85  //- Momentum predictor
87 
88  //- Number of outer correctors
90 
91  //- Number of PISO-like correctors
92  label nCorr_;
93 
94  //- Number of non-orthogonal correctors
96 
97  //- Cumulative continuity error
98  scalar cumulativeContErr_;
99 
100  //- Small delta
102 
103  //- Film thickness above which Courant number calculation in valid
104  scalar deltaCoLimit_;
105 
106 
107  // Thermo properties
108 
109  // Fields
110 
111  //- Density / [kg/m3]
113 
114  //- Dynamic viscosity / [Pa.s]
116 
117  //- Surface tension / [m/s2]
119 
120 
121  // Fields
122 
123  //- Film thickness / [m]
125 
126  //- Film coverage indicator, 1 = covered, 0 = uncovered / []
128 
129  //- Velocity - mean / [m/s]
131 
132  //- Velocity - surface / [m/s]
134 
135  //- Velocity - wall / [m/s]
137 
138  //- Film thickness*density (helper field) / [kg/m2]
140 
141  //- Mass flux (includes film thickness) / [kg.m/s]
143 
144 
145  // Transfer fields
146 
147  //- Film mass available for transfer to the primary region
149 
150  //- Film mass available for transfer to cloud
152 
153  //- Parcel diameters originating from film to cloud
155 
156 
157  // Source term fields
158 
159  // Film region - registered to the film region mesh
160  // Note: need boundary value mapped from primary region, and then
161  // pushed into the patch internal field
162 
163  //- Momementum / [kg/m/s2]
165 
166  //- Pressure / [Pa]
168 
169  //- Mass / [kg/m2/s]
171 
172 
173  // Primary region - registered to the primary region mesh
174  // Internal use only - not read-in
175 
176  //- Momementum / [kg/m/s2]
178 
179  //- Pressure / [Pa]
181 
182  //- Mass / [kg/m2/s]
184 
185 
186  // Fields mapped from primary region - registered to the film region
187  // Note: need both boundary AND patch internal fields to be mapped
188 
189  //- Velocity / [m/s]
191 
192  //- Pressure / [Pa]
194 
195  //- Density / [kg/m3]
197 
198  //- Viscosity / [Pa.s]
200 
201 
202  // Sub-models
203 
204  //- Film thermo model
206 
207  //- Available mass for transfer via sub-models
209 
210  //- Cloud injection
212 
213  //- Turbulence model
215 
216  //- List of film forces
218 
219 
220  // Checks
221 
222  //- Cumulative mass added via sources [kg]
223  scalar addedMassTotal_;
224 
225 
226  // Protected member functions
227 
228  //- Read control parameters from dictionary
229  virtual bool read();
230 
231  //- Correct the thermo fields
232  virtual void correctThermoFields();
233 
234  //- Reset source term fields
235  virtual void resetPrimaryRegionSourceTerms();
236 
237  //- Transfer thermo fields from the primary region to the film region
238  virtual void transferPrimaryRegionThermoFields();
239 
240  //- Transfer source fields from the primary region to the film region
241  virtual void transferPrimaryRegionSourceFields();
242 
243  //- Explicit pressure source contribution
244  virtual tmp<volScalarField> pu();
245 
246  //- Implicit pressure source coefficient
247  virtual tmp<volScalarField> pp();
248 
249  //- Correct film coverage field
250  virtual void correctAlpha();
251 
252  //- Update the film sub-models
253  virtual void updateSubmodels();
254 
255  //- Continuity check
256  virtual void continuityCheck();
257 
258  //- Update film surface velocities
259  virtual void updateSurfaceVelocities();
260 
261  //- Constrain a film region master/slave boundaries of a field to a
262  // given value
263  template<class Type>
264  void constrainFilmField
265  (
266  Type& field,
267  const typename Type::cmptType& value
268  );
269 
270 
271  // Equations
272 
273  //- Solve continuity equation
274  virtual void solveContinuity();
275 
276  //- Solve for film velocity
278  (
279  const volScalarField& pu,
280  const volScalarField& pp
281  );
282 
283  //- Solve coupled velocity-thickness equations
284  virtual void solveThickness
285  (
286  const volScalarField& pu,
287  const volScalarField& pp,
288  const fvVectorMatrix& UEqn
289  );
290 
291 
292 public:
293 
294  //- Runtime type information
295  TypeName("kinematicSingleLayer");
296 
297 
298  // Constructors
299 
300  //- Construct from components
302  (
303  const word& modelType,
304  const fvMesh& mesh,
305  const dimensionedVector& g,
306  const word& regionType,
307  const bool readFields = true
308  );
309 
310 
311  //- Destructor
312  virtual ~kinematicSingleLayer();
313 
314 
315  // Member Functions
316 
317  // Solution parameters
318 
319  //- Courant number evaluation
320  virtual scalar CourantNumber() const;
321 
322  //- Return the momentum predictor
323  inline const Switch& momentumPredictor() const;
324 
325  //- Return the number of outer correctors
326  inline label nOuterCorr() const;
327 
328  //- Return the number of PISO correctors
329  inline label nCorr() const;
330 
331  //- Return the number of non-orthogonal correctors
332  inline label nNonOrthCorr() const;
333 
334  //- Return small delta
335  inline const dimensionedScalar& deltaSmall() const;
336 
337 
338  // Thermo properties
339 
340  //- Return const access to the dynamic viscosity / [Pa.s]
341  inline const volScalarField& mu() const;
342 
343  //- Return const access to the surface tension / [m/s2]
344  inline const volScalarField& sigma() const;
345 
346 
347  // Fields
348 
349  //- Return const access to the film thickness / [m]
350  inline const volScalarField& delta() const;
351 
352  //- Return the film coverage, 1 = covered, 0 = uncovered / []
353  inline const volScalarField& alpha() const;
354 
355  //- Return the film velocity [m/s]
356  virtual const volVectorField& U() const;
357 
358  //- Return the film surface velocity [m/s]
359  virtual const volVectorField& Us() const;
360 
361  //- Return the film wall velocity [m/s]
362  virtual const volVectorField& Uw() const;
363 
364  //- Return the film thickness*density (helper field) [kg/m3]
365  virtual const volScalarField& deltaRho() const;
366 
367  //- Return the film flux [kg.m/s]
368  virtual const surfaceScalarField& phi() const;
369 
370  //- Return the film density [kg/m3]
371  virtual const volScalarField& rho() const;
372 
373  //- Return the film mean temperature [K]
374  virtual const volScalarField& T() const;
375 
376  //- Return the film surface temperature [K]
377  virtual const volScalarField& Ts() const;
378 
379  //- Return the film wall temperature [K]
380  virtual const volScalarField& Tw() const;
381 
382  //- Return the film specific heat capacity [J/kg/K]
383  virtual const volScalarField& Cp() const;
384 
385  //- Return the film thermal conductivity [W/m/K]
386  virtual const volScalarField& kappa() const;
387 
388 
389  // Transfer fields - to the primary region
390 
391  //- Return mass transfer source - Eulerian phase only
392  virtual tmp<volScalarField> primaryMassTrans() const;
393 
394  //- Return the film mass available for transfer to cloud
395  virtual const volScalarField& cloudMassTrans() const;
396 
397  //- Return the parcel diameters originating from film to cloud
398  virtual const volScalarField& cloudDiameterTrans() const;
399 
400 
401  // External helper functions
402 
403  //- External hook to add sources to the film
404  virtual void addSources
405  (
406  const label patchI, // patchI on primary region
407  const label faceI, // faceI of patchI
408  const scalar massSource, // [kg]
409  const vector& momentumSource, // [kg.m/s] (tang'l momentum)
410  const scalar pressureSource, // [kg.m/s] (normal momentum)
411  const scalar energySource = 0 // [J]
412  );
413 
414 
415  // Source fields (read/write access)
416 
417  // Primary region
418 
419  //- Momementum / [kg/m/s2]
420  inline volVectorField& USpPrimary();
421 
422  //- Pressure / [Pa]
423  inline volScalarField& pSpPrimary();
424 
425  //- Mass / [kg/m2/s]
426  inline volScalarField& rhoSpPrimary();
427 
428 
429  // Film region
430 
431  //- Momentum / [kg/m/s2]
432  inline volVectorField& USp();
433 
434  //- Pressure / [Pa]
435  inline volScalarField& pSp();
436 
437  //- Mass / [kg/m2/s]
438  inline volScalarField& rhoSp();
439 
440  //- Momentum / [kg/m/s2]
441  inline const volVectorField& USp() const;
442 
443  //- Pressure / [Pa]
444  inline const volScalarField& pSp() const;
445 
446  //- Mass / [kg/m2/s]
447  inline const volScalarField& rhoSp() const;
448 
449 
450  // Fields mapped from primary region
451 
452  //- Velocity / [m/s]
453  inline const volVectorField& UPrimary() const;
454 
455  //- Pressure / [Pa]
456  inline const volScalarField& pPrimary() const;
457 
458  //- Density / [kg/m3]
459  inline const volScalarField& rhoPrimary() const;
460 
461  //- Viscosity / [Pa.s]
462  inline const volScalarField& muPrimary() const;
463 
464 
465  // Sub-models
466 
467  //- Film thermo
468  inline const filmThermoModel& filmThermo() const;
469 
470  //- Injection
471  inline injectionModelList& injection();
472 
473  //- Turbulence
474  inline const filmTurbulenceModel& turbulence() const;
475 
476 
477  // Helper functions
478 
479  //- Return the current film mass
480  inline tmp<volScalarField> mass() const;
481 
482  //- Return the net film mass available over the next integration
483  inline tmp<volScalarField> netMass() const;
484 
485  //- Return the change in film mass due to sources/sinks
486  inline tmp<volScalarField> deltaMass() const;
487 
488  //- Return the gravity normal-to-patch component contribution
489  inline tmp<volScalarField> gNorm() const;
490 
491  //- Return the gravity normal-to-patch component contribution
492  // Clipped so that only non-zero if g & nHat_ < 0
493  inline tmp<volScalarField> gNormClipped() const;
494 
495  //- Return the gravity tangential component contributions
496  inline tmp<volVectorField> gTan() const;
497 
498  //- Return the gravity tangential component contributions for patchI
499  inline tmp<vectorField> gTan(const label patchI) const;
500 
501 
502  // Evolution
503 
504  //- Pre-evolve film hook
505  virtual void preEvolveRegion();
506 
507  //- Evolve the film equations
508  virtual void evolveRegion();
509 
510  //- Post-evolve film hook
511  virtual void postEvolveRegion();
512 
513 
514  // Source fields
515 
516  // Mapped into primary region
517 
518  //- Return total mass source - Eulerian phase only
520 
521  //- Return mass source for specie i - Eulerian phase only
523  (
524  const label i
525  ) const;
526 
527  //- Return enthalpy source - Eulerian phase only
528  virtual tmp<DimensionedField<scalar, volMesh> > Sh() const;
529 
530 
531  // I-O
532 
533  //- Provide some feedback
534  virtual void info();
535 };
536 
537 
538 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
539 
540 } // End namespace surfaceFilmModels
541 } // End namespace regionModels
542 } // End namespace Foam
543 
544 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
545 
546 #ifdef NoRepository
548 #endif
549 
550 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
551 
552 #include "kinematicSingleLayerI.H"
553 
554 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
555 
556 #endif
557 
558 // ************************************************************************* //
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::T
virtual const volScalarField & T() const
Return the film mean temperature [K].
Definition: kinematicSingleLayer.C:1006
volFields.H
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::read
virtual bool read()
Read control parameters from dictionary.
Definition: kinematicSingleLayer.C:55
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::correctThermoFields
virtual void correctThermoFields()
Correct the thermo fields.
Definition: kinematicSingleLayer.C:74
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoSpPrimary_
volScalarField rhoSpPrimary_
Mass / [kg/m2/s].
Definition: kinematicSingleLayer.H:182
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::deltaCoLimit_
scalar deltaCoLimit_
Film thickness above which Courant number calculation in valid.
Definition: kinematicSingleLayer.H:103
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::~kinematicSingleLayer
virtual ~kinematicSingleLayer()
Destructor.
Definition: kinematicSingleLayer.C:831
injectionModelList.H
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:60
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::transferPrimaryRegionThermoFields
virtual void transferPrimaryRegionThermoFields()
Transfer thermo fields from the primary region to the film region.
Definition: kinematicSingleLayer.C:95
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Us_
volVectorField Us_
Velocity - surface / [m/s].
Definition: kinematicSingleLayer.H:132
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Srho
virtual tmp< DimensionedField< scalar, volMesh > > Srho() const
Return total mass source - Eulerian phase only.
Definition: kinematicSingleLayer.C:1109
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::USpPrimary_
volVectorField USpPrimary_
Momementum / [kg/m/s2].
Definition: kinematicSingleLayer.H:176
kinematicSingleLayerTemplates.C
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::USp_
volVectorField USp_
Momementum / [kg/m/s2].
Definition: kinematicSingleLayer.H:163
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::turbulence
const filmTurbulenceModel & turbulence() const
Turbulence.
Definition: kinematicSingleLayerI.H:185
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Us
virtual const volVectorField & Us() const
Return the film surface velocity [m/s].
Definition: kinematicSingleLayer.C:976
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::injection_
injectionModelList injection_
Cloud injection.
Definition: kinematicSingleLayer.H:210
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::preEvolveRegion
virtual void preEvolveRegion()
Pre-evolve film hook.
Definition: kinematicSingleLayer.C:863
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::primaryMassTrans_
volScalarField primaryMassTrans_
Film mass available for transfer to the primary region.
Definition: kinematicSingleLayer.H:147
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Ts
virtual const volScalarField & Ts() const
Return the film surface temperature [K].
Definition: kinematicSingleLayer.C:1015
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::turbulence_
autoPtr< filmTurbulenceModel > turbulence_
Turbulence model.
Definition: kinematicSingleLayer.H:213
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::filmThermo_
autoPtr< filmThermoModel > filmThermo_
Film thermo model.
Definition: kinematicSingleLayer.H:204
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::U
virtual const volVectorField & U() const
Return the film velocity [m/s].
Definition: kinematicSingleLayer.C:970
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::deltaSmall_
const dimensionedScalar deltaSmall_
Small delta.
Definition: kinematicSingleLayer.H:100
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::evolveRegion
virtual void evolveRegion()
Evolve the film equations.
Definition: kinematicSingleLayer.C:890
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::solveContinuity
virtual void solveContinuity()
Solve continuity equation.
Definition: kinematicSingleLayer.C:258
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::operator=
void operator=(const kinematicSingleLayer &)
Disallow default bitwise assignment.
Foam::regionModels::surfaceFilmModels::filmTurbulenceModel
Base class for film turbulence models.
Definition: filmTurbulenceModel.H:56
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::nCorr
label nCorr() const
Return the number of PISO correctors.
Definition: kinematicSingleLayerI.H:53
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer
Definition: kinematicSingleLayer.H:63
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::solveThickness
virtual void solveThickness(const volScalarField &pu, const volScalarField &pp, const fvVectorMatrix &UEqn)
Solve coupled velocity-thickness equations.
Definition: kinematicSingleLayer.C:351
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rho_
volScalarField rho_
Density / [kg/m3].
Definition: kinematicSingleLayer.H:111
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::mu
const volScalarField & mu() const
Return const access to the dynamic viscosity / [Pa.s].
Definition: kinematicSingleLayerI.H:71
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::sigma_
volScalarField sigma_
Surface tension / [m/s2].
Definition: kinematicSingleLayer.H:117
surfaceFields.H
Foam::surfaceFields.
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::continuityCheck
virtual void continuityCheck()
Continuity check.
Definition: kinematicSingleLayer.C:225
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::UPrimary
const volVectorField & UPrimary() const
Velocity / [m/s].
Definition: kinematicSingleLayerI.H:149
Foam::readFields
This function object reads fields from the time directories and adds them to the mesh database for fu...
Definition: readFields.H:104
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::availableMass_
scalarField availableMass_
Available mass for transfer via sub-models.
Definition: kinematicSingleLayer.H:207
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pPrimary_
volScalarField pPrimary_
Pressure / [Pa].
Definition: kinematicSingleLayer.H:192
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::deltaSmall
const dimensionedScalar & deltaSmall() const
Return small delta.
Definition: kinematicSingleLayerI.H:65
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::muPrimary_
volScalarField muPrimary_
Viscosity / [Pa.s].
Definition: kinematicSingleLayer.H:198
Foam::regionModels::surfaceFilmModels::surfaceFilmModel::g
const dimensionedVector & g() const
Return the accleration due to gravity.
Definition: surfaceFilmModelI.H:39
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::gNormClipped
tmp< volScalarField > gNormClipped() const
Return the gravity normal-to-patch component contribution.
Definition: kinematicSingleLayerI.H:235
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoPrimary
const volScalarField & rhoPrimary() const
Density / [kg/m3].
Definition: kinematicSingleLayerI.H:161
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::injection
injectionModelList & injection()
Injection.
Definition: kinematicSingleLayerI.H:179
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::cloudDiameterTrans_
volScalarField cloudDiameterTrans_
Parcel diameters originating from film to cloud.
Definition: kinematicSingleLayer.H:153
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::deltaRho_
volScalarField deltaRho_
Film thickness*density (helper field) / [kg/m2].
Definition: kinematicSingleLayer.H:138
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::cloudMassTrans
virtual const volScalarField & cloudMassTrans() const
Return the film mass available for transfer to cloud.
Definition: kinematicSingleLayer.C:1073
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::addedMassTotal_
scalar addedMassTotal_
Cumulative mass added via sources [kg].
Definition: kinematicSingleLayer.H:222
Foam::label
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
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::transferPrimaryRegionSourceFields
virtual void transferPrimaryRegionSourceFields()
Transfer source fields from the primary region to the film region.
Definition: kinematicSingleLayer.C:112
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::gNorm
tmp< volScalarField > gNorm() const
Return the gravity normal-to-patch component contribution.
Definition: kinematicSingleLayerI.H:212
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::delta_
volScalarField delta_
Film thickness / [m].
Definition: kinematicSingleLayer.H:123
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::deltaMass
tmp< volScalarField > deltaMass() const
Return the change in film mass due to sources/sinks.
Definition: kinematicSingleLayerI.H:206
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::USpPrimary
volVectorField & USpPrimary()
Momementum / [kg/m/s2].
Definition: kinematicSingleLayerI.H:95
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::filmThermo
const filmThermoModel & filmThermo() const
Film thermo.
Definition: kinematicSingleLayerI.H:173
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::kappa
virtual const volScalarField & kappa() const
Return the film thermal conductivity [W/m/K].
Definition: kinematicSingleLayer.C:1042
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pp
virtual tmp< volScalarField > pp()
Implicit pressure source coefficient.
Definition: kinematicSingleLayer.C:181
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::cloudDiameterTrans
virtual const volScalarField & cloudDiameterTrans() const
Return the parcel diameters originating from film to cloud.
Definition: kinematicSingleLayer.C:1079
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::gTan
tmp< volVectorField > gTan() const
Return the gravity tangential component contributions.
Definition: kinematicSingleLayerI.H:261
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::momentumPredictor
const Switch & momentumPredictor() const
Return the momentum predictor.
Definition: kinematicSingleLayerI.H:41
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::info
virtual void info()
Provide some feedback.
Definition: kinematicSingleLayer.C:1085
UEqn
tmp< fvVectorMatrix > UEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
filmTurbulenceModel.H
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Uw_
volVectorField Uw_
Velocity - wall / [m/s].
Definition: kinematicSingleLayer.H:135
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::primaryMassTrans
virtual tmp< volScalarField > primaryMassTrans() const
Return mass transfer source - Eulerian phase only.
Definition: kinematicSingleLayer.C:1051
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::resetPrimaryRegionSourceTerms
virtual void resetPrimaryRegionSourceTerms()
Reset source term fields.
Definition: kinematicSingleLayer.C:82
surfaceFilmModel.H
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pPrimary
const volScalarField & pPrimary() const
Pressure / [Pa].
Definition: kinematicSingleLayerI.H:155
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::solveMomentum
virtual tmp< fvVectorMatrix > solveMomentum(const volScalarField &pu, const volScalarField &pp)
Solve for film velocity.
Definition: kinematicSingleLayer.C:293
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::USp
volVectorField & USp()
Momentum / [kg/m/s2].
Definition: kinematicSingleLayerI.H:113
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::phi
virtual const surfaceScalarField & phi() const
Return the film flux [kg.m/s].
Definition: kinematicSingleLayer.C:994
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::forces_
forceList forces_
List of film forces.
Definition: kinematicSingleLayer.H:216
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::muPrimary
const volScalarField & muPrimary() const
Viscosity / [Pa.s].
Definition: kinematicSingleLayerI.H:167
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::correctAlpha
virtual void correctAlpha()
Correct film coverage field.
Definition: kinematicSingleLayer.C:201
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::momentumPredictor_
Switch momentumPredictor_
Momentum predictor.
Definition: kinematicSingleLayer.H:85
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::addSources
virtual void addSources(const label patchI, const label faceI, const scalar massSource, const vector &momentumSource, const scalar pressureSource, const scalar energySource=0)
External hook to add sources to the film.
Definition: kinematicSingleLayer.C:838
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::U_
volVectorField U_
Velocity - mean / [m/s].
Definition: kinematicSingleLayer.H:129
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::nNonOrthCorr
label nNonOrthCorr() const
Return the number of non-orthogonal correctors.
Definition: kinematicSingleLayerI.H:59
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::nCorr_
label nCorr_
Number of PISO-like correctors.
Definition: kinematicSingleLayer.H:91
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::mass
tmp< volScalarField > mass() const
Return the current film mass.
Definition: kinematicSingleLayerI.H:191
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pu
virtual tmp< volScalarField > pu()
Explicit pressure source contribution.
Definition: kinematicSingleLayer.C:159
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::UPrimary_
volVectorField UPrimary_
Velocity / [m/s].
Definition: kinematicSingleLayer.H:189
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pSp_
volScalarField pSp_
Pressure / [Pa].
Definition: kinematicSingleLayer.H:166
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::cumulativeContErr_
scalar cumulativeContErr_
Cumulative continuity error.
Definition: kinematicSingleLayer.H:97
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::phi_
surfaceScalarField phi_
Mass flux (includes film thickness) / [kg.m/s].
Definition: kinematicSingleLayer.H:141
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::sigma
const volScalarField & sigma() const
Return const access to the surface tension / [m/s2].
Definition: kinematicSingleLayerI.H:77
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pSp
volScalarField & pSp()
Pressure / [Pa].
Definition: kinematicSingleLayerI.H:119
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::constrainFilmField
void constrainFilmField(Type &field, const typename Type::cmptType &value)
Constrain a film region master/slave boundaries of a field to a.
Definition: kinematicSingleLayerTemplates.C:41
Foam::regionModels::surfaceFilmModels::injectionModelList
Definition: injectionModelList.H:54
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoSpPrimary
volScalarField & rhoSpPrimary()
Mass / [kg/m2/s].
Definition: kinematicSingleLayerI.H:107
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
kinematicSingleLayerI.H
Foam::regionModels::surfaceFilmModels::filmThermoModel
Definition: filmThermoModel.H:54
filmThermoModel
Base class for film thermo models.
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::alpha_
volScalarField alpha_
Film coverage indicator, 1 = covered, 0 = uncovered / [].
Definition: kinematicSingleLayer.H:126
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pSpPrimary_
volScalarField pSpPrimary_
Pressure / [Pa].
Definition: kinematicSingleLayer.H:179
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::mu_
volScalarField mu_
Dynamic viscosity / [Pa.s].
Definition: kinematicSingleLayer.H:114
Foam::Vector< scalar >
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::nNonOrthCorr_
label nNonOrthCorr_
Number of non-orthogonal correctors.
Definition: kinematicSingleLayer.H:94
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Uw
virtual const volVectorField & Uw() const
Return the film wall velocity [m/s].
Definition: kinematicSingleLayer.C:982
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoPrimary_
volScalarField rhoPrimary_
Density / [kg/m3].
Definition: kinematicSingleLayer.H:195
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::delta
const volScalarField & delta() const
Return const access to the film thickness / [m].
Definition: kinematicSingleLayerI.H:83
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::updateSurfaceVelocities
virtual void updateSurfaceVelocities()
Update film surface velocities.
Definition: kinematicSingleLayer.C:275
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::alpha
const volScalarField & alpha() const
Return the film coverage, 1 = covered, 0 = uncovered / [].
Definition: kinematicSingleLayerI.H:89
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::netMass
tmp< volScalarField > netMass() const
Return the net film mass available over the next integration.
Definition: kinematicSingleLayerI.H:197
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Cp
virtual const volScalarField & Cp() const
Return the film specific heat capacity [J/kg/K].
Definition: kinematicSingleLayer.C:1033
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::kinematicSingleLayer
kinematicSingleLayer(const kinematicSingleLayer &)
Disallow default bitwise copy construct.
forceList.H
Foam::regionModels::surfaceFilmModels::surfaceFilmModel
Base class for surface film models.
Definition: surfaceFilmModel.H:61
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Tw
virtual const volScalarField & Tw() const
Return the film wall temperature [K].
Definition: kinematicSingleLayer.C:1024
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::postEvolveRegion
virtual void postEvolveRegion()
Post-evolve film hook.
Definition: kinematicSingleLayer.C:927
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rho
virtual const volScalarField & rho() const
Return the film density [kg/m3].
Definition: kinematicSingleLayer.C:1000
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::deltaRho
virtual const volScalarField & deltaRho() const
Return the film thickness*density (helper field) [kg/m3].
Definition: kinematicSingleLayer.C:988
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::updateSubmodels
virtual void updateSubmodels()
Update the film sub-models.
Definition: kinematicSingleLayer.C:207
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::TypeName
TypeName("kinematicSingleLayer")
Runtime type information.
Foam::regionModels::surfaceFilmModels::forceList
Definition: forceList.H:53
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoSp_
volScalarField rhoSp_
Mass / [kg/m2/s].
Definition: kinematicSingleLayer.H:169
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::rhoSp
volScalarField & rhoSp()
Mass / [kg/m2/s].
Definition: kinematicSingleLayerI.H:125
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::cloudMassTrans_
volScalarField cloudMassTrans_
Film mass available for transfer to cloud.
Definition: kinematicSingleLayer.H:150
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::Sh
virtual tmp< DimensionedField< scalar, volMesh > > Sh() const
Return enthalpy source - Eulerian phase only.
Definition: kinematicSingleLayer.C:1156
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::pSpPrimary
volScalarField & pSpPrimary()
Pressure / [Pa].
Definition: kinematicSingleLayerI.H:101
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::nOuterCorr_
label nOuterCorr_
Number of outer correctors.
Definition: kinematicSingleLayer.H:88
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::CourantNumber
virtual scalar CourantNumber() const
Courant number evaluation.
Definition: kinematicSingleLayer.C:939
Foam::regionModels::surfaceFilmModels::kinematicSingleLayer::nOuterCorr
label nOuterCorr() const
Return the number of outer correctors.
Definition: kinematicSingleLayerI.H:47