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-2018 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 
28 // * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
29 
30 template<class CloudType>
32 (
33  const dictionary& dict,
34  CloudType& cloud
35 )
36 :
37  PatchInteractionModel<CloudType>(dict, cloud, typeName),
38  patchData_(cloud.mesh(), this->coeffDict()),
39  nEscape_(patchData_.size(), 0),
40  massEscape_(patchData_.size(), 0.0),
41  nStick_(patchData_.size(), 0),
42  massStick_(patchData_.size(), 0.0),
43  writeFields_(this->coeffDict().lookupOrDefault("writeFields", false)),
44  massEscapePtr_(nullptr),
45  massStickPtr_(nullptr)
46 {
47  if (writeFields_)
48  {
49  word massEscapeName(this->owner().name() + ":massEscape");
50  word massStickName(this->owner().name() + ":massStick");
51  Info<< " Interaction fields will be written to " << massEscapeName
52  << " and " << massStickName << endl;
53 
54  (void)massEscape();
55  (void)massStick();
56  }
57  else
58  {
59  Info<< " Interaction fields will not be written" << endl;
60  }
61 
62  // check that interactions are valid/specified
63  forAll(patchData_, patchi)
64  {
65  const word& interactionTypeName =
66  patchData_[patchi].interactionTypeName();
68  this->wordToInteractionType(interactionTypeName);
69 
71  {
72  const word& patchName = patchData_[patchi].patchName();
74  << "Unknown patch interaction type "
75  << interactionTypeName << " for patch " << patchName
76  << ". Valid selections are:"
78  << nl << exit(FatalError);
79  }
80  }
81 }
82 
83 
84 template<class CloudType>
86 (
88 )
89 :
91  patchData_(pim.patchData_),
92  nEscape_(pim.nEscape_),
93  massEscape_(pim.massEscape_),
94  nStick_(pim.nStick_),
95  massStick_(pim.massStick_),
96  writeFields_(pim.writeFields_),
97  massEscapePtr_(nullptr),
98  massStickPtr_(nullptr)
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
103 
104 template<class CloudType>
106 {}
107 
108 
109 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
110 
111 template<class CloudType>
113 {
114  if (!massEscapePtr_.valid())
115  {
116  const fvMesh& mesh = this->owner().mesh();
117 
118  massEscapePtr_.reset
119  (
120  new volScalarField
121  (
122  IOobject
123  (
124  this->owner().name() + ":massEscape",
125  mesh.time().timeName(),
126  mesh,
127  IOobject::READ_IF_PRESENT,
128  IOobject::AUTO_WRITE
129  ),
130  mesh,
132  )
133  );
134  }
135 
136  return massEscapePtr_();
137 }
138 
139 
140 template<class CloudType>
142 {
143  if (!massStickPtr_.valid())
144  {
145  const fvMesh& mesh = this->owner().mesh();
146 
147  massStickPtr_.reset
148  (
149  new volScalarField
150  (
151  IOobject
152  (
153  this->owner().name() + ":massStick",
154  mesh.time().timeName(),
155  mesh,
156  IOobject::READ_IF_PRESENT,
157  IOobject::AUTO_WRITE
158  ),
159  mesh,
161  )
162  );
163  }
164 
165  return massStickPtr_();
166 }
167 
168 
169 template<class CloudType>
171 (
172  typename CloudType::parcelType& p,
173  const polyPatch& pp,
174  bool& keepParticle
175 )
176 {
177  label patchi = patchData_.applyToPatch(pp.index());
178 
179  if (patchi >= 0)
180  {
181  vector& U = p.U();
182  bool& active = p.active();
183 
185  this->wordToInteractionType
186  (
187  patchData_[patchi].interactionTypeName()
188  );
189 
190  switch (it)
191  {
193  {
194  return false;
195  }
197  {
198  scalar dm = p.mass()*p.nParticle();
199 
200  keepParticle = false;
201  active = false;
202  U = Zero;
203  nEscape_[patchi]++;
204  massEscape_[patchi] += dm;
205  if (writeFields_)
206  {
207  label pI = pp.index();
208  label fI = pp.whichFace(p.face());
209  massEscape().boundaryFieldRef()[pI][fI] += dm;
210  }
211  break;
212  }
214  {
215  scalar dm = p.mass()*p.nParticle();
216 
217  keepParticle = true;
218  active = false;
219  U = Zero;
220  nStick_[patchi]++;
221  massStick_[patchi] += dm;
222  if (writeFields_)
223  {
224  label pI = pp.index();
225  label fI = pp.whichFace(p.face());
226  massStick().boundaryFieldRef()[pI][fI] += dm;
227  }
228  break;
229  }
231  {
232  keepParticle = true;
233  active = true;
234 
235  vector nw;
236  vector Up;
237 
238  this->owner().patchData(p, pp, nw, Up);
239 
240  // Calculate motion relative to patch velocity
241  U -= Up;
242 
243  scalar Un = U & nw;
244  vector Ut = U - Un*nw;
245 
246  if (Un > 0)
247  {
248  U -= (1.0 + patchData_[patchi].e())*Un*nw;
249  }
250 
251  U -= patchData_[patchi].mu()*Ut;
252 
253  // Return velocity to global space
254  U += Up;
255 
256  break;
257  }
258  default:
259  {
261  << "Unknown interaction type "
262  << patchData_[patchi].interactionTypeName()
263  << "(" << it << ") for patch "
264  << patchData_[patchi].patchName()
265  << ". Valid selections are:" << this->interactionTypeNames_
266  << endl << abort(FatalError);
267  }
268  }
269 
270  return true;
271  }
272 
273  return false;
274 }
275 
276 
277 template<class CloudType>
279 {
280  // retrieve any stored data
281  labelList npe0(patchData_.size(), 0);
282  this->getModelProperty("nEscape", npe0);
283 
284  scalarList mpe0(patchData_.size(), 0.0);
285  this->getModelProperty("massEscape", mpe0);
286 
287  labelList nps0(patchData_.size(), 0);
288  this->getModelProperty("nStick", nps0);
289 
290  scalarList mps0(patchData_.size(), 0.0);
291  this->getModelProperty("massStick", mps0);
292 
293  // accumulate current data
294  labelList npe(nEscape_);
295  Pstream::listCombineGather(npe, plusEqOp<label>());
296  npe = npe + npe0;
297 
298  scalarList mpe(massEscape_);
299  Pstream::listCombineGather(mpe, plusEqOp<scalar>());
300  mpe = mpe + mpe0;
301 
302  labelList nps(nStick_);
303  Pstream::listCombineGather(nps, plusEqOp<label>());
304  nps = nps + nps0;
305 
306  scalarList mps(massStick_);
307  Pstream::listCombineGather(mps, plusEqOp<scalar>());
308  mps = mps + mps0;
309 
310 
311  forAll(patchData_, i)
312  {
313  os << " Parcel fate (number, mass) : patch "
314  << patchData_[i].patchName() << nl
315  << " - escape = " << npe[i]
316  << ", " << mpe[i] << nl
317  << " - stick = " << nps[i]
318  << ", " << mps[i] << nl;
319  }
320 
321  if (this->writeTime())
322  {
323  this->setModelProperty("nEscape", npe);
324  nEscape_ = 0;
325 
326  this->setModelProperty("massEscape", mpe);
327  massEscape_ = 0.0;
328 
329  this->setModelProperty("nStick", nps);
330  nStick_ = 0;
331 
332  this->setModelProperty("massStick", mps);
333  massStick_ = 0.0;
334  }
335 }
336 
337 
338 // ************************************************************************* //
Patch interaction specified on a patch-by-patch basis.
volScalarField & massEscape()
Return access to the massEscape field.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
virtual void info(Ostream &os)
Write patch interaction info to stream.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual ~LocalInteraction()
Destructor.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:622
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
Templated patch interaction model class.
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
LocalInteraction(const dictionary &dict, CloudType &owner)
Construct from dictionary.
static const zero Zero
Definition: zero.H:97
errorManip< error > abort(error &err)
Definition: errorManip.H:131
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:215
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
volScalarField & massStick()
Return access to the massStick field.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label patchi
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label index() const
Return the index of this patch in the boundaryMesh.
messageStream Info
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:69
virtual bool correct(typename CloudType::parcelType &p, const polyPatch &pp, bool &keepParticle)
Apply velocity correction.
label nw
Definition: createFields.H:25
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:383