QwGeant4
QweakSimTrackInformation Class Reference

Class with additional track information like Q2. More...

#include <QweakSimTrackInformation.hh>

Inherits G4VUserTrackInformation.

Public Member Functions

 QweakSimTrackInformation ()
 
 QweakSimTrackInformation (const G4Track *aTrack)
 
 QweakSimTrackInformation (const QweakSimTrackInformation *aTrackInfo)
 
virtual ~QweakSimTrackInformation ()
 
void * operator new (size_t)
 
void operator delete (void *aTrackInfo)
 
int operator== (const QweakSimTrackInformation &right) const
 
void Initialize ()
 
void Print () const
 
void PrintPrimaryTrackInfo () const
 
void PrintSourceTrackInfo () const
 
void PrintCerenkovImpactTrackInfo () const
 
void StoreParticleDefinition (G4ParticleDefinition *pdef)
 
G4ParticleDefinition * GetParticleDefinitionAtIndex (G4int ind)
 
G4int GetParticleHistoryLength ()
 
void StoreParentEnergy (G4double eng)
 
G4double GetParentEnergyAtIndex (G4int ind)
 
void StoreCerenkovHitEnergy (G4int ind, G4double eng)
 
G4double GetCerenkovHitEnergyAtIndex (G4int ind)
 
void StoreCreatorProcess (G4String proc)
 
G4String GetCreatorProcessAtIndex (G4int ind)
 
void StoreOriginVertex (G4ThreeVector vert)
 
G4ThreeVector GetOriginVertex (G4int ind)
 
void StorePrimaryTrackID (G4int trackid)
 
G4int GetPrimaryTrackID () const
 
void StorePrimaryParticle (G4ParticleDefinition *pdef)
 
G4ParticleDefinition * GetPrimaryParticle () const
 
void StorePrimaryPosition (G4ThreeVector xyz)
 
G4ThreeVector GetPrimaryPosition () const
 
void StorePrimaryMomentum (G4ThreeVector pxyz)
 
G4ThreeVector GetPrimaryMomentum () const
 
void StorePrimaryEnergy (G4double ekin)
 
G4double GetPrimaryEnergy () const
 
void StorePrimaryKineticEnergy (G4double ekin)
 
G4double GetPrimaryKineticEnergy () const
 
void StorePrimaryTime (G4double otime)
 
G4double GetPrimaryTime () const
 
void SetTrackingStatus (G4int i)
 
G4int GetTrackingStatus () const
 
void SetTrackIsPrimaryStatus (G4int i)
 
G4int GetTrackIsPrimaryStatus () const
 
void SetSourceTrackInformation (const G4Track *aTrack)
 
G4int GetSourceTrackID () const
 
G4ParticleDefinition * GetSourceParticle () const
 
G4ThreeVector GetSourcePosition () const
 
G4ThreeVector GetSourceMomentum () const
 
G4double GetSourceEnergy () const
 
G4double GetSourceTime () const
 
G4ThreeVector GetSourceOriginPosition () const
 
G4double GetSourceOriginEnergy () const
 
void SetImpactTrackInformationForCerenkov (const G4Track *aTrack)
 
G4int GetImpactTrackID () const
 
G4ParticleDefinition * GetImpactParticle () const
 
G4ThreeVector GetImpactPosition () const
 
G4ThreeVector GetImpactMomentum () const
 
G4double GetImpactEnergy () const
 
G4double GetImpactTime () const
 
G4ThreeVector GetImpactOriginPosition () const
 
G4double GetImpactOriginEnergy () const
 
void AddTrackInfoToCerenkovTrackHistory (const G4Track *aTrack)
 

Private Attributes

G4int primaryTrackID
 
G4ParticleDefinition * particleDefinition
 
G4ThreeVector primaryPosition
 
G4ThreeVector primaryMomentum
 
G4double primaryEnergy
 
G4double primaryKineticEnergy
 
G4double primaryTime
 
std::vector
< QweakSimTrackHistory * > 
theCerenkovTrackHistory
 
std::vector< G4ThreeVector > OriginVertex
 
std::vector< G4double > CerenkovHitEnergy
 
std::vector
< G4ParticleDefinition * > 
ParticleHistory
 
std::vector< G4String > ParticleCreatorProcess
 
std::vector< G4double > ParentEnergy
 
G4int trackingStatus
 
G4int trackIsPrimaryStatus
 
G4int sourceTrackID
 
G4ParticleDefinition * sourceDefinition
 
G4ThreeVector sourcePosition
 
G4ThreeVector sourceMomentum
 
G4double sourceEnergy
 
G4double sourceTime
 
G4ThreeVector sourceOriginPosition
 
G4double sourceOriginEnergy
 
G4int cerenkovImpactTrackID
 
