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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "ParticleErosion.H"
30 
31 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
32 
33 template<class CloudType>
35 {
36  if (QPtr_)
37  {
38  QPtr_->primitiveFieldRef() = 0.0;
39  }
40  else
41  {
42  const fvMesh& mesh = this->owner().mesh();
43 
44  QPtr_.reset
45  (
46  new volScalarField
47  (
48  IOobject
49  (
50  this->owner().name() + "Q",
51  mesh.time().timeName(),
52  mesh,
53  IOobject::READ_IF_PRESENT,
54  IOobject::NO_WRITE
55  ),
56  mesh,
58  )
59  );
60  }
61 }
62 
63 
64 template<class CloudType>
66 (
67  const label globalPatchi
68 ) const
69 {
70  forAll(patchIDs_, i)
71  {
72  if (patchIDs_[i] == globalPatchi)
73  {
74  return i;
75  }
76  }
77 
78  return -1;
79 }
80 
81 
82 template<class CloudType>
84 {
85  if (QPtr_)
86  {
87  QPtr_->write();
88  }
89  else
90  {
92  << "QPtr not valid" << abort(FatalError);
93  }
94 }
95 
96 
97 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
98 
99 template<class CloudType>
101 (
102  const dictionary& dict,
103  CloudType& owner,
104  const word& modelName
105 )
106 :
107  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
108  QPtr_(nullptr),
109  patchIDs_(),
110  p_(this->coeffDict().getScalar("p")),
111  psi_(this->coeffDict().template getOrDefault<scalar>("psi", 2.0)),
112  K_(this->coeffDict().template getOrDefault<scalar>("K", 2.0))
113 {
114  const wordList allPatchNames(owner.mesh().boundaryMesh().names());
115  const wordRes patchNames
116  (
117  this->coeffDict().template get<wordRes>("patches")
118  );
119 
120  labelHashSet uniqIds;
121  for (const wordRe& re : patchNames)
122  {
123  labelList ids = findStrings(re, allPatchNames);
124 
125  if (ids.empty())
126  {
128  << "Cannot find any patch names matching " << re
129  << endl;
130  }
131 
132  uniqIds.insert(ids);
133  }
134 
135  patchIDs_ = uniqIds.sortedToc();
136 
137  // Trigger creation of the Q field
138  resetQ();
139 }
140 
141 
142 template<class CloudType>
144 (
146 )
147 :
149  QPtr_(nullptr),
150  patchIDs_(pe.patchIDs_),
151  p_(pe.p_),
152  psi_(pe.psi_),
153  K_(pe.K_)
154 {}
155 
156 
157 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158 
159 template<class CloudType>
161 (
162  const typename parcelType::trackingData& td
163 )
164 {
165  resetQ();
166 }
167 
168 
169 template<class CloudType>
171 (
172  const parcelType& p,
173  const polyPatch& pp,
174  bool&
175 )
176 {
177  const label patchi = pp.index();
178 
179  const label localPatchi = applyToPatch(patchi);
180 
181  if (localPatchi != -1)
182  {
183  vector nw;
184  vector Up;
185 
186  // patch-normal direction
187  this->owner().patchData(p, pp, nw, Up);
188 
189  // particle velocity relative to patch
190  const vector& U = p.U() - Up;
191 
192  // quick reject if particle travelling away from the patch
193  if ((nw & U) < 0)
194  {
195  return;
196  }
197 
198  const scalar magU = mag(U);
199  const vector Udir = U/magU;
200 
201  // determine impact angle, alpha
202  const scalar alpha = mathematical::piByTwo - acos(nw & Udir);
203 
204  const scalar coeff = p.nParticle()*p.mass()*sqr(magU)/(p_*psi_*K_);
205 
206  const label patchFacei = pp.whichFace(p.face());
207  scalar& Q = QPtr_->boundaryFieldRef()[patchi][patchFacei];
208  if (tan(alpha) < K_/6.0)
209  {
210  Q += coeff*(sin(2.0*alpha) - 6.0/K_*sqr(sin(alpha)));
211  }
212  else
213  {
214  Q += coeff*(K_*sqr(cos(alpha))/6.0);
215  }
216  }
217 }
218 
219 
220 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::tan
dimensionedScalar tan(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:266
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::ParticleErosion::applyToPatch
label applyToPatch(const label globalPatchi) const
Returns local patchi if patch is in patchIds_ list.
Definition: ParticleErosion.C:66
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::ParticleErosion::write
virtual void write()
Write post-processing info.
Definition: ParticleErosion.C:83
Foam::HashSet< label, Hash< label > >
Foam::wordRe
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition: wordRe.H:72
ParticleErosion.H
Foam::polyBoundaryMesh::names
wordList names() const
Return a list of patch names.
Definition: polyBoundaryMesh.C:593
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::DSMCCloud::mesh
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
patchNames
wordList patchNames(nPatches)
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::constant::atomic::re
const dimensionedScalar re
Classical electron radius: default SI units: [m].
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam::findStrings
labelList findStrings(const regExp &matcher, const UList< StringType > &input, const bool invert=false)
Return list indices for strings matching the regular expression.
Definition: stringListOps.H:76
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::HashTable::sortedToc
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:136
Foam::ParticleErosion::postPatch
virtual void postPatch(const parcelType &p, const polyPatch &pp, bool &keepParticle)
Post-patch hook.
Definition: ParticleErosion.C:171
U
U
Definition: pEqn.H:72
Foam::constant::mathematical::piByTwo
constexpr scalar piByTwo(0.5 *M_PI)
Foam::ParticleErosion::preEvolve
virtual void preEvolve(const typename parcelType::trackingData &td)
Pre-evolve hook.
Definition: ParticleErosion.C:161
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::polyPatch::whichFace
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:400
Foam::Vector< scalar >
Foam::List< word >
Foam::CloudFunctionObject
Templated cloud function object base class.
Definition: CloudFunctionObject.H:62
Foam::acos
dimensionedScalar acos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:268
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::ParticleErosion::resetQ
void resetQ()
Create|read|reset the Q field.
Definition: ParticleErosion.C:34
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:181
Foam::primitiveMesh::reset
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
Definition: primitiveMesh.C:208
Foam::ParticleErosion
Creates particle erosion field, Q.
Definition: ParticleErosion.H:56
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::ParticleErosion::ParticleErosion
ParticleErosion(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: ParticleErosion.C:101
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
Foam::patchIdentifier::index
label index() const
The index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:158
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265