kinematicParcel.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 "kinematicParcel.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 {
43  makeCloudFvModel(kinematicParcel);
44 }
45 namespace functionObjects
46 {
47  makeCloudFunctionObject(kinematicParcel);
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
53 
55 {
56  cloud::initialise(predict);
57  coupled::initialise(predict);
58 }
59 
60 
62 {
65 }
66 
67 
70 (
71  const LagrangianSubMesh& subMesh
72 ) const
73 {
74  const LagrangianSubScalarSubField& v = this->v.ref(subMesh);
75  const LagrangianSubVectorSubField& U = this->U.ref(subMesh);
76 
77  return
79  ? (Lagrangianc::Ddt(v, U) - Lagrangianc::Ddt(v)*U)/v
81 }
82 
83 
85 {
86  const bool dUdt = tracking == trackingType::parabolic;
87 
88  const LagrangianSubMesh subMesh = this->mesh().subNone();
89 
90  LagrangianSubScalarSubField& number = this->number.ref(subMesh);
91  LagrangianSubScalarSubField& v = this->v.ref(subMesh);
92  LagrangianSubVectorSubField& U = this->U.ref(subMesh);
93 
94  bool result = false;
95 
96  if (LagrangianModels().addsSupToField(word::null))
97  {
98  result = Lagrangianm::initDdt(dimless, number) || result;
99  }
100 
101  if (LagrangianModels().addsSupToField(v))
102  {
103  result = Lagrangianm::initDdt(dimless, v, dUdt) || result;
104 
105  /*
106  if (context != contextType::functionObject)
107  {
108  result = Lagrangianm::initDdt(dimVolume, onec()) || result;
109  }
110  */
111  }
112 
113  result = Lagrangianm::initDdt(dimVolume, U, dUdt) || result;
114 
115  if (context != contextType::functionObject)
116  {
117  result = Lagrangianm::initDdt(dimVolume, Uc(subMesh)) || result;
118  }
119 
120  return result;
121 }
122 
123 
125 (
126  const LagrangianSubScalarField& deltaT,
127  const bool final
128 )
129 {
130  const LagrangianSubMesh& subMesh = deltaT.mesh();
131 
132  LagrangianSubScalarSubField& number = this->number.ref(subMesh);
133  LagrangianSubScalarSubField& v = this->v.ref(subMesh);
134  LagrangianSubVectorSubField& U = this->U.ref(subMesh);
135 
136  // Solve the number equation if a model provides a number source
137  if (LagrangianModels().addsSupToField(word::null))
138  {
139  LagrangianEqn<scalar> numberEqn
140  (
141  Lagrangianm::Ddt(deltaT, number)
142  ==
143  number*LagrangianModels().source(deltaT)
144  );
145 
146  numberEqn.solve(final);
147  }
148 
149  // Solve the volume equation if a model provides a volume source
150  if (LagrangianModels().addsSupToField(v))
151  {
153  (
154  Lagrangianm::Ddt(deltaT, v)
155  ==
156  LagrangianModels().source(deltaT, v)
157  );
158 
159  vEqn.solve(final);
160 
161  // Correct the diameter
163 
164  /*
165  if (context != contextType::functionObject && final)
166  {
167  LagrangianEqn<scalar> vcEqn
168  (
169  Lagrangianm::noDdt(deltaT, dimVolume, onec())
170  ==
171  LagrangianModels().sourceProxy(deltaT, v, onec())
172  );
173 
174  carrierEqn(one()) += number*vcEqn;
175  }
176  */
177  }
178 
179  // Solve the velocity equation
180  {
182  (
183  Lagrangianm::Ddt(deltaT, v, U)
184  ==
185  LagrangianModels().source(deltaT, v, U)
186  );
187 
188  UEqn.solve(final);
189 
190  if (context != contextType::functionObject && final)
191  {
193  (
194  Lagrangianm::noDdt(deltaT, dimVolume, Uc(subMesh))
195  ==
196  LagrangianModels().sourceProxy(deltaT, v, U, Uc(subMesh))
197  );
198 
199  carrierEqn(Uc) += number*UcEqn;
200  }
201  }
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
206 
208 (
209  const polyMesh& mesh,
210  const word& name,
211  const contextType context,
214 )
215 :
216  cloud(mesh, name, context, readOption, writeOption),
217  grouped(static_cast<const cloud&>(*this)),
218  spherical(*this, *this),
219  coupledToIncompressibleFluid(static_cast<const cloud&>(*this)),
220  sphericalCoupled(*this, *this, *this)
221 {
223 }
224 
225 
226 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
227 
229 {}
230 
231 
232 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
233 
235 {
236  // Pre-solve operations ...
239 
240  cloud::solve();
241 
242  // Post-solve operations ...
243 }
244 
245 
246 // ************************************************************************* //
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...
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 an incompressible fluid.
void updateNuc()
Update the carrier kinematic viscosity.
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
Basic cloud with spherical, constant density, particles, grouped into parcels.
virtual ~kinematicParcel()
Destructor.
virtual bool reCalculateModified()
Do we need to re-calculate particles that are modified?
virtual tmp< LagrangianSubVectorField > dUdt(const LagrangianSubMesh &) const
Return the acceleration with which to do second-order tracking.
virtual void initialise(const bool predict)
Initialisation hook.
kinematicParcel(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.
virtual void partition()
Partition hook.
virtual void calculate(const LagrangianSubScalarField &deltaT, const bool final)
Update the cloud properties.
virtual void solve()
Solve the cloud's evolution over the current time-step.
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.
labelList fv(nPoints)