ACTARSim
ACTAR TPC Simulation Reference Guide
ActarSimSciRingSD.cc
Go to the documentation of this file.
1 // - AUTHOR: Hector Alvarez-Pol 04/2008
2 /******************************************************************
3  * Copyright (C) 2005-2016, Hector Alvarez-Pol *
4  * All rights reserved. *
5  * *
6  * License according to GNU LESSER GPL (see lgpl-3.0.txt). *
7  * For the list of contributors see CREDITS. *
8  ******************************************************************/
9 //////////////////////////////////////////////////////////////////
10 /// \class ActarSimSciRingSD
11 /// SD for the Scintillators
12 /////////////////////////////////////////////////////////////////
13 
14 #include "ActarSimSciRingSD.hh"
15 
16 #include "G4HCofThisEvent.hh"
17 #include "G4Step.hh"
18 #include "G4ThreeVector.hh"
19 #include "G4SDManager.hh"
20 #include "G4ios.hh"
21 #include "G4UnitsTable.hh"
22 
23 #include "G4VProcess.hh"
24 #include "G4VPhysicalVolume.hh"
25 #include "G4TouchableHistory.hh"
26 #include "G4VTouchable.hh"
27 
28 //////////////////////////////////////////////////////////////////
29 /// Constructor, just naming the Hit collection
31  :G4VSensitiveDetector(name){
32  G4String HCname;
33  collectionName.insert(HCname="SciRingCollection");
34 }
35 
36 //////////////////////////////////////////////////////////////////
37 /// Destructor, nothing to do
39 }
40 
41 //////////////////////////////////////////////////////////////////
42 /// Initializing the ActarSimSciGeantHitsCollection object
43 /// Invoked automatically at the beggining of each event
44 void ActarSimSciRingSD::Initialize(G4HCofThisEvent* HCE){
46  (SensitiveDetectorName,collectionName[0]);
47  static G4int HCID = -1;
48  if(HCID<0)
49  HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
50 
51  HCE->AddHitsCollection( HCID, hitsCollection );
52 }
53 
54 //////////////////////////////////////////////////////////////////
55 /// Filling the ActarSimSciGeantHit information with the step info
56 /// Invoked by G4SteppingManager for each step
57 G4bool ActarSimSciRingSD::ProcessHits(G4Step* aStep,G4TouchableHistory*){
58  G4double edep = aStep->GetTotalEnergyDeposit();
59 
60  if(edep==0.) return false;
61 
63 
64  newHit->SetEdep(edep);
65 
66  newHit->SetPos(aStep->GetPostStepPoint()->GetPosition());
67  //The elements should be taken of the PostStep, the real parameters of the
68  //interaction. The PreStep keeps record of where the particle comes from
69  newHit->SetPrePos(aStep->GetPreStepPoint()->GetPosition());
70 
71  //for the local translation, x = R-1 (x_m - T)
72  G4ThreeVector theLocalPos =
73  (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetObjectRotationValue().inverse()) *
74  (aStep->GetPostStepPoint()->GetPosition()-aStep->GetPostStepPoint()->GetPhysicalVolume()->GetObjectTranslation());
75  newHit->SetLocalPos(theLocalPos);
76 
77  theLocalPos = //just using the same variable...
78  (aStep->GetPreStepPoint()->GetPhysicalVolume()->GetObjectRotationValue().inverse())*
79  (aStep->GetPreStepPoint()->GetPosition()-aStep->GetPreStepPoint()->GetPhysicalVolume()->GetObjectTranslation());
80  newHit->SetLocalPrePos(theLocalPos);
81 
82  newHit->SetDetName(aStep->GetTrack()->GetVolume()->GetName());
83  newHit->SetPreDetName(aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName());
84  newHit->SetPostDetName(aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName());
85  newHit->SetDetID(aStep->GetTrack()->GetVolume()->GetCopyNo());
86 
87  newHit->SetToF(aStep->GetPostStepPoint()->GetGlobalTime());
88 
89  newHit->SetTrackID(aStep->GetTrack()->GetTrackID());
90  newHit->SetParentID(aStep->GetTrack()->GetParentID());
91  newHit->SetParticleID(aStep->GetTrack()->GetDefinition()->GetPDGEncoding());
92  newHit->SetParticleCharge(aStep->GetTrack()->GetDefinition()->GetPDGCharge());
93  newHit->SetParticleMass(aStep->GetTrack()->GetDefinition()->GetPDGMass());
94 
95  hitsCollection->insert(newHit);
96 
97  // newHit cannot be deleted here !
98  // It should be deleted after the end of the event
99 
100  return true;
101 }
102 
103 //////////////////////////////////////////////////////////////////
104 /// Just prints and draws the event hits (class ActarSimSciGeantHit)
105 /// The recollection of the hits energy deposition in the plastic
106 /// is done in the ActarSimROOTAnalysis::EndOfEventAction()
107 void ActarSimSciRingSD::EndOfEvent(G4HCofThisEvent*){
108  G4int NbHits = hitsCollection->entries();
109  if (verboseLevel>0) {
110  G4cout << "Hits Collection: in this event they are " << NbHits
111  << " (GEANT-like) hits in the SciRing: " << G4endl;
112  for (G4int i=0;i<NbHits;i++) (*hitsCollection)[i]->Print();
113  }
114  //for (G4int i=0;i<NbHits;i++) (*hitsCollection)[i]->Draw();
115 }
void SetLocalPrePos(G4ThreeVector xyz)
void SetPrePos(G4ThreeVector xyz)
void EndOfEvent(G4HCofThisEvent *)
~ActarSimSciRingSD()
Destructor, nothing to do.
ActarSimSciRingSD(G4String)
Constructor, just naming the Hit collection.
void SetPos(G4ThreeVector xyz)
ActarSimSciRingGeantHitsCollection * hitsCollection
Geant step-like hits collect.
void SetParticleCharge(G4double charge)
void SetPostDetName(G4String Name)
G4bool ProcessHits(G4Step *, G4TouchableHistory *)
G4THitsCollection< ActarSimSciRingGeantHit > ActarSimSciRingGeantHitsCollection
void SetLocalPos(G4ThreeVector xyz)
void SetParticleMass(G4double mass)
void SetPreDetName(G4String Name)
void Initialize(G4HCofThisEvent *)