VTK-m  2.0
Particle.h
Go to the documentation of this file.
1 //============================================================================
2 // Copyright (c) Kitware, Inc.
3 // All rights reserved.
4 // See LICENSE.txt for details.
5 //
6 // This software is distributed WITHOUT ANY WARRANTY; without even
7 // the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
8 // PURPOSE. See the above copyright notice for more information.
9 //============================================================================
10 #ifndef vtk_m_Particle_h
11 #define vtk_m_Particle_h
12 
13 #include <ostream>
14 #include <vtkm/Bitset.h>
15 #include <vtkm/VecVariable.h>
16 #include <vtkm/VectorAnalysis.h>
18 
19 namespace vtkm
20 {
21 
22 //Bit field describing the status:
23 class ParticleStatus : public vtkm::Bitset<vtkm::UInt8>
24 {
25 public:
27  {
28  this->SetOk();
29  this->ClearTerminate();
30  }
31 
32  VTKM_EXEC_CONT void SetOk() { this->set(this->SUCCESS_BIT); }
33  VTKM_EXEC_CONT bool CheckOk() const { return this->test(this->SUCCESS_BIT); }
34 
35  VTKM_EXEC_CONT void SetFail() { this->reset(this->SUCCESS_BIT); }
36  VTKM_EXEC_CONT bool CheckFail() const { return !this->test(this->SUCCESS_BIT); }
37 
38  VTKM_EXEC_CONT void SetTerminate() { this->set(this->TERMINATE_BIT); }
40  VTKM_EXEC_CONT bool CheckTerminate() const { return this->test(this->TERMINATE_BIT); }
41 
44  VTKM_EXEC_CONT bool CheckSpatialBounds() const { return this->test(this->SPATIAL_BOUNDS_BIT); }
45 
48  VTKM_EXEC_CONT bool CheckTemporalBounds() const { return this->test(this->TEMPORAL_BOUNDS_BIT); }
49 
52  VTKM_EXEC_CONT bool CheckTookAnySteps() const { return this->test(this->TOOK_ANY_STEPS_BIT); }
53 
56  VTKM_EXEC_CONT bool CheckInGhostCell() const { return this->test(this->IN_GHOST_CELL_BIT); }
57 
60  VTKM_EXEC_CONT bool CheckZeroVelocity() const { return this->test(this->ZERO_VELOCITY); }
61 
63  {
64  return this->CheckOk() && !this->CheckTerminate() && !this->CheckSpatialBounds() &&
65  !this->CheckTemporalBounds() && !this->CheckInGhostCell() && !this->CheckZeroVelocity();
66  }
67 
68 private:
69  static constexpr vtkm::Id SUCCESS_BIT = 0;
70  static constexpr vtkm::Id TERMINATE_BIT = 1;
71  static constexpr vtkm::Id SPATIAL_BOUNDS_BIT = 2;
72  static constexpr vtkm::Id TEMPORAL_BOUNDS_BIT = 3;
73  static constexpr vtkm::Id TOOK_ANY_STEPS_BIT = 4;
74  static constexpr vtkm::Id IN_GHOST_CELL_BIT = 5;
75  static constexpr vtkm::Id ZERO_VELOCITY = 6;
76 };
77 
78 inline VTKM_CONT std::ostream& operator<<(std::ostream& s, const vtkm::ParticleStatus& status)
79 {
80  s << "[ok= " << status.CheckOk();
81  s << " term= " << status.CheckTerminate();
82  s << " spat= " << status.CheckSpatialBounds();
83  s << " temp= " << status.CheckTemporalBounds();
84  s << " ghst= " << status.CheckInGhostCell();
85  s << " zvel= " << status.CheckZeroVelocity();
86  s << "]";
87  return s;
88 }
89 
90 class Particle
91 {
92 public:
94  Particle() {}
95 
98  const vtkm::Id& id,
99  const vtkm::Id& numSteps = 0,
100  const vtkm::ParticleStatus& status = vtkm::ParticleStatus(),
101  const vtkm::FloatDefault& time = 0)
102  : Position(p)
103  , ID(id)
104  , NumSteps(numSteps)
105  , Status(status)
106  , Time(time)
107  {
108  }
109 
112  : Position(p.Position)
113  , ID(p.ID)
114  , NumSteps(p.NumSteps)
115  , Status(p.Status)
116  , Time(p.Time)
117  {
118  }
119 
120  vtkm::Particle& operator=(const vtkm::Particle&) = default;
121 
123  {
124  // This must not be defaulted, since defaulted virtual destructors are
125  // troublesome with CUDA __host__ __device__ markup.
126  }
127 
128  VTKM_EXEC_CONT const vtkm::Vec3f& GetPosition() const { return this->Position; }
129  VTKM_EXEC_CONT void SetPosition(const vtkm::Vec3f& position) { this->Position = position; }
130 
131  VTKM_EXEC_CONT vtkm::Id GetID() const { return this->ID; }
132  VTKM_EXEC_CONT void SetID(vtkm::Id id) { this->ID = id; }
133 
135  VTKM_EXEC_CONT void SetNumberOfSteps(vtkm::Id numSteps) { this->NumSteps = numSteps; }
136 
139  VTKM_EXEC_CONT void SetStatus(vtkm::ParticleStatus status) { this->Status = status; }
140 
141  VTKM_EXEC_CONT vtkm::FloatDefault GetTime() const { return this->Time; }
142  VTKM_EXEC_CONT void SetTime(vtkm::FloatDefault time) { this->Time = time; }
143 
146  const vtkm::FloatDefault& vtkmNotUsed(length))
147  {
148  // Velocity is evaluated from the Velocity field
149  // and is not influenced by the particle
150  VTKM_ASSERT(vectors.GetNumberOfComponents() > 0);
151  return vectors[0];
152  }
153 
156  {
157  (void)deltaT; // unused for a general particle advection case
158  return this->Position;
159  }
160 
161  inline VTKM_CONT friend std::ostream& operator<<(std::ostream& out, const vtkm::Particle& p)
162  {
163  out << "v(" << p.Time << ") = " << p.Position << ", ID: " << p.ID
164  << ", NumSteps: " << p.NumSteps << ", Status: " << p.Status;
165  return out;
166  }
167 
168 private:
170  vtkm::Id ID = -1;
174 
175 public:
176  static size_t Sizeof()
177  {
178  constexpr std::size_t sz = sizeof(vtkm::Vec3f) // Pos
179  + sizeof(vtkm::Id) // ID
180  + sizeof(vtkm::Id) // NumSteps
181  + sizeof(vtkm::UInt8) // Status
182  + sizeof(vtkm::FloatDefault); // Time
183 
184  return sz;
185  }
186 };
187 
189 {
190 public:
193 
195  ChargedParticle(const vtkm::Vec3f& position,
196  const vtkm::Id& id,
197  const vtkm::Float64& mass,
198  const vtkm::Float64& charge,
199  const vtkm::Float64& weighting,
200  const vtkm::Vec3f& momentum,
201  const vtkm::Id& numSteps = 0,
202  const vtkm::ParticleStatus& status = vtkm::ParticleStatus(),
203  const vtkm::FloatDefault& time = 0)
204  : Position(position)
205  , ID(id)
206  , NumSteps(numSteps)
207  , Status(status)
208  , Time(time)
209  , Mass(mass)
210  , Charge(charge)
211  , Weighting(weighting)
212  , Momentum(momentum)
213  {
214  }
215 
216  VTKM_EXEC_CONT const vtkm::Vec3f& GetPosition() const { return this->Position; }
217  VTKM_EXEC_CONT void SetPosition(const vtkm::Vec3f& position) { this->Position = position; }
218 
219  VTKM_EXEC_CONT vtkm::Id GetID() const { return this->ID; }
220  VTKM_EXEC_CONT void SetID(vtkm::Id id) { this->ID = id; }
221 
223  VTKM_EXEC_CONT void SetNumberOfSteps(vtkm::Id numSteps) { this->NumSteps = numSteps; }
224 
227  VTKM_EXEC_CONT void SetStatus(vtkm::ParticleStatus status) { this->Status = status; }
228 
229  VTKM_EXEC_CONT vtkm::FloatDefault GetTime() const { return this->Time; }
230  VTKM_EXEC_CONT void SetTime(vtkm::FloatDefault time) { this->Time = time; }
231 
233  vtkm::Float64 Gamma(vtkm::Vec3f momentum, bool reciprocal = false) const
234  {
236  const vtkm::Float64 fMom2 = vtkm::MagnitudeSquared(momentum);
237  const vtkm::Float64 m2 = this->Mass * this->Mass;
238  const vtkm::Float64 m2_c2_reci = 1.0 / (m2 * c2);
239  if (reciprocal)
240  return vtkm::RSqrt(1.0 + fMom2 * m2_c2_reci);
241  else
242  return vtkm::Sqrt(1.0 + fMom2 * m2_c2_reci);
243  }
244 
247  const vtkm::FloatDefault& length)
248  {
249  VTKM_ASSERT(vectors.GetNumberOfComponents() == 2);
250 
251  // Suppress unused warning
252  (void)this->Weighting;
253 
254  vtkm::Vec3f eField = vectors[0];
255  vtkm::Vec3f bField = vectors[1];
256 
257  const vtkm::Float64 QoM = this->Charge / this->Mass;
258  const vtkm::Vec3f mom_minus = this->Momentum + (0.5 * this->Charge * eField * length);
259 
260  // Get reciprocal of Gamma
261  vtkm::Vec3f gamma_reci = static_cast<vtkm::FloatDefault>(this->Gamma(mom_minus, true));
262  const vtkm::Vec3f t = 0.5 * QoM * length * bField * gamma_reci;
263  const vtkm::Vec3f s = 2.0f * t * (1.0 / (1.0 + vtkm::Magnitude(t)));
264  const vtkm::Vec3f mom_prime = mom_minus + vtkm::Cross(mom_minus, t);
265  const vtkm::Vec3f mom_plus = mom_minus + vtkm::Cross(mom_prime, s);
266 
267  const vtkm::Vec3f mom_new = mom_plus + 0.5 * this->Charge * eField * length;
268  //TODO : Sould this be a const method?
269  // If yes, need a better way to update momentum
270  this->Momentum = mom_new;
271 
272  // momentum = velocity * mass * gamma;
273  // --> velocity = momentum / (mass * gamma)
274  // --> velocity = ( momentum / mass ) * gamma_reci
275  vtkm::Vec3f velocity = (mom_new / this->Mass) * this->Gamma(mom_new, true);
276  return velocity;
277  }
278 
281  {
282  // Translation is in -ve Z direction,
283  // this needs to be a parameter.
284  auto translation = this->NumSteps * deltaT * SPEED_OF_LIGHT * vtkm::Vec3f{ 0., 0., -1.0 };
285  return this->Position + translation;
286  }
287 
288  inline VTKM_CONT friend std::ostream& operator<<(std::ostream& out,
289  const vtkm::ChargedParticle& p)
290  {
291  out << "v(" << p.Time << ") = " << p.Position << ", ID: " << p.ID
292  << ", NumSteps: " << p.NumSteps << ", Status: " << p.Status;
293  return out;
294  }
295 
296 private:
298  vtkm::Id ID = -1;
307  static_cast<vtkm::FloatDefault>(2.99792458e8);
308 
309  friend struct mangled_diy_namespace::Serialization<vtkm::ChargedParticle>;
310 
311 public:
312  static size_t Sizeof()
313  {
314  constexpr std::size_t sz = sizeof(vtkm::Vec3f) // Pos
315  + sizeof(vtkm::Id) // ID
316  + sizeof(vtkm::Id) // NumSteps
317  + sizeof(vtkm::UInt8) // Status
318  + sizeof(vtkm::FloatDefault) // Time
319  + sizeof(vtkm::Float64) //Mass
320  + sizeof(vtkm::Float64) //Charge
321  + sizeof(vtkm::Float64) //Weighting
322  + sizeof(vtkm::Vec3f); //Momentum
323 
324  return sz;
325  }
326 };
327 
328 } //namespace vtkm
329 
330 
332 {
333 template <>
334 struct Serialization<vtkm::Particle>
335 {
336 public:
337  static VTKM_CONT void save(BinaryBuffer& bb, const vtkm::Particle& p)
338  {
339  vtkmdiy::save(bb, p.GetPosition());
340  vtkmdiy::save(bb, p.GetID());
341  vtkmdiy::save(bb, p.GetNumberOfSteps());
342  vtkmdiy::save(bb, p.GetStatus());
343  vtkmdiy::save(bb, p.GetTime());
344  }
345 
346  static VTKM_CONT void load(BinaryBuffer& bb, vtkm::Particle& p)
347  {
348  vtkm::Vec3f pos;
349  vtkmdiy::load(bb, pos);
350  p.SetPosition(pos);
351 
352  vtkm::Id id;
353  vtkmdiy::load(bb, id);
354  p.SetID(id);
355 
356  vtkm::Id numSteps;
357  vtkmdiy::load(bb, numSteps);
358  p.SetNumberOfSteps(numSteps);
359 
360  vtkm::ParticleStatus status;
361  vtkmdiy::load(bb, status);
362  p.SetStatus(status);
363 
364  vtkm::FloatDefault time;
365  vtkmdiy::load(bb, time);
366  p.SetTime(time);
367  }
368 };
369 
370 template <>
371 struct Serialization<vtkm::ChargedParticle>
372 {
373 public:
374  static VTKM_CONT void save(BinaryBuffer& bb, const vtkm::ChargedParticle& e)
375  {
376  vtkmdiy::save(bb, e.Position);
377  vtkmdiy::save(bb, e.ID);
378  vtkmdiy::save(bb, e.NumSteps);
379  vtkmdiy::save(bb, e.Status);
380  vtkmdiy::save(bb, e.Time);
381  vtkmdiy::save(bb, e.Mass);
382  vtkmdiy::save(bb, e.Charge);
383  vtkmdiy::save(bb, e.Weighting);
384  vtkmdiy::save(bb, e.Momentum);
385  }
386 
387  static VTKM_CONT void load(BinaryBuffer& bb, vtkm::ChargedParticle& e)
388  {
389  vtkmdiy::load(bb, e.Position);
390  vtkmdiy::load(bb, e.ID);
391  vtkmdiy::load(bb, e.NumSteps);
392  vtkmdiy::load(bb, e.Status);
393  vtkmdiy::load(bb, e.Time);
394  vtkmdiy::load(bb, e.Mass);
395  vtkmdiy::load(bb, e.Charge);
396  vtkmdiy::load(bb, e.Weighting);
397  vtkmdiy::load(bb, e.Momentum);
398  }
399 };
400 }
401 
402 #endif // vtk_m_Particle_h
vtkm::Bitset< vtkm::UInt8 >::reset
VTKM_EXEC_CONT void reset(vtkm::Id bitIndex)
Definition: Bitset.h:43
vtkm::Bitset
A bitmap to serve different needs.
Definition: Bitset.h:28
Bitset.h
vtkm::ChargedParticle::Gamma
VTKM_EXEC_CONT vtkm::Float64 Gamma(vtkm::Vec3f momentum, bool reciprocal=false) const
Definition: Particle.h:233
vtkm::ChargedParticle::GetEvaluationPosition
VTKM_EXEC_CONT vtkm::Vec3f GetEvaluationPosition(const vtkm::FloatDefault &deltaT) const
Definition: Particle.h:280
vtkm::Sqrt
VTKM_EXEC_CONT vtkm::Float32 Sqrt(vtkm::Float32 x)
Compute the square root of x.
Definition: Math.h:958
vtkm::ParticleStatus::ClearInGhostCell
VTKM_EXEC_CONT void ClearInGhostCell()
Definition: Particle.h:55
vtkm::ChargedParticle::GetStatus
VTKM_EXEC_CONT vtkm::ParticleStatus & GetStatus()
Definition: Particle.h:226
vtkm::ChargedParticle::Weighting
vtkm::Float64 Weighting
Definition: Particle.h:304
mangled_diy_namespace::Serialization< vtkm::ChargedParticle >::load
static VTKM_CONT void load(BinaryBuffer &bb, vtkm::ChargedParticle &e)
Definition: Particle.h:387
vtkm::ChargedParticle::operator<<
VTKM_CONT friend std::ostream & operator<<(std::ostream &out, const vtkm::ChargedParticle &p)
Definition: Particle.h:288
vtkm::MagnitudeSquared
VTKM_EXEC_CONT detail::FloatingPointReturnType< T >::Type MagnitudeSquared(const T &x)
Returns the square of the magnitude of a vector.
Definition: VectorAnalysis.h:64
vtkm::ParticleStatus::CheckTerminate
VTKM_EXEC_CONT bool CheckTerminate() const
Definition: Particle.h:40
vtkm::Particle::Status
vtkm::ParticleStatus Status
Definition: Particle.h:172
vtkm
Groups connected points that have the same field value.
Definition: Atomic.h:19
vtkm::ChargedParticle::Velocity
VTKM_EXEC_CONT vtkm::Vec3f Velocity(const vtkm::VecVariable< vtkm::Vec3f, 2 > &vectors, const vtkm::FloatDefault &length)
Definition: Particle.h:246
vtkm::ParticleStatus::ClearSpatialBounds
VTKM_EXEC_CONT void ClearSpatialBounds()
Definition: Particle.h:43
vtkm::ChargedParticle::ID
vtkm::Id ID
Definition: Particle.h:298
vtkm::Particle::Particle
VTKM_EXEC_CONT Particle(const vtkm::Vec3f &p, const vtkm::Id &id, const vtkm::Id &numSteps=0, const vtkm::ParticleStatus &status=vtkm::ParticleStatus(), const vtkm::FloatDefault &time=0)
Definition: Particle.h:97
VTKM_ASSERT
#define VTKM_ASSERT(condition)
Definition: Assert.h:43
VTKM_EXEC_CONT
#define VTKM_EXEC_CONT
Definition: ExportMacros.h:52
vtkm::ChargedParticle::GetStatus
VTKM_EXEC_CONT vtkm::ParticleStatus GetStatus() const
Definition: Particle.h:225
vtkm::Particle
Definition: Particle.h:90
vtkm::Particle::GetEvaluationPosition
VTKM_EXEC_CONT vtkm::Vec3f GetEvaluationPosition(const vtkm::FloatDefault &deltaT) const
Definition: Particle.h:155
vtkm::Particle::SetTime
VTKM_EXEC_CONT void SetTime(vtkm::FloatDefault time)
Definition: Particle.h:142
vtkm::Bitset< vtkm::UInt8 >::test
VTKM_EXEC_CONT bool test(vtkm::Id bitIndex) const
Definition: Bitset.h:53
vtkm::Particle::operator=
vtkm::Particle & operator=(const vtkm::Particle &)=default
vtkm::Magnitude
VTKM_EXEC_CONT detail::FloatingPointReturnType< T >::Type Magnitude(const T &x)
Returns the magnitude of a vector.
Definition: VectorAnalysis.h:100
vtkm::ChargedParticle::Status
vtkm::ParticleStatus Status
Definition: Particle.h:300
vtkm::ChargedParticle::SetPosition
VTKM_EXEC_CONT void SetPosition(const vtkm::Vec3f &position)
Definition: Particle.h:217
vtkm::ChargedParticle::Momentum
vtkm::Vec3f Momentum
Definition: Particle.h:305
vtkm::ParticleStatus::IN_GHOST_CELL_BIT
static constexpr vtkm::Id IN_GHOST_CELL_BIT
Definition: Particle.h:74
vtkm::ChargedParticle::SetID
VTKM_EXEC_CONT void SetID(vtkm::Id id)
Definition: Particle.h:220
vtkm::ChargedParticle::Position
vtkm::Vec3f Position
Definition: Particle.h:297
vtkm::ParticleStatus::TOOK_ANY_STEPS_BIT
static constexpr vtkm::Id TOOK_ANY_STEPS_BIT
Definition: Particle.h:73
vtkm::ParticleStatus::ParticleStatus
VTKM_EXEC_CONT ParticleStatus()
Definition: Particle.h:26
vtkm::ParticleStatus::CheckSpatialBounds
VTKM_EXEC_CONT bool CheckSpatialBounds() const
Definition: Particle.h:44
vtkm::ParticleStatus::TEMPORAL_BOUNDS_BIT
static constexpr vtkm::Id TEMPORAL_BOUNDS_BIT
Definition: Particle.h:72
vtkm::ChargedParticle::Mass
vtkm::Float64 Mass
Definition: Particle.h:302
mangled_diy_namespace
Definition: Particle.h:331
vtkm::Id
vtkm::Int32 Id
Represents an ID (index into arrays).
Definition: Types.h:191
vtkm::exec::arg::load
VTKM_SUPPRESS_EXEC_WARNINGS VTKM_EXEC T load(const U &u, vtkm::Id v)
Definition: FetchTagArrayDirectIn.h:36
vtkm::ParticleStatus::SetTerminate
VTKM_EXEC_CONT void SetTerminate()
Definition: Particle.h:38
vtkm::ParticleStatus::ClearZeroVelocity
VTKM_EXEC_CONT void ClearZeroVelocity()
Definition: Particle.h:59
vtkm::ParticleStatus::CheckInGhostCell
VTKM_EXEC_CONT bool CheckInGhostCell() const
Definition: Particle.h:56
vtkm::ParticleStatus::CheckOk
VTKM_EXEC_CONT bool CheckOk() const
Definition: Particle.h:33
VectorAnalysis.h
vtkm::ChargedParticle::GetPosition
const VTKM_EXEC_CONT vtkm::Vec3f & GetPosition() const
Definition: Particle.h:216
vtkm::ChargedParticle::Charge
vtkm::Float64 Charge
Definition: Particle.h:303
vtkm::VecVariable::GetNumberOfComponents
VTKM_EXEC_CONT vtkm::IdComponent GetNumberOfComponents() const
Definition: VecVariable.h:53
vtkm::ChargedParticle
Definition: Particle.h:188
vtkm::ParticleStatus::ClearTerminate
VTKM_EXEC_CONT void ClearTerminate()
Definition: Particle.h:39
vtkm::Particle::SetID
VTKM_EXEC_CONT void SetID(vtkm::Id id)
Definition: Particle.h:132
vtkm::Particle::GetStatus
VTKM_EXEC_CONT vtkm::ParticleStatus & GetStatus()
Definition: Particle.h:138
vtkm::ParticleStatus::ZERO_VELOCITY
static constexpr vtkm::Id ZERO_VELOCITY
Definition: Particle.h:75
vtkm::ChargedParticle::SetTime
VTKM_EXEC_CONT void SetTime(vtkm::FloatDefault time)
Definition: Particle.h:230
vtkm::ChargedParticle::Sizeof
static size_t Sizeof()
Definition: Particle.h:312
vtkm::VecVariable
A short variable-length array with maximum length.
Definition: VecVariable.h:30
vtkm::ParticleStatus::SetFail
VTKM_EXEC_CONT void SetFail()
Definition: Particle.h:35
vtkm::Cross
VTKM_EXEC_CONT vtkm::Vec< typename detail::FloatingPointReturnType< T >::Type, 3 > Cross(const vtkm::Vec< T, 3 > &x, const vtkm::Vec< T, 3 > &y)
Find the cross product of two vectors.
Definition: VectorAnalysis.h:177
Serialization.h
vtkm::ChargedParticle::SetNumberOfSteps
VTKM_EXEC_CONT void SetNumberOfSteps(vtkm::Id numSteps)
Definition: Particle.h:223
mangled_diy_namespace::Serialization< vtkm::Particle >::load
static VTKM_CONT void load(BinaryBuffer &bb, vtkm::Particle &p)
Definition: Particle.h:346
vtkm::Particle::GetTime
VTKM_EXEC_CONT vtkm::FloatDefault GetTime() const
Definition: Particle.h:141
vtkm::ParticleStatus::CheckFail
VTKM_EXEC_CONT bool CheckFail() const
Definition: Particle.h:36
vtkm::Particle::Sizeof
static size_t Sizeof()
Definition: Particle.h:176
vtkm::operator<<
VTKM_CONT std::ostream & operator<<(std::ostream &stream, const vtkm::Bounds &bounds)
Helper function for printing bounds during testing.
Definition: Bounds.h:237
vtkm::ChargedParticle::SPEED_OF_LIGHT
constexpr static vtkm::FloatDefault SPEED_OF_LIGHT
Definition: Particle.h:306
VTKM_CONT
#define VTKM_CONT
Definition: ExportMacros.h:57
vtkm::Particle::~Particle
VTKM_EXEC_CONT ~Particle() noexcept
Definition: Particle.h:122
vtkm::ParticleStatus::SetInGhostCell
VTKM_EXEC_CONT void SetInGhostCell()
Definition: Particle.h:54
vtkm::Particle::Velocity
VTKM_EXEC_CONT vtkm::Vec3f Velocity(const vtkm::VecVariable< vtkm::Vec3f, 2 > &vectors, const vtkm::FloatDefault &vtkmNotUsed(length))
Definition: Particle.h:145
vtkm::UInt8
uint8_t UInt8
Definition: Types.h:157
vtkm::ParticleStatus
Definition: Particle.h:23
vtkm::Particle::GetStatus
VTKM_EXEC_CONT vtkm::ParticleStatus GetStatus() const
Definition: Particle.h:137
vtkmNotUsed
#define vtkmNotUsed(parameter_name)
Simple macro to identify a parameter as unused.
Definition: ExportMacros.h:128
vtkm::ChargedParticle::ChargedParticle
VTKM_EXEC_CONT ChargedParticle(const vtkm::Vec3f &position, const vtkm::Id &id, const vtkm::Float64 &mass, const vtkm::Float64 &charge, const vtkm::Float64 &weighting, const vtkm::Vec3f &momentum, const vtkm::Id &numSteps=0, const vtkm::ParticleStatus &status=vtkm::ParticleStatus(), const vtkm::FloatDefault &time=0)
Definition: Particle.h:195
vtkm::ParticleStatus::SUCCESS_BIT
static constexpr vtkm::Id SUCCESS_BIT
Definition: Particle.h:69
vtkm::Particle::operator<<
VTKM_CONT friend std::ostream & operator<<(std::ostream &out, const vtkm::Particle &p)
Definition: Particle.h:161
vtkm::Particle::Position
vtkm::Vec3f Position
Definition: Particle.h:169
vtkm::Particle::Time
vtkm::FloatDefault Time
Definition: Particle.h:173
mangled_diy_namespace::Serialization< vtkm::Particle >::save
static VTKM_CONT void save(BinaryBuffer &bb, const vtkm::Particle &p)
Definition: Particle.h:337
vtkm::Particle::SetNumberOfSteps
VTKM_EXEC_CONT void SetNumberOfSteps(vtkm::Id numSteps)
Definition: Particle.h:135
vtkm::ChargedParticle::GetNumberOfSteps
VTKM_EXEC_CONT vtkm::Id GetNumberOfSteps() const
Definition: Particle.h:222
vtkm::ChargedParticle::NumSteps
vtkm::Id NumSteps
Definition: Particle.h:299
vtkm::ChargedParticle::GetTime
VTKM_EXEC_CONT vtkm::FloatDefault GetTime() const
Definition: Particle.h:229
vtkm::Vec3f
vtkm::Vec< vtkm::FloatDefault, 3 > Vec3f
Vec3f corresponds to a 3-dimensional vector of floating point values.
Definition: Types.h:1014
vtkm::ParticleStatus::CheckZeroVelocity
VTKM_EXEC_CONT bool CheckZeroVelocity() const
Definition: Particle.h:60
vtkm::Particle::SetStatus
VTKM_EXEC_CONT void SetStatus(vtkm::ParticleStatus status)
Definition: Particle.h:139
vtkm::ParticleStatus::ClearTookAnySteps
VTKM_EXEC_CONT void ClearTookAnySteps()
Definition: Particle.h:51
vtkm::ParticleStatus::CheckTemporalBounds
VTKM_EXEC_CONT bool CheckTemporalBounds() const
Definition: Particle.h:48
vtkm::Vec< vtkm::FloatDefault, 3 >
vtkm::FloatDefault
vtkm::Float32 FloatDefault
The floating point type to use when no other precision is specified.
Definition: Types.h:198
VecVariable.h
vtkm::Particle::ID
vtkm::Id ID
Definition: Particle.h:170
vtkm::ParticleStatus::SPATIAL_BOUNDS_BIT
static constexpr vtkm::Id SPATIAL_BOUNDS_BIT
Definition: Particle.h:71
vtkm::ParticleStatus::TERMINATE_BIT
static constexpr vtkm::Id TERMINATE_BIT
Definition: Particle.h:70
vtkm::Particle::Particle
VTKM_EXEC_CONT Particle(const vtkm::Particle &p)
Definition: Particle.h:111
vtkm::Particle::GetID
VTKM_EXEC_CONT vtkm::Id GetID() const
Definition: Particle.h:131
vtkm::ChargedParticle::Time
vtkm::FloatDefault Time
Definition: Particle.h:301
vtkm::Bitset< vtkm::UInt8 >::set
VTKM_EXEC_CONT void set(vtkm::Id bitIndex)
Definition: Bitset.h:30
vtkm::Particle::GetPosition
const VTKM_EXEC_CONT vtkm::Vec3f & GetPosition() const
Definition: Particle.h:128
vtkm::Particle::GetNumberOfSteps
VTKM_EXEC_CONT vtkm::Id GetNumberOfSteps() const
Definition: Particle.h:134
vtkm::Float64
double Float64
Definition: Types.h:155
vtkm::ParticleStatus::SetTemporalBounds
VTKM_EXEC_CONT void SetTemporalBounds()
Definition: Particle.h:46
vtkm::ChargedParticle::SetStatus
VTKM_EXEC_CONT void SetStatus(vtkm::ParticleStatus status)
Definition: Particle.h:227
vtkm::Particle::Particle
VTKM_EXEC_CONT Particle()
Definition: Particle.h:94
vtkm::ParticleStatus::CanContinue
VTKM_EXEC_CONT bool CanContinue() const
Definition: Particle.h:62
vtkm::ParticleStatus::SetSpatialBounds
VTKM_EXEC_CONT void SetSpatialBounds()
Definition: Particle.h:42
vtkm::ParticleStatus::SetOk
VTKM_EXEC_CONT void SetOk()
Definition: Particle.h:32
vtkm::ParticleStatus::CheckTookAnySteps
VTKM_EXEC_CONT bool CheckTookAnySteps() const
Definition: Particle.h:52
vtkm::Particle::SetPosition
VTKM_EXEC_CONT void SetPosition(const vtkm::Vec3f &position)
Definition: Particle.h:129
vtkm::ChargedParticle::ChargedParticle
VTKM_EXEC_CONT ChargedParticle()
Definition: Particle.h:192
vtkm::ParticleStatus::ClearTemporalBounds
VTKM_EXEC_CONT void ClearTemporalBounds()
Definition: Particle.h:47
vtkm::ParticleStatus::SetZeroVelocity
VTKM_EXEC_CONT void SetZeroVelocity()
Definition: Particle.h:58
vtkm::ParticleStatus::SetTookAnySteps
VTKM_EXEC_CONT void SetTookAnySteps()
Definition: Particle.h:50
vtkm::Particle::NumSteps
vtkm::Id NumSteps
Definition: Particle.h:171
mangled_diy_namespace::Serialization< vtkm::ChargedParticle >::save
static VTKM_CONT void save(BinaryBuffer &bb, const vtkm::ChargedParticle &e)
Definition: Particle.h:374
vtkm::ChargedParticle::GetID
VTKM_EXEC_CONT vtkm::Id GetID() const
Definition: Particle.h:219