Geant4 9.6.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//
27// $Id$
28//
29//---------------------------------------------------------------
30//
31// G4Trajectory.hh
32//
33// class description:
34// This class represents the trajectory of a particle tracked.
35// It includes information of
36// 1) List of trajectory points which compose the trajectory,
37// 2) static information of particle which generated the
38// trajectory,
39// 3) trackID and parent particle ID of the trajectory,
40//
41// Contact:
42// Questions and comments to this code should be sent to
43// Katsuya Amako (e-mail: [email protected])
44// Makoto Asai (e-mail: [email protected])
45// Takashi Sasaki (e-mail: [email protected])
46//
47// ---------------------------------------------------------------
48
49class G4Trajectory;
50
51#ifndef G4Trajectory_h
52#define G4Trajectory_h 1
53
54#include "G4VTrajectory.hh"
55#include "G4Allocator.hh"
56#include <stdlib.h> // Include from 'system'
57#include "G4ios.hh" // Include from 'system'
58#include <vector> // G4RWTValOrderedVector
59#include "globals.hh" // Include from 'global'
60#include "G4ParticleDefinition.hh" // Include from 'particle+matter'
61#include "G4TrajectoryPoint.hh" // Include from 'tracking'
62#include "G4Track.hh"
63#include "G4Step.hh"
64
65class G4Polyline; // Forward declaration.
66
67typedef std::vector<G4VTrajectoryPoint*> TrajectoryPointContainer;
68///////////////////
70///////////////////
71{
72
73//--------
74public: // with description
75//--------
76
77// Constructor/Destrcutor
78
80
81 G4Trajectory(const G4Track* aTrack);
83 virtual ~G4Trajectory();
84
85// Operators
86 inline void* operator new(size_t);
87 inline void operator delete(void*);
88 inline int operator == (const G4Trajectory& right) const
89 {return (this==&right);}
90
91// Get/Set functions
92 inline G4int GetTrackID() const
93 { return fTrackID; }
94 inline G4int GetParentID() const
95 { return fParentID; }
97 { return ParticleName; }
98 inline G4double GetCharge() const
99 { return PDGCharge; }
100 inline G4int GetPDGEncoding() const
101 { return PDGEncoding; }
103 { return initialKineticEnergy; }
105 { return initialMomentum; }
106
107// Other member functions
108 virtual void ShowTrajectory(std::ostream& os=G4cout) const;
109 //virtual void DrawTrajectory() const;
110 virtual void DrawTrajectory(G4int i_mode = 0) const;
111 virtual void AppendStep(const G4Step* aStep);
112 virtual int GetPointEntries() const { return positionRecord->size(); }
114 { return (*positionRecord)[i]; }
115 virtual void MergeTrajectory(G4VTrajectory* secondTrajectory);
116
118
119 virtual const std::map<G4String,G4AttDef>* GetAttDefs() const;
120 virtual std::vector<G4AttValue>* CreateAttValues() const;
121
122//---------
123 private:
124//---------
125
126 TrajectoryPointContainer* positionRecord;
127 G4int fTrackID;
128 G4int fParentID;
129 G4int PDGEncoding;
130 G4double PDGCharge;
131 G4String ParticleName;
132 G4double initialKineticEnergy;
133 G4ThreeVector initialMomentum;
134
135};
136
137#if defined G4TRACKING_ALLOC_EXPORT
139#else
141#endif
142
143inline void* G4Trajectory::operator new(size_t)
144{
145 void* aTrajectory;
146 aTrajectory = (void*)aTrajectoryAllocator.MallocSingle();
147 return aTrajectory;
148}
149
150inline void G4Trajectory::operator delete(void* aTrajectory)
151{
152 aTrajectoryAllocator.FreeSingle((G4Trajectory*)aTrajectory);
153}
154
155#endif
156
157
158
159
160
161
162
163
164
165
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
std::vector< G4VTrajectoryPoint * > TrajectoryPointContainer
Definition: G4Trajectory.hh:67
G4DLLIMPORT G4Allocator< G4Trajectory > aTrajectoryAllocator
Definition: G4Trajectory.cc:55
double G4double
Definition: G4Types.hh:64
#define G4DLLIMPORT
Definition: G4Types.hh:56
#define G4DLLEXPORT
Definition: G4Types.hh:55
int G4int
Definition: G4Types.hh:66
G4DLLIMPORT std::ostream G4cout
Definition: G4Step.hh:78
G4ParticleDefinition * GetParticleDefinition()
G4int GetTrackID() const
Definition: G4Trajectory.hh:92
virtual G4VTrajectoryPoint * GetPoint(G4int i) const
virtual void DrawTrajectory(G4int i_mode=0) const
G4double GetInitialKineticEnergy() const
int operator==(const G4Trajectory &right) const
Definition: G4Trajectory.hh:88
virtual void ShowTrajectory(std::ostream &os=G4cout) const
G4String GetParticleName() const
Definition: G4Trajectory.hh:96
virtual void AppendStep(const G4Step *aStep)
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)
G4double GetCharge() const
Definition: G4Trajectory.hh:98
G4int GetParentID() const
Definition: G4Trajectory.hh:94
G4int GetPDGEncoding() const
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual int GetPointEntries() const
virtual ~G4Trajectory()
Definition: G4Trajectory.cc:96
G4ThreeVector GetInitialMomentum() const