CMS 3D CMS Logo

CastorShowerLibraryMaker.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: ShowerLibraryProducer
4 // Class : CastorShowerLibraryMaker
5 //
6 // Implementation:
7 // <Notes on implementation>
8 //
9 // Original Author: P. Katsas
10 // Created: 02/2007
11 //
12 // Adapted by W. Carvalho , 02/2009
13 //
15 
18 
24 
25 // Classes for shower library Root file
28 
37 
40 
42 
43 #include "G4RunManager.hh"
44 #include "G4SDManager.hh"
45 #include "G4Step.hh"
46 #include "G4Track.hh"
47 #include "G4Event.hh"
48 #include "G4PrimaryVertex.hh"
49 #include "G4VProcess.hh"
50 #include "G4HCofThisEvent.hh"
51 #include "G4UserEventAction.hh"
52 
53 #include <CLHEP/Random/Randomize.h>
54 #include <CLHEP/Units/SystemOfUnits.h>
55 #include <CLHEP/Units/PhysicalConstants.h>
56 #include <CLHEP/Units/SystemOfUnits.h>
57 #include <CLHEP/Units/GlobalPhysicalConstants.h>
58 
59 #include "TROOT.h"
60 #include "TFile.h"
61 #include "TH1.h"
62 #include "TH2.h"
63 #include "TProfile.h"
64 #include "TNtuple.h"
65 #include "TRandom.h"
66 #include "TLorentzVector.h"
67 #include "TUnixSystem.h"
68 #include "TSystem.h"
69 #include "TMath.h"
70 #include "TF1.h"
71 
72 #include <cassert>
73 #include <cmath>
74 #include <cstdlib>
75 #include <iomanip>
76 #include <iostream>
77 #include <map>
78 #include <set>
79 #include <sstream>
80 #include <string>
81 #include <vector>
82 
83 typedef std::vector<std::vector<std::vector<std::vector<CastorShowerEvent> > > > SLBin3D; // bin in energy, eta and phi
84 
86  public Observer<const BeginOfJob*>,
87  public Observer<const BeginOfRun*>,
88  public Observer<const EndOfRun*>,
89  public Observer<const BeginOfEvent*>,
90  public Observer<const EndOfEvent*>,
91  public Observer<const G4Step*> {
92 public:
94  ~CastorShowerLibraryMaker() override;
95 
96 private:
97  typedef int ebin;
98  typedef int etabin;
99  typedef int phibin;
100  // private structures
101  struct ShowerLib {
103  SLBin3D SLCollection; // the showers
104  std::vector<double> SLEnergyBins;
105  std::vector<double> SLEtaBins;
106  std::vector<double> SLPhiBins;
107  unsigned int nEvtPerBinE;
108  unsigned int nEvtPerBinEta;
109  unsigned int nEvtPerBinPhi;
110  std::vector<int> nEvtInBinE;
111  std::vector<std::vector<int> > nEvtInBinEta;
112  std::vector<std::vector<std::vector<int> > > nEvtInBinPhi;
113  };
114 
115  // observer classes
116  void update(const BeginOfJob* run) override;
117  void update(const BeginOfRun* run) override;
118  void update(const EndOfRun* run) override;
119  void update(const BeginOfEvent* evt) override;
120  void update(const EndOfEvent* evt) override;
121  void update(const G4Step* step) override;
122 
123 private:
124  void Finish();
125 
126  // Job general parameters
129 
130  unsigned int NPGParticle; // number of particles requested to Particle Gun
131  std::vector<int> PGParticleIDs; //p. gun particle IDs
132  bool DoHadSL; // true if hadronic SL should be produced
133  bool DoEmSL; // true if electromag. SL should be produced
134  bool InsideCastor; // true if particle step inside CASTOR
135  bool DeActivatePhysicsProcess; //cfg parameter: True if phys. proc. should be off from IP to Castor
136  std::vector<G4PrimaryParticle*> thePrims; // list of primaries for this event
137 
138  // Pointers for user defined class objects to be stored to Root file
145  ShowerLib* SLShowerptr; // pointer to the current shower collection (above)
146  std::map<int, std::set<int> > MapOfSecondaries; // map to hold all secondaries ID keyed by
147  // the PDG code of the primary
148 
149  std::map<int, G4ThreeVector> PrimaryMomentum;
150  std::map<int, G4ThreeVector> PrimaryPosition;
151  double MaxEta; // limits the eta region, the lower limit is given by the SL bins
152  double MaxPhi; // limits the phi region, the lower limit is given by the SL bins
153  // private methods
154  int FindEnergyBin(double e);
155  int FindEtaBin(double eta);
156  int FindPhiBin(double phi);
157  bool SLacceptEvent(int, int, int);
158  bool IsSLReady();
159  void GetKinematics(G4PrimaryParticle*, double& px, double& py, double& pz, double& pInit, double& eta, double& phi);
160  void GetKinematics(int, double& px, double& py, double& pz, double& pInit, double& eta, double& phi);
161 
162  std::vector<G4PrimaryParticle*> GetPrimary(const G4Event*);
164  void InitSLHolder(ShowerLib&);
165 
166  void printSLstatus(int, int, int);
167  int& SLnEvtInBinE(int ebin);
168  int& SLnEvtInBinEta(int ebin, int etabin);
169  int& SLnEvtInBinPhi(int ebin, int etabin, int phibin);
170  bool SLisEBinFilled(int ebin);
171  bool SLisEtaBinFilled(int ebin, int etabin);
172  bool SLisPhiBinFilled(int ebin, int etabin, int phibin);
173  void KillSecondaries(const G4Step* step);
174  void GetMissingEnergy(CaloG4HitCollection*, double&, double&);
175 
176  // Root pointers
177  TFile* theFile;
178  TTree* theTree;
179 
181  int stepIndex; // ignore, please
182 };
183 
185  : NPGParticle(0),
186  DoHadSL(false),
187  DoEmSL(false),
188  DeActivatePhysicsProcess(false),
189  emShower(nullptr),
190  hadShower(nullptr) {
191  MapOfSecondaries.clear();
192  hadInfo = nullptr;
193  emInfo = nullptr;
194  edm::ParameterSet p_SLM = p.getParameter<edm::ParameterSet>("CastorShowerLibraryMaker");
195  verbosity = p_SLM.getParameter<int>("Verbosity");
196  eventNtFileName = p_SLM.getParameter<std::string>("EventNtupleFileName");
197  hadSLHolder.nEvtPerBinPhi = p_SLM.getParameter<int>("nhadEvents");
198  emSLHolder.nEvtPerBinPhi = p_SLM.getParameter<int>("nemEvents");
199  hadSLHolder.SLEnergyBins = p_SLM.getParameter<std::vector<double> >("SLhadEnergyBins");
200  hadSLHolder.SLEtaBins = p_SLM.getParameter<std::vector<double> >("SLhadEtaBins");
201  hadSLHolder.SLPhiBins = p_SLM.getParameter<std::vector<double> >("SLhadPhiBins");
202  emSLHolder.SLEnergyBins = p_SLM.getParameter<std::vector<double> >("SLemEnergyBins");
203  emSLHolder.SLEtaBins = p_SLM.getParameter<std::vector<double> >("SLemEtaBins");
204  emSLHolder.SLPhiBins = p_SLM.getParameter<std::vector<double> >("SLemPhiBins");
205  PGParticleIDs = p_SLM.getParameter<std::vector<int> >("PartID");
206  DeActivatePhysicsProcess = p_SLM.getParameter<bool>("DeActivatePhysicsProcess");
207  MaxPhi = p_SLM.getParameter<double>("SLMaxPhi");
208  MaxEta = p_SLM.getParameter<double>("SLMaxEta");
209  //
210  NPGParticle = PGParticleIDs.size();
211  for (unsigned int i = 0; i < PGParticleIDs.size(); i++) {
212  switch (std::abs(PGParticleIDs.at(i))) {
213  case 11:
214  case 22:
215  DoEmSL = true;
216  break;
217  default:
218  DoHadSL = true;
219  }
220  }
225 
226  edm::LogVerbatim("HcalSim") << "============================================================================";
227  edm::LogVerbatim("HcalSim") << "CastorShowerLibraryMaker:: Initialized as observer";
228  edm::LogVerbatim("HcalSim") << " Event Ntuple will be created";
229  edm::LogVerbatim("HcalSim") << " Event Ntuple file: " << eventNtFileName;
230  edm::LogVerbatim("HcalSim") << " Number of Hadronic events in E bins: " << hadSLHolder.nEvtPerBinE;
231  edm::LogVerbatim("HcalSim") << " Number of Hadronic events in Eta bins: " << hadSLHolder.nEvtPerBinEta;
232  edm::LogVerbatim("HcalSim") << " Number of Hadronic events in Phi bins: " << hadSLHolder.nEvtPerBinPhi;
233  edm::LogVerbatim("HcalSim") << " Number of Electromag. events in E bins: " << emSLHolder.nEvtPerBinE;
234  edm::LogVerbatim("HcalSim") << " Number of Electromag. events in Eta bins: " << emSLHolder.nEvtPerBinEta;
235  edm::LogVerbatim("HcalSim") << " Number of Electromag. events in Phi bins: " << emSLHolder.nEvtPerBinPhi;
236  edm::LogVerbatim("HcalSim") << "============================================================================\n";
237 
238  // Initializing the SL collections
241 }
243  int nBinsE, nBinsEta, nBinsPhi, nEvtPerBinPhi;
244  nBinsE = showerholder.SLEnergyBins.size();
245  nBinsEta = showerholder.SLEtaBins.size();
246  nBinsPhi = showerholder.SLPhiBins.size();
247  nEvtPerBinPhi = showerholder.nEvtPerBinPhi;
248  //
249  // Info
250  //
251  showerholder.SLInfo.Energy.setNEvts(nEvtPerBinPhi * nBinsPhi * nBinsEta * nBinsE);
252  showerholder.SLInfo.Energy.setNEvtPerBin(nEvtPerBinPhi * nBinsPhi * nBinsEta);
253  showerholder.SLInfo.Energy.setNBins(nBinsE);
254  showerholder.SLInfo.Energy.setBin(showerholder.SLEnergyBins);
255  //
256  showerholder.SLInfo.Eta.setNEvts(nEvtPerBinPhi * nBinsPhi * nBinsEta);
257  showerholder.SLInfo.Eta.setNEvtPerBin(nEvtPerBinPhi * nBinsPhi);
258  showerholder.SLInfo.Eta.setNBins(nBinsEta);
259  showerholder.SLInfo.Eta.setBin(showerholder.SLEtaBins);
260  //
261  showerholder.SLInfo.Phi.setNEvts(nEvtPerBinPhi * nBinsPhi);
262  showerholder.SLInfo.Phi.setNEvtPerBin(nEvtPerBinPhi);
263  showerholder.SLInfo.Phi.setNBins(nBinsPhi);
264  showerholder.SLInfo.Phi.setBin(showerholder.SLPhiBins);
265  //
266  // Shower
267  showerholder.SLCollection.assign(nBinsE, std::vector<std::vector<std::vector<CastorShowerEvent> > >());
268  showerholder.nEvtInBinE.assign(nBinsE, 0);
269  showerholder.nEvtInBinEta.assign(nBinsE, std::vector<int>(0));
270  showerholder.nEvtInBinPhi.assign(nBinsE, std::vector<std::vector<int> >());
271  for (int i = 0; i < nBinsE; i++) {
272  showerholder.SLCollection.at(i).assign(nBinsEta, std::vector<std::vector<CastorShowerEvent> >());
273  showerholder.nEvtInBinEta.at(i).assign(nBinsEta, 0);
274  showerholder.nEvtInBinPhi.at(i).assign(nBinsEta, std::vector<int>(0));
275  for (int j = 0; j < nBinsEta; j++) {
276  showerholder.SLCollection.at(i).at(j).assign(nBinsPhi, std::vector<CastorShowerEvent>());
277  showerholder.nEvtInBinPhi.at(i).at(j).assign(nBinsPhi, 0);
278  for (int k = 0; k < nBinsPhi; k++)
279  showerholder.SLCollection.at(i).at(j).at(k).assign(nEvtPerBinPhi, CastorShowerEvent());
280  }
281  }
282 }
283 
284 //===============================================================================================
285 
287  Finish();
288 
289  edm::LogVerbatim("HcalSim") << "CastorShowerLibraryMaker: End of process";
290 }
291 
292 //=================================================================== per EVENT
294  edm::LogVerbatim("HcalSim") << " CastorShowerLibraryMaker::Starting new job ";
295 }
296 
297 //==================================================================== per RUN
299  edm::LogVerbatim("HcalSim") << "\nCastorShowerLibraryMaker: Starting Run";
300 
301  edm::LogVerbatim("HcalSim") << "CastorShowerLibraryMaker: output event root file created";
302 
303  TString eventfilename = eventNtFileName;
304  theFile = new TFile(eventfilename, "RECREATE");
305  theTree = new TTree("CastorCherenkovPhotons", "Cherenkov Photons");
306 
307  Int_t split = 1;
308  Int_t bsize = 64000;
310  emShower = new CastorShowerEvent();
313  // Create Branchs
314  theTree->Branch("emShowerLibInfo.", "CastorShowerLibraryInfo", &emInfo, bsize, split);
315  theTree->Branch("emParticles.", "CastorShowerEvent", &emShower, bsize, split);
316  theTree->Branch("hadShowerLibInfo.", "CastorShowerLibraryInfo", &hadInfo, bsize, split);
317  theTree->Branch("hadParticles.", "CastorShowerEvent", &hadShower, bsize, split);
318 
319  // set the Info for electromagnetic shower
320  // set the energy bins info
325  // set the eta bins info
330  // set the eta bins info
335  // The same for the hadronic shower
336  // set the energy bins info
341  // set the eta bins info
346  // set the eta bins info
351  // int flag = theTree->GetBranch("CastorShowerLibInfo")->Fill();
352  // Loop on all leaves of this branch to fill Basket buffer.
353  // The function returns the number of bytes committed to the memory basket.
354  // If a write error occurs, the number of bytes returned is -1.
355  // If no data are written, because e.g. the branch is disabled,
356  // the number of bytes returned is 0.
357  // if(flag==-1) {
358  // edm::LogVerbatim("CastorAnalyzer") << " WARNING: Error writing to Branch \"CastorShowerLibInfo\" \n" ;
359  // } else
360  // if(flag==0) {
361  // edm::LogVerbatim("CastorAnalyzer") << " WARNING: No data written to Branch \"CastorShowerLibInfo\" \n" ;
362  // }
363 
364  // Initialize "accounting" variables
365 
366  eventIndex = 0;
367 }
368 
369 //=================================================================== per EVENT
371  eventIndex++;
372  stepIndex = 0;
373  InsideCastor = false;
374  PrimaryMomentum.clear();
375  PrimaryPosition.clear();
376  int NAccepted = 0;
377  // reset the pointers to the shower objects
378  SLShowerptr = nullptr;
379  MapOfSecondaries.clear();
380  thePrims.clear();
381  G4EventManager* e_mgr = G4EventManager::GetEventManager();
382  if (IsSLReady()) {
383  printSLstatus(-1, -1, -1);
384  update((EndOfRun*)nullptr);
385  return;
386  }
387 
388  thePrims = GetPrimary((*evt)());
389  for (unsigned int i = 0; i < thePrims.size(); i++) {
390  G4PrimaryParticle* thePrim = thePrims.at(i);
391  int particleType = thePrim->GetPDGcode();
392 
393  std::string SLType("");
394  if (particleType == 11) {
396  SLType = "Electromagnetic";
397  } else {
399  SLType = "Hadronic";
400  }
401  double px = 0., py = 0., pz = 0., pInit = 0., eta = 0., phi = 0.;
402  GetKinematics(thePrim, px, py, pz, pInit, eta, phi);
403  int ebin = FindEnergyBin(pInit);
404  int etabin = FindEtaBin(eta);
405  int phibin = FindPhiBin(phi);
406  if (verbosity)
407  edm::LogVerbatim("HcalSim") << "\n========================================================================"
408  << "\nBeginOfEvent: E : " << pInit << "\t bin : " << ebin
409  << "\n Eta : " << eta << "\t bin : " << etabin
410  << "\n Phi : " << phi << "\t bin : " << phibin
411  << "\n========================================================================";
412 
413  if (ebin < 0 || etabin < 0 || phibin < 0)
414  continue;
415  bool accept = false;
416  if (!(SLacceptEvent(ebin, etabin, phibin))) {
417  /*
418 // To increase the chance of a particle arriving at CASTOR inside a not full bin,
419 // check if there is available phase space in the neighboring bins
420  unsigned int ebin_min = std::max(0,ebin-3);
421  unsigned int eta_bin_min = std::max(0,etabin-2);
422  unsigned int eta_bin_max = std::min(etabin,etabin+2);
423  unsigned int phi_bin_min = std::max(0,phibin-2);
424  unsigned int phi_bin_max = std::min(phibin,phibin+2);
425  for(unsigned int i_ebin=ebin_min;i_ebin<=(unsigned int)ebin;i_ebin++) {
426  for (unsigned int i_etabin=eta_bin_min;i_etabin<=eta_bin_max;i_etabin++) {
427  for (unsigned int i_phibin=phi_bin_min;i_phibin<=phi_bin_max;i_phibin++) {
428  if (SLacceptEvent((int)i_ebin,(int)i_etabin,(int)i_phibin)) {accept=true;break;}
429  }
430  if (accept) break;
431  }
432  if (accept) break;
433  }
434 */
435  if (!accept)
436  edm::LogVerbatim("CastorShowerLibraryMaker")
437  << "Event not accepted for ebin=" << ebin << ",etabin=" << etabin << ",phibin=" << phibin;
438  } else {
439  accept = true;
440  }
441  if (accept)
442  NAccepted++;
443  }
444 
445  if (NAccepted == 0) {
446  const_cast<G4Event*>((*evt)())->SetEventAborted();
447  const_cast<G4Event*>((*evt)())->KeepTheEvent((G4bool) false);
448  e_mgr->AbortCurrentEvent();
449  }
450  SLShowerptr = nullptr;
451  //
452  edm::LogVerbatim("HcalSim") << "CastorShowerLibraryMaker: Processing Event Number: " << eventIndex;
453 }
454 
455 //=================================================================== per STEP
456 void CastorShowerLibraryMaker::update(const G4Step* aStep) {
457  static thread_local int CurrentPrimary = 0;
458  G4Track* trk = aStep->GetTrack();
459  int pvec_size;
460  if (trk->GetCurrentStepNumber() == 1) {
461  if (trk->GetParentID() == 0) {
462  CurrentPrimary = (int)trk->GetDynamicParticle()->GetPDGcode();
463  if (CurrentPrimary == 0)
464  SimG4Exception("CastorShowerLibraryMaker::update(G4Step) -> Primary particle undefined");
465  InsideCastor = false;
466  // Deactivate the physics process
468  G4ProcessManager* p_mgr = trk->GetDefinition()->GetProcessManager();
469  G4ProcessVector* pvec = p_mgr->GetProcessList();
470  pvec_size = pvec->size();
471  for (int i = 0; i < pvec_size; i++) {
472  G4VProcess* proc = (*pvec)(i);
473  if (proc->GetProcessName() != "Transportation" && proc->GetProcessName() != "Decay") {
474  edm::LogVerbatim("HcalSim") << "DeActivating process: " << proc->GetProcessName();
475  p_mgr->SetProcessActivation(proc, false);
476  }
477  }
478  }
479  // move track to z of CASTOR
480  G4ThreeVector pos;
481  pos.setZ(-14390);
482  double t = std::abs((pos.z() - trk->GetPosition().z())) / trk->GetVelocity();
483  double r = (pos.z() - trk->GetPosition().z()) / trk->GetMomentum().cosTheta();
484  pos.setX(r * sin(trk->GetMomentum().theta()) * cos(trk->GetMomentum().phi()) + trk->GetPosition().x());
485  pos.setY(r * sin(trk->GetMomentum().theta()) * sin(trk->GetMomentum().phi()) + trk->GetPosition().y());
486  trk->SetPosition(pos);
487  trk->SetGlobalTime(trk->GetGlobalTime() + t);
488  trk->AddTrackLength(r);
489  } else if (!InsideCastor) {
490  edm::LogVerbatim("HcalSim") << "CastorShowerLibraryMaker::update(G4Step) -> Killing spurious track";
491  trk->SetTrackStatus(fKillTrackAndSecondaries);
492  return;
493  }
494  MapOfSecondaries[CurrentPrimary].insert((int)trk->GetTrackID());
495  }
496  // Checks if primary already inside CASTOR
497  std::string CurVolume = trk->GetVolume()->GetName();
498  if (!InsideCastor && (
499  //CurVolume=="C3EF"||CurVolume=="C4EF"||CurVolume=="CAEL"||
500  //CurVolume=="CAHL"||CurVolume=="C3HF"||CurVolume=="C4HF")) {
501  //CurVolume=="CastorB"||
502  CurVolume == "CAST")) {
503  //CurVolume=="CAIR")) {
504  InsideCastor = true;
505  // Activate the physics process
506  if (trk->GetParentID() == 0 && DeActivatePhysicsProcess) {
507  G4ProcessManager* p_mgr = trk->GetDefinition()->GetProcessManager();
508  G4ProcessVector* pvec = p_mgr->GetProcessList();
509  pvec_size = pvec->size();
510  for (int i = 0; i < pvec_size; i++) {
511  G4VProcess* proc = (*pvec)(i);
512  if (proc->GetProcessName() != "Transportation" && proc->GetProcessName() != "Decay") {
513  edm::LogVerbatim("HcalSim") << " Activating process: " << proc->GetProcessName();
514  p_mgr->SetProcessActivation(proc, true);
515  }
516  }
517  }
518  //PrimaryMomentum[CurrentPrimary]=aStep->GetPreStepPoint()->GetMomentum();
519  // check fiducial eta and phi
520  if (trk->GetMomentum().phi() > MaxPhi || trk->GetMomentum().eta() > MaxEta) {
521  trk->SetTrackStatus(fKillTrackAndSecondaries);
522  InsideCastor = false;
523  return;
524  }
525  PrimaryMomentum[CurrentPrimary] = trk->GetMomentum();
526  PrimaryPosition[CurrentPrimary] = trk->GetPosition();
527  KillSecondaries(aStep);
528  return;
529  }
530  // Kill the secondaries if they have been produced before entering castor
531  if (CurrentPrimary != 0 && trk->GetParentID() == 0 && !InsideCastor) {
532  KillSecondaries(aStep);
533  if (verbosity) {
534  double pre_phi = aStep->GetPreStepPoint()->GetMomentum().phi();
535  double cur_phi = trk->GetMomentum().phi();
536  if (pre_phi != cur_phi) {
537  edm::LogVerbatim("HcalSim") << "Primary track phi : " << pre_phi << " changed in current step: " << cur_phi
538  << " by processes:";
539  const G4VProcess* proc = aStep->GetPreStepPoint()->GetProcessDefinedStep();
540  edm::LogVerbatim("HcalSim") << " " << proc->GetProcessName() << " In volume " << CurVolume;
541  }
542  }
543  }
544 
545  //==============================================
546  /*
547 */
548  /*
549  if(aStep->IsFirstStepInVolume()) {
550  edm::LogVerbatim("CastorShowerLibraryMaker") << "CastorShowerLibraryMaker::update(const G4Step * aStep):"
551  << "\n IsFirstStepInVolume , "
552  << "time = " << aStep->GetTrack()->GetGlobalTime() ;
553  }
554  stepIndex++;
555 */
556 }
557 
558 //================= End of EVENT ===============
560  // check if the job is done!
561  if ((*evt)()->IsAborted()) {
562  edm::LogVerbatim("HcalSim") << "\n========================================================================"
563  << "\nEndOfEvent: EVENT ABORTED"
564  << "\n========================================================================";
565  return;
566  }
567  //DynamicRangeFlatRandomEGunProducer* pgun = edm::DynamicRangeFlatRandomEGunKernel::get_instance();
568  //edm::LogVerbatim("HcalSim") << pgun->EGunMaxE();
569  /*
570  edm::LogVerbatim("HcalSim") << "Minimum energy in Particle Gun : " << pgun->EGunMinE() << "\nMaximum energy in Particle Gun : " << pgun->EGunMaxE();
571 */
572  if (verbosity)
573  edm::LogVerbatim("HcalSim") << "CastorShowerLibraryMaker: End of Event: " << eventIndex;
574  // Get the pointer to the primary particle
575  if (thePrims.empty()) {
576  edm::LogVerbatim("CastorShowerLibraryMaker") << "No valid primary particle found. Skipping event";
577  return;
578  }
579  // access to the G4 hit collections
580  G4HCofThisEvent* allHC = (*evt)()->GetHCofThisEvent();
581  int CAFIid = G4SDManager::GetSDMpointer()->GetCollectionID("CastorFI");
582  CaloG4HitCollection* theCAFI = (CaloG4HitCollection*)allHC->GetHC(CAFIid);
583  if (verbosity)
584  edm::LogVerbatim("CastorShowerLibraryMaker") << " update(*evt) --> accessed all HC ";
585  edm::LogVerbatim("CastorShowerLibraryMaker") << "Found " << theCAFI->entries() << " hits in G4HitCollection";
586  if (theCAFI->entries() == 0) {
587  edm::LogVerbatim("CastorShowerLibraryMaker") << "\n Empty G4HitCollection";
588  return;
589  }
590 
591  // Loop over primaries
592  int NEvtAccepted = 0;
593  int NHitInEvent = 0;
594  for (unsigned int i = 0; i < thePrims.size(); i++) {
595  G4PrimaryParticle* thePrim = thePrims.at(i);
596  if (!thePrim) {
597  edm::LogVerbatim("CastorShowerLibraryMaker") << "nullptr Pointer to the primary";
598  continue;
599  }
600  // Check primary particle type
601  int particleType = thePrim->GetPDGcode();
602 
603  // set the pointer to the shower collection
604  std::string SLType("");
605  if (particleType == 11) {
607  SLType = "Electromagnetic";
608  } else {
610  SLType = "Hadronic";
611  }
612  edm::LogVerbatim("CastorShowerLibraryMaker") << "\n Primary (thePrim) trackID is " << thePrim->GetTrackID() << "\n";
613 
614  // Obtain primary particle's initial momentum (pInit)
615  double px = 0., py = 0., pz = 0., pInit = 0., eta = 0., phi = 0.;
616  GetKinematics(particleType, px, py, pz, pInit, eta, phi);
617  // Check if current event falls into any bin
618  // first: energy
619  if (pInit == 0) {
620  edm::LogVerbatim("CastorShowerLibraryMaker") << "Primary did not hit CASTOR";
621  continue;
622  }
623  int ebin = FindEnergyBin(pInit);
624  int etabin = FindEtaBin(eta);
625  int phibin = FindPhiBin(phi);
626  edm::LogVerbatim("HcalSim") << SLType;
628  if (!SLacceptEvent(ebin, etabin, phibin)) {
629  edm::LogVerbatim("CastorShowerLibraryMaker")
630  << "Event not accepted for ebin=" << ebin << ",etabin=" << etabin << ",phibin=" << phibin << "(" << pInit
631  << "," << eta << "," << phi << ")";
632  continue;
633  }
634  //
635  // event passed. Fill the vector accordingly
636  //
637  // Look for the Hit Collection
638  edm::LogVerbatim("CastorShowerLibraryMaker")
639  << "\n CastorShowerLibraryMaker::update(EndOfEvent * evt) - event #" << (*evt)()->GetEventID();
640 
641  /*
642  edm::LogVerbatim("HcalSim") << "Number of collections : " << allHC->GetNumberOfCollections();
643  for(int ii = 0;ii<allHC->GetNumberOfCollections();ii++)
644  edm::LogVerbatim("HcalSim") << "Name of collection " << ii << " : " << allHC->GetHC(ii)->GetName();
645 */
646 
647  CastorShowerEvent* shower = nullptr;
648  int cur_evt_idx = SLShowerptr->nEvtInBinPhi.at(ebin).at(etabin).at(phibin);
649  shower = &(SLShowerptr->SLCollection.at(ebin).at(etabin).at(phibin).at(cur_evt_idx));
650 
651  // Get Hit information
652  if (FillShowerEvent(theCAFI, shower, particleType)) {
653  // Primary particle information
654  /*
655  edm::LogVerbatim("CastorShowerLibraryMaker") << "New SL event: Primary = " << particleType << "; Energy = " << pInit << "; Eta = " << eta << "; Phi = " << phi << "; Nhits = " << shower->getNhit();
656 */
657  shower->setPrimE(pInit);
658  shower->setPrimEta(eta);
659  shower->setPrimPhi(phi);
663  SLnEvtInBinE(ebin)++;
666  NHitInEvent += shower->getNhit();
667  NEvtAccepted++;
668  } else {
669  shower->Clear();
670  }
671  }
672  // Check for unassociated energy
673  int thecafi_entries = theCAFI->entries();
674  if (NEvtAccepted == int(thePrims.size()) && thecafi_entries != NHitInEvent) {
675  edm::LogWarning("HcalSim") << "WARNING: Inconsistent Number of Hits -> Hits in collection: " << theCAFI->entries()
676  << " Hits in the showers: " << NHitInEvent;
677  double miss_energy = 0;
678  double tot_energy = 0;
679  GetMissingEnergy(theCAFI, miss_energy, tot_energy);
680  if (miss_energy > 0) {
681  edm::LogVerbatim("HcalSim") << "Total missing energy: " << miss_energy
682  << " for an incident energy: " << tot_energy;
683  }
684  }
685 
686  /*
687  for (int i=emSLHolder.SLEnergyBins.size()-1;i>0;i--) {
688  if (emSLHolder.nEvtInBinE.at(i)==(int)emSLHolder.nEvtPerBinE) {
689  std::ostringstream out;
690  out << emSLHolder.SLEnergyBins.at(i);
691  edm::LogVerbatim("HcalSim") << "Bin Limit: " << out.str();
692  setenv("CASTOR_SL_PG_MAXE",out.str().c_str(),1);
693  }
694  break;
695  }
696 */
697  //int iEvt = (*evt)()->GetEventID();
698  //double xint;
699  /*
700  if (modf(log10(iEvt),&xint)==0)
701  edm::LogVerbatim("HcalSim") << " CastorShowerLibraryMaker Event " << iEvt;
702 */
703  // edm::LogVerbatim("HcalSim") << "\n===>>> Done writing user histograms ";
704 }
705 
706 //========================= End of RUN ======================
708  // Fill the tree with the collected objects
709  if (!IsSLReady())
710  SimG4Exception("\n\nShower Library NOT READY.\n\n");
711 
712  unsigned int ibine, ibineta, ibinphi, ievt; // indexes for em shower
713  unsigned int jbine, jbineta, jbinphi, jevt; // indexes for had shower
714 
715  ibine = ibineta = ibinphi = ievt = jbine = jbineta = jbinphi = jevt = 0;
716 
717  int nEvtInTree = 0;
718  int nEMevt = emSLHolder.nEvtPerBinE * emSLHolder.SLEnergyBins.size();
719  int nHadevt = hadSLHolder.nEvtPerBinE * hadSLHolder.SLEnergyBins.size();
720  int maxEvtInTree = std::max(nEMevt, nHadevt);
721 
724 
725  while (nEvtInTree < maxEvtInTree) {
726  if (emShower)
727  emShower->Clear();
728  if (hadShower)
729  hadShower->Clear();
730  while (ibine < emSLHolder.SLEnergyBins.size() && nEMevt > 0) {
731  emShower = &(emSLHolder.SLCollection.at(ibine).at(ibineta).at(ibinphi).at(ievt));
732  ievt++;
733  if (ievt == emSLHolder.nEvtPerBinPhi) {
734  ievt = 0;
735  ibinphi++;
736  }
737  if (ibinphi == emSLHolder.SLPhiBins.size()) {
738  ibinphi = 0;
739  ibineta++;
740  }
741  if (ibineta == emSLHolder.SLEtaBins.size()) {
742  ibineta = 0;
743  ibine++;
744  }
745  break;
746  }
747  while (jbine < hadSLHolder.SLEnergyBins.size() && nHadevt > 0) {
748  hadShower = &(hadSLHolder.SLCollection.at(jbine).at(jbineta).at(jbinphi).at(jevt));
749  jevt++;
750  if (jevt == hadSLHolder.nEvtPerBinPhi) {
751  jevt = 0;
752  jbinphi++;
753  }
754  if (jbinphi == hadSLHolder.SLPhiBins.size()) {
755  jbinphi = 0;
756  jbineta++;
757  }
758  if (jbineta == hadSLHolder.SLEtaBins.size()) {
759  jbineta = 0;
760  jbine++;
761  }
762  break;
763  }
764  theTree->Fill();
765  nEvtInTree++;
766  if (nEvtInTree == 1) {
767  theTree->SetBranchStatus("emShowerLibInfo.", false);
768  theTree->SetBranchStatus("hadShowerLibInfo.", false);
769  }
770  }
771  // check if run is nullptr and exit
772  if (run == nullptr)
773  throw SimG4Exception("\n\nNumber of needed trigger events reached in CastorShowerLibraryMaker\n\n");
774 }
775 
776 //============================================================
778  // if (doNTcastorevent) {
779 
780  theFile->cd();
781  theTree->Write("", TObject::kOverwrite);
782  edm::LogVerbatim("HcalSim") << "CastorShowerLibraryMaker: Ntuple event written";
783  theFile->Close();
784  edm::LogVerbatim("HcalSim") << "CastorShowerLibraryMaker: Event file closed";
785 
786  // Delete pointers to objects, now that TTree has been written and TFile closed
787  // delete info;
788  // delete emShower;
789  // delete hadShower;
790  // }
791 }
793  //
794  // returns the integer index of the energy bin, taken from SLenergies vector
795  // returns -1 if ouside valid range
796  //
797  if (!SLShowerptr) {
798  edm::LogVerbatim("CastorShowerLibraryMaker") << "\n\nFindEnergyBin can be called only after BeginOfEvent\n\n";
799  throw SimG4Exception("\n\nnullptr Pointer to the shower library.\n\n");
800  }
801  const std::vector<double>& SLenergies = SLShowerptr->SLEnergyBins;
802  if (energy >= SLenergies.back())
803  return SLenergies.size() - 1;
804 
805  unsigned int i = 0;
806  for (; i < SLenergies.size() - 1; i++)
807  if (energy >= SLenergies.at(i) && energy < SLenergies.at(i + 1))
808  return (int)i;
809 
810  // now i points to the last but 1 bin
811  if (energy >= SLenergies.at(i))
812  return (int)i;
813  // energy outside bin range
814  return -1;
815 }
817  //
818  // returns the integer index of the eta bin, taken from SLetas vector
819  // returns -1 if ouside valid range
820  //
821  if (!SLShowerptr) {
822  edm::LogVerbatim("CastorShowerLibraryMaker") << "\n\nFindEtaBin can be called only after BeginOfEvent\n\n";
823  throw SimG4Exception("\n\nnullptr Pointer to the shower library.\n\n");
824  }
825  const std::vector<double>& SLetas = SLShowerptr->SLEtaBins;
826  if (eta >= SLetas.back())
827  return SLetas.size() - 1;
828  unsigned int i = 0;
829  for (; i < SLetas.size() - 1; i++)
830  if (eta >= SLetas.at(i) && eta < SLetas.at(i + 1))
831  return (int)i;
832  // now i points to the last but 1 bin
833  if (eta >= SLetas.at(i))
834  return (int)i;
835  // eta outside bin range
836  return -1;
837 }
839  //
840  // returns the integer index of the phi bin, taken from SLphis vector
841  // returns -1 if ouside valid range
842  //
843  // needs protection in case phi is outside range -pi,pi
844  //
845  if (!SLShowerptr) {
846  edm::LogVerbatim("CastorShowerLibraryMaker") << "\n\nFindPhiBin can be called only after BeginOfEvent\n\n";
847  throw SimG4Exception("\n\nnullptr Pointer to the shower library.\n\n");
848  }
849  const std::vector<double>& SLphis = SLShowerptr->SLPhiBins;
850  if (phi >= SLphis.back())
851  return SLphis.size() - 1;
852  unsigned int i = 0;
853  for (; i < SLphis.size() - 1; i++)
854  if (phi >= SLphis.at(i) && phi < SLphis.at(i + 1))
855  return (int)i;
856  // now i points to the last but 1 bin
857  if (phi >= SLphis.at(i))
858  return (int)i;
859  // phi outside bin range
860  return -1;
861 }
863  // at this point, the pointer to the shower library should be nullptr
864  if (SLShowerptr) {
865  edm::LogVerbatim("CastorShowerLibraryMaker") << "\n\nIsSLReady must be called when a new event starts.\n\n";
866  throw SimG4Exception("\n\nNOT nullptr Pointer to the shower library.\n\n");
867  }
868  // it is enough to check if all the energy bin is filled
869  if (DoEmSL) {
871  for (unsigned int i = 0; i < SLShowerptr->SLEnergyBins.size(); i++) {
872  if (!SLisEBinFilled(i)) {
873  SLShowerptr = nullptr;
874  return false;
875  }
876  }
877  }
878  if (DoHadSL) {
880  for (unsigned int i = 0; i < SLShowerptr->SLEnergyBins.size(); i++) {
881  if (!SLisEBinFilled(i)) {
882  SLShowerptr = nullptr;
883  return false;
884  }
885  }
886  }
887  SLShowerptr = nullptr;
888  return true;
889 }
891  int thePrim, double& px, double& py, double& pz, double& pInit, double& eta, double& phi) {
892  if (thePrim == 0)
893  return;
894  if (PrimaryMomentum.find(thePrim) == PrimaryMomentum.end())
895  return;
896  px = PrimaryMomentum[thePrim].x() / GeV;
897  py = PrimaryMomentum[thePrim].y() / GeV;
898  pz = PrimaryMomentum[thePrim].z() / GeV;
899  pInit = PrimaryMomentum[thePrim].mag() / GeV;
900  if (pInit == 0)
901  return;
902  double costheta = pz / pInit;
903  double theta = acos(std::min(std::max(costheta, double(-1.)), double(1.)));
904  eta = -log(tan(theta / 2.0));
905  phi = (px == 0 && py == 0) ? 0 : atan2(py, px); // the recommended way of calculating phi
906  phi = PrimaryMomentum[thePrim].phi();
907 }
909  G4PrimaryParticle* thePrim, double& px, double& py, double& pz, double& pInit, double& eta, double& phi) {
910  px = py = pz = phi = eta = 0.0;
911  if (thePrim == nullptr)
912  return;
913  px = thePrim->GetMomentum().x() / GeV;
914  py = thePrim->GetMomentum().y() / GeV;
915  pz = thePrim->GetMomentum().z() / GeV;
916  pInit = thePrim->GetMomentum().mag() / GeV;
917  //pInit = sqrt(pow(px,2.)+pow(py,2.)+pow(pz,2.));
918  if (pInit == 0)
919  return;
920  double costheta = pz / pInit;
921  double theta = acos(std::min(std::max(costheta, double(-1.)), double(1.)));
922  eta = -log(tan(theta / 2.0));
923  phi = (px == 0 && py == 0) ? 0 : atan2(py, px); // the recommended way of calculating phi
924  phi = thePrim->GetMomentum().phi();
925  //if (px!=0) phi=atan(py/px);
926 }
927 std::vector<G4PrimaryParticle*> CastorShowerLibraryMaker::GetPrimary(const G4Event* evt) {
928  // Find Primary info:
929  int trackID = 0;
930  std::vector<G4PrimaryParticle*> thePrims;
931  G4PrimaryParticle* thePrim = nullptr;
932  G4int nvertex = evt->GetNumberOfPrimaryVertex();
933  edm::LogVerbatim("CastorShowerLibraryMaker") << "Event has " << nvertex << " vertex";
934  if (nvertex != 1) {
935  edm::LogVerbatim("CastorShowerLibraryMaker") << "CastorShowerLibraryMaker::GetPrimary ERROR: no vertex";
936  return thePrims;
937  }
938 
939  for (int i = 0; i < nvertex; i++) {
940  G4PrimaryVertex* avertex = evt->GetPrimaryVertex(i);
941  if (avertex == nullptr) {
942  edm::LogVerbatim("CastorShowerLibraryMaker")
943  << "CastorShowerLibraryMaker::GetPrimary ERROR: pointer to vertex = 0";
944  continue;
945  }
946  unsigned int npart = avertex->GetNumberOfParticle();
947  if (npart != NPGParticle)
948  continue;
949  for (unsigned int j = 0; j < npart; j++) {
950  unsigned int k = 0;
951  //int test_pID = 0;
952  trackID = j;
953  thePrim = avertex->GetPrimary(trackID);
954  while (k < NPGParticle && PGParticleIDs.at(k++) != thePrim->GetPDGcode()) {
955  ;
956  };
957  if (k > NPGParticle)
958  continue; // ID not found in the requested particles
959  thePrims.push_back(thePrim);
960  }
961  }
962  return thePrims;
963 }
965  if (!SLShowerptr) {
966  edm::LogVerbatim("CastorShowerLibraryInfo") << "nullptr shower pointer. Printing both";
967  edm::LogVerbatim("HcalSim") << "Electromagnetic";
969  this->printSLstatus(ebin, etabin, phibin);
970  edm::LogVerbatim("HcalSim") << "Hadronic";
972  this->printSLstatus(ebin, etabin, phibin);
973  SLShowerptr = nullptr;
974  return;
975  }
976  int nBinsE = SLShowerptr->SLEnergyBins.size();
977  int nBinsEta = SLShowerptr->SLEtaBins.size();
978  int nBinsPhi = SLShowerptr->SLPhiBins.size();
979  std::vector<double> SLenergies = SLShowerptr->SLEnergyBins;
980  std::ostringstream st1;
981  for (int n = 0; n < 11 + (nBinsEta * nBinsPhi); n++)
982  st1 << "=";
983  edm::LogVerbatim("HcalSim") << st1.str();
984  for (int i = 0; i < nBinsE; i++) {
985  std::ostringstream st1;
986  st1 << "E bin " << std::setw(6) << SLenergies.at(i) << " : ";
987  for (int j = 0; j < nBinsEta; j++) {
988  for (int k = 0; k < nBinsPhi; k++) {
989  (SLisPhiBinFilled(i, j, k)) ? st1 << "1" : st1 << "-";
990  }
991  if (j < nBinsEta - 1)
992  st1 << "|";
993  }
994  st1 << " (" << SLnEvtInBinE(i) << " events)";
995  edm::LogVerbatim("HcalSim") << st1.str();
996  if (ebin != i)
997  continue;
998  std::ostringstream st2;
999  st2 << " ";
1000  for (int j = 0; j < nBinsEta; j++) {
1001  for (int k = 0; k < nBinsPhi; k++) {
1002  (ebin == i && etabin == j && phibin == k) ? st2 << "^" : st2 << " ";
1003  }
1004  if (j < nBinsEta - 1)
1005  st2 << " ";
1006  }
1007  edm::LogVerbatim("HcalSim") << st2.str();
1008  }
1009  std::ostringstream st2;
1010  for (int n = 0; n < 11 + (nBinsEta * nBinsPhi); n++)
1011  st2 << "=";
1012  edm::LogVerbatim("HcalSim") << st2.str();
1013 }
1015  if (SLShowerptr == nullptr) {
1016  edm::LogVerbatim("CastorShowerLibraryMaker::SLacceptEvent:") << "Error. nullptr pointer to CastorShowerEvent";
1017  return false;
1018  }
1019  if (ebin < 0 || ebin >= int(SLShowerptr->SLEnergyBins.size()))
1020  return false;
1021  if (SLisEBinFilled(ebin))
1022  return false;
1023 
1024  if (etabin < 0 || etabin >= int(SLShowerptr->SLEtaBins.size()))
1025  return false;
1027  return false;
1028 
1029  if (phibin < 0 || phibin >= int(SLShowerptr->SLPhiBins.size()))
1030  return false;
1032  return false;
1033  return true;
1034 }
1036  unsigned int volumeID = 0;
1037  double en_in_fi = 0.;
1038  //double totalEnergy = 0;
1039 
1040  int nentries = theCAFI->entries();
1041 
1042  // Compute Total Energy in CastorFI volume
1043  /*
1044  for(int ihit = 0; ihit < nentries; ihit++) {
1045  CaloG4Hit* aHit = (*theCAFI)[ihit];
1046  totalEnergy += aHit->getEnergyDeposit();
1047  }
1048 */
1049  if (!shower) {
1050  edm::LogVerbatim("CastorShowerLibraryMaker") << "Error. nullptr pointer to CastorShowerEvent";
1051  return false;
1052  }
1053 
1054  CastorNumberingScheme* theCastorNumScheme = new CastorNumberingScheme();
1055  // Hit position
1058  int nHits;
1059  nHits = 0;
1060  for (int ihit = 0; ihit < nentries; ihit++) {
1061  CaloG4Hit* aHit = (*theCAFI)[ihit];
1062  int hit_particleID = aHit->getTrackID();
1063  if (MapOfSecondaries[ipart].find(hit_particleID) == MapOfSecondaries[ipart].end()) {
1064  if (verbosity)
1065  edm::LogVerbatim("CastorShowerLibraryMaker") << "Skipping hit from trackID " << hit_particleID;
1066  continue;
1067  }
1068  volumeID = aHit->getUnitID();
1069  double hitEnergy = aHit->getEnergyDeposit();
1070  en_in_fi += aHit->getEnergyDeposit();
1071  float time = aHit->getTimeSlice();
1072  int zside, sector, zmodule;
1073  theCastorNumScheme->unpackIndex(volumeID, zside, sector, zmodule);
1074  entry = aHit->getEntry();
1075  position = aHit->getPosition();
1076  if (verbosity)
1077  edm::LogVerbatim("CastorShowerLibraryMaker") << "\n side , sector , module = " << zside << " , " << sector
1078  << " , " << zmodule << "\n nphotons = " << hitEnergy;
1079 
1080  if (verbosity)
1081  edm::LogVerbatim("CastorShowerLibraryMaker")
1082  << "\n packIndex = " << theCastorNumScheme->packIndex(zside, sector, zmodule);
1083 
1084  if (verbosity && time > 100.) {
1085  edm::LogVerbatim("CastorShowerLibraryMaker")
1086  << "\n nentries = " << nentries << "\n time[" << ihit << "] = " << time << "\n trackID[" << ihit
1087  << "] = " << aHit->getTrackID() << "\n volumeID[" << ihit << "] = " << volumeID << "\n nphotons[" << ihit
1088  << "] = " << hitEnergy << "\n side, sector, module = " << zside << ", " << sector << ", " << zmodule
1089  << "\n packIndex " << theCastorNumScheme->packIndex(zside, sector, zmodule) << "\n X,Y,Z = " << entry.x()
1090  << "," << entry.y() << "," << entry.z();
1091  }
1092  if (verbosity)
1093  edm::LogVerbatim("CastorShowerLibraryMaker") << "\n Incident Energy = " << aHit->getIncidentEnergy() << " \n";
1094 
1095  // CaloG4Hit information
1096  shower->setDetID(volumeID);
1097  shower->setHitPosition(position);
1098  shower->setNphotons(hitEnergy);
1099  shower->setTime(time);
1100  nHits++;
1101  }
1102  // Write number of hits to CastorShowerEvent instance
1103  if (nHits == 0) {
1104  edm::LogVerbatim("CastorShowerLibraryMaker") << "No hits found for this track (trackID=" << ipart << ").";
1105  if (theCastorNumScheme)
1106  delete theCastorNumScheme;
1107  return false;
1108  }
1109  shower->setNhit(nHits);
1110 
1111  // update the event counters
1112  if (theCastorNumScheme)
1113  delete theCastorNumScheme;
1114  return true;
1115 }
1117  if (!SLShowerptr) {
1118  edm::LogVerbatim("CastorShowerLibraryMaker") << "\n\nSLnEvtInBinE can be called only after BeginOfEvent\n\n";
1119  throw SimG4Exception("\n\nnullptr Pointer to the shower library.");
1120  }
1121  return SLShowerptr->nEvtInBinE.at(ebin);
1122 }
1123 
1125  if (!SLShowerptr) {
1126  edm::LogVerbatim("CastorShowerLibraryMaker") << "\n\nSLnEvtInBinEta can be called only after BeginOfEvent\n\n";
1127  throw SimG4Exception("\n\nnullptr Pointer to the shower library.");
1128  }
1129  return SLShowerptr->nEvtInBinEta.at(ebin).at(etabin);
1130 }
1131 
1133  if (!SLShowerptr) {
1134  edm::LogVerbatim("CastorShowerLibraryMaker") << "\n\nSLnEvtInBinPhi can be called only after BeginOfEvent\n\n";
1135  throw SimG4Exception("\n\nnullptr Pointer to the shower library.");
1136  }
1137  return SLShowerptr->nEvtInBinPhi.at(ebin).at(etabin).at(phibin);
1138 }
1140  if (!SLShowerptr) {
1141  edm::LogVerbatim("CastorShowerLibraryMaker") << "\n\nSLisEBinFilled can be called only after BeginOfEvent\n\n";
1142  throw SimG4Exception("\n\nnullptr Pointer to the shower library.");
1143  }
1145  return false;
1146  return true;
1147 }
1149  if (!SLShowerptr) {
1150  edm::LogVerbatim("CastorShowerLibraryMaker") << "\n\nSLisEtaBinFilled can be called only after BeginOfEvent\n\n";
1151  throw SimG4Exception("\n\nnullptr Pointer to the shower library.");
1152  }
1154  return false;
1155  return true;
1156 }
1158  if (!SLShowerptr) {
1159  edm::LogVerbatim("CastorShowerLibraryMaker") << "\n\nSLisPhiBinFilled can be called only after BeginOfEvent\n\n";
1160  throw SimG4Exception("\n\nnullptr Pointer to the shower library.");
1161  }
1163  return false;
1164  return true;
1165 }
1167  const G4TrackVector* p_sec = aStep->GetSecondary();
1168  for (int i = 0; i < int(p_sec->size()); i++) {
1169  edm::LogVerbatim("HcalSim") << "Killing track ID " << p_sec->at(i)->GetTrackID()
1170  << " and its secondaries Produced by Process "
1171  << p_sec->at(i)->GetCreatorProcess()->GetProcessName() << " in the volume "
1172  << aStep->GetTrack()->GetVolume()->GetName();
1173  p_sec->at(i)->SetTrackStatus(fKillTrackAndSecondaries);
1174  }
1175 }
1176 
1177 void CastorShowerLibraryMaker::GetMissingEnergy(CaloG4HitCollection* theCAFI, double& miss_energy, double& tot_energy) {
1178  // Get the total deposited energy and the one from hit not associated to a primary
1179  miss_energy = 0;
1180  tot_energy = 0;
1181  int nhits = theCAFI->entries();
1182  for (int ihit = 0; ihit < nhits; ihit++) {
1183  int id = (*theCAFI)[ihit]->getTrackID();
1184  tot_energy += (*theCAFI)[ihit]->getEnergyDeposit();
1185  int hit_prim = 0;
1186  for (unsigned int i = 0; i < thePrims.size(); i++) {
1187  int particleType = thePrims.at(i)->GetPDGcode();
1189  hit_prim = particleType;
1190  }
1191  if (hit_prim == 0) {
1192  edm::LogVerbatim("HcalSim") << "Track ID " << id << " produced a hit which is not associated with a primary.";
1193  miss_energy += (*theCAFI)[ihit]->getEnergyDeposit();
1194  }
1195  }
1196 }
1197 
1200 
void Clear(Option_t *option="") override
Log< level::Info, true > LogVerbatim
int getTrackID() const
Definition: CaloG4Hit.h:64
void setNBins(unsigned int n)
#define DEFINE_SIMWATCHER(type)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
void setTime(float t)
void setPrimEta(float eta)
void setPrimZ(float z)
void setDetID(unsigned int id)
int & SLnEvtInBinEta(int ebin, int etabin)
void setPrimPhi(float phi)
void setBin(double val)
std::vector< G4PrimaryParticle * > thePrims
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
bool FillShowerEvent(CaloG4HitCollection *, CastorShowerEvent *, int)
std::map< int, G4ThreeVector > PrimaryPosition
void KillSecondaries(const G4Step *step)
double npart
Definition: HydjetWrapper.h:48
int zside(DetId const &)
std::vector< std::vector< std::vector< std::vector< CastorShowerEvent > > > > SLBin3D
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
bool accept(const edm::Event &event, const edm::TriggerResults &triggerTable, const std::string &triggerPath)
Definition: TopDQMHelpers.h:31
bool SLisEtaBinFilled(int ebin, int etabin)
static void unpackIndex(const uint32_t &idx, int &z, int &sector, int &zmodule)
void setNEvtPerBin(unsigned int n)
std::map< int, std::set< int > > MapOfSecondaries
std::vector< G4PrimaryParticle * > GetPrimary(const G4Event *)
CastorShowerLibraryInfo * hadInfo
math::XYZPoint getPosition() const
Definition: CaloG4Hit.h:52
math::XYZPoint getEntry() const
Definition: CaloG4Hit.h:46
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< std::vector< std::vector< int > > > nEvtInBinPhi
CastorShowerLibraryInfo * emInfo
CastorShowerLibraryMaker(const edm::ParameterSet &p)
void setPrimE(float e)
void setPrimY(float y)
void GetMissingEnergy(CaloG4HitCollection *, double &, double &)
void GetKinematics(G4PrimaryParticle *, double &px, double &py, double &pz, double &pInit, double &eta, double &phi)
void setPrimX(float x)
int & SLnEvtInBinPhi(int ebin, int etabin, int phibin)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< std::vector< int > > nEvtInBinEta
double getEnergyDeposit() const
Definition: CaloG4Hit.h:79
uint32_t getUnitID() const
Definition: CaloG4Hit.h:66
void setNEvts(unsigned int n)
G4THitsCollection< CaloG4Hit > CaloG4HitCollection
static int position[264][3]
Definition: ReadPGInfo.cc:289
bool SLisPhiBinFilled(int ebin, int etabin, int phibin)
std::map< int, G4ThreeVector > PrimaryMomentum
double getTimeSlice() const
Definition: CaloG4Hit.h:67
step
Definition: StallMonitor.cc:83
void update(const BeginOfJob *run) override
This routine will be called when the appropriate signal arrives.
Log< level::Warning, false > LogWarning
static uint32_t packIndex(int z, int sector, int zmodule)
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
void setNphotons(float np)
double getIncidentEnergy() const
Definition: CaloG4Hit.h:61
Geom::Theta< T > theta() const
void setNhit(unsigned int i)
void setHitPosition(const Point &p)
unsigned int getNhit()