G4ParticleDefinition * cerenkovImpactDefinition
 
G4ThreeVector cerenkovImpactPosition
 
G4ThreeVector cerenkovImpactMomentum
 
G4double cerenkovImpactEnergy
 
G4double cerenkovImpactTime
 
G4ThreeVector cerenkovImpactOriginPosition
 
G4double cerenkovImpactOriginEnergy
 

Detailed Description

Class with additional track information like Q2.

Placeholder for a long explaination

Definition at line 53 of file QweakSimTrackInformation.hh.

Constructor & Destructor Documentation

QweakSimTrackInformation::QweakSimTrackInformation ( )

Definition at line 36 of file QweakSimTrackInformation.cc.

References Initialize().

37 {
38  //G4cout << G4endl << "###### Calling QweakSimTrackInformation::QweakSimTrackInformation() " << G4endl << G4endl;
39 
40  Initialize();
41 
42  //G4cout << G4endl << "###### Leaving QweakSimTrackInformation::QweakSimTrackInformation() " << G4endl << G4endl;
43 }

+ Here is the call graph for this function:

QweakSimTrackInformation::QweakSimTrackInformation ( const G4Track *  aTrack)

Definition at line 51 of file QweakSimTrackInformation.cc.

References CerenkovHitEnergy, cerenkovImpactDefinition, cerenkovImpactEnergy, cerenkovImpactMomentum, cerenkovImpactPosition, cerenkovImpactTime, cerenkovImpactTrackID, OriginVertex, ParentEnergy, ParticleCreatorProcess, particleDefinition, ParticleHistory, primaryEnergy, primaryKineticEnergy, primaryMomentum, primaryPosition, primaryTime, primaryTrackID, sourceDefinition, sourceEnergy, sourceMomentum, sourcePosition, sourceTime, sourceTrackID, trackingStatus, and trackIsPrimaryStatus.

52 {
53 
54  primaryTrackID = aTrack->GetTrackID();
55  particleDefinition = aTrack->GetDefinition();
56 
57  G4double tmpE = -1.0*MeV;
58  CerenkovHitEnergy.push_back(tmpE);
59  ParticleHistory.push_back(aTrack->GetDefinition());
60  ParentEnergy.push_back(aTrack->GetTotalEnergy());
61  ParticleCreatorProcess.push_back("Primary");
62  OriginVertex.push_back(aTrack->GetPosition());
63 
64  primaryPosition = aTrack->GetPosition();
65  primaryMomentum = aTrack->GetMomentum();
66  primaryEnergy = aTrack->GetTotalEnergy();
67  primaryKineticEnergy = aTrack->GetKineticEnergy();
68  primaryTime = aTrack->GetGlobalTime();
69 
70  trackingStatus = 0;
72 
73  sourceTrackID = -1;
74  sourceDefinition = 0;
75  sourcePosition = G4ThreeVector(0.,0.,0.);
76  sourceMomentum = G4ThreeVector(0.,0.,0.);
77  sourceEnergy = 0.;
78  sourceTime = 0.;
79 
82  cerenkovImpactPosition = G4ThreeVector(0.,0.,0.);
83  cerenkovImpactMomentum = G4ThreeVector(0.,0.,0.);
85  cerenkovImpactTime = 0.;
86 }
std::vector< G4ThreeVector > OriginVertex
G4ParticleDefinition * cerenkovImpactDefinition
G4ParticleDefinition * particleDefinition
std::vector< G4ParticleDefinition * > ParticleHistory
G4ParticleDefinition * sourceDefinition
std::vector< G4String > ParticleCreatorProcess
std::vector< G4double > CerenkovHitEnergy
std::vector< G4double > ParentEnergy
QweakSimTrackInformation::QweakSimTrackInformation ( const QweakSimTrackInformation aTrackInfo)

Definition at line 90 of file QweakSimTrackInformation.cc.

References CerenkovHitEnergy, cerenkovImpactDefinition, cerenkovImpactEnergy, cerenkovImpactMomentum, cerenkovImpactOriginEnergy, cerenkovImpactOriginPosition, cerenkovImpactPosition, cerenkovImpactTime, cerenkovImpactTrackID, OriginVertex, ParentEnergy, ParticleCreatorProcess, particleDefinition, ParticleHistory, primaryEnergy, primaryMomentum, primaryPosition, primaryTime, primaryTrackID, sourceDefinition, sourceEnergy, sourceMomentum, sourceOriginEnergy, sourceOriginPosition, sourcePosition, sourceTime, sourceTrackID, trackingStatus, and trackIsPrimaryStatus.

