ParticleErosion.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 "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_(nullptr),
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 creation of the Q field
101  preEvolve();
102 }
103 
104 
105 template<class CloudType>
107 (
109 )
110 :
112  QPtr_(nullptr),
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,
154  )
155  );
156  }
157 }
158 
159 
160 template<class CloudType>
162 (
163  const parcelType& p,
164  const polyPatch& pp,
165  bool&
166 )
167 {
168  const label patchi = pp.index();
169  const label localPatchi = applyToPatch(patchi);
170 
171  if (localPatchi != -1)
172  {
173  // Get patch data
174  vector nw, Up;
175  this->owner().patchData(p, pp, nw, Up);
176 
177  // Particle velocity relative to patch
178  const vector& U = p.U() - Up;
179 
180  // Quick rejection if the particle is travelling away from the patch
181  if ((nw & U) < 0)
182  {
183  return;
184  }
185 
186  const scalar magU = mag(U);
187  const vector UHat = U/magU;
188 
189  // Impact angle
190  const scalar alpha = mathematical::pi/2 - acos(nw & UHat);
191 
192  // Get the face value to accumulate into
193  const label patchFacei = pp.whichFace(p.face());
194  scalar& Q = QPtr_->boundaryFieldRef()[patchi][patchFacei];
195 
196  // Finnie's model
197  const scalar coeff = p.nParticle()*p.mass()*sqr(magU)/(p_*psi_*K_);
198  if (tan(alpha) < K_/6)
199  {
200  Q += coeff*(sin(2*alpha) - 6/K_*sqr(sin(alpha)));
201  }
202  else
203  {
204  Q += coeff*(K_*sqr(cos(alpha))/6);
205  }
206  }
207 }
208 
209 
210 // ************************************************************************* //
#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:434
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:313
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
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
ParticleErosion(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Function object to create a field of eroded volume, Q, on a specified list of patches. The volume is calculated by the model of Finnie et al. The implementation follows the description given by the review of Yadav et al.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:446
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:626
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:108
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.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
dimensionedScalar sin(const dimensionedScalar &ds)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
label patchi
U
Definition: pEqn.H:72
#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
label index() const
Return the index of this patch in the boundaryMesh.
dimensioned< scalar > mag(const dimensioned< Type > &)
List< Key > toc() const
Return the table of contents.
Definition: HashTable.C:202
label applyToPatch(const label globalPatchi) const
Returns local patchi if patch is in patchIds_ list.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
virtual void postPatch(const parcelType &p, const polyPatch &pp, bool &keepParticle)
Post-patch hook.
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)
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
Templated cloud function object base class.
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:380