CMS 3D CMS Logo

MaterialBudgetForward.cc
Go to the documentation of this file.
2 
7 
11 
12 #include "G4LogicalVolumeStore.hh"
13 #include "G4Step.hh"
14 #include "G4Track.hh"
15 
16 #include <iostream>
17 //#define DebugLog
19 
21 
22  edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("MaterialBudgetForward");
23  detTypes = m_p.getParameter<std::vector<std::string> >("DetectorTypes");
24  constituents = m_p.getParameter<std::vector<int> >("Constituents");
25  stackOrder = m_p.getParameter<std::vector<int> >("StackOrder");
26  detNames = m_p.getParameter<std::vector<std::string> >("DetectorNames");
27  detLevels = m_p.getParameter<std::vector<int> >("DetectorLevels");
28  etaRegions = m_p.getParameter<std::vector<double> >("EtaBoundaries");
29  regionTypes = m_p.getParameter<std::vector<int> >("RegionTypes");
30  boundaries = m_p.getParameter<std::vector<double> >("Boundaries");
31  edm::LogInfo("MaterialBudget") << "MaterialBudgetForward initialized for "
32  << detTypes.size() << " detector types";
33  unsigned int nc = 0;
34  for (unsigned int ii=0; ii<detTypes.size(); ++ii) {
35  edm::LogInfo("MaterialBudget") << "Type[" << ii << "] : " << detTypes[ii]
36  << " with " << constituents[ii] <<", order "
37  << stackOrder[ii] << " constituents --> ";
38  for (int kk=0; kk<constituents[ii]; ++kk) {
39  std::string name = "Unknown"; int level = -1;
40  if (nc < detNames.size()) {
41  name = detNames[nc]; level = detLevels[nc]; ++nc;
42  }
43  edm::LogInfo("MaterialBudget") << " Constituent[" << kk << "]: "
44  << name << " at level " << level;
45  }
46  }
47  edm::LogInfo("MaterialBudget") << "MaterialBudgetForward Stop condition for "
48  << etaRegions.size() << " eta regions";
49  for (unsigned int ii=0; ii<etaRegions.size(); ++ii) {
50  edm::LogInfo("MaterialBudget") << "Region[" << ii << "] : Eta < "
51  << etaRegions[ii] << " boundary type "
52  << regionTypes[ii] << " limit "
53  << boundaries[ii];
54  }
55  book(m_p);
56 }
57 
59 }
60 
62 
63  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
64  std::vector<G4LogicalVolume *>::const_iterator lvcite;
65 
66  unsigned int kount=detNames.size();
67  for (unsigned int ii=0; ii<kount; ++ii)
68  logVolumes.push_back(nullptr);
69 
70  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++) {
71  for (unsigned int ii=0; ii<detNames.size(); ++ii) {
72  if (strcmp(detNames[ii].c_str(),(*lvcite)->GetName().c_str()) == 0) {
73  logVolumes[ii] = (*lvcite);
74  kount--;
75  break;
76  }
77  }
78  if (kount <= 0) break;
79  }
80  edm::LogInfo("MaterialBudget") << "MaterialBudgetForward::Finds "
81  << detNames.size()-kount << " out of "
82  << detNames.size() << " LV addresses";
83  for (unsigned int ii=0; ii<detNames.size(); ++ii) {
84  std::string name("Unknown");
85  if (logVolumes[ii] != nullptr) name = logVolumes[ii]->GetName();
86  edm::LogInfo("MaterialBudget") << "LV[" << ii << "] : " << detNames[ii]
87  << " Address " << logVolumes[ii] << " | "
88  << name;
89  }
90 
91  for (unsigned int ii=0; ii<(detTypes.size()+1); ++ii) {
92  stepLen.push_back(0); radLen.push_back(0); intLen.push_back(0);
93  }
94  stackOrder.push_back(0);
95 }
96 
98 
99  for (unsigned int ii=0; ii<(detTypes.size()+1); ++ii) {
100  stepLen[ii] = 0; radLen[ii] = 0; intLen[ii] = 0;
101  }
102 
103  const G4Track * aTrack = (*trk)(); // recover G4 pointer if wanted
104  const G4ThreeVector& mom = aTrack->GetMomentum() ;
105  if (mom.theta() != 0 ) {
106  eta = mom.eta();
107  } else {
108  eta = -99;
109  }
110  phi = mom.phi();
111  stepT = 0;
112 
113 #ifdef DebugLog
114  double theEnergy = aTrack->GetTotalEnergy();
115  int theID = (int)(aTrack->GetDefinition()->GetPDGEncoding());
116  edm::LogInfo("MaterialBudget") << "MaterialBudgetHcalHistos: Track "
117  << aTrack->GetTrackID() << " Code " << theID
118  << " Energy " <<theEnergy/CLHEP::GeV
119  << " GeV; Eta " << eta << " Phi "
120  << phi/CLHEP::deg << " PT "
121  << mom.perp()/CLHEP::GeV << " GeV *****";
122 #endif
123 }
124 
125 void MaterialBudgetForward::update(const G4Step* aStep) {
126 
127  //---------- each step
128  G4Material * material = aStep->GetPreStepPoint()->GetMaterial();
129  double step = aStep->GetStepLength();
130  double radl = material->GetRadlen();
131  double intl = material->GetNuclearInterLength();
132 
133  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
134  unsigned int indx = detTypes.size();
135  unsigned int nc = 0;
136  for (unsigned int ii=0; ii<detTypes.size(); ++ii) {
137  for (int kk=0; kk<constituents[ii]; ++kk) {
138  if (detLevels[nc+kk] <= (int)((touch->GetHistoryDepth())+1)) {
139  int jj = (int)((touch->GetHistoryDepth())+1)-detLevels[nc+kk];
140  if ((touch->GetVolume(jj)->GetLogicalVolume()) == logVolumes[nc+kk]) {
141  indx = ii;
142  break;
143  }
144  }
145  }
146  nc += (unsigned int)(constituents[ii]);
147  if (indx == ii) break;
148  }
149 
150  stepT += step;
151  stepLen[indx] = stepT;
152  radLen[indx] += (step/radl);
153  intLen[indx] += (step/intl);
154 #ifdef DebugLog
155  edm::LogInfo("MaterialBudget") << "MaterialBudgetForward::Step in "
156  << touch->GetVolume(0)->GetLogicalVolume()->GetName()
157  << " Index " << indx <<" Step " << step
158  << " RadL " << step/radl << " IntL "
159  << step/intl;
160 #endif
161  //----- Stop tracking after selected position
162  if (stopAfter(aStep)) {
163  G4Track* track = aStep->GetTrack();
164  track->SetTrackStatus( fStopAndKill );
165  }
166 }
167 
168 
170 
171  for (unsigned int ii=0; ii<detTypes.size(); ++ii) {
172  for (unsigned int jj=0; jj<=detTypes.size(); ++jj) {
173  if (stackOrder[jj] == (int)(ii+1)) {
174  for (unsigned int kk=0; kk<=detTypes.size(); ++kk) {
175  if (stackOrder[kk] == (int)(ii)) {
176  radLen[jj] += radLen[kk];
177  intLen[jj] += intLen[kk];
178 #ifdef DebugLog
179  edm::LogInfo("MaterialBudget") << "MaterialBudgetForward::Add "
180  << kk << ":" << stackOrder[kk]
181  << " to " << jj <<":"
182  << stackOrder[jj] <<" RadL "
183  << radLen[kk] << " : " << radLen[jj]
184  << " IntL " << intLen[kk] << " : "
185  << intLen[jj] <<" StepL "
186  << stepLen[kk]<< " : "
187  << stepLen[jj];
188 #endif
189  break;
190  }
191  }
192  break;
193  }
194  }
195  }
196 
197  for (unsigned int ii=0; ii<=detTypes.size(); ++ii) {
198  me100[ii]->Fill(eta, radLen[ii]);
199  me200[ii]->Fill(eta, intLen[ii]);
200  me300[ii]->Fill(eta, stepLen[ii]);
201  me400[ii]->Fill(eta);
202  me500[ii]->Fill(eta, phi, radLen[ii]);
203  me600[ii]->Fill(eta, phi, intLen[ii]);
204  me700[ii]->Fill(eta, phi, stepLen[ii]);
205  me800[ii]->Fill(eta, phi);
206 
207  std::string name("Unknown");
208  if (ii < detTypes.size()) name = detTypes[ii];
209 #ifdef DebugLog
210  edm::LogInfo("MaterialBudget") << "MaterialBudgetForward::Volume[" << ii
211  << "]: " << name << " eta " << eta
212  << " == Step " << stepLen[ii] << " RadL "
213  << radLen[ii] << " IntL " << intLen[ii];
214 #endif
215  }
216 }
217 
219 
220  // Book histograms
222 
223  if ( !tfile.isAvailable() )
224  throw cms::Exception("BadConfig") << "TFileService unavailable: "
225  << "please add it to config file";
226 
227  int binEta = m_p.getUntrackedParameter<int>("NBinEta", 320);
228  int binPhi = m_p.getUntrackedParameter<int>("NBinPhi", 180);
229  double minEta = m_p.getUntrackedParameter<double>("MinEta",-8.0);
230  double maxEta = m_p.getUntrackedParameter<double>("MaxEta", 8.0);
231  double maxPhi = CLHEP::pi;
232  edm::LogInfo("MaterialBudget") << "MaterialBudgetForward: Booking user "
233  << "histos === with " << binEta << " bins "
234  << "in eta from " << minEta << " to "
235  << maxEta << " and " << binPhi << " bins "
236  << "in phi from " << -maxPhi << " to "
237  << maxPhi;
238 
239  char name[20], title[80];
240  std::string named;
241  for (int i=0; i<std::min((int)(detTypes.size()+1),maxSet); i++) {
242  if (i >= (int)(detTypes.size())) named = "Unknown";
243  else named = detTypes[i];
244  sprintf(name, "%d", i+100);
245  sprintf(title, "MB(X0) prof Eta in %s", named.c_str());
246  me100[i] = tfile->make<TProfile>(name, title, binEta, minEta, maxEta);
247  sprintf(name, "%d", i+200);
248  sprintf(title, "MB(L0) prof Eta in %s", named.c_str());
249  me200[i] = tfile->make<TProfile>(name, title, binEta, minEta, maxEta);
250  sprintf(name, "%d", i+300);
251  sprintf(title, "MB(Step) prof Eta in %s", named.c_str());
252  me300[i] = tfile->make<TProfile>(name, title, binEta, minEta, maxEta);
253  sprintf(name, "%d", i+400);
254  sprintf(title, "Eta in %s", named.c_str());
255  me400[i] = tfile->make<TH1F>(name, title, binEta, minEta, maxEta);
256  sprintf(name, "%d", i+500);
257  sprintf(title, "MB(X0) prof Eta Phi in %s", named.c_str());
258  me500[i] = tfile->make<TProfile2D>(name, title, binEta/2, minEta, maxEta,
259  binPhi/2, -maxPhi, maxPhi);
260  sprintf(name, "%d", i+600);
261  sprintf(title, "MB(L0) prof Eta Phi in %s", named.c_str());
262  me600[i]= tfile->make<TProfile2D>(name, title, binEta/2, minEta, maxEta,
263  binPhi/2, -maxPhi, maxPhi);
264  sprintf(name, "%d", i+700);
265  sprintf(title, "MB(Step) prof Eta Phi in %s", named.c_str());
266  me700[i]= tfile->make<TProfile2D>(name, title, binEta/2, minEta, maxEta,
267  binPhi/2, -maxPhi, maxPhi);
268  sprintf(name, "%d", i+800);
269  sprintf(title, "Eta vs Phi in %s", named.c_str());
270  me800[i]= tfile->make<TH2F>(name, title, binEta/2, minEta, maxEta,
271  binPhi/2, -maxPhi, maxPhi);
272  }
273 
274  edm::LogInfo("MaterialBudget") << "MaterialBudgetForward: Booking user "
275  << "histos done ===";
276 
277 }
278 
279 bool MaterialBudgetForward::stopAfter(const G4Step* aStep) {
280 
281  G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
282  if (aStep->GetPostStepPoint() != nullptr)
283  hitPoint = aStep->GetPostStepPoint()->GetPosition();
284  double rr = hitPoint.perp();
285  double zz = std::abs(hitPoint.z());
286 
287  bool flag(false);
288  for (unsigned int ii=0; ii<etaRegions.size(); ++ii) {
289 #ifdef DebugLog
290  edm::LogInfo("MaterialBudget") << " MaterialBudgetForward::Eta " << eta
291  << " in Region[" << ii << "] with "
292  << etaRegions[ii] << " type "
293  << regionTypes[ii] << "|" << boundaries[ii];
294 #endif
295  if (fabs(eta) < etaRegions[ii]) {
296  if (regionTypes[ii] == 0) {
297  if (rr >= boundaries[ii]-0.001) flag = true;
298  } else {
299  if (zz >= boundaries[ii]-0.001) flag = true;
300  }
301 #ifdef DebugLog
302  if (flag)
303  edm::LogInfo("MaterialBudget") <<" MaterialBudgetForward::Stop after R = "
304  << rr << " and Z = " << zz;
305 #endif
306  break;
307  }
308  }
309 #ifdef DebugLog
310  edm::LogInfo("MaterialBudget") <<" MaterialBudgetForward:: R = " << rr
311  << " and Z = " << zz << " Flag " << flag;
312 #endif
313  return flag;
314 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
TProfile2D * me700[maxSet]
const double GeV
Definition: MathUtil.h:16
TProfile * me100[maxSet]
MaterialBudgetForward(const edm::ParameterSet &)
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
double maxEta
TProfile2D * me600[maxSet]
std::vector< double > etaRegions
std::vector< int > detLevels
std::vector< double > radLen
std::vector< int > stackOrder
const Double_t pi
void update(const BeginOfRun *) override
This routine will be called when the appropriate signal arrives.
TProfile * me200[maxSet]
std::vector< double > intLen
bool isAvailable() const
Definition: Service.h:46
std::vector< double > stepLen
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< G4LogicalVolume * > logVolumes
T min(T a, T b)
Definition: MathUtil.h:58
std::vector< std::string > detNames
ii
Definition: cuy.py:590
std::vector< int > constituents
bool stopAfter(const G4Step *)
void book(const edm::ParameterSet &)
std::vector< int > regionTypes
TProfile * me300[maxSet]
std::vector< double > boundaries
TProfile2D * me500[maxSet]
step
std::vector< std::string > detTypes