91 {
92  //G4cout << G4endl << "###### Calling QweakSimTrackInformation::QweakSimTrackInformation(const QweakSimTrackInformation* aTrackInfo) " << G4endl << G4endl;
93 
94 
95  // Return TrackInfo which can not derived from G4Track
96  // Must be filled/provided manually
97 // primaryQ2 = aTrackInfo->primaryQ2;
98 // crossSection = aTrackInfo->crossSection;
99 // crossSectionWeight = aTrackInfo->crossSectionWeight;
100 // primaryEventNumber = aTrackInfo->primaryEventNumber;
101  //----------
102  primaryTrackID = aTrackInfo->primaryTrackID;
104  primaryPosition = aTrackInfo->primaryPosition;
105  primaryMomentum = aTrackInfo->primaryMomentum;
106  primaryEnergy = aTrackInfo->primaryEnergy;
107  primaryTime = aTrackInfo->primaryTime;
108  //---------
109  sourceTrackID = aTrackInfo->sourceTrackID;
110  sourceDefinition = aTrackInfo->sourceDefinition;
111  sourcePosition = aTrackInfo->sourcePosition;
112  sourceMomentum = aTrackInfo->sourceMomentum;
113  sourceEnergy = aTrackInfo->sourceEnergy;
114  sourceTime = aTrackInfo->sourceTime;
117  //---------
126  //---------
127  trackingStatus = aTrackInfo->trackingStatus;
129  //---------
130 
131  ParticleHistory.resize(aTrackInfo->ParticleHistory.size());
132  ParentEnergy.resize(aTrackInfo->ParticleHistory.size());
133  CerenkovHitEnergy.resize(aTrackInfo->ParticleHistory.size());
134  ParticleCreatorProcess.resize(aTrackInfo->ParticleHistory.size());
135  OriginVertex.resize(aTrackInfo->ParticleHistory.size());
136  for(size_t i = 0; i < aTrackInfo->ParticleHistory.size(); i++)
137  {
138  ParticleHistory[i] = aTrackInfo->ParticleHistory[i];
139  ParentEnergy[i] = aTrackInfo->ParentEnergy[i];
140  CerenkovHitEnergy[i] = aTrackInfo->CerenkovHitEnergy[i];
142  OriginVertex[i] = aTrackInfo->OriginVertex[i];
143  }
144 
145 }
std::vector< G4ThreeVector > OriginVertex
G4ParticleDefinition * cerenkovImpactDefinition
G4ParticleDefinition * particleDefinition
std::vector< G4ParticleDefinition * > ParticleHistory
G4ParticleDefinition * sourceDefinition
std::vector< G4String > ParticleCreatorProcess
std::vector< G4double > CerenkovHitEnergy
std::vector< G4double > ParentEnergy
QweakSimTrackInformation::~QweakSimTrackInformation ( )
virtual

Definition at line 46 of file QweakSimTrackInformation.cc.

47 {
48  //G4cout << G4endl << "###### Calling/Leaving QweakSimTrackInformation::~QweakSimTrackInformation() " << G4endl << G4endl;
49 }

Member Function Documentation

void QweakSimTrackInformation::AddTrackInfoToCerenkovTrackHistory ( const G4Track *  aTrack)

Definition at line 175 of file QweakSimTrackInformation.cc.

References QweakSimTrackHistory::AddTrackInfo(), and theCerenkovTrackHistory.

176 {
177 
178  QweakSimTrackHistory* theTrackHistory = new QweakSimTrackHistory();
179 
180  theTrackHistory->AddTrackInfo(aTrack);
181 
182  theCerenkovTrackHistory.push_back(theTrackHistory);
183 }
void AddTrackInfo(const G4Track *aTrack)
std::vector< QweakSimTrackHistory * > theCerenkovTrackHistory

+ Here is the call graph for this function:

G4double QweakSimTrackInformation::GetCerenkovHitEnergyAtIndex ( G4int  ind)

Definition at line 304 of file QweakSimTrackInformation.cc.

References CerenkovHitEnergy.

305 {
306  if (CerenkovHitEnergy.size() == 0) return -1.0;
307  if (ind < 0 || ind >= (G4int) CerenkovHitEnergy.size()) return CerenkovHitEnergy[0];
308  return CerenkovHitEnergy[ind];
309 }
std::vector< G4double > CerenkovHitEnergy
G4String QweakSimTrackInformation::GetCreatorProcessAtIndex ( G4int  ind)

Definition at line 325 of file QweakSimTrackInformation.cc.

References ParticleCreatorProcess.

326 {
327  if (ParticleCreatorProcess.size() == 0) return "None";
328  if (ind > (G4int) ParticleCreatorProcess.size() - 1) return "None";
329  return ParticleCreatorProcess[ind];
330 }
std::vector< G4String > ParticleCreatorProcess
G4double QweakSimTrackInformation::GetImpactEnergy ( ) const
inline

Definition at line 215 of file QweakSimTrackInformation.hh.

