LBNETrajectory.hh
Go to the documentation of this file.
1 //
2 // LBNETrajectory.hh
3 //
4 
5 #ifndef LBNETrajectory_h
6 #define LBNETrajectory_h 1
7 
8 #include "G4VTrajectory.hh"
9 #include "G4Allocator.hh"
10 #include <stdlib.h>
11 #include "G4ThreeVector.hh"
12 #include "G4ios.hh"
13 #include <vector>
14 #include <iostream>
15 #include "globals.hh"
16 #include "G4ParticleDefinition.hh"
17 #include "G4TrajectoryPoint.hh"
18 #include "G4Track.hh"
19 #include "G4Step.hh"
20 #include "G4UserSteppingAction.hh"
21 #include "G4Colour.hh"
22 
23 class LBNETrajectory;
24 class G4Polyline;
25 
26 typedef std::vector<G4VTrajectoryPoint*> LBNETrajectoryPointContainer;
27 typedef std::vector<G4ThreeVector> LBNETrajectoryMomentumContainer;
28 typedef std::vector<G4String> LBNETrajectoryVolumeName;
29 typedef std::vector<G4String> LBNETrajectoryMaterialName;
30 typedef std::vector<G4double> DVec;
31 
32 class LBNETrajectory : public G4VTrajectory
33 {
34  public:
36  LBNETrajectory(const G4Track* aTrack);
39  virtual ~LBNETrajectory();
40 
41  inline void* operator new(size_t);
42  inline void operator delete(void*);
43  inline int operator == (const LBNETrajectory& right) const
44  {return (this==&right);}
45 
46  inline G4int GetTrackID() const
47  { return fTrackID; }
48  inline G4int GetParentID() const
49  { return fParentID; }
50  inline G4String GetParticleName() const
51  { return fParticleName; }
52  inline G4double GetCharge() const
53  { return fPDGCharge; }
54  inline G4double GetMass() const
55  { return fParticleMass; }
56  inline G4int GetPDGEncoding() const
57  { return fPDGEncoding; }
58  inline const G4ThreeVector& GetVertexPosition() const
59  { return fVertexPosition; }
60  virtual int GetPointEntries() const
61  { return fPositionRecord->size(); }
62  virtual G4VTrajectoryPoint* GetPoint(G4int i) const
63  { return (*fPositionRecord)[i]; }
64  virtual G4ThreeVector GetMomentum(G4int i) const
65  { return (*fMomentumRecord)[i]; }
66  virtual G4String GetPreStepVolumeName(G4int i) const
67  { return (*fPreStepVolume)[i]; }
68  inline G4ThreeVector GetInitialMomentum() const
69  { return fMomentum; }
70  virtual G4int GetTgen() const
71  { return fTgen;}
72  inline G4int GetDecayCode() const
73  { return fDecayCode;}
74  virtual G4double GetNImpWt() const
75  { return fNImpWt;}
76  virtual G4double GetStepLength(G4int i) const
77  {return (*fStepLength)[i];}
78 
79  inline G4int GetMaterialNumber1rst() const
80  {return fMaterialNumber1rst;}
81  inline void SetMaterialNumber1rst(G4double matNum)
82  {fMaterialNumber1rst = matNum;}
83 
84  inline G4int GetMaterialNumberLast() const
85  {return fMaterialNumberLast;}
86  inline void SetMaterialNumberLast(G4double matNum)
87  {fMaterialNumberLast = matNum;}
88 
89  inline G4String GetVolName1rst() const
90  {return fVolName1rst;}
91 
92  inline G4double GetTimeStart() const
93  {return fTimeStart;}
94 
95  inline G4ThreeVector GetPolarization() const
96  { return fPolarization;}
97 
98  inline G4int GetPDGNucleus() const {
99  return fPDGNucleus; }
100 
101  inline G4String GetProcessName() const
102  { return fProcessName; }
103 
104  inline G4String GetMaterialName1rst() const
105  { return fMaterialName1rst; }
106  virtual G4String GetMaterialName(G4int i)const
107  {return (*fMaterialName)[i];}
108  inline const G4ThreeVector GetParentMomentumAtThisProduction() const {
110  }
111  virtual void ShowTrajectory() const;
112  virtual void ShowTrajectory(std::ostream& o) const;
113  virtual void DrawTrajectory() const;
114 
115  virtual const std::map<G4String,G4AttDef>* GetAttDefs() const;
116  virtual std::vector<G4AttValue>* CreateAttValues() const;
117 
118  virtual void AppendStep(const G4Step* aStep);
119  virtual void MergeTrajectory(G4VTrajectory* secondTrajectory);
120 
121  G4ParticleDefinition* GetParticleDefinition();
122 
123  private:
126  G4int fDecayCode;
127  G4int fEventNo;
128  G4int fTgen;
129  G4double fNImpWt;
130  G4int fTrackID;
131  G4int fParentID;
132  G4ParticleDefinition* fParticleDefinition;
133  G4String fParticleName;
134  G4double fPDGCharge;
136  G4ThreeVector fMomentum;
137  G4ThreeVector fVertexPosition;
138  G4double fParticleMass;
141 // April 2014 : added to support Dk2Nu Ntuple output
142  G4int fMaterialNumber1rst; // at the first point
143  G4String fVolName1rst; // at the first point
144  G4int fMaterialNumberLast; // at the last recorded point
145  G4double fTimeStart;
146  G4ThreeVector fPolarization;
147  G4int fPDGNucleus; // The nucleus that caused the scatter..
148  // if I can find it..
149  G4String fProcessName; // Process name at creation time.
150  G4ThreeVector fParentMomentumAtThisProduction; //for ppfx...A.Bashyal
151  LBNETrajectoryMaterialName* fMaterialName; //for ppf....A.Bashyal
152  G4String fMaterialName1rst; // First material name where the corresponding track
153  // has been created.
154 
155 };
156 extern G4Allocator<LBNETrajectory> myTrajectoryAllocator;
157 
158 inline void* LBNETrajectory::operator new(size_t)
159 {
160  void* aTrajectory;
161  aTrajectory = (void*)myTrajectoryAllocator.MallocSingle();
162  return aTrajectory;
163 }
164 
165 inline void LBNETrajectory::operator delete(void* aTrajectory)
166 {
167  myTrajectoryAllocator.FreeSingle((LBNETrajectory*)aTrajectory);
168 }
169 #endif
170 
171 
172 
G4String GetVolName1rst() const
G4double GetTimeStart() const
std::vector< G4double > DVec
G4ThreeVector GetPolarization() const
constexpr auto const & right(const_AssnsIter< L, R, D, Dir > const &a, const_AssnsIter< L, R, D, Dir > const &b)
Definition: AssnsIter.h:104
virtual std::vector< G4AttValue > * CreateAttValues() const
virtual void MergeTrajectory(G4VTrajectory *secondTrajectory)
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
G4ParticleDefinition * GetParticleDefinition()
virtual G4String GetMaterialName(G4int i) const
G4double GetMass() const
LBNETrajectoryVolumeName * fPreStepVolume
virtual ~LBNETrajectory()
G4int GetDecayCode() const
G4int GetTrackID() const
std::vector< G4VTrajectoryPoint * > LBNETrajectoryPointContainer
G4int GetPDGNucleus() const
G4ParticleDefinition * fParticleDefinition
virtual void ShowTrajectory() const
void SetMaterialNumber1rst(G4double matNum)
virtual G4double GetNImpWt() const
virtual void DrawTrajectory() const
virtual G4String GetPreStepVolumeName(G4int i) const
const G4ThreeVector GetParentMomentumAtThisProduction() const
G4int GetMaterialNumberLast() const
LBNETrajectoryMaterialName * fMaterialName
virtual G4VTrajectoryPoint * GetPoint(G4int i) const
LBNETrajectoryMomentumContainer * fMomentumRecord
virtual G4int GetTgen() const
G4ThreeVector GetInitialMomentum() const
G4String fMaterialName1rst
G4String fVolName1rst
G4double fParticleMass
void SetMaterialNumberLast(G4double matNum)
G4String GetProcessName() const
G4int GetMaterialNumber1rst() const
G4ThreeVector fVertexPosition
G4ThreeVector fPolarization
std::vector< G4String > LBNETrajectoryVolumeName
std::vector< G4String > LBNETrajectoryMaterialName
G4int GetPDGEncoding() const
G4int GetParentID() const
G4double GetCharge() const
virtual void AppendStep(const G4Step *aStep)
G4String GetMaterialName1rst() const
G4String GetParticleName() const
LBNETrajectoryPointContainer * fPositionRecord
G4ThreeVector fMomentum
const G4ThreeVector & GetVertexPosition() const
G4ThreeVector fParentMomentumAtThisProduction
virtual int GetPointEntries() const
G4String fParticleName
virtual G4ThreeVector GetMomentum(G4int i) const
virtual G4double GetStepLength(G4int i) const
int operator==(const LBNETrajectory &right) const
std::vector< G4ThreeVector > LBNETrajectoryMomentumContainer
G4String fProcessName
G4Allocator< LBNETrajectory > myTrajectoryAllocator