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  std::string pName = fp.fullPath();
10  if (pName.find(".") == 0)
11  pName.erase(0, 2);
12  theDataFile = pName;
13  readUserData();
14  produces<HFShowerPhotonCollection> ("emParticles");
15  produces<HFShowerPhotonCollection> ("hadParticles");
16  produces<HFShowerLibraryEventInfo> ("EventInfo");
17 
18 }
19 
21 
22 }
23 
25  std::auto_ptr<HFShowerPhotonCollection> product_em(new HFShowerPhotonCollection);
26  std::auto_ptr<HFShowerPhotonCollection> product_had(new HFShowerPhotonCollection);
27 
28  //EventInfo
29  std::auto_ptr<HFShowerLibraryEventInfo> product_evtInfo(new HFShowerLibraryEventInfo);
30  float hfShowerLibV = 1.1;
31  float phyListV = 3.6;
32  std::vector<double> en;
33  double momBin[12] = { 10., 15., 20., 35., 50., 80., 100., 150., 250., 350., 500., 1000. };
34  for (int i = 0; i < 12; ++i)
35  en.push_back(momBin[i]);
36  HFShowerLibraryEventInfo evtInfo(60000, 12, 5000, hfShowerLibV, phyListV, en);
37  *product_evtInfo = evtInfo;
38  iEvent.put(product_evtInfo, "EventInfo");
39 
40 
41  //shower photons
44 
45  int n = theFileHandle.size();
46  for (int i = 0; i < n; ++i) {
47  std::string fn = theFileHandle[i].name;
48  std::string particle = theFileHandle[i].id;
49  // int momBin = theFileHandle[i].momentum;
50  TFile* theFile = new TFile(fn.c_str(), "READ");
51  TTree* theTree = (TTree*) theFile->FindObjectAny("CherenkovPhotons");
52  int nphot = 0;
53  float x[10000];
54  float y[10000];
55  float z[10000];
56  float t[10000];
57  float lambda[10000];
58  int fiberId[10000];
59  for (int kk = 0; kk < 10000; ++kk) {
60  x[kk] = 0.;
61  y[kk] = 0.;
62  z[kk] = 0.;
63  t[kk] = 0.;
64  lambda[kk] = 0.;
65  fiberId[kk] = 0;
66  }
67  theTree->SetBranchAddress("nphot", &nphot);
68  theTree->SetBranchAddress("x", &x);
69  theTree->SetBranchAddress("y", &y);
70  theTree->SetBranchAddress("z", &z);
71  theTree->SetBranchAddress("t", &t);
72  theTree->SetBranchAddress("lambda", &lambda);
73  theTree->SetBranchAddress("fiberId", &fiberId);
74  int nentry = int(theTree->GetEntries());
75  if (particle == "electron") {
76  for (int iev = 0; iev < nentry; iev++) {
77  for (int iph = 0; iph < nphot; ++iph) {
78  HFShowerPhoton::Point pos(x[iph], y[iph], z[iph]);
79  HFShowerPhoton aPhoton(pos, t[iph], lambda[iph]);
80  emColl.push_back(aPhoton);
81  }
82  }
83  }
84  if (particle == "pion") {
85  }
86  theFile->Close();
87  if (theFile)
88  delete theFile;
89  if (theTree)
90  delete theTree;
91  }
92  *product_em = emColl;
93  *product_had = hadColl;
94 
95  //fillEvent(product_em,product_had);
96  //fillEvent(product_em, product_had);
97  iEvent.put(product_em, "emParticles");
98  iEvent.put(product_had, "hadParticles");
99 }
100 
102 }
103 //void HcalForwardLibWriter::fillEvent(HFShowerPhotonCollection& em, HFShowerPhotonCollection& had) {
104 //
105 // HFShowerPhotonCollection emColl;
106 // HFShowerPhotonCollection hadColl;
107 //
108 // int n = theFileHandle.size();
109 // for (int i = 0; i < n; ++i) {
110 // std::string fn = theFileHandle[i].name;
111 // std::string particle = theFileHandle[i].id;
112 // int momBin = theFileHandle[i].momentum;
113 // TFile* theFile = new TFile(fn.c_str(), "READ");
114 // TTree* theTree = (TTree*) theFile->FindObjectAny("CherenkovPhotons");
115 // int nphot = 0;
116 // float x[10000];
117 // float y[10000];
118 // float z[10000];
119 // float t[10000];
120 // float lambda[10000];
121 // int fiberId[10000];
122 // for (int kk = 0; kk < 10000; ++kk) {
123 // x[kk] = 0.;
124 // y[kk] = 0.;
125 // z[kk] = 0.;
126 // t[kk] = 0.;
127 // lambda[kk] = 0.;
128 // fiberId[kk] = 0;
129 // }
130 // theTree->SetBranchAddress("nphot", &nphot);
131 // theTree->SetBranchAddress("x", &x);
132 // theTree->SetBranchAddress("y", &y);
133 // theTree->SetBranchAddress("z", &z);
134 // theTree->SetBranchAddress("t", &t);
135 // theTree->SetBranchAddress("lambda", &lambda);
136 // theTree->SetBranchAddress("fiberId", &fiberId);
137 // int nentry = int(theTree->GetEntries());
138 // if (particle == "electron") {
139 // for (int iev = 0; iev < nentry; iev++) {
140 // for (int iph = 0; iph < nphot; ++iph) {
141 // HFShowerPhoton::Point pos(x[iph], y[iph], z[iph]);
142 // HFShowerPhoton aPhoton(pos, t[iph], lambda[iph]);
143 // emColl.push_back(aPhoton);
144 // }
145 // }
146 // }
147 // if (particle == "pion") {
148 // }
149 // theFile->Close();
150 // if (theFile)
151 // delete theFile;
152 // if (theTree)
153 // delete theTree;
154 // }
155 // *em = emColl;
156 // *had = hadColl;
157 //}
158 
160 }
162  std::ifstream input(theDataFile.c_str());
163  if (input.fail()) {
164  return 0;
165  }
166  std::string theFileName, thePID;
167  int mom;
168  int k = 0;
169  while (!input.eof()) {
170  input >> theFileName >> thePID >> mom;
171  if (!input.fail()) {
172  FileHandle aFile;
173  aFile.name = theFileName;
174  aFile.id = thePID;
175  aFile.momentum = mom;
176  theFileHandle.push_back(aFile);
177  ++k;
178  } else {
179  input.clear();
180  }
181  input.ignore(999, '\n');
182  }
183  return k;
184 }
185 //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
double double double z
int iEvent
Definition: GenABIO.cc:243
std::vector< HFShowerPhoton > HFShowerPhotonCollection
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
int k[5][pyjets_maxn]
HcalForwardLibWriter(const edm::ParameterSet &)
std::vector< FileHandle > theFileHandle
virtual void produce(edm::Event &, const edm::EventSetup &)
std::string fullPath() const
Definition: FileInPath.cc:171
x
Definition: VDTMath.h:216