All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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  Virtual abstract base class for parcel clouds. Inserted by ParcelCloudBase
29  into the base of the cloud template hierarchy and adds virtualisation of
30  most methods defined by the clouds.
31 
32  Note that the "evolve" method is not virtualised here. Due to the way in
33  which TrackCloudType and trackingData templating work, it is not possible
34  to virtualise this method directly. Instead it has to be wrapped. That is
35  achieved by the parcelCloud and ParcelCloud classes.
36 
37 SourceFiles
38  parcelCloudBase.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef parcelCloudBase_H
43 #define parcelCloudBase_H
44 
45 #include "volFields.H"
46 #include "fvMatricesFwd.H"
47 
48 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 class parcelCloudBase
56 {
57 public:
58 
59  //- Runtime type information
60  TypeName("parcelCloudBase");
61 
62 
63  // Constructors
64 
65  //- Null constructor
67 
68  //- Disallow default bitwise copy construction
69  parcelCloudBase(const parcelCloudBase&) = delete;
70 
71 
72  //- Destructor
73  virtual ~parcelCloudBase();
74 
75 
76  // Member Functions
77 
78  // Check
79 
80  //- Number of parcels
81  virtual label nParcels() const = 0;
82 
83  //- Total mass in system
84  virtual scalar massInSystem() const = 0;
85 
86 
87  // Fields
88 
89  //- Volume swept rate of parcels per cell
90  virtual const tmp<volScalarField> vDotSweep() const = 0;
91 
92  //- Return the particle volume fraction field
93  virtual const tmp<volScalarField> theta() const = 0;
94 
95  //- Return the particle mass fraction field
96  virtual const tmp<volScalarField> alpha() const = 0;
97 
98  //- Return the particle effective density field
99  virtual const tmp<volScalarField> rhoEff() const = 0;
100 
101 
102  // Sources
103 
104  // Momentum
105 
106  //- Return momentum source term [kg m/s^2]
107  virtual tmp<fvVectorMatrix> SU
108  (
109  const volVectorField& U
110  ) const = 0;
111 
112  //- Momentum transfer [kg m/s]
113  virtual tmp<volVectorField::Internal> UTrans() const = 0;
114 
115  //- Momentum transfer coefficient [kg]
116  virtual tmp<volScalarField::Internal> UCoeff() const = 0;
117 
118 
119  // Energy
120 
121  //- Return sensible enthalpy source term [J/s]
122  virtual tmp<fvScalarMatrix> Sh
123  (
124  const volScalarField& hs
125  ) const = 0;
126 
127  //- Sensible enthalpy transfer [J]
128  virtual tmp<volScalarField::Internal> hsTrans() const = 0;
129 
130  //- Sensible enthalpy transfer coefficient [J/K]
131  virtual tmp<volScalarField::Internal> hsCoeff() const = 0;
132 
133  //- Return equivalent particulate emission [kg/m/s^3]
134  virtual tmp<volScalarField> Ep() const = 0;
135 
136  //- Return equivalent particulate absorption [1/m]
137  virtual tmp<volScalarField> ap() const = 0;
138 
139  //- Return equivalent particulate scattering factor [1/m]
140  virtual tmp<volScalarField> sigmap() const = 0;
141 
142 
143  // Mass
144 
145  //- Return mass source term for specie [kg/s]
146  virtual tmp<fvScalarMatrix> SYi
147  (
148  const label i,
149  const volScalarField& Yi
150  ) const = 0;
151 
152  //- Return total mass source term [kg/s]
153  virtual tmp<fvScalarMatrix> Srho
154  (
155  const volScalarField& rho
156  ) const = 0;
157 
158  //- Return total mass source [kg/m^3/s]
159  virtual tmp<volScalarField::Internal> Srho() const = 0;
160 
161 
162  // I-O
163 
164  //- Print cloud information
165  virtual void info() = 0;
166 
167 
168  // Member Operators
169 
170  //- Disallow default bitwise assignment
171  void operator=(const parcelCloudBase&) = delete;
172 };
173 
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 } // End namespace Foam
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 #endif
182 
183 // ************************************************************************* //
virtual tmp< volScalarField::Internal > Srho() const =0
Return total mass source [kg/m^3/s].
virtual tmp< volScalarField::Internal > hsTrans() const =0
Sensible enthalpy transfer [J].
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
void operator=(const parcelCloudBase &)=delete
Disallow default bitwise assignment.
virtual tmp< volScalarField > Ep() const =0
Return equivalent particulate emission [kg/m/s^3].
virtual const tmp< volScalarField > theta() const =0
Return the particle volume fraction field.
virtual tmp< fvScalarMatrix > SYi(const label i, const volScalarField &Yi) const =0
Return mass source term for specie [kg/s].
virtual tmp< volScalarField::Internal > hsCoeff() const =0
Sensible enthalpy transfer coefficient [J/K].
virtual const tmp< volScalarField > vDotSweep() const =0
Volume swept rate of parcels per cell.
virtual tmp< fvVectorMatrix > SU(const volVectorField &U) const =0
Return momentum source term [kg m/s^2].
virtual scalar massInSystem() const =0
Total mass in system.
virtual const tmp< volScalarField > rhoEff() const =0
Return the particle effective density field.
parcelCloudBase()
Null constructor.
TypeName("parcelCloudBase")
Runtime type information.
virtual ~parcelCloudBase()
Destructor.
Forward declarations of fvMatrix specialisations.
U
Definition: pEqn.H:72
virtual const tmp< volScalarField > alpha() const =0
Return the particle mass fraction field.
Virtual abstract base class for parcel clouds. Inserted by ParcelCloudBase into the base of the cloud...
virtual tmp< fvScalarMatrix > Sh(const volScalarField &hs) const =0
Return sensible enthalpy source term [J/s].
virtual tmp< volScalarField > ap() const =0
Return equivalent particulate absorption [1/m].
virtual tmp< volScalarField::Internal > UCoeff() const =0
Momentum transfer coefficient [kg].
A class for managing temporary objects.
Definition: PtrList.H:53
virtual label nParcels() const =0
Number of parcels.
virtual tmp< volVectorField::Internal > UTrans() const =0
Momentum transfer [kg m/s].
virtual void info()=0
Print cloud information.
Namespace for OpenFOAM.
virtual tmp< volScalarField > sigmap() const =0
Return equivalent particulate scattering factor [1/m].