QwGeant4
QweakSimLumi_DetectorSD Class Reference

#include <QweakSimLumi_DetectorSD.hh>

Inherits G4VSensitiveDetector.

Public Member Functions

 QweakSimLumi_DetectorSD (G4String name)
 
 ~QweakSimLumi_DetectorSD ()
 
void Initialize (G4HCofThisEvent *HCE)
 
G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROhist)
 
void EndOfEvent (G4HCofThisEvent *HCE)
 

Private Attributes

QweakSimLumi_DetectorHitsCollectionLumi_DetectorHitsCollection
 
G4int Lumi_CollectionID
 

Detailed Description

Definition at line 7 of file QweakSimLumi_DetectorSD.hh.

Constructor & Destructor Documentation

QweakSimLumi_DetectorSD::QweakSimLumi_DetectorSD ( G4String  name)

Definition at line 6 of file QweakSimLumi_DetectorSD.cc.

References Lumi_CollectionID.

6  :G4VSensitiveDetector(name) {
7  collectionName.insert("LumiCollection");
9 }
QweakSimLumi_DetectorSD::~QweakSimLumi_DetectorSD ( )

Definition at line 11 of file QweakSimLumi_DetectorSD.cc.

11  {
12  //delete Lumi_HitCollection;
13 }

Member Function Documentation

void QweakSimLumi_DetectorSD::EndOfEvent ( G4HCofThisEvent *  HCE)

Definition at line 115 of file QweakSimLumi_DetectorSD.cc.

115  {
116  // placeholder
117 }
void QweakSimLumi_DetectorSD::Initialize ( G4HCofThisEvent *  HCE)

Definition at line 15 of file QweakSimLumi_DetectorSD.cc.

References Lumi_CollectionID, and Lumi_DetectorHitsCollection.

15  {
16  Lumi_DetectorHitsCollection = new QweakSimLumi_DetectorHitsCollection(SensitiveDetectorName, collectionName[0]);
17  Lumi_CollectionID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
18  HCE->AddHitsCollection(Lumi_CollectionID, Lumi_DetectorHitsCollection);
19 }
QweakSimLumi_DetectorHitsCollection * Lumi_DetectorHitsCollection
G4THitsCollection< QweakSimLumi_DetectorHit > QweakSimLumi_DetectorHitsCollection
G4bool QweakSimLumi_DetectorSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *  ROhist 
)

Definition at line 21 of file QweakSimLumi_DetectorSD.cc.

References Lumi_DetectorHitsCollection.

21  {
22  // Kill chargeless particles that don't deposit much energy
23  G4double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge();
24  if(charge==0. && aStep->GetTrack()->GetTotalEnergy()/MeV < 0.1) {
25  return false;
26  }
27 
28  // Kill optical photons (don't need to save these tracks
29  if( aStep->GetTrack()->GetDefinition()->GetParticleName() == "opticalphoton") {
30  return false;
31  }
32  // begin experimental
33  if(aStep->GetPreStepPoint()->GetStepStatus() != fGeomBoundary)
34  return false;
35  // end experimental
36  G4TouchableHandle theTouchable = aStep->GetPreStepPoint()->GetTouchableHandle();
37  G4int motherCopyNo = theTouchable->GetVolume(1)->GetCopyNo();
38 
39  G4double trackID = aStep->GetTrack()->GetTrackID();
40 
41  G4ParticleDefinition* fpParticleDefinition = aStep->GetTrack()->GetDefinition();
42  G4String particleName = fpParticleDefinition->GetParticleName();
43 
44  G4int pdgEncoding = fpParticleDefinition->GetPDGEncoding();
45 
46  G4double globalTime = aStep->GetPreStepPoint()->GetGlobalTime();
47 
48  G4ThreeVector worldPos = aStep->GetPreStepPoint()->GetPosition();
49  G4ThreeVector localPos = theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldPos);
50  G4ThreeVector originVertexPosition = aStep->GetTrack()->GetVertexPosition();
51 
52  G4ThreeVector worldMomentum = aStep->GetPreStepPoint()->GetMomentum();
53  G4ThreeVector localMomentum = theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldMomentum);
54  G4ThreeVector currentMomentumDirection = aStep->GetTrack()->GetMomentumDirection();
55  G4ThreeVector originVertexMomentumDirection = aStep->GetTrack()->GetVertexMomentumDirection();
56 
57  G4double originVertexKineticEnergy = aStep->GetTrack()->GetVertexKineticEnergy();
58  G4double originVertexTotalEnergy = aStep->GetTrack()->GetVertexKineticEnergy();//GetVertexTotalEnergy();
59 
60  G4double currentKineticEnergy = aStep->GetTrack()->GetKineticEnergy();
61  G4double currentTotalEnergy = aStep->GetTrack()->GetTotalEnergy();
62 
63  G4double depositedEnergy = aStep->GetTotalEnergyDeposit();
64 
65  QweakSimLumi_DetectorHit* aHit = new QweakSimLumi_DetectorHit(motherCopyNo);
66 
67  aHit -> StoreTrackID(trackID);
68 
69  aHit -> StoreParticleName(particleName);
70  aHit -> StoreParticleType(pdgEncoding);
71 
72  aHit -> StoreGlobalTime(globalTime);
73 
74  //aHit -> StoreHasBeenHit(hasBeenHit);
75  //aHit -> StoreEdgeEventFlag(edgeEventFlag);
76  //aHit -> SotreNbOfHits(nbOfHits);
77 
78  aHit -> StoreWorldPosition(worldPos);
79  aHit -> StoreLocalPosition(localPos);
80  //aHit -> StoreLocalExitPosition(localExitPos);
81  aHit -> StoreOriginVertexPosition(originVertexPosition);
82 
83  aHit -> StoreWorldMomentum(worldMomentum);
84  aHit -> StoreLocalMomentum(localMomentum);
85  aHit -> StoreMomentumDirection(currentMomentumDirection);
86  aHit -> StoreOriginVertexMomentumDirection(originVertexMomentumDirection);
87 
88  //aHit -> StoreOriginVertexThetaAngle; //calculated in EventAction
89  //aHit -> StoreOriginVertexPhiAngle;
90 
91  aHit -> StoreOriginVertexKineticEnergy(originVertexKineticEnergy);
92  aHit -> StoreOriginVertexTotalEnergy(originVertexTotalEnergy); // kinetic energy
93 
94  aHit -> StoreKineticEnergy(currentKineticEnergy);
95  aHit -> StoreTotalEnergy(currentTotalEnergy);
96 
97  aHit -> StoreHitDepositedEnergy(depositedEnergy);
98 
99  if ( !(aHit -> GetLogVolume()) )
100  {
101  //--- store translation and rotation matrix of the drift cell
102  //--- for the sake of drawing the hit
103  aHit -> StoreLogVolume(theTouchable->GetVolume()->GetLogicalVolume());
104  G4AffineTransform aTrans = theTouchable->GetHistory()->GetTopTransform();
105  aTrans.Invert();
106  aHit -> StoreCellRotation(aTrans.NetRotation());
107  aHit -> StoreCellPosition(aTrans.NetTranslation());
108  }
109 
110  Lumi_DetectorHitsCollection->insert(aHit);
111 
112  return true;
113 }
QweakSimLumi_DetectorHitsCollection * Lumi_DetectorHitsCollection

Field Documentation

G4int QweakSimLumi_DetectorSD::Lumi_CollectionID
private

Definition at line 18 of file QweakSimLumi_DetectorSD.hh.

Referenced by Initialize(), and QweakSimLumi_DetectorSD().

QweakSimLumi_DetectorHitsCollection* QweakSimLumi_DetectorSD::Lumi_DetectorHitsCollection
private

Definition at line 17 of file QweakSimLumi_DetectorSD.hh.

Referenced by Initialize(), and ProcessHits().


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