pressureGradientForce.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 "pressureGradientForce.H"
29 #include "coupledToFluid.H"
30 #include "shaped.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace Lagrangian
37 {
40  (
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::Lagrangian::pressureGradientForce::addUSup
52 (
54  LagrangianEqn<vector>& eqn
55 ) const
56 {
57  const LagrangianSubMesh& subMesh = U.mesh();
58 
59  const clouds::shaped& sCloud = cloud<clouds::shaped>();
60  const clouds::coupled& cCloud = cloud<clouds::coupled>();
61 
62  const LagrangianSubScalarField& v = sCloud.v(subMesh);
63 
64  const LagrangianSubScalarField mcByMOrMc
65  (
66  isCloud<clouds::coupledToIncompressibleFluid>()
67  ? v/cloud<clouds::coupledToIncompressibleFluid>().rhoByRhoc
68  : v*cloud<clouds::coupledToFluid>().rhoc(subMesh)
69  );
70 
71  eqn.Su += mcByMOrMc*cCloud.DUDtc(subMesh);
72 }
73 
74 
75 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
76 
78 (
79  const word& name,
80  const LagrangianMesh& mesh,
81  const dictionary& modelDict,
82  const dictionary& stateDict
83 )
84 :
86  cloudLagrangianModel(static_cast<const LagrangianModel&>(*this))
87 {
88  cloud<clouds::coupled>().DUDtc.psi();
89 }
90 
91 
92 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
93 
95 {
96  return wordList(1, cloud().U.name());
97 }
98 
99 
101 (
102  const LagrangianSubScalarField& deltaT,
105 ) const
106 {
107  assertCloud<clouds::coupledToIncompressibleFluid>();
108 
109  addUSup(U, eqn);
110 }
111 
112 
114 (
115  const LagrangianSubScalarField& deltaT,
116  const LagrangianSubScalarSubField& vOrM,
119 ) const
120 {
121  assertCloud<clouds::coupledToIncompressibleFluid, clouds::coupledToFluid>();
122 
123  addUSup(U, eqn);
124 }
125 
126 
127 // ************************************************************************* //
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
Class containing Lagrangian geometry and topology.
Base class for Lagrangian models.
virtual wordList addSupFields() const
Return the name of the velocity field.
pressureGradientForce(const word &name, const LagrangianMesh &mesh, const dictionary &modelDict, const dictionary &stateDict)
Construct from components.
virtual void addSup(const LagrangianSubScalarField &deltaT, const LagrangianSubVectorSubField &U, LagrangianEqn< vector > &eqn) const
Add a source term to the velocity equation.
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
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)
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
LagrangianSubSubField< vector > LagrangianSubVectorSubField
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
LagrangianSubField< scalar > LagrangianSubScalarField