QwGeant4
QweakSimPMTOnly_DetectorSD Class Reference

#include <QweakSimPMTOnly_DetectorSD.hh>

Inherits G4VSensitiveDetector.

Public Member Functions

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

Private Attributes

QweakSimPMTOnly_DetectorHitsCollectionPMTOnly_DetectorHitsCollection
 
G4int PMTOnly_CollectionID
 

Detailed Description

Definition at line 19 of file QweakSimPMTOnly_DetectorSD.hh.

Constructor & Destructor Documentation

QweakSimPMTOnly_DetectorSD::QweakSimPMTOnly_DetectorSD ( G4String  name)

Definition at line 19 of file QweakSimPMTOnly_DetectorSD.cc.

References PMTOnly_CollectionID.

20 :G4VSensitiveDetector(name)
21 {
22  collectionName.insert("PMTOnlyCollection");
24 }
QweakSimPMTOnly_DetectorSD::~QweakSimPMTOnly_DetectorSD ( )

Definition at line 30 of file QweakSimPMTOnly_DetectorSD.cc.

31 {
32  //delete PMTOnly_DetectorHitsCollection;
33 }

Member Function Documentation

void QweakSimPMTOnly_DetectorSD::EndOfEvent ( G4HCofThisEvent *  HCE)

Definition at line 144 of file QweakSimPMTOnly_DetectorSD.cc.

145 {
146  //G4int NbDCHits = DC_hitsCollection->entries();
147 
148  //G4cout << "\n-------->Hits Collection: in this event they are " << NbDCHits
149  // << " hits in the Drift Cells : " << G4endl;
150  //for (G4int i=0;i<NbDCHits;i++) (
151  // *DC_hitsCollection)[i]->Print();
152 
153 }
void QweakSimPMTOnly_DetectorSD::Initialize ( G4HCofThisEvent *  HCE)

Definition at line 39 of file QweakSimPMTOnly_DetectorSD.cc.

References PMTOnly_CollectionID, and PMTOnly_DetectorHitsCollection.

40 {
41  PMTOnly_DetectorHitsCollection = new QweakSimPMTOnly_DetectorHitsCollection(SensitiveDetectorName,collectionName[0]);
42  //if(PMTOnly_CollectionID<0)
43  PMTOnly_CollectionID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
44  HCE->AddHitsCollection(PMTOnly_CollectionID , PMTOnly_DetectorHitsCollection);
45 
46 }
QweakSimPMTOnly_DetectorHitsCollection * PMTOnly_DetectorHitsCollection
G4THitsCollection< QweakSimPMTOnly_DetectorHit > QweakSimPMTOnly_DetectorHitsCollection
G4bool QweakSimPMTOnly_DetectorSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *  ROhist 
)

Definition at line 51 of file QweakSimPMTOnly_DetectorSD.cc.

References PMTOnly_DetectorHitsCollection.

52 {
53  //--- Not Entering Geometry
54  if (aStep->GetPreStepPoint()->GetStepStatus() != fGeomBoundary ) return false;
55 
56  G4TouchableHandle theTouchable = aStep->GetPreStepPoint()->GetTouchableHandle();
57 
58  G4int motherCopyNo = theTouchable->GetVolume(1)->GetCopyNo(); // one Mother Volume
59 
60  //G4double parentID = aStep->GetTrack()->GetParentID();
61  G4double trackID = aStep->GetTrack()->GetTrackID();
62 
63  G4ParticleDefinition* fpParticleDefinition = aStep->GetTrack()->GetDefinition();
64  G4String particleName = fpParticleDefinition->GetParticleName();
65  //G4double pdgCharge = fpParticleDefinition->GetPDGCharge();
66  G4int pdgEncoding = fpParticleDefinition->GetPDGEncoding();
67 
68  //GlobalTimeOfHit
69  G4double globalTime = aStep->GetPreStepPoint()->GetGlobalTime();
70 
71  G4ThreeVector worldPos = aStep->GetPreStepPoint()->GetPosition();
72  G4ThreeVector localPos = theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldPos);
73 
74  // Local Origin (where did the hit come from)
75  G4ThreeVector originVertexPosition = aStep->GetTrack()->GetVertexPosition();
76 
77  G4ThreeVector worldMomentum = aStep->GetPreStepPoint()->GetMomentum();
78  G4ThreeVector localMomentum = theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldMomentum);
79  G4ThreeVector currentMomentumDirection = aStep->GetTrack()->GetMomentumDirection();
80  G4ThreeVector originVertexMomentumDirection = aStep->GetTrack()->GetVertexMomentumDirection();
81 
82  G4double originVertexKineticEnergy = aStep->GetTrack()->GetVertexKineticEnergy();
83  G4double originVertexTotalEnergy = aStep->GetTrack()->GetVertexKineticEnergy();//GetVertexTotalEnergy();
84 
85  G4double currentKineticEnergy = aStep->GetTrack()->GetKineticEnergy();
86  G4double currentTotalEnergy = aStep->GetTrack()->GetTotalEnergy();
87 
88  //--- PMTOnly deposited energy
89  G4double depositedEnergy = aStep->GetTotalEnergyDeposit();
90 
91 
92  //----------------- Store Hit information
93 
95 
96  aHit -> StoreTrackID(trackID);
97 
98  aHit -> StoreParticleName(particleName);
99  aHit -> StoreParticleType(pdgEncoding);
100 
101  aHit -> StoreGlobalTime(globalTime);
102 
103  aHit -> StoreWorldPosition(worldPos);
104  aHit -> StoreLocalPosition(localPos);
105  aHit -> StoreOriginVertexPosition(originVertexPosition);
106 
107  aHit -> StoreWorldMomentum(worldMomentum);
108  aHit -> StoreLocalMomentum(localMomentum);
109  aHit -> StoreMomentumDirection(currentMomentumDirection);
110  aHit -> StoreOriginVertexMomentumDirection(originVertexMomentumDirection);
111 
112  aHit -> StoreOriginVertexKineticEnergy(originVertexKineticEnergy);
113  aHit -> StoreOriginVertexTotalEnergy(originVertexTotalEnergy); // kinetic energy
114 
115  aHit -> StoreKineticEnergy(currentKineticEnergy);
116  aHit -> StoreTotalEnergy(currentTotalEnergy);
117 
118  //--- PMTOnly deposited energy
119  aHit -> StoreDepositedEnergy(depositedEnergy);
120 
121 
122  //----------------- check if it is first touch
123 
124  if ( !(aHit -> GetLogVolume()) )
125  {
126  //--- store translation and rotation matrix of the drift cell
127  //--- for the sake of drawing the hit
128  aHit -> StoreLogVolume(theTouchable->GetVolume()->GetLogicalVolume());
129  G4AffineTransform aTrans = theTouchable->GetHistory()->GetTopTransform();
130  aTrans.Invert();
131  aHit -> StoreCellRotation(aTrans.NetRotation());
132  aHit -> StoreCellPosition(aTrans.NetTranslation());
133  }
134 
135  PMTOnly_DetectorHitsCollection->insert(aHit);
136 
137  return true;
138 }
QweakSimPMTOnly_DetectorHitsCollection * PMTOnly_DetectorHitsCollection

Field Documentation

G4int QweakSimPMTOnly_DetectorSD::PMTOnly_CollectionID
private

Definition at line 33 of file QweakSimPMTOnly_DetectorSD.hh.

Referenced by Initialize(), and QweakSimPMTOnly_DetectorSD().

QweakSimPMTOnly_DetectorHitsCollection* QweakSimPMTOnly_DetectorSD::PMTOnly_DetectorHitsCollection
private

Definition at line 32 of file QweakSimPMTOnly_DetectorSD.hh.

Referenced by Initialize(), and ProcessHits().


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