References cerenkovImpactEnergy.

G4ThreeVector QweakSimTrackInformation::GetImpactMomentum ( ) const
inline

Definition at line 214 of file QweakSimTrackInformation.hh.

References cerenkovImpactMomentum.

G4double QweakSimTrackInformation::GetImpactOriginEnergy ( ) const
inline
G4ThreeVector QweakSimTrackInformation::GetImpactOriginPosition ( ) const
inline

Definition at line 217 of file QweakSimTrackInformation.hh.

References cerenkovImpactOriginPosition.

G4ParticleDefinition* QweakSimTrackInformation::GetImpactParticle ( ) const
inline

Definition at line 212 of file QweakSimTrackInformation.hh.

References cerenkovImpactDefinition.

212 {return cerenkovImpactDefinition;}
G4ParticleDefinition * cerenkovImpactDefinition
G4ThreeVector QweakSimTrackInformation::GetImpactPosition ( ) const
inline

Definition at line 213 of file QweakSimTrackInformation.hh.

References cerenkovImpactPosition.

G4double QweakSimTrackInformation::GetImpactTime ( ) const
inline

Definition at line 216 of file QweakSimTrackInformation.hh.

References cerenkovImpactTime.

G4int QweakSimTrackInformation::GetImpactTrackID ( ) const
inline

Definition at line 211 of file QweakSimTrackInformation.hh.

References cerenkovImpactTrackID.

G4ThreeVector QweakSimTrackInformation::GetOriginVertex ( G4int  ind)

Definition at line 332 of file QweakSimTrackInformation.cc.

References OriginVertex.

333 {
334  G4ThreeVector vec(-10000*cm,-10000*cm,-10000*cm);
335 
336  if (OriginVertex.size() == 0) return vec;
337  if (ind > (G4int) OriginVertex.size() - 1) return vec;
338  return OriginVertex[ind];
339 }
std::vector< G4ThreeVector > OriginVertex
G4double QweakSimTrackInformation::GetParentEnergyAtIndex ( G4int  ind)

Definition at line 318 of file QweakSimTrackInformation.cc.

References ParentEnergy.

319 {
320  if (ParentEnergy.size() == 0) return -1.0;
321  if (ind > (G4int) ParentEnergy.size() - 1) return ParentEnergy[0];
322  return ParentEnergy[ind];
323 }
std::vector< G4double > ParentEnergy
G4ParticleDefinition * QweakSimTrackInformation::GetParticleDefinitionAtIndex ( G4int  ind)

Definition at line 311 of file QweakSimTrackInformation.cc.

References ParticleHistory.

312 {
313  if (ParticleHistory.size() == 0) return NULL;
314  if (ind > (G4int) ParticleHistory.size() - 1) return ParticleHistory[0];
315  return ParticleHistory[ind];
316 }
std::vector< G4ParticleDefinition * > ParticleHistory
G4int QweakSimTrackInformation::GetParticleHistoryLength ( )
inline

Definition at line 139 of file QweakSimTrackInformation.hh.

References ParticleHistory.

Referenced by StoreCerenkovHitEnergy().

139 {return ParticleHistory.size();};
std::vector< G4ParticleDefinition * > ParticleHistory

+ Here is the caller graph for this function:

G4double QweakSimTrackInformation::GetPrimaryEnergy ( ) const
inline

Definition at line 168 of file QweakSimTrackInformation.hh.

References primaryEnergy.

G4double QweakSimTrackInformation::GetPrimaryKineticEnergy ( ) const
inline

Definition at line 171 of file QweakSimTrackInformation.hh.

References primaryKineticEnergy.

G4ThreeVector QweakSimTrackInformation::GetPrimaryMomentum ( ) const
inline

Definition at line 165 of file QweakSimTrackInformation.hh.

References primaryMomentum.

165 {return primaryMomentum;}
G4ParticleDefinition* QweakSimTrackInformation::GetPrimaryParticle ( ) const
inline

Definition at line 159 of file QweakSimTrackInformation.hh.

References particleDefinition.

159 {return particleDefinition;}
G4ParticleDefinition * particleDefinition
G4ThreeVector QweakSimTrackInformation::GetPrimaryPosition ( ) const
inline

Definition at line 162 of file QweakSimTrackInformation.hh.

References primaryPosition.

162 {return primaryPosition;}
G4double QweakSimTrackInformation::GetPrimaryTime ( ) const
inline

Definition at line 174 of file QweakSimTrackInformation.hh.

References primaryTime.

G4int QweakSimTrackInformation::GetPrimaryTrackID ( ) const
inline

Definition at line 156 of file QweakSimTrackInformation.hh.

References primaryTrackID.

G4double QweakSimTrackInformation::GetSourceEnergy ( ) const
inline

