ACTARSim
ACTAR TPC Simulation Reference Guide
ActarSimGasSD.cc
Go to the documentation of this file.
1 // - AUTHOR: Hector Alvarez-Pol 04/2006
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 ActarSimGasSD
11 /// SD for the gas volume in the detector
12 /////////////////////////////////////////////////////////////////
13 
14 #include "ActarSimGasSD.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 "G4VPhysicalVolume.hh"
24 #include "G4TouchableHistory.hh"
25 #include "G4VTouchable.hh"
26 
27 #include "G4PhysicalConstants.hh"
28 #include "G4SystemOfUnits.hh"
29 
30 //////////////////////////////////////////////////////////////////
31 /// Constructor
32 ActarSimGasSD::ActarSimGasSD(G4String name):G4VSensitiveDetector(name){
33  G4String HCname;
34  collectionName.insert(HCname="gasCollection");
35 }
36 
37 //////////////////////////////////////////////////////////////////
38 /// Destructor
40 }
41 
42 //////////////////////////////////////////////////////////////////
43 /// Initializing the ActarSimCalGeantHitsCollection object
44 /// Invoked automatically at the beggining of each event
45 void ActarSimGasSD::Initialize(G4HCofThisEvent* HCE){
46 
48  (SensitiveDetectorName,collectionName[0]);
49  static G4int HCID = -1;
50  if(HCID<0)
51  HCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
52 
53  HCE->AddHitsCollection( HCID, hitsCollection );
54 }
55 
56 //////////////////////////////////////////////////////////////////
57 /// Filling the ActarSimCalGeantHit information with the step info.
58 /// Invoked by G4SteppingManager for each step
59 G4bool ActarSimGasSD::ProcessHits(G4Step* aStep,G4TouchableHistory*){
60 
61  //G4double edep = aStep->GetTotalEnergyDeposit()/MeV;
62  G4double edep = -aStep->GetDeltaEnergy()/MeV;
63 
64  if(edep==0.) return false;
65 
67 
68  newHit->SetTrackID(aStep->GetTrack()->GetTrackID());
69  newHit->SetParentID(aStep->GetTrack()->GetParentID());
70  //newHit->SetStep(aStep);
71 
72  newHit->SetEdep(edep);
73  newHit->SetParticleCharge(aStep->GetTrack()->GetDefinition()->GetPDGCharge());
74  newHit->SetParticleMass(aStep->GetTrack()->GetDefinition()->GetPDGMass());
75  newHit->SetParticleID(aStep->GetTrack()->GetDefinition()->GetPDGEncoding());
76  newHit->SetPrePos(aStep->GetPreStepPoint()->GetPosition()/mm);
77  newHit->SetPostPos(aStep->GetPostStepPoint()->GetPosition()/mm);
78  newHit->SetPreToF(aStep->GetPreStepPoint()->GetGlobalTime()/ns);
79  newHit->SetPostToF(aStep->GetPostStepPoint()->GetGlobalTime()/ns);
80  newHit->SetStepLength(aStep->GetStepLength()/mm);
81  newHit->SetStepEnergy(aStep->GetTrack()->GetKineticEnergy()/MeV);
82 
83  newHit->SetDetName(aStep->GetTrack()->GetVolume()->GetName());
84  newHit->SetDetID(aStep->GetTrack()->GetVolume()->GetCopyNo());
85 
86  //newHit->Print();
87 
88  hitsCollection->insert(newHit);
89 
90  // newHit cannot be deleted here !
91  // It should be deleted after the end of the event
92 
93  return true;
94 }
95 
96 //////////////////////////////////////////////////////////////////
97 /// Just prints and draws the event hits (class ActarSimGasGeantHit).
98 /// The recollection of the hits energy deposition in a crystal
99 /// is done in the ActarSimAnalysis::EndOfEventAction()
100 void ActarSimGasSD::EndOfEvent(G4HCofThisEvent*){
101 
102  G4int NbHits = hitsCollection->entries();
103  if (verboseLevel>0) {
104  G4cout << "Hits Collection: in this event they are " << NbHits
105  << " (GEANT-like) hits in the gas volume: " << G4endl;
106  for (G4int i=0;i<NbHits;i++) (*hitsCollection)[i]->Print();
107  }
108  //for (G4int i=0;i<NbHits;i++) (*hitsCollection)[i]->Draw();
109 }
void SetParentID(G4int track)
ActarSimGasSD(G4String)
Constructor.
void SetStepLength(G4double len)
void SetPostToF(G4double Time)
ActarSimGasGeantHitsCollection * hitsCollection
Geant step-like hits collect.
void SetEdep(G4double de)
G4THitsCollection< ActarSimGasGeantHit > ActarSimGasGeantHitsCollection
void SetTrackID(G4int track)
~ActarSimGasSD()
Destructor.
void SetParticleCharge(G4double pc)
void SetDetName(G4String Name)
void SetPrePos(G4ThreeVector xyz)
void SetStepEnergy(G4double sten)
void SetPostPos(G4ThreeVector xyz)
void SetParticleMass(G4double pm)
void EndOfEvent(G4HCofThisEvent *)
void SetParticleID(G4int pi)
void SetPreToF(G4double Time)
void Initialize(G4HCofThisEvent *)
G4bool ProcessHits(G4Step *, G4TouchableHistory *)