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 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,
35 trackingData& td,
36 const scalar trackTime
37)
38{
39 typename TrackCloudType::particleType& p =
40 static_cast<typename TrackCloudType::particleType&>(*this);
41
42
43 // Check position is inside tet
44 //checkInside();
45
46 td.switchProcessor = false;
47 td.keepParticle = true;
48
49 scalar tEnd = (1.0 - stepFraction())*trackTime;
50 scalar maxDt = mesh().bounds().mag();
51
52 while
53 (
55 && !td.switchProcessor
56 && lifeTime_ > 0
57 )
58 {
59 // set the lagrangian time-step
60 scalar dt = maxDt;
61
62 --lifeTime_;
63
64 // Get sampled velocity and fields. Store if position changed.
65 vector U = p.sample(td);
66
67 // !user parameter!
68 if (dt < SMALL)
69 {
70 // Force removal
71 lifeTime_ = 0;
72 break;
73 }
74
75
76 if (td.trackLength_ < GREAT)
77 {
78 dt = td.trackLength_;
79 }
80
81
82 scalar fraction = trackToEdge(cloud, td, localPosition_ + dt*U);
83 dt *= fraction;
84
85 tEnd -= dt;
86 stepFraction() = 1.0 - tEnd/trackTime;
87
88
89 if (tEnd <= ROOTVSMALL)
90 {
91 // Force removal
92 lifeTime_ = 0;
93 }
94 }
95
96
97 if (!td.keepParticle || lifeTime_ == 0)
98 {
99 if (lifeTime_ == 0)
100 {
101 if (debug)
102 {
103 Pout<< "wallBoundedStreamLineParticle :"
104 << " Removing stagnant particle:"
106 << " sampled positions:" << sampledPositions_.size()
107 << endl;
108 }
109 td.keepParticle = false;
110 }
111 else
112 {
113 // Normal exit. Store last position and fields
114 sample(td);
115
116 if (debug)
117 {
118 Pout<< "wallBoundedStreamLineParticle : Removing particle:"
120 << " sampled positions:" << sampledPositions_.size()
121 << endl;
122 }
123 }
124
125 // Transfer particle data into trackingData.
126 {
127 //td.allPositions_.append(sampledPositions_);
128 td.allPositions_.append(vectorList());
129 vectorList& top = td.allPositions_.last();
131 }
132
134 {
135 //td.allScalars_[i].append(sampledScalars_[i]);
136 td.allScalars_[i].append(scalarList());
137 scalarList& top = td.allScalars_[i].last();
139 }
141 {
142 //td.allVectors_[i].append(sampledVectors_[i]);
143 td.allVectors_[i].append(vectorList());
144 vectorList& top = td.allVectors_[i].last();
146 }
147 }
148
149 return td.keepParticle;
150}
151
152
153// ************************************************************************* //
Minimal example by using system/controlDict.functions:
void transfer(List< T > &list)
Definition: List.C:447
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:133
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
virtual void move()=0
bool switchProcessor
Flag to switch processor.
Definition: particle.H:102
bool keepParticle
Flag to indicate whether to keep particle (false = delete)
Definition: particle.H:105
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition: particleI.H:137
scalar stepFraction() const noexcept
Return the fraction of time-step completed.
Definition: particleI.H:197
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:462
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.
List< DynamicList< scalar > > sampledScalars_
Sampled scalars.
U
Definition: pEqn.H:72
volScalarField & p
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
List< vector > vectorList
A List of vectors.
Definition: vectorList.H:64
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333