StandardWallInteraction.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) 2015-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
30 
31 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
32 
33 template<class CloudType>
35 {
37 
38  forAll(nEscape_, patchi)
39  {
40  const word& patchName = mesh_.boundary()[patchi].name();
41 
42  forAll(nEscape_[patchi], injectori)
43  {
44  const word suffix = Foam::name(injectori);
45  this->writeTabbed(os, patchName + "_nEscape_" + suffix);
46  this->writeTabbed(os, patchName + "_massEscape_" + suffix);
47  this->writeTabbed(os, patchName + "_nStick_" + suffix);
48  this->writeTabbed(os, patchName + "_massStick_" + suffix);
49  }
50  }
51 }
52 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
56 template<class CloudType>
58 (
59  const dictionary& dict,
61 )
62 :
64  mesh_(cloud.mesh()),
65  interactionType_
66  (
67  this->wordToInteractionType(this->coeffDict().getWord("type"))
68  ),
69  e_(0.0),
70  mu_(0.0),
71  nEscape_(mesh_.boundaryMesh().nNonProcessor()),
72  massEscape_(nEscape_.size()),
73  nStick_(nEscape_.size()),
74  massStick_(nEscape_.size()),
75  injIdToIndex_()
76 {
77  const bool outputByInjectorId =
78  this->coeffDict().getOrDefault("outputByInjectorId", false);
79 
80  switch (interactionType_)
81  {
83  {
84  const word interactionTypeName(this->coeffDict().getWord("type"));
85 
87  << "Unknown interaction result type "
88  << interactionTypeName
89  << ". Valid selections are:" << this->interactionTypeNames_
90  << endl << exit(FatalError);
91 
92  break;
93  }
95  {
96  e_ = this->coeffDict().getOrDefault("e", 1.0);
97  mu_ = this->coeffDict().getOrDefault("mu", 0.0);
98  break;
99  }
100  default:
101  {}
102  }
103 
104  // Determine the number of injectors and the injector mapping
105  label nInjectors = 0;
106  if (outputByInjectorId)
107  {
108  for (const auto& inj : cloud.injectors())
109  {
110  injIdToIndex_.insert(inj.injectorID(), nInjectors++);
111  }
112  }
113 
114  // The normal case, and safety if injector mapping was somehow null.
115  if (injIdToIndex_.empty())
116  {
117  nInjectors = 1;
118  }
119 
120  forAll(nEscape_, patchi)
121  {
122  nEscape_[patchi].setSize(nInjectors, Zero);
123  massEscape_[patchi].setSize(nInjectors, Zero);
124  nStick_[patchi].setSize(nInjectors, Zero);
125  massStick_[patchi].setSize(nInjectors, Zero);
126  }
127 }
128 
129 
130 template<class CloudType>
132 (
133  const StandardWallInteraction<CloudType>& pim
134 )
135 :
136  PatchInteractionModel<CloudType>(pim),
137  mesh_(pim.mesh_),
138  interactionType_(pim.interactionType_),
139  e_(pim.e_),
140  mu_(pim.mu_),
141  nEscape_(pim.nEscape_),
142  massEscape_(pim.massEscape_),
143  nStick_(pim.nStick_),
144  massStick_(pim.massStick_),
145  injIdToIndex_(pim.injIdToIndex_)
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
151 template<class CloudType>
153 (
154  typename CloudType::parcelType& p,
155  const polyPatch& pp,
156  bool& keepParticle
157 )
158 {
159  vector& U = p.U();
160 
161  if (isA<wallPolyPatch>(pp))
162  {
163  // Location for storing the stats.
164  const label idx =
165  (
166  injIdToIndex_.size()
167  ? injIdToIndex_.lookup(p.typeId(), 0)
168  : 0
169  );
170 
171  switch (interactionType_)
172  {
174  {
175  return false;
176  }
178  {
179  keepParticle = false;
180  p.active(false);
181  U = Zero;
182 
183  const scalar dm = p.nParticle()*p.mass();
184 
185  nEscape_[pp.index()][idx]++;
186  massEscape_[pp.index()][idx] += dm;
187  break;
188  }
190  {
191  keepParticle = true;
192  p.active(false);
193  U = Zero;
194 
195  const scalar dm = p.nParticle()*p.mass();
196 
197  nStick_[pp.index()][idx]++;
198  massStick_[pp.index()][idx] += dm;
199  break;
200  }
202  {
203  keepParticle = true;
204  p.active(true);
205 
206  vector nw;
207  vector Up;
208 
209  this->owner().patchData(p, pp, nw, Up);
210 
211  // Calculate motion relative to patch velocity
212  U -= Up;
213 
214  if (mag(Up) > 0 && mag(U) < this->Urmax())
215  {
217  << "Particle U the same as patch "
218  << " The particle has been removed" << nl << endl;
219 
220  keepParticle = false;
221  p.active(false);
222  U = Zero;
223  break;
224  }
225 
226  scalar Un = U & nw;
227  vector Ut = U - Un*nw;
228 
229  if (Un > 0)
230  {
231  U -= (1.0 + e_)*Un*nw;
232  }
233 
234  U -= mu_*Ut;
235 
236  // Return velocity to global space
237  U += Up;
238 
239 
240  break;
241  }
242  default:
243  {
245  << "Unknown interaction type "
246  << this->interactionTypeToWord(interactionType_)
247  << "(" << interactionType_ << ")" << endl
248  << abort(FatalError);
249  }
250  }
251 
252  return true;
253  }
254 
255  return false;
256 }
257 
258 
259 template<class CloudType>
261 {
263 
264  labelListList npe0(nEscape_.size());
265  scalarListList mpe0(nEscape_.size());
266  labelListList nps0(nEscape_.size());
267  scalarListList mps0(nEscape_.size());
268 
269  forAll(nEscape_, patchi)
270  {
271  label lsd = nEscape_[patchi].size();
272  npe0[patchi].setSize(lsd, Zero);
273  mpe0[patchi].setSize(lsd, Zero);
274  nps0[patchi].setSize(lsd, Zero);
275  mps0[patchi].setSize(lsd, Zero);
276  }
277 
278  this->getModelProperty("nEscape", npe0);
279  this->getModelProperty("massEscape", mpe0);
280  this->getModelProperty("nStick", nps0);
281  this->getModelProperty("massStick", mps0);
282 
283  // Accumulate current data
284  labelListList npe(nEscape_);
285 
286  forAll(npe, i)
287  {
288  Pstream::listCombineGather(npe[i], plusEqOp<label>());
289  npe[i] = npe[i] + npe0[i];
290  }
291 
292  scalarListList mpe(massEscape_);
293  forAll(mpe, i)
294  {
295  Pstream::listCombineGather(mpe[i], plusEqOp<scalar>());
296  mpe[i] = mpe[i] + mpe0[i];
297  }
298 
299  labelListList nps(nStick_);
300  forAll(nps, i)
301  {
302  Pstream::listCombineGather(nps[i], plusEqOp<label>());
303  nps[i] = nps[i] + nps0[i];
304  }
305 
306  scalarListList mps(massStick_);
307  forAll(nps, i)
308  {
309  Pstream::listCombineGather(mps[i], plusEqOp<scalar>());
310  mps[i] = mps[i] + mps0[i];
311  }
312 
313  if (injIdToIndex_.size())
314  {
315  // Since injIdToIndex_ is a one-to-one mapping (starting as zero),
316  // can simply invert it.
317  labelList indexToInjector(injIdToIndex_.size());
318  forAllConstIters(injIdToIndex_, iter)
319  {
320  indexToInjector[iter.val()] = iter.key();
321  }
322 
323  forAll(npe, patchi)
324  {
325  forAll(mpe[patchi], indexi)
326  {
327  const word& patchName = mesh_.boundary()[patchi].name() ;
328 
329  Log_<< " Parcel fate: patch " << patchName
330  << " (number, mass)" << nl
331  << " - escape (injector " << indexToInjector[indexi]
332  << ") = " << npe[patchi][indexi]
333  << ", " << mpe[patchi][indexi] << nl
334  << " - stick (injector " << indexToInjector[indexi]
335  << ") = " << nps[patchi][indexi]
336  << ", " << mps[patchi][indexi] << nl;
337 
338  this->file()
339  << tab << npe[patchi][indexi] << tab << mpe[patchi][indexi]
340  << tab << nps[patchi][indexi] << tab << mps[patchi][indexi];
341  }
342  }
343 
344  this->file() << endl;
345  }
346  else
347  {
348  forAll(npe, patchi)
349  {
350  const word& patchName = mesh_.boundary()[patchi].name();
351 
352  Log_<< " Parcel fate: patch (number, mass) "
353  << patchName << nl
354  << " - escape = "
355  << npe[patchi][0] << ", " << mpe[patchi][0] << nl
356  << " - stick = "
357  << nps[patchi][0] << ", " << mps[patchi][0] << nl;
358 
359  this->file()
360  << tab << npe[patchi][0] << tab << mpe[patchi][0]
361  << tab << nps[patchi][0] << tab << mps[patchi][0];
362  }
363 
364  this->file() << endl;
365  }
366 
367  if (this->writeTime())
368  {
369  this->setModelProperty("nEscape", npe);
370  this->setModelProperty("massEscape", mpe);
371  this->setModelProperty("nStick", nps);
372  this->setModelProperty("massStick", mps);
373 
374  nEscape_ = Zero;
375  massEscape_ = Zero;
376  nStick_ = Zero;
377  massStick_ = Zero;
378  }
379 }
380 
381 
382 // ************************************************************************* //
dictionary dict
DSMCCloud< dsmcParcel > CloudType
virtual void info()
Write patch interaction info.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
scalar mu_
The unity minus the restitution coefficient.
virtual void writeFileHeader(Ostream &os)
Output file header information.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
List< List< label > > nStick_
Number of parcels stuck to patches.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
virtual bool correct(typename CloudType::parcelType &p, const polyPatch &pp, bool &keepParticle)
Apply velocity correction.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
List< List< scalar > > massStick_
Mass of parcels stuck to patches.
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
virtual void info()
Write patch interaction info.
Templated patch interaction model class.
PatchInteractionModel< CloudType >::interactionType interactionType_
Interaction type.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:152
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
#define Log_
Report write to Foam::Info if the class log switch is true.
List< scalarList > scalarListList
List of scalarList.
Definition: scalarList.H:35
A class for handling words, derived from Foam::string.
Definition: word.H:63
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
Vector< scalar > vector
Definition: vector.H:57
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:122
errorManip< error > abort(error &err)
Definition: errorManip.H:139
scalar e_
Elasticity coefficient.
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:290
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
bool empty() const noexcept
True if the hash table is empty.
Definition: HashTable.H:337
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
List< List< scalar > > massEscape_
Mass of parcels escaped.
Map< label > injIdToIndex_
InjectorId to index map, when outputting escaped/stick/...
U
Definition: pEqn.H:72
#define WarningInFunction
Report a warning using Foam::Warning.
StandardWallInteraction(const dictionary &dict, CloudType &cloud)
Construct from dictionary.
List< List< label > > nEscape_
Number of parcels escaped.
static void listCombineGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127