Loading...
Searching...
No Matches
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 Copyright (C) 2022-2024 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
27Class
28 Foam::streamLineParticle
29
30Description
31 Particle class that samples fields as it passes through.
32 Used in streamline calculation.
33
34SourceFiles
35 streamLineParticle.C
36
37\*---------------------------------------------------------------------------*/
38
39#ifndef Foam_streamLineParticle_H
40#define Foam_streamLineParticle_H
41
42#include "particle.H"
43#include "autoPtr.H"
44#include "interpolation.H"
45#include "vectorList.H"
46
47// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48
49namespace Foam
50{
51
52// Forward Declarations
55
58/*---------------------------------------------------------------------------*\
59 Class streamLineParticle Declaration
60\*---------------------------------------------------------------------------*/
61
64 public particle
65{
66public:
67
68 class trackingData
69 :
71 {
72 public:
74 // Public Data
79 const label nSubCycle_;
80 const scalar trackLength_;
81
85
86
87 // Constructors
88
89 //- Construct from components
91 (
93 const PtrList<interpolation<scalar>>& vsInterp,
94 const PtrList<interpolation<vector>>& vvInterp,
95 const interpolation<vector>& UInterp,
96 const label nSubCycle,
97 const scalar trackLength,
98 DynamicList<List<point>>& allPositions,
99 List<DynamicList<scalarList>>& allScalars,
100 List<DynamicList<vectorList>>& allVectors
101 )
102 :
104 vsInterp_(vsInterp),
105 vvInterp_(vvInterp),
106 UInterp_(UInterp),
107 nSubCycle_(nSubCycle),
108 trackLength_(trackLength),
109 allPositions_(allPositions),
110 allScalars_(allScalars),
111 allVectors_(allVectors)
112 {}
113 };
114
115
116private:
117
118 // Private Data
119
120 //- Whether particle transports with +U or -U
121 bool trackForward_;
122
123 //- Lifetime of particle. Particle dies when reaches 0.
124 label lifeTime_;
125
126 //- Sampled positions
127 DynamicList<point> sampledPositions_;
128
129 //- Sampled scalars
130 List<DynamicList<scalar>> sampledScalars_;
131
132 //- Sampled vectors
133 List<DynamicList<vector>> sampledVectors_;
134
135
136 // Private Member Functions
137
138 //- Interpolate quantities; return interpolated velocity.
139 vector interpolateFields
140 (
141 const trackingData& td,
142 const barycentric& tetCoords,
143 const tetIndices& tetIs
144 );
145
146
147public:
148
149 // Constructors
150
151 //- Construct from components
153 (
154 const polyMesh& mesh,
155 const vector& position,
156 const label celli,
157 const bool trackForward,
158 const label lifeTime
159 );
160
161 //- Construct from Istream
163 (
164 const polyMesh& mesh,
165 Istream& is,
166 bool readFields = true,
167 bool newFormat = true
168 );
169
170 //- Construct copy
172
173 //- Return a clone
174 virtual autoPtr<particle> clone() const
175 {
176 return particle::Clone(*this);
177 }
178
179 //- Factory class to read-construct particles used for parallel transfer
180 class iNew
181 {
182 const polyMesh& mesh_;
183
184 public:
185
186 iNew(const polyMesh& mesh)
187 :
188 mesh_(mesh)
189 {}
190
195 new streamLineParticle(mesh_, is, true)
196 );
197 }
198 };
200
201 // Member Functions
202
203 // Tracking
204
205 //- Track all particles to their end point
206 bool move(streamLineParticleCloud&, trackingData&, const scalar);
207
208 //- Overridable function to handle the particle hitting a patch
209 // Executed before other patch-hitting functions
211
212 //- Overridable function to handle the particle hitting a wedge
214
215 //- Overridable function to handle the particle hitting a
216 // symmetry plane
218
219 //- Overridable function to handle the particle hitting a
220 // symmetry patch
222
223 //- Overridable function to handle the particle hitting a cyclic
225
226 //- Overridable function to handle the particle hitting a
227 // cyclicAMIPatch
229 (
232 const vector& direction
233 );
234
235 //- Overridable function to handle the particle hitting a
236 // cyclicACMIPatch
238 (
241 const vector& direction
242 );
243
244 //- Overridable function to handle the particle hitting a
245 //- processorPatch
247
248 //- Overridable function to handle the particle hitting a wallPatch
250
251
252 // I-O
253
254 //- Read
256
257 //- Write
258 static void writeFields(const Cloud<streamLineParticle>&);
259
260
261 // Ostream Operator
262
264};
265
266
267// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
268
269} // End namespace Foam
270
271// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
272
273#endif
274
275// ************************************************************************* //
Base cloud calls templated on particle type.
Definition Cloud.H:64
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
Abstract base class for volume field interpolation.
vector position() const
Return current particle position.
Definition particleI.H:283
static autoPtr< particle > Clone(const Derived &p)
Clone a particle.
Definition particle.H:552
const polyMesh & mesh() const noexcept
Return the mesh database.
Definition particleI.H:110
particle(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from components.
Definition particle.C:507
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A Cloud of streamLine particles.
autoPtr< streamLineParticle > operator()(Istream &is) const
const PtrList< interpolation< vector > > & vvInterp_
const interpolation< vector > & UInterp_
trackingData(streamLineParticleCloud &cloud, const PtrList< interpolation< scalar > > &vsInterp, const PtrList< interpolation< vector > > &vvInterp, const interpolation< vector > &UInterp, 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.
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.
virtual autoPtr< particle > clone() const
Return a clone.
void hitCyclicAMIPatch(streamLineParticleCloud &, trackingData &, const vector &direction)
Overridable function to handle the particle hitting a.
static void readFields(Cloud< streamLineParticle > &)
Read.
streamLineParticle(const polyMesh &mesh, const vector &position, const label celli, const bool trackForward, const label lifeTime)
Construct from components.
void hitWedgePatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a wedge.
friend Ostream & operator<<(Ostream &, const streamLineParticle &)
void hitProcessorPatch(streamLineParticleCloud &, trackingData &)
Overridable function to handle the particle hitting a processorPatch.
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.
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition tetIndices.H:79
volScalarField & p
dynamicFvMesh & mesh
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for OpenFOAM.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
uint8_t direction
Definition direction.H:49
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
Vector< scalar > vector
Definition vector.H:57