Loading...
Searching...
No Matches
LocalInteraction.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) 2015-2025 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 "LocalInteraction.H"
30
31// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * //
32
33template<class CloudType>
35{
37
38 forAll(nEscape_, patchi)
39 {
40 const word& patchName = patchData_[patchi].patchName();
41
42 forAll(nEscape_[patchi], injectori)
43 {
44 const word suffix = Foam::name(injectori);
45 this->writeTabbed(os, patchName + "_nEscape_" + suffix);
46 this->writeTabbed(os, patchName + "_massEscape_" + suffix);
47 this->writeTabbed(os, patchName + "_nStick_" + suffix);
48 this->writeTabbed(os, patchName + "_massStick_" + suffix);
49 }
50 }
51}
52
53
54template<class CloudType>
56(
57 const dictionary& dict,
58 CloudType& cloud
59)
60:
61 PatchInteractionModel<CloudType>(dict, cloud, typeName),
62 patchData_(cloud.mesh(), this->coeffDict()),
63 nEscape_(patchData_.size()),
64 massEscape_(nEscape_.size()),
65 nStick_(nEscape_.size()),
66 massStick_(nEscape_.size()),
67 writeFields_(this->coeffDict().getOrDefault("writeFields", false)),
68 injIdToIndex_(),
69 massEscapePtr_(nullptr),
70 massStickPtr_(nullptr)
71{
72 const bool outputByInjectorId
73 = this->coeffDict().getOrDefault("outputByInjectorId", false);
74
75 if (writeFields_)
76 {
77 Info<< " Interaction fields will be written to "
78 << IOobject::scopedName(this->owner().name(), "massEscape")
79 << " and "
80 << IOobject::scopedName(this->owner().name(), "massStick") << endl;
81
82 (void)massEscape();
83 (void)massStick();
84 }
85 else
86 {
87 Info<< " Interaction fields will not be written" << endl;
88 }
89
90 // Determine the number of injectors and the injector mapping
91 label nInjectors = 0;
92 if (outputByInjectorId)
93 {
94 for (const auto& inj : cloud.injectors())
95 {
96 injIdToIndex_.insert(inj.injectorID(), nInjectors++);
97 }
98 }
99
100 // The normal case, and safety if injector mapping was somehow null.
101 if (!nInjectors)
102 {
103 nInjectors = 1;
104 }
105
106 // Check that interactions are valid/specified
107 forAll(patchData_, patchi)
108 {
109 const word& interactionTypeName =
110 patchData_[patchi].interactionTypeName();
112 this->wordToInteractionType(interactionTypeName);
113
115 {
116 const word& patchName = patchData_[patchi].patchName();
118 << "Unknown patch interaction type "
119 << interactionTypeName << " for patch " << patchName
120 << ". Valid selections are:"
122 << nl << exit(FatalError);
123 }
124
125 nEscape_[patchi].setSize(nInjectors, Zero);
126 massEscape_[patchi].setSize(nInjectors, Zero);
127 nStick_[patchi].setSize(nInjectors, Zero);
128 massStick_[patchi].setSize(nInjectors, Zero);
129 }
130}
131
132
133template<class CloudType>
135(
136 const LocalInteraction<CloudType>& pim
137)
138:
139 PatchInteractionModel<CloudType>(pim),
140 patchData_(pim.patchData_),
141 nEscape_(pim.nEscape_),
142 massEscape_(pim.massEscape_),
143 nStick_(pim.nStick_),
144 massStick_(pim.massStick_),
145 writeFields_(pim.writeFields_),
146 injIdToIndex_(pim.injIdToIndex_),
147 massEscapePtr_(nullptr),
148 massStickPtr_(nullptr)
149{}
150
151
152// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
153
154template<class CloudType>
156{
157 if (!massEscapePtr_)
158 {
159 const fvMesh& mesh = this->owner().mesh();
160
161 massEscapePtr_.reset
162 (
164 (
166 (
167 IOobject::scopedName(this->owner().name(), "massEscape"),
168 mesh.time().timeName(),
169 mesh.thisDb(),
173 ),
174 mesh,
176 )
177 );
179
180 return *massEscapePtr_;
181}
182
183
184template<class CloudType>
186{
187 if (!massStickPtr_)
188 {
189 const fvMesh& mesh = this->owner().mesh();
190
191 massStickPtr_.reset
192 (
194 (
195 IOobject
196 (
197 IOobject::scopedName(this->owner().name(), "massStick"),
198 mesh.time().timeName(),
199 mesh.thisDb(),
203 ),
204 mesh,
206 )
207 );
209
210 return *massStickPtr_;
211}
212
213
214template<class CloudType>
216(
217 typename CloudType::parcelType& p,
218 const polyPatch& pp,
219 bool& keepParticle
220)
221{
222 const label patchi = patchData_.applyToPatch(pp.index());
223
224 if (patchi >= 0)
225 {
226 vector& U = p.U();
227
228 // Location for storing the stats.
229 const label idx =
230 (
231 injIdToIndex_.size()
232 ? injIdToIndex_.lookup(p.typeId(), 0)
233 : 0
234 );
235
237 this->wordToInteractionType
238 (
239 patchData_[patchi].interactionTypeName()
240 );
241
242 switch (it)
243 {
245 {
246 return false;
247 }
249 {
250 keepParticle = false;
251 p.active(false);
252 U = Zero;
253
254 const scalar dm = p.mass()*p.nParticle();
255
256 nEscape_[patchi][idx]++;
257 massEscape_[patchi][idx] += dm;
258
259 if (writeFields_)
260 {
261 const label pI = pp.index();
262 const label fI = pp.whichFace(p.face());
263 massEscape().boundaryFieldRef()[pI][fI] += dm;
264 }
265 break;
266 }
268 {
269 keepParticle = true;
270 p.active(false);
271 U = Zero;
272
273 const scalar dm = p.mass()*p.nParticle();
274
275 nStick_[patchi][idx]++;
276 massStick_[patchi][idx] += dm;
277
278 if (writeFields_)
279 {
280 const label pI = pp.index();
281 const label fI = pp.whichFace(p.face());
282 massStick().boundaryFieldRef()[pI][fI] += dm;
283 }
284 break;
285 }
287 {
288 keepParticle = true;
289 p.active(true);
290
291 vector nw;
292 vector Up;
293
294 this->owner().patchData(p, pp, nw, Up);
295
296 // Calculate motion relative to patch velocity
297 U -= Up;
298
299 if (mag(Up) > 0 && mag(U) < this->Urmax())
300 {
302 << "Particle U the same as patch "
303 << " The particle has been removed" << nl << endl;
304
305 keepParticle = false;
306 p.active(false);
307 U = Zero;
308 break;
309 }
310
311 scalar Un = U & nw;
312 vector Ut = U - Un*nw;
313
314 if (Un > 0)
315 {
316 U -= (1.0 + patchData_[patchi].e())*Un*nw;
317 }
318
319 U -= patchData_[patchi].mu()*Ut;
320
321 // Return velocity to global space
322 U += Up;
323
324 break;
325 }
326 default:
327 {
329 << "Unknown interaction type "
330 << patchData_[patchi].interactionTypeName()
331 << "(" << it << ") for patch "
332 << patchData_[patchi].patchName()
333 << ". Valid selections are:" << this->interactionTypeNames_
334 << endl << abort(FatalError);
335 }
336 }
337
338 return true;
340
341 return false;
342}
343
344
345template<class CloudType>
347{
349
350 // retrieve any stored data
351 labelListList npe0(patchData_.size());
352 scalarListList mpe0(patchData_.size());
353 labelListList nps0(patchData_.size());
354 scalarListList mps0(patchData_.size());
355
356 forAll(patchData_, patchi)
357 {
358 label lsd = nEscape_[patchi].size();
359 npe0[patchi].setSize(lsd, Zero);
360 mpe0[patchi].setSize(lsd, Zero);
361 nps0[patchi].setSize(lsd, Zero);
362 mps0[patchi].setSize(lsd, Zero);
363 }
364
365
366 this->getModelProperty("nEscape", npe0);
367 this->getModelProperty("massEscape", mpe0);
368 this->getModelProperty("nStick", nps0);
369 this->getModelProperty("massStick", mps0);
370
371 // accumulate current data
372 labelListList npe(nEscape_);
373 forAll(npe, i)
374 {
376 npe[i] = npe[i] + npe0[i];
377 }
378
379 scalarListList mpe(massEscape_);
380 forAll(mpe, i)
381 {
383 mpe[i] = mpe[i] + mpe0[i];
384 }
385
386 labelListList nps(nStick_);
387 forAll(nps, i)
388 {
390 nps[i] = nps[i] + nps0[i];
391 }
392
393 scalarListList mps(massStick_);
394 forAll(nps, i)
395 {
397 mps[i] = mps[i] + mps0[i];
398 }
399
400 if (injIdToIndex_.size())
401 {
402 // Since injIdToIndex_ is a one-to-one mapping (starting at zero),
403 // can simply invert it.
404 labelList indexToInjector(injIdToIndex_.size());
405 forAllConstIters(injIdToIndex_, iter)
406 {
407 indexToInjector[iter.val()] = iter.key();
408 }
409
410 forAll(patchData_, patchi)
411 {
412 forAll(mpe[patchi], indexi)
413 {
414 const word& patchName = patchData_[patchi].patchName();
415
416 Log_<< " Parcel fate: patch " << patchName
417 << " (number, mass)" << nl
418 << " - escape (injector " << indexToInjector[indexi]
419 << " ) = " << npe[patchi][indexi]
420 << ", " << mpe[patchi][indexi] << nl
421 << " - stick (injector " << indexToInjector[indexi]
422 << " ) = " << nps[patchi][indexi]
423 << ", " << mps[patchi][indexi] << nl;
424 }
425 }
426 }
427 else
428 {
429 forAll(patchData_, patchi)
430 {
431 const word& patchName = patchData_[patchi].patchName();
432
433 Log_<< " Parcel fate: patch " << patchName
434 << " (number, mass)" << nl
435 << " - escape = "
436 << npe[patchi][0] << ", " << mpe[patchi][0] << nl
437 << " - stick = "
438 << nps[patchi][0] << ", " << mps[patchi][0] << nl;
439 }
440 }
441
442 forAll(npe, patchi)
443 {
444 forAll(npe[patchi], injectori)
445 {
446 this->file()
447 << tab << npe[patchi][injectori]
448 << tab << mpe[patchi][injectori]
449 << tab << nps[patchi][injectori]
450 << tab << mps[patchi][injectori];
451 }
452 }
453
454 this->file() << endl;
455
456 if (this->writeTime())
457 {
458 this->setModelProperty("nEscape", npe);
459 this->setModelProperty("massEscape", mpe);
460 this->setModelProperty("nStick", nps);
461 this->setModelProperty("massStick", mps);
462
463 nEscape_ = Zero;
464 massEscape_ = Zero;
465 nStick_ = Zero;
466 massStick_ = Zero;
467 }
468}
469
470
471// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const CloudType & owner() const
Return const access to the owner cloud.
virtual bool writeTime() const
Flag to indicate when to write a property.
@ REGISTER
Request registration (bool: true).
@ READ_IF_PRESENT
Reading is optional [identical to LAZY_READ].
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
void setSize(label n)
Alias for resize().
Definition List.H:536
LocalInteraction(const dictionary &dict, CloudType &owner)
Construct from dictionary.
volScalarField & massStick()
Return access to the massStick field.
volScalarField & massEscape()
Return access to the massEscape field.
virtual bool correct(typename CloudType::parcelType &p, const polyPatch &pp, bool &keepParticle)
Apply velocity correction.
virtual void info()
Write patch interaction info.
virtual void writeFileHeader(Ostream &os)
Output file header information.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
const scalar & Urmax() const
Return Urmax.
PatchInteractionModel(CloudType &owner)
Construct null from owner.
virtual void info()
Write patch interaction info.
static interactionType wordToInteractionType(const word &itWord)
Convert word to interaction result.
virtual void writeFileHeader(Ostream &os)
Output file header information.
static void listGather(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather (reduce) list elements, applying bop to each list element.
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition writeFile.C:334
virtual OFstream & file()
Return access to the file (if only 1).
Definition writeFile.C:270
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
bool getModelProperty(const word &entryName, Type &value) const
Retrieve generic property from the sub-model.
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
void setModelProperty(const word &entryName, const Type &value)
Add generic property to the sub-model.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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
OBJstream os(runTime.globalPath()/outputName)
auto & name
#define Log_
Report write to Foam::Info if the class log switch is true.
#define WarningInFunction
Report a warning using Foam::Warning.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
DSMCCloud< dsmcParcel > CloudType
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
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...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235