CMS 3D CMS Logo

CheckSecondary.cc
Go to the documentation of this file.
4 
7 
11 
13 
14 #include "CLHEP/Units/GlobalPhysicalConstants.h"
15 #include "CLHEP/Units/GlobalSystemOfUnits.h"
16 #include "G4HCofThisEvent.hh"
17 #include "G4Step.hh"
18 #include "G4Track.hh"
19 #include "G4VProcess.hh"
20 
21 #include <cmath>
22 #include <iomanip>
23 #include <iostream>
24 
26  : treatSecondary(nullptr),
27  typeEnumerator(nullptr),
28  nsec(nullptr),
29  procids(nullptr),
30  px(nullptr),
31  py(nullptr),
32  pz(nullptr),
33  mass(nullptr),
34  deltae(nullptr),
35  procs(nullptr),
36  file(nullptr),
37  tree(nullptr) {
38  edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("CheckSecondary");
39  std::string saveFile = m_p.getUntrackedParameter<std::string>("SaveInFile", "None");
40  treatSecondary = new TreatSecondary(m_p);
42 
43  nsec = new std::vector<int>();
44  px = new std::vector<double>();
45  py = new std::vector<double>();
46  pz = new std::vector<double>();
47  mass = new std::vector<double>();
48  deltae = new std::vector<double>();
49  procids = new std::vector<int>();
50  procs = new std::vector<std::string>();
51 
52  if (saveFile != "None") {
53  saveToTree = true;
55  edm::LogInfo("CheckSecondary") << "Instantiate CheckSecondary with first"
56  << " hadronic interaction information"
57  << " to be saved in file " << saveFile;
58  } else {
59  saveToTree = false;
60  edm::LogInfo("CheckSecondary") << "Instantiate CheckSecondary with first"
61  << " hadronic interaction information"
62  << " not saved";
63  }
64 }
65 
67  if (saveToTree)
68  endTree();
69  if (nsec)
70  delete nsec;
71  if (px)
72  delete px;
73  if (py)
74  delete py;
75  if (pz)
76  delete pz;
77  if (mass)
78  delete mass;
79  if (deltae)
80  delete deltae;
81  if (procs)
82  delete procs;
83  if (procids)
84  delete procids;
85  if (typeEnumerator)
86  delete typeEnumerator;
87  if (treatSecondary)
88  delete treatSecondary;
89 }
90 
92  file = new TFile(fileName.c_str(), "RECREATE");
93  file->cd();
94 
95  TTree *t1 = new TTree("T1", "Secondary Particle Information");
96  t1->Branch("SecondaryPx", "std::vector<double>", &px);
97  t1->Branch("SecondaryPy", "std::vector<double>", &py);
98  t1->Branch("SecondaryPz", "std::vector<double>", &pz);
99  t1->Branch("SecondaryMass", "std::vector<double>", &mass);
100  t1->Branch("NumberSecondaries", "std::vector<int>", &nsec);
101  t1->Branch("DeltaEinInteract", "std::vector<double>", &deltae);
102  t1->Branch("ProcessID", "std::vector<int>", &procids);
103  t1->Branch("ProcessNames", "std::vector<std::string>", &procs);
104  return t1;
105 }
106 
108  edm::LogInfo("CheckSecondary") << "Save the Secondary Tree " << tree->GetName() << " (" << tree << ") in file "
109  << file->GetName() << " (" << file << ")";
110  file->cd();
111  tree->Write();
112  file->Close();
113  delete file;
114 }
115 
117  int iev = (*evt)()->GetEventID();
118  LogDebug("CheckSecondary") << "CheckSecondary::=====> Begin of event = " << iev;
119 
120  (*nsec).clear();
121  (*procs).clear();
122  (*procids).clear();
123  (*deltae).clear();
124  (*px).clear();
125  (*py).clear();
126  (*pz).clear();
127  (*mass).clear();
128 }
129 
131  const G4Track *thTk = (*trk)();
132  treatSecondary->initTrack(thTk);
133  if (thTk->GetParentID() <= 0)
134  storeIt = true;
135  else
136  storeIt = false;
137  nHad = 0;
138  LogDebug("CheckSecondary") << "CheckSecondary:: Track " << thTk->GetTrackID() << " Parent " << thTk->GetParentID()
139  << " Flag " << storeIt;
140 }
141 
142 void CheckSecondary::update(const G4Step *aStep) {
144  int procID;
145  bool hadrInt;
146  double deltaE;
147  std::vector<int> charge;
148  std::vector<math::XYZTLorentzVector> tracks = treatSecondary->tracks(aStep, name, procID, hadrInt, deltaE, charge);
149  if (storeIt && hadrInt) {
150  double pInit = (aStep->GetPreStepPoint()->GetMomentum()).mag();
151  double pEnd = (aStep->GetPostStepPoint()->GetMomentum()).mag();
152  nHad++;
153  int sec = (int)(tracks.size());
154  LogDebug("CheckSecondary") << "CheckSecondary:: Hadronic Interaction " << nHad << " of type " << name << " ID "
155  << procID << " with " << sec << " secondaries "
156  << " and Momentum (MeV/c) at start " << pInit << " and at end " << pEnd;
157  (*nsec).push_back(sec);
158  (*procs).push_back(name);
159  (*procids).push_back(procID);
160  (*deltae).push_back(deltaE);
161  if (nHad == 1) {
162  for (int i = 0; i < sec; i++) {
163  double m = tracks[i].M();
164  if (charge[i] < 0)
165  m = -m;
166  (*px).push_back(tracks[i].Px());
167  (*py).push_back(tracks[i].Py());
168  (*pz).push_back(tracks[i].Pz());
169  (*mass).push_back(m);
170  }
171  }
172  }
173 }
174 
176  LogDebug("CheckSecondary") << "CheckSecondary::EndofEvent =====> Event " << (*evt)()->GetEventID() << " with "
177  << (*nsec).size() << " hadronic collisions with"
178  << " secondaries produced in each step";
179  for (unsigned int i = 0; i < (*nsec).size(); i++)
180  LogDebug("CheckSecondary") << " " << (*nsec)[i] << " from " << (*procs)[i] << " ID " << (*procids)[i] << " ("
181  << typeEnumerator->processG4Name((*procids)[i]) << ") deltaE = " << (*deltae)[i]
182  << " MeV";
183  LogDebug("CheckSecondary") << "And " << (*mass).size() << " secondaries "
184  << "produced in the first interactions";
185  for (unsigned int i = 0; i < (*mass).size(); i++)
186  LogDebug("CheckSecondary") << "Secondary " << i << " (" << (*px)[i] << ", " << (*py)[i] << ", " << (*pz)[i] << ", "
187  << (*mass)[i] << ")";
188 
189  if (saveToTree)
190  tree->Fill();
191 }
std::vector< std::string > * procs
~CheckSecondary() override
std::vector< double > * pz
std::vector< double > * deltae
TkSoAView< TrackerTraits > HitToTuple< TrackerTraits > const *__restrict__ int32_t int32_t int iev
std::vector< math::XYZTLorentzVector > tracks(const G4Step *step, std::string &procName, int &procID, bool &intr, double &deltaE, std::vector< int > &charges)
T getUntrackedParameter(std::string const &, T const &) const
std::vector< int > * nsec
Log< level::Info, false > LogInfo
G4ProcessTypeEnumerator * typeEnumerator
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
TTree * bookTree(std::string)
void initTrack(const G4Track *trk)
std::string processG4Name(int) const
std::vector< double > * py
void update(const BeginOfEvent *evt) override
This routine will be called when the appropriate signal arrives.
CheckSecondary(const edm::ParameterSet &p)
std::vector< int > * procids
Definition: tree.py:1
std::vector< double > * px
TreatSecondary * treatSecondary
std::vector< double > * mass
#define LogDebug(id)