CMS 3D CMS Logo

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