Definition at line 203 of file QweakSimTrackInformation.hh.

References sourceEnergy.

G4ThreeVector QweakSimTrackInformation::GetSourceMomentum ( ) const
inline

Definition at line 202 of file QweakSimTrackInformation.hh.

References sourceMomentum.

202 {return sourceMomentum;}
G4double QweakSimTrackInformation::GetSourceOriginEnergy ( ) const
inline

Definition at line 206 of file QweakSimTrackInformation.hh.

References sourceOriginEnergy.

G4ThreeVector QweakSimTrackInformation::GetSourceOriginPosition ( ) const
inline

Definition at line 205 of file QweakSimTrackInformation.hh.

References sourceOriginPosition.

G4ParticleDefinition* QweakSimTrackInformation::GetSourceParticle ( ) const
inline

Definition at line 200 of file QweakSimTrackInformation.hh.

References sourceDefinition.

200 {return sourceDefinition;}
G4ParticleDefinition * sourceDefinition
G4ThreeVector QweakSimTrackInformation::GetSourcePosition ( ) const
inline

Definition at line 201 of file QweakSimTrackInformation.hh.

References sourcePosition.

201 {return sourcePosition;}
G4double QweakSimTrackInformation::GetSourceTime ( ) const
inline

Definition at line 204 of file QweakSimTrackInformation.hh.

References sourceTime.

G4int QweakSimTrackInformation::GetSourceTrackID ( ) const
inline

Definition at line 199 of file QweakSimTrackInformation.hh.

References sourceTrackID.

G4int QweakSimTrackInformation::GetTrackingStatus ( ) const
inline

Definition at line 191 of file QweakSimTrackInformation.hh.

References trackingStatus.

G4int QweakSimTrackInformation::GetTrackIsPrimaryStatus ( ) const
inline

Definition at line 194 of file QweakSimTrackInformation.hh.

References trackIsPrimaryStatus.

void QweakSimTrackInformation::Initialize ( )

Definition at line 240 of file QweakSimTrackInformation.cc.

References CerenkovHitEnergy, cerenkovImpactDefinition, cerenkovImpactEnergy, cerenkovImpactMomentum, cerenkovImpactOriginEnergy, cerenkovImpactOriginPosition, cerenkovImpactPosition, cerenkovImpactTime, cerenkovImpactTrackID, OriginVertex, ParentEnergy, ParticleCreatorProcess, particleDefinition, ParticleHistory, primaryEnergy, primaryMomentum, primaryPosition, primaryTime, primaryTrackID, sourceDefinition, sourceEnergy, sourceMomentum, sourceOriginEnergy, sourceOriginPosition, sourcePosition, sourceTime, sourceTrackID, theCerenkovTrackHistory, trackingStatus, and trackIsPrimaryStatus.

Referenced by QweakSimTrackInformation().

241 {
242  //G4cout << G4endl << "###### Calling QweakSimTrackInformation::Initialize() " << G4endl << G4endl;
243 
244  primaryTrackID = 0;
245  particleDefinition = 0;
246  primaryPosition = G4ThreeVector(0.0,0.0,0.0);
247  primaryMomentum = G4ThreeVector(0.0,0.0,0.0);
248  primaryEnergy = 0.0;
249  primaryTime = 0.0;
250 
251  ParticleHistory.resize(0);
252  ParentEnergy.resize(0);
253  CerenkovHitEnergy.resize(0);
254  ParticleCreatorProcess.resize(0);
255  OriginVertex.resize(0);
256 
257 // primaryQ2 = 0.0;
258 // crossSection = 0.0;
259 // crossSectionWeight = 0.0;
260 // primaryEventNumber = 0;
261 
262  trackingStatus = 0;
264 
265  sourceTrackID = -1;
266  sourceTrackID = -1;
267  sourceDefinition = 0;
268  sourcePosition = G4ThreeVector(0.,0.,0.);
269  sourceMomentum = G4ThreeVector(0.,0.,0.);
270  sourceEnergy = 0.;
271  sourceTime = 0.;
272  sourceOriginPosition = G4ThreeVector(0.,0.,0.);
273  sourceOriginEnergy = 0.;
274 
278  cerenkovImpactPosition = G4ThreeVector(0.,0.,0.);
279  cerenkovImpactMomentum = G4ThreeVector(0.,0.,0.);
281  cerenkovImpactTime = 0.;
282  cerenkovImpactOriginPosition = G4ThreeVector(0.,0.,0.);
284 
285 
286 
287  theCerenkovTrackHistory.clear();
288  theCerenkovTrackHistory.resize(0);
289 
290  //G4cout << G4endl << "###### Leaving QweakSimTrackInformation::Initialize() " << G4endl << G4endl;
291 
292 }
std::vector< G4ThreeVector > OriginVertex
G4ParticleDefinition * cerenkovImpactDefinition
G4ParticleDefinition * particleDefinition
std::vector< G4ParticleDefinition * > ParticleHistory
G4ParticleDefinition * sourceDefinition
std::vector< G4String > ParticleCreatorProcess
std::vector< G4double > CerenkovHitEnergy
std::vector< QweakSimTrackHistory * > theCerenkovTrackHistory
std::vector< G4double > ParentEnergy

