LocalInteraction.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 
29 #include "LocalInteraction.H"
30 
31 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
32 
33 template<class CloudType>
35 {
37 
38  forAll(nEscape_, patchi)
39  {
40  const word& patchName = patchData_[patchi].patchName();
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 template<class CloudType>
56 (
57  const dictionary& dict,
58  CloudType& cloud
59 )
60 :
61  PatchInteractionModel<CloudType>(dict, cloud, typeName),
62  patchData_(cloud.mesh(), this->coeffDict()),
63  nEscape_(patchData_.size()),
64  massEscape_(nEscape_.size()),
65  nStick_(nEscape_.size()),
66  massStick_(nEscape_.size()),
67  writeFields_(this->coeffDict().getOrDefault("writeFields", false)),
68  injIdToIndex_(),
69  massEscapePtr_(nullptr),
70  massStickPtr_(nullptr)
71 {
72  const bool outputByInjectorId
73  = this->coeffDict().getOrDefault("outputByInjectorId", false);
74 
75  if (writeFields_)
76  {
77  Info<< " Interaction fields will be written to "
78  << this->owner().name() << ":massEscape"
79  << " and "
80  << this->owner().name() << ":massStick" << endl;
81 
82  (void)massEscape();
83  (void)massStick();
84  }
85  else
86  {
87  Info<< " Interaction fields will not be written" << endl;
88  }
89 
90  // Determine the number of injectors and the injector mapping
91  label nInjectors = 0;
92  if (outputByInjectorId)
93  {
94  for (const auto& inj : cloud.injectors())
95  {
96  injIdToIndex_.insert(inj.injectorID(), nInjectors++);
97  }
98  }
99 
100  // The normal case, and safety if injector mapping was somehow null.
101  if (!nInjectors)
102  {
103  nInjectors = 1;
104  }
105 
106  // Check that interactions are valid/specified
107  forAll(patchData_, patchi)
108  {
109  const word& interactionTypeName =
110  patchData_[patchi].interactionTypeName();
112  this->wordToInteractionType(interactionTypeName);
113 
115  {
116  const word& patchName = patchData_[patchi].patchName();
118  << "Unknown patch interaction type "
119  << interactionTypeName << " for patch " << patchName
120  << ". Valid selections are:"
122  << nl << exit(FatalError);
123  }
124 
125  nEscape_[patchi].setSize(nInjectors, Zero);
126  massEscape_[patchi].setSize(nInjectors, Zero);
127  nStick_[patchi].setSize(nInjectors, Zero);
128  massStick_[patchi].setSize(nInjectors, Zero);
129  }
130 }
131 
132 
133 template<class CloudType>
135 (
136  const LocalInteraction<CloudType>& pim
137 )
138 :
139  PatchInteractionModel<CloudType>(pim),
140  patchData_(pim.patchData_),
141  nEscape_(pim.nEscape_),
142  massEscape_(pim.massEscape_),
143  nStick_(pim.nStick_),
144  massStick_(pim.massStick_),
145  writeFields_(pim.writeFields_),
146  injIdToIndex_(pim.injIdToIndex_),
147  massEscapePtr_(nullptr),
148  massStickPtr_(nullptr)
149 {}
150 
151 
152 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
153 
154 template<class CloudType>
156 {
157  if (!massEscapePtr_)
158  {
159  const fvMesh& mesh = this->owner().mesh();
160 
161  massEscapePtr_.reset
162  (
163  new volScalarField
164  (
165  IOobject
166  (
167  this->owner().name() + ":massEscape",
168  mesh.time().timeName(),
169  mesh,
172  ),
173  mesh,
175  )
176  );
177  }
178 
179  return *massEscapePtr_;
180 }
181 
182 
183 template<class CloudType>
185 {
186  if (!massStickPtr_)
187  {
188  const fvMesh& mesh = this->owner().mesh();
189 
190  massStickPtr_.reset
191  (
192  new volScalarField
193  (
194  IOobject
195  (
196  this->owner().name() + ":massStick",
197  mesh.time().timeName(),
198  mesh,
201  ),
202  mesh,
204  )
205  );
206  }
208  return *massStickPtr_;
209 }
210 
211 
212 template<class CloudType>
214 (
215  typename CloudType::parcelType& p,
216  const polyPatch& pp,
217  bool& keepParticle
218 )
219 {
220  const label patchi = patchData_.applyToPatch(pp.index());
221 
222  if (patchi >= 0)
223  {
224  vector& U = p.U();
225 
226  // Location for storing the stats.
227  const label idx =
228  (
229  injIdToIndex_.size()
230  ? injIdToIndex_.lookup(p.typeId(), 0)
231  : 0
232  );
233 
235  this->wordToInteractionType
236  (
237  patchData_[patchi].interactionTypeName()
238  );
239 
240  switch (it)
241  {
243  {
244  return false;
245  }
247  {
248  keepParticle = false;
249  p.active(false);
250  U = Zero;
251 
252  const scalar dm = p.mass()*p.nParticle();
253 
254  nEscape_[patchi][idx]++;
255  massEscape_[patchi][idx] += dm;
256 
257  if (writeFields_)
258  {
259  const label pI = pp.index();
260  const label fI = pp.whichFace(p.face());
261  massEscape().boundaryFieldRef()[pI][fI] += dm;
262  }
263  break;
264  }
266  {
267  keepParticle = true;
268  p.active(false);
269  U = Zero;
270 
271  const scalar dm = p.mass()*p.nParticle();
272 
273  nStick_[patchi][idx]++;
274  massStick_[patchi][idx] += dm;
275 
276  if (writeFields_)
277  {
278  const label pI = pp.index();
279  const label fI = pp.whichFace(p.face());
280  massStick().boundaryFieldRef()[pI][fI] += dm;
281  }
282  break;
283  }
285  {
286  keepParticle = true;
287  p.active(true);
288 
289  vector nw;
290  vector Up;
291 
292  this->owner().patchData(p, pp, nw, Up);
293 
294  // Calculate motion relative to patch velocity
295  U -= Up;
296 
297  if (mag(Up) > 0 && mag(U) < this->Urmax())
298  {
300  << "Particle U the same as patch "
301  << " The particle has been removed" << nl << endl;
302 
303  keepParticle = false;
304  p.active(false);
305  U = Zero;
306  break;
307  }
308 
309  scalar Un = U & nw;
310  vector Ut = U - Un*nw;
311 
312  if (Un > 0)
313  {
314  U -= (1.0 + patchData_[patchi].e())*Un*nw;
315  }
316 
317  U -= patchData_[patchi].mu()*Ut;
318 
319  // Return velocity to global space
320  U += Up;
321 
322  break;
323  }
324  default:
325  {
327  << "Unknown interaction type "
328  << patchData_[patchi].interactionTypeName()
329  << "(" << it << ") for patch "
330  << patchData_[patchi].patchName()
331  << ". Valid selections are:" << this->interactionTypeNames_
332  << endl << abort(FatalError);
333  }
334  }
335 
336  return true;
337  }
338 
339  return false;
340 }
341 
342 
343 template<class CloudType>
345 {
347 
348  // retrieve any stored data
349  labelListList npe0(patchData_.size());
350  scalarListList mpe0(patchData_.size());
351  labelListList nps0(patchData_.size());
352  scalarListList mps0(patchData_.size());
353 
354  forAll(patchData_, patchi)
355  {
356  label lsd = nEscape_[patchi].size();
357  npe0[patchi].setSize(lsd, Zero);
358  mpe0[patchi].setSize(lsd, Zero);
359  nps0[patchi].setSize(lsd, Zero);
360  mps0[patchi].setSize(lsd, Zero);
361  }
362 
363 
364  this->getModelProperty("nEscape", npe0);
365  this->getModelProperty("massEscape", mpe0);
366  this->getModelProperty("nStick", nps0);
367  this->getModelProperty("massStick", mps0);
368 
369  // accumulate current data
370  labelListList npe(nEscape_);
371  forAll(npe, i)
372  {
373  Pstream::listCombineGather(npe[i], plusEqOp<label>());
374  npe[i] = npe[i] + npe0[i];
375  }
376 
377  scalarListList mpe(massEscape_);
378  forAll(mpe, i)
379  {
380  Pstream::listCombineGather(mpe[i], plusEqOp<scalar>());
381  mpe[i] = mpe[i] + mpe0[i];
382  }
383 
384  labelListList nps(nStick_);
385  forAll(nps, i)
386  {
387  Pstream::listCombineGather(nps[i], plusEqOp<label>());
388  nps[i] = nps[i] + nps0[i];
389  }
390 
391  scalarListList mps(massStick_);
392  forAll(nps, i)
393  {
394  Pstream::listCombineGather(mps[i], plusEqOp<scalar>());
395  mps[i] = mps[i] + mps0[i];
396  }
397 
398  if (injIdToIndex_.size())
399  {
400  // Since injIdToIndex_ is a one-to-one mapping (starting at zero),
401  // can simply invert it.
402  labelList indexToInjector(injIdToIndex_.size());
403  forAllConstIters(injIdToIndex_, iter)
404  {
405  indexToInjector[iter.val()] = iter.key();
406  }
407 
408  forAll(patchData_, patchi)
409  {
410  forAll(mpe[patchi], indexi)
411  {
412  const word& patchName = patchData_[patchi].patchName();
413 
414  Log_<< " Parcel fate: patch " << patchName
415  << " (number, mass)" << nl
416  << " - escape (injector " << indexToInjector[indexi]
417  << " ) = " << npe[patchi][indexi]
418  << ", " << mpe[patchi][indexi] << nl
419  << " - stick (injector " << indexToInjector[indexi]
420  << " ) = " << nps[patchi][indexi]
421  << ", " << mps[patchi][indexi] << nl;
422  }
423  }
424  }
425  else
426  {
427  forAll(patchData_, patchi)
428  {
429  const word& patchName = patchData_[patchi].patchName();
430 
431  Log_<< " Parcel fate: patch " << patchName
432  << " (number, mass)" << nl
433  << " - escape = "
434  << npe[patchi][0] << ", " << mpe[patchi][0] << nl
435  << " - stick = "
436  << nps[patchi][0] << ", " << mps[patchi][0] << nl;
437  }
438  }
439 
440  forAll(npe, patchi)
441  {
442  forAll(npe[patchi], injectori)
443  {
444  this->file()
445  << tab << npe[patchi][injectori]
446  << tab << mpe[patchi][injectori]
447  << tab << nps[patchi][injectori]
448  << tab << mps[patchi][injectori];
449  }
450  }
451 
452  this->file() << endl;
453 
454  if (this->writeTime())
455  {
456  this->setModelProperty("nEscape", npe);
457  this->setModelProperty("massEscape", mpe);
458  this->setModelProperty("nStick", nps);
459  this->setModelProperty("massStick", mps);
460 
461  nEscape_ = Zero;
462  massEscape_ = Zero;
463  nStick_ = Zero;
464  massStick_ = Zero;
465  }
466 }
467 
468 
469 // ************************************************************************* //
volScalarField & massEscape()
Return access to the massEscape field.
dictionary dict
DSMCCloud< dsmcParcel > CloudType
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
constexpr char tab
The tab &#39;\t&#39; character(0x09)
Definition: Ostream.H:49
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
virtual void info()
Write patch interaction info.
Templated patch interaction model class.
const CloudType & owner() const
Return const access to the owner cloud.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition: HashTableI.H:152
virtual void info()
Write patch interaction info.
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
Reading is optional [identical to LAZY_READ].
LocalInteraction(const dictionary &dict, CloudType &owner)
Construct from dictionary.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:122
errorManip< error > abort(error &err)
Definition: errorManip.H:139
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
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
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
OBJstream os(runTime.globalPath()/outputName)
volScalarField & massStick()
Return access to the massStick field.
U
Definition: pEqn.H:72
#define WarningInFunction
Report a warning using Foam::Warning.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Automatically write from objectRegistry::writeObject()
virtual void writeFileHeader(Ostream &os)
Output file header information.
messageStream Info
Information stream (stdout output on master, null elsewhere)
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...
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
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
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
static interactionType wordToInteractionType(const word &itWord)
Convert word to interaction result.
virtual bool correct(typename CloudType::parcelType &p, const polyPatch &pp, bool &keepParticle)
Apply velocity correction.
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