Loading...
Searching...
No Matches
wallBoundedStreamLineParticleTemplates.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-2016 OpenFOAM Foundation
9 Copyright (C) 2017-2023 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// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
30
31template<class TrackCloudType>
33(
34 TrackCloudType& cloud,
36 const scalar trackTime
37)
38{
39 auto& p = static_cast<typename TrackCloudType::particleType&>(*this);
40
41 // Check position is inside tet
42 //checkInside();
43
44 td.switchProcessor = false;
45 td.keepParticle = true;
46
47 scalar tEnd = (1.0 - stepFraction())*trackTime;
48 scalar maxDt = mesh().bounds().mag();
49
50 while (td.keepParticle && !td.switchProcessor && lifeTime_ > 0)
51 {
52 // Set the lagrangian time-step
53 scalar dt = maxDt;
54
55 --lifeTime_;
56
57 // Get sampled velocity and fields. Store if position changed.
58 vector U = p.sample(td);
59
60 // !user parameter!
61 if (dt < SMALL)
62 {
63 // Force removal
64 lifeTime_ = 0;
65 break;
66 }
67
68
69 if (td.trackLength_ < GREAT)
70 {
71 dt = td.trackLength_;
72 }
73
74
75 scalar fraction = trackToEdge(cloud, td, localPosition_ + dt*U);
76 dt *= fraction;
77
78 tEnd -= dt;
79 stepFraction() = 1.0 - tEnd/trackTime;
80
81
82 if (tEnd <= ROOTVSMALL)
83 {
84 // Force removal
85 lifeTime_ = 0;
86 }
87 }
88
89
90 if (!td.keepParticle || lifeTime_ == 0)
91 {
92 if (lifeTime_ == 0)
93 {
94 if (debug)
95 {
96 Pout<< "wallBoundedStreamLineParticle :"
97 << " Removing stagnant particle:"
99 << " sampled positions:" << sampledPositions_.size()
100 << endl;
101 }
102 td.keepParticle = false;
103 }
104 else
105 {
106 // Normal exit. Store last position and fields
107 sample(td);
108
109 if (debug)
110 {
111 Pout<< "wallBoundedStreamLineParticle : Removing particle:"
113 << " sampled positions:" << sampledPositions_.size()
114 << endl;
115 }
116 }
117
118 // Transfer particle data into trackingData.
119 {
120 td.allPositions_.emplace_back().transfer(sampledPositions_);
121 }
122
124 {
125 td.allScalars_[i].emplace_back().transfer(sampledScalars_[i]);
126 }
128 {
129 td.allVectors_[i].emplace_back().transfer(sampledVectors_[i]);
130 }
131 }
132
133 return td.keepParticle;
134}
135
136
137// ************************************************************************* //
scalar mag() const
The magnitude/length of the bounding box diagonal.
Definition boundBoxI.H:198
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition particleI.H:110
scalar stepFraction() const noexcept
Return the fraction of time-step completed.
Definition particleI.H:170
const boundBox & bounds() const noexcept
Return mesh bounding box.
Definition polyMesh.H:617
scalar trackToEdge(TrackCloudType &cloud, trackingData &td, const vector &endPosition)
Equivalent of trackToFace.
point localPosition_
Particle position is updated locally as opposed to via track.
Class used to pass tracking data to the trackToEdge function.
label lifeTime_
Lifetime of particle. Particle dies when reaches 0.
DynamicList< point > sampledPositions_
Sampled positions.
List< DynamicList< vector > > sampledVectors_
Sampled vectors.
bool move(TrackCloudType &cloud, trackingData &td, const scalar trackTime)
Track all particles to their end point.
List< DynamicList< scalar > > sampledScalars_
Sampled scalars.
U
Definition pEqn.H:72
volScalarField & p
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for handling debugging switches.
Definition debug.C:45
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Vector< scalar > vector
Definition vector.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299