CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HCalSD.cc
Go to the documentation of this file.
1 // File: HCalSD.cc
3 // Description: Sensitive Detector class for calorimeters
5 
17 
20 
21 #include "G4LogicalVolumeStore.hh"
22 #include "G4LogicalVolume.hh"
23 #include "G4Step.hh"
24 #include "G4Track.hh"
25 #include "G4ParticleTable.hh"
26 #include "G4VProcess.hh"
27 
28 #include "G4SystemOfUnits.hh"
29 #include "G4PhysicalConstants.hh"
30 
31 #include <iostream>
32 #include <fstream>
33 #include <iomanip>
34 
35 //#define DebugLog
36 //#define plotDebug
37 
38 HCalSD::HCalSD(G4String name, const DDCompactView & cpv,
39  const SensitiveDetectorCatalog & clg,
40  edm::ParameterSet const & p, const SimTrackManager* manager) :
41  CaloSD(name, cpv, clg, p, manager,
42  p.getParameter<edm::ParameterSet>("HCalSD").getParameter<int>("TimeSliceUnit"),
43  p.getParameter<edm::ParameterSet>("HCalSD").getParameter<bool>("IgnoreTrackID")),
44  numberingFromDDD(0), numberingScheme(0), showerLibrary(0), hfshower(0),
45  showerParam(0), showerPMT(0), showerBundle(0), m_HEDarkening(0),
46  m_HFDarkening(0) {
47 
48  //static SimpleConfigurable<bool> on1(false, "HCalSD:UseBirkLaw");
49  //static SimpleConfigurable<double> bk1(0.013, "HCalSD:BirkC1");
50  //static SimpleConfigurable<double> bk2(0.0568,"HCalSD:BirkC2");
51  //static SimpleConfigurable<double> bk3(1.75, "HCalSD:BirkC3");
52  // Values from NIM 80 (1970) 239-244: as implemented in Geant3
53  //static SimpleConfigurable<bool> on2(true,"HCalSD:UseShowerLibrary");
55  useBirk = m_HC.getParameter<bool>("UseBirkLaw");
56  birk1 = m_HC.getParameter<double>("BirkC1")*(g/(MeV*cm2));
57  birk2 = m_HC.getParameter<double>("BirkC2");
58  birk3 = m_HC.getParameter<double>("BirkC3");
59  useShowerLibrary = m_HC.getParameter<bool>("UseShowerLibrary");
60  useParam = m_HC.getParameter<bool>("UseParametrize");
61  testNumber = m_HC.getParameter<bool>("TestNumberingScheme");
62  usePMTHit = m_HC.getParameter<bool>("UsePMTHits");
63  betaThr = m_HC.getParameter<double>("BetaThreshold");
64  eminHitHB = m_HC.getParameter<double>("EminHitHB")*MeV;
65  eminHitHE = m_HC.getParameter<double>("EminHitHE")*MeV;
66  eminHitHO = m_HC.getParameter<double>("EminHitHO")*MeV;
67  eminHitHF = m_HC.getParameter<double>("EminHitHF")*MeV;
68  useFibreBundle = m_HC.getParameter<bool>("UseFibreBundleHits");
69  deliveredLumi = m_HC.getParameter<double>("DelivLuminosity");
70  bool ageingFlagHE= m_HC.getParameter<bool>("HEDarkening");
71  bool ageingFlagHF= m_HC.getParameter<bool>("HFDarkening");
72  useHF = m_HC.getUntrackedParameter<bool>("UseHF",true);
73  bool forTBH2 = m_HC.getUntrackedParameter<bool>("ForTBH2",false);
74  useLayerWt = m_HC.getUntrackedParameter<bool>("UseLayerWt",false);
75  std::string file = m_HC.getUntrackedParameter<std::string>("WtFile","None");
76  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShower");
77  applyFidCut = m_HF.getParameter<bool>("ApplyFiducialCut");
78 
79 #ifdef DebugLog
80  LogDebug("HcalSim") << "***************************************************"
81  << "\n"
82  << "* *"
83  << "\n"
84  << "* Constructing a HCalSD with name " << name << "\n"
85  << "* *"
86  << "\n"
87  << "***************************************************";
88 #endif
89  edm::LogInfo("HcalSim") << "HCalSD:: Use of HF code is set to " << useHF
90  << "\nUse of shower parametrization set to "
91  << useParam << "\nUse of shower library is set to "
92  << useShowerLibrary << "\nUse PMT Hit is set to "
93  << usePMTHit << " with beta Threshold "<< betaThr
94  << "\nUSe of FibreBundle Hit set to "<<useFibreBundle
95  << "\n Use of Birks law is set to "
96  << useBirk << " with three constants kB = "
97  << birk1 << ", C1 = " << birk2 << ", C2 = " << birk3;
98  edm::LogInfo("HcalSim") << "HCalSD:: Suppression Flag " << suppressHeavy
99  << " protons below " << kmaxProton << " MeV,"
100  << " neutrons below " << kmaxNeutron << " MeV and"
101  << " ions below " << kmaxIon << " MeV\n"
102  << " Threshold for storing hits in HB: "
103  << eminHitHB << " HE: " << eminHitHE << " HO: "
104  << eminHitHO << " HF: " << eminHitHF << "\n"
105  << "Delivered luminosity for Darkening "
106  << deliveredLumi << " Flag (HE) " << ageingFlagHE
107  << " Flag (HF) " << ageingFlagHF << "\n"
108  << "Application of Fiducial Cut " << applyFidCut;
109 
110  numberingFromDDD = new HcalNumberingFromDDD(name, cpv);
111  HcalNumberingScheme* scheme;
112  if (testNumber || forTBH2)
113  scheme = dynamic_cast<HcalNumberingScheme*>(new HcalTestNumberingScheme(forTBH2));
114  else
115  scheme = new HcalNumberingScheme();
116  setNumberingScheme(scheme);
117 
118  const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
119  std::vector<G4LogicalVolume *>::const_iterator lvcite;
120  G4LogicalVolume* lv;
121  std::string attribute, value;
122  if (useHF) {
123  if (useParam) {
124  showerParam = new HFShowerParam(name, cpv, p);
125  } else {
126  if (useShowerLibrary) showerLibrary = new HFShowerLibrary(name, cpv, p);
127  hfshower = new HFShower(name, cpv, p, 0);
128  }
129 
130  // HF volume names
131  attribute = "Volume";
132  value = "HF";
133  DDSpecificsFilter filter0;
134  DDValue ddv0(attribute, value, 0);
135  filter0.setCriteria(ddv0, DDSpecificsFilter::equals);
136  DDFilteredView fv0(cpv);
137  fv0.addFilter(filter0);
138  hfNames = getNames(fv0);
139  fv0.firstChild();
140  DDsvalues_type sv0(fv0.mergedSpecifics());
141  std::vector<double> temp = getDDDArray("Levels",sv0);
142  edm::LogInfo("HcalSim") << "HCalSD: Names to be tested for " << attribute
143  << " = " << value << " has " << hfNames.size()
144  << " elements";
145  for (unsigned int i=0; i < hfNames.size(); ++i) {
146  G4String namv = hfNames[i];
147  lv = 0;
148  for(lvcite=lvs->begin(); lvcite!=lvs->end(); lvcite++)
149  if((*lvcite)->GetName()==namv) {
150  lv = (*lvcite);
151  break;
152  }
153  hfLV.push_back(lv);
154  int level = static_cast<int>(temp[i]);
155  hfLevels.push_back(level);
156  edm::LogInfo("HcalSim") << "HCalSD: HF[" << i << "] = " << hfNames[i]
157  << " LV " << hfLV[i] << " at level "
158  << hfLevels[i];
159  }
160 
161  // HF Fibre volume names
162  value = "HFFibre";
163  DDSpecificsFilter filter1;
164  DDValue ddv1(attribute,value,0);
165  filter1.setCriteria(ddv1, DDSpecificsFilter::equals);
166  DDFilteredView fv1(cpv);
167  fv1.addFilter(filter1);
168  fibreNames = getNames(fv1);
169  edm::LogInfo("HcalSim") << "HCalSD: Names to be tested for " << attribute
170  << " = " << value << ":";
171  for (unsigned int i=0; i<fibreNames.size(); ++i) {
172  G4String namv = fibreNames[i];
173  lv = 0;
174  for (lvcite = lvs->begin(); lvcite != lvs->end(); ++lvcite) {
175  if ((*lvcite)->GetName() == namv) {
176  lv = (*lvcite);
177  break;
178  }
179  }
180  fibreLV.push_back(lv);
181  edm::LogInfo("HcalSim") << "HCalSD: (" << i << ") " << fibreNames[i]
182  << " LV " << fibreLV[i];
183  }
184 
185  // HF PMT volume names
186  value = "HFPMT";
187  DDSpecificsFilter filter3;
188  DDValue ddv3(attribute,value,0);
190  DDFilteredView fv3(cpv);
191  fv3.addFilter(filter3);
192  std::vector<G4String> pmtNames = getNames(fv3);
193  edm::LogInfo("HcalSim") << "HCalSD: Names to be tested for " << attribute
194  << " = " << value << " have " << pmtNames.size()
195  << " entries";
196  for (unsigned int i=0; i<pmtNames.size(); ++i) {
197  G4String namv = pmtNames[i];
198  lv = 0;
199  for (lvcite = lvs->begin(); lvcite != lvs->end(); ++lvcite)
200  if ((*lvcite)->GetName() == namv) {
201  lv = (*lvcite);
202  break;
203  }
204  pmtLV.push_back(lv);
205  edm::LogInfo("HcalSim") << "HCalSD: (" << i << ") " << pmtNames[i]
206  << " LV " << pmtLV[i];
207  }
208  if (pmtNames.size() > 0) showerPMT = new HFShowerPMT (name, cpv, p);
209 
210  // HF Fibre bundles
211  value = "HFFibreBundleStraight";
212  DDSpecificsFilter filter4;
213  DDValue ddv4(attribute,value,0);
215  DDFilteredView fv4(cpv);
216  fv4.addFilter(filter4);
217  std::vector<G4String> fibreNames = getNames(fv4);
218  edm::LogInfo("HcalSim") << "HCalSD: Names to be tested for " << attribute
219  << " = " << value << " have " << fibreNames.size()
220  << " entries";
221  for (unsigned int i=0; i<fibreNames.size(); ++i) {
222  G4String namv = fibreNames[i];
223  lv = 0;
224  for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
225  if ((*lvcite)->GetName() == namv) {
226  lv = (*lvcite);
227  break;
228  }
229  fibre1LV.push_back(lv);
230  edm::LogInfo("HcalSim") << "HCalSD: (" << i << ") " << fibreNames[i]
231  << " LV " << fibre1LV[i];
232  }
233 
234  // Geometry parameters for HF
235  value = "HFFibreBundleConical";
236  DDSpecificsFilter filter5;
237  DDValue ddv5(attribute,value,0);
239  DDFilteredView fv5(cpv);
240  fv5.addFilter(filter5);
241  fibreNames = getNames(fv5);
242  edm::LogInfo("HcalSim") << "HCalSD: Names to be tested for " << attribute
243  << " = " << value << " have " << fibreNames.size()
244  << " entries";
245  for (unsigned int i=0; i<fibreNames.size(); ++i) {
246  G4String namv = fibreNames[i];
247  lv = 0;
248  for (lvcite = lvs->begin(); lvcite != lvs->end(); ++lvcite)
249  if ((*lvcite)->GetName() == namv) {
250  lv = (*lvcite);
251  break;
252  }
253  fibre2LV.push_back(lv);
254  edm::LogInfo("HcalSim") << "HCalSD: (" << i << ") " << fibreNames[i]
255  << " LV " << fibre2LV[i];
256  }
257  if (fibre1LV.size() > 0 || fibre2LV.size() > 0)
258  showerBundle = new HFShowerFibreBundle (name, cpv, p);
259 
260  attribute = "ReadOutName";
261  value = name;
262  DDSpecificsFilter filter6;
263  DDValue ddv6(attribute,value,0);
265  DDFilteredView fv6(cpv);
266  fv6.addFilter(filter6);
267  if (fv6.firstChild()) {
269  //Special Geometry parameters
270  gpar = getDDDArray("gparHF",sv);
271  edm::LogInfo("HFShower") << "HFShowerParam: " << gpar.size()
272  << " gpar (cm)";
273  for (unsigned int ig=0; ig<gpar.size(); ig++)
274  edm::LogInfo("HFShower") << "HFShowerParam: gpar[" << ig << "] = "
275  << gpar[ig]/cm << " cm";
276  } else {
277  edm::LogWarning("HFShower") << "HFShowerParam: cannot get filtered "
278  << " view for " << attribute << " matching "
279  << name;
280  }
281  }
282 
283  //Material list for HB/HE/HO sensitive detectors
284  attribute = "ReadOutName";
285  DDSpecificsFilter filter2;
286  DDValue ddv2(attribute,name,0);
288  DDFilteredView fv2(cpv);
289  fv2.addFilter(filter2);
290  bool dodet = fv2.firstChild();
291 
293  //Layer0 Weight
294  layer0wt = getDDDArray("Layer0Wt",sv);
295  edm::LogInfo("HcalSim") << "HCalSD: " << layer0wt.size() << " Layer0Wt";
296  for (unsigned int it=0; it<layer0wt.size(); ++it)
297  edm::LogInfo("HcalSim") << "HCalSD: [" << it << "] " << layer0wt[it];
298 
299  const G4MaterialTable * matTab = G4Material::GetMaterialTable();
300  std::vector<G4Material*>::const_iterator matite;
301  while (dodet) {
302  const DDLogicalPart & log = fv2.logicalPart();
303  G4String namx = log.name().name();
304  if (!isItHF(namx) && !isItFibre(namx)) {
305  bool notIn = true;
306  for (unsigned int i=0; i<matNames.size(); ++i) {
307  if (!strcmp(matNames[i].c_str(),log.material().name().name().c_str())){
308  notIn = false;
309  break;
310  }
311  }
312  if (notIn) {
313  namx = log.material().name().name();
314  matNames.push_back(namx);
315  G4Material* mat = nullptr;
316  for (matite = matTab->begin(); matite != matTab->end(); ++matite) {
317  if ((*matite)->GetName() == namx) {
318  mat = (*matite);
319  break;
320  }
321  }
322  materials.push_back(mat);
323  }
324  }
325  dodet = fv2.next();
326  }
327 
328  edm::LogInfo("HcalSim") << "HCalSD: Material names for " << attribute
329  << " = " << name << ":";
330  for (unsigned int i=0; i<matNames.size(); ++i)
331  edm::LogInfo("HcalSim") << "HCalSD: (" << i << ") " << matNames[i]
332  << " pointer " << materials[i];
333 
334  mumPDG = mupPDG = 0;
335 
336  if (useLayerWt) readWeightFromFile(file);
337 
338  for (int i=0; i<9; ++i) hit_[i] = time_[i]= dist_[i] = 0;
339  hzvem = hzvhad = 0;
340 
341  if (ageingFlagHE) m_HEDarkening = new HEDarkening();
342  if (ageingFlagHF) m_HFDarkening = new HFDarkening();
343 #ifdef plotDebug
345 
346  if ( tfile.isAvailable() ) {
347  static const char * const labels[] = {"HB", "HE", "HO", "HF Absorber", "HF PMT",
348  "HF Absorber Long", "HF Absorber Short",
349  "HF PMT Long", "HF PMT Short"};
350  TFileDirectory hcDir = tfile->mkdir("ProfileFromHCalSD");
351  char name[20], title[60];
352  for (int i=0; i<9; ++i) {
353  sprintf (title, "Hit energy in %s", labels[i]);
354  sprintf (name, "HCalSDHit%d", i);
355  hit_[i] = hcDir.make<TH1F>(name, title, 2000, 0., 2000.);
356  sprintf (title, "Energy (MeV)");
357  hit_[i]->GetXaxis()->SetTitle(title);
358  hit_[i]->GetYaxis()->SetTitle("Hits");
359  sprintf (title, "Time of the hit in %s", labels[i]);
360  sprintf (name, "HCalSDTime%d", i);
361  time_[i] = hcDir.make<TH1F>(name, title, 2000, 0., 2000.);
362  sprintf (title, "Time (ns)");
363  time_[i]->GetXaxis()->SetTitle(title);
364  time_[i]->GetYaxis()->SetTitle("Hits");
365  sprintf (title, "Longitudinal profile in %s", labels[i]);
366  sprintf (name, "HCalSDDist%d", i);
367  dist_[i] = hcDir.make<TH1F>(name, title, 2000, 0., 2000.);
368  sprintf (title, "Distance (mm)");
369  dist_[i]->GetXaxis()->SetTitle(title);
370  dist_[i]->GetYaxis()->SetTitle("Hits");
371  }
372  if (useHF && (!useParam)) {
373  hzvem = hcDir.make<TH1F>("hzvem", "Longitudinal Profile (EM Part)",330,0.0,1650.0);
374  hzvem->GetXaxis()->SetTitle("Longitudinal Profile (EM Part)");
375  hzvhad = hcDir.make<TH1F>("hzvhad","Longitudinal Profile (Had Part)",330,0.0,1650.0);
376  hzvhad->GetXaxis()->SetTitle("Longitudinal Profile (Hadronic Part)");
377  }
378  }
379 #endif
380 }
381 
383 
385  if (numberingScheme) delete numberingScheme;
386  if (showerLibrary) delete showerLibrary;
387  if (hfshower) delete hfshower;
388  if (showerParam) delete showerParam;
389  if (showerPMT) delete showerPMT;
390  if (showerBundle) delete showerBundle;
391  if (m_HEDarkening) delete m_HEDarkening;
392  if (m_HFDarkening) delete m_HFDarkening;
393 }
394 
395 bool HCalSD::ProcessHits(G4Step * aStep, G4TouchableHistory * ) {
396 
397  NaNTrap( aStep ) ;
398 
399  if (aStep == NULL) {
400  return true;
401  } else {
402  G4LogicalVolume* lv =
403  aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
404  G4String nameVolume = lv->GetName();
405  if (isItHF(aStep)) {
406  G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
407  double weight(1.0);
408  if (m_HFDarkening != 0) {
409  G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
410  double r = hitPoint.perp()/CLHEP::cm;
411  double z = std::abs(hitPoint.z())/CLHEP::cm;
412  float dose_acquired = 0.;
413  if (z>=1100 && z <= 1300) {
414  int hfZLayer = (int)((z - 1100)/20);
415  if (hfZLayer > 9) hfZLayer = 9;
416  float normalized_lumi = m_HFDarkening->int_lumi(deliveredLumi);
417  for (int i = hfZLayer; i <= 9; ++i) {
418  dose_acquired = m_HFDarkening->dose(i,r);
419  weight *= m_HFDarkening->degradation(normalized_lumi*dose_acquired);
420  }
421  }
422 #ifdef DebugLog
423  LogDebug("HcalSim") << "HCalSD: HFLumiDarkening at r = " << r
424  << ", z = " << z << " Dose " << dose_acquired
425  << " weight " << weight;
426 #endif
427  }
428  if (useParam) {
429 #ifdef DebugLog
430  LogDebug("HcalSim") << "HCalSD: " << getNumberOfHits()
431  << " hits from parametrization in " << nameVolume
432  << " for Track " << aStep->GetTrack()->GetTrackID()
433  <<" (" << aStep->GetTrack()->GetDefinition()->GetParticleName()
434  <<")";
435 #endif
436  getFromParam(aStep, weight);
437 #ifdef DebugLog
438  LogDebug("HcalSim") << "HCalSD: " << getNumberOfHits()
439  << " hits afterParamS*";
440 #endif
441  } else {
442  bool notaMuon = true;
443  if (parCode == mupPDG || parCode == mumPDG ) notaMuon = false;
444  if (useShowerLibrary && notaMuon) {
445 #ifdef DebugLog
446  LogDebug("HcalSim") << "HCalSD: Starts shower library from "
447  << nameVolume << " for Track "
448  << aStep->GetTrack()->GetTrackID() << " ("
449  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
450 #endif
451  getFromLibrary(aStep, weight);
452  } else if (isItFibre(lv)) {
453 #ifdef DebugLog
454  LogDebug("HcalSim") << "HCalSD: Hit at Fibre in " << nameVolume
455  << " for Track "
456  << aStep->GetTrack()->GetTrackID() <<" ("
457  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
458 #endif
459  hitForFibre(aStep, weight);
460  }
461  }
462  } else if (isItPMT(lv)) {
463 #ifdef DebugLog
464  LogDebug("HcalSim") << "HCalSD: Hit from PMT parametrization from "
465  << nameVolume << " for Track "
466  << aStep->GetTrack()->GetTrackID() << " ("
467  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
468 #endif
469  if (usePMTHit && showerPMT) getHitPMT(aStep);
470  } else if (isItStraightBundle(lv) || isItConicalBundle(lv)) {
471 #ifdef DebugLog
472  LogDebug("HcalSim") << "HCalSD: Hit from FibreBundle from "
473  << nameVolume << " for Track "
474  << aStep->GetTrack()->GetTrackID() << " ("
475  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
476 #endif
479  } else {
480 #ifdef DebugLog
481  LogDebug("HcalSim") << "HCalSD: Hit from standard path from "
482  << nameVolume << " for Track "
483  << aStep->GetTrack()->GetTrackID() << " ("
484  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
485 #endif
486  if (getStepInfo(aStep)) {
487 #ifdef plotDebug
488  if (edepositEM+edepositHAD > 0)
489  plotProfile(aStep, aStep->GetPreStepPoint()->GetPosition(),
490  edepositEM+edepositHAD,aStep->GetPostStepPoint()->GetGlobalTime(),0);
491 #endif
492  if (hitExists() == false && edepositEM+edepositHAD>0.) currentHit = createNewHit();
493  }
494  }
495  return true;
496  }
497 }
498 
499 double HCalSD::getEnergyDeposit(G4Step* aStep) {
500  double destep = aStep->GetTotalEnergyDeposit();
501  double weight = 1;
502  G4Track* theTrack = aStep->GetTrack();
503 
504  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
505  int depth = (touch->GetReplicaNumber(0))%10;
506  int det = ((touch->GetReplicaNumber(1))/1000)-3;
507  if (depth==0 && (det==0 || det==1)) weight = layer0wt[det];
508  if (useLayerWt) {
509  int lay = (touch->GetReplicaNumber(0)/10)%100 + 1;
510  G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
511  weight = layerWeight(det+3, hitPoint, depth, lay);
512  }
513 
514  if (m_HEDarkening !=0 && det == 1) {
515  int lay = (touch->GetReplicaNumber(0)/10)%100 + 1;
516  int ieta;
517  uint32_t detid = setDetUnitId(aStep);
518  if(testNumber) {
519  int det, z, depth, eta, phi, lay;
520  HcalTestNumbering::unpackHcalIndex(detid,det,z,depth,eta,phi,lay);
521  ieta = eta;
522  }
523  else ieta = HcalDetId(detid).ietaAbs();
524 #ifdef DebugLog
525  edm::LogInfo("HcalSimDark") << "HCalSD:HE_Darkening >>> ieta: "<< ieta //<< " vs. ietaAbs(): " << HcalDetId(detid).ietaAbs()
526  << " lay: " << lay-2;
527 #endif
528  float dweight = m_HEDarkening->degradation(deliveredLumi,ieta,lay-2);//NB:diff. layer count
529  weight *= dweight;
530 #ifdef DebugLog
531  edm::LogInfo("HcalSimDark") << "HCalSD: >>> Lumi: " << deliveredLumi
532  << " coefficient = " << dweight;
533 #endif
534  }
535 
536  if (suppressHeavy) {
537  TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
538  if (trkInfo) {
539  int pdg = theTrack->GetDefinition()->GetPDGEncoding();
540  if (!(trkInfo->isPrimary())) { // Only secondary particles
541  double ke = theTrack->GetKineticEnergy()/MeV;
542  if ( pdg/1000000000 == 1 && (pdg/10000)%100 > 0 &&
543  (pdg/10)%100 > 0 && ke <kmaxIon ) weight = 0;
544  if ((pdg == 2212) && (ke < kmaxProton)) weight = 0;
545  if ((pdg == 2112) && (ke < kmaxNeutron)) weight = 0;
546 #ifdef DebugLog
547  if (weight == 0)
548  edm::LogInfo("HcalSim") << "HCalSD:Ignore Track "
549  << theTrack->GetTrackID() << " Type "
550  << theTrack->GetDefinition()->GetParticleName()
551  << " Kinetic Energy " << ke << " MeV";
552 #endif
553  }
554  }
555  }
556 #ifdef DebugLog
557  double weight0 = weight;
558 #endif
559  if (useBirk) {
560  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
561  if (isItScintillator(mat)) weight *= getAttenuation(aStep, birk1, birk2, birk3);
562  }
563  double wt1 = getResponseWt(theTrack);
564  double wt2 = theTrack->GetWeight();
565  /*
566  if(wt2 != 1.0) {
567  std::cout << "HCalSD: Detector " << det+3 << " Depth " << depth
568  << " weight= " << weight << " wt1= " << wt1
569  << " wt2= " << wt2 << std::endl;
570  std::cout << theTrack->GetDefinition()->GetParticleName()
571  << " " << theTrack->GetKineticEnergy()
572  << " Id=" << theTrack->GetTrackID()
573  << " IdP=" << theTrack->GetParentID();
574  const G4VProcess* pr = theTrack->GetCreatorProcess();
575  if(pr) std::cout << " from " << pr->GetProcessName();
576  std::cout << std::endl;
577  }
578  */
579 #ifdef DebugLog
580  edm::LogInfo("HcalSim") << "HCalSD: Detector " << det+3 << " Depth " << depth
581  << " weight " << weight0 << " " << weight << " " << wt1
582  << " " << wt2;
583 #endif
584  double edep = weight*wt1*destep;
585  if(wt2 > 0.0) { edep *= wt2; }
586  return edep;
587 }
588 
589 uint32_t HCalSD::setDetUnitId(G4Step * aStep) {
590 
591  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
592  const G4VTouchable* touch = preStepPoint->GetTouchable();
593  G4ThreeVector hitPoint = preStepPoint->GetPosition();
594 
595  int depth = (touch->GetReplicaNumber(0))%10 + 1;
596  int lay = (touch->GetReplicaNumber(0)/10)%100 + 1;
597  int det = (touch->GetReplicaNumber(1))/1000;
598 
599  return setDetUnitId (det, hitPoint, depth, lay);
600 }
601 
603  if (scheme != 0) {
604  edm::LogInfo("HcalSim") << "HCalSD: updates numbering scheme for " << GetName();
605  if (numberingScheme) delete numberingScheme;
606  numberingScheme = scheme;
607  }
608 }
609 
611  G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
612  G4String particleName;
613  mumPDG = theParticleTable->FindParticle(particleName="mu-")->GetPDGEncoding();
614  mupPDG = theParticleTable->FindParticle(particleName="mu+")->GetPDGEncoding();
615 #ifdef DebugLog
616  LogDebug("HcalSim") << "HCalSD: Particle code for mu- = " << mumPDG
617  << " for mu+ = " << mupPDG;
618 #endif
619  if (showerLibrary) showerLibrary->initRun(theParticleTable);
620  if (showerParam) showerParam->initRun(theParticleTable);
621  if (hfshower) hfshower->initRun(theParticleTable);
622 }
623 
624 bool HCalSD::filterHit(CaloG4Hit* aHit, double time) {
625  double threshold=0;
626  DetId theId(aHit->getUnitID());
627  switch (theId.subdetId()) {
628  case HcalBarrel:
629  threshold = eminHitHB; break;
630  case HcalEndcap:
631  threshold = eminHitHE; break;
632  case HcalOuter:
633  threshold = eminHitHO; break;
634  case HcalForward:
635  threshold = eminHitHF; break;
636  default:
637  break;
638  }
639  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > threshold));
640 }
641 
642 
643 uint32_t HCalSD::setDetUnitId (int det, const G4ThreeVector& pos, int depth, int lay=1) {
644  uint32_t id = 0;
645  if (numberingFromDDD) {
646  //get the ID's as eta, phi, depth, ... indices
647  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, pos, depth, lay);
648  //get the ID
649  if (numberingScheme) id = numberingScheme->getUnitID(tmp);
650  }
651  return id;
652 }
653 
654 std::vector<double> HCalSD::getDDDArray(const std::string & str,
655  const DDsvalues_type & sv) {
656 #ifdef DebugLog
657  LogDebug("HcalSim") << "HCalSD:getDDDArray called for " << str;
658 #endif
659  DDValue value(str);
660  if (DDfetch(&sv,value)) {
661 #ifdef DebugLog
662  LogDebug("HcalSim") << value;
663 #endif
664  const std::vector<double> & fvec = value.doubles();
665  int nval = fvec.size();
666  if (nval < 1) {
667  edm::LogError("HcalSim") << "HCalSD : # of " << str << " bins " << nval
668  << " < 2 ==> illegal";
669  throw cms::Exception("Unknown", "HCalSD") << "nval < 2 for array " << str << "\n";
670  }
671 
672  return fvec;
673  } else {
674  edm::LogError("HcalSim") << "HCalSD : cannot get array " << str;
675  throw cms::Exception("Unknown", "HCalSD") << "cannot get array " << str << "\n";
676  }
677 }
678 
679 std::vector<G4String> HCalSD::getNames(DDFilteredView& fv) {
680 
681  std::vector<G4String> tmp;
682  bool dodet = fv.firstChild();
683  while (dodet) {
684  const DDLogicalPart & log = fv.logicalPart();
685  bool ok = true;
686 
687  for (unsigned int i=0; i<tmp.size(); ++i) {
688  if (!strcmp(tmp[i].c_str(), log.name().name().c_str())) {
689  ok = false;
690  break;
691  }
692  }
693  if (ok) tmp.push_back(log.name().name());
694  dodet = fv.next();
695  }
696  return tmp;
697 }
698 
699 bool HCalSD::isItHF(G4Step * aStep) {
700  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
701  int levels = (touch->GetHistoryDepth()) + 1;
702  for (unsigned int it=0; it < hfNames.size(); ++it) {
703  if (levels >= hfLevels[it]) {
704  G4LogicalVolume* lv = touch->GetVolume(levels-hfLevels[it])->GetLogicalVolume();
705  if (lv == hfLV[it]) return true;
706  }
707  }
708  return false;
709 }
710 
711 bool HCalSD::isItHF (G4String name) {
712  std::vector<G4String>::const_iterator it = hfNames.begin();
713  for (; it != hfNames.end(); ++it) if (name == *it) return true;
714  return false;
715 }
716 
717 bool HCalSD::isItFibre (G4LogicalVolume* lv) {
718  std::vector<G4LogicalVolume*>::const_iterator ite = fibreLV.begin();
719  for (; ite != fibreLV.end(); ++ite) if (lv == *ite) return true;
720  return false;
721 }
722 
723 bool HCalSD::isItFibre (G4String name) {
724  std::vector<G4String>::const_iterator it = fibreNames.begin();
725  for (; it != fibreNames.end(); ++it) if (name == *it) return true;
726  return false;
727 }
728 
729 bool HCalSD::isItPMT (G4LogicalVolume* lv) {
730  std::vector<G4LogicalVolume*>::const_iterator ite = pmtLV.begin();
731  for (; ite != pmtLV.end(); ++ite) if (lv == *ite) return true;
732  return false;
733 }
734 
735 bool HCalSD::isItStraightBundle (G4LogicalVolume* lv) {
736  std::vector<G4LogicalVolume*>::const_iterator ite = fibre1LV.begin();
737  for (; ite != fibre1LV.end(); ++ite) if (lv == *ite) return true;
738  return false;
739 }
740 
741 bool HCalSD::isItConicalBundle (G4LogicalVolume* lv) {
742  std::vector<G4LogicalVolume*>::const_iterator ite = fibre2LV.begin();
743  for (; ite != fibre2LV.end(); ++ite) if (lv == *ite) return true;
744  return false;
745 }
746 
747 bool HCalSD::isItScintillator (G4Material* mat) {
748  std::vector<G4Material*>::const_iterator ite = materials.begin();
749  for (; ite != materials.end(); ++ite) if (mat == *ite) return true;
750  return false;
751 }
752 
753 bool HCalSD::isItinFidVolume (G4ThreeVector& hitPoint) {
754  bool flag = true;
755  if (applyFidCut) {
756  int npmt = HFFibreFiducial:: PMTNumber(hitPoint);
757 #ifdef DebugLog
758  edm::LogInfo("HcalSim") << "HCalSD::isItinFidVolume:#PMT= " << npmt
759  << " for hit point " << hitPoint;
760 #endif
761  if (npmt <= 0) flag = false;
762  }
763 #ifdef DebugLog
764  edm::LogInfo("HcalSim") << "HCalSD::isItinFidVolume: point " << hitPoint
765  << " return flag " << flag;
766 #endif
767  return flag;
768 }
769 
770 void HCalSD::getFromLibrary (G4Step* aStep, double weight) {
771  preStepPoint = aStep->GetPreStepPoint();
772  theTrack = aStep->GetTrack();
773  int det = 5;
774  bool ok;
775 
776  std::vector<HFShowerLibrary::Hit> hits = showerLibrary->getHits(aStep, ok, weight, false);
777 
778  double etrack = preStepPoint->GetKineticEnergy();
779  int primaryID = setTrackID(aStep);
780 
781  // Reset entry point for new primary
782  posGlobal = preStepPoint->GetPosition();
783  resetForNewPrimary(posGlobal, etrack);
784 
785  G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding();
786  if (particleCode==emPDG || particleCode==epPDG || particleCode==gammaPDG) {
787  edepositEM = 1.*GeV;
788  edepositHAD = 0.;
789  } else {
790  edepositEM = 0.;
791  edepositHAD = 1.*GeV;
792  }
793 #ifdef DebugLog
794  edm::LogInfo("HcalSim") << "HCalSD::getFromLibrary " <<hits.size()
795  << " hits for " << GetName() << " of " << primaryID
796  << " with " << theTrack->GetDefinition()->GetParticleName()
797  << " of " << preStepPoint->GetKineticEnergy()/GeV << " GeV";
798 #endif
799  for (unsigned int i=0; i<hits.size(); ++i) {
800  G4ThreeVector hitPoint = hits[i].position;
801  if (isItinFidVolume (hitPoint)) {
802  int depth = hits[i].depth;
803  double time = hits[i].time;
804  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
805  currentID.setID(unitID, time, primaryID, 0);
806 #ifdef plotDebug
807  plotProfile(aStep, hitPoint, 1.0*GeV, time, depth);
808  bool emType = false;
809  if (particleCode==emPDG || particleCode==epPDG || particleCode==gammaPDG)
810  emType = true;
811  plotHF(hitPoint,emType);
812 #endif
813 
814  // check if it is in the same unit and timeslice as the previous one
815  if (currentID == previousID) {
817  } else {
818  if (!checkHit()) currentHit = createNewHit();
819  }
820  }
821  }
822 
823  //Now kill the current track
824  if (ok) {
825  theTrack->SetTrackStatus(fStopAndKill);
826  G4TrackVector tv = *(aStep->GetSecondary());
827  for (unsigned int kk=0; kk<tv.size(); ++kk)
828  if (tv[kk]->GetVolume() == preStepPoint->GetPhysicalVolume())
829  tv[kk]->SetTrackStatus(fStopAndKill);
830  }
831 }
832 
833 void HCalSD::hitForFibre (G4Step* aStep, double weight) { // if not ParamShower
834 
835  preStepPoint = aStep->GetPreStepPoint();
836  theTrack = aStep->GetTrack();
837  int primaryID = setTrackID(aStep);
838 
839  int det = 5;
840  std::vector<HFShower::Hit> hits = hfshower->getHits(aStep, weight);
841 
842  G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding();
843  if (particleCode==emPDG || particleCode==epPDG || particleCode==gammaPDG) {
844  edepositEM = 1.*GeV;
845  edepositHAD = 0.;
846  } else {
847  edepositEM = 0.;
848  edepositHAD = 1.*GeV;
849  }
850 
851 #ifdef DebugLog
852  edm::LogInfo("HcalSim") << "HCalSD::hitForFibre " << hits.size()
853  << " hits for " << GetName() << " of " << primaryID
854  << " with " << theTrack->GetDefinition()->GetParticleName()
855  << " of " << preStepPoint->GetKineticEnergy()/GeV
856  << " GeV in detector type " << det;
857 #endif
858  if (hits.size() > 0) {
859  for (unsigned int i=0; i<hits.size(); ++i) {
860  G4ThreeVector hitPoint = hits[i].position;
861  if (isItinFidVolume (hitPoint)) {
862  int depth = hits[i].depth;
863  double time = hits[i].time;
864  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
865  currentID.setID(unitID, time, primaryID, 0);
866 #ifdef plotDebug
867  plotProfile(aStep, hitPoint, edepositEM, time, depth);
868  bool emType = false;
869  if (particleCode==emPDG || particleCode==epPDG || particleCode==gammaPDG)
870  emType = true;
871  plotHF(hitPoint,emType);
872 #endif
873  // check if it is in the same unit and timeslice as the previous one
874  if (currentID == previousID) {
876  } else {
877  posGlobal = preStepPoint->GetPosition();
878  if (!checkHit()) currentHit = createNewHit();
879  }
880  }
881  }
882  }
883 }
884 
885 void HCalSD::getFromParam (G4Step* aStep, double weight) {
886  std::vector<HFShowerParam::Hit> hits = showerParam->getHits(aStep, weight);
887  int nHit = static_cast<int>(hits.size());
888 
889  if (nHit > 0) {
890  preStepPoint = aStep->GetPreStepPoint();
891  int primaryID = setTrackID(aStep);
892 
893  int det = 5;
894 #ifdef DebugLog
895  edm::LogInfo("HcalSim") << "HCalSD::getFromParam " << nHit << " hits for "
896  << GetName() << " of " << primaryID << " with "
897  << aStep->GetTrack()->GetDefinition()->GetParticleName()
898  << " of " << preStepPoint->GetKineticEnergy()/GeV
899  << " GeV in detector type " << det;
900 #endif
901  for (int i = 0; i<nHit; ++i) {
902  G4ThreeVector hitPoint = hits[i].position;
903  int depth = hits[i].depth;
904  double time = hits[i].time;
905  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
906  currentID.setID(unitID, time, primaryID, 0);
907  edepositEM = hits[i].edep*GeV;
908  edepositHAD = 0.;
909 #ifdef plotDebug
910  plotProfile(aStep, hitPoint, edepositEM, time, depth);
911 #endif
912 
913  // check if it is in the same unit and timeslice as the previous one
914  if (currentID == previousID) {
916  } else {
917  posGlobal = preStepPoint->GetPosition();
918  if (!checkHit()) currentHit = createNewHit();
919  }
920  }
921  }
922 }
923 
924 void HCalSD::getHitPMT (G4Step * aStep) {
925 
926  preStepPoint = aStep->GetPreStepPoint();
927  theTrack = aStep->GetTrack();
928  double edep = showerPMT->getHits(aStep);
929 
930  if (edep >= 0) {
931  double etrack = preStepPoint->GetKineticEnergy();
932  int primaryID = 0;
933  if (etrack >= energyCut) {
934  primaryID = theTrack->GetTrackID();
935  } else {
936  primaryID = theTrack->GetParentID();
937  if (primaryID == 0) primaryID = theTrack->GetTrackID();
938  }
939  // Reset entry point for new primary
940  posGlobal = preStepPoint->GetPosition();
941  resetForNewPrimary(posGlobal, etrack);
942 
943  //
944  int det = static_cast<int>(HcalForward);
945  G4ThreeVector hitPoint = preStepPoint->GetPosition();
946  double rr = (hitPoint.x()*hitPoint.x() + hitPoint.y()*hitPoint.y());
947  double phi = (rr == 0. ? 0. :atan2(hitPoint.y(),hitPoint.x()));
948  double etaR = showerPMT->getRadius();
949  int depth = 1;
950  if (etaR < 0) {
951  depth = 2;
952  etaR =-etaR;
953  }
954  if (hitPoint.z() < 0) etaR =-etaR;
955 #ifdef DebugLog
956  edm::LogInfo("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR "
957  << etaR << " phi " << phi/deg << " depth " <<depth;
958 #endif
959  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
960  uint32_t unitID = 0;
961  if (numberingFromDDD) {
963  depth,1);
964  if (numberingScheme) unitID = numberingScheme->getUnitID(tmp);
965  }
966  currentID.setID(unitID, time, primaryID, 1);
967 
968  edepositHAD = aStep->GetTotalEnergyDeposit();
969  edepositEM =-edepositHAD + (edep*GeV);
970 #ifdef plotDebug
971  plotProfile(aStep, hitPoint, edep*GeV, time, depth);
972 #endif
973 #ifdef DebugLog
974  double beta = preStepPoint->GetBeta();
975  LogDebug("HcalSim") << "HCalSD::getHitPMT 1 hit for " << GetName()
976  << " of " << primaryID << " with "
977  << theTrack->GetDefinition()->GetParticleName()
978  << " of " << preStepPoint->GetKineticEnergy()/GeV
979  << " GeV with velocity " << beta << " UnitID "
980  << std::hex << unitID << std::dec;
981 #endif
982  // check if it is in the same unit and timeslice as the previous one
983  if (currentID == previousID) {
985  } else {
986  if (!checkHit()) currentHit = createNewHit();
987  }
988  }
989 }
990 
991 void HCalSD::getHitFibreBundle (G4Step* aStep, bool type) {
992  preStepPoint = aStep->GetPreStepPoint();
993  theTrack = aStep->GetTrack();
994  double edep = showerBundle->getHits(aStep, type);
995 
996  if (edep >= 0) {
997  double etrack = preStepPoint->GetKineticEnergy();
998  int primaryID = 0;
999  if (etrack >= energyCut) {
1000  primaryID = theTrack->GetTrackID();
1001  } else {
1002  primaryID = theTrack->GetParentID();
1003  if (primaryID == 0) primaryID = theTrack->GetTrackID();
1004  }
1005  // Reset entry point for new primary
1006  posGlobal = preStepPoint->GetPosition();
1007  resetForNewPrimary(posGlobal, etrack);
1008 
1009  //
1010  int det = static_cast<int>(HcalForward);
1011  G4ThreeVector hitPoint = preStepPoint->GetPosition();
1012  double rr = (hitPoint.x()*hitPoint.x() + hitPoint.y()*hitPoint.y());
1013  double phi = (rr == 0. ? 0. :atan2(hitPoint.y(),hitPoint.x()));
1014  double etaR = showerBundle->getRadius();
1015  int depth = 1;
1016  if (etaR < 0) {
1017  depth = 2;
1018  etaR =-etaR;
1019  }
1020  if (hitPoint.z() < 0) etaR =-etaR;
1021 #ifdef DebugLog
1022  LogDebug("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR "
1023  << etaR << " phi " << phi/deg << " depth " <<depth;
1024 #endif
1025  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
1026  uint32_t unitID = 0;
1027  if (numberingFromDDD) {
1028  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det,etaR,phi,depth,1);
1029  if (numberingScheme) unitID = numberingScheme->getUnitID(tmp);
1030  }
1031  if (type) currentID.setID(unitID, time, primaryID, 3);
1032  else currentID.setID(unitID, time, primaryID, 2);
1033 
1034  edepositHAD = aStep->GetTotalEnergyDeposit();
1035  edepositEM =-edepositHAD + (edep*GeV);
1036 #ifdef plotDebug
1037  plotProfile(aStep, hitPoint, edep*GeV, time, depth);
1038 #endif
1039 #ifdef DebugLog
1040  double beta = preStepPoint->GetBeta();
1041  LogDebug("HcalSim") << "HCalSD::getHitFibreBundle 1 hit for " << GetName()
1042  << " of " << primaryID << " with "
1043  << theTrack->GetDefinition()->GetParticleName()
1044  << " of " << preStepPoint->GetKineticEnergy()/GeV
1045  << " GeV with velocity " << beta << " UnitID "
1046  << std::hex << unitID << std::dec;
1047 #endif
1048  // check if it is in the same unit and timeslice as the previous one
1050  else if (!checkHit()) currentHit = createNewHit();
1051  } // non-zero energy deposit
1052 }
1053 
1054 int HCalSD::setTrackID (G4Step* aStep) {
1055  theTrack = aStep->GetTrack();
1056 
1057  double etrack = preStepPoint->GetKineticEnergy();
1058  TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
1059  int primaryID = trkInfo->getIDonCaloSurface();
1060  if (primaryID == 0) {
1061 #ifdef DebugLog
1062  edm::LogInfo("HcalSim") << "HCalSD: Problem with primaryID **** set by "
1063  << "force to TkID **** " <<theTrack->GetTrackID();
1064 #endif
1065  primaryID = theTrack->GetTrackID();
1066  }
1067 
1068  if (primaryID != previousID.trackID())
1069  resetForNewPrimary(preStepPoint->GetPosition(), etrack);
1070 
1071  return primaryID;
1072 }
1073 
1075 
1076  std::ifstream infile;
1077  int entry=0;
1078  infile.open(fName.c_str(), std::ios::in);
1079  if (infile) {
1080  int det, zside, etaR, phi, lay;
1081  double wt;
1082  while (infile >> det >> zside >> etaR >> phi >> lay >> wt) {
1083  uint32_t id = HcalTestNumbering::packHcalIndex(det,zside,1,etaR,phi,lay);
1084  layerWeights.insert(std::pair<uint32_t,double>(id,wt));
1085  ++entry;
1086 #ifdef DebugLog
1087  edm::LogInfo("HcalSim") << "HCalSD::readWeightFromFile:Entry " << entry
1088  << " ID " << std::hex << id << std::dec << " ("
1089  << det << "/" << zside << "/1/" << etaR << "/"
1090  << phi << "/" << lay << ") Weight " << wt;
1091 #endif
1092  }
1093  infile.close();
1094  }
1095  edm::LogInfo("HcalSim") << "HCalSD::readWeightFromFile: reads " << entry
1096  << " weights from " << fName;
1097  if (entry <= 0) useLayerWt = false;
1098 }
1099 
1100 double HCalSD::layerWeight(int det, const G4ThreeVector& pos, int depth, int lay) {
1101 
1102  double wt = 1.;
1103  if (numberingFromDDD) {
1104  //get the ID's as eta, phi, depth, ... indices
1106  depth, lay);
1107  uint32_t id = HcalTestNumbering::packHcalIndex(tmp.subdet, tmp.zside, 1,
1108  tmp.etaR, tmp.phis,tmp.lay);
1109  std::map<uint32_t,double>::const_iterator ite = layerWeights.find(id);
1110  if (ite != layerWeights.end()) wt = ite->second;
1111 #ifdef DebugLog
1112  edm::LogInfo("HcalSim") << "HCalSD::layerWeight: ID " << std::hex << id
1113  << std::dec << " (" << tmp.subdet << "/"
1114  << tmp.zside << "/1/" << tmp.etaR << "/"
1115  << tmp.phis << "/" << tmp.lay << ") Weight " <<wt;
1116 #endif
1117  }
1118  return wt;
1119 }
1120 
1121 void HCalSD::plotProfile(G4Step* aStep,const G4ThreeVector& global, double edep,
1122  double time, int id) {
1123 
1124  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
1125  static G4String modName[8] = {"HEModule", "HVQF" , "HBModule", "MBAT",
1126  "MBBT" , "MBBTC", "MBBT_R1P", "MBBT_R1M"};
1127  G4ThreeVector local;
1128  bool found=false;
1129  double depth=-2000;
1130  int idx = 4;
1131  for (int n=0; n<touch->GetHistoryDepth(); ++n) {
1132  G4String name = touch->GetVolume(n)->GetName();
1133 #ifdef DebugLog
1134  LogDebug("HcalSim") << "plotProfile Depth " << n << " Name " << name;
1135 #endif
1136  for (unsigned int ii=0; ii<8; ++ii) {
1137  if (name == modName[ii]) {
1138  found = true;
1139  int dn = touch->GetHistoryDepth() - n;
1140  local = touch->GetHistory()->GetTransform(dn).TransformPoint(global);
1141  if (ii == 0) {depth = local.z() - 4006.5; idx = 1;}
1142  else if (ii == 1) {depth = local.z() + 825.0; idx = 3;}
1143  else if (ii == 2) {depth = local.x() - 1775.; idx = 0;}
1144  else {depth = local.y() + 15.; idx = 2;}
1145  break;
1146  }
1147  }
1148  if (found) break;
1149  }
1150  if (!found) depth = std::abs(global.z()) - 11500;
1151 #ifdef DebugLog
1152  LogDebug("HcalSim") << "plotProfile Found " << found << " Global " << global
1153  << " Local " << local << " depth " << depth << " ID "
1154  << id << " EDEP " << edep << " Time " << time;
1155 #endif
1156  if (hit_[idx] != 0) hit_[idx]->Fill(edep);
1157  if (time_[idx] != 0) time_[idx]->Fill(time,edep);
1158  if (dist_[idx] != 0) dist_[idx]->Fill(depth,edep);
1159  int jd = 2*idx + id - 7;
1160  if (jd >= 0 && jd < 4) {
1161  jd += 5;
1162  if (hit_[jd] != 0) hit_[jd]->Fill(edep);
1163  if (time_[jd] != 0) time_[jd]->Fill(time,edep);
1164  if (dist_[jd] != 0) dist_[jd]->Fill(depth,edep);
1165  }
1166 }
1167 
1168 void HCalSD::plotHF(G4ThreeVector& hitPoint, bool emType) {
1169  double zv = std::abs(hitPoint.z()) - gpar[4];
1170  if (emType) {
1171  if (hzvem != 0) hzvem->Fill(zv);
1172  } else {
1173  if (hzvhad != 0) hzvhad->Fill(zv);
1174  }
1175 }
#define LogDebug(id)
float edepositEM
Definition: CaloSD.h:121
const double beta
double energyCut
Definition: CaloSD.h:123
HCalSD(G4String, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HCalSD.cc:38
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
std::vector< double > layer0wt
Definition: HCalSD.h:92
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
bool useParam
Definition: HCalSD.h:88
double eminHitHE
Definition: HCalSD.h:89
HFShowerParam * showerParam
Definition: HCalSD.h:81
TH1F * time_[9]
Definition: HCalSD.h:98
G4int emPDG
Definition: CaloSD.h:136
double getRadius()
Definition: HFShowerPMT.cc:150
float int_lumi(float intlumi)
Definition: HFDarkening.cc:147
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
const std::vector< double > & doubles() const
a reference to the double-valued values stored in the given instance of DDValue
Definition: DDValue.cc:139
const double GeV
Definition: MathUtil.h:16
double kmaxNeutron
Definition: CaloSD.h:134
void addFilter(const DDFilter &, log_op op=AND)
void updateHit(CaloG4Hit *)
Definition: CaloSD.cc:440
void hitForFibre(G4Step *step, double weight)
Definition: HCalSD.cc:833
const N & name() const
Definition: DDBase.h:82
void getFromParam(G4Step *step, double weight)
Definition: HCalSD.cc:885
int getIDonCaloSurface() const
bool useLayerWt
Definition: HCalSD.h:86
Definition: CaloSD.h:42
void initRun(G4ParticleTable *)
Definition: HFShower.cc:478
bool isItStraightBundle(G4LogicalVolume *)
Definition: HCalSD.cc:735
virtual bool ProcessHits(G4Step *, G4TouchableHistory *)
Definition: HCalSD.cc:395
std::vector< double > gpar
Definition: HCalSD.h:92
bool useFibreBundle
Definition: HCalSD.h:86
bool isItinFidVolume(G4ThreeVector &)
Definition: HCalSD.cc:753
double betaThr
Definition: HCalSD.h:87
double deliveredLumi
Definition: HCalSD.h:90
std::vector< Hit > getHits(G4Step *aStep, bool &ok, double weight, bool onlyLong=false)
virtual uint32_t getUnitID(const HcalNumberingFromDDD::HcalID &id)
HFShowerFibreBundle * showerBundle
Definition: HCalSD.h:83
std::vector< G4LogicalVolume * > fibre1LV
Definition: HCalSD.h:96
std::vector< G4LogicalVolume * > fibreLV
Definition: HCalSD.h:96
double eminHitHB
Definition: HCalSD.h:89
bool useShowerLibrary
Definition: HCalSD.h:88
std::vector< double > getDDDArray(const std::string &, const DDsvalues_type &)
Definition: HCalSD.cc:654
double getHits(G4Step *aStep, bool type)
#define NULL
Definition: scimark2.h:8
virtual uint32_t setDetUnitId(G4Step *step)
Definition: HCalSD.cc:589
double birk2
Definition: HCalSD.h:87
int setTrackID(G4Step *step)
Definition: HCalSD.cc:1054
void getFromLibrary(G4Step *step, double weight)
Definition: HCalSD.cc:770
bool usePMTHit
Definition: HCalSD.h:86
G4ThreeVector posGlobal
Definition: CaloSD.h:113
int zside(DetId const &)
T eta() const
int ii
Definition: cuy.py:588
void setNumberingScheme(HcalNumberingScheme *)
Definition: HCalSD.cc:602
type of data representation of DDCompactView
Definition: DDCompactView.h:77
double birk1
Definition: HCalSD.h:87
double kmaxProton
Definition: CaloSD.h:134
bool DDfetch(const DDsvalues_type *, DDValue &)
helper for retrieving DDValues from DDsvalues_type *.
Definition: DDsvalues.cc:80
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
float float float z
static int PMTNumber(const G4ThreeVector &pe_effect)
bool isItHF(G4Step *)
Definition: HCalSD.cc:699
G4bool checkHit()
Definition: CaloSD.cc:331
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
TH1F * hzvhad
Definition: HCalSD.h:98
double eminHitHF
Definition: HCalSD.h:89
const double MeV
bool isItConicalBundle(G4LogicalVolume *)
Definition: HCalSD.cc:741
static uint32_t packHcalIndex(int det, int z, int depth, int eta, int phi, int lay)
TH1F * hit_[9]
Definition: HCalSD.h:98
void resetForNewPrimary(const G4ThreeVector &, double)
Definition: CaloSD.cc:454
virtual void initRun()
Definition: HCalSD.cc:610
bool testNumber
Definition: HCalSD.h:86
std::vector< G4LogicalVolume * > pmtLV
Definition: HCalSD.h:96
std::vector< G4LogicalVolume * > hfLV
Definition: HCalSD.h:96
double kmaxIon
Definition: CaloSD.h:134
bool suppressHeavy
Definition: CaloSD.h:132
bool isItScintillator(G4Material *)
Definition: HCalSD.cc:747
bool next()
set current node to the next node in the filtered tree
HFShower * hfshower
Definition: HCalSD.h:80
bool useBirk
Definition: HCalSD.h:86
float edepositHAD
Definition: CaloSD.h:121
double birk3
Definition: HCalSD.h:87
int trackID() const
Definition: CaloHitID.h:25
virtual double getEnergyDeposit(G4Step *)
Definition: HCalSD.cc:499
G4int epPDG
Definition: CaloSD.h:136
void plotHF(G4ThreeVector &pos, bool emType)
Definition: HCalSD.cc:1168
std::vector< std::pair< unsigned int, DDValue > > DDsvalues_type
std::maps an index to a DDValue. The index corresponds to the index assigned to the name of the std::...
Definition: DDsvalues.h:19
float degradation(float mrad)
Definition: HFDarkening.cc:143
G4int gammaPDG
Definition: CaloSD.h:136
bool isAvailable() const
Definition: Service.h:46
bool applyFidCut
Definition: HCalSD.h:88
std::vector< Hit > getHits(G4Step *aStep, double weight)
Definition: HFShower.cc:71
TH1F * dist_[9]
Definition: HCalSD.h:98
virtual ~HCalSD()
Definition: HCalSD.cc:382
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static void unpackHcalIndex(const uint32_t &idx, int &det, int &z, int &depth, int &eta, int &phi, int &lay)
A DDLogicalPart aggregates information concerning material, solid and sensitveness ...
Definition: DDLogicalPart.h:88
HcalNumberingScheme * numberingScheme
Definition: HCalSD.h:78
void getHitPMT(G4Step *step)
Definition: HCalSD.cc:924
T * make(const Args &...args) const
make new ROOT object
int getNumberOfHits()
Definition: CaloSD.cc:361
HFDarkening * m_HFDarkening
Definition: HCalSD.h:85
CaloHitID previousID
Definition: CaloSD.h:117
double getAttenuation(G4Step *aStep, double birk1, double birk2, double birk3)
Definition: CaloSD.cc:465
CaloG4Hit * currentHit
Definition: CaloSD.h:128
std::vector< int > hfLevels
Definition: HCalSD.h:93
virtual std::vector< std::string > getNames()
void plotProfile(G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:1121
HEDarkening * m_HEDarkening
Definition: HCalSD.h:84
void NaNTrap(G4Step *step)
G4Track * theTrack
Definition: CaloSD.h:118
virtual G4bool getStepInfo(G4Step *aStep)
Definition: CaloSD.cc:250
bool isItPMT(G4LogicalVolume *)
Definition: HCalSD.cc:729
G4int mupPDG
Definition: HCalSD.h:91
double eminHitHO
Definition: HCalSD.h:89
double tmaxHit
Definition: CaloSD.h:123
int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.h:34
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:44
TH1F * hzvem
Definition: HCalSD.h:98
G4StepPoint * preStepPoint
Definition: CaloSD.h:120
Definition: DetId.h:18
static const G4LogicalVolume * GetVolume(const std::string &name)
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
CaloHitID currentID
Definition: CaloSD.h:117
std::vector< Hit > getHits(G4Step *aStep, double weight)
void readWeightFromFile(std::string)
Definition: HCalSD.cc:1074
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
void getHitFibreBundle(G4Step *step, bool type)
Definition: HCalSD.cc:991
bool isPrimary() const
DDsvalues_type mergedSpecifics() const
double getHits(G4Step *aStep)
Definition: HFShowerPMT.cc:104
void initRun(G4ParticleTable *theParticleTable)
HcalNumberingFromDDD * numberingFromDDD
Definition: HCalSD.h:77
double layerWeight(int, const G4ThreeVector &, int, int)
Definition: HCalSD.cc:1100
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
float dose(int layer, float radius)
Definition: HFDarkening.cc:92
double getResponseWt(G4Track *)
Definition: CaloSD.cc:607
void initRun(G4ParticleTable *)
std::vector< G4String > fibreNames
Definition: HCalSD.h:94
G4bool hitExists()
Definition: CaloSD.cc:310
bool firstChild()
set the current node to the first child ...
std::vector< G4Material * > materials
Definition: HCalSD.h:95
virtual bool filterHit(CaloG4Hit *, double)
Definition: HCalSD.cc:624
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
Definition: DDFilter.cc:285
std::vector< G4String > hfNames
Definition: HCalSD.h:94
G4int mumPDG
Definition: HCalSD.h:91
std::map< uint32_t, double > layerWeights
Definition: HCalSD.h:97
tuple level
Definition: testEve_cfg.py:34
HFShowerPMT * showerPMT
Definition: HCalSD.h:82
uint32_t getUnitID() const
Definition: CaloG4Hit.h:69
int weight
Definition: histoStyle.py:50
HcalID unitID(int det, const CLHEP::Hep3Vector &pos, int depth, int lay=-1) const
float degradation(float intlumi, int ieta, int lay)
Definition: HEDarkening.cc:52
HFShowerLibrary * showerLibrary
Definition: HCalSD.h:79
std::vector< G4LogicalVolume * > fibre2LV
Definition: HCalSD.h:96
const std::string & name() const
Returns the name.
Definition: DDName.cc:87
const DDMaterial & material(void) const
Returns a reference object of the material this LogicalPart is made of.
CaloG4Hit * createNewHit()
Definition: CaloSD.cc:363
std::vector< G4String > matNames
Definition: HCalSD.h:94
bool isItFibre(G4LogicalVolume *)
Definition: HCalSD.cc:717
double getEnergyDeposit() const
Definition: CaloG4Hit.h:81
bool useHF
Definition: HCalSD.h:88
Definition: DDAxes.h:10
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
Definition: DDFilter.h:37