CMS 3D CMS Logo

HcalForwardLibWriter.cc
Go to the documentation of this file.
3 
5  edm::ParameterSet theParms = iConfig.getParameter<edm::ParameterSet>("hcalForwardLibWriterParameters");
6  edm::FileInPath fp = theParms.getParameter<edm::FileInPath>("FileName");
7  nbins = theParms.getParameter<int>("Nbins");
8  nshowers = theParms.getParameter<int>("Nshowers");
9  bsize = theParms.getParameter<int>("BufSize");
10  splitlevel = theParms.getParameter<int>("SplitLevel");
11  compressionAlgo = theParms.getParameter<int>("CompressionAlgo");
12  compressionLevel = theParms.getParameter<int>("CompressionLevel");
13 
14  std::string pName = fp.fullPath();
15  if (pName.find('.') == 0)
16  pName.erase(0, 2);
17  theDataFile = pName;
18  readUserData();
19 
21  fs->file().cd();
22  fs->file().SetCompressionAlgorithm(compressionAlgo);
23  fs->file().SetCompressionLevel(compressionLevel);
24 
25  LibTree = new TTree("HFSimHits", "HFSimHits");
26 
27  //https://root.cern/root/html534/TTree.html
28  // TBranch*Branch(const char* name, const char* classname, void** obj, Int_t bufsize = 32000, Int_t splitlevel = 99)
29  emBranch = LibTree->Branch("emParticles", "HFShowerPhotons-emParticles", &emColl, bsize, splitlevel);
30  hadBranch = LibTree->Branch("hadParticles", "HFShowerPhotons-hadParticles", &hadColl, bsize, splitlevel);
31 }
32 
35  desc.add<edm::FileInPath>("FileName", edm::FileInPath("SimG4CMS/ShowerLibraryProducer/data/fileList.txt"));
36  desc.add<int>("Nbins", 16);
37  desc.add<int>("Nshowers", 10000);
38  desc.add<int>("BufSize", 1);
39  desc.add<int>("SplitLevel", 2);
40  desc.add<int>("CompressionAlgo", 4);
41  desc.add<int>("CompressionLevel", 4);
42  descriptions.add("hcalForwardLibWriterParameters", desc);
43 }
44 
46  // Event info
47  std::vector<double> en;
48  double momBin[16] = {2, 3, 5, 7, 10, 15, 20, 30, 50, 75, 100, 150, 250, 350, 500, 1000};
49  en.reserve(nbins);
50  for (int i = 0; i < nbins; ++i)
51  en.push_back(momBin[i]);
52 
53  int nem = 0;
54  int nhad = 0;
55 
56  //shower photons
57  int n = theFileHandle.size();
58  // cycle over files ++++++++++++++++++++++++++++++++++++++++++++++++++++
59  for (int i = 0; i < n; ++i) {
60  std::string fn = theFileHandle[i].name;
61  std::string particle = theFileHandle[i].id;
62 
63  // edm::LogVerbatim("HcalSim") << "*** Input file " << i << " " << fn;
64 
65  TFile* theFile = new TFile(fn.c_str(), "READ");
66  TTree* theTree = (TTree*)gDirectory->Get("g4SimHits/CherenkovPhotons");
67  int nphot = 0;
68 
69  const int size = 10000;
70  if (nshowers > size) {
71  edm::LogError("HcalForwardLibWriter") << "Too big Nshowers number";
72  return;
73  }
74 
75  float x[size] = {0.};
76  float y[size] = {0.};
77  float z[size] = {0.};
78  float t[size] = {0.};
79  float lambda[size] = {0.};
80  int fiberId[size] = {0};
81  float primZ; // added
82 
83  theTree->SetBranchAddress("nphot", &nphot);
84  theTree->SetBranchAddress("x", &x);
85  theTree->SetBranchAddress("y", &y);
86  theTree->SetBranchAddress("z", &z);
87  theTree->SetBranchAddress("t", &t);
88  theTree->SetBranchAddress("lambda", &lambda);
89  theTree->SetBranchAddress("fiberId", &fiberId);
90  theTree->SetBranchAddress("primZ", &primZ); // added
91  int nentries = int(theTree->GetEntries());
92  int ngood = 0;
93  // cycle over showers ====================================================
94  for (int iev = 0; iev < nentries; iev++) {
95  if (primZ < 990.)
96  continue; // exclude showers with interactions in front of HF (1m of air)
97  ngood++;
98  if (ngood > nshowers)
99  continue;
100  if (particle == "electron") {
101  emColl.clear();
102  } else {
103  hadColl.clear();
104  }
105  float nphot_long = 0;
106  float nphot_short = 0;
107  // cycle over photons in shower -------------------------------------------
108  for (int iph = 0; iph < nphot; ++iph) {
109  if (fiberId[iph] == 1) {
110  nphot_long++;
111  } else {
112  nphot_short++;
113  z[iph] = -z[iph];
114  }
115 
116  HFShowerPhoton::Point pos(x[iph], y[iph], z[iph]);
117  HFShowerPhoton aPhoton(pos, t[iph], lambda[iph]);
118  if (particle == "electron") {
119  emColl.push_back(aPhoton);
120  } else {
121  hadColl.push_back(aPhoton);
122  }
123  }
124  // end of cycle over photons in shower -------------------------------------------
125 
126  if (particle == "electron") {
127  LibTree->SetEntries(nem + 1);
128  emBranch->Fill();
129  nem++;
130  emColl.clear();
131  } else {
132  LibTree->SetEntries(nhad + 1);
133  nhad++;
134  hadBranch->Fill();
135  hadColl.clear();
136  }
137  }
138  // end of cycle over showers ====================================================
139  theFile->Close();
140  }
141  // end of cycle over files ++++++++++++++++++++++++++++++++++++++++++++++++++++
142 }
143 
145 
148  fs->file().cd();
149  LibTree->Write();
150  LibTree->Print();
151 }
152 
154  std::ifstream input(theDataFile.c_str());
155  if (input.fail()) {
156  return 0;
157  }
158  std::string theFileName, thePID;
159  int mom;
160  int k = 0;
161  while (!input.eof()) {
162  input >> theFileName >> thePID >> mom;
163  if (!input.fail()) {
164  FileHandle aFile;
165  aFile.name = theFileName;
166  aFile.id = thePID;
167  aFile.momentum = mom;
168  theFileHandle.push_back(aFile);
169  ++k;
170  } else {
171  input.clear();
172  }
173  input.ignore(999, '\n');
174  }
175  return k;
176 }
177 
178 //define this as a plug-in
size
Write out results.
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
math::XYZPointF Point
point in the space
TkSoAView< TrackerTraits > HitToTuple< TrackerTraits > const *__restrict__ int32_t int32_t int iev
Log< level::Error, false > LogError
static std::string const input
Definition: EdmProvDump.cc:50
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
HcalForwardLibWriter(const edm::ParameterSet &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void analyze(const edm::Event &, const edm::EventSetup &) override
std::vector< FileHandle > theFileHandle
HFShowerPhotonCollection hadColl
HFShowerPhotonCollection emColl