All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2021 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 a solver in order
30  to support the coupled simulation of multiple clouds.
31 
32 SourceFiles
33  parcelCloudList.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef parcelCloudList_H
38 #define parcelCloudList_H
39 
40 #include "parcelCloud.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 class parcelCloudList
50 :
51  public PtrList<parcelCloud>
52 {
53 private:
54 
55  // Private data
56 
57  //- Reference to the mesh
58  const fvMesh& mesh_;
59 
60 
61  // Private member functions
62 
63  //- Initialise the pointer list of clouds
64  template <class ... Args>
65  void initialise(const Args& ... args);
66 
67 
68 public:
69 
70  // Static Data Members
71 
72  //- The name of the clouds
73  static const word cloudsName;
74 
75 
76  // Constructors
77 
78  //- Construct with given carrier fields
80  (
81  const volScalarField& rho,
82  const volVectorField& U,
83  const volScalarField& mu,
84  const dimensionedVector& g
85  );
86 
87  //- Construct with given carrier fields and thermo
89  (
90  const volScalarField& rho,
91  const volVectorField& U,
92  const dimensionedVector& g,
93  const fluidThermo& carrierThermo
94  );
95 
96  //- Disallow default bitwise copy construction
97  parcelCloudList(const parcelCloudList&) = delete;
98 
99 
100  //- Destructor
102 
103 
104  // Member Functions
105 
106  // Fields
107 
108  //- Return the particle volume fraction field
109  const tmp<volScalarField> theta() const;
110 
111 
112  // Sources
113 
114  // Momentum
115 
116  //- Return momentum source term [kg m/s^2]
117  tmp<fvVectorMatrix> SU(const volVectorField& U) const;
118 
119  //- Momentum transfer [kg m/s]
121 
122  //- Momentum transfer coefficient [kg]
124 
125 
126  // Energy
127 
128  //- Return sensible enthalpy source term [J/s]
129  tmp<fvScalarMatrix> Sh(const volScalarField& hs) const;
130 
131  //- Sensible enthalpy transfer [J]
133 
134  //- Sensible enthalpy transfer coefficient [J/K]
136 
137  //- Return equivalent particulate emission [kg/m/s^3]
138  tmp<volScalarField> Ep() const;
139 
140  //- Return equivalent particulate absorption [1/m]
141  tmp<volScalarField> ap() const;
142 
143  //- Return equivalent particulate scattering factor [1/m]
144  tmp<volScalarField> sigmap() const;
145 
146 
147  // Mass
148 
149  //- Return mass source term for specie [kg/s]
151  (
152  const label speciei,
153  const volScalarField& Yi
154  ) const;
155 
156  //- Return total mass source term [kg/s]
158 
159  //- Return total mass source [kg/m^3/s]
161 
162 
163  // I-O
164 
165  //- Print cloud information
166  void info();
167 
168 
169  // Evolution
170 
171  //- Evolve the cloud
172  void evolve();
173 
174 
175  //- Call this before a topology change. Stores the particles global
176  // positions in the database for use during mapping.
177  void storeGlobalPositions();
178 
179 
180  // Member Operators
181 
182  //- Disallow default bitwise assignment
183  void operator=(const parcelCloudList&) = delete;
184 };
185 
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 } // End namespace Foam
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 #endif
194 
195 // ************************************************************************* //
const tmp< volScalarField > theta() const
Return the particle volume fraction field.
tmp< fvScalarMatrix > SYi(const label speciei, const volScalarField &Yi) const
Return mass source term for specie [kg/s].
void operator=(const parcelCloudList &)=delete
Disallow default bitwise assignment.
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
tmp< volScalarField::Internal > hsCoeff() const
Sensible enthalpy transfer coefficient [J/K].
parcelCloudList(const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g)
Construct with given carrier fields.
tmp< volScalarField > sigmap() const
Return equivalent particulate scattering factor [1/m].
~parcelCloudList()
Destructor.
List of parcel clouds, with the same interface as an individual parcel cloud. This is the object that...
tmp< fvVectorMatrix > SU(const volVectorField &U) const
Return momentum source term [kg m/s^2].
tmp< volVectorField::Internal > UTrans() const
Momentum transfer [kg m/s].
void info()
Print cloud information.
A class for handling words, derived from string.
Definition: word.H:59
tmp< volScalarField::Internal > UCoeff() const
Momentum transfer coefficient [kg].
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
void storeGlobalPositions()
Call this before a topology change. Stores the particles global.
tmp< volScalarField::Internal > Srho() const
Return total mass source [kg/m^3/s].
const dimensionedScalar mu
Atomic mass unit.
tmp< volScalarField > ap() const
Return equivalent particulate absorption [1/m].
tmp< volScalarField::Internal > hsTrans() const
Sensible enthalpy transfer [J].
U
Definition: pEqn.H:72
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void evolve()
Evolve the cloud.
static const word cloudsName
The name of the clouds.
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedVector & g
Foam::argList args(argc, argv)
tmp< volScalarField > Ep() const
Return equivalent particulate emission [kg/m/s^3].
tmp< fvScalarMatrix > Sh(const volScalarField &hs) const
Return sensible enthalpy source term [J/s].
Namespace for OpenFOAM.