1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
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
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/>.
28 Templated base class for kinematic cloud
30 - cloud function objects
32 - particle forces, e.g.
41 - patch interaction model
48 \*---------------------------------------------------------------------------*/
50 #ifndef KinematicCloud_H
51 #define KinematicCloud_H
55 #include "kinematicCloud.H"
56 #include "IOdictionary.H"
58 #include "cachedRandom.H"
60 #include "volFields.H"
61 #include "fvMatrices.H"
62 #include "IntegrationSchemesFwd.H"
63 #include "cloudSolution.H"
65 #include "ParticleForceList.H"
66 #include "CloudFunctionObjectList.H"
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 // Forward declaration of classes
75 template<class CloudType>
76 class DispersionModel;
78 template<class CloudType>
81 template<class CloudType>
82 class PatchInteractionModel;
84 template<class CloudType>
85 class SurfaceFilmModel;
88 /*---------------------------------------------------------------------------*\
89 Class KinematicCloud Declaration
90 \*---------------------------------------------------------------------------*/
92 template<class CloudType>
102 //- Type of cloud this cloud was instantiated for
103 typedef CloudType cloudType;
105 //- Type of parcel the cloud was instantiated for
106 typedef typename CloudType::particleType parcelType;
108 //- Convenience typedef for this cloud type
109 typedef KinematicCloud<CloudType> kinematicCloudType;
111 //- Force models type
112 typedef ParticleForceList<KinematicCloud<CloudType> > forceType;
114 //- Function object type
115 typedef CloudFunctionObjectList<KinematicCloud<CloudType> >
123 //- Cloud copy pointer
124 autoPtr<KinematicCloud<CloudType> > cloudCopyPtr_;
127 // Private Member Functions
129 //- Disallow default bitwise copy construct
130 KinematicCloud(const KinematicCloud&);
132 //- Disallow default bitwise assignment
133 void operator=(const KinematicCloud&);
140 //- References to the mesh and time databases
143 //- Dictionary of particle properties
144 IOdictionary particleProperties_;
146 //- Solution properties
147 cloudSolution solution_;
149 //- Parcel constant properties
150 typename parcelType::constantProperties constProps_;
152 //- Sub-models dictionary
153 const dictionary subModelProperties_;
155 //- Random number generator - used by some injection routines
156 cachedRandom rndGen_;
158 //- Cell occupancy information for each parcel, (demand driven)
159 autoPtr<List<DynamicList<parcelType*> > > cellOccupancyPtr_;
162 // References to the carrier gas fields
165 const volScalarField& rho_;
168 const volVectorField& U_;
170 //- Dynamic viscosity [Pa.s]
171 const volScalarField& mu_;
174 // Environmental properties
177 const dimensionedVector& g_;
179 //- Averaged ambient domain pressure
183 //- Optional particle forces
186 //- Optional cloud function objects
187 functionType functions_;
190 // References to the cloud sub-models
193 autoPtr<DispersionModel<KinematicCloud<CloudType> > >
197 autoPtr<InjectionModel<KinematicCloud<CloudType> > >
200 //- Patch interaction model
201 autoPtr<PatchInteractionModel<KinematicCloud<CloudType> > >
202 patchInteractionModel_;
204 //- Surface film model
205 autoPtr<SurfaceFilmModel<KinematicCloud<CloudType> > >
209 // Reference to the particle integration schemes
211 //- Velocity integration
212 autoPtr<vectorIntegrationScheme> UIntegrator_;
218 autoPtr<DimensionedField<vector, volMesh> > UTrans_;
220 //- Coefficient for carrier phase U equation
221 autoPtr<DimensionedField<scalar, volMesh> > UCoeff_;
226 //- Set cloud sub-models
230 // Cloud evolution functions
232 //- Solve the cloud - calls all evolution functions
233 template<class TrackData>
234 void solve(TrackData& td);
236 //- Build the cellOccupancy
237 void buildCellOccupancy();
239 //- Update (i.e. build) the cellOccupancy if it has
241 void updateCellOccupancy();
244 template<class TrackData>
245 void evolveCloud(TrackData& td);
250 //- Reset state of cloud
251 void cloudReset(KinematicCloud<CloudType>& c);
258 //- Construct given carrier gas fields
261 const word& cloudName,
262 const volScalarField& rho,
263 const volVectorField& U,
264 const volScalarField& mu,
265 const dimensionedVector& g,
266 bool readFields = true
269 //- Copy constructor with new name
272 KinematicCloud<CloudType>& c,
276 //- Copy constructor with new name - creates bare cloud
281 const KinematicCloud<CloudType>& c
284 //- Construct and return clone based on (this) with new name
285 virtual autoPtr<Cloud<parcelType> > clone(const word& name)
287 return autoPtr<Cloud<parcelType> >
289 new KinematicCloud(*this, name)
293 //- Construct and return bare clone based on (this) with new name
294 virtual autoPtr<Cloud<parcelType> > cloneBare(const word& name) const
296 return autoPtr<Cloud<parcelType> >
298 new KinematicCloud(this->mesh(), name, *this)
304 virtual ~KinematicCloud();
311 //- Return a reference to the cloud copy
312 inline const KinematicCloud& cloudCopy() const;
314 //- Switch to specify if particles of the cloud can return
315 // non-zero wall distance values - true for kinematic parcels
316 virtual bool hasWallImpactDistance() const;
319 // References to the mesh and databases
321 //- Return reference to the mesh
322 inline const fvMesh& mesh() const;
324 //- Return particle properties dictionary
325 inline const IOdictionary& particleProperties() const;
327 //- Return const access to the solution properties
328 inline const cloudSolution& solution() const;
330 //- Return access to the solution properties
331 inline cloudSolution& solution();
333 //- Return the constant properties
334 inline const typename parcelType::constantProperties&
337 //- Return reference to the sub-models dictionary
338 inline const dictionary& subModelProperties() const;
343 //- Return reference to the random object
344 inline cachedRandom& rndGen();
346 //- Return the cell occupancy information for each
347 // parcel, non-const access, the caller is
348 // responsible for updating it for its own purposes
349 // if particles are removed or created.
350 inline List<DynamicList<parcelType*> >& cellOccupancy();
353 // References to the carrier gas fields
355 //- Return carrier gas velocity
356 inline const volVectorField& U() const;
358 //- Return carrier gas density
359 inline const volScalarField& rho() const;
361 //- Return carrier gas dynamic viscosity
362 inline const volScalarField& mu() const;
365 // Environmental properties
368 inline const dimensionedVector& g() const;
370 //- Return const-access to the ambient pressure
371 inline scalar pAmbient() const;
373 //- Return reference to the ambient pressure
374 inline scalar& pAmbient();
377 //- Optional particle forces
378 // inline const typename parcelType::forceType& forces() const;
379 inline const forceType& forces() const;
381 //- Optional cloud function objects
382 inline functionType& functions();
387 //- Return const-access to the dispersion model
388 inline const DispersionModel<KinematicCloud<CloudType> >&
391 //- Return reference to the dispersion model
392 inline DispersionModel<KinematicCloud<CloudType> >&
395 //- Return const access to the injection model
396 inline const InjectionModel<KinematicCloud<CloudType> >&
399 //- Return reference to the injection model
400 inline InjectionModel<KinematicCloud<CloudType> >&
403 //- Return const-access to the patch interaction model
404 inline const PatchInteractionModel<KinematicCloud<CloudType> >&
405 patchInteraction() const;
407 //- Return reference to the patch interaction model
408 inline PatchInteractionModel<KinematicCloud<CloudType> >&
411 //- Return const-access to the surface film model
412 inline const SurfaceFilmModel<KinematicCloud<CloudType> >&
415 //- Return reference to the surface film model
416 inline SurfaceFilmModel<KinematicCloud<CloudType> >&
420 // Integration schemes
422 //-Return reference to velocity integration
423 inline const vectorIntegrationScheme& UIntegrator() const;
430 //- Return reference to momentum source
431 inline DimensionedField<vector, volMesh>& UTrans();
433 //- Return const reference to momentum source
434 inline const DimensionedField<vector, volMesh>&
437 //- Return coefficient for carrier phase U equation
438 inline DimensionedField<scalar, volMesh>& UCoeff();
440 //- Return const coefficient for carrier phase U equation
441 inline const DimensionedField<scalar, volMesh>&
444 //- Return tmp momentum source term
445 inline tmp<fvVectorMatrix> SU(volVectorField& U) const;
450 //- Total mass injected
451 inline scalar massInjected() const;
453 //- Total mass in system
454 inline scalar massInSystem() const;
456 //- Total linear momentum of the system
457 inline vector linearMomentumOfSystem() const;
459 //- Total linear kinetic energy in the system
460 inline scalar linearKineticEnergyOfSystem() const;
462 //- Total rotational kinetic energy in the system
463 inline scalar rotationalKineticEnergyOfSystem() const;
465 //- Penetration for percentage of the current total mass
466 inline scalar penetration(const scalar& prc) const;
468 //- Mean diameter Dij
469 inline scalar Dij(const label i, const label j) const;
474 //- Return the particle volume fraction field
475 // Note: for particles belonging to this cloud only
476 inline const tmp<volScalarField> theta() const;
478 //- Return the particle mass fraction field
479 // Note: for particles belonging to this cloud only
480 inline const tmp<volScalarField> alpha() const;
482 //- Return the particle effective density field
483 // Note: for particles belonging to this cloud only
484 inline const tmp<volScalarField> rhoEff() const;
487 // Cloud evolution functions
489 //- Set parcel thermo properties
490 void setParcelThermoProperties
493 const scalar lagrangianDt
496 //- Check parcel properties
497 void checkParcelProperties
500 const scalar lagrangianDt,
501 const bool fullyDescribed
504 //- Store the current cloud state
507 //- Reset the current cloud to the previously stored state
510 //- Reset the cloud source terms
511 void resetSourceTerms();
517 DimensionedField<Type, volMesh>& field,
518 const DimensionedField<Type, volMesh>& field0,
526 DimensionedField<Type, volMesh>& field,
530 //- Apply relaxation to (steady state) cloud sources
531 void relaxSources(const KinematicCloud<CloudType>& cloudOldTime);
533 //- Apply scaling to (transient) cloud sources
543 template<class TrackData>
544 void motion(TrackData& td);
546 //- Print cloud information
551 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
553 } // End namespace Foam
555 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
557 #include "KinematicCloudI.H"
559 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
562 # include "KinematicCloud.C"
565 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
569 // ************************************************************************* //