ACTARSim
ACTAR TPC Simulation Reference Guide
ActarSimSilSD.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 ActarSimSilSD
11 /// SD for the Silicons
12 /////////////////////////////////////////////////////////////////
13 
14 #include "ActarSimSilSD.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 #include "G4PhysicalConstants.hh"
29 #include "G4SystemOfUnits.hh"
30 
31 //////////////////////////////////////////////////////////////////
32 /// Constructor, just naming the Hit collection
34  :G4VSensitiveDetector(name){
35  G4String HCname;
36  collectionName.insert(HCname="SilCollection");
37 }
38 
39 //////////////////////////////////////////////////////////////////
40 /// Destructor, nothing to do
42 }
43 
44 //////////////////////////////////////////////////////////////////
45 /// Initializing the ActarSimSilGeantHitsCollection object
46 /// Invoked automatically at the beggining of each event
47 void ActarSimSilSD::Initialize(G4HCofThisEvent* HCE){
49  (SensitiveDetectorName,collectionName[0]);
50  static G4int HCID = -1;
51  if(HCID<0)
52  HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
53 
54  HCE->AddHitsCollection( HCID, hitsCollection );
55 }
56 
57 //////////////////////////////////////////////////////////////////
58 /// Filling the ActarSimSilGeantHit information with the step info
59 /// Invoked by G4SteppingManager for each step
60 G4bool ActarSimSilSD::ProcessHits(G4Step* aStep,G4TouchableHistory*){
61  //G4double edep = aStep->GetTotalEnergyDeposit();
62  G4double edep = -aStep->GetDeltaEnergy()/MeV;
63 
64  if(edep==0.) return false;
65 
67 
68  newHit->SetEdep(edep/MeV);
69  newHit->SetEBeforeSil(aStep->GetPreStepPoint()->GetKineticEnergy()/MeV);
70  newHit->SetEAfterSil(aStep->GetPostStepPoint()->GetKineticEnergy()/MeV);
71 
72  newHit->SetPos(aStep->GetPostStepPoint()->GetPosition());
73  //The elements should be taken of the PostStep, the real parameters of the
74  //interaction. The PreStep keeps record of where the particle comes from
75  newHit->SetPrePos(aStep->GetPreStepPoint()->GetPosition());
76 
77  //for the local translation, x = R-1 (x_m - T)
78  G4ThreeVector theLocalPos =
79  (aStep->GetPostStepPoint()->GetPhysicalVolume()->GetObjectRotationValue().inverse()) *
80  (aStep->GetPostStepPoint()->GetPosition()-aStep->GetPostStepPoint()->GetPhysicalVolume()->GetObjectTranslation());
81  newHit->SetLocalPos(theLocalPos);
82 
83  theLocalPos = //just using the same variable...
84  (aStep->GetPreStepPoint()->GetPhysicalVolume()->GetObjectRotationValue().inverse())*
85  (aStep->GetPreStepPoint()->GetPosition()-aStep->GetPreStepPoint()->GetPhysicalVolume()->GetObjectTranslation());
86  newHit->SetLocalPrePos(theLocalPos);
87 
88  newHit->SetDetName(aStep->GetTrack()->GetVolume()->GetName());
89  newHit->SetPreDetName(aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName());
90  newHit->SetPostDetName(aStep->GetPostStepPoint()->GetPhysicalVolume()->GetName());
91  newHit->SetDetID(aStep->GetTrack()->GetVolume()->GetCopyNo());
92 
93  newHit->SetToF(aStep->GetPostStepPoint()->GetGlobalTime());
94 
95  newHit->SetTrackID(aStep->GetTrack()->GetTrackID());
96  newHit->SetParentID(aStep->GetTrack()->GetParentID());
97  newHit->SetParticleID(aStep->GetTrack()->GetDefinition()->GetPDGEncoding());
98  newHit->SetParticleCharge(aStep->GetTrack()->GetDefinition()->GetPDGCharge());
99  newHit->SetParticleMass(aStep->GetTrack()->GetDefinition()->GetPDGMass());
100 
101  hitsCollection->insert(newHit);
102 
103  // newHit cannot be deleted here !
104  // It should be deleted after the end of the event
105 
106  return true;
107 }
108 
109 //////////////////////////////////////////////////////////////////
110 /// Just prints and draws the event hits (class ActarSimSilGeantHit)
111 /// The recollection of the hits energy deposition in the plastic
112 /// is done in the ActarSimROOTAnalysis::EndOfEventAction()
113 void ActarSimSilSD::EndOfEvent(G4HCofThisEvent*){
114  G4int NbHits = hitsCollection->entries();
115  if (verboseLevel>0) {
116  G4cout << "Hits Collection: in this event they are " << NbHits
117  << " (GEANT-like) hits in the Sil: " << G4endl;
118  for (G4int i=0;i<NbHits;i++) (*hitsCollection)[i]->Print();
119  }
120  //for (G4int i=0;i<NbHits;i++) (*hitsCollection)[i]->Draw();
121 }
void SetPrePos(G4ThreeVector xyz)
void SetPostDetName(G4String Name)
void SetParticleID(G4int ID)
ActarSimSilSD(G4String)
Constructor, just naming the Hit collection.
void SetDetName(G4String Name)
void SetLocalPos(G4ThreeVector xyz)
void SetLocalPrePos(G4ThreeVector xyz)
void SetEBeforeSil(G4double eb)
void SetEAfterSil(G4double ea)
void SetParticleMass(G4double mass)
void SetPreDetName(G4String Name)
void SetParticleCharge(G4double charge)
void EndOfEvent(G4HCofThisEvent *)
void SetToF(G4double Time)
G4THitsCollection< ActarSimSilGeantHit > ActarSimSilGeantHitsCollection
void Initialize(G4HCofThisEvent *)
void SetParentID(G4int id)
void SetPos(G4ThreeVector xyz)
~ActarSimSilSD()
Destructor, nothing to do.
G4bool ProcessHits(G4Step *, G4TouchableHistory *)
void SetEdep(G4double de)
ActarSimSilGeantHitsCollection * hitsCollection
Geant step-like hits collect.
void SetTrackID(G4int track)