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 | 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 "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_(NULL),
45  massStickPtr_(NULL)
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_(NULL),
98  massStickPtr_(NULL)
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,
131  dimensionedScalar("zero", dimMass, 0.0)
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,
160  dimensionedScalar("zero", dimMass, 0.0)
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  const scalar trackFraction,
176  const tetIndices& tetIs
177 )
178 {
179  label patchi = patchData_.applyToPatch(pp.index());
180 
181  if (patchi >= 0)
182  {
183  vector& U = p.U();
184  bool& active = p.active();
185 
187  this->wordToInteractionType
188  (
189  patchData_[patchi].interactionTypeName()
190  );
191 
192  switch (it)
193  {
195  {
196  scalar dm = p.mass()*p.nParticle();
197 
198  keepParticle = false;
199  active = false;
200  U = Zero;
201  nEscape_[patchi]++;
202  massEscape_[patchi] += dm;
203  if (writeFields_)
204  {
205  label pI = pp.index();
206  label fI = pp.whichFace(p.face());
207  massEscape().boundaryFieldRef()[pI][fI] += dm;
208  }
209  break;
210  }
212  {
213  scalar dm = p.mass()*p.nParticle();
214 
215  keepParticle = true;
216  active = false;
217  U = Zero;
218  nStick_[patchi]++;
219  massStick_[patchi] += dm;
220  if (writeFields_)
221  {
222  label pI = pp.index();
223  label fI = pp.whichFace(p.face());
224  massStick().boundaryFieldRef()[pI][fI] += dm;
225  }
226  break;
227  }
229  {
230  keepParticle = true;
231  active = true;
232 
233  vector nw;
234  vector Up;
235 
236  this->owner().patchData(p, pp, trackFraction, tetIs, nw, Up);
237 
238  // Calculate motion relative to patch velocity
239  U -= Up;
240 
241  scalar Un = U & nw;
242  vector Ut = U - Un*nw;
243 
244  if (Un > 0)
245  {
246  U -= (1.0 + patchData_[patchi].e())*Un*nw;
247  }
248 
249  U -= patchData_[patchi].mu()*Ut;
250 
251  // Return velocity to global space
252  U += Up;
253 
254  break;
255  }
256  default:
257  {
259  << "Unknown interaction type "
260  << patchData_[patchi].interactionTypeName()
261  << "(" << it << ") for patch "
262  << patchData_[patchi].patchName()
263  << ". Valid selections are:" << this->interactionTypeNames_
264  << endl << abort(FatalError);
265  }
266  }
267 
268  return true;
269  }
270 
271  return false;
272 }
273 
274 
275 template<class CloudType>
277 {
278  // retrieve any stored data
279  labelList npe0(patchData_.size(), 0);
280  this->getModelProperty("nEscape", npe0);
281 
282  scalarList mpe0(patchData_.size(), 0.0);
283  this->getModelProperty("massEscape", mpe0);
284 
285  labelList nps0(patchData_.size(), 0);
286  this->getModelProperty("nStick", nps0);
287 
288  scalarList mps0(patchData_.size(), 0.0);
289  this->getModelProperty("massStick", mps0);
290 
291  // accumulate current data
292  labelList npe(nEscape_);
293  Pstream::listCombineGather(npe, plusEqOp<label>());
294  npe = npe + npe0;
295 
296  scalarList mpe(massEscape_);
297  Pstream::listCombineGather(mpe, plusEqOp<scalar>());
298  mpe = mpe + mpe0;
299 
300  labelList nps(nStick_);
301  Pstream::listCombineGather(nps, plusEqOp<label>());
302  nps = nps + nps0;
303 
304  scalarList mps(massStick_);
305  Pstream::listCombineGather(mps, plusEqOp<scalar>());
306  mps = mps + mps0;
307 
308 
309  forAll(patchData_, i)
310  {
311  os << " Parcel fate (number, mass) : patch "
312  << patchData_[i].patchName() << nl
313  << " - escape = " << npe[i]
314  << ", " << mpe[i] << nl
315  << " - stick = " << nps[i]
316  << ", " << mps[i] << nl;
317  }
318 
319  if (this->writeTime())
320  {
321  this->setModelProperty("nEscape", npe);
322  nEscape_ = 0;
323 
324  this->setModelProperty("massEscape", mpe);
325  massEscape_ = 0.0;
326 
327  this->setModelProperty("nStick", nps);
328  nStick_ = 0;
329 
330  this->setModelProperty("massStick", mps);
331  massStick_ = 0.0;
332  }
333 }
334 
335 
336 // ************************************************************************* //
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:428
U
Definition: pEqn.H:83
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
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:377
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:137
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:253
virtual ~LocalInteraction()
Destructor.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
Templated patch interaction model class.
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
virtual bool correct(typename CloudType::parcelType &p, const polyPatch &pp, bool &keepParticle, const scalar trackFraction, const tetIndices &tetIs)
Apply velocity correction.
LocalInteraction(const dictionary &dict, CloudType &owner)
Construct from dictionary.
static const zero Zero
Definition: zero.H:91
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
errorManip< error > abort(error &err)
Definition: errorManip.H:131
ParcelType parcelType
Type of parcel the cloud was instantiated for.
Definition: DSMCCloud.H:217
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
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
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
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:91
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
label index() const
Return the index of this patch in the boundaryMesh.
label nw
Definition: createFields.H:25
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243