QwGeant4
QweakSimVDC_DriftCellBackSD Class Reference

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

#include <QweakSimVDC_DriftCellBackSD.hh>

Inherits G4VSensitiveDetector.

+ Collaboration diagram for QweakSimVDC_DriftCellBackSD:

Public Member Functions

 QweakSimVDC_DriftCellBackSD (G4String name)
 
 ~QweakSimVDC_DriftCellBackSD ()
 
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 63 of file QweakSimVDC_DriftCellBackSD.hh.

Constructor & Destructor Documentation

QweakSimVDC_DriftCellBackSD::QweakSimVDC_DriftCellBackSD ( G4String  name)

Definition at line 54 of file QweakSimVDC_DriftCellBackSD.cc.

References DC_ID.

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

Definition at line 67 of file QweakSimVDC_DriftCellBackSD.cc.

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

Member Function Documentation

void QweakSimVDC_DriftCellBackSD::EndOfEvent ( G4HCofThisEvent *  HCE)

Definition at line 275 of file QweakSimVDC_DriftCellBackSD.cc.

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

Definition at line 74 of file QweakSimVDC_DriftCellBackSD.cc.

References DC_hitsCollection, and DC_ID.

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

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

+ Here is the call graph for this function:

static void QweakSimVDC_DriftCellBackSD::SetNumberOfDriftCellsPerPlane ( G4int  dc_npp)
inlinestatic

Definition at line 74 of file QweakSimVDC_DriftCellBackSD.hh.

References DCNumberPerPlane.

Referenced by QweakSimVDC::SetVDC_DriftCellGeometryUpdate().

+ Here is the caller graph for this function:

static void QweakSimVDC_DriftCellBackSD::StoreDCFullThickness ( G4double  dc_ft)
inlinestatic

Definition at line 76 of file QweakSimVDC_DriftCellBackSD.hh.

References DCFullThickness.

Referenced by QweakSimVDC::SetVDC_DriftCellFullThickness().

76 { DCFullThickness = dc_ft; }

+ Here is the caller graph for this function:

static void QweakSimVDC_DriftCellBackSD::StoreDCUPlaneWireAngle ( G4double  dc_ua)
inlinestatic

Definition at line 77 of file QweakSimVDC_DriftCellBackSD.hh.

References DCUPlaneWireAngle.

Referenced by QweakSimVDC::SetVDC_DriftCellFrontWireAngle().

77 { DCUPlaneWireAngle = dc_ua; }

+ Here is the caller graph for this function:

static void QweakSimVDC_DriftCellBackSD::StoreDCVPlaneWireAngle ( G4double  dc_va)
inlinestatic

Definition at line 78 of file QweakSimVDC_DriftCellBackSD.hh.

References DCVPlaneWireAngle.

Referenced by QweakSimVDC::SetVDC_DriftCellBackWireAngle().

78 { DCVPlaneWireAngle = dc_va; }

+ Here is the caller graph for this function:

static void QweakSimVDC_DriftCellBackSD::StoreDCWidthOnFrame ( G4double  dc_w)
inlinestatic

Definition at line 75 of file QweakSimVDC_DriftCellBackSD.hh.

References DCWidthOnFrame.

Referenced by QweakSimVDC::SetVDC_DriftCellFullWidthOnFrame().

75 { DCWidthOnFrame = dc_w; }

+ Here is the caller graph for this function:

Field Documentation

QweakSimVDC_DriftCellHitsCollection* QweakSimVDC_DriftCellBackSD::DC_hitsCollection
private

Definition at line 82 of file QweakSimVDC_DriftCellBackSD.hh.

Referenced by Initialize(), and ProcessHits().

G4int QweakSimVDC_DriftCellBackSD::DC_ID
private

Definition at line 84 of file QweakSimVDC_DriftCellBackSD.hh.

Referenced by Initialize(), and QweakSimVDC_DriftCellBackSD().

G4double QweakSimVDC_DriftCellBackSD::DCFullThickness = 26.0*mm
staticprivate

Definition at line 88 of file QweakSimVDC_DriftCellBackSD.hh.

Referenced by ProcessHits(), and StoreDCFullThickness().

G4int QweakSimVDC_DriftCellBackSD::DCNumberPerPlane = 401
staticprivate

Definition at line 86 of file QweakSimVDC_DriftCellBackSD.hh.

Referenced by SetNumberOfDriftCellsPerPlane().

G4double QweakSimVDC_DriftCellBackSD::DCUPlaneWireAngle = 60.0*degree
staticprivate

Definition at line 89 of file QweakSimVDC_DriftCellBackSD.hh.

Referenced by ProcessHits(), and StoreDCUPlaneWireAngle().

G4double QweakSimVDC_DriftCellBackSD::DCVPlaneWireAngle = -60.0*degree
staticprivate

Definition at line 90 of file QweakSimVDC_DriftCellBackSD.hh.

Referenced by ProcessHits(), and StoreDCVPlaneWireAngle().

G4double QweakSimVDC_DriftCellBackSD::DCWidthOnFrame = 11.0*mm
staticprivate

Definition at line 87 of file QweakSimVDC_DriftCellBackSD.hh.

Referenced by ProcessHits(), and StoreDCWidthOnFrame().

QweakSimVDC* QweakSimVDC_DriftCellBackSD::pQweakSimVDCSetup
private

Definition at line 92 of file QweakSimVDC_DriftCellBackSD.hh.


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