gravity.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 "gravity.H"
30 #include "coupledToFluid.H"
31 #include "massive.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace Lagrangian
38 {
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  const word& name,
50  const LagrangianMesh& mesh,
51  const dictionary& modelDict,
52  const dictionary& stateDict
53 )
54 :
56  cloudLagrangianModel(static_cast<const LagrangianModel&>(*this)),
58 {}
59 
60 
61 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
62 
64 {
65  return wordList(1, cloud().U.name());
66 }
67 
68 
70 (
71  const LagrangianSubScalarField& deltaT,
74 ) const
75 {
76  if
77  (
78  isCloud<clouds::coupled>()
79  && eqn.isPsi(cloud<clouds::coupled>().Uc(U.mesh()))
80  )
81  {
82  // Gravity does not contribute to the carrier momentum equation
83  }
84  else if (isCloud<clouds::coupledToIncompressibleFluid>())
85  {
86  const clouds::coupledToIncompressibleFluid& ctifCloud =
88 
89  const dimensionedScalar& rhoByRhoc = ctifCloud.rhoByRhoc;
90 
91  eqn.Su += g*(1 - 1/rhoByRhoc);
92  }
93  else
94  {
95  eqn.Su += g;
96  }
97 }
98 
99 
101 (
102  const LagrangianSubScalarField& deltaT,
103  const LagrangianSubScalarSubField& vOrM,
106 ) const
107 {
108  if
109  (
110  isCloud<clouds::coupled>()
111  && eqn.isPsi(cloud<clouds::coupled>().Uc(U.mesh()))
112  )
113  {
114  // Gravity does not contribute to the carrier momentum equation
115  }
116  else if (isCloud<clouds::coupledToIncompressibleFluid>())
117  {
118  const clouds::coupledToIncompressibleFluid& ctifCloud =
120 
121  const dimensionedScalar& rhoByRhoc = ctifCloud.rhoByRhoc;
122 
123  eqn.Su += vOrM*g*(1 - 1/rhoByRhoc);
124  }
125  else if (isCloud<clouds::coupledToFluid>() && isCloud<clouds::massive>())
126  {
127  const clouds::massive& mCloud = cloud<clouds::massive>();
128  const clouds::coupledToFluid& ctfCloud =
130 
131  const LagrangianSubScalarSubField rho(mCloud.rho(U.mesh()));
132  const LagrangianSubScalarField& rhoc = ctfCloud.rhoc(U.mesh());
133 
134  eqn.Su += vOrM*g*(1 - rhoc/rho);
135  }
136  else if (isCloud<clouds::coupledToFluid>() && !isCloud<clouds::massive>())
137  {
139  }
140  else
141  {
142  eqn.Su += vOrM*g;
143  }
144 }
145 
146 
147 // ************************************************************************* //
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...
This class stores the coefficients of a Lagrangian equation, and facilitates solving that equation an...
Definition: LagrangianEqn.H:56
LagrangianCoeff< Type, false > Su
Explicit coefficient.
bool isPsi(const LagrangianSubField< Type, PrimitiveField > &psi) const
Return whether the given field is that of the equation.
Class containing Lagrangian geometry and topology.
Base class for Lagrangian models.
virtual wordList addSupFields() const
Return the name of the velocity field.
Definition: gravity.C:63
virtual void addSup(const LagrangianSubScalarField &deltaT, const LagrangianSubVectorSubField &U, LagrangianEqn< vector > &eqn) const
Add a source term to the velocity equation.
Definition: gravity.C:70
gravity(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict, const dictionary &stateDict)
Construct from components.
Definition: gravity.C:48
Base class for Lagrangian models that refer to a cloud. Not a Lagrangian model in itself....
Base class for clouds. Provides a basic evolution algorithm, models, and a database for caching deriv...
Definition: cloud.H:63
Base class for clouds which are coupled to a fluid.
const CarrierField< scalar > & rhoc
Carrier density.
Base class for clouds which are coupled to an incompressible fluid.
const dimensionedScalar rhoByRhoc
Cloud/carrier density ratio.
Base class for clouds with particles with mass.
Definition: massive.H:51
CloudStateField< scalar > & rho
Density.
Definition: massive.H:66
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
U
Definition: pEqn.H:72
addToRunTimeSelectionTable(LagrangianModel, constantCoefficientVirtualMass, dictionary)
defineTypeNameAndDebug(constantCoefficientVirtualMass, 0)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const UniformDimensionedField< Type > & lookupUniformDimensionedField(const objectRegistry &db, const word &name)