Berkeley Nuclear Data Software
SimulationDataStructures.h
Go to the documentation of this file.
1 //Constructed by Josh Brown.
2 //brown.ja@berkeley.edu
3 #ifndef _SIMULATION_DATA_STRUCTURES_H_
4 #define _SIMULATION_DATA_STRUCTURES_H_
5 
6 
7 #ifdef OSLY_PRESENT
8 #include "OScintLightYield.h"
9 #endif
10 
11 
12 #ifndef __ROOTCLING__
13 #include "G4VHit.hh"
14 #include "G4THitsCollection.hh"
15 #include "G4Allocator.hh"
16 #include "G4ThreeVector.hh"
17 #include "G4Step.hh"
18 #else
19 #endif
20 
21 #include <string>
22 #include <vector>
23 #include <iostream>
24 #include <cmath>
25 namespace NSDG4
26 {
30 {
31 public:
32  float m_x;
33  float m_y;
34  float m_z;
35  float m_kE;
36  double m_t;
37  // float m_w;
38 };
39 
40 
41 class stepInfo
42 {
43 public:
44  float m_x;
45  float m_y;
46  float m_z;
47  float m_kE;
48  float m_ED;
49  double m_t;
50  int m_lvID;
51  double inline getDeltaL(stepInfo a_rhs);
52 };
53 
54 
56 {
57 public:
59  int m_lvID;
61  double m_amplitude;
65  double m_time;
66 };
67 // #ifndef __ROOTCLING__
68 // class TrackInfo : public G4VHit
69 
70 // #else
71 class TrackInfo
72 // #endif
73 {
74 public:
76  TrackInfo();
77 
78  //setters
79  void setTrackID(int a_trackID);
80  void setParentTrackID(int a_parentTrackID);
81  void setlvID(int a_lvID);
82  void setLvName(std::string a_lvName);
83  void setParticleName(std::string a_particleName);
84  void setCreatorProcess(std::string a_creatorProcess);
85  void setVertexInfo(const VertexInfo& a_vertexInfo);
86  void setPostSteps(const std::vector<stepInfo>& a_postSteps);
87 
88  //getters
89  int getTrackID() const;
90  int getParentTrackID() const;
91  int getlvID() const;
92  std::string getLvName() const;
93  std::string getParticleName() const;
94  std::string getCreatorProcess() const;
95  VertexInfo getVertexInfo() const;
96  std::vector<stepInfo> getPostSteps() const;
97  std::ostream& printInfo(std::ostream& a_stream=std::cout);
98 
99 
100 
103 #ifdef OSLY_PRESENT
104  double getLight(OScintLightYield& a_ly);
105 #endif
107  double getTotalED();
109  double getTime();
115  bool fullyStopped(double a_fraction=0.97);
122  std::vector<float> getDTheta();
123 
127  float getMeanED(float& a_x,
128  float& a_y,
129  float& a_z
130  );
131 
132 //stuff hidden from root used in geant
133 #ifndef __ROOTCLING__
134  void initialize(G4Step* a_step, int a_lvID=0);
135  void addStep(G4Step* a_step, int a_lvID=0);
136 #endif
137 
138 private:
140  int m_trackID;
142  int m_parentTrackID;
144  int m_lvID;
146  std::string m_lvName;
148  std::string m_particleName;
151  std::string m_creatorProcess;
154  VertexInfo m_vertexInfo;
158  std::vector<stepInfo> m_postSteps;
159 };
160 
161 
163 {
164 public:
165  std::ostream& printInfo(std::ostream& a_stream=std::cout);
166  std::vector<TrackInfo> m_tracks;
167  //for the particle generated at the beginning of the event
168  double m_primEnergy;
169  double m_primTheta;
170  double m_primPhi;
174  std::vector<ScintEventInfo> getScintEventInfo(std::string a_lvName
175  );
176  //calculations
179  bool hasED();
182  bool hasParticles(std::vector<std::string> a_partList,
183  std::string a_lvName =""
184  );
185 
186  //this function gets returns a vector of the tracks associated with a
187  //specific logical volume
188  std::vector<TrackInfo> getTracks(std::string a_lvName);
189  std::vector<TrackInfo> getTracks(std::vector<std::string> a_lvList);
190  std::vector<double> getSourcePosition();
191 
192 };
193 
194 #ifndef __ROOTCLING__
195 class TrackInfoHit : public G4VHit
196 {
197 public:
199  inline void* operator new(size_t);
200  inline void operator delete(void*);
201 };
202 
203 typedef G4THitsCollection<TrackInfoHit> TrackInfoHitCollection;
204 
205 extern G4ThreadLocal G4Allocator<TrackInfoHit>* TrackInfoHitAllocator;
206 
207 inline void* TrackInfoHit::operator new(size_t)
208 {
210  TrackInfoHitAllocator = new G4Allocator<TrackInfoHit>;
211  return (void *) TrackInfoHitAllocator->MallocSingle();
212 }
213 
214 inline void TrackInfoHit::operator delete(void *hit)
215 {
216  TrackInfoHitAllocator->FreeSingle((TrackInfoHit*) hit);
217 }
218 #endif
219 
221 {
222  double delX = m_x - a_rhs.m_x;
223  double delY = m_y - a_rhs.m_y;
224  double delZ = m_z - a_rhs.m_z;
225  return std::sqrt(delX*delX+delY*delY+delZ*delZ);
226 }
227 
228 
230 {
231 public:
238 };
239 
240 }
241 
242 #endif
Definition: SimulationDataStructures.h:163
std::vector< double > getSourcePosition()
Definition: SimulationDataStructures.cpp:452
bool hasED()
Definition: SimulationDataStructures.cpp:383
std::vector< TrackInfo > getTracks(std::string a_lvName)
Definition: SimulationDataStructures.cpp:425
double m_primTheta
Definition: SimulationDataStructures.h:169
bool hasParticles(std::vector< std::string > a_partList, std::string a_lvName="")
Definition: SimulationDataStructures.cpp:399
double m_primEnergy
Definition: SimulationDataStructures.h:168
std::ostream & printInfo(std::ostream &a_stream=std::cout)
Definition: SimulationDataStructures.cpp:367
std::vector< TrackInfo > m_tracks
Definition: SimulationDataStructures.h:166
std::vector< ScintEventInfo > getScintEventInfo(std::string a_lvName)
Definition: SimulationDataStructures.cpp:469
double m_primPhi
Definition: SimulationDataStructures.h:170
Definition: SimulationDataStructures.h:230
double m_initialEnergy
start energy
Definition: SimulationDataStructures.h:233
int m_numberOfBounces
number of bounces
Definition: SimulationDataStructures.h:235
bool m_isDetected
detected
Definition: SimulationDataStructures.h:237
Definition: SimulationDataStructures.h:56
double m_amplitude
light or energy depending on build conditions
Definition: SimulationDataStructures.h:61
double m_initPartEnergy
the inital event particles energy
Definition: SimulationDataStructures.h:63
int m_lvID
logical vulume copy number
Definition: SimulationDataStructures.h:59
double m_time
the time the energy was deposited
Definition: SimulationDataStructures.h:65
Definition: SimulationDataStructures.h:196
Definition: SimulationDataStructures.h:73
TrackInfo()
default constructor
Definition: SimulationDataStructures.cpp:108
void setParentTrackID(int a_parentTrackID)
Definition: SimulationDataStructures.cpp:126
std::string getParticleName() const
Definition: SimulationDataStructures.cpp:170
void initialize(G4Step *a_step, int a_lvID=0)
Definition: SimulationDataStructures.cpp:21
std::vector< stepInfo > getPostSteps() const
Definition: SimulationDataStructures.cpp:182
void setVertexInfo(const VertexInfo &a_vertexInfo)
Definition: SimulationDataStructures.cpp:146
void setPostSteps(const std::vector< stepInfo > &a_postSteps)
Definition: SimulationDataStructures.cpp:150
void setlvID(int a_lvID)
Definition: SimulationDataStructures.cpp:130
std::ostream & printInfo(std::ostream &a_stream=std::cout)
Definition: SimulationDataStructures.cpp:352
int getlvID() const
Definition: SimulationDataStructures.cpp:162
void addStep(G4Step *a_step, int a_lvID=0)
Definition: SimulationDataStructures.cpp:77
void setLvName(std::string a_lvName)
Definition: SimulationDataStructures.cpp:134
int getParentTrackID() const
Definition: SimulationDataStructures.cpp:158
void setParticleName(std::string a_particleName)
Definition: SimulationDataStructures.cpp:138
void setCreatorProcess(std::string a_creatorProcess)
Definition: SimulationDataStructures.cpp:142
std::string getLvName() const
Definition: SimulationDataStructures.cpp:166
std::vector< float > getDTheta()
Definition: SimulationDataStructures.cpp:266
int getTrackID() const
Definition: SimulationDataStructures.cpp:154
bool fullyStopped(double a_fraction=0.97)
Definition: SimulationDataStructures.cpp:250
std::string getCreatorProcess() const
Definition: SimulationDataStructures.cpp:174
void setTrackID(int a_trackID)
Definition: SimulationDataStructures.cpp:122
double getTotalED()
returns the total energy deposited by the track
Definition: SimulationDataStructures.cpp:231
VertexInfo getVertexInfo() const
Definition: SimulationDataStructures.cpp:178
double getTime()
returns the time of the start of the track
Definition: SimulationDataStructures.cpp:242
float getMeanED(float &a_x, float &a_y, float &a_z)
Definition: SimulationDataStructures.cpp:316
Definition: SimulationDataStructures.h:30
float m_x
Definition: SimulationDataStructures.h:32
float m_kE
Definition: SimulationDataStructures.h:35
float m_z
Definition: SimulationDataStructures.h:34
float m_y
Definition: SimulationDataStructures.h:33
double m_t
Definition: SimulationDataStructures.h:36
Definition: SimulationDataStructures.h:42
int m_lvID
Definition: SimulationDataStructures.h:50
double getDeltaL(stepInfo a_rhs)
Definition: SimulationDataStructures.h:220
float m_z
Definition: SimulationDataStructures.h:46
float m_y
Definition: SimulationDataStructures.h:45
double m_t
Definition: SimulationDataStructures.h:49
float m_x
Definition: SimulationDataStructures.h:44
float m_kE
Definition: SimulationDataStructures.h:47
float m_ED
Definition: SimulationDataStructures.h:48
Definition: OScintLightYield.h:18
Definition: AbsLYAna.h:7
G4THitsCollection< TrackInfoHit > TrackInfoHitCollection
Definition: SimulationDataStructures.h:203
G4ThreadLocal G4Allocator< TrackInfoHit > * TrackInfoHitAllocator