Loading...
Searching...
No Matches
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\*---------------------------------------------------------------------------*/
29#include "ParticleErosion.H"
30#include "wordRes.H"
31
32// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
33
34template<class CloudType>
36{
37 if (QPtr_)
38 {
39 QPtr_->primitiveFieldRef() = 0.0;
40 }
41 else
42 {
43 const fvMesh& mesh = this->owner().mesh();
44
45 QPtr_.reset
46 (
48 (
50 (
51 this->owner().name() + "Q",
52 mesh.time().timeName(),
53 mesh,
56 ),
57 mesh,
59 )
60 );
61 }
62}
63
64
65template<class CloudType>
67(
68 const label globalPatchi
69) const
70{
71 return patchIDs_.find(globalPatchi);
72}
73
74
75template<class CloudType>
77{
78 if (QPtr_)
79 {
80 QPtr_->write();
81 }
82 else
83 {
85 << "QPtr not valid" << abort(FatalError);
86 }
87}
88
89
90// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
91
92template<class CloudType>
94(
95 const dictionary& dict,
96 CloudType& owner,
97 const word& modelName
98)
99:
100 CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
101 QPtr_(nullptr),
102 patchIDs_(),
103 p_(this->coeffDict().getScalar("p")),
104 psi_(this->coeffDict().template getOrDefault<scalar>("psi", 2.0)),
105 K_(this->coeffDict().template getOrDefault<scalar>("K", 2.0))
106{
107 const wordList allPatchNames(owner.mesh().boundaryMesh().names());
108 const wordRes patchNames
109 (
110 this->coeffDict().template get<wordRes>("patches")
111 );
112
113 labelHashSet uniqIds;
114 for (const wordRe& select : patchNames)
115 {
116 labelList ids = wordRes::matching(select, allPatchNames);
117
118 if (ids.empty())
119 {
121 << "Cannot find any patch names matching "
122 << select << nl;
123 }
124
125 uniqIds.insert(ids);
126 }
127
128 patchIDs_ = uniqIds.sortedToc();
130 // Trigger creation of the Q field
131 resetQ();
132}
133
134
135template<class CloudType>
137(
139)
140:
142 QPtr_(nullptr),
143 patchIDs_(pe.patchIDs_),
144 p_(pe.p_),
145 psi_(pe.psi_),
146 K_(pe.K_)
147{}
148
149
150// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
151
152template<class CloudType>
154(
155 const typename parcelType::trackingData& td
157{
158 resetQ();
159}
160
161
162template<class CloudType>
164(
165 const parcelType& p,
166 const polyPatch& pp,
167 const typename parcelType::trackingData& td
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, nw, Up);
181
182 // particle velocity relative 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 true;
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::piByTwo - 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 return true;
212}
213
214
215// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Templated cloud function object base class.
CloudFunctionObject(CloudType &owner)
Construct null from owner.
const CloudType & owner() const
Return const access to the owner cloud.
const fvMesh & mesh() const
Return reference to the mesh.
Definition DSMCCloudI.H:37
particle::trackingData trackingData
Definition DSMCParcel.H:155
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition HashTable.C:156
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Creates particle erosion field, Q.
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.
virtual bool postPatch(const parcelType &p, const polyPatch &pp, const typename parcelType::trackingData &td)
Post-patch hook.
ParticleErosion(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
wordList names() const
Return a list of patch names.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
const word & modelName() const
Return const access to the name of the sub-model.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A wordRe is a Foam::word, but can contain a regular expression for matching words or strings.
Definition wordRe.H:81
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
static labelList matching(const wordRe &select, const UList< StringType > &input, const bool invert=false)
Determine the list indices for all matches.
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
auto & name
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
#define WarningInFunction
Report a warning using Foam::Warning.
List< word > wordList
List of word.
Definition fileName.H:60
DSMCCloud< dsmcParcel > CloudType
dimensionedScalar tan(const dimensionedScalar &ds)
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionedScalar sin(const dimensionedScalar &ds)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
wordList patchNames(nPatches)
volScalarField & alpha
dictionary dict