Geant4 11.2.2
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 "G4Allocator.hh"
47#include "G4ParticleDefinition.hh" // Include from 'particle+matter'
48#include "G4Step.hh"
49#include "G4Track.hh"
50#include "G4TrajectoryPoint.hh" // Include from 'tracking'
51#include "G4VTrajectory.hh"
52#include "G4ios.hh" // Include from 'system'
53#include "globals.hh" // Include from 'global'
54
55#include "trkgdefs.hh"
56#include <stdlib.h> // Include from 'system'
57
58#include <vector>
59
60class G4Polyline;
61
63{
64 using G4TrajectoryPointContainer = std::vector<G4VTrajectoryPoint*>;
65
66 public:
67 // Constructors/Destructor
68
69 G4Trajectory() = default;
70 G4Trajectory(const G4Track* aTrack);
72 ~G4Trajectory() override;
73
74 // Operators
75
76 inline void* operator new(size_t);
77 inline void operator delete(void*);
78 inline G4bool operator==(const G4Trajectory& r) const;
79
80 // Get/Set functions
81
82 inline G4int GetTrackID() const override { return fTrackID; }
83 inline G4int GetParentID() const override { return fParentID; }
84 inline G4String GetParticleName() const override { return ParticleName; }
85 inline G4double GetCharge() const override { return PDGCharge; }
86 inline G4int GetPDGEncoding() const override { return PDGEncoding; }
87 inline G4double GetInitialKineticEnergy() const { return initialKineticEnergy; }
88 inline G4ThreeVector GetInitialMomentum() const override { return initialMomentum; }
89
90 // Other member functions
91
92 void ShowTrajectory(std::ostream& os = G4cout) const override;
93 void DrawTrajectory() const override;
94 void AppendStep(const G4Step* aStep) override;
95 G4int GetPointEntries() const override { return G4int(positionRecord->size()); }
96 G4VTrajectoryPoint* GetPoint(G4int i) const override { return (*positionRecord)[i]; }
97 void MergeTrajectory(G4VTrajectory* secondTrajectory) override;
98
100
101 const std::map<G4String, G4AttDef>* GetAttDefs() const override;
102 std::vector<G4AttValue>* CreateAttValues() const override;
103
104 private:
105 G4TrajectoryPointContainer* positionRecord = nullptr;
106 G4int fTrackID = 0;
107 G4int fParentID = 0;
108 G4int PDGEncoding = 0;
109 G4double PDGCharge = 0.0;
110 G4String ParticleName = "";
111 G4double initialKineticEnergy = 0.0;
112 G4ThreeVector initialMomentum;
113};
114
116
117inline void* G4Trajectory::operator new(size_t)
118{
119 if (aTrajectoryAllocator() == nullptr) {
121 }
122 return (void*)aTrajectoryAllocator()->MallocSingle();
123}
124
125inline void G4Trajectory::operator delete(void* aTrajectory)
126{
127 aTrajectoryAllocator()->FreeSingle((G4Trajectory*)aTrajectory);
128}
129
130inline G4bool G4Trajectory::operator==(const G4Trajectory& r) const { return (this == &r); }
131
132#endif
G4TRACKING_DLL G4Allocator< G4Trajectory > *& aTrajectoryAllocator()
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
G4String GetParticleName() const override
void DrawTrajectory() const override
G4ParticleDefinition * GetParticleDefinition()
void ShowTrajectory(std::ostream &os=G4cout) const override
G4double GetInitialKineticEnergy() const
G4int GetPointEntries() const override
void MergeTrajectory(G4VTrajectory *secondTrajectory) override
G4int GetPDGEncoding() const override
G4VTrajectoryPoint * GetPoint(G4int i) const override
const std::map< G4String, G4AttDef > * GetAttDefs() const override
std::vector< G4AttValue > * CreateAttValues() const override
G4bool operator==(const G4Trajectory &r) const
void AppendStep(const G4Step *aStep) override
G4Trajectory()=default
~G4Trajectory() override
G4double GetCharge() const override
G4ThreeVector GetInitialMomentum() const override
G4int GetParentID() const override
G4int GetTrackID() const override
#define G4TRACKING_DLL
Definition trkgdefs.hh:45