ParticleErosion.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-2013 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 "ParticleErosion.H"
27 
28 // * * * * * * * * * * * * * Protectd Member Functions * * * * * * * * * * * //
29 
30 template<class CloudType>
32 (
33  const label globalPatchI
34 ) const
35 {
36  forAll(patchIDs_, i)
37  {
38  if (patchIDs_[i] == globalPatchI)
39  {
40  return i;
41  }
42  }
43 
44  return -1;
45 }
46 
47 
48 template<class CloudType>
50 {
51  if (QPtr_.valid())
52  {
53  QPtr_->write();
54  }
55  else
56  {
57  FatalErrorIn("void Foam::ParticleErosion<CloudType>::write()")
58  << "QPtr not valid" << abort(FatalError);
59  }
60 }
61 
62 
63 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
64 
65 template<class CloudType>
67 (
68  const dictionary& dict,
69  CloudType& owner,
70  const word& modelName
71 )
72 :
73  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
74  QPtr_(NULL),
75  patchIDs_(),
76  p_(readScalar(this->coeffDict().lookup("p"))),
77  psi_(this->coeffDict().template lookupOrDefault<scalar>("psi", 2.0)),
78  K_(this->coeffDict().template lookupOrDefault<scalar>("K", 2.0))
79 {
80  const wordList allPatchNames = owner.mesh().boundaryMesh().names();
81  wordList patchName(this->coeffDict().lookup("patches"));
82 
83  labelHashSet uniquePatchIDs;
84  forAllReverse(patchName, i)
85  {
86  labelList patchIDs = findStrings(patchName[i], allPatchNames);
87 
88  if (patchIDs.empty())
89  {
90  WarningIn
91  (
92  "Foam::ParticleErosion<CloudType>::ParticleErosion"
93  "("
94  "const dictionary&, "
95  "CloudType& "
96  ")"
97  ) << "Cannot find any patch names matching " << patchName[i]
98  << endl;
99  }
100 
101  uniquePatchIDs.insert(patchIDs);
102  }
103 
104  patchIDs_ = uniquePatchIDs.toc();
105 
106  // trigger ther creation of the Q field
107  preEvolve();
108 }
109 
110 
111 template<class CloudType>
113 (
115 )
116 :
118  QPtr_(NULL),
119  patchIDs_(pe.patchIDs_),
120  p_(pe.p_),
121  psi_(pe.psi_),
122  K_(pe.K_)
123 {}
124 
125 
126 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
127 
128 template<class CloudType>
130 {}
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
135 template<class CloudType>
137 {
138  if (QPtr_.valid())
139  {
140  QPtr_->internalField() = 0.0;
141  }
142  else
143  {
144  const fvMesh& mesh = this->owner().mesh();
145 
146  QPtr_.reset
147  (
148  new volScalarField
149  (
150  IOobject
151  (
152  this->owner().name() + "Q",
153  mesh.time().timeName(),
154  mesh,
155  IOobject::READ_IF_PRESENT,
156  IOobject::NO_WRITE
157  ),
158  mesh,
159  dimensionedScalar("zero", dimVolume, 0.0)
160  )
161  );
162  }
163 }
164 
165 
166 template<class CloudType>
168 (
169  const parcelType& p,
170  const polyPatch& pp,
171  const scalar trackFraction,
172  const tetIndices& tetIs,
173  bool&
174 )
175 {
176  const label patchI = pp.index();
177 
178  const label localPatchI = applyToPatch(patchI);
179 
180  if (localPatchI != -1)
181  {
182  vector nw;
183  vector Up;
184 
185  // patch-normal direction
186  this->owner().patchData(p, pp, trackFraction, tetIs, nw, Up);
187 
188  // particle velocity reletive to patch
189  const vector& U = p.U() - Up;
190 
191  // quick reject if particle travelling away from the patch
192  if ((nw & U) < 0)
193  {
194  return;
195  }
196 
197  const scalar magU = mag(U);
198  const vector Udir = U/magU;
199 
200  // determine impact angle, alpha
201  const scalar alpha = mathematical::pi/2.0 - acos(nw & Udir);
202 
203  const scalar coeff = p.nParticle()*p.mass()*sqr(magU)/(p_*psi_*K_);
204 
205  const label patchFaceI = pp.whichFace(p.face());
206  scalar& Q = QPtr_->boundaryField()[patchI][patchFaceI];
207  if (tan(alpha) < K_/6.0)
208  {
209  Q += coeff*(sin(2.0*alpha) - 6.0/K_*sqr(sin(alpha)));
210  }
211  else
212  {
213  Q += coeff*(K_*sqr(cos(alpha))/6.0);
214  }
215  }
216 }
217 
218 
219 // ************************************************************************* //
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual void postPatch(const parcelType &p, const polyPatch &pp, const scalar trackFraction, const tetIndices &tetIs, bool &keepParticle)
Post-patch hook.
#define readScalar
Definition: doubleScalar.C:38
dimensioned< scalar > mag(const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
bool empty() const
Return true if the UList is empty (ie, size() is zero).
Definition: UListI.H:313
label index() const
Return the index of this patch in the boundaryMesh.
A class for handling words, derived from string.
Definition: word.H:59
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
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
label nw
Definition: createFields.H:25
dictionary dict
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
stressControl lookup("compactNormalStress") >> compactNormalStress
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
#define WarningIn(functionName)
Report a warning using Foam::Warning.
label applyToPatch(const label globalPatchI) const
Returns local patchI if patch is in patchIds_ list.
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
#define forAll(list, i)
Definition: UList.H:421
dimensionedScalar cos(const dimensionedScalar &ds)
virtual ~ParticleErosion()
Destructor.
dimensionedScalar acos(const dimensionedScalar &ds)
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:201
virtual void write()
Write post-processing info.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensionedScalar tan(const dimensionedScalar &ds)
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
#define forAllReverse(list, i)
Definition: UList.H:424
error FatalError
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
virtual void preEvolve()
Pre-evolve hook.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:68
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Templated cloud function object base class.
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:73
This function object calculates and outputs the second invariant of the velocity gradient tensor [1/s...
Definition: Q.H:66
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
ParticleErosion(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
U
Definition: pEqn.H:82
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:377
Creates particle erosion field, Q.
dimensionedScalar sin(const dimensionedScalar &ds)
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116