ParcelCloudBase.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::ParcelCloudBase
26 
27 Description
28  Base template for parcel clouds. Inserts the parcelCloudBase virtualisation
29  layer into the class. Also defines default zero-return source methods to
30  enable the less functional clouds to be used in more complex situations.
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #ifndef ParcelCloudBase_H
35 #define ParcelCloudBase_H
36 
37 #include "Cloud.H"
38 #include "parcelCloudBase.H"
39 #include "fvMatrices.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 class fluidThermo;
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 template<class ParticleType>
51 class ParcelCloudBase
52 :
53  public Cloud<ParticleType>,
54  virtual public parcelCloudBase
55 {
56 protected:
57 
58  // Protected data
59 
60  //- Reference to the mesh
61  const fvMesh& mesh_;
62 
63 
64 public:
65 
66  // Constructors
67 
68  //- Construct given carrier fields
70  (
71  const word& cloudName,
72  const volScalarField& rho,
73  const volVectorField& U,
74  const volScalarField& mu,
75  const dimensionedVector& g,
76  const bool readFields = true
77  )
78  :
79  Cloud<ParticleType>(rho.mesh(), cloudName, false),
80  mesh_(rho.mesh())
81  {}
82 
83  //- Construct given carrier fields and thermo
85  (
86  const word& cloudName,
87  const volScalarField& rho,
88  const volVectorField& U,
89  const dimensionedVector& g,
90  const fluidThermo& carrierThermo,
91  const bool readFields = true
92  )
93  :
94  Cloud<ParticleType>(rho.mesh(), cloudName, false),
95  mesh_(rho.mesh())
96  {}
97 
98  //- Copy constructor with new name
100  (
102  const word& name
103  )
104  :
105  Cloud<ParticleType>(c.mesh_, name, c),
106  mesh_(c.mesh_)
107  {}
108 
109  //- Copy constructor with new name - creates bare cloud
111  (
112  const fvMesh& mesh,
113  const word& name,
115  )
116  :
117  Cloud<ParticleType>(mesh, name, IDLList<ParticleType>()),
118  mesh_(mesh)
119  {}
120 
121 
122  // Member Functions
123 
124  // Access
125 
126  //- Return reference to the mesh
127  inline const fvMesh& mesh() const
128  {
129  return mesh_;
130  }
131 
132 
133  // Sources
134 
135  // Momentum
136 
137  //- Return momentum source term [kg m/s^2]
138  tmp<fvVectorMatrix> SU(const volVectorField& U) const
139  {
140  return tmp<fvVectorMatrix>
141  (
143  );
144  }
145 
146  //- Momentum transfer [kg m/s]
148  {
150  (
151  this->name() + ":UTrans",
152  this->mesh(),
154  );
155  }
156 
157  //- Momentum transfer coefficient [kg]
159  {
161  (
162  this->name() + ":UCoeffs",
163  this->mesh(),
165  );
166  }
167 
168 
169  // Energy
170 
171  //- Return sensible enthalpy source term [J/s]
172  tmp<fvScalarMatrix> Sh(const volScalarField& hs) const
173  {
174  return tmp<fvScalarMatrix>
175  (
177  );
178  }
179 
180  //- Sensible enthalpy transfer [J]
182  {
184  (
185  this->name() + ":hsTrans",
186  this->mesh(),
188  );
189  }
190 
191  //- Sensible enthalpy transfer coefficient [J/T]
193  {
195  (
196  this->name() + ":hsCoeffs",
197  this->mesh(),
199  );
200  }
201 
202  //- Return equivalent particulate emission [kg/m/s^3]
203  tmp<volScalarField> Ep() const
204  {
205  return volScalarField::New
206  (
207  this->name() + ":radiation:Ep",
208  this->mesh(),
210  );
211  }
212 
213  //- Return equivalent particulate absorption [1/m]
214  tmp<volScalarField> ap() const
215  {
216  return volScalarField::New
217  (
218  this->name() + ":radiation:ap",
219  this->mesh(),
221  );
222  }
223 
224  //- Return equivalent particulate scattering factor [1/m]
226  {
227  return volScalarField::New
228  (
229  this->name() + ":radiation:sigmap",
230  this->mesh(),
232  );
233  }
234 
235 
236  // Mass
237 
238  //- Return mass source term for specie [kg/s]
240  (
241  const label i,
242  const volScalarField& Yi
243  ) const
244  {
245  return tmp<fvScalarMatrix>
246  (
248  );
249  }
250 
251  //- Return total mass source term [kg/s]
252  tmp<fvScalarMatrix> Srho(const volScalarField& rho) const
253  {
254  return tmp<fvScalarMatrix>
255  (
256  new fvScalarMatrix(rho, dimMass/dimTime)
257  );
258  }
259 
260  //- Return total mass source [kg/m^3/s]
262  {
264  (
266  (
267  this->name() + ":Srho",
268  this->mesh(),
270  )
271  );
272  }
273 };
274 
275 
276 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
277 
278 } // End namespace Foam
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
282 #endif
283 
284 // ************************************************************************* //
Template class for intrusive linked lists.
Definition: ILList.H:50
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
const word & name() const
Return name.
Definition: IOobject.H:315
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
tmp< fvScalarMatrix > SYi(const label i, const volScalarField &Yi) const
Return mass source term for specie [kg/s].
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject *> &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
const fvMesh & mesh_
Reference to the mesh.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
tmp< volScalarField::Internal > hsCoeff() const
Sensible enthalpy transfer coefficient [J/T].
ParcelCloudBase(const word &cloudName, const volScalarField &rho, const volVectorField &U, const volScalarField &mu, const dimensionedVector &g, const bool readFields=true)
Construct given carrier fields.
const dimensionSet dimless
tmp< volScalarField > sigmap() const
Return equivalent particulate scattering factor [1/m].
const dimensionedScalar c
Speed of light in a vacuum.
const dimensionSet dimLength
tmp< volScalarField > ap() const
Return equivalent particulate absorption [1/m].
const dimensionSet dimTime
const dimensionSet dimAcceleration
A class for handling words, derived from string.
Definition: word.H:59
tmp< fvScalarMatrix > Sh(const volScalarField &hs) const
Return sensible enthalpy source term [J/s].
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
const dimensionSet dimDensity
tmp< volScalarField::Internal > Srho() const
Return total mass source [kg/m^3/s].
static const zero Zero
Definition: zero.H:97
tmp< volVectorField::Internal > UTrans() const
Momentum transfer [kg m/s].
tmp< fvVectorMatrix > SU(const volVectorField &U) const
Return momentum source term [kg m/s^2].
Base cloud calls templated on particle type.
Definition: Cloud.H:52
const dimensionSet dimVelocity
tmp< volScalarField::Internal > UCoeff() const
Momentum transfer coefficient [kg].
const dimensionSet dimEnergy
const dimensionSet dimMass
dimensionedScalar pow3(const dimensionedScalar &ds)
Virtual abstract base class for parcel clouds. Inserted by ParcelCloudBase into the base of the cloud...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const fvMesh & mesh() const
Return reference to the mesh.
tmp< volScalarField > Ep() const
Return equivalent particulate emission [kg/m/s^3].
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
A special matrix type and solver, designed for finite volume solutions of scalar equations.
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< volScalarField::Internal > hsTrans() const
Sensible enthalpy transfer [J].
Base template for parcel clouds. Inserts the parcelCloudBase virtualisation layer into the class...
const dimensionSet dimTemperature
Namespace for OpenFOAM.