test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 
12  std::string pName = fp.fullPath();
13  if (pName.find(".") == 0)
14  pName.erase(0, 2);
15  theDataFile = pName;
16  readUserData();
17 
18  int bsize = 64000;
19  fs->file().cd();
20  LibTree = new TTree("HFSimHits", "HFSimHits");
21  LibTree->Branch("emParticles", "HFShowerPhotons-emParticles", &emColl, bsize);
22  LibTree->Branch("hadParticles", "HFShowerPhotons-hadParticles", &hadColl, bsize);
23 }
24 
26 }
27 
29 
30 // Event info
31  std::vector<double> en;
32  double momBin[16] = {2,3,5,7,10,15,20,30,50,75,100,150,250,350,500,1000};
33  for (int i = 0; i < nbins; ++i) en.push_back(momBin[i]);
34 
35  //shower photons
36  int n = theFileHandle.size();
37 // cycle over files ++++++++++++++++++++++++++++++++++++++++++++++++++++
38  for (int i = 0; i < n; ++i) {
39  std::string fn = theFileHandle[i].name;
40  std::string particle = theFileHandle[i].id;
41  TFile* theFile = new TFile(fn.c_str(), "READ");
42  TTree *theTree = (TTree*)gDirectory->Get("g4SimHits/CherenkovPhotons");
43  int nphot = 0;
44  float x[10000];
45  float y[10000];
46  float z[10000];
47  float t[10000];
48  float lambda[10000];
49  int fiberId[10000];
50  for (int kk = 0; kk < 10000; ++kk) {
51  x[kk] = 0.;
52  y[kk] = 0.;
53  z[kk] = 0.;
54  t[kk] = 0.;
55  lambda[kk] = 0.;
56  fiberId[kk] = 0;
57  }
58  theTree->SetBranchAddress("nphot", &nphot);
59  theTree->SetBranchAddress("x", &x);
60  theTree->SetBranchAddress("y", &y);
61  theTree->SetBranchAddress("z", &z);
62  theTree->SetBranchAddress("t", &t);
63  theTree->SetBranchAddress("lambda", &lambda);
64  theTree->SetBranchAddress("fiberId", &fiberId);
65  int nentries = int(theTree->GetEntries());
66  if(nentries>5000) nentries=5000;
67  int nbytes = 0;
68 // cycle over showers ====================================================
69  for (int iev = 0; iev < nentries; iev++) {
70  nbytes += theTree->GetEntry(iev);
71  if (particle == "electron") {
72  emColl.clear();
73  } else {
74  hadColl.clear();
75  }
76  float nphot_long=0;
77  float nphot_short=0;
78 // cycle over photons in shower -------------------------------------------
79  for (int iph = 0; iph < nphot; ++iph) {
80 
81  if (fiberId[iph]==1) {
82  nphot_long++;
83  } else {
84  nphot_short++;
85  z[iph]= -z[iph];
86  }
87 
88  HFShowerPhoton::Point pos(x[iph], y[iph], z[iph]);
89  HFShowerPhoton aPhoton(pos, t[iph], lambda[iph]);
90  if (particle == "electron") {
91  emColl.push_back(aPhoton);
92  } else {
93  hadColl.push_back(aPhoton);
94  }
95  }
96 // end of cycle over photons in shower -------------------------------------------
97 
98 // if(iev>0) LibTree->SetBranchStatus("HFShowerLibInfo",0);
99  if (particle == "electron") {
100  LibTree->SetBranchStatus("hadParticles",0);
101  } else {
102  LibTree->SetBranchStatus("emParticles",0);
103  }
104  LibTree->Fill();
105  if (particle == "electron") {
106  emColl.clear();
107  } else {
108  hadColl.clear();
109  }
110  }
111 // end of cycle over showers ====================================================
112  theFile->Close();
113  }
114 // end of cycle over files ++++++++++++++++++++++++++++++++++++++++++++++++++++
115 }
116 
118 }
119 
121  fs->file().cd();
122  LibTree->Write();
123  LibTree->Print();
124 }
125 
127  std::ifstream input(theDataFile.c_str());
128  if (input.fail()) {
129  return 0;
130  }
131  std::string theFileName, thePID;
132  int mom;
133  int k = 0;
134  while (!input.eof()) {
135  input >> theFileName >> thePID >> mom;
136  if (!input.fail()) {
137  FileHandle aFile;
138  aFile.name = theFileName;
139  aFile.id = thePID;
140  aFile.momentum = mom;
141  theFileHandle.push_back(aFile);
142  ++k;
143  } else {
144  input.clear();
145  }
146  input.ignore(999, '\n');
147  }
148  return k;
149 }
150 //define this as a plug-in
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
math::XYZPointF Point
point in the space
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::Service< TFileService > fs
static std::string const input
Definition: EdmProvDump.cc:44
int iEvent
Definition: GenABIO.cc:230
virtual void analyze(const edm::Event &, const edm::EventSetup &)
TFile & file() const
return opened TFile
Definition: TFileService.h:37
HcalForwardLibWriter(const edm::ParameterSet &)
std::vector< FileHandle > theFileHandle
HFShowerPhotonCollection hadColl
HFShowerPhotonCollection emColl