ParticleForceList.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "ParticleForceList.H"
30#include "entry.H"
31
32// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33
34template<class CloudType>
36(
37 CloudType& owner,
38 const fvMesh& mesh
39)
40:
42 owner_(owner),
43 mesh_(mesh),
44 dict_(),
45 calcCoupled_(true),
46 calcNonCoupled_(true)
47{}
48
49
50template<class CloudType>
52(
53 CloudType& owner,
54 const fvMesh& mesh,
55 const dictionary& dict,
56 const bool readFields
57)
58:
60 owner_(owner),
61 mesh_(mesh),
62 dict_(dict),
63 calcCoupled_(true),
64 calcNonCoupled_(true)
65{
66 if (readFields)
67 {
68 Info<< "Constructing particle forces" << endl;
69
70 this->resize(dict.size());
71
72 label count = 0;
73 for (const entry& dEntry : dict)
74 {
75 const word& modelName = dEntry.keyword();
76 if (dEntry.isDict())
77 {
78 this->set
79 (
80 count,
82 (
84 mesh,
85 dEntry.dict(),
86 modelName
87 )
88 );
89 }
90 else
91 {
92 this->set
93 (
94 count,
96 (
97 owner,
98 mesh,
99 dict,
100 modelName
101 )
102 );
103 }
104 ++count;
105 }
106
107 if (!count)
108 {
109 Info<< " none" << endl;
110 }
111 }
112}
113
114
115template<class CloudType>
117(
118 const ParticleForceList& pf
119)
120:
122 owner_(pf.owner_),
123 mesh_(pf.mesh_),
124 dict_(pf.dict_)
126
127
128// * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
129
130template<class CloudType>
132{}
133
134
135// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136
137template<class CloudType>
140 forAll(*this, i)
141 {
142 this->operator[](i).cacheFields(store);
143 }
144}
145
146
147template<class CloudType>
149(
150 const typename CloudType::parcelType& p,
151 const typename CloudType::parcelType::trackingData& td,
152 const scalar dt,
153 const scalar mass,
154 const scalar Re,
155 const scalar muc
156) const
157{
158 forceSuSp value(Zero);
159
160 if (calcCoupled_)
161 {
162 forAll(*this, i)
163 {
164 value += this->operator[](i).calcCoupled(p, td, dt, mass, Re, muc);
165 }
166 }
167
168 return value;
169}
170
171
172template<class CloudType>
174(
175 const typename CloudType::parcelType& p,
176 const typename CloudType::parcelType::trackingData& td,
177 const scalar dt,
178 const scalar mass,
179 const scalar Re,
180 const scalar muc
181) const
182{
183 forceSuSp value(Zero);
184
185 if (calcNonCoupled_)
186 {
187 forAll(*this, i)
188 {
189 value +=
190 this->operator[](i).calcNonCoupled(p, td, dt, mass, Re, muc);
191 }
192 }
193
194 return value;
195}
196
197
198template<class CloudType>
200(
201 const typename CloudType::parcelType& p,
202 const typename CloudType::parcelType::trackingData& td,
203 const scalar mass
204) const
205{
206 scalar massEff = mass;
207 forAll(*this, i)
208 {
209 massEff += this->operator[](i).massAdd(p, td, mass);
210 }
211
212 return massEff;
213}
214
215
216// ************************************************************************* //
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
List of particle forces.
virtual forceSuSp calcNonCoupled(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar dt, const scalar mass, const scalar Re, const scalar muc) const
Calculate the non-coupled force.
virtual scalar massEff(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar mass) const
Return the effective mass.
virtual void cacheFields(const bool store)
Cache fields.
const dictionary & dict() const
Return the forces dictionary.
const CloudType & owner() const
Return const access to the cloud owner.
virtual ~ParticleForceList()
Destructor.
const fvMesh & mesh() const
Return the mesh database.
virtual forceSuSp calcCoupled(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar dt, const scalar mass, const scalar Re, const scalar muc) const
Calculate the coupled force.
Abstract base class for particle forces.
Definition: ParticleForce.H:60
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
const ParticleForce< CloudType > * set(const label i) const
Definition: PtrList.H:138
void resize(const label newLen)
Adjust size of PtrList.
Definition: PtrList.C:103
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:70
Class used to pass data into container.
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:67
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
dynamicFvMesh & mesh
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
messageStream Info
Information stream (stdout output on master, null elsewhere)
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333