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 | Copyright (C) 2011-2016 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 "ParticleForceList.H"
27 #include "entry.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class CloudType>
33 (
34  CloudType& owner,
35  const fvMesh& mesh
36 )
37 :
39  owner_(owner),
40  mesh_(mesh),
41  dict_(dictionary::null),
42  calcCoupled_(true),
43  calcNonCoupled_(true)
44 {}
45 
46 
47 template<class CloudType>
49 (
50  CloudType& owner,
51  const fvMesh& mesh,
52  const dictionary& dict,
53  const bool readFields
54 )
55 :
57  owner_(owner),
58  mesh_(mesh),
59  dict_(dict),
60  calcCoupled_(true),
61  calcNonCoupled_(true)
62 {
63  if (readFields)
64  {
65  wordList modelNames(dict.toc());
66 
67  Info<< "Constructing particle forces" << endl;
68 
69  if (modelNames.size() > 0)
70  {
71  this->setSize(modelNames.size());
72 
73  label i = 0;
74  forAllConstIter(IDLList<entry>, dict, iter)
75  {
76  const word& model = iter().keyword();
77  if (iter().isDict())
78  {
79  this->set
80  (
81  i++,
83  (
84  owner,
85  mesh,
86  iter().dict(),
87  model
88  )
89  );
90  }
91  else
92  {
93  this->set
94  (
95  i++,
97  (
98  owner,
99  mesh,
100  dict,
101  model
102  )
103  );
104  }
105  }
106  }
107  else
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 scalar dt,
152  const scalar mass,
153  const scalar Re,
154  const scalar muc
155 ) const
156 {
157  forceSuSp value(Zero, 0.0);
158 
159  if (calcCoupled_)
160  {
161  forAll(*this, i)
162  {
163  value += this->operator[](i).calcCoupled(p, dt, mass, Re, muc);
164  }
165  }
166 
167  return value;
168 }
169 
170 
171 template<class CloudType>
173 (
174  const typename CloudType::parcelType& p,
175  const scalar dt,
176  const scalar mass,
177  const scalar Re,
178  const scalar muc
179 ) const
180 {
181  forceSuSp value(Zero, 0.0);
182 
183  if (calcNonCoupled_)
184  {
185  forAll(*this, i)
186  {
187  value += this->operator[](i).calcNonCoupled(p, dt, mass, Re, muc);
188  }
189  }
190 
191  return value;
192 }
193 
194 
195 template<class CloudType>
197 (
198  const typename CloudType::parcelType& p,
199  const scalar mass
200 ) const
201 {
202  scalar massEff = mass;
203  forAll(*this, i)
204  {
205  massEff += this->operator[](i).massAdd(p, mass);
206  }
207 
208  return massEff;
209 }
210 
211 
212 // ************************************************************************* //
Template class for intrusive linked lists.
Definition: ILList.H:50
dictionary dict
virtual void cacheFields(const bool store)
Cache fields.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual scalar massEff(const typename CloudType::parcelType &p, const scalar mass) const
Return the effective mass.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual forceSuSp calcNonCoupled(const typename CloudType::parcelType &p, const scalar dt, const scalar mass, const scalar Re, const scalar muc) const
Calculate the non-coupled force.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
ParticleForceList(CloudType &owner, const fvMesh &mesh)
Null constructor.
wordList toc() const
Return the table of contents.
Definition: dictionary.C:776
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Helper container for force Su and Sp terms.
Definition: forceSuSp.H:61
points setSize(newPointi)
A class for handling words, derived from string.
Definition: word.H:59
virtual forceSuSp calcCoupled(const typename CloudType::parcelType &p, const scalar dt, const scalar mass, const scalar Re, const scalar muc) const
Calculate the coupled force.
static const zero Zero
Definition: zero.H:91
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:218
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
messageStream Info
virtual ~ParticleForceList()
Destructor.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
List of particle forces.