Loading...
Searching...
No Matches
KinematicCloud.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) 2016-2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::KinematicCloud
29
30Group
31 grpLagrangianIntermediateClouds
32
33Description
34 Templated base class for kinematic cloud
35
36 - cloud function objects
37
38 - particle forces, e.g.
39 - buoyancy
40 - drag
41 - pressure gradient
42 - ...
43
44 - sub-models:
45 - dispersion model
46 - injection model
47 - patch interaction model
48 - stochastic collision model
49 - surface film model
50
51SourceFiles
52 KinematicCloudI.H
53 KinematicCloud.C
54
55\*---------------------------------------------------------------------------*/
56
57#ifndef KinematicCloud_H
58#define KinematicCloud_H
59
60#include "particle.H"
61#include "Cloud.H"
62#include "kinematicCloud.H"
63#include "IOdictionary.H"
64#include "autoPtr.H"
65#include "Random.H"
66#include "fvMesh.H"
67#include "volFields.H"
68#include "fvMatrices.H"
69#include "cloudSolution.H"
70
71#include "ParticleForceList.H"
73
74// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
75
76namespace Foam
77{
78
79// Forward declaration of classes
80
82
83template<class CloudType>
85
86template<class CloudType>
87class DispersionModel;
88
89template<class CloudType>
91
92template<class CloudType>
94
95template<class CloudType>
97
98template<class CloudType>
99class PackingModel;
100
101template<class CloudType>
102class DampingModel;
103
104template<class CloudType>
105class IsotropyModel;
106
107
108/*---------------------------------------------------------------------------*\
109 Class KinematicCloud Declaration
110\*---------------------------------------------------------------------------*/
111
112template<class CloudType>
113class KinematicCloud
114:
115 public CloudType,
116 public kinematicCloud
117{
118public:
119
120 // Public typedefs
121
122 //- Type of cloud this cloud was instantiated for
123 typedef CloudType cloudType;
124
125 //- Type of parcel the cloud was instantiated for
126 typedef typename CloudType::particleType parcelType;
127
128 //- Convenience typedef for this cloud type
129 typedef KinematicCloud<CloudType> kinematicCloudType;
131 //- Force models type
133
134 //- Function object type
137
138
139private:
140
141 // Private Data
142
143 //- Cloud copy pointer
145
146
147 // Private Member Functions
148
149 //- No copy construct
150 KinematicCloud(const KinematicCloud&) = delete;
151
152 //- No copy assignment
153 void operator=(const KinematicCloud&) = delete;
154
155
156protected:
157
158 // Protected Data
159
160 //- References to the mesh and time databases
161 const fvMesh& mesh_;
162
163 //- Dictionary of particle properties
165
166 //- Dictionary of output properties
168
169 //- Solution properties
171
172 //- Parcel constant properties
173 typename parcelType::constantProperties constProps_;
175 //- Sub-models dictionary
177
178 //- Random number generator - used by some injection routines
180
181 //- Cell occupancy information for each parcel, (demand driven)
183
184 //- Cell length scale
186
187
188 // References to the carrier gas fields
190 //- Density [kg/m3]
191 const volScalarField& rho_;
192
193 //- Velocity [m/s]
195
196 //- Dynamic viscosity [Pa.s]
197 const volScalarField& mu_;
198
200 // Environmental properties
201
202 //- Gravity
203 const dimensionedVector& g_;
205 //- Averaged ambient domain pressure
206 scalar pAmbient_;
207
208
209 //- Optional particle forces
211
212 //- Optional cloud function objects
215
216 // References to the cloud sub-models
217
218 //- Injector models
220
221 //- Dispersion model
224
225 //- Patch interaction model
228
229 //- Stochastic collision model
233 //- Surface film model
236
237 //- Packing model
241 //- Damping model
244
245 //- Exchange model
248
249
250 // Reference to the particle integration schemes
252 //- Velocity integration
254
255
256 // Sources
257
258 //- Mass for kinematic cloud
260
261 //- Momentum
263
264 //- Coefficient for carrier phase U equation
266
267
268 // Initialisation
269
270 //- Set cloud sub-models
271 void setModels();
272
273
274 // Cloud evolution functions
275
276 //- Solve the cloud - calls all evolution functions
277 template<class TrackCloudType>
278 void solve
279 (
280 TrackCloudType& cloud,
281 typename parcelType::trackingData& td
282 );
283
284 //- Build the cellOccupancy
285 void buildCellOccupancy();
286
287 //- Update (i.e. build) the cellOccupancy if it has
288 // already been used
289 void updateCellOccupancy();
290
291 //- Evolve the cloud
292 template<class TrackCloudType>
293 void evolveCloud
295 TrackCloudType& cloud,
296 typename parcelType::trackingData& td
297 );
298
299 //- Post-evolve
300 void postEvolve(const typename parcelType::trackingData& td);
301
302 //- Reset state of cloud
303 void cloudReset(KinematicCloud<CloudType>& c);
304
305
306public:
307
308 // Public Data
309
310 //- Flag to write log into Info
311 bool log;
312
313
314 // Constructors
315
316 //- Construct given carrier gas fields
317 KinematicCloud
318 (
319 const word& cloudName,
320 const volScalarField& rho,
321 const volVectorField& U,
323 const dimensionedVector& g,
324 bool readFields = true
325 );
326
327 //- Copy constructor with new name
328 KinematicCloud
329 (
330 KinematicCloud<CloudType>& c,
331 const word& name
332 );
333
334 //- Copy constructor with new name - creates bare cloud
335 KinematicCloud
336 (
337 const fvMesh& mesh,
338 const word& name,
339 const KinematicCloud<CloudType>& c
340 );
341
342 //- Construct and return clone based on (this) with new name
344 {
346 (
347 new KinematicCloud(*this, name)
348 );
350
351 //- Construct and return bare clone based on (this) with new name
352 virtual autoPtr<Cloud<parcelType>> cloneBare(const word& name) const
353 {
355 (
356 new KinematicCloud(this->mesh(), name, *this)
357 );
359
360
361 //- Destructor
362 virtual ~KinematicCloud() = default;
363
364
365 // Member Functions
366
367 // Access
368
369 //- Return a reference to the cloud copy
370 inline const KinematicCloud& cloudCopy() const;
372
373 // References to the mesh and databases
374
375 //- Return reference to the mesh
376 inline const fvMesh& mesh() const;
377
378 //- Return particle properties dictionary
379 inline const IOdictionary& particleProperties() const;
381 //- Return output properties dictionary
382 inline const IOdictionary& outputProperties() const;
383
384 //- Return non-const access to the output properties dictionary
386
387 //- Return const access to the solution properties
388 inline const cloudSolution& solution() const;
389
390 //- Return access to the solution properties
391 inline cloudSolution& solution();
392
393 //- Return the constant properties
394 inline const typename parcelType::constantProperties&
395 constProps() const;
396
397 //- Return access to the constant properties
398 inline typename parcelType::constantProperties& constProps();
399
400 //- Return reference to the sub-models dictionary
401 inline const dictionary& subModelProperties() const;
402
404 // Cloud data
405
406 //- Return reference to the random object
407 inline Random& rndGen() const;
408
409 //- Return the cell occupancy information for each
410 // parcel, non-const access, the caller is
411 // responsible for updating it for its own purposes
412 // if particles are removed or created.
414
415 //- Return the cell length scale
416 inline const scalarField& cellLengthScale() const;
417
418
419 // References to the carrier gas fields
420
421 //- Return carrier gas velocity
422 inline const volVectorField& U() const;
423
424 //- Return carrier gas density
425 inline const volScalarField& rho() const;
426
427 //- Return carrier gas dynamic viscosity
428 inline const volScalarField& mu() const;
429
430
431 // Environmental properties
432
433 //- Gravity
434 inline const dimensionedVector& g() const;
436 //- Return const-access to the ambient pressure
437 inline scalar pAmbient() const;
438
439 //- Return reference to the ambient pressure
440 inline scalar& pAmbient();
441
442
443 //- Optional particle forces
444 inline const forceType& forces() const;
445
446 //- Return the optional particle forces
447 inline forceType& forces();
448
449 //- Optional cloud function objects
450 inline functionType& functions();
451
452
453 // Sub-models
454
455 //- Return const access to the injection model
457 injectors() const;
459 //- Return reference to the injection model
461 injectors();
462
463 //- Return const-access to the dispersion model
465 dispersion() const;
466
467 //- Return reference to the dispersion model
469 dispersion();
470
471 //- Return const-access to the patch interaction model
473 patchInteraction() const;
474
475 //- Return reference to the patch interaction model
478
479 //- Return const-access to the stochastic collision model
480 inline const
482 stochasticCollision() const;
483
484 //- Return reference to the stochastic collision model
487
488 //- Return const-access to the surface film model
490 surfaceFilm() const;
492 //- Return reference to the surface film model
494 surfaceFilm();
495
497 //- Return const access to the packing model
499 packingModel() const;
500
501 //- Return a reference to the packing model
503 packingModel();
504
505 //- Return const access to the damping model
508
509 //- Return a reference to the damping model
511 dampingModel();
513 //- Return const access to the isotropy model
515 isotropyModel() const;
516
517 //- Return a reference to the isotropy model
520
521
522 // Integration schemes
523
524 //-Return reference to velocity integration
525 inline const integrationScheme& UIntegrator() const;
526
527
528 // Sources
529
530 //- Transfer the effect of parcel to the carrier phase
531 inline void transferToCarrier
532 (
533 const parcelType& p,
534 const typename parcelType::trackingData& td
535 );
536
537
538 // Momentum
540 //- Return reference to mass for kinematic source
542
543 //- Return const reference to mass for kinematic source
544 inline const volScalarField::Internal& rhokTrans() const;
545
546 //- Return reference to momentum source
548
549 //- Return const reference to momentum source
550 inline const volVectorField::Internal& UTrans() const;
551
552 //- Return coefficient for carrier phase U equation
554
555 //- Return const coefficient for carrier phase U equation
556 inline const volScalarField::Internal& UCoeff() const;
558 //- Return tmp mass source for kinematic
560
561 //- Return tmp momentum source term (compressible)
563 (
565 bool incompressible = false
566 ) const;
567
568
569 // Check
571 //- Total number of parcels
572 virtual label nParcels() const
573 {
574 return CloudType::nParcels();
576
577 //- Total mass in system
578 inline scalar massInSystem() const;
579
580 //- Total linear momentum of the system
582
583 //- Average particle per parcel
584 inline scalar totalParticlePerParcel() const;
585
586 //- Total linear kinetic energy in the system
587 inline scalar linearKineticEnergyOfSystem() const;
588
589 //- Total rotational kinetic energy in the system
590 inline scalar rotationalKineticEnergyOfSystem() const;
592 //- Mean diameter Dij
593 inline scalar Dij(const label i, const label j) const;
594
595 //- Max diameter
596 inline scalar Dmax() const;
597
598
599 // Fields
601 //- Volume swept rate of parcels per cell
602 inline const tmp<volScalarField> vDotSweep() const;
603
604 //- Return the particle volume fraction field
605 // Note: for particles belonging to this cloud only
606 inline const tmp<volScalarField> theta() const;
607
608 //- Return the particle mass fraction field
609 // Note: for particles belonging to this cloud only
610 inline const tmp<volScalarField> alpha() const;
611
612 //- Return the particle effective density field
613 // Note: for particles belonging to this cloud only
614 inline const tmp<volScalarField> rhoEff() const;
615
616
617 // Cloud evolution functions
619 //- Set parcel thermo properties
621 (
622 parcelType& parcel,
623 const scalar lagrangianDt
624 );
625
626 //- Check parcel properties
628 (
629 parcelType& parcel,
630 const scalar lagrangianDt,
631 const bool fullyDescribed
632 );
633
634 //- Store the current cloud state
635 void storeState();
636
637 //- Reset the current cloud to the previously stored state
638 void restoreState();
639
640 //- Reset the cloud source terms
641 void resetSourceTerms();
642
643 //- Relax field
644 template<class Type>
645 void relax
646 (
649 const word& name
650 ) const;
651
652 //- Scale field
653 template<class Type>
654 void scale
657 const word& name
658 ) const;
659
660 //- Apply relaxation to (steady state) cloud sources
661 void relaxSources(const KinematicCloud<CloudType>& cloudOldTime);
663 //- Apply scaling to (transient) cloud sources
664 void scaleSources();
665
666 //- Pre-evolve
667 void preEvolve
669 const typename parcelType::trackingData& td
670 );
671
672 //- Evolve the cloud
673 void evolve();
675 //- Particle motion
676 template<class TrackCloudType>
677 void motion
678 (
679 TrackCloudType& cloud,
680 typename parcelType::trackingData& td
681 );
682
683 //- Calculate the patch normal and velocity to interact with,
684 // accounting for patch motion if required.
685 void patchData
687 const parcelType& p,
688 const polyPatch& pp,
689 vector& normal,
690 vector& Up
691 ) const;
693
694 // Mapping
695
696 //- Update mesh
697 void updateMesh();
698
699 //- Remap the cells of particles corresponding to the
700 // mesh topology change with a default tracking data object
701 virtual void autoMap(const mapPolyMesh&);
702
703
704 // I-O
705
706 //- Print cloud information
707 void info();
709 //- Read particle fields from objects in the obr registry
710 virtual void readObjects(const objectRegistry& obr);
711
712 //- Write particle fields as objects into the obr registry
713 virtual void writeObjects(objectRegistry& obr) const;
714};
715
716
717// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
718
719} // End namespace Foam
721// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
722
723#include "KinematicCloudI.H"
724
725// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
726
727#ifdef NoRepository
728 #include "KinematicCloud.C"
729#endif
731// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
732
733#endif
734
735// ************************************************************************* //
const word cloudName(propsDict.get< word >("cloud"))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
List of cloud function objects.
virtual label nParcels() const
Return the number of particles in the cloud.
Definition Cloud.H:206
ParticleType particleType
Definition Cloud.H:130
Base class for collisional damping models.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Base class for dispersion modelling.
DimensionedField< scalar, volMesh > Internal
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
autoPtr< IOobject > clone() const
Clone.
Definition IOobject.H:641
List of injection models.
Base class for collisional return-to-isotropy models.
Templated base class for kinematic cloud.
DispersionModel< KinematicCloud< CloudType > > & dispersion()
Return reference to the dispersion model.
DampingModel< KinematicCloud< CloudType > > & dampingModel()
Return a reference to the damping model.
const parcelType::constantProperties & constProps() const
Return the constant properties.
scalar & pAmbient()
Return reference to the ambient pressure.
const dictionary & subModelProperties() const
Return reference to the sub-models dictionary.
SurfaceFilmModel< KinematicCloud< CloudType > > & surfaceFilm()
Return reference to the surface film model.
autoPtr< volScalarField::Internal > rhokTrans_
autoPtr< DampingModel< KinematicCloud< CloudType > > > dampingModel_
void setModels()
Set cloud sub-models.
scalar massInSystem() const
Total mass in system.
scalar rotationalKineticEnergyOfSystem() const
Total rotational kinetic energy in the system.
scalar totalParticlePerParcel() const
Average particle per parcel.
virtual label nParcels() const
Total number of parcels.
InjectionModelList< KinematicCloud< CloudType > > & injectors()
Return reference to the injection model.
void relaxSources(const KinematicCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
PatchInteractionModel< KinematicCloud< CloudType > > & patchInteraction()
Return reference to the patch interaction model.
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
void scale(DimensionedField< Type, volMesh > &field, const word &name) const
Scale field.
tmp< volScalarField::Internal > Srhok() const
Return tmp mass source for kinematic.
const tmp< volScalarField > theta() const
Return the particle volume fraction field.
void storeState()
Store the current cloud state.
autoPtr< SurfaceFilmModel< KinematicCloud< CloudType > > > surfaceFilmModel_
const tmp< volScalarField > alpha() const
Return the particle mass fraction field.
void patchData(const parcelType &p, const polyPatch &pp, vector &normal, vector &Up) const
Calculate the patch normal and velocity to interact with,.
const scalarField & cellLengthScale() const
Return the cell length scale.
void setParcelThermoProperties(parcelType &parcel, const scalar lagrangianDt)
Set parcel thermo properties.
const SurfaceFilmModel< KinematicCloud< CloudType > > & surfaceFilm() const
Return const-access to the surface film model.
virtual void readObjects(const objectRegistry &obr)
Read particle fields from objects in the obr registry.
IsotropyModel< KinematicCloud< CloudType > > & isotropyModel()
Return a reference to the isotropy model.
void cloudReset(KinematicCloud< CloudType > &c)
Reset state of cloud.
const PatchInteractionModel< KinematicCloud< CloudType > > & patchInteraction() const
Return const-access to the patch interaction model.
void checkParcelProperties(parcelType &parcel, const scalar lagrangianDt, const bool fullyDescribed)
Check parcel properties.
KinematicCloud(const word &cloudName, const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g, bool readFields=true)
Construct given carrier gas fields.
vector linearMomentumOfSystem() const
Total linear momentum of the system.
StochasticCollisionModel< KinematicCloud< CloudType > > & stochasticCollision()
Return reference to the stochastic collision model.
scalar pAmbient() const
Return const-access to the ambient pressure.
void relax(DimensionedField< Type, volMesh > &field, const DimensionedField< Type, volMesh > &field0, const word &name) const
Relax field.
const KinematicCloud & cloudCopy() const
Return a reference to the cloud copy.
InjectionModelList< KinematicCloud< CloudType > > injectors_
PackingModel< KinematicCloud< CloudType > > & packingModel()
Return a reference to the packing model.
autoPtr< IsotropyModel< KinematicCloud< CloudType > > > isotropyModel_
autoPtr< DispersionModel< KinematicCloud< CloudType > > > dispersionModel_
void scaleSources()
Apply scaling to (transient) cloud sources.
autoPtr< PatchInteractionModel< KinematicCloud< CloudType > > > patchInteractionModel_
const integrationScheme & UIntegrator() const
Return reference to velocity integration.
List< DynamicList< parcelType * > > & cellOccupancy()
Return the cell occupancy information for each.
functionType & functions()
Optional cloud function objects.
Cloud< basicKinematicParcel >::particleType parcelType
const tmp< volScalarField > rhoEff() const
Return the particle effective density field.
const cloudSolution & solution() const
Return const access to the solution properties.
void updateCellOccupancy()
Update (i.e. build) the cellOccupancy if it has.
const tmp< volScalarField > vDotSweep() const
Volume swept rate of parcels per cell.
const PackingModel< KinematicCloud< CloudType > > & packingModel() const
Return const access to the packing model.
virtual void autoMap(const mapPolyMesh &)
Remap the cells of particles corresponding to the.
virtual autoPtr< Cloud< parcelType > > clone(const word &name)
Construct and return clone based on (this) with new name.
const DispersionModel< KinematicCloud< CloudType > > & dispersion() const
Return const-access to the dispersion model.
virtual autoPtr< Cloud< parcelType > > cloneBare(const word &name) const
Construct and return bare clone based on (this) with new name.
void evolve()
Evolve the cloud.
forceType & forces()
Return the optional particle forces.
void postEvolve(const typename parcelType::trackingData &td)
Post-evolve.
volScalarField::Internal & UCoeff()
Return coefficient for carrier phase U equation.
const StochasticCollisionModel< KinematicCloud< CloudType > > & stochasticCollision() const
Return const-access to the stochastic collision model.
void evolveCloud(TrackCloudType &cloud, typename parcelType::trackingData &td)
Evolve the cloud.
scalar Dmax() const
Max diameter.
autoPtr< List< DynamicList< parcelType * > > > cellOccupancyPtr_
const DampingModel< KinematicCloud< CloudType > > & dampingModel() const
Return const access to the damping model.
CloudFunctionObjectList< KinematicCloud< CloudType > > functionType
scalar Dij(const label i, const label j) const
Mean diameter Dij.
parcelType::constantProperties & constProps()
Return access to the constant properties.
const IOdictionary & particleProperties() const
Return particle properties dictionary.
volVectorField::Internal & UTrans()
Return reference to momentum source.
void preEvolve(const typename parcelType::trackingData &td)
Pre-evolve.
void info()
Print cloud information.
void restoreState()
Reset the current cloud to the previously stored state.
tmp< fvVectorMatrix > SU(volVectorField &U, bool incompressible=false) const
Return tmp momentum source term (compressible).
virtual ~KinematicCloud()=default
Destructor.
void transferToCarrier(const parcelType &p, const typename parcelType::trackingData &td)
Transfer the effect of parcel to the carrier phase.
const volScalarField::Internal & rhokTrans() const
Return const reference to mass for kinematic source.
const InjectionModelList< KinematicCloud< CloudType > > & injectors() const
Return const access to the injection model.
Random & rndGen() const
Return reference to the random object.
void resetSourceTerms()
Reset the cloud source terms.
volScalarField::Internal & rhokTrans()
Return reference to mass for kinematic source.
const forceType & forces() const
Optional particle forces.
ParticleForceList< KinematicCloud< CloudType > > forceType
KinematicCloud(const fvMesh &mesh, const word &name, const KinematicCloud< CloudType > &c)
Copy constructor with new name - creates bare cloud.
autoPtr< PackingModel< KinematicCloud< CloudType > > > packingModel_
const IOdictionary & outputProperties() const
Return output properties dictionary.
const IsotropyModel< KinematicCloud< CloudType > > & isotropyModel() const
Return const access to the isotropy model.
KinematicCloud(KinematicCloud< CloudType > &c, const word &name)
Copy constructor with new name.
void buildCellOccupancy()
Build the cellOccupancy.
IOdictionary & outputProperties()
Return non-const access to the output properties dictionary.
void updateMesh()
Update mesh.
const volScalarField::Internal & UCoeff() const
Return const coefficient for carrier phase U equation.
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
autoPtr< StochasticCollisionModel< KinematicCloud< CloudType > > > stochasticCollisionModel_
const volVectorField::Internal & UTrans() const
Return const reference to momentum source.
cloudSolution & solution()
Return access to the solution properties.
void solve(TrackCloudType &cloud, typename parcelType::trackingData &td)
Solve the cloud - calls all evolution functions.
virtual void writeObjects(objectRegistry &obr) const
Write particle fields as objects into the obr registry.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
Base class for packing models.
List of particle forces.
Templated patch interaction model class.
Random number generator.
Definition Random.H:56
Templated stochastic collision model class.
Templated wall surface film model class.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Stores all relevant solution info for cloud.
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Base for a set of schemes which integrate simple ODEs which arise from semi-implcit rate expressions.
kinematicCloud()=default
Null constructor.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Registry of regIOobjects.
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
UEqn relax()
rDeltaTY field()
A special matrix type and solver, designed for finite volume solutions of scalar equations.
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for OpenFOAM.
DSMCCloud< dsmcParcel > CloudType
GeometricField< vector, fvPatchField, volMesh > volVectorField
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.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Vector< scalar > vector
Definition vector.H:57
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
CEqn solve()