BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / src / regionModels / surfaceFilmModels / kinematicSingleLayer / kinematicSingleLayer.H
blobb084f885d48905bd801862d658605745658a203f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 Class
25     Foam::kinematicSingleLayer
27 Description
28     Kinematic form of single-cell layer surface film model
30 SourceFiles
31     kinematicSingleLayer.C
33 \*---------------------------------------------------------------------------*/
35 #ifndef kinematicSingleLayer_H
36 #define kinematicSingleLayer_H
38 #include "surfaceFilmModel.H"
39 #include "fvMesh.H"
40 #include "volFields.H"
41 #include "surfaceFields.H"
42 #include "fvMatrices.H"
44 #include "injectionModelList.H"
45 #include "forceList.H"
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 namespace Foam
51 namespace regionModels
53 namespace surfaceFilmModels
56 /*---------------------------------------------------------------------------*\
57                    Class kinematicSingleLayer Declaration
58 \*---------------------------------------------------------------------------*/
60 class kinematicSingleLayer
62     public surfaceFilmModel
64 private:
66     // Private member functions
68         //- Disallow default bitwise copy construct
69         kinematicSingleLayer(const kinematicSingleLayer&);
71         //- Disallow default bitwise assignment
72         void operator=(const kinematicSingleLayer&);
75 protected:
77     // Protected data
79         // Solution parameters
81             //- Momentum predictor
82             Switch momentumPredictor_;
84             //- Number of outer correctors
85             label nOuterCorr_;
87             //- Number of PISO-like correctors
88             label nCorr_;
90             //- Number of non-orthogonal correctors
91             label nNonOrthCorr_;
93             //- Cumulative continuity error
94             scalar cumulativeContErr_;
97         // Thermo properties
99             // Fields
101                 //- Density / [kg/m3]
102                 volScalarField rho_;
104                 //- Dynamic viscosity / [Pa.s]
105                 volScalarField mu_;
107                 //- Surface tension / [m/s2]
108                 volScalarField sigma_;
111         // Fields
113             //- Film thickness / [m]
114             volScalarField delta_;
116             //- Velocity - mean / [m/s]
117             volVectorField U_;
119             //- Velocity - surface / [m/s]
120             volVectorField Us_;
122             //- Velocity - wall / [m/s]
123             volVectorField Uw_;
125             //- Film thickness*density (helper field) / [kg/m2]
126             volScalarField deltaRho_;
128             //- Mass flux (includes film thickness) / [kg.m/s]
129             surfaceScalarField phi_;
132             // Transfer fields
134                 //- Film mass available for transfer to the primary region
135                 volScalarField primaryMassTrans_;
137                 //- Film mass available for transfer to cloud
138                 volScalarField cloudMassTrans_;
140                 //- Parcel diameters originating from film to cloud
141                 volScalarField cloudDiameterTrans_;
144         // Source term fields
146             // Film region - registered to the film region mesh
147             // Note: need boundary value mapped from primary region, and then
148             // pushed into the patch internal field
150                 //- Momementum / [kg/m/s2]
151                 volVectorField USp_;
153                 //- Pressure / [Pa]
154                 volScalarField pSp_;
156                 //- Mass / [kg/m2/s]
157                 volScalarField rhoSp_;
160             // Primary region - registered to the primary region mesh
161             // Internal use only - not read-in
163                 //- Momementum / [kg/m/s2]
164                 volVectorField USpPrimary_;
166                 //- Pressure / [Pa]
167                 volScalarField pSpPrimary_;
169                 //- Mass / [kg/m2/s]
170                 volScalarField rhoSpPrimary_;
173         // Fields mapped from primary region - registered to the film region
174         // Note: need both boundary AND patch internal fields to be mapped
176             //- Velocity / [m/s]
177             volVectorField UPrimary_;
179             //- Pressure / [Pa]
180             volScalarField pPrimary_;
182             //- Density / [kg/m3]
183             volScalarField rhoPrimary_;
185             //- Viscosity / [Pa.s]
186             volScalarField muPrimary_;
189         // Sub-models
191             //- Available mass for transfer via sub-models
192             scalarField availableMass_;
194             //- Cloud injection
195             injectionModelList injection_;
197             //- List of film forces
198             forceList forces_;
201        // Checks
203            //- Cumulative mass added via sources [kg]
204            scalar addedMassTotal_;
207     // Protected member functions
209         //- Read control parameters from dictionary
210         virtual bool read();
212         //- Correct the thermo fields
213         virtual void correctThermoFields();
215         //- Reset source term fields
216         virtual void resetPrimaryRegionSourceTerms();
218         //- Transfer thermo fields from the primary region to the film region
219         virtual void transferPrimaryRegionThermoFields();
221         //- Transfer source fields from the primary region to the film region
222         virtual void transferPrimaryRegionSourceFields();
224         // Explicit pressure source contribution
225         virtual tmp<volScalarField> pu();
227         // Implicit pressure source coefficient
228         virtual tmp<volScalarField> pp();
230         //- Update the film sub-models
231         virtual void updateSubmodels();
233         //- Continuity check
234         virtual void continuityCheck();
236         //- Update film surface velocities
237         virtual void updateSurfaceVelocities();
239         //- Constrain a film region master/slave boundaries of a field to a
240         //  given value
241         template<class Type>
242         void constrainFilmField
243         (
244             Type& field,
245             const typename Type::cmptType& value
246         );
249         // Equations
251             //- Solve continuity equation
252             virtual void solveContinuity();
254             //- Solve for film velocity
255             virtual tmp<fvVectorMatrix> solveMomentum
256             (
257                 const volScalarField& pu,
258                 const volScalarField& pp
259             );
261             //- Solve coupled velocity-thickness equations
262             virtual void solveThickness
263             (
264                 const volScalarField& pu,
265                 const volScalarField& pp,
266                 const fvVectorMatrix& UEqn
267             );
270 public:
272     //- Runtime type information
273     TypeName("kinematicSingleLayer");
276     // Constructors
278         //- Construct from components
279         kinematicSingleLayer
280         (
281             const word& modelType,
282             const fvMesh& mesh,
283             const dimensionedVector& g,
284             const bool readFields = true
285         );
288     //- Destructor
289     virtual ~kinematicSingleLayer();
292     // Member Functions
294         // Solution parameters
296             //- Courant number evaluation
297             virtual scalar CourantNumber() const;
299             //- Return the momentum predictor
300             inline const Switch& momentumPredictor() const;
302             //- Return the number of outer correctors
303             inline label nOuterCorr() const;
305             //- Return the number of PISO correctors
306             inline label nCorr() const;
308             //- Return the number of non-orthogonal correctors
309             inline label nNonOrthCorr() const;
312         // Thermo properties
314             //- Return const access to the dynamic viscosity / [Pa.s]
315             inline const volScalarField& mu() const;
317             //- Return const access to the surface tension / [m/s2]
318             inline const volScalarField& sigma() const;
321         // Fields
323             //- Return const access to the film thickness / [m]
324             inline const volScalarField& delta() const;
326             //- Return the film velocity [m/s]
327             virtual const volVectorField& U() const;
329             //- Return the film surface velocity [m/s]
330             virtual const volVectorField& Us() const;
332             //- Return the film wall velocity [m/s]
333             virtual const volVectorField& Uw() const;
335             //- Return the film flux [kg.m/s]
336             virtual const surfaceScalarField& phi() const;
338             //- Return the film density [kg/m3]
339             virtual const volScalarField& rho() const;
341             //- Return the film mean temperature [K]
342             virtual const volScalarField& T() const;
344             //- Return the film surface temperature [K]
345             virtual const volScalarField& Ts() const;
347             //- Return the film wall temperature [K]
348             virtual const volScalarField& Tw() const;
350             //- Return the film specific heat capacity [J/kg/K]
351             virtual const volScalarField& Cp() const;
353             //- Return the film thermal conductivity [W/m/K]
354             virtual const volScalarField& kappa() const;
357             // Transfer fields - to the primary region
359                 //- Return mass transfer source - Eulerian phase only
360                 virtual tmp<volScalarField> primaryMassTrans() const;
362                 //- Return the film mass available for transfer to cloud
363                 virtual const volScalarField& cloudMassTrans() const;
365                 //- Return the parcel diameters originating from film to cloud
366                 virtual const volScalarField& cloudDiameterTrans() const;
369         // External helper functions
371             //- External hook to add sources to the film
372             virtual void addSources
373             (
374                 const label patchI,            // patchI on primary region
375                 const label faceI,             // faceI of patchI
376                 const scalar massSource,       // [kg]
377                 const vector& momentumSource,  // [kg.m/s] (tang'l momentum)
378                 const scalar pressureSource,   // [kg.m/s] (normal momentum)
379                 const scalar energySource = 0  // [J]
380             );
383          // Source fields (read/write access)
385             // Primary region
387                 //- Momementum / [kg/m/s2]
388                 inline volVectorField& USpPrimary();
390                 //- Pressure / [Pa]
391                 inline volScalarField& pSpPrimary();
393                 //- Mass / [kg/m2/s]
394                 inline volScalarField& rhoSpPrimary();
397             // Film region
399                 //- Momentum / [kg/m/s2]
400                 inline volVectorField& USp();
402                 //- Pressure / [Pa]
403                 inline volScalarField& pSp();
405                 //- Mass / [kg/m2/s]
406                 inline volScalarField& rhoSp();
408                 //- Momentum / [kg/m/s2]
409                 inline const volVectorField& USp() const;
411                 //- Pressure / [Pa]
412                 inline const volScalarField& pSp() const;
414                 //- Mass / [kg/m2/s]
415                 inline const volScalarField& rhoSp() const;
418         // Fields mapped from primary region
420             //- Velocity / [m/s]
421             inline const volVectorField& UPrimary() const;
423             //- Pressure / [Pa]
424             inline const volScalarField& pPrimary() const;
426             //- Density / [kg/m3]
427             inline const volScalarField& rhoPrimary() const;
429             //- Viscosity / [Pa.s]
430             inline const volScalarField& muPrimary() const;
433         // Sub-models
435             //- Injection
436             inline injectionModelList& injection();
439         // Helper functions
441             //- Return the current film mass
442             inline tmp<volScalarField> mass() const;
444             //- Return the net film mass available over the next integration
445             inline tmp<volScalarField> netMass() const;
447             //- Return the gravity normal-to-patch component contribution
448             inline tmp<volScalarField> gNorm() const;
450             //- Return the gravity normal-to-patch component contribution
451             //  Clipped so that only non-zero if g & nHat_ < 0
452             inline tmp<volScalarField> gNormClipped() const;
454             //- Return the gravity tangential component contributions
455             inline tmp<volVectorField> gTan() const;
458         // Evolution
460             //- Pre-evolve film hook
461             virtual void preEvolveRegion();
463             //- Evolve the film equations
464             virtual void evolveRegion();
467         // Source fields
469             // Mapped into primary region
471                 //- Return total mass source - Eulerian phase only
472                 virtual tmp<DimensionedField<scalar, volMesh> > Srho() const;
474                 //- Return mass source for specie i - Eulerian phase only
475                 virtual tmp<DimensionedField<scalar, volMesh> > Srho
476                 (
477                     const label i
478                 ) const;
480                 //- Return enthalpy source - Eulerian phase only
481                 virtual tmp<DimensionedField<scalar, volMesh> > Sh() const;
484         // I-O
486             //- Provide some feedback
487             virtual void info() const;
491 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
493 } // End namespace surfaceFilmModels
494 } // End namespace regionModels
495 } // End namespace Foam
497 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
499 #ifdef NoRepository
500 #   include "kinematicSingleLayerTemplates.C"
501 #endif
503 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
505 #include "kinematicSingleLayerI.H"
507 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
509 #endif
511 // ************************************************************************* //