CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HFChamberSD.cc
Go to the documentation of this file.
3 
4 #include "G4VPhysicalVolume.hh"
5 #include "G4PVPlacement.hh"
6 #include "G4HCofThisEvent.hh"
7 #include "G4TouchableHistory.hh"
8 #include "G4Track.hh"
9 #include "G4Step.hh"
10 #include "G4VSolid.hh"
11 #include "G4DynamicParticle.hh"
12 #include "G4ParticleDefinition.hh"
13 #include "G4SDManager.hh"
14 #include "G4ios.hh"
15 
16 #include "G4PhysicalConstants.hh"
17 #include "G4SystemOfUnits.hh"
18 
21  const SimTrackManager* manager) :
22  SensitiveCaloDetector(name, cpv, clg, p), theName(name),
23  m_trackManager(manager), theHCID(-1), theHC(0), theNSteps(0) {
24 
25  collectionName.insert(name);
26  LogDebug("FiberSim") << "***************************************************"
27  << "\n"
28  << "* *"
29  << "\n"
30  << "* Constructing a HFChamberSD with name " << GetName()
31  << "\n"
32  << "* *"
33  << "\n"
34  << "***************************************************";
35  //
36  // Now attach the right detectors (LogicalVolumes) to me
37  //
38  std::vector<std::string> lvNames = clg.logicalNames(name);
39  this->Register();
40  for (std::vector<std::string>::iterator it=lvNames.begin();
41  it !=lvNames.end(); it++){
42  this->AssignSD(*it);
43  LogDebug("FiberSim") << "HFChamberSD : Assigns SD to LV " << (*it);
44  }
45 }
46 
48  if (theHC) delete theHC;
49 }
50 
51 void HFChamberSD::Initialize(G4HCofThisEvent * HCE) {
52 
53  LogDebug("FiberSim") << "HFChamberSD : Initialize called for " << GetName();
54  theHC = new HFShowerG4HitsCollection(GetName(), collectionName[0]);
55  if (theHCID<0)
56  theHCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
57  HCE->AddHitsCollection(theHCID, theHC);
58 
59 }
60 
61 G4bool HFChamberSD::ProcessHits(G4Step * aStep, G4TouchableHistory*) {
62 
63  //do not process hits other than primary particle hits:
64  double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge();
65  int trackID = aStep->GetTrack()->GetTrackID();
66  if(charge == 0. || trackID != 1 ||aStep->GetTrack()->GetParentID() != 0 || aStep->GetTrack()->GetCreatorProcess() != NULL) return false;
67  ++theNSteps;
68  //if(theNSteps>1)return false;
69 
70  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
71  const G4VTouchable* touch = preStepPoint->GetTouchable();
72  int detID = setDetUnitId(aStep);
73 
74  double edep = aStep->GetTotalEnergyDeposit();
75  double time = (preStepPoint->GetGlobalTime())/ns;
76 
77  G4ThreeVector globalPos = preStepPoint->GetPosition();
78  G4ThreeVector localPos = touch->GetHistory()->GetTopTransform().TransformPoint(globalPos);
79  const G4DynamicParticle* particle = aStep->GetTrack()->GetDynamicParticle();
80  G4ThreeVector momDir = particle->GetMomentumDirection();
81 
82  HFShowerG4Hit *aHit = new HFShowerG4Hit(detID, trackID, edep, time);
83  aHit->setLocalPos(localPos);
84  aHit->setGlobalPos(globalPos);
85  aHit->setPrimMomDir(momDir);
86 
87  LogDebug("FiberSim") << "HFChamberSD: Hit created in ("
88  << touch->GetVolume(0)->GetLogicalVolume()->GetName()
89  << ") " << " ID " << detID << " Track " << trackID
90  << " Edep: " << edep/MeV << " MeV; Time: " << time
91  << " ns; Position (local) " << localPos << " (global ) "
92  << globalPos << " direction " << momDir;
93 
94  theHC->insert(aHit);
95  return true;
96 }
97 
98 void HFChamberSD::EndOfEvent(G4HCofThisEvent * HCE) {
99 
100  LogDebug("FiberSim") << "HFChamberSD: Sees" << theHC->entries() << " hits";
101  clear();
102 }
103 
105  theNSteps = 0;
106 }
107 
109 
111 
113 
114 uint32_t HFChamberSD::setDetUnitId(G4Step* aStep) {
115  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
116  return (touch->GetReplicaNumber(0));
117 }
118 
#define LogDebug(id)
std::vector< PCaloHit > PCaloHitContainer
std::vector< std::string > logicalNames(std::string &readoutName)
virtual void clearHits()
Definition: HFChamberSD.cc:112
void setLocalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:47
#define NULL
Definition: scimark2.h:8
virtual uint32_t setDetUnitId(G4Step *)
Definition: HFChamberSD.cc:114
virtual void fillHits(edm::PCaloHitContainer &, std::string)
Definition: HFChamberSD.cc:119
type of data representation of DDCompactView
Definition: DDCompactView.h:77
double charge(const std::vector< uint8_t > &Ampls)
virtual void PrintAll()
Definition: HFChamberSD.cc:110
const double MeV
HFShowerG4HitsCollection * theHC
Definition: HFChamberSD.h:46
std::string const collectionName[nCollections]
Definition: Collections.h:45
virtual void EndOfEvent(G4HCofThisEvent *HCE)
Definition: HFChamberSD.cc:98
virtual void DrawAll()
Definition: HFChamberSD.cc:108
void setPrimMomDir(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:49
G4int theHCID
Definition: HFChamberSD.h:45
HFChamberSD(std::string, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HFChamberSD.cc:19
virtual void AssignSD(std::string &vname)
void setGlobalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:48
virtual void clear()
Definition: HFChamberSD.cc:104
virtual void Initialize(G4HCofThisEvent *HCE)
Definition: HFChamberSD.cc:51
G4THitsCollection< HFShowerG4Hit > HFShowerG4HitsCollection
Definition: HFShowerG4Hit.h:60
virtual ~HFChamberSD()
Definition: HFChamberSD.cc:47
virtual G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *ROhist)
Definition: HFChamberSD.cc:61