CMS 3D CMS Logo

HGCalTBMB.cc
Go to the documentation of this file.
2 
7 
10 
11 #include "G4LogicalVolumeStore.hh"
12 #include "G4Step.hh"
13 #include "G4Track.hh"
14 
15 #include <iostream>
16 
17 //#define EDM_ML_DEBUG
18 
20 
21  edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("HGCalTBMB");
22  listNames_ = m_p.getParameter<std::vector<std::string> >("DetectorNames");
23  stopName_ = m_p.getParameter<std::string>("StopName");
24  stopZ_ = m_p.getParameter<double>("MaximumZ");
25  nList_ = listNames_.size();
26  edm::LogInfo("HGCSim") << "HGCalTBMB initialized for " << nList_ <<" volumes\n";
27  for (unsigned int k=0; k<nList_; ++k)
28  edm::LogInfo("HGCSim") << " [" << k << "] " << listNames_[k] << std::endl;
29  edm::LogInfo("HGCSim") << "Stop after " << stopZ_ << " or reaching volume "
30  << stopName_ << std::endl;
31 
33  if ( !tfile.isAvailable() )
34  throw cms::Exception("BadConfig") << "TFileService unavailable: "
35  << "please add it to config file";
36  char name[20], title[80];
37  TH1D* hist;
38  for (unsigned int i=0; i<=nList_; i++) {
39  std::string named = (i == nList_) ? "Total" : listNames_[i];
40  sprintf(name, "RadL%d", i);
41  sprintf(title, "MB(X0) for (%s)", named.c_str());
42  hist = tfile->make<TH1D>(name,title,1000,0.0,100.0);
43  hist->Sumw2(true);
44  me100_.push_back(hist);
45  sprintf(name, "IntL%d", i);
46  sprintf(title, "MB(L0) for (%s)", named.c_str());
47  hist = tfile->make<TH1D>(name,title,1000,0.0,10.0);
48  hist->Sumw2(true);
49  me200_.push_back(hist);
50  sprintf(name, "StepL%d", i);
51  sprintf(title, "MB(Step) for (%s)", named.c_str());
52  hist = tfile->make<TH1D>(name,title,1000,0.0,50000.0);
53  hist->Sumw2(true);
54  me300_.push_back(hist);
55  }
56  edm::LogInfo("HGCSim") << "HGCalTBMB: Booking user histos done ===";
57 }
58 
60 
61 
62 void HGCalTBMB::update(const BeginOfTrack* trk) {
63 
64  radLen_ = std::vector<double>(nList_+1,0);
65  intLen_ = std::vector<double>(nList_+1,0);
66  stepLen_ = std::vector<double>(nList_+1,0);
67 
68 #ifdef EDM_ML_DEBUG
69  const G4Track * aTrack = (*trk)(); // recover G4 pointer if wanted
70  const G4ThreeVector& mom = aTrack->GetMomentum() ;
71  double theEnergy = aTrack->GetTotalEnergy();
72  int theID = (int)(aTrack->GetDefinition()->GetPDGEncoding());
73  std::cout << "MaterialBudgetHcalHistos: Track " << aTrack->GetTrackID()
74  << " Code " << theID << " Energy " <<theEnergy/CLHEP::GeV
75  << " GeV; Momentum " << mom << std::endl;
76 #endif
77 }
78 
79 void HGCalTBMB::update(const G4Step* aStep) {
80 
81  G4Material * material = aStep->GetPreStepPoint()->GetMaterial();
82  double step = aStep->GetStepLength();
83  double radl = material->GetRadlen();
84  double intl = material->GetNuclearInterLength();
85 
86  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
87  int indx = findVolume(touch,false);
88 
89  if (indx >= 0) {
90  stepLen_[indx] += step;
91  radLen_[indx] += (step/radl);
92  intLen_[indx] += (step/intl);
93  }
94  stepLen_[nList_] += step;
95  radLen_[nList_] += (step/radl);
96  intLen_[nList_] += (step/intl);
97 #ifdef EDM_ML_DEBUG
98  std::cout << "HGCalTBMB::Step in "
99  << touch->GetVolume(0)->GetLogicalVolume()->GetName()
100  << " Index " << indx <<" Step " << step << " RadL " << step/radl
101  << " IntL " << step/intl << std::endl;
102 #endif
103 
104  if (stopAfter(aStep)) {
105  G4Track* track = aStep->GetTrack();
106  track->SetTrackStatus(fStopAndKill);
107  }
108 }
109 
110 void HGCalTBMB::update(const EndOfTrack* trk) {
111 
112  for (unsigned int ii=0; ii<=nList_; ++ii) {
113  me100_[ii]->Fill(radLen_[ii]);
114  me200_[ii]->Fill(intLen_[ii]);
115  me300_[ii]->Fill(stepLen_[ii]);
116 #ifdef EDM_ML_DEBUG
117  std::string name("Total");
118  if (ii < nList_) name = listNames_[ii];
119  std::cout << "HGCalTBMB::Volume[" << ii << "]: " << name << " == Step "
120  << stepLen_[ii] << " RadL " << radLen_[ii] << " IntL "
121  << intLen_[ii] << std::endl;
122 #endif
123  }
124 }
125 
126 bool HGCalTBMB::stopAfter(const G4Step* aStep) {
127 
128  bool flag(false);
129  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
130  G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
131  if (aStep->GetPostStepPoint() != nullptr)
132  hitPoint = aStep->GetPostStepPoint()->GetPosition();
133  double zz = hitPoint.z();
134 
135  if ((findVolume(touch,true) == 0) || (zz > stopZ_)) flag = true;
136 #ifdef EDM_ML_DEBUG
137  std::cout << " HGCalTBMB::Name " << touch->GetVolume(0)->GetName() << " z "
138  << zz << " Flag" << flag << std::endl;
139 #endif
140  return flag;
141 }
142 
143 int HGCalTBMB::findVolume(const G4VTouchable* touch, bool stop) const {
144 
145  int ivol =-1;
146  int level = (touch->GetHistoryDepth())+1;
147  for (int ii = 0; ii < level; ii++) {
148  std::string name = touch->GetVolume(ii)->GetName();
149  if (stop) {
150  if (strcmp(name.c_str(),stopName_.c_str()) == 0) ivol = 0;
151  } else {
152  for (unsigned int k=0; k<nList_; ++k) {
153  if (strcmp(name.c_str(),listNames_[k].c_str()) == 0) {
154  ivol = k; break;
155  }
156  }
157  }
158  if (ivol >= 0) break;
159  }
160  return ivol;
161 }
162 
T getParameter(std::string const &) const
const double GeV
Definition: MathUtil.h:16
~HGCalTBMB() override
Definition: HGCalTBMB.cc:59
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
int findVolume(const G4VTouchable *touch, bool stop) const
Definition: HGCalTBMB.cc:143
bool stopAfter(const G4Step *)
Definition: HGCalTBMB.cc:126
double stopZ_
Definition: HGCalTBMB.h:41
std::vector< double > intLen_
Definition: HGCalTBMB.h:43
std::string stopName_
Definition: HGCalTBMB.h:40
std::vector< TH1D * > me100_
Definition: HGCalTBMB.h:44
bool isAvailable() const
Definition: Service.h:46
unsigned int nList_
Definition: HGCalTBMB.h:42
ii
Definition: cuy.py:590
std::vector< TH1D * > me300_
Definition: HGCalTBMB.h:44
int k[5][pyjets_maxn]
void update(const BeginOfTrack *) override
This routine will be called when the appropriate signal arrives.
Definition: HGCalTBMB.cc:62
HGCalTBMB(const edm::ParameterSet &)
Definition: HGCalTBMB.cc:19
std::vector< TH1D * > me200_
Definition: HGCalTBMB.h:44
std::vector< double > radLen_
Definition: HGCalTBMB.h:43
step
std::vector< std::string > listNames_
Definition: HGCalTBMB.h:39
std::vector< double > stepLen_
Definition: HGCalTBMB.h:43