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-------------------------------------------------------------------------------
11License
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
33template<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 (
47 (
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
64template<class CloudType>
66(
67 const label globalPatchi
68) const
69{
70 return patchIDs_.find(globalPatchi);
71}
72
73
74template<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
91template<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
134template<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
151template<class CloudType>
153(
154 const typename parcelType::trackingData& td
155)
156{
157 resetQ();
158}
159
160
161template<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// ************************************************************************* //
Templated cloud function object base class.
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
const fvMesh & mesh() const
Return reference to the mesh.
Definition: DSMCCloudI.H:44
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Creates particle erosion field, Q.
virtual void postPatch(const parcelType &p, const polyPatch &pp, bool &keepParticle)
Post-patch hook.
label applyToPatch(const label globalPatchi) const
Returns local patchi if patch is in patchIds_ list.
virtual void write()
Write post-processing info.
virtual void preEvolve(const typename parcelType::trackingData &td)
Pre-evolve hook.
void resetQ()
Create|read|reset the Q field.
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Class used to pass data into container.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
label index() const noexcept
The index of this patch in the boundaryMesh.
wordList names() const
Return a list of patch names.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:451
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:131
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition: wordRe.H:83
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
#define WarningInFunction
Report a warning using Foam::Warning.
dimensionedScalar tan(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sin(const dimensionedScalar &ds)
labelList findMatchingStrings(const UnaryMatchPredicate &matcher, const UList< StringType > &input, const bool invert=false)
Extract list indices for all matches.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
error FatalError
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
wordList patchNames(nPatches)
volScalarField & alpha
dictionary dict