QwGeant4
QweakSimVDC_DriftCellFrontSD Class Reference

Checks for a valid VDC U-DriftCell and/or V-DriftCell hit and stores the hit information. More...

#include <QweakSimVDC_DriftCellFrontSD.hh>

Inherits G4VSensitiveDetector.

+ Collaboration diagram for QweakSimVDC_DriftCellFrontSD:

Public Member Functions

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

Static Public Member Functions

static void SetNumberOfDriftCellsPerPlane (G4int dc_npp)
 
static void StoreDCWidthOnFrame (G4double dc_w)
 
static void StoreDCFullThickness (G4double dc_ft)
 
static void StoreDCUPlaneWireAngle (G4double dc_ua)
 
static void StoreDCVPlaneWireAngle (G4double dc_va)
 

Private Attributes

QweakSimVDC_DriftCellHitsCollectionDC_hitsCollection
 
G4int DC_ID
 
QweakSimVDCpQweakSimVDCSetup
 

Static Private Attributes

static G4int DCNumberPerPlane = 401
 
static G4double DCWidthOnFrame = 11.0*mm
 
static G4double DCFullThickness = 26.0*mm
 
static G4double DCUPlaneWireAngle = 60.0*degree
 
static G4double DCVPlaneWireAngle = -60.0*degree
 

Detailed Description

Checks for a valid VDC U-DriftCell and/or V-DriftCell hit and stores the hit information.

Placeholder for a long explaination

Definition at line 61 of file QweakSimVDC_DriftCellFrontSD.hh.

Constructor & Destructor Documentation

QweakSimVDC_DriftCellFrontSD::QweakSimVDC_DriftCellFrontSD ( G4String  name)

Definition at line 55 of file QweakSimVDC_DriftCellFrontSD.cc.

References DC_ID.

56 :G4VSensitiveDetector(name)
57 {
58  //G4cout << G4endl << "###### Calling QweakSimVDC_DriftCellFrontSD::QweakSimVDC_DriftCellFrontSD() " << G4endl << G4endl;
59 
60  collectionName.insert("DriftCellFrontCollection");
61 
62  DC_ID = -1;
63 
64  //G4cout << G4endl << "###### Leaving QweakSimVDC_DriftCellFrontSD::QweakSimVDC_DriftCellFrontSD() " << G4endl << G4endl;
65 }
QweakSimVDC_DriftCellFrontSD::~QweakSimVDC_DriftCellFrontSD ( )

Definition at line 68 of file QweakSimVDC_DriftCellFrontSD.cc.

69 {
70  //G4cout << G4endl << "###### Calling QweakSimVDC_DriftCellFrontSD::~QweakSimVDC_DriftCellFrontSD() " << G4endl << G4endl;
71 
72 }

Member Function Documentation

void QweakSimVDC_DriftCellFrontSD::EndOfEvent ( G4HCofThisEvent *  HCE)

Definition at line 276 of file QweakSimVDC_DriftCellFrontSD.cc.

277 {
278  //G4cout << G4endl << "###### Calling QweakSimVDC_DriftCellFrontSD::EndOfEvent() " << G4endl << G4endl;
279 
280 // G4int NbDCHits = DC_hitsCollection->entries();
281 
282 // G4cout << "\n-------->Hits Collection: in this event they are " << NbDCHits
283 // << " hits in the Drift Cells : " << G4endl;
284 // for (G4int i=0;i<NbDCHits;i++) (*DC_hitsCollection)[i]->Print();
285 
286 
287 
288  //G4cout << G4endl << "###### Leaving QweakSimVDC_DriftCellFrontSD::EndOfEvent() " << G4endl << G4endl;
289 
290 }
void QweakSimVDC_DriftCellFrontSD::Initialize ( G4HCofThisEvent *  HCE)

Definition at line 75 of file QweakSimVDC_DriftCellFrontSD.cc.

References DC_hitsCollection, and DC_ID.

76 {
77  //G4cout << G4endl << "###### Calling QweakSimVDC_DriftCellFrontSD::Initialize() " << G4endl << G4endl;
78 
79  DC_hitsCollection = new QweakSimVDC_DriftCellHitsCollection(SensitiveDetectorName,collectionName[0]);
80 
81  // for(int i=0;i<DCNumberPerPlane;i++)
82  //{ DC_hitsCollection->insert(new QweakSimVDC_DriftCellHit); }
83 
84  if(DC_ID<0)
85  { DC_ID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]); }
86  //{ DC_ID = G4SDManager::GetSDMpointer()->GetCollectionID(DC_hitsCollection); }
87  HCE->AddHitsCollection( DC_ID, DC_hitsCollection );
88 
89  //G4cout << G4endl << "###### Leaving QweakSimVDC_DriftCellFrontSD::Initialize() " << G4endl << G4endl;
90 }
G4THitsCollection< QweakSimVDC_DriftCellHit > QweakSimVDC_DriftCellHitsCollection
QweakSimVDC_DriftCellHitsCollection * DC_hitsCollection
G4bool QweakSimVDC_DriftCellFrontSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *  ROhist 
)

Definition at line 93 of file QweakSimVDC_DriftCellFrontSD.cc.

References DC_hitsCollection, DCFullThickness, DCUPlaneWireAngle, DCVPlaneWireAngle, DCWidthOnFrame, QweakSimVDC_DriftCellHit::GetLogV(), QweakSimVDC_DriftCellHit::StoreCellPos(), QweakSimVDC_DriftCellHit::StoreCellRot(), QweakSimVDC_DriftCellHit::StoreDCFullThickness(), QweakSimVDC_DriftCellHit::StoreDCUPlaneWireAngle(), QweakSimVDC_DriftCellHit::StoreDCVPlaneWireAngle(), QweakSimVDC_DriftCellHit::StoreDCWidthOnFrame(), QweakSimVDC_DriftCellHit::StoreDriftCellPlaneID(), QweakSimVDC_DriftCellHit::StoreKineticEnergy(), QweakSimVDC_DriftCellHit::StoreLocalPos(), QweakSimVDC_DriftCellHit::StoreLogV(), QweakSimVDC_DriftCellHit::StoreMomentumDirection(), QweakSimVDC_DriftCellHit::StoreOriginVertexKineticEnergy(), QweakSimVDC_DriftCellHit::StoreOriginVertexMomentumDirection(), QweakSimVDC_DriftCellHit::StoreOriginVertexPosition(), QweakSimVDC_DriftCellHit::StoreTime(), QweakSimVDC_DriftCellHit::StoreTotalEnergy(), and QweakSimVDC_DriftCellHit::StoreWorldPos().

94 {
95  //G4cout << G4endl << "###### Calling QweakSimVDC_DriftCellFrontSD::ProcessHits() " << G4endl << G4endl;
96 
97  //G4double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge();
98 
99  // G4cout << "%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% Charge =" << charge << G4endl;
100  //
101  // dismiss photons
102  // if(charge==0.) {
103 
104  // G4cout << "==================================================================" << G4endl;
105  // G4cout << ">>>> Neutral Particle detected (e.g. Photon) <<<<<<" << G4endl;
106  // G4cout << ">>>> Up to now we don't count HITs for neutrals, therefore: <<<<<<" << G4endl;
107  // G4cout << ">>>> Aborting QweakSimVDC_DriftCellFrontSD::ProcessHits <<<<<<" << G4endl;
108  // G4cout << "==================================================================" << G4endl;
109 
110  // G4cout << G4endl << "###### Leaving QweakSimVDC_DriftCellFrontSD::ProcessHits() " << G4endl << G4endl;
111  // return false;
112  // }
113 
114 
115 
116 
117  //=====================================================================================
118  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
119  //G4StepPoint* postStepPoint = aStep->GetPreStepPoint();
120  G4TouchableHistory* theTouchable = (G4TouchableHistory*)(preStepPoint->GetTouchable());
121  G4VPhysicalVolume* physVol = theTouchable->GetVolume();
122  //====================================================================================
123 
124 
125 // G4cout << "==================================================================" << G4endl;
126 // G4cout << ">>>>>> We are in volume :" << physVol->GetName() << "<<<<<<" << G4endl;
127 // G4cout << "==================================================================" << G4endl;
128 
129 
130  if( physVol->GetName().compare("VDC_DriftCellFront_Physical")==0 )
131  {
132  // G4cout << "==================================================================" << G4endl;
133  // G4cout << ">>>>>>> Particle crossing : VDC_DriftCellFront_Physical <<<<<<<<<" << G4endl;
134  // G4cout << "==================================================================" << G4endl;
135 
136  //const G4VProcess* process = aStep->GetPostStepPoint()->GetProcessDefinedStep();
137  //G4cout << "@@@@@@@@@@@@@> Physical Processes in VDC_DriftCellFront_Physical = "<< process->GetProcessName() << G4endl;
138 
139  }
140  else
141  {
142  // G4cout << "==================================================================" << G4endl;
143  // G4cout << ">>>> ERROR: We are not within a physical DriftCell Volume <<<<<<" << G4endl;
144  // G4cout << ">>>> Aborting QweakSimVDC_DriftCellFrontSD::ProcessHits <<<<<<" << G4endl;
145  // G4cout << "==================================================================" << G4endl;
146 
147  return false;
148  }
149 
150  //================================================================================================
151  // here we are only interested if the particle is crossing a boundary
152  // otherwise we get additional hits due to physical processes within the sensitive volume
153  //
154  // If you want to identify the first step in a volume: pick fGeomBoundary status in preStepPoint
155  // If you want to identify a step out of a volume : pick fGeomBoundary status in postStepPoint
156  //=================================================================================================
157  //
158  if (preStepPoint->GetStepStatus() != fGeomBoundary)
159  {
160  // G4cout << "=======================================================================" << G4endl;
161  // G4cout << ">>>> QweakSimVDC_DriftCellFrontSD: We are NOT crossing a boundary <<<<<<" << G4endl;
162  // G4cout << ">>>> Aborting QweakSimVDC_DriftCellFrontSD::ProcessHits <<<<<<" << G4endl;
163  // G4cout << "=======================================================================" << G4endl;
164  //
165  // G4cout << G4endl << "###### Leaving QweakSimVDC_DriftCellFrontSD::ProcessHits() " << G4endl << G4endl;
166 
167  return false;
168  }
169 
170 
171  G4double currentKineticEnergy = aStep->GetTrack()->GetKineticEnergy();
172  G4double currentTotalEnergy = aStep->GetTrack()->GetTotalEnergy();
173  G4ThreeVector currentMomentumDirection = aStep->GetTrack()->GetMomentumDirection();
174 
175  //G4int trackID = aStep->GetTrack()->GetTrackID();
176  //G4cout << "====> Track ID : " << trackID << G4endl;
177 
178  //G4int parentID = aStep->GetTrack()->GetParentID();
179  //G4cout << "====> Parent ID : " << parentID << G4endl;
180 
181  // Vertex: where this track was generated
182  G4ThreeVector originVertexPosition = aStep->GetTrack()->GetVertexPosition();
183  G4ThreeVector originVertexMomentumDirection = aStep->GetTrack()->GetVertexMomentumDirection();
184  G4double originVertexKineticEnergy = aStep->GetTrack()->GetVertexKineticEnergy();
185 
186 // G4cout << "====> originVertexPosition_X [cm] : " << originVertexPosition.x()/cm << G4endl;
187 // G4cout << "====> originVertexPosition_Y [cm] : " << originVertexPosition.y()/cm << G4endl;
188 // G4cout << "====> originVertexPosition_Z [cm] : " << originVertexPosition.z()/cm << G4endl;
189 
190 // G4cout << "====> originVertexMomentumDirection_X : " << originVertexMomentumDirection.x() << G4endl;
191 // G4cout << "====> originVertexMomentumDirection_Y : " << originVertexMomentumDirection.y() << G4endl;
192 // G4cout << "====> originVertexMomentumDirection_Z : " << originVertexMomentumDirection.z() << G4endl;
193 
194 // G4cout << "====> originVertexKineticEnergy [Mev] : " << originVertexKineticEnergy/MeV << G4endl;
195 
196 
197  //G4int MotherCopyNo = theTouchable->GetVolume(1)->GetCopyNo(); // one Mother Volume
198  //G4int DriftCellCopyNo = theTouchable->GetVolume()->GetCopyNo(); // but several Driftcells per MV
199  G4int DriftCellReplicaNo = theTouchable->GetReplicaNumber(); // but several Driftcells per MV
200 
201 // G4cout << "%%%%%%%%%%%%%%%%%%%%%%% VDC_DriftCellFront MV CopyNumber :" << MotherCopyNo << G4endl;
202 // G4cout << "%%%%%%%%%%%%%%%%%%%%%%% VDC_DriftCellFront CopyNumber :" << DriftCellCopyNo << G4endl;
203 // G4cout << "%%%%%%%%%%%%%%%%%%%%%%% VDC_DriftCellFront ReplicaNumber :" << DriftCellReplicaNo << G4endl;
204 
205 
206  G4ThreeVector worldPos = preStepPoint->GetPosition();
207 // G4cout << "%%%%%%%%%%%%%%%%%%%%%%% VDC_DriftCellFront WorldPos.X [cm]:" << worldPos.x()/cm << G4endl;
208 // G4cout << "%%%%%%%%%%%%%%%%%%%%%%% VDC_DriftCellFront WorldPos.Y [cm]:" << worldPos.y()/cm << G4endl;
209 // G4cout << "%%%%%%%%%%%%%%%%%%%%%%% VDC_DriftCellFront WorldPos.Z [cm]:" << worldPos.z()/cm << G4endl;
210 
211  G4ThreeVector localPos = theTouchable->GetHistory()->GetTopTransform().TransformPoint(worldPos);
212 // G4cout << "%%%%%%%%%%%%%%%%%%%%%%% VDC_DriftCellFront LocalPos.X [cm]:" << localPos.x()/cm << G4endl;
213 // G4cout << "%%%%%%%%%%%%%%%%%%%%%%% VDC_DriftCellFront LocalPos.Y [cm]:" << localPos.y()/cm << G4endl;
214 // G4cout << "%%%%%%%%%%%%%%%%%%%%%%% VDC_DriftCellFront LocalPos.Z [cm]:" << localPos.z()/cm << G4endl;
215 
216  QweakSimVDC_DriftCellHit* aHit = new QweakSimVDC_DriftCellHit(DriftCellReplicaNo);
217  //QweakSimVDC_DriftCellHit* aHit = (*DC_hitsCollection)[DriftCellReplicaNo];
218 
219 // G4cout << "%%%%%%%%%%%%%%%%%%%%%%% got QweakSimVDC_DriftCellHit* aHit" << G4endl;
220 
221  aHit->StoreWorldPos(worldPos);
222  aHit->StoreLocalPos(localPos);
223  aHit->StoreTime(preStepPoint->GetGlobalTime());
224 
225  aHit->StoreOriginVertexPosition(originVertexPosition);
226  aHit->StoreOriginVertexKineticEnergy(originVertexKineticEnergy);
227  aHit->StoreOriginVertexMomentumDirection(originVertexMomentumDirection);
228 
229  aHit->StoreMomentumDirection(currentMomentumDirection);
230  aHit->StoreKineticEnergy(currentKineticEnergy);
231  aHit->StoreTotalEnergy(currentTotalEnergy);
232 
233 
234  G4double myDCUPlaneWireAngle = 90.0*degree - DCUPlaneWireAngle;
235  G4double myDCVPlaneWireAngle = -90.0*degree - DCVPlaneWireAngle;
236 
237 // G4cout << " DCWidthOnFrame [mm] = " << DCWidthOnFrame/mm << G4endl;
238 // G4cout << " DCFullThickness [mm] = " << DCFullThickness/mm << G4endl;
239 // G4cout << " DCUPlaneWireAngle [deg] = " << myDCUPlaneWireAngle/degree << G4endl;
240 // G4cout << " DCVPlaneWireAngle [deg] = " << myDCVPlaneWireAngle/degree << G4endl;
241 
244  aHit->StoreDCUPlaneWireAngle(myDCUPlaneWireAngle);
245  aHit->StoreDCVPlaneWireAngle(myDCVPlaneWireAngle);
246 
247 
248 
249  // check if it is first touch
250  if(!(aHit->GetLogV()))
251  {
252  // store translation and rotation matrix of the drift cell
253  // for the sake of drawing the hit
254  aHit->StoreLogV(physVol->GetLogicalVolume());
255  G4AffineTransform aTrans = theTouchable->GetHistory()->GetTopTransform();
256  aTrans.Invert();
257  aHit->StoreCellRot(aTrans.NetRotation());
258  aHit->StoreCellPos(aTrans.NetTranslation());
259  }
260 
261  if( physVol->GetName().compare("VDC_DriftCellFront_Physical") == 0 ) { aHit->StoreDriftCellPlaneID(0);}
262  if( physVol->GetName().compare("VDC_DriftCellBack_Physical") == 0 ) { aHit->StoreDriftCellPlaneID(1);}
263 
264 // G4cout << "%%%%%%%%%%%%%%%%%%%%%%% before inserting the hit" << G4endl;
265 
266  DC_hitsCollection->insert(aHit);
267 
268 // G4cout << G4endl << "###### Leaving QweakSimVDC_DriftCellFrontSD::ProcessHits() " << G4endl << G4endl;
269 
270  return true;
271 }
void StoreCellPos(G4ThreeVector xyz)
void StoreOriginVertexKineticEnergy(G4double ekin)
void StoreDCUPlaneWireAngle(G4double dca)
void StoreDriftCellPlaneID(G4int dcplane_id)
QweakSimVDC_DriftCellHitsCollection * DC_hitsCollection
Region 3 Vertical Drift Chamber Drift Cell Hit.
void StoreKineticEnergy(G4double ekin)
void StoreMomentumDirection(G4ThreeVector pxyz)
void StoreDCVPlaneWireAngle(G4double dca)
void StoreLogV(G4LogicalVolume *val)
void StoreWorldPos(G4ThreeVector xyz)
void StoreOriginVertexPosition(G4ThreeVector xyz)
void StoreOriginVertexMomentumDirection(G4ThreeVector pxyz)
void StoreLocalPos(G4ThreeVector xyz)
const G4LogicalVolume * GetLogV() const
void StoreCellRot(G4RotationMatrix rmat)

+ Here is the call graph for this function:

static void QweakSimVDC_DriftCellFrontSD::SetNumberOfDriftCellsPerPlane ( G4int  dc_npp)
inlinestatic

Definition at line 72 of file QweakSimVDC_DriftCellFrontSD.hh.

References DCNumberPerPlane.

Referenced by QweakSimVDC::SetVDC_DriftCellGeometryUpdate().

+ Here is the caller graph for this function:

static void QweakSimVDC_DriftCellFrontSD::StoreDCFullThickness ( G4double  dc_ft)
inlinestatic

Definition at line 74 of file QweakSimVDC_DriftCellFrontSD.hh.

References DCFullThickness.

Referenced by QweakSimVDC::SetVDC_DriftCellFullThickness().

74 { DCFullThickness = dc_ft; }

+ Here is the caller graph for this function:

static void QweakSimVDC_DriftCellFrontSD::StoreDCUPlaneWireAngle ( G4double  dc_ua)
inlinestatic

Definition at line 75 of file QweakSimVDC_DriftCellFrontSD.hh.

References DCUPlaneWireAngle.

Referenced by QweakSimVDC::SetVDC_DriftCellFrontWireAngle().

+ Here is the caller graph for this function:

static void QweakSimVDC_DriftCellFrontSD::StoreDCVPlaneWireAngle ( G4double  dc_va)
inlinestatic

Definition at line 76 of file QweakSimVDC_DriftCellFrontSD.hh.

References DCVPlaneWireAngle.

Referenced by QweakSimVDC::SetVDC_DriftCellBackWireAngle().

+ Here is the caller graph for this function:

static void QweakSimVDC_DriftCellFrontSD::StoreDCWidthOnFrame ( G4double  dc_w)
inlinestatic

Definition at line 73 of file QweakSimVDC_DriftCellFrontSD.hh.

References DCWidthOnFrame.

Referenced by QweakSimVDC::SetVDC_DriftCellFullWidthOnFrame().

+ Here is the caller graph for this function:

Field Documentation

QweakSimVDC_DriftCellHitsCollection* QweakSimVDC_DriftCellFrontSD::DC_hitsCollection
private

Definition at line 80 of file QweakSimVDC_DriftCellFrontSD.hh.

Referenced by Initialize(), and ProcessHits().

G4int QweakSimVDC_DriftCellFrontSD::DC_ID
private

Definition at line 82 of file QweakSimVDC_DriftCellFrontSD.hh.

Referenced by Initialize(), and QweakSimVDC_DriftCellFrontSD().

G4double QweakSimVDC_DriftCellFrontSD::DCFullThickness = 26.0*mm
staticprivate

Definition at line 86 of file QweakSimVDC_DriftCellFrontSD.hh.

Referenced by ProcessHits(), and StoreDCFullThickness().

G4int QweakSimVDC_DriftCellFrontSD::DCNumberPerPlane = 401
staticprivate

Definition at line 84 of file QweakSimVDC_DriftCellFrontSD.hh.

Referenced by SetNumberOfDriftCellsPerPlane().

G4double QweakSimVDC_DriftCellFrontSD::DCUPlaneWireAngle = 60.0*degree
staticprivate

Definition at line 87 of file QweakSimVDC_DriftCellFrontSD.hh.

Referenced by ProcessHits(), and StoreDCUPlaneWireAngle().

G4double QweakSimVDC_DriftCellFrontSD::DCVPlaneWireAngle = -60.0*degree
staticprivate

Definition at line 88 of file QweakSimVDC_DriftCellFrontSD.hh.

Referenced by ProcessHits(), and StoreDCVPlaneWireAngle().

G4double QweakSimVDC_DriftCellFrontSD::DCWidthOnFrame = 11.0*mm
staticprivate

Definition at line 85 of file QweakSimVDC_DriftCellFrontSD.hh.

Referenced by ProcessHits(), and StoreDCWidthOnFrame().

QweakSimVDC* QweakSimVDC_DriftCellFrontSD::pQweakSimVDCSetup
private

Definition at line 90 of file QweakSimVDC_DriftCellFrontSD.hh.


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