QwGeant4
QweakSimHDC_WirePlaneSD Class Reference

Checks for a valid HDC WirePlane hit and stores the hit information (Up to now I only activated only 1 out of 6 sensitive wireplanes per HDC) More...

#include <QweakSimHDC_WirePlaneSD.hh>

Inherits G4VSensitiveDetector.

Public Member Functions

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

Private Attributes

QweakSimHDC_WirePlane_HitsCollectionhitsCollection
 
G4int HCID
 

Detailed Description

Checks for a valid HDC WirePlane hit and stores the hit information (Up to now I only activated only 1 out of 6 sensitive wireplanes per HDC)

Placeholder for a long explaination

Definition at line 60 of file QweakSimHDC_WirePlaneSD.hh.

Constructor & Destructor Documentation

QweakSimHDC_WirePlaneSD::QweakSimHDC_WirePlaneSD ( G4String  name)

Definition at line 47 of file QweakSimHDC_WirePlaneSD.cc.

References HCID.

48 :G4VSensitiveDetector(name)
49 {
50 // G4cout << G4endl << "###### Calling QweakSimHDC_WirePlaneSD::QweakSimHDC_WirePlaneSD() " << G4endl << G4endl;
51 
52  collectionName.insert("HDCWirePlaneCollection"); //collectionName[0]
53  HCID = -1;
54 
55 // G4cout << G4endl << "###### Leaving QweakSimHDC_WirePlaneSD::QweakSimHDC_WirePlaneSD() " << G4endl << G4endl;
56 }
QweakSimHDC_WirePlaneSD::~QweakSimHDC_WirePlaneSD ( )

Definition at line 59 of file QweakSimHDC_WirePlaneSD.cc.

60 {
61  //G4cout << G4endl << "###### Calling QweakSimHDC_WirePlaneSD::~QweakSimHDC_WirePlaneSD() " << G4endl << G4endl;
62  //delete hitsCollection;
63  //G4cout << G4endl << "###### Leaving QweakSimHDC_WirePlaneSD::~QweakSimHDC_WirePlaneSD() " << G4endl << G4endl;
64 }

Member Function Documentation

void QweakSimHDC_WirePlaneSD::EndOfEvent ( G4HCofThisEvent *  HCE)

Definition at line 243 of file QweakSimHDC_WirePlaneSD.cc.

244 {
245  //G4cout << G4endl << "###### Calling QweakSimHDC_WirePlaneSD::EndOfEvent() " << G4endl << G4endl;
246 
247 // G4int NbHits = hitsCollection->entries();
248 // G4cout << "\n-------->Hits Collection: in this event they are " << NbHits
249 // << " hits in the Wire Planes: " << G4endl;
250 // for (G4int i=0;i<NbHits;i++) (*hitsCollection)[i]->Print();
251 
252  //G4cout << G4endl << "###### Leaving QweakSimHDC_WirePlaneSD::EndOfEvent() " << G4endl << G4endl;
253 
254 }
void QweakSimHDC_WirePlaneSD::Initialize ( G4HCofThisEvent *  HCE)

Definition at line 67 of file QweakSimHDC_WirePlaneSD.cc.

References HCID, and hitsCollection.

68 {
69  //G4cout << G4endl << "###### Calling QweakSimHDC_WirePlaneSD::Initialize() " << G4endl << G4endl;
70 
71  hitsCollection = new QweakSimHDC_WirePlane_HitsCollection(SensitiveDetectorName,collectionName[0]);
72 
73  HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
74  // HCID = G4SDManager::GetSDMpointer()->GetCollectionID(hitsCollection);
75 
76  HCE->AddHitsCollection( HCID, hitsCollection );
77 
78  //G4cout << G4endl << "###### Leaving QweakSimHDC_WirePlaneSD::Initialize() " << G4endl << G4endl;
79 }
QweakSimHDC_WirePlane_HitsCollection * hitsCollection
G4THitsCollection< QweakSimHDC_WirePlaneHit > QweakSimHDC_WirePlane_HitsCollection
G4bool QweakSimHDC_WirePlaneSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *  ROhist 
)

