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-2020 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
81class integrationScheme;
82
83template<class CloudType>
84class InjectionModelList;
85
86template<class CloudType>
87class DispersionModel;
88
89template<class CloudType>
90class PatchInteractionModel;
91
92template<class CloudType>
93class SurfaceFilmModel;
94
95template<class CloudType>
96class StochasticCollisionModel;
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
130
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
174
175 //- Sub-models dictionary
177
178 //- Random number generator - used by some injection routines
179 mutable Random rndGen_;
180
181 //- Cell occupancy information for each parcel, (demand driven)
183
184 //- Cell length scale
186
187
188 // References to the carrier gas fields
189
190 //- Density [kg/m3]
191 const volScalarField& rho_;
192
193 //- Velocity [m/s]
194 const volVectorField& U_;
195
196 //- Dynamic viscosity [Pa.s]
197 const volScalarField& mu_;
198
199
200 // Environmental properties
201
202 //- Gravity
203 const dimensionedVector& g_;
204
205 //- Averaged ambient domain pressure
206 scalar pAmbient_;
207
208
209 //- Optional particle forces
211
212 //- Optional cloud function objects
214
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
232
233 //- Surface film model
236
237 //- Packing model
240
241 //- Damping model
244
245 //- Exchange model
248
249
250 // Reference to the particle integration schemes
251
252 //- Velocity integration
254
255
256 // Sources
257
258 //- Momentum
260
261 //- Coefficient for carrier phase U equation
263
264
265 // Initialisation
266
267 //- Set cloud sub-models
268 void setModels();
269
270
271 // Cloud evolution functions
272
273 //- Solve the cloud - calls all evolution functions
274 template<class TrackCloudType>
275 void solve
276 (
277 TrackCloudType& cloud,
278 typename parcelType::trackingData& td
279 );
280
281 //- Build the cellOccupancy
282 void buildCellOccupancy();
283
284 //- Update (i.e. build) the cellOccupancy if it has
285 // already been used
286 void updateCellOccupancy();
287
288 //- Evolve the cloud
289 template<class TrackCloudType>
290 void evolveCloud
291 (
292 TrackCloudType& cloud,
293 typename parcelType::trackingData& td
294 );
295
296 //- Post-evolve
297 void postEvolve(const typename parcelType::trackingData& td);
298
299 //- Reset state of cloud
301
302
303public:
304
305 // Constructors
306
307 //- Construct given carrier gas fields
309 (
310 const word& cloudName,
311 const volScalarField& rho,
312 const volVectorField& U,
313 const volScalarField& mu,
314 const dimensionedVector& g,
315 bool readFields = true
316 );
317
318 //- Copy constructor with new name
320 (
322 const word& name
323 );
324
325 //- Copy constructor with new name - creates bare cloud
327 (
328 const fvMesh& mesh,
329 const word& name,
331 );
332
333 //- Construct and return clone based on (this) with new name
335 {
337 (
338 new KinematicCloud(*this, name)
339 );
340 }
341
342 //- Construct and return bare clone based on (this) with new name
343 virtual autoPtr<Cloud<parcelType>> cloneBare(const word& name) const
344 {
346 (
347 new KinematicCloud(this->mesh(), name, *this)
348 );
349 }
350
351
352 //- Destructor
353 virtual ~KinematicCloud() = default;
354
355
356 // Member Functions
357
358 // Access
359
360 //- Return a reference to the cloud copy
361 inline const KinematicCloud& cloudCopy() const;
362
363
364 // References to the mesh and databases
365
366 //- Return reference to the mesh
367 inline const fvMesh& mesh() const;
368
369 //- Return particle properties dictionary
370 inline const IOdictionary& particleProperties() const;
371
372 //- Return output properties dictionary
373 inline const IOdictionary& outputProperties() const;
374
375 //- Return non-const access to the output properties dictionary
377
378 //- Return const access to the solution properties
379 inline const cloudSolution& solution() const;
380
381 //- Return access to the solution properties
382 inline cloudSolution& solution();
383
384 //- Return the constant properties
385 inline const typename parcelType::constantProperties&
386 constProps() const;
387
388 //- Return access to the constant properties
390
391 //- Return reference to the sub-models dictionary
392 inline const dictionary& subModelProperties() const;
393
394
395 // Cloud data
396
397 //- Return reference to the random object
398 inline Random& rndGen() const;
399
400 //- Return the cell occupancy information for each
401 // parcel, non-const access, the caller is
402 // responsible for updating it for its own purposes
403 // if particles are removed or created.
405
406 //- Return the cell length scale
407 inline const scalarField& cellLengthScale() const;
408
409
410 // References to the carrier gas fields
411
412 //- Return carrier gas velocity
413 inline const volVectorField& U() const;
414
415 //- Return carrier gas density
416 inline const volScalarField& rho() const;
417
418 //- Return carrier gas dynamic viscosity
419 inline const volScalarField& mu() const;
420
421
422 // Environmental properties
423
424 //- Gravity
425 inline const dimensionedVector& g() const;
426
427 //- Return const-access to the ambient pressure
428 inline scalar pAmbient() const;
429
430 //- Return reference to the ambient pressure
431 inline scalar& pAmbient();
432
433
434 //- Optional particle forces
435 inline const forceType& forces() const;
436
437 //- Return the optional particle forces
438 inline forceType& forces();
439
440 //- Optional cloud function objects
441 inline functionType& functions();
442
443
444 // Sub-models
445
446 //- Return const access to the injection model
448 injectors() const;
449
450 //- Return reference to the injection model
452 injectors();
453
454 //- Return const-access to the dispersion model
456 dispersion() const;
457
458 //- Return reference to the dispersion model
460 dispersion();
461
462 //- Return const-access to the patch interaction model
464 patchInteraction() const;
465
466 //- Return reference to the patch interaction model
469
470 //- Return const-access to the stochastic collision model
471 inline const
473 stochasticCollision() const;
474
475 //- Return reference to the stochastic collision model
478
479 //- Return const-access to the surface film model
481 surfaceFilm() const;
482
483 //- Return reference to the surface film model
485 surfaceFilm();
486
487
488 //- Return const access to the packing model
490 packingModel() const;
491
492 //- Return a reference to the packing model
494 packingModel();
495
496 //- Return const access to the damping model
498 dampingModel() const;
499
500 //- Return a reference to the damping model
502 dampingModel();
503
504 //- Return const access to the isotropy model
506 isotropyModel() const;
507
508 //- Return a reference to the isotropy model
511
512
513 // Integration schemes
514
515 //-Return reference to velocity integration
516 inline const integrationScheme& UIntegrator() const;
517
518
519 // Sources
520
521 // Momentum
522
523 //- Return reference to momentum source
525
526 //- Return const reference to momentum source
527 inline const volVectorField::Internal&
528 UTrans() const;
529
530 //- Return coefficient for carrier phase U equation
532
533 //- Return const coefficient for carrier phase U equation
534 inline const volScalarField::Internal&
535 UCoeff() const;
536
537 //- Return tmp momentum source term (compressible)
539 (
541 bool incompressible = false
542 ) const;
543
544
545 // Check
546
547 //- Total number of parcels
548 virtual label nParcels() const
549 {
550 return CloudType::nParcels();
551 }
552
553 //- Total mass in system
554 inline scalar massInSystem() const;
555
556 //- Total linear momentum of the system
557 inline vector linearMomentumOfSystem() const;
558
559 //- Average particle per parcel
560 inline scalar totalParticlePerParcel() const;
561
562 //- Total linear kinetic energy in the system
563 inline scalar linearKineticEnergyOfSystem() const;
564
565 //- Total rotational kinetic energy in the system
566 inline scalar rotationalKineticEnergyOfSystem() const;
567
568 //- Mean diameter Dij
569 inline scalar Dij(const label i, const label j) const;
570
571 //- Max diameter
572 inline scalar Dmax() const;
573
574
575 // Fields
576
577 //- Volume swept rate of parcels per cell
578 inline const tmp<volScalarField> vDotSweep() const;
579
580 //- Return the particle volume fraction field
581 // Note: for particles belonging to this cloud only
582 inline const tmp<volScalarField> theta() const;
583
584 //- Return the particle mass fraction field
585 // Note: for particles belonging to this cloud only
586 inline const tmp<volScalarField> alpha() const;
587
588 //- Return the particle effective density field
589 // Note: for particles belonging to this cloud only
590 inline const tmp<volScalarField> rhoEff() const;
591
592
593 // Cloud evolution functions
594
595 //- Set parcel thermo properties
597 (
598 parcelType& parcel,
599 const scalar lagrangianDt
600 );
601
602 //- Check parcel properties
604 (
605 parcelType& parcel,
606 const scalar lagrangianDt,
607 const bool fullyDescribed
608 );
609
610 //- Store the current cloud state
611 void storeState();
612
613 //- Reset the current cloud to the previously stored state
614 void restoreState();
615
616 //- Reset the cloud source terms
617 void resetSourceTerms();
618
619 //- Relax field
620 template<class Type>
621 void relax
622 (
625 const word& name
626 ) const;
627
628 //- Scale field
629 template<class Type>
630 void scale
631 (
633 const word& name
634 ) const;
635
636 //- Apply relaxation to (steady state) cloud sources
637 void relaxSources(const KinematicCloud<CloudType>& cloudOldTime);
638
639 //- Apply scaling to (transient) cloud sources
640 void scaleSources();
641
642 //- Pre-evolve
643 void preEvolve
644 (
645 const typename parcelType::trackingData& td
646 );
647
648 //- Evolve the cloud
649 void evolve();
650
651 //- Particle motion
652 template<class TrackCloudType>
653 void motion
654 (
655 TrackCloudType& cloud,
656 typename parcelType::trackingData& td
657 );
658
659 //- Calculate the patch normal and velocity to interact with,
660 // accounting for patch motion if required.
661 void patchData
662 (
663 const parcelType& p,
664 const polyPatch& pp,
665 vector& normal,
666 vector& Up
667 ) const;
668
669
670 // Mapping
671
672 //- Update mesh
673 void updateMesh();
674
675 //- Remap the cells of particles corresponding to the
676 // mesh topology change with a default tracking data object
677 virtual void autoMap(const mapPolyMesh&);
678
679
680 // I-O
681
682 //- Print cloud information
683 void info();
684
685 //- Read particle fields from objects in the obr registry
686 virtual void readObjects(const objectRegistry& obr);
687
688 //- Write particle fields as objects into the obr registry
689 virtual void writeObjects(objectRegistry& obr) const;
690};
691
692
693// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
694
695} // End namespace Foam
696
697// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
698
699#include "KinematicCloudI.H"
700
701// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
702
703#ifdef NoRepository
704 #include "KinematicCloud.C"
705#endif
706
707// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
708
709#endif
710
711// ************************************************************************* //
List of cloud function objects.
virtual label nParcels() const
Return the number of particles in the cloud.
Definition: Cloud.H:164
ParticleType particleType
Definition: Cloud.H:114
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
const word & cloudName() const
Return the cloud type.
Definition: DSMCCloudI.H:37
Class to hold DSMC particle constant properties.
Definition: DSMCParcel.H:82
Base class for collisional damping models.
Definition: DampingModel.H:66
Base class for dispersion modelling.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:459
List of injection models.
Base class for collisional return-to-isotropy models.
Definition: IsotropyModel.H:66
Templated base class for kinematic cloud.
const parcelType::constantProperties & constProps() const
Return the constant properties.
const dictionary & subModelProperties() const
Return reference to the sub-models dictionary.
cloudSolution solution_
Solution properties.
autoPtr< DampingModel< KinematicCloud< CloudType > > > dampingModel_
Damping model.
autoPtr< volVectorField::Internal > UTrans_
Momentum.
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.
void relaxSources(const KinematicCloud< CloudType > &cloudOldTime)
Apply relaxation to (steady state) cloud sources.
void motion(TrackCloudType &cloud, typename parcelType::trackingData &td)
Particle motion.
forceType forces_
Optional particle forces.
void scale(DimensionedField< Type, volMesh > &field, const word &name) const
Scale field.
const tmp< volScalarField > theta() const
Return the particle volume fraction field.
const fvMesh & mesh_
References to the mesh and time databases.
void storeState()
Store the current cloud state.
autoPtr< SurfaceFilmModel< KinematicCloud< CloudType > > > surfaceFilmModel_
Surface film model.
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.
scalar pAmbient_
Averaged ambient domain pressure.
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.
scalarField cellLengthScale_
Cell length scale.
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.
autoPtr< integrationScheme > UIntegrator_
Velocity integration.
vector linearMomentumOfSystem() const
Total linear momentum of the system.
const volVectorField & U() const
Return carrier gas velocity.
scalar pAmbient() const
Return const-access to the ambient pressure.
const KinematicCloud & cloudCopy() const
Return a reference to the cloud copy.
InjectionModelList< KinematicCloud< CloudType > > injectors_
Injector models.
const volScalarField & rho() const
Return carrier gas density.
autoPtr< IsotropyModel< KinematicCloud< CloudType > > > isotropyModel_
Exchange model.
autoPtr< DispersionModel< KinematicCloud< CloudType > > > dispersionModel_
Dispersion model.
void scaleSources()
Apply scaling to (transient) cloud sources.
KinematicCloud< CloudType > kinematicCloudType
Convenience typedef for this cloud type.
autoPtr< PatchInteractionModel< KinematicCloud< CloudType > > > patchInteractionModel_
Patch interaction model.
const volVectorField & U_
Velocity [m/s].
const integrationScheme & UIntegrator() const
Return reference to velocity integration.
List< DynamicList< parcelType * > > & cellOccupancy()
Return the cell occupancy information for each.
const volScalarField & mu_
Dynamic viscosity [Pa.s].
const dimensionedVector & g_
Gravity.
functionType & functions()
Optional cloud function objects.
CloudType::particleType parcelType
Type of parcel the cloud was instantiated for.
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 dictionary subModelProperties_
Sub-models dictionary.
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.
CloudType cloudType
Type of cloud this cloud was instantiated for.
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.
void postEvolve(const typename parcelType::trackingData &td)
Post-evolve.
volScalarField::Internal & UCoeff()
Return coefficient for carrier phase U equation.
functionType functions_
Optional cloud function objects.
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.
parcelType::constantProperties constProps_
Parcel constant properties.
scalar Dmax() const
Max diameter.
autoPtr< List< DynamicList< parcelType * > > > cellOccupancyPtr_
Cell occupancy information for each parcel, (demand driven)
const DampingModel< KinematicCloud< CloudType > > & dampingModel() const
Return const access to the damping model.
const volScalarField & rho_
Density [kg/m3].
CloudFunctionObjectList< KinematicCloud< CloudType > > functionType
Function object type.
scalar Dij(const label i, const label j) const
Mean diameter Dij.
IOdictionary particleProperties_
Dictionary of particle 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)
const dimensionedVector & g() const
Gravity.
virtual ~KinematicCloud()=default
Destructor.
const fvMesh & mesh() const
Return reference to the mesh.
Random rndGen_
Random number generator - used by some injection routines.
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.
const forceType & forces() const
Optional particle forces.
ParticleForceList< KinematicCloud< CloudType > > forceType
Force models type.
autoPtr< PackingModel< KinematicCloud< CloudType > > > packingModel_
Packing model.
const IOdictionary & outputProperties() const
Return output properties dictionary.
const IsotropyModel< KinematicCloud< CloudType > > & isotropyModel() const
Return const access to the isotropy model.
void buildCellOccupancy()
Build the cellOccupancy.
autoPtr< volScalarField::Internal > UCoeff_
Coefficient for carrier phase U equation.
void updateMesh()
Update mesh.
IOdictionary outputProperties_
Dictionary of output properties.
scalar linearKineticEnergyOfSystem() const
Total linear kinetic energy in the system.
const volScalarField & mu() const
Return carrier gas dynamic viscosity.
autoPtr< StochasticCollisionModel< KinematicCloud< CloudType > > > stochasticCollisionModel_
Stochastic collision model.
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:77
Base class for packing models.
Definition: PackingModel.H:71
List of particle forces.
Templated patch interaction model class.
Random number generator.
Definition: Random.H:60
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:66
Stores all relevant solution info for cloud.
Definition: cloudSolution.H:54
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Class used to pass data into container.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Base for a set of schemes which integrate simple ODEs which arise from semi-implcit rate expressions.
Virtual abstract base class for templated KinematicCloud.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:162
Registry of regIOobjects.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
UEqn relax()
rDeltaTY field()
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Namespace for OpenFOAM.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
CEqn solve()