SaffmanMeiLift.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 "SaffmanMeiLift.H"
29 #include "coupledToFluid.H"
30 #include "sphericalCoupled.H"
31 #include "LagrangianmSp.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace Lagrangian
38 {
41  (
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
53 Foam::Lagrangian::SaffmanMeiLift::calcL
54 (
55  const LagrangianModelRef& model,
56  const LagrangianSubMesh& subMesh
57 ) const
58 {
59  using namespace constant::mathematical;
60 
61  const clouds::spherical& sCloud = cloud<clouds::spherical>();
62  const clouds::coupled& cCloud = cloud<clouds::coupled>();
63  const clouds::sphericalCoupled& scCloud = cloud<clouds::sphericalCoupled>();
64 
65  const LagrangianSubVectorSubField U(cloud().U(model, subMesh));
66  const LagrangianSubScalarSubField d(sCloud.d(model, subMesh));
67  const LagrangianSubScalarField& v = sCloud.v(model, subMesh);
68  const LagrangianSubScalarField& nuc = cCloud.nuc(model, subMesh);
69  const LagrangianSubScalarField& Re = scCloud.Re(model, subMesh);
70  const LagrangianSubVectorField& curlUc = cCloud.curlUc(model, subMesh);
71 
72  const LagrangianSubScalarField mcByMOrMc
73  (
74  isCloud<clouds::coupledToIncompressibleFluid>()
75  ? v/cloud<clouds::coupledToIncompressibleFluid>().rhoByRhoc
76  : v*cloud<clouds::coupledToFluid>().rhoc(model, subMesh)
77  );
78 
79  const LagrangianSubScalarField Rew(mag(curlUc)*sqr(d)/nuc);
80  const LagrangianSubScalarField beta(Rew/(Re + rootVSmall)/2);
81  const LagrangianSubScalarField alpha(0.3314*sqrt(beta));
82  const LagrangianSubScalarField Cld
83  (
84  neg(Re - 40)*6.46*((1 - alpha)*exp(-0.1*Re) + alpha)
85  + pos(Re - 40)*6.46*0.0524*sqrt(beta*Re)
86  );
87 
88  return
90  (
91  "L:" + Foam::name(subMesh.group()),
92  -mcByMOrMc*3/(twoPi*sqrt(Rew + rootVSmall))*Cld*(*curlUc)
93  );
94 }
95 
96 
97 void Foam::Lagrangian::SaffmanMeiLift::addUSup
98 (
100  LagrangianEqn<vector>& eqn
101 ) const
102 {
103  const LagrangianSubTensorField& L = this->L(U.mesh());
104 
105  const LagrangianSubVectorField& Uc = cloud<clouds::coupled>().Uc(U.mesh());
106 
107  if (eqn.isPsi(U))
108  {
109  eqn.Su += L & Uc;
110  eqn -= Lagrangianm::Sp(L, U);
111  }
112  else if (eqn.isPsi(Uc))
113  {
114  eqn += Lagrangianm::Sp(L, Uc);
115  eqn.Su -= L & U;
116  }
117  else
118  {
119  eqn.Su += L & (Uc - U);
120  }
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
125 
127 (
128  const word& name,
129  const LagrangianMesh& mesh,
130  const dictionary& modelDict,
131  const dictionary& stateDict
132 )
133 :
135  cloudLagrangianModel(static_cast<const LagrangianModel&>(*this)),
136  L
137  (
138  cloud().derivedField<tensor>
139  (
140  *this,
141  &SaffmanMeiLift::calcL
142  )
143  )
144 {
145  cloud<clouds::coupled>().curlUc.psi();
146 }
147 
148 
149 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
150 
152 {
153  return wordList(1, cloud().U.name());
154 }
155 
156 
158 (
159  const LagrangianSubScalarField& deltaT,
162 ) const
163 {
164  assertCloud<clouds::coupledToIncompressibleFluid>();
165 
166  addUSup(U, eqn);
167 }
168 
169 
171 (
172  const LagrangianSubScalarField& deltaT,
173  const LagrangianSubScalarSubField& vOrM,
176 ) const
177 {
178  assertCloud<clouds::coupledToIncompressibleFluid, clouds::coupledToFluid>();
179 
180  addUSup(U, eqn);
181 }
182 
183 
184 // ************************************************************************* //
Functions for calculating sources in a Lagrangian equation.
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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,.
This class stores the coefficients of a Lagrangian equation, and facilitates solving that equation an...
Definition: LagrangianEqn.H:56
Class containing Lagrangian geometry and topology.
Base class for Lagrangian models.
const LagrangianMesh & mesh() const
The mesh.
const word & name() const
The source name.
virtual wordList addSupFields() const
Return the name of the velocity field.
virtual void addSup(const LagrangianSubScalarField &deltaT, const LagrangianSubVectorSubField &U, LagrangianEqn< vector > &eqn) const
Add a source term to the velocity equation.
SaffmanMeiLift(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict, const dictionary &stateDict)
Construct from components.
const CloudDerivedField< tensor > & L
Lift force.
Base class for Lagrangian models that refer to a cloud. Not a Lagrangian model in itself....
const Cloud & cloud() const
Get a reference to the cloud.
Base class for clouds. Provides a basic evolution algorithm, models, and a database for caching deriv...
Definition: cloud.H:63
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
static dictionary stateDict(const word &name, const objectRegistry &db)
Construct and return the state dictionary for reading.
Definition: stateModel.C:137
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(LagrangianModel, constantCoefficientVirtualMass, dictionary)
defineTypeNameAndDebug(constantCoefficientVirtualMass, 0)
tmp< LagrangianEqn< Type > > Sp(const LagrangianSubField< SpType > &Sp, const LagrangianSubSubField< Type > &psi)
const scalar twoPi(2 *pi)
Namespace for OpenFOAM.
dimensionedScalar pos(const dimensionedScalar &ds)
List< word > wordList
A List of words.
Definition: fileName.H:54
dimensionedScalar exp(const dimensionedScalar &ds)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
LagrangianSubSubField< vector > LagrangianSubVectorSubField
dimensionedScalar neg(const dimensionedScalar &ds)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
LagrangianSubSubField< scalar > LagrangianSubScalarSubField
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
LagrangianSubField< tensor > LagrangianSubTensorField
LagrangianSubField< scalar > LagrangianSubScalarField
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
LagrangianSubField< vector > LagrangianSubVectorField
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97