19#ifndef BDSTRAJECTORY_H
20#define BDSTRAJECTORY_H
21#include "BDSTrajectoryOptions.hh"
22#include "BDSTrajectoryPoint.hh"
23#include "G4Trajectory.hh"
30class G4TrajectoryContainer;
31class G4VTrajectoryPoint;
33typedef std::vector<BDSTrajectoryPoint*> BDSTrajectoryPointsContainer;
55 inline void*
operator new(size_t);
56 inline void operator delete(
void*);
57 inline int operator == (
const BDSTrajectory& right)
const {
return (
this==&right);}
86 inline void SetTrajIndex(G4int trajIndexIn) {trajIndex = trajIndexIn;}
87 inline G4int GetTrajIndex()
const {
return trajIndex;}
91 inline void SetParentIndex(G4int parentIndexIn) {parentIndex = parentIndexIn;}
92 inline G4int GetParentIndex()
const {
return parentIndex;}
96 inline G4int GetParentStepIndex()
const {
return parentStepIndex;}
100 inline void SetDepth(G4int depthIn) {depth = depthIn;}
105 inline G4int GetCreatorProcessType()
const {
return creatorProcessType;}
106 inline G4int GetCreatorProcessSubType()
const {
return creatorProcessSubType;}
118 G4int creatorProcessType;
119 G4int creatorProcessSubType;
122 G4bool suppressTransportationAndNotInteractive;
127 G4int parentStepIndex;
136extern G4Allocator<BDSTrajectory> bdsTrajectoryAllocator;
138inline void* BDSTrajectory::operator
new(size_t)
141 aTrajectory = (
void*)bdsTrajectoryAllocator.MallocSingle();
145inline void BDSTrajectory::operator
delete(
void* aTrajectory)
146{bdsTrajectoryAllocator.FreeSingle((
BDSTrajectory*)aTrajectory);}
A Point in a trajectory with extra information.
Trajectory information from track including last scatter etc.
virtual int GetPointEntries() const
Get number of trajectory points in this trajectory.
BDSTrajectoryPoint * FirstInteraction() const
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)
Merge another trajectory into this one.
void SetTrajIndex(G4int trajIndexIn)
void SetParentStepIndex(G4int parentStepIndexIn)
The index of the step along the parent trajectory from which this one was created.
virtual G4VTrajectoryPoint * GetPoint(G4int i) const
Access a point - use this class's container.
BDSTrajectory(BDSTrajectory &)=delete
copy constructor is not needed
friend std::ostream & operator<<(std::ostream &out, BDSTrajectory const &t)
Output stream.
void SetParent(BDSTrajectory *parentIn)
Record the parent trajectory.
void CleanPoint(BDSTrajectoryPoint *point) const
virtual void AppendStep(const G4Step *aStep)
virtual G4bool IsPrimary() const
Method to identify which one is a primary. Overridden in derived class.
void SetParentIndex(G4int parentIndexIn)
BDSTrajectoryPointsContainer * fpBDSPointsContainer
G4int GetDepth() const
Depth in the tree. Will be filled later once all trajectories are created and sorted.