+ Here is the caller graph for this function:

void QweakSimTrackInformation::operator delete ( void *  aTrackInfo)
inline

Definition at line 240 of file QweakSimTrackInformation.hh.

References aTrackInformationAllocator.

241 { aTrackInformationAllocator.FreeSingle((QweakSimTrackInformation*)aTrackInfo);}
Class with additional track information like Q2.
G4Allocator< QweakSimTrackInformation > aTrackInformationAllocator
void * QweakSimTrackInformation::operator new ( size_t  )
inline

Definition at line 231 of file QweakSimTrackInformation.hh.

References aTrackInformationAllocator.

232 {
233  void* aTrackInfo;
234  aTrackInfo = (void*) aTrackInformationAllocator.MallocSingle();
235  return aTrackInfo;
236 }
G4Allocator< QweakSimTrackInformation > aTrackInformationAllocator
int QweakSimTrackInformation::operator== ( const QweakSimTrackInformation right) const
inline

Definition at line 65 of file QweakSimTrackInformation.hh.

66  {return (this==&right);}
void QweakSimTrackInformation::Print ( ) const
inline

Definition at line 130 of file QweakSimTrackInformation.hh.

130 {;}
void QweakSimTrackInformation::PrintCerenkovImpactTrackInfo ( ) const

Definition at line 223 of file QweakSimTrackInformation.cc.

References cerenkovImpactDefinition, cerenkovImpactEnergy, cerenkovImpactOriginEnergy, cerenkovImpactOriginPosition, cerenkovImpactPosition, and cerenkovImpactTrackID.

224 {
225 // G4cout << "###### Calling QweakSimTrackInformation::PrintCerenkovImpactTrackInfo() " << G4endl;
226 
227  G4cout << "########################################################################" << G4endl;
228  G4cout << "Cerenkov Impact Track ID = " << cerenkovImpactTrackID << G4endl;
229  G4cout << "Cerenkov Impact Particle Name = " << cerenkovImpactDefinition->GetParticleName() << G4endl;
230  G4cout << "Cerenkov Impact Particle Energy [GeV] = " << cerenkovImpactEnergy/GeV << G4endl;
231  G4cout << "Cerenkov Impact Particle Position [mm] = " << cerenkovImpactPosition << G4endl;
232  G4cout << "Cerenkov Impact Origin Position [mm] = " << cerenkovImpactOriginPosition << G4endl;
233  G4cout << "Cerenkov Impact Origin Energy [GeV] = " << cerenkovImpactOriginEnergy << G4endl;
234  G4cout << "########################################################################" << G4endl;
235 
236 // G4cout << "###### Leaving QweakSimTrackInformation::PrintCerenkovImpactTrackInfo() " << G4endl;
237 }
G4ParticleDefinition * cerenkovImpactDefinition
void QweakSimTrackInformation::PrintPrimaryTrackInfo ( ) const

Definition at line 186 of file QweakSimTrackInformation.cc.

References particleDefinition, primaryEnergy, primaryPosition, and primaryTrackID.

187 {
188 // G4cout << "###### Calling QweakSimTrackInformation::PrintPrimaryTrackInfo() " << G4endl;
189 
190  G4cout << "########################################################################" << G4endl;
191  G4cout << "Primary Track ID = " << primaryTrackID << G4endl;
192  G4cout << "Primary Particle Name = " << particleDefinition->GetParticleName() << G4endl;
193  G4cout << "Primary Particle Energy [GeV] = " << primaryEnergy/GeV << G4endl;
194  G4cout << "Primary Start Position [mm] = " << primaryPosition << G4endl;
195  G4cout << "------------------------------------------------------------------------" << G4endl;
196 // G4cout << "Primary track with Q2 = " << primaryQ2 << G4endl;
197 // G4cout << "Primary track with CS = " << crossSection << G4endl;
198 // G4cout << "Primary track with CSW = " << crossSectionWeight << G4endl;
199 // G4cout << "Primary track event number = " << primaryEventNumber << G4endl;
200  G4cout << "########################################################################" << G4endl;
201 
202 // G4cout << "###### Leaving QweakSimTrackInformation::PrintPrimaryTrackInfo() " << G4endl;
203 }
G4ParticleDefinition * particleDefinition
void QweakSimTrackInformation::PrintSourceTrackInfo ( ) const

Definition at line 206 of file QweakSimTrackInformation.cc.

References sourceDefinition, sourceEnergy, sourceOriginEnergy, sourceOriginPosition, sourcePosition, and sourceTrackID.

207 {
208 // G4cout << "###### Calling QweakSimTrackInformation::PrintSourceTrackInfo() " << G4endl;
209 
210  G4cout << "########################################################################" << G4endl;
211  G4cout << "Source Track ID = " << sourceTrackID << G4endl;
212  G4cout << "Source Particle Name = " << sourceDefinition->GetParticleName() << G4endl;
213  G4cout << "Source Particle Energy [GeV] = " << sourceEnergy/GeV << G4endl;
214  G4cout << "Source Particle Position [mm] = " << sourcePosition << G4endl;
215  G4cout << "Source Origin Position [mm] = " << sourceOriginPosition << G4endl;
216  G4cout << "Source Origin Energy [GeV] = " << sourceOriginEnergy << G4endl;
217  G4cout << "########################################################################" << G4endl;
218 
219 // G4cout << "###### Leaving QweakSimTrackInformation::PrintSourceTrackInfo() " << G4endl;
220 }
G4ParticleDefinition * sourceDefinition
void QweakSimTrackInformation::SetImpactTrackInformationForCerenkov ( const G4Track *  aTrack)

Definition at line 162 of file QweakSimTrackInformation.cc.

References cerenkovImpactDefinition, cerenkovImpactEnergy, cerenkovImpactMomentum, cerenkovImpactOriginEnergy, cerenkovImpactOriginPosition, cerenkovImpactPosition, cerenkovImpactTime, and cerenkovImpactTrackID.

163 {
164  cerenkovImpactTrackID = aTrack->GetTrackID();
165  cerenkovImpactDefinition = aTrack->GetDefinition();
166  cerenkovImpactPosition = aTrack->GetPosition();
167  cerenkovImpactMomentum = aTrack->GetMomentum();
168  cerenkovImpactEnergy = aTrack->GetTotalEnergy();
169  cerenkovImpactTime = aTrack->GetGlobalTime();
170  cerenkovImpactOriginPosition = aTrack->GetVertexPosition();
171  cerenkovImpactOriginEnergy = aTrack->GetVertexKineticEnergy();
172 }
G4ParticleDefinition * cerenkovImpactDefinition
void QweakSimTrackInformation::SetSourceTrackInformation ( const G4Track *  aTrack)

Definition at line 148 of file QweakSimTrackInformation.cc.

References sourceDefinition, sourceEnergy, sourceMomentum, sourceOriginEnergy, sourceOriginPosition, sourcePosition, sourceTime, and sourceTrackID.

Referenced by QweakSimTrackingAction::PreUserTrackingAction().

149 {
150  sourceTrackID = aTrack->GetTrackID();
151  sourceDefinition = aTrack->GetDefinition();
152  sourcePosition = aTrack->GetPosition();
153  sourceMomentum = aTrack->GetMomentum();
154  sourceEnergy = aTrack->GetTotalEnergy();
155  sourceTime = aTrack->GetGlobalTime();
156  sourceOriginPosition = aTrack->GetVertexPosition();
157  sourceOriginEnergy = aTrack->GetVertexKineticEnergy();
158 
159 }
G4ParticleDefinition * sourceDefinition

+ Here is the caller graph for this function:

void QweakSimTrackInformation::SetTrackingStatus ( G4int  i)
inline

Definition at line 190 of file QweakSimTrackInformation.hh.

References trackingStatus.

void QweakSimTrackInformation::SetTrackIsPrimaryStatus ( G4int  i)
inline
void QweakSimTrackInformation::StoreCerenkovHitEnergy ( G4int  ind,
G4double  eng 
)

Definition at line 294 of file QweakSimTrackInformation.cc.

References CerenkovHitEnergy, and GetParticleHistoryLength().

295 {
296  if(ind < 0 || ind >= GetParticleHistoryLength()) {
297  CerenkovHitEnergy.push_back(eng);
298  return;
299  }
300 
301  CerenkovHitEnergy[ind] = eng;
302 }
std::vector< G4double > CerenkovHitEnergy

+ Here is the call graph for this function:

void QweakSimTrackInformation::StoreCreatorProcess ( G4String  proc)
inline

Definition at line 147 of file QweakSimTrackInformation.hh.

References ParticleCreatorProcess.

147 {ParticleCreatorProcess.push_back(proc);};
std::vector< G4String > ParticleCreatorProcess
void QweakSimTrackInformation::StoreOriginVertex ( G4ThreeVector  vert)
inline

Definition at line 150 of file QweakSimTrackInformation.hh.

References OriginVertex.

150 {OriginVertex.push_back(vert);};
std::vector< G4ThreeVector > OriginVertex
void QweakSimTrackInformation::StoreParentEnergy ( G4double  eng)
inline

Definition at line 141 of file QweakSimTrackInformation.hh.

References ParentEnergy.

141 {ParentEnergy.push_back(eng);};
std::vector< G4double > ParentEnergy
void QweakSimTrackInformation::StoreParticleDefinition ( G4ParticleDefinition *  pdef)
inline

Definition at line 137 of file QweakSimTrackInformation.hh.

References ParticleHistory.

137 {ParticleHistory.push_back(pdef);};
std::vector< G4ParticleDefinition * > ParticleHistory
void QweakSimTrackInformation::StorePrimaryEnergy ( G4double  ekin)
inline

Definition at line 167 of file QweakSimTrackInformation.hh.

References primaryEnergy.

void QweakSimTrackInformation::StorePrimaryKineticEnergy ( G4double  ekin)
inline

Definition at line 170 of file QweakSimTrackInformation.hh.

References primaryKineticEnergy.

void QweakSimTrackInformation::StorePrimaryMomentum ( G4ThreeVector  pxyz)
inline

Definition at line 164 of file QweakSimTrackInformation.hh.

References primaryMomentum.

164 {primaryMomentum = pxyz;}
void QweakSimTrackInformation::StorePrimaryParticle ( G4ParticleDefinition *  pdef)
inline

Definition at line 158 of file QweakSimTrackInformation.hh.

References particleDefinition.

158 {particleDefinition = pdef;}
G4ParticleDefinition * particleDefinition
void QweakSimTrackInformation::StorePrimaryPosition ( G4ThreeVector  xyz)
inline

Definition at line 161 of file QweakSimTrackInformation.hh.

References primaryPosition.

161 {primaryPosition = xyz;}
void QweakSimTrackInformation::StorePrimaryTime ( G4double  otime)
inline

Definition at line 173 of file QweakSimTrackInformation.hh.

References primaryTime.

void QweakSimTrackInformation::StorePrimaryTrackID ( G4int  trackid)
inline

Definition at line 155 of file QweakSimTrackInformation.hh.

References primaryTrackID.

Field Documentation

std::vector<G4double> QweakSimTrackInformation::CerenkovHitEnergy
private
G4ParticleDefinition* QweakSimTrackInformation::cerenkovImpactDefinition
private
G4double QweakSimTrackInformation::cerenkovImpactEnergy
private
G4ThreeVector QweakSimTrackInformation::cerenkovImpactMomentum
private
G4double QweakSimTrackInformation::cerenkovImpactOriginEnergy
private
G4ThreeVector QweakSimTrackInformation::cerenkovImpactOriginPosition
private
G4ThreeVector QweakSimTrackInformation::cerenkovImpactPosition
private
G4double QweakSimTrackInformation::cerenkovImpactTime
private
G4int QweakSimTrackInformation::cerenkovImpactTrackID
private
std::vector<G4ThreeVector> QweakSimTrackInformation::OriginVertex
private
std::vector<G4double> QweakSimTrackInformation::ParentEnergy
private
std::vector<G4String> QweakSimTrackInformation::ParticleCreatorProcess
private
G4ParticleDefinition* QweakSimTrackInformation::particleDefinition
private
std::vector<G4ParticleDefinition*> QweakSimTrackInformation::ParticleHistory
private
G4double QweakSimTrackInformation::primaryEnergy
private
G4double QweakSimTrackInformation::primaryKineticEnergy
private
G4ThreeVector QweakSimTrackInformation::primaryMomentum
private
G4ThreeVector QweakSimTrackInformation::primaryPosition
private
G4double QweakSimTrackInformation::primaryTime
private
G4int QweakSimTrackInformation::primaryTrackID
private
G4ParticleDefinition* QweakSimTrackInformation::sourceDefinition
private
G4double QweakSimTrackInformation::sourceEnergy
private
G4ThreeVector QweakSimTrackInformation::sourceMomentum
private
G4double QweakSimTrackInformation::sourceOriginEnergy
private
G4ThreeVector QweakSimTrackInformation::sourceOriginPosition
private
G4ThreeVector QweakSimTrackInformation::sourcePosition
private
G4double QweakSimTrackInformation::sourceTime
private
G4int QweakSimTrackInformation::sourceTrackID
private
std::vector<QweakSimTrackHistory*> QweakSimTrackInformation::theCerenkovTrackHistory
private

Definition at line 79 of file QweakSimTrackInformation.hh.

Referenced by AddTrackInfoToCerenkovTrackHistory(), and Initialize().

G4int QweakSimTrackInformation::trackingStatus
private
G4int QweakSimTrackInformation::trackIsPrimaryStatus
private

The documentation for this class was generated from the following files: