BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
|
Trajectory information from track including last scatter etc. More...
#include <BDSTrajectory.hh>
Public Member Functions | |
BDSTrajectory (const G4Track *aTrack, G4bool interactiveIn, const BDS::TrajectoryOptions &storageOptionsIn) | |
BDSTrajectory (BDSTrajectory &)=delete | |
copy constructor is not needed | |
void * | operator new (size_t) |
void | operator delete (void *) |
int | operator== (const BDSTrajectory &right) const |
virtual void | AppendStep (const G4Step *aStep) |
void | AppendStep (const BDSTrajectoryPoint *pointIn) |
void | CleanPoint (BDSTrajectoryPoint *point) const |
virtual void | MergeTrajectory (G4VTrajectory *secondTrajectory) |
Merge another trajectory into this one. | |
virtual G4VTrajectoryPoint * | GetPoint (G4int i) const |
Access a point - use this class's container. | |
virtual int | GetPointEntries () const |
Get number of trajectory points in this trajectory. | |
virtual G4bool | IsPrimary () const |
Method to identify which one is a primary. Overridden in derived class. | |
void | SetTrajIndex (G4int trajIndexIn) |
G4int | GetTrajIndex () const |
void | SetParentIndex (G4int parentIndexIn) |
G4int | GetParentIndex () const |
void | SetParentStepIndex (G4int parentStepIndexIn) |
The index of the step along the parent trajectory from which this one was created. | |
G4int | GetParentStepIndex () const |
G4int | GetDepth () const |
Depth in the tree. Will be filled later once all trajectories are created and sorted. | |
void | SetDepth (G4int depthIn) |
void | SetParent (BDSTrajectory *parentIn) |
Record the parent trajectory. | |
BDSTrajectory * | GetParent () const |
G4int | GetCreatorProcessType () const |
G4int | GetCreatorProcessSubType () const |
BDSTrajectoryPoint * | FirstInteraction () const |
BDSTrajectoryPoint * | LastInteraction () const |
Protected Attributes | |
G4int | creatorProcessType |
G4int | creatorProcessSubType |
G4double | weight |
G4bool | interactive |
G4bool | suppressTransportationAndNotInteractive |
const BDS::TrajectoryOptions | storageOptions |
BDSTrajectory * | parent |
G4int | trajIndex |
G4int | parentIndex |
G4int | parentStepIndex |
G4int | depth |
BDSTrajectoryPointsContainer * | fpBDSPointsContainer |
Friends | |
std::ostream & | operator<< (std::ostream &out, BDSTrajectory const &t) |
Output stream. | |
Trajectory information from track including last scatter etc.
BDSTrajectory stores BDSTrajectoryPoints
Definition at line 43 of file BDSTrajectory.hh.
BDSTrajectory::BDSTrajectory | ( | const G4Track * | aTrack, |
G4bool | interactiveIn, | ||
const BDS::TrajectoryOptions & | storageOptionsIn | ||
) |
Definition at line 35 of file BDSTrajectory.cc.
|
virtual |
Definition at line 70 of file BDSTrajectory.cc.
void BDSTrajectory::AppendStep | ( | const BDSTrajectoryPoint * | pointIn | ) |
Append a step point. Use a pre-made BDSTrajectoryPoint to save creating it again, which involves coordinate transforms.
Definition at line 79 of file BDSTrajectory.cc.
References CleanPoint(), fpBDSPointsContainer, BDSTrajectoryPoint::GetMaterial(), and BDSTrajectoryPoint::NotTransportationLimitedStep().
|
virtual |
Append a step point to this trajectory. This is required for the trajectory points to show up in the visualisation correctly.
Reimplemented in BDSTrajectoryPrimary.
Definition at line 112 of file BDSTrajectory.cc.
References fpBDSPointsContainer.
Referenced by BDSTrajectoryPrimary::AppendStep().
void BDSTrajectory::CleanPoint | ( | BDSTrajectoryPoint * | point | ) | const |
Apply the options to store a trajectory extra delete them if not storing them. Used to compensate when copying a fuller primary trajectory point that will always have extra information, but may not be needed when appending to the primary trajectory.
Definition at line 102 of file BDSTrajectory.cc.
Referenced by AppendStep().
BDSTrajectoryPoint * BDSTrajectory::FirstInteraction | ( | ) | const |
Find the first point in a trajectory where the post step process isn't fTransportation AND the post step process isn't fGeneral in combination with the post step process subtype isn't step_limiter.
Definition at line 166 of file BDSTrajectory.cc.
References GetPoint(), GetPointEntries(), and BDSTrajectoryPoint::IsScatteringPoint().
|
inline |
Definition at line 106 of file BDSTrajectory.hh.
|
inline |
Definition at line 105 of file BDSTrajectory.hh.
|
inline |
Depth in the tree. Will be filled later once all trajectories are created and sorted.
Definition at line 99 of file BDSTrajectory.hh.
|
inline |
Definition at line 104 of file BDSTrajectory.hh.
|
inline |
Definition at line 92 of file BDSTrajectory.hh.
|
inline |
Definition at line 96 of file BDSTrajectory.hh.
|
inlinevirtual |
Access a point - use this class's container.
Definition at line 76 of file BDSTrajectory.hh.
References fpBDSPointsContainer.
Referenced by BDSOutputROOTEventTrajectory::FillIndividualTrajectory(), FirstInteraction(), and BDSEventAction::IdentifyTrajectoriesForStorage().
|
inlinevirtual |
Get number of trajectory points in this trajectory.
Definition at line 79 of file BDSTrajectory.hh.
References fpBDSPointsContainer.
Referenced by FirstInteraction(), BDSEventAction::IdentifyTrajectoriesForStorage(), and MergeTrajectory().
|
inline |
Definition at line 87 of file BDSTrajectory.hh.
|
inlinevirtual |
Method to identify which one is a primary. Overridden in derived class.
Reimplemented in BDSTrajectoryPrimary.
Definition at line 82 of file BDSTrajectory.hh.
BDSTrajectoryPoint * BDSTrajectory::LastInteraction | ( | ) | const |
Definition at line 181 of file BDSTrajectory.cc.
|
virtual |
Merge another trajectory into this one.
Reimplemented in BDSTrajectoryPrimary.
Definition at line 152 of file BDSTrajectory.cc.
References fpBDSPointsContainer, and GetPointEntries().
Referenced by BDSTrajectoryPrimary::MergeTrajectory().
|
inline |
Definition at line 145 of file BDSTrajectory.hh.
|
inline |
Definition at line 138 of file BDSTrajectory.hh.
|
inline |
Definition at line 57 of file BDSTrajectory.hh.
|
inline |
Definition at line 100 of file BDSTrajectory.hh.
|
inline |
Record the parent trajectory.
Definition at line 103 of file BDSTrajectory.hh.
Referenced by BDSEventAction::IdentifyTrajectoriesForStorage().
|
inline |
Record the TrajIndex (i.e. output index) of the trajectory of the parent trajectory for this one.
Definition at line 91 of file BDSTrajectory.hh.
|
inline |
The index of the step along the parent trajectory from which this one was created.
Definition at line 95 of file BDSTrajectory.hh.
|
inline |
The index of the trajectory assigned in the output from the reduced set of indices. This is why it will not be the same as the track ID.
Definition at line 86 of file BDSTrajectory.hh.
|
friend |
Output stream.
Definition at line 196 of file BDSTrajectory.cc.
|
protected |
Definition at line 119 of file BDSTrajectory.hh.
|
protected |
Definition at line 118 of file BDSTrajectory.hh.
|
protected |
Definition at line 128 of file BDSTrajectory.hh.
|
protected |
Container of all points. This is really a vector so all memory is dynamically allocated and there's no need to make this dynamically allocated itself a la all Geant4 examples.
Definition at line 133 of file BDSTrajectory.hh.
Referenced by AppendStep(), GetPoint(), GetPointEntries(), and MergeTrajectory().
|
protected |
Definition at line 121 of file BDSTrajectory.hh.
|
protected |
Definition at line 124 of file BDSTrajectory.hh.
|
protected |
Definition at line 126 of file BDSTrajectory.hh.
|
protected |
Definition at line 127 of file BDSTrajectory.hh.
|
protected |
Definition at line 123 of file BDSTrajectory.hh.
|
protected |
Definition at line 122 of file BDSTrajectory.hh.
|
protected |
Definition at line 125 of file BDSTrajectory.hh.
|
protected |
Definition at line 120 of file BDSTrajectory.hh.