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-2024 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 lagrangian::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  lagrangian::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  lagrangian::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  lagrangian::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  lagrangian::Cloud<ParticleType>
118  (
119  mesh,
120  name,
121  IDLList<ParticleType>()
122  ),
123  mesh_(mesh)
124  {}
125 
126 
127  // Member Functions
128 
129  // Access
130 
131  //- Return reference to the mesh
132  inline const fvMesh& mesh() const
133  {
134  return mesh_;
135  }
136 
137 
138  // Sources
139 
140  // Momentum
141 
142  //- Return momentum source term [kg m/s^2]
143  tmp<fvVectorMatrix> SU(const volVectorField& U) const
144  {
145  return tmp<fvVectorMatrix>
146  (
148  );
149  }
150 
151  //- Momentum transfer [kg m/s]
153  {
155  (
156  this->name() + ":UTrans",
157  this->mesh(),
159  );
160  }
161 
162  //- Momentum transfer coefficient [kg]
164  {
166  (
167  this->name() + ":UCoeffs",
168  this->mesh(),
170  );
171  }
172 
173 
174  // Energy
175 
176  //- Return sensible enthalpy source term [J/s]
178  {
179  return tmp<fvScalarMatrix>
180  (
182  );
183  }
184 
185  //- Sensible enthalpy transfer [J]
187  {
189  (
190  this->name() + ":hsTrans",
191  this->mesh(),
193  );
194  }
195 
196  //- Sensible enthalpy transfer coefficient [J/T]
198  {
200  (
201  this->name() + ":hsCoeffs",
202  this->mesh(),
204  );
205  }
206 
207  //- Return equivalent particulate emission [kg/m/s^3]
208  tmp<volScalarField> Ep() const
209  {
210  return volScalarField::New
211  (
212  this->name() + ":radiation:Ep",
213  this->mesh(),
215  );
216  }
217 
218  //- Return equivalent particulate absorption [1/m]
219  tmp<volScalarField> ap() const
220  {
221  return volScalarField::New
222  (
223  this->name() + ":radiation:ap",
224  this->mesh(),
226  );
227  }
228 
229  //- Return equivalent particulate scattering factor [1/m]
231  {
232  return volScalarField::New
233  (
234  this->name() + ":radiation:sigmap",
235  this->mesh(),
237  );
238  }
239 
240 
241  // Mass
242 
243  //- Return mass source term for specie [kg/s]
245  (
246  const label i,
247  const volScalarField& Yi
248  ) const
249  {
250  return tmp<fvScalarMatrix>
251  (
253  );
254  }
255 
256  //- Return total mass source term [kg/s]
258  {
259  return tmp<fvScalarMatrix>
260  (
262  );
263  }
264 
265  //- Return total mass source [kg/m^3/s]
267  {
269  (
271  (
272  this->name() + ":Srho",
273  this->mesh(),
275  )
276  );
277  }
278 };
279 
280 
281 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
282 
283 } // End namespace Foam
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 #endif
288 
289 // ************************************************************************* //
scalar hs(const scalar p, const scalar T) const
Definition: EtoHthermo.H:11
Base cloud calls templated on particle type.
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
Template class for intrusive linked lists.
Definition: ILList.H:67
const word & name() const
Return name.
Definition: IOobject.H:307
Base template for parcel clouds. Inserts the parcelCloudBase virtualisation layer into the class....
tmp< volScalarField > sigmap() const
Return equivalent particulate scattering factor [1/m].
tmp< volScalarField::Internal > UCoeff() const
Momentum transfer coefficient [kg].
tmp< fvVectorMatrix > SU(const volVectorField &U) const
Return momentum source term [kg m/s^2].
const fvMesh & mesh_
Reference to the mesh.
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 fvMesh & mesh() const
Return reference to the mesh.
tmp< volVectorField::Internal > UTrans() const
Momentum transfer [kg m/s].
tmp< volScalarField > Ep() const
Return equivalent particulate emission [kg/m/s^3].
tmp< volScalarField::Internal > hsCoeff() const
Sensible enthalpy transfer coefficient [J/T].
tmp< volScalarField::Internal > Srho() const
Return total mass source [kg/m^3/s].
tmp< volScalarField > ap() const
Return equivalent particulate absorption [1/m].
tmp< fvScalarMatrix > SYi(const label i, const volScalarField &Yi) const
Return mass source term for specie [kg/s].
tmp< fvScalarMatrix > Sh(const volScalarField &hs) const
Return sensible enthalpy source term [J/s].
tmp< volScalarField::Internal > hsTrans() const
Sensible enthalpy transfer [J].
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Virtual abstract base class for parcel clouds. Inserted by ParcelCloudBase into the base of the cloud...
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
A special matrix type and solver, designed for finite volume solutions of scalar equations.
U
Definition: pEqn.H:72
const dimensionedScalar mu
Atomic mass unit.
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
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
const dimensionSet dimEnergy
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
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 dimensionSet dimless
const dimensionSet dimLength
const dimensionSet dimTemperature
const dimensionSet dimAcceleration
const dimensionSet dimTime
const dimensionSet dimDensity
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
const dimensionSet dimMass
const dimensionSet dimVelocity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const word cloudName(propsDict.lookup("cloudName"))