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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "ParticleForceList.H"
29 #include "entry.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class CloudType>
35 (
36  CloudType& owner,
37  const fvMesh& mesh
38 )
39 :
41  owner_(owner),
42  mesh_(mesh),
43  dict_(dictionary::null),
44  calcCoupled_(true),
45  calcNonCoupled_(true)
46 {}
47 
48 
49 template<class CloudType>
51 (
52  CloudType& owner,
53  const fvMesh& mesh,
54  const dictionary& dict,
55  const bool readFields
56 )
57 :
59  owner_(owner),
60  mesh_(mesh),
61  dict_(dict),
62  calcCoupled_(true),
63  calcNonCoupled_(true)
64 {
65  if (readFields)
66  {
67  Info<< "Constructing particle forces" << endl;
68 
69  this->resize(dict.size());
70 
71  label count = 0;
72  for (const entry& dEntry : dict)
73  {
74  const word& model = dEntry.keyword();
75  if (dEntry.isDict())
76  {
77  this->set
78  (
79  count,
81  (
82  owner,
83  mesh,
84  dEntry.dict(),
85  model
86  )
87  );
88  }
89  else
90  {
91  this->set
92  (
93  count,
95  (
96  owner,
97  mesh,
98  dict,
99  model
100  )
101  );
102  }
103 
104  ++count;
105  }
106 
107  if (!count)
108  {
109  Info<< " none" << endl;
110  }
111  }
112 }
113 
114 
115 template<class CloudType>
117 (
118  const ParticleForceList& pf
119 )
120 :
122  owner_(pf.owner_),
123  mesh_(pf.mesh_),
124  dict_(pf.dict_)
125 {}
126 
127 
128 // * * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * //
129 
130 template<class CloudType>
132 {}
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
137 template<class CloudType>
139 {
140  forAll(*this, i)
141  {
142  this->operator[](i).cacheFields(store);
143  }
144 }
145 
146 
147 template<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 
172 template<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 
198 template<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 // ************************************************************************* //
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:67
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::BitOps::set
void set(List< bool > &bools, const labelRange &range)
Set the specified range 'on' in a boolList.
Definition: BitOps.C:37
Foam::ParticleForceList::cacheFields
virtual void cacheFields(const bool store)
Cache fields.
Definition: ParticleForceList.C:138
Foam::ParticleForceList::calcNonCoupled
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.
Definition: ParticleForceList.C:174
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::ParticleForceList::~ParticleForceList
virtual ~ParticleForceList()
Destructor.
Definition: ParticleForceList.C:131
entry.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
resize
patchWriters resize(patchIds.size())
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
ParticleForceList.H
Foam::forceSuSp
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:64
Foam::ParticleForce
Abstract base class for particle forces.
Definition: ParticleForce.H:59
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:62
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::ParticleForceList::ParticleForceList
ParticleForceList(CloudType &owner, const fvMesh &mesh)
Null constructor.
Definition: ParticleForceList.C:35
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam::readFields
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
Definition: ReadFieldsTemplates.C:312
Foam::ParticleForceList::calcCoupled
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.
Definition: ParticleForceList.C:149
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
Foam::ParticleForceList::massEff
virtual scalar massEff(const typename CloudType::parcelType &p, const typename CloudType::parcelType::trackingData &td, const scalar mass) const
Return the effective mass.
Definition: ParticleForceList.C:200
Foam::Re
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
Foam::DSMCCloud::parcelType
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:220
Foam::ParticleForceList
List of particle forces.
Definition: ParticleForceList.H:53