LocalInteraction.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 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 "LocalInteraction.H"
27 #include "wordAndDictionary.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 
35 :
36  public Tuple2<wordRe, dictionary>
37 {
39 
41 
43 };
44 
45 
47 {
48  wd.first() = wordRe(is);
49  dictionary d(is);
50  wd.second().transfer(d);
51  return is;
52 }
53 
54 
56 {
57  return os << wd.first() << token::SPACE << wd.second();
58 }
59 
60 
62 {}
63 
64 
66 {
67  is >> *this;
68 }
69 
70 }
71 
72 
73 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
74 
75 template<class CloudType>
77 (
78  const dictionary& dict,
80 )
81 :
83  patchInteractionTypes_
84  (
85  this->owner().mesh().boundaryMesh().size(),
87  ),
88  patchEs_(this->owner().mesh().boundaryMesh().size(), NaN),
89  patchMus_(this->owner().mesh().boundaryMesh().size(), NaN),
90  nEscape_(this->owner().mesh().boundaryMesh().size(), 0),
91  massEscape_(this->owner().mesh().boundaryMesh().size(), scalar(0)),
92  nStick_(this->owner().mesh().boundaryMesh().size(), 0),
93  massStick_(this->owner().mesh().boundaryMesh().size(), scalar(0)),
94  writeFields_(this->coeffDict().lookupOrDefault("writeFields", false)),
95  massEscapePtr_(nullptr),
96  massStickPtr_(nullptr)
97 {
98  const polyBoundaryMesh& patches = this->owner().mesh().boundaryMesh();
99 
100  if (writeFields_)
101  {
102  const word massEscapeName(this->owner().name() + ":massEscape");
103  const word massStickName(this->owner().name() + ":massStick");
104 
105  Info<< " Interaction fields will be written to " << massEscapeName
106  << " and " << massStickName << endl;
107 
108  (void)massEscape();
109  (void)massStick();
110  }
111  else
112  {
113  Info<< " Interaction fields will not be written" << endl;
114  }
115 
116  // Get the patch-settings dictionaries
117  dictionary patchesDict;
118  if (this->coeffDict().isDict("patches"))
119  {
120  patchesDict = this->coeffDict().subDict("patches");
121  }
122  else
123  {
124  const List<wordReAndDictionary> patchNameAndDicts
125  (
126  this->coeffDict().lookup("patches")
127  );
128 
129  forAll(patchNameAndDicts, dicti)
130  {
131  patchesDict.set
132  (
133  keyType(string(patchNameAndDicts[dicti].first())),
134  patchNameAndDicts[dicti].second()
135  );
136  }
137  }
138 
139  // Read the patch settings
140  wordList unspecifiedNonConstraintPatches;
142  {
143  const word& patchName = patches[patchi].name();
144 
145  const bool havePatchDict = patchesDict.found(patchName);
146 
147  const bool patchIsConstraint =
149 
150  // No settings for constrained patch. No model.
151  if (!havePatchDict && patchIsConstraint)
152  {
153  patchInteractionTypes_[patchi] =
155  continue;
156  }
157 
158  // No settings for non-constrained patch. Error.
159  if (!havePatchDict && !patchIsConstraint)
160  {
161  unspecifiedNonConstraintPatches.append(patches[patchi].name());
162  continue;
163  }
164 
165  const dictionary& patchDict = patchesDict.subDict(patchName);
166 
167  // Settings for constrained patch. Ignored unless "patchType" is
168  // correctly specified.
169  if (havePatchDict && patchIsConstraint)
170  {
171  if (!patchDict.found("patchType"))
172  {
173  patchInteractionTypes_[patchi] =
175  continue;
176  }
177 
178  const word patchType = patchDict.lookup<word>("patchType");
179 
180  if (patchType != patches[patchi].type())
181  {
183  << "Type " << patchType
184  << " specified for patch " << patchName
185  << " does not match the patch type "
186  << patches[patchi].type() << exit(FatalError);
187  }
188  }
189 
190  // Read and set the interaction model
191  const word itName = patchDict.lookup<word>("type");
192  const interactionType it = this->wordToInteractionType(itName);
193 
195  {
197  << "Unknown patch interaction type "
198  << itName << " for patch " << patchName
199  << ". Valid types are:"
201  << nl << exit(FatalError);
202  }
203 
204  patchInteractionTypes_[patchi] = it;
205 
207  {
208  patchEs_[patchi] = patchDict.lookupOrDefault<scalar>("e", 1);
209  patchMus_[patchi] = patchDict.lookupOrDefault<scalar>("mu", 0);
210  }
211  }
212 
213  // Error if interactions are unspecified for non-constraint patches
214  if (!unspecifiedNonConstraintPatches.empty())
215  {
217  << "No interaction type was specified for non-constraint patches: "
218  << unspecifiedNonConstraintPatches
219  << exit(FatalError);
220  }
221 }
222 
223 
224 template<class CloudType>
226 (
227  const LocalInteraction<CloudType>& pim
228 )
229 :
231  patchInteractionTypes_(pim.patchInteractionTypes_),
232  patchEs_(pim.patchEs_),
233  patchMus_(pim.patchMus_),
234  nEscape_(pim.nEscape_),
235  massEscape_(pim.massEscape_),
236  nStick_(pim.nStick_),
237  massStick_(pim.massStick_),
238  writeFields_(pim.writeFields_),
239  massEscapePtr_(nullptr),
240  massStickPtr_(nullptr)
241 {}
242 
243 
244 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
245 
246 template<class CloudType>
248 {}
249 
250 
251 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
252 
253 template<class CloudType>
255 {
256  if (!massEscapePtr_.valid())
257  {
258  const fvMesh& mesh = this->owner().mesh();
259 
260  massEscapePtr_.reset
261  (
262  new volScalarField
263  (
264  IOobject
265  (
266  this->owner().name() + ":massEscape",
267  mesh.time().name(),
268  mesh,
271  ),
272  mesh,
274  )
275  );
276  }
277 
278  return massEscapePtr_();
279 }
280 
281 
282 template<class CloudType>
284 {
285  if (!massStickPtr_.valid())
286  {
287  const fvMesh& mesh = this->owner().mesh();
288 
289  massStickPtr_.reset
290  (
291  new volScalarField
292  (
293  IOobject
294  (
295  this->owner().name() + ":massStick",
296  mesh.time().name(),
297  mesh,
300  ),
301  mesh,
303  )
304  );
305  }
306 
307  return massStickPtr_();
308 }
309 
310 
311 template<class CloudType>
313 (
314  typename CloudType::parcelType& p,
315  const polyPatch& pp,
316  bool& keepParticle
317 )
318 {
319  const label patchi = pp.index();
320 
321  if (isA<processorPolyPatch>(pp))
322  {
323  return false;
324  }
325 
326  switch (patchInteractionTypes_[patchi])
327  {
329  {
330  return false;
331  }
332 
334  {
335  const scalar dm = p.mass()*p.nParticle();
336 
337  keepParticle = false;
338  p.moving() = false;
339  p.U() = Zero;
340  nEscape_[patchi] ++;
341  massEscape_[patchi] += dm;
342 
343  if (writeFields_)
344  {
345  const label patchFacei = pp.whichFace(p.face());
346  massEscape().boundaryFieldRef()[patchi][patchFacei] += dm;
347  }
348 
349  return true;
350  }
351 
353  {
354  const scalar dm = p.mass()*p.nParticle();
355 
356  keepParticle = true;
357  p.moving() = false;
358  p.U() = Zero;
359  nStick_[patchi] ++;
360  massStick_[patchi] += dm;
361 
362  if (writeFields_)
363  {
364  const label patchFacei = pp.whichFace(p.face());
365  massStick().boundaryFieldRef()[patchi][patchFacei] += dm;
366  }
367 
368  return true;
369  }
370 
372  {
373  keepParticle = true;
374  p.moving() = true;
375 
376  vector nw, Up;
377  this->owner().patchData(p, pp, nw, Up);
378 
379  // Make motion relative to patch velocity
380  p.U() -= Up;
381 
382  const scalar Un = p.U() & nw;
383  const vector Ut = p.U() - Un*nw;
384 
385  if (Un > 0)
386  {
387  p.U() -= (1 + patchEs_[patchi])*Un*nw;
388  }
389 
390  p.U() -= patchMus_[patchi]*Ut;
391 
392  // Return velocity to global space
393  p.U() += Up;
394 
395  return true;
396  }
397 
398  default:
399  {
400  return false;
401  }
402  }
403 
404  return false;
405 }
406 
407 
408 template<class CloudType>
410 {
411  const polyBoundaryMesh& patches = this->owner().mesh().boundaryMesh();
412 
413  // Determine the number of non-processor patches
415  for (; isA<processorPolyPatch>(patches[nPatches - 1]); nPatches --);
416 
417  // Retrieve any stored data
418  labelList npe0(nPatches, 0);
419  this->getModelProperty("nEscape", npe0);
420 
421  scalarList mpe0(nPatches, scalar(0));
422  this->getModelProperty("massEscape", mpe0);
423 
424  labelList nps0(nPatches, 0);
425  this->getModelProperty("nStick", nps0);
426 
427  scalarList mps0(nPatches, scalar(0));
428  this->getModelProperty("massStick", mps0);
429 
430  // Accumulate current data
431  labelList npe(SubList<label>(nEscape_, nPatches));
433  npe = npe + npe0;
434 
435  scalarList mpe(SubList<scalar>(massEscape_, nPatches));
437  mpe = mpe + mpe0;
438 
439  labelList nps(SubList<label>(nStick_, nPatches));
441  nps = nps + nps0;
442 
443  scalarList mps(SubList<scalar>(massStick_, nPatches));
445  mps = mps + mps0;
446 
447  for (label patchi = 0; patchi < nPatches; ++ patchi)
448  {
449  if
450  (
451  patchInteractionTypes_[patchi]
453  )
454  {
455  os << " Parcel fate (number, mass) : patch "
456  << this->owner().mesh().boundaryMesh()[patchi].name() << nl
457  << " - escape = " << npe[patchi]
458  << ", " << mpe[patchi] << nl
459  << " - stick = " << nps[patchi]
460  << ", " << mps[patchi] << nl;
461  }
462  }
463 
464  if (this->writeTime())
465  {
466  this->setModelProperty("nEscape", npe);
467  nEscape_ = 0;
468 
469  this->setModelProperty("massEscape", mpe);
470  massEscape_ = scalar(0);
471 
472  this->setModelProperty("nStick", nps);
473  nStick_ = 0;
474 
475  this->setModelProperty("massStick", mps);
476  massStick_ = scalar(0);
477  }
478 }
479 
480 
481 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:225
const fvMesh & mesh() const
Return references to the mesh.
Definition: DSMCCloudI.H:41
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:294
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
Patch interaction specified on a patch-by-patch basis.
virtual void info(Ostream &os)
Write patch interaction info to stream.
LocalInteraction(const dictionary &dict, CloudType &owner)
Construct from dictionary.
volScalarField & massStick()
Return access to the massStick field.
volScalarField & massEscape()
Return access to the massEscape field.
virtual ~LocalInteraction()
Destructor.
virtual bool correct(typename CloudType::parcelType &p, const polyPatch &pp, bool &keepParticle)
Apply velocity correction.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Templated patch interaction model class.
static interactionType wordToInteractionType(const word &itWord)
Convert word to interaction result.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
A List obtained as a section of another List.
Definition: SubList.H:56
A 2-tuple for storing two objects of different types.
Definition: Tuple2.H:66
const Type2 & second() const
Return second.
Definition: Tuple2.H:131
const Type1 & first() const
Return first.
Definition: Tuple2.H:119
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
void transfer(dictionary &)
Transfer the contents of the argument and annul the argument.
Definition: dictionary.C:1353
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:1152
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:797
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
const polyMesh & mesh() const
Return reference to polyMesh.
Definition: fvMesh.H:441
A class for handling keywords in dictionaries.
Definition: keyType.H:69
label index() const
Return the index of this patch in the boundaryMesh.
Foam::polyBoundaryMesh.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:360
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:238
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:128
A wordRe is a word, but can also have a regular expression for matching words.
Definition: wordRe.H:77
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
const fvPatchList & patches
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
labelList second(const UList< labelPair > &p)
Definition: patchToPatch.C:49
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
Istream & operator>>(Istream &, pistonPointEdgeData &)
const dimensionSet dimMass
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
error FatalError
static const char nl
Definition: Ostream.H:266
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
label nPatches
Definition: readKivaGrid.H:396
dictionary dict
volScalarField & p