Geant4-11
G4VTrajectory.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// G4VTrajectory
27//
28// Class description:
29//
30// This class is the abstract base class representing a trajectory of
31// a particle being tracked.
32// Its concrete class includes information of:
33// 1) List of trajectory points composing the trajectory;
34// 2) Static information of the particle which generated the 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: Katsuya.Amako@kek.jp)
40// Makoto Asai (e-mail: asai@slac.stanford.edu)
41// Takashi Sasaki (e-mail: Takashi.Sasaki@kek.jp)
42// --------------------------------------------------------------------
43#ifndef G4VTrajectory_hh
44#define G4VTrajectory_hh 1
45
46#include <vector>
47#include <map>
48
49#include "globals.hh"
50#include "G4ThreeVector.hh"
51
52class G4Step;
54class G4AttDef;
55class G4AttValue;
56
58{
59 public:
60
62 virtual ~G4VTrajectory();
63 // Constructor/Destrcutor
64
65 G4bool operator == (const G4VTrajectory& right) const;
66 // Equality operator
67
68 virtual G4int GetTrackID() const = 0;
69 virtual G4int GetParentID() const = 0;
70 virtual G4String GetParticleName() const = 0;
71 virtual G4double GetCharge() const = 0;
72 // Accessors
73
74 virtual G4int GetPDGEncoding() const = 0;
75 // Charge is that of G4DynamicParticle
76 virtual G4ThreeVector GetInitialMomentum() const = 0;
77 // Zero will be returned if the particle does not have PDG code.
78 // Momentum at the origin of the track in global coordinate system
79
80 virtual G4int GetPointEntries() const = 0;
81 // Returns the number of trajectory points
82 virtual G4VTrajectoryPoint* GetPoint(G4int i) const = 0;
83 // Returns i-th trajectory point
84 virtual void ShowTrajectory(std::ostream& os=G4cout) const;
85 // Converts attributes in trajectory (and trajectory point if
86 // needed) to ostream. A default implementation in this base class
87 // may be used or may be overridden in the concrete class. Note:
88 // the user needs to follow with new-line or end-of-string,
89 // depending on the nature of os
90 virtual void DrawTrajectory() const;
91 // Draw the trajectory. A default implementation in this base
92 // class may be used or may be overridden in the concrete class
93 virtual const std::map<G4String,G4AttDef>* GetAttDefs() const
94 { return nullptr; }
95 // If implemented by a derived class, returns a pointer to a map of
96 // attribute definitions for the attribute values below. The user
97 // must test the validity of this pointer. See G4Trajectory for an
98 // example of a concrete implementation of this method
99 virtual std::vector<G4AttValue>* CreateAttValues() const
100 { return nullptr; }
101 // If implemented by a derived class, returns a pointer to a list
102 // of attribute values suitable, e.g., for picking. Each must
103 // refer to an attribute definition in the above map; its name is
104 // the key. The user must test the validity of this pointer (it
105 // must be non-zero and conform to the G4AttDefs, which may be
106 // checked with G4AttCheck) and delete the list after use. See
107 // G4Trajectory for an example of a concrete implementation of this
108 // method and G4VTrajectory::ShowTrajectory for an example of its use.
109
110 virtual void AppendStep(const G4Step* aStep) = 0;
111 virtual void MergeTrajectory(G4VTrajectory* secondTrajectory) = 0;
112 // Methods invoked exclusively by G4TrackingManager
113};
114
115#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
Definition: G4Step.hh:62
virtual G4VTrajectoryPoint * GetPoint(G4int i) const =0
virtual G4int GetPointEntries() const =0
virtual void ShowTrajectory(std::ostream &os=G4cout) const
virtual G4String GetParticleName() const =0
virtual G4int GetTrackID() const =0
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual G4ThreeVector GetInitialMomentum() const =0
virtual void AppendStep(const G4Step *aStep)=0
G4bool operator==(const G4VTrajectory &right) const
virtual G4double GetCharge() const =0
virtual G4int GetParentID() const =0
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)=0
virtual void DrawTrajectory() const
virtual G4int GetPDGEncoding() const =0
virtual ~G4VTrajectory()
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const