parcel.C
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) 2025 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 \*---------------------------------------------------------------------------*/
25 
26 #include "parcel.H"
27 #include "cloud_fvModel.H"
28 #include "cloud_functionObject.H"
29 #include "LagrangiancDdt.H"
30 #include "LagrangianmDdt.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace clouds
37 {
40 }
41 namespace fv
42 {
44 }
45 namespace functionObjects
46 {
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
53 
54 void Foam::clouds::parcel::initialise(const bool predict)
55 {
56  cloud::initialise(predict);
57  coupled::initialise(predict);
58 }
59 
60 
62 {
65 }
66 
67 
69 (
70  const LagrangianSubMesh& subMesh
71 ) const
72 {
73  const LagrangianSubScalarSubField& m = this->m.ref(subMesh);
74  const LagrangianSubVectorSubField& U = this->U.ref(subMesh);
75 
76  return
78  ? (Lagrangianc::Ddt(m, U) - Lagrangianc::Ddt(m)*U)/m
80 }
81 
82 
84 {
85  const bool dUdt = tracking == trackingType::parabolic;
86 
87  const LagrangianSubMesh subMesh = this->mesh().subNone();
88 
89  LagrangianSubScalarSubField& number = this->number.ref(subMesh);
90  LagrangianSubScalarSubField& m = this->m.ref(subMesh);
91  LagrangianSubVectorSubField& U = this->U.ref(subMesh);
92 
93  bool result = false;
94 
95  if (LagrangianModels().addsSupToField(word::null))
96  {
97  result = Lagrangianm::initDdt(dimless, number) || result;
98  }
99 
100  if (LagrangianModels().addsSupToField(m))
101  {
102  result = Lagrangianm::initDdt(dimless, m, dUdt) || result;
103 
104  if (context != contextType::functionObject)
105  {
106  result = Lagrangianm::initDdt(dimVolume, rhoc(subMesh)) || result;
107  }
108  }
109 
110  result = Lagrangianm::initDdt(dimVolume, U, dUdt) || result;
111 
112  if (context != contextType::functionObject)
113  {
114  result = Lagrangianm::initDdt(dimVolume, Uc(subMesh)) || result;
115  }
116 
117  return result;
118 }
119 
120 
122 (
123  const LagrangianSubScalarField& deltaT,
124  const bool final
125 )
126 {
127  const LagrangianSubMesh& subMesh = deltaT.mesh();
128 
129  LagrangianSubScalarSubField& number = this->number.ref(subMesh);
130  LagrangianSubScalarSubField& m = this->m.ref(subMesh);
131  const LagrangianSubScalarSubField rho(subMesh.sub(this->rho));
132  LagrangianSubVectorSubField& U = this->U.ref(subMesh);
133 
134  // Solve the number equation if a model provides a number source
135  if (LagrangianModels().addsSupToField(number))
136  {
137  LagrangianEqn<scalar> numberEqn
138  (
139  Lagrangianm::Ddt(deltaT, number)
140  ==
141  number*LagrangianModels().source(deltaT)
142  );
143 
144  numberEqn.solve(final);
145  }
146 
147  // Solve the mass equation if a model provides a mass source
148  if (LagrangianModels().addsSupToField(m))
149  {
151  (
152  Lagrangianm::Ddt(deltaT, number, m)
153  ==
154  number*LagrangianModels().source(deltaT, m)
155  );
156 
157  mEqn.solve(final);
158 
159  // Correct the diameter, assuming the density remains constant
161 
162  if (context != contextType::functionObject && final)
163  {
165  (
166  Lagrangianm::noDdt(deltaT, dimVolume, rhoc(subMesh))
167  ==
168  LagrangianModels().sourceProxy(deltaT, m, rhoc(subMesh))
169  );
170 
171  carrierEqn(rhoc) += number*mcEqn;
172  }
173  }
174 
175  // Solve the velocity equation
176  {
178  (
179  Lagrangianm::Ddt(deltaT, m, U)
180  ==
181  LagrangianModels().source(deltaT, m, U)
182  );
183 
184  UEqn.solve(final);
185 
186  if (context != contextType::functionObject && final)
187  {
189  (
190  Lagrangianm::noDdt(deltaT, dimMass, Uc(subMesh))
191  ==
192  LagrangianModels().sourceProxy(deltaT, m, U, Uc(subMesh))
193  );
194 
195  carrierEqn(Uc) += number*UcEqn;
196  }
197  }
198 }
199 
200 
201 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
202 
204 (
205  const polyMesh& mesh,
206  const word& name,
207  const contextType context,
210 )
211 :
212  cloud(mesh, name, context, readOption, writeOption),
213  grouped(static_cast<const cloud&>(*this)),
214  spherical(*this, *this),
215  massive(*this, *this),
216  coupledToFluid(static_cast<const cloud&>(*this)),
217  sphericalCoupled(*this, *this, *this)
218 {
220 }
221 
222 
223 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
224 
226 {}
227 
228 
229 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
230 
232 {
233  // Pre-solve operations ...
235 
236  cloud::solve();
237 
238  // Post-solve operations ...
239 }
240 
241 
242 // ************************************************************************* //
Functions for calculating the time derivative for a Lagrangian equation.
Functions for calculating the time derivative for a Lagrangian equation.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const
Return mesh.
readOption
Enumeration defining the read options.
Definition: IOobject.H:117
writeOption
Enumeration defining the write options.
Definition: IOobject.H:126
This class stores the coefficients of a Lagrangian equation, and facilitates solving that equation an...
Definition: LagrangianEqn.H:56
void solve(const bool final)
Solve.
List of Lagrangian models, constructed as a (Lagrangian) mesh object. Provides similar functions to t...
bool addsSupToField(const word &) const
Return true if the LagrangianModels adds a source term to the.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
SubList< Type > sub(const List< Type > &list) const
Return a sub-list corresponding to this sub-mesh.
Base class for clouds. Provides a basic evolution algorithm, models, and a database for caching deriv...
Definition: cloud.H:63
contextType
Context in which this cloud is used.
Definition: cloud.H:222
virtual void initialise(const bool predict)
Initialisation hook.
Definition: cloud.C:476
virtual void partition()
Partition hook.
Definition: cloud.C:484
virtual void solve()
Solve the cloud's evolution over the current time-step.
Definition: cloud.C:634
Base class for clouds which are coupled to a fluid.
void initialise(const bool predict)
Initialisation hook.
Definition: coupled.C:188
void clearCarrierEqns()
Clear the carrier equations.
Definition: coupled.C:171
void partition()
Partition hook.
Definition: coupled.C:194
Base class for clouds in which particles are grouped into parcels.
Definition: grouped.H:51
Base class for clouds with particles with mass.
Definition: massive.H:51
Basic cloud with spherical, variable density, particles, grouped into parcels.
Definition: parcel.H:72
virtual bool reCalculateModified()
Do we need to re-calculate particles that are modified?
Definition: parcel.C:83
virtual tmp< LagrangianSubVectorField > dUdt(const LagrangianSubMesh &) const
Return the acceleration with which to do second-order tracking.
Definition: parcel.C:69
virtual void initialise(const bool predict)
Initialisation hook.
Definition: parcel.C:54
parcel(const polyMesh &mesh, const word &name, const contextType context, const IOobject::readOption readOption=IOobject::READ_IF_PRESENT, const IOobject::writeOption writeOption=IOobject::AUTO_WRITE)
Construct from a mesh and a name.
Definition: parcel.C:204
virtual ~parcel()
Destructor.
Definition: parcel.C:225
virtual void partition()
Partition hook.
Definition: parcel.C:61
virtual void calculate(const LagrangianSubScalarField &deltaT, const bool final)
Update the cloud properties.
Definition: parcel.C:122
virtual void solve()
Solve the cloud's evolution over the current time-step.
Definition: parcel.C:231
Base class for clouds of spherical particles which are coupled to a fluid.
Base class for clouds with spherical particles.
Definition: spherical.H:53
void correct(const LagrangianSubScalarSubField &v)
Correct the shape to match the given volume.
Definition: spherical.C:79
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
fvVectorMatrix & UEqn
Definition: UEqn.H:11
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
U
Definition: pEqn.H:72
tmp< LagrangianSubField< Type > > Ddt(const LagrangianSubSubField< Type > &psi)
tmp< LagrangianEqn< Type > > Ddt(const LagrangianSubScalarField &deltaT, LagrangianSubSubField< Type > &psi)
bool initDdt(const dimensionSet &mDims, const LagrangianSubSubField< Type > &psi, const bool instantaneousDdt=false)
tmp< LagrangianEqn< Type > > noDdt(const LagrangianSubScalarField &deltaT, const dimensionSet &mDims, const LagrangianSubSubField< Type > &psi)
defineTypeNameAndDebug(coupled, 0)
addToRunTimeSelectionTable(cloud, kinematicParcel, polyMesh)
makeCloudFunctionObject(kinematicParcel)
makeCloudFvModel(kinematicParcel)
Namespace for OpenFOAM.
const dimensionSet dimless
const dimensionSet dimVolume
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
tmp< DimensionedField< Type, GeoMesh, SubField > > toSubField(const DimensionedField< Type, GeoMesh, Field > &)
Return a temporary sub-field from a reference to a field.
labelList fv(nPoints)