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