CMS 3D CMS Logo

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