Definition at line 82 of file QweakSimHDC_WirePlaneSD.cc.

References hitsCollection, QweakSimHDC_WirePlaneHit::StoreGlobalTime(), QweakSimHDC_WirePlaneHit::StoreHDCID(), QweakSimHDC_WirePlaneHit::StoreKineticEnergy(), QweakSimHDC_WirePlaneHit::StoreLocalMomentum(), QweakSimHDC_WirePlaneHit::StoreLocalPosition(), QweakSimHDC_WirePlaneHit::StoreOriginVertexKineticEnergy(), QweakSimHDC_WirePlaneHit::StoreOriginVertexMomentumDirection(), QweakSimHDC_WirePlaneHit::StoreOriginVertexPosition(), QweakSimHDC_WirePlaneHit::StoreParticleName(), QweakSimHDC_WirePlaneHit::StoreParticleType(), QweakSimHDC_WirePlaneHit::StoreTotalEnergy(), QweakSimHDC_WirePlaneHit::StoreWirePlaneID(), QweakSimHDC_WirePlaneHit::StoreWorldMomentum(), and QweakSimHDC_WirePlaneHit::StoreWorldPosition().

83 {
84  //G4cout << G4endl << "###### Calling QweakSimHDC_WirePlaneSD::ProcessHits() " << G4endl << G4endl;
85 
86 
87  G4double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge();
88  G4String particlename = aStep->GetTrack()->GetDefinition()->GetParticleName();
89  G4int particletype = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
90 
91 // G4cout << "%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% Charge =" << charge << G4endl;
92 
93  // reject non-charged particle hits
94  if (fabs(charge)<0.1)
95  return false;
96 
97  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
98  G4StepPoint* postStepPoint = aStep->GetPreStepPoint();
99  G4TouchableHistory* theTouchable = (G4TouchableHistory*)(preStepPoint->GetTouchable());
100  G4VPhysicalVolume* physVol = theTouchable->GetVolume();
101 
102  // we are only interested if the particle is crossing a boundary (glass->vaccuum)
103  if (postStepPoint->GetStepStatus() != fGeomBoundary) {
104 // G4cout << "=======================================================================" << G4endl;
105 // G4cout << ">>>> QweakSimHDC_WirePlaneSD: We are NOT crossing a boundary <<<<<<" << G4endl;
106 // G4cout << ">>>> Aborting QweakSimHDC_WirePlaneSD::ProcessHits <<<<<<" << G4endl;
107 // G4cout << "=======================================================================" << G4endl;
108 
109  return false;
110  }
111 
112 // G4cout << "================== We are in volume :" << physVol->GetName() << G4endl << G4endl;
113 
114 
115  if(strcmp(physVol->GetName(),"HDC_WirePlane_Physical")==0) {
116 
117 // G4cout << "=============================================================" << G4endl;
118 // G4cout << ">>>>>>> Particle crossing : HDC_WirePlane_Physical <<<<<<<<<" << G4endl;
119 // G4cout << "=============================================================" << G4endl;
120 
121 // G4int trackID = aStep->GetTrack()->GetTrackID();
122 // G4cout << "====> Track ID : " << trackID << G4endl;
123 
124 // G4int parentID = aStep->GetTrack()->GetParentID();
125 // G4cout << "====> Parent ID : " << parentID << G4endl;
126 
127 
128 
129 // G4int MotherCopyNo = theTouchable->GetVolume(1)->GetCopyNo(); // Several MotherVolumes
130  G4int WirePlaneCopyNo = theTouchable->GetVolume()->GetCopyNo(); // but only one WirePlane per MV
131 // G4int MotherReplicaNo = theTouchable->GetReplicaNumber(1); // Several MotherVolumes
132 
133  G4int MotherCopyNo2 = theTouchable->GetVolume(2)->GetCopyNo(); // Several MotherVolumes
134 // G4int MotherReplicaNo2 = theTouchable->GetReplicaNumber(2); // Several MotherVolumes
135 
136 
137 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane MV CopyNumber :" << MotherCopyNo << G4endl;
138 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane CopyNumber :" << WirePlaneCopyNo << G4endl;
139 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane ReplicaNumber :" << MotherReplicaNo << G4endl;
140 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane MV CopyNumber2 :" << MotherCopyNo2 << G4endl;
141 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane ReplicaNumber2 :" << MotherReplicaNo2 << G4endl;
142 
143 
144  G4ThreeVector worldPos = preStepPoint->GetPosition();
145 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane WorldPos.X [cm]:" << worldPos.x()/cm << G4endl;
146 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane WorldPos.Y [cm]:" << worldPos.y()/cm << G4endl;
147 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane WorldPos.Z [cm]:" << worldPos.z()/cm << G4endl;
148 
149  G4ThreeVector localPos = theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldPos);
150 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane LocalPos.X [cm]:" << localPos.x()/cm << G4endl;
151 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane LocalPos.Y [cm]:" << localPos.y()/cm << G4endl;
152 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane LocalPos.Z [cm]:" << localPos.z()/cm << G4endl;
153 
154  G4ThreeVector worldMomentum = preStepPoint->GetMomentum();
155 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane WorldMomentum.X [GeV]:" << worldMomentum.x()/GeV << G4endl;
156 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane WorldMomentum.Y [GeV]:" << worldMomentum.y()/GeV << G4endl;
157 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane WorldMomentum.Z [GeV]:" << worldMomentum.z()/GeV << G4endl;
158 
159  G4ThreeVector localMomentum = theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldMomentum);
160 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane LocalMomentum.X [GeV]:" << localMomentum.x()/GeV << G4endl;
161 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane LocalMomentum.Y [GeV]:" << localMomentum.y()/GeV << G4endl;
162 // G4cout << "%%%%%%%%%%%%%%%%%%% WirePlane LocalMomentum.Z [GeV]:" << localMomentum.z()/GeV << G4endl;
163 
164 
165 
166  // Vertex: where this track was generated
167  G4ThreeVector originVertexPosition = aStep->GetTrack()->GetVertexPosition();
168  G4ThreeVector originVertexMomentumDirection = aStep->GetTrack()->GetVertexMomentumDirection();
169  G4double originVertexKineticEnergy = aStep->GetTrack()->GetVertexKineticEnergy();
170 
171 // G4cout << "====> originVertexPosition_X [cm] : " << originVertexPosition.x()/cm << G4endl;
172 // G4cout << "====> originVertexPosition_Y [cm] : " << originVertexPosition.y()/cm << G4endl;
173 // G4cout << "====> originVertexPosition_Z [cm] : " << originVertexPosition.z()/cm << G4endl;
174 
175 // G4cout << "====> originVertexMomentumDirection_X : " << originVertexMomentumDirection.x() << G4endl;
176 // G4cout << "====> originVertexMomentumDirection_Y : " << originVertexMomentumDirection.y() << G4endl;
177 // G4cout << "====> originVertexMomentumDirection_Z : " << originVertexMomentumDirection.z() << G4endl;
178 
179 // G4cout << "====> originVertexKineticEnergy [Mev] : " << originVertexKineticEnergy/MeV << G4endl;
180 
181 
182  G4double currentKineticEnergy = preStepPoint->GetKineticEnergy();
183  G4double currentTotalEnergy = preStepPoint->GetTotalEnergy();
184 
185  // get User Track Info
186  //QweakSimTrackInformation* info = (QweakSimTrackInformation*) (aStep->GetTrack()->GetUserInformation());
187  //info->Print();
188 
189 // G4double primaryQ2 = info->GetPrimaryQ2();
190 // G4double crossSection = info->GetCrossSection();
191 // G4double crossSectionWeight = info->GetCrossSectionWeight();
192 // G4int primaryEventNumber = info->GetPrimaryEventNumber();
193 
194 // G4cout << " @@@@@@@@@@@@@@@@@@@@@@@ Original Track Q2 : " << primaryQ2 << G4endl;
195 // G4cout << " @@@@@@@@@@@@@@@@@@@@@@@ Original Track Event Number : " << primaryEventNumber << G4endl;
196 // G4cout << " @@@@@@@@@@@@@@@@@@@@@@@ Original CrossSection Weight (ub) : " << crossSectionWeight << G4endl << G4endl;
197 
198 
199 
200  //G4cout << " =====> Storing HDC hit information into aHit" << G4endl;
201 
202  // QweakSimHDCHit* aHit = new QweakSimHDCHit(MotherCopyNo); // there is only one plane per motherVolume
203  QweakSimHDC_WirePlaneHit* aHit = new QweakSimHDC_WirePlaneHit(); // there is only one plane per motherVolume
204 
205  aHit->StoreHDCID(MotherCopyNo2); // 0: Front HDC , 1: Back HDC
206  aHit->StoreWirePlaneID(WirePlaneCopyNo); // [0,5] : WirePlane #0 to #5 (6 wire plane in total per HDC)
207 
208  aHit->StoreWorldPosition(worldPos);
209  aHit->StoreLocalPosition(localPos);
210  aHit->StoreGlobalTime(preStepPoint->GetGlobalTime());
211 
212  aHit->StoreOriginVertexPosition(originVertexPosition);
213  aHit->StoreOriginVertexMomentumDirection(originVertexMomentumDirection);
214  aHit->StoreOriginVertexKineticEnergy(originVertexKineticEnergy);
215 
216 // aHit->StorePrimaryQ2(primaryQ2);
217 // aHit->StoreCrossSection(crossSection);
218 // aHit->StoreCrossSectionWeight(crossSectionWeight);
219 // aHit->StorePrimaryEventNumber(primaryEventNumber);
220 
221  aHit->StoreWorldMomentum(worldMomentum);
222  aHit->StoreLocalMomentum(localMomentum);
223 
224  aHit->StoreParticleName(particlename);
225  aHit->StoreParticleType(particletype);
226 
227  aHit->StoreKineticEnergy(currentKineticEnergy);
228  aHit->StoreTotalEnergy(currentTotalEnergy);
229 
230  // now store all the hit information into the hit collection
231  hitsCollection->insert(aHit);
232 
233  } // end of if(strcmp(physVol->GetName(),"HDC_WirePlane_Physical")==0)
234 
235 
236  //G4cout << G4endl << "###### Leaving QweakSimHDC_WirePlaneSD::ProcessHits() " << G4endl << G4endl;
237 
238  return true;
239 }
void StoreKineticEnergy(G4double ekin)
void StoreWorldPosition(G4ThreeVector gxyz)
QweakSimHDC_WirePlane_HitsCollection * hitsCollection
void StoreOriginVertexPosition(G4ThreeVector xyz)
void StoreWorldMomentum(G4ThreeVector gpxyz)
void StoreWirePlaneID(G4int wireplane_ID)
Handles hits of the HDC wire planes.
void StoreOriginVertexMomentumDirection(G4ThreeVector pxyz)
void StoreLocalPosition(G4ThreeVector lxyz)
void StoreLocalMomentum(G4ThreeVector lpxyz)
void StoreOriginVertexKineticEnergy(G4double ekin)

+ Here is the call graph for this function:

Field Documentation

G4int QweakSimHDC_WirePlaneSD::HCID
private

Definition at line 74 of file QweakSimHDC_WirePlaneSD.hh.

Referenced by Initialize(), and QweakSimHDC_WirePlaneSD().

QweakSimHDC_WirePlane_HitsCollection* QweakSimHDC_WirePlaneSD::hitsCollection
private

Definition at line 72 of file QweakSimHDC_WirePlaneSD.hh.

Referenced by Initialize(), and ProcessHits().


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