Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Trajectory.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4Trajectory
27//
28// Class description:
29//
30// This class represents the trajectory of a particle being tracked.
31// It includes information of:
32// 1) List of trajectory points which compose the trajectory;
33// 2) Static information of the particle which generated the
34// trajectory;
35// 3) Track ID and parent particle ID of the trajectory.
36
37// Contact:
38// Questions and comments to this code should be sent to
39// Katsuya Amako (e-mail: [email protected])
40// Makoto Asai (e-mail: [email protected])
41// Takashi Sasaki (e-mail: [email protected])
42// --------------------------------------------------------------------
43#ifndef G4Trajectory_hh
44#define G4Trajectory_hh 1
45
46#include <stdlib.h> // Include from 'system'
47#include <vector>
48
49#include "trkgdefs.hh"
50#include "G4VTrajectory.hh"
51#include "G4Allocator.hh"
52#include "G4ios.hh" // Include from 'system'
53#include "globals.hh" // Include from 'global'
54#include "G4ParticleDefinition.hh" // Include from 'particle+matter'
55#include "G4TrajectoryPoint.hh" // Include from 'tracking'
56#include "G4Track.hh"
57#include "G4Step.hh"
58
59class G4Polyline;
60
62{
63
64 using G4TrajectoryPointContainer = std::vector<G4VTrajectoryPoint*>;
65
66 public:
67
68 // Constructors/Destructor
69
71 G4Trajectory(const G4Track* aTrack);
73 virtual ~G4Trajectory();
74
75 // Operators
76
77 inline void* operator new(size_t);
78 inline void operator delete(void*);
79 inline G4int operator == (const G4Trajectory& r) const;
80
81 // Get/Set functions
82
83 inline G4int GetTrackID() const
84 { return fTrackID; }
85 inline G4int GetParentID() const
86 { return fParentID; }
88 { return ParticleName; }
89 inline G4double GetCharge() const
90 { return PDGCharge; }
91 inline G4int GetPDGEncoding() const
92 { return PDGEncoding; }
94 { return initialKineticEnergy; }
96 { return initialMomentum; }
97
98 // Other member functions
99
100 virtual void ShowTrajectory(std::ostream& os=G4cout) const;
101 virtual void DrawTrajectory() const;
102 virtual void AppendStep(const G4Step* aStep);
103 virtual G4int GetPointEntries() const
104 { return G4int(positionRecord->size()); }
106 { return (*positionRecord)[i]; }
107 virtual void MergeTrajectory(G4VTrajectory* secondTrajectory);
108
110
111 virtual const std::map<G4String,G4AttDef>* GetAttDefs() const;
112 virtual std::vector<G4AttValue>* CreateAttValues() const;
113
114 private:
115
116 G4TrajectoryPointContainer* positionRecord = nullptr;
117 G4int fTrackID = 0;
118 G4int fParentID = 0;
119 G4int PDGEncoding = 0;
120 G4double PDGCharge = 0.0;
121 G4String ParticleName = "";
122 G4double initialKineticEnergy = 0.0;
123 G4ThreeVector initialMomentum;
124};
125
127
128inline void* G4Trajectory::operator new(size_t)
129{
130 if (aTrajectoryAllocator() == nullptr)
131 {
133 }
134 return (void*)aTrajectoryAllocator()->MallocSingle();
135}
136
137inline void G4Trajectory::operator delete(void* aTrajectory)
138{
139 aTrajectoryAllocator()->FreeSingle((G4Trajectory*)aTrajectory);
140}
141
143{
144 return (this==&r);
145}
146
147#endif
G4TRACKING_DLL G4Allocator< G4Trajectory > *& aTrajectoryAllocator()
Definition: G4Trajectory.cc:49
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
Definition: G4Step.hh:62
G4ParticleDefinition * GetParticleDefinition()
G4int GetTrackID() const
Definition: G4Trajectory.hh:83
virtual G4int GetPointEntries() const
virtual G4VTrajectoryPoint * GetPoint(G4int i) const
G4double GetInitialKineticEnergy() const
Definition: G4Trajectory.hh:93
virtual void ShowTrajectory(std::ostream &os=G4cout) const
G4String GetParticleName() const
Definition: G4Trajectory.hh:87
virtual void AppendStep(const G4Step *aStep)
virtual void DrawTrajectory() const
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)
G4double GetCharge() const
Definition: G4Trajectory.hh:89
G4int operator==(const G4Trajectory &r) const
G4int GetParentID() const
Definition: G4Trajectory.hh:85
G4int GetPDGEncoding() const
Definition: G4Trajectory.hh:91
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual ~G4Trajectory()
Definition: G4Trajectory.cc:96
G4ThreeVector GetInitialMomentum() const
Definition: G4Trajectory.hh:95
#define G4TRACKING_DLL
Definition: trkgdefs.hh:45