parcelCloudList.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2020-2022 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::parcelCloudList
26 
27 Description
28  List of parcel clouds, with the same interface as an individual parcel
29  cloud. This is the object that should be constructed by an fvModel, or any
30  system that can call this class' mesh change functions. A solver should
31  *not* construct this object, as that would not provide a mechanism for the
32  mesh change functions to be executed. A solver should construct a
33  parcelClouds object instead.
34 
35 SourceFiles
36  parcelCloudList.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef parcelCloudList_H
41 #define parcelCloudList_H
42 
43 #include "parcelCloud.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 class parcelCloudList
53 :
54  public PtrList<parcelCloud>
55 {
56 private:
57 
58  // Private Static Member Functions
59 
60  //- Get the cloud names for this case
61  static wordList cloudNames(const objectRegistry& db);
62 
63 
64  // Private data
65 
66  //- Reference to the mesh
67  const fvMesh& mesh_;
68 
69 
70 public:
71 
72  // Static Data Members
73 
74  //- The default cloud name
75  static const word defaultCloudName;
76 
77  //- The default cloud names (i.e., a list of length one with the
78  // defaultCloudName as its value)
79  static const wordList defaultCloudNames;
80 
81  //- The name of the clouds file in which multiple cloud names are
82  // specified
83  static const word cloudNamesName;
84 
85 
86  // Constructors
87 
88  //- Construct specified clouds with given carrier fields
90  (
91  const wordList& cloudNames,
92  const volScalarField& rho,
93  const volVectorField& U,
94  const volScalarField& mu,
95  const dimensionedVector& g
96  );
97 
98  //- Construct specified clouds with given carrier fields and thermo
100  (
101  const wordList& cloudNames,
102  const volScalarField& rho,
103  const volVectorField& U,
104  const dimensionedVector& g,
105  const fluidThermo& carrierThermo
106  );
107 
108  //- Construct detected clouds with given carrier fields
110  (
111  const volScalarField& rho,
112  const volVectorField& U,
113  const volScalarField& mu,
114  const dimensionedVector& g
115  );
116 
117  //- Construct detected clouds with given carrier fields and thermo
119  (
120  const volScalarField& rho,
121  const volVectorField& U,
122  const dimensionedVector& g,
123  const fluidThermo& carrierThermo
124  );
125 
126  //- Disallow default bitwise copy construction
127  parcelCloudList(const parcelCloudList&) = delete;
128 
129 
130  //- Destructor
132 
133 
134  // Member Functions
135 
136  // Fields
137 
138  //- Return the particle volume fraction field
139  const tmp<volScalarField> theta() const;
140 
141 
142  // Sources
143 
144  // Momentum
145 
146  //- Return momentum source term [kg m/s^2]
147  tmp<fvVectorMatrix> SU(const volVectorField& U) const;
148 
149  //- Momentum transfer [kg m/s]
151 
152  //- Momentum transfer coefficient [kg]
154 
155 
156  // Energy
157 
158  //- Return sensible enthalpy source term [J/s]
159  tmp<fvScalarMatrix> Sh(const volScalarField& hs) const;
160 
161  //- Sensible enthalpy transfer [J]
163 
164  //- Sensible enthalpy transfer coefficient [J/K]
166 
167  //- Return equivalent particulate emission [kg/m/s^3]
168  tmp<volScalarField> Ep() const;
169 
170  //- Return equivalent particulate absorption [1/m]
171  tmp<volScalarField> ap() const;
172 
173  //- Return equivalent particulate scattering factor [1/m]
174  tmp<volScalarField> sigmap() const;
175 
176 
177  // Mass
178 
179  //- Return mass source term for specie [kg/s]
181  (
182  const label speciei,
183  const volScalarField& Yi
184  ) const;
185 
186  //- Return total mass source term [kg/s]
188 
189  //- Return total mass source [kg/m^3/s]
191 
192 
193  // I-O
194 
195  //- Print cloud information
196  void info();
197 
198 
199  // Evolution
200 
201  //- Evolve the cloud
202  void evolve();
203 
204 
205  // Mesh changes
206 
207  //- Call this before a topology change. Stores the particles global
208  // positions in the database for use during mapping.
209  void storeGlobalPositions();
210 
211  //- Update topology using the given map
212  void topoChange(const polyTopoChangeMap&);
213 
214  //- Update from another mesh using the given map
215  void mapMesh(const polyMeshMap&);
216 
217  //- Redistribute or update using the given distribution map
218  void distribute(const polyDistributionMap&);
219 
220 
221  // Member Operators
222 
223  //- Disallow default bitwise assignment
224  void operator=(const parcelCloudList&) = delete;
225 };
226 
227 
228 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
229 
230 } // End namespace Foam
231 
232 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
233 
234 #endif
235 
236 // ************************************************************************* //
Generic GeometricField class.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:57
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Registry of regIOobjects.
List of parcel clouds, with the same interface as an individual parcel cloud. This is the object that...
static const wordList defaultCloudNames
The default cloud names (i.e., a list of length one with the.
tmp< fvScalarMatrix > Sh(const volScalarField &hs) const
Return sensible enthalpy source term [J/s].
tmp< volScalarField::Internal > UCoeff() const
Momentum transfer coefficient [kg].
tmp< volVectorField::Internal > UTrans() const
Momentum transfer [kg m/s].
const tmp< volScalarField > theta() const
Return the particle volume fraction field.
~parcelCloudList()
Destructor.
void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
tmp< fvScalarMatrix > SYi(const label speciei, const volScalarField &Yi) const
Return mass source term for specie [kg/s].
void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
tmp< volScalarField > sigmap() const
Return equivalent particulate scattering factor [1/m].
void evolve()
Evolve the cloud.
parcelCloudList(const wordList &cloudNames, const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g)
Construct specified clouds with given carrier fields.
void operator=(const parcelCloudList &)=delete
Disallow default bitwise assignment.
tmp< volScalarField::Internal > hsCoeff() const
Sensible enthalpy transfer coefficient [J/K].
static const word defaultCloudName
The default cloud name.
void info()
Print cloud information.
tmp< fvVectorMatrix > SU(const volVectorField &U) const
Return momentum source term [kg m/s^2].
tmp< volScalarField::Internal > Srho() const
Return total mass source [kg/m^3/s].
void storeGlobalPositions()
Call this before a topology change. Stores the particles global.
tmp< volScalarField > Ep() const
Return equivalent particulate emission [kg/m/s^3].
tmp< volScalarField::Internal > hsTrans() const
Sensible enthalpy transfer [J].
static const word cloudNamesName
The name of the clouds file in which multiple cloud names are.
tmp< volScalarField > ap() const
Return equivalent particulate absorption [1/m].
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
const dimensionedScalar mu
Atomic mass unit.
Namespace for OpenFOAM.
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