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