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-2021 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  return patchIDs_.find(globalPatchi);
71 }
72 
73 
74 template<class CloudType>
76 {
77  if (QPtr_)
78  {
79  QPtr_->write();
80  }
81  else
82  {
84  << "QPtr not valid" << abort(FatalError);
85  }
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
91 template<class CloudType>
93 (
94  const dictionary& dict,
95  CloudType& owner,
96  const word& modelName
97 )
98 :
99  CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
100  QPtr_(nullptr),
101  patchIDs_(),
102  p_(this->coeffDict().getScalar("p")),
103  psi_(this->coeffDict().template getOrDefault<scalar>("psi", 2.0)),
104  K_(this->coeffDict().template getOrDefault<scalar>("K", 2.0))
105 {
106  const wordList allPatchNames(owner.mesh().boundaryMesh().names());
107  const wordRes patchNames
108  (
109  this->coeffDict().template get<wordRes>("patches")
110  );
111 
112  labelHashSet uniqIds;
113  for (const wordRe& re : patchNames)
114  {
115  labelList ids = findMatchingStrings(re, allPatchNames);
116 
117  if (ids.empty())
118  {
120  << "Cannot find any patch names matching " << re
121  << endl;
122  }
123 
124  uniqIds.insert(ids);
125  }
126 
127  patchIDs_ = uniqIds.sortedToc();
128 
129  // Trigger creation of the Q field
130  resetQ();
131 }
132 
133 
134 template<class CloudType>
136 (
138 )
139 :
141  QPtr_(nullptr),
142  patchIDs_(pe.patchIDs_),
143  p_(pe.p_),
144  psi_(pe.psi_),
145  K_(pe.K_)
146 {}
147 
148 
149 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150 
151 template<class CloudType>
153 (
154  const typename parcelType::trackingData& td
155 )
156 {
157  resetQ();
158 }
159 
160 
161 template<class CloudType>
163 (
164  const parcelType& p,
165  const polyPatch& pp,
166  bool&
167 )
168 {
169  const label patchi = pp.index();
170 
171  const label localPatchi = applyToPatch(patchi);
172 
173  if (localPatchi != -1)
174  {
175  vector nw;
176  vector Up;
177 
178  // patch-normal direction
179  this->owner().patchData(p, pp, nw, Up);
180 
181  // particle velocity relative to patch
182  const vector& U = p.U() - Up;
183 
184  // quick reject if particle travelling away from the patch
185  if ((nw & U) < 0)
186  {
187  return;
188  }
189 
190  const scalar magU = mag(U);
191  const vector Udir = U/magU;
192 
193  // determine impact angle, alpha
194  const scalar alpha = mathematical::piByTwo - acos(nw & Udir);
195 
196  const scalar coeff = p.nParticle()*p.mass()*sqr(magU)/(p_*psi_*K_);
197 
198  const label patchFacei = pp.whichFace(p.face());
199  scalar& Q = QPtr_->boundaryFieldRef()[patchi][patchFacei];
200  if (tan(alpha) < K_/6.0)
201  {
202  Q += coeff*(sin(2.0*alpha) - 6.0/K_*sqr(sin(alpha)));
203  }
204  else
205  {
206  Q += coeff*(K_*sqr(cos(alpha))/6.0);
207  }
208  }
209 }
210 
211 
212 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
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:65
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:369
Foam::ParticleErosion::write
virtual void write()
Write post-processing info.
Definition: ParticleErosion.C:75
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:80
ParticleErosion.H
Foam::polyBoundaryMesh::names
wordList names() const
Return a list of patch names.
Definition: polyBoundaryMesh.C:555
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::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
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:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::patchIdentifier::index
label index() const noexcept
The index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:147
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::ParticleErosion::postPatch
virtual void postPatch(const parcelType &p, const polyPatch &pp, bool &keepParticle)
Post-patch hook.
Definition: ParticleErosion.C:163
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:153
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
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:448
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:191
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:207
Foam::ParticleErosion
Creates particle erosion field, Q.
Definition: ParticleErosion.H:56
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::ParticleErosion::ParticleErosion
ParticleErosion(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: ParticleErosion.C:93
Foam::findMatchingStrings
labelList findMatchingStrings(const UnaryMatchPredicate &matcher, const UList< StringType > &input, const bool invert=false)
Extract list indices for all matches.
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265