streamLineParticle.H
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-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::streamLineParticle
28
29Description
30 Particle class that samples fields as it passes through. Used in streamline
31 calculation.
32
33SourceFiles
34 streamLineParticle.C
35
36\*---------------------------------------------------------------------------*/
37
38#ifndef streamLineParticle_H
39#define streamLineParticle_H
40
41#include "particle.H"
42#include "autoPtr.H"
43#include "interpolation.H"
44#include "vectorList.H"
45
46// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47
48namespace Foam
49{
50
51// Forward declaration of friend functions and operators
52
53class streamLineParticle;
54class streamLineParticleCloud;
55
56Ostream& operator<<(Ostream&, const streamLineParticle&);
57
58/*---------------------------------------------------------------------------*\
59 Class streamLineParticle Declaration
60\*---------------------------------------------------------------------------*/
63:
64 public particle
65{
66public:
68 class trackingData
69 :
71 {
72 public:
73
74 // Public data
80 const label UIndex_;
82 const label nSubCycle_;
84 const scalar trackLength_;
91
92
93 // Constructors
94
95 //- Construct from components
97 (
99 const PtrList<interpolation<scalar>>& vsInterp,
100 const PtrList<interpolation<vector>>& vvInterp,
101 const label UIndex,
102 const label nSubCycle,
103 const scalar trackLength,
104 DynamicList<List<point>>& allPositions,
105 List<DynamicList<scalarList>>& allScalars,
106 List<DynamicList<vectorList>>& allVectors
107 )
108 :
110 vsInterp_(vsInterp),
111 vvInterp_(vvInterp),
112 UIndex_(UIndex),
113 nSubCycle_(nSubCycle),
114 trackLength_(trackLength),
115 allPositions_(allPositions),
116 allScalars_(allScalars),
117 allVectors_(allVectors)
118 {}
119 };
120
121
122private:
123
124 // Private data
125
126 //- Whether particle transports with +U or -U
127 bool trackForward_;
128
129 //- Lifetime of particle. Particle dies when reaches 0.
130 label lifeTime_;
131
132 //- Sampled positions
133 DynamicList<point> sampledPositions_;
134
135 //- Sampled scalars
136 List<DynamicList<scalar>> sampledScalars_;
137
138 //- Sampled vectors
139 List<DynamicList<vector>> sampledVectors_;
140
141
142 // Private Member Functions
143
144 //- Interpolate all quantities; return interpolated velocity.
145 vector interpolateFields
146 (
147 const trackingData&,
148 const point&,
149 const label celli,
150 const label facei
151 );
152
153
154public:
155
156 // Constructors
157
158 //- Construct from components
160 (
161 const polyMesh& c,
162 const vector& position,
163 const label celli,
164 const bool trackForward,
165 const label lifeTime
166 );
167
168 //- Construct from Istream
170 (
171 const polyMesh& c,
172 Istream& is,
173 bool readFields = true,
174 bool newFormat = true
175 );
176
177 //- Construct copy
179
180 //- Construct and return a clone
182 {
183 return autoPtr<particle>(new streamLineParticle(*this));
184 }
185
186 //- Factory class to read-construct particles used for parallel transfer
187 class iNew
188 {
189 const polyMesh& mesh_;
190
191 public:
193 iNew(const polyMesh& mesh)
194 :
195 mesh_(mesh)
196 {}
199 {
201 (
202 new streamLineParticle(mesh_, is, true)
203 );
204 }
205 };
206
207
208 // Member Functions
209
210 // Tracking
211
212 //- Track all particles to their end point
213 bool move(streamLineParticleCloud&, trackingData&, const scalar);
214
215 //- Overridable function to handle the particle hitting a patch
216 // Executed before other patch-hitting functions
217 bool hitPatch(streamLineParticleCloud&, trackingData&);
218
219 //- Overridable function to handle the particle hitting a wedge
220 void hitWedgePatch(streamLineParticleCloud&, trackingData&);
221
222 //- Overridable function to handle the particle hitting a
223 // symmetry plane
225
226 //- Overridable function to handle the particle hitting a
227 // symmetry patch
228 void hitSymmetryPatch(streamLineParticleCloud&, trackingData&);
229
230 //- Overridable function to handle the particle hitting a cyclic
231 void hitCyclicPatch(streamLineParticleCloud&, trackingData&);
232
233 //- Overridable function to handle the particle hitting a
234 // cyclicAMIPatch
236 (
238 trackingData&,
239 const vector& direction
240 );
241
242 //- Overridable function to handle the particle hitting a
243 // cyclicACMIPatch
245 (
247 trackingData&,
248 const vector& direction
249 );
250
251 //- Overridable function to handle the particle hitting a
252 //- processorPatch
253 void hitProcessorPatch(streamLineParticleCloud&, trackingData&);
254
255 //- Overridable function to handle the particle hitting a wallPatch
256 void hitWallPatch(streamLineParticleCloud&, trackingData&);
257
258
259 // I-O
260
261 //- Read
263
264 //- Write
265 static void writeFields(const Cloud<streamLineParticle>&);
266
267
268 // Ostream Operator
271};
272
273
274// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275
276} // End namespace Foam
277
278// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
279
280#endif
281
282// ************************************************************************* //
Base cloud calls templated on particle type.
Definition: Cloud.H:68
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:64
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
Abstract base class for volume field interpolation.
Definition: interpolation.H:60
Base particle class.
Definition: particle.H:79
vector position() const
Return current particle position.
Definition: particleI.H:314
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition: particleI.H:137
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A Cloud of streamLine particles.
Factory class to read-construct particles used for parallel transfer.
autoPtr< streamLineParticle > operator()(Istream &is) const
const PtrList< interpolation< vector > > & vvInterp_
trackingData(streamLineParticleCloud &cloud, const PtrList< interpolation< scalar > > &vsInterp, const PtrList< interpolation< vector > > &vvInterp, const label UIndex, const label nSubCycle, const scalar trackLength, DynamicList< List< point > > &allPositions, List< DynamicList< scalarList > > &allScalars, List< DynamicList< vectorList > > &allVectors)
Construct from components.
DynamicList< vectorList > & allPositions_
const PtrList< interpolation< scalar > > & vsInterp_
List< DynamicList< vectorList > > & allVectors_
List< DynamicList< scalarList > > & allScalars_
Particle class that samples fields as it passes through. Used in streamline calculation.
static void writeFields(const Cloud< streamLineParticle > &)
Write.
autoPtr< particle > clone() const
Construct and return a clone.
void hitCyclicACMIPatch(streamLineParticleCloud &, trackingData &, const vector &direction)
Overridable function to handle the particle hitting a.
void hitSymmetryPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a.
void hitSymmetryPlanePatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a.
void hitCyclicPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a cyclic.
void hitCyclicAMIPatch(streamLineParticleCloud &, trackingData &, const vector &direction)
Overridable function to handle the particle hitting a.
static void readFields(Cloud< streamLineParticle > &)
Read.
void hitWedgePatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a wedge.
friend Ostream & operator<<(Ostream &, const streamLineParticle &)
void hitProcessorPatch(streamLineParticleCloud &, trackingData &)
bool move(streamLineParticleCloud &, trackingData &, const scalar)
Track all particles to their end point.
bool hitPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a patch.
void hitWallPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a wallPatch.
volScalarField & p
Namespace for OpenFOAM.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces)
Definition: boundaryPatch.C:83
uint8_t direction
Definition: direction.H:56