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-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 "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  {
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  {
91  << "Cannot find any patch names matching " << patchName[i]
92  << endl;
93  }
94 
95  uniquePatchIDs.insert(patchIDs);
96  }
97 
98  patchIDs_ = uniquePatchIDs.toc();
99 
100  // trigger ther creation of the Q field
101  preEvolve();
102 }
103 
104 
105 template<class CloudType>
107 (
109 )
110 :
112  QPtr_(NULL),
113  patchIDs_(pe.patchIDs_),
114  p_(pe.p_),
115  psi_(pe.psi_),
116  K_(pe.K_)
117 {}
118 
119 
120 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
121 
122 template<class CloudType>
124 {}
125 
126 
127 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
128 
129 template<class CloudType>
131 {
132  if (QPtr_.valid())
133  {
134  QPtr_->primitiveFieldRef() = 0.0;
135  }
136  else
137  {
138  const fvMesh& mesh = this->owner().mesh();
139 
140  QPtr_.reset
141  (
142  new volScalarField
143  (
144  IOobject
145  (
146  this->owner().name() + "Q",
147  mesh.time().timeName(),
148  mesh,
149  IOobject::READ_IF_PRESENT,
150  IOobject::NO_WRITE
151  ),
152  mesh,
153  dimensionedScalar("zero", dimVolume, 0.0)
154  )
155  );
156  }
157 }
158 
159 
160 template<class CloudType>
162 (
163  const parcelType& p,
164  const polyPatch& pp,
165  const scalar trackFraction,
166  const tetIndices& tetIs,
167  bool&
168 )
169 {
170  const label patchi = pp.index();
171 
172  const label localPatchi = applyToPatch(patchi);
173 
174  if (localPatchi != -1)
175  {
176  vector nw;
177  vector Up;
178 
179  // patch-normal direction
180  this->owner().patchData(p, pp, trackFraction, tetIs, nw, Up);
181 
182  // particle velocity reletive to patch
183  const vector& U = p.U() - Up;
184 
185  // quick reject if particle travelling away from the patch
186  if ((nw & U) < 0)
187  {
188  return;
189  }
190 
191  const scalar magU = mag(U);
192  const vector Udir = U/magU;
193 
194  // determine impact angle, alpha
195  const scalar alpha = mathematical::pi/2.0 - acos(nw & Udir);
196 
197  const scalar coeff = p.nParticle()*p.mass()*sqr(magU)/(p_*psi_*K_);
198 
199  const label patchFacei = pp.whichFace(p.face());
200  scalar& Q = QPtr_->boundaryFieldRef()[patchi][patchFacei];
201  if (tan(alpha) < K_/6.0)
202  {
203  Q += coeff*(sin(2.0*alpha) - 6.0/K_*sqr(sin(alpha)));
204  }
205  else
206  {
207  Q += coeff*(K_*sqr(cos(alpha))/6.0);
208  }
209  }
210 }
211 
212 
213 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
dictionary dict
dimensionedScalar acos(const dimensionedScalar &ds)
#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
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:313
ParticleErosion(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Creates particle erosion field, Q.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:440
virtual void preEvolve()
Pre-evolve hook.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
virtual ~ParticleErosion()
Destructor.
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
stressControl lookup("compactNormalStress") >> compactNormalStress
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
A class for handling words, derived from string.
Definition: word.H:59
virtual void write()
Write post-processing info.
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
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:198
dimensionedScalar sin(const dimensionedScalar &ds)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label patchi
#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
dimensioned< scalar > mag(const dimensioned< Type > &)
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
dimensionedScalar tan(const dimensionedScalar &ds)
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
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
Templated cloud function object base class.
label index() const
Return the index of this patch in the boundaryMesh.
label applyToPatch(const label globalPatchi) const
Returns local patchi if patch is in patchIds_ list.
label nw
Definition: createFields.H:25
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
virtual void postPatch(const parcelType &p, const polyPatch &pp, const scalar trackFraction, const tetIndices &tetIs, bool &keepParticle)
Post-patch hook.