CMS 3D CMS Logo

EcalPnGraphs.cc
Go to the documentation of this file.
1 
14 
17 
20 
21 #include <iostream>
22 #include <vector>
23 #include <map>
24 
25 #include "TFile.h"
26 #include "TGraph.h"
27 
28 //=============================================================================
30  //=============================================================================
31 
32  digiProducer_ = ps.getParameter<std::string>("digiProducer");
33  fileName = ps.getUntrackedParameter<std::string>("fileName", std::string("toto"));
34 
35  first_Pn = 0;
36 
37  listPns = ps.getUntrackedParameter<std::vector<int> >("listPns", std::vector<int>());
38  numPn = ps.getUntrackedParameter<int>("numPn");
39 
40  std::vector<int> listDefaults;
41  listDefaults.push_back(-1);
42  feds_ = ps.getUntrackedParameter<std::vector<int> >("requestedFeds", listDefaults);
43  bool fedIsGiven = false;
44 
45  std::vector<std::string> ebDefaults;
46  ebDefaults.push_back("none");
47  ebs_ = ps.getUntrackedParameter<std::vector<std::string> >("requestedEbs", ebDefaults);
48 
49  //FEDs and EBs
50  if (feds_[0] != -1) {
51  edm::LogInfo("EcalPnGraphs") << "FED id is given! Goining to beginJob! ";
52  fedIsGiven = true;
53  } else {
54  feds_.clear();
55  if (ebs_[0] != "none") {
56  //EB id is given and convert to FED id
57  fedMap = new EcalFedMap();
58  for (std::vector<std::string>::const_iterator ebItr = ebs_.begin(); ebItr != ebs_.end(); ++ebItr) {
59  feds_.push_back(fedMap->getFedFromSlice(*ebItr));
60  }
61  delete fedMap;
62  } else {
63  //Select all FEDs in the Event
64  for (int i = 601; i < 655; ++i) {
65  feds_.push_back(i);
66  }
67  }
68  }
69 
70  // consistency checks checks
71  inputIsOk = true;
72  //check with FEDs
73  if (fedIsGiven) {
74  std::vector<int>::iterator fedIter;
75  for (fedIter = feds_.begin(); fedIter != feds_.end(); ++fedIter) {
76  if ((*fedIter) < 601 || (*fedIter) > 654) {
77  std::cout << "[EcalPnGraphs] pn number : " << (*fedIter) << " found in listFeds. "
78  << " Valid range is [601-654]. Returning." << std::endl;
79  inputIsOk = false;
80  return;
81  }
82  }
83  }
84 
85  //Check with Pns
86  if (listPns[0] != -1) {
87  std::vector<int>::iterator intIter;
88  for (intIter = listPns.begin(); intIter != listPns.end(); intIter++) {
89  if (((*intIter) < 1) || (10 < (*intIter))) {
90  std::cout << "[EcalPnGraphs] pn number : " << (*intIter) << " found in listPns. "
91  << " Valid range is 1-10. Returning." << std::endl;
92  inputIsOk = false;
93  return;
94  }
95  if (!first_Pn)
96  first_Pn = (*intIter);
97  }
98  } else {
99  listPns.clear();
100  listPns.push_back(5);
101  listPns.push_back(6);
102  }
103 
104  // setting the abcissa array once for all
105  for (int i = 0; i < 50; i++)
106  abscissa[i] = i;
107 
108  // local event counter (in general different from LV1)
109  eventCounter = 0;
110 }
111 
112 //=============================================================================
114  //=============================================================================
115  //delete *;
116 }
117 
118 //=============================================================================
120  //=============================================================================
121  edm::LogInfo("EcalPhGraphs") << "entering beginJob! ";
122 }
123 
124 //=============================================================================
126  //=============================================================================
127 
128  eventCounter++;
129  if (!inputIsOk)
130  return;
131 
132  // retrieving crystal PN diodes from Event
134  try {
135  e.getByLabel(digiProducer_, pn_digis);
136  } catch (cms::Exception& ex) {
137  edm::LogError("EcalPnGraphs") << "PNs were not found!";
138  }
139 
140  // getting the list of all the Pns which will be dumped on TGraph
141  // - listPns is the list as given by the user
142  // -numPn is the number of Pns (centered at Pn from listPns) for which graphs are required
143  std::vector<int>::iterator pn_it;
144  for (pn_it = listPns.begin(); pn_it != listPns.end(); pn_it++) {
145  int ipn = (*pn_it);
146  int hpn = numPn;
147 
148  for (int u = (-hpn); u <= hpn; u++) {
149  int ipn_c = ipn + u;
150  if (ipn_c < 1 || ipn_c > 10)
151  continue;
152  std::vector<int>::iterator notInList = find(listAllPns.begin(), listAllPns.end(), ipn_c);
153  if (notInList == listAllPns.end()) {
154  listAllPns.push_back(ipn_c);
155  }
156  }
157  }
158 
159  //Loop over PN digis
160  for (EcalPnDiodeDigiCollection::const_iterator pnItr = pn_digis->begin(); pnItr != pn_digis->end(); ++pnItr) {
161  //Get PNid of a digi
162  int ipn = (*pnItr).id().iPnId();
163  //Get DCC id where the digi is from
164  int ieb = EcalPnDiodeDetId((*pnItr).id()).iDCCId();
165 
166  //Make sure that these are PnDigis from the requested FEDid
167  int FEDid = ieb + 600;
168 
169  std::vector<int>::iterator fedIter = find(feds_.begin(), feds_.end(), FEDid);
170 
171  if (fedIter == feds_.end()) {
172  edm::LogWarning("EcalPnGraphs") << "For Event " << eventCounter
173  << " PnDigis are not found from requested SM!. returning...";
174  return;
175  }
176  // selecting desired Pns only
177  std::vector<int>::iterator iPnIter;
178  iPnIter = find(listAllPns.begin(), listAllPns.end(), ipn);
179  if (iPnIter == listAllPns.end())
180  continue;
181 
182  for (int i = 0; i < (*pnItr).size() && i < 50; ++i) {
183  ordinate[i] = (*pnItr).sample(i).adc();
184  }
185  //make grapn of ph digis
186  TGraph oneGraph(50, abscissa, ordinate);
188  title = "Graph_ev" + intToString(eventCounter) + "_FED" + intToString(FEDid) + "_ipn" + intToString(ipn);
189  oneGraph.SetTitle(title.c_str());
190  oneGraph.SetName(title.c_str());
191  graphs.push_back(oneGraph);
192 
193  } // loop over Pn digis
194 }
195 
197  //
198  // outputs the number into the string stream and then flushes
199  // the buffer (makes sure the output is put into the stream)
200  //
201  std::ostringstream myStream; //creates an ostringstream object
202  myStream << num << std::flush;
203 
204  return (myStream.str()); //returns the string form of the stringstream object
205 }
206 
207 //=============================================================================
209  //=============================================================================
210  fileName += (std::string("_Pn") + intToString(first_Pn));
211  fileName += ".graph.root";
212 
213  root_file = new TFile(fileName.c_str(), "RECREATE");
214  std::vector<TGraph>::iterator gr_it;
215  for (gr_it = graphs.begin(); gr_it != graphs.end(); gr_it++)
216  (*gr_it).Write();
217  root_file->Close();
218 
219  edm::LogInfo("EcalPnGraphs") << "DONE!.... ";
220 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
int abscissa[50]
Definition: EcalPnGraphs.h:61
int getFedFromSlice(std::string)
Definition: EcalFedMap.cc:96
std::vector< T >::const_iterator const_iterator
std::string intToString(int num)
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
std::vector< int > feds_
Definition: EcalPnGraphs.h:41
T getUntrackedParameter(std::string const &, T const &) const
std::string fileName
Definition: EcalPnGraphs.h:52
void analyze(const edm::Event &e, const edm::EventSetup &c) override
void beginJob() override
EcalFedMap * fedMap
Definition: EcalPnGraphs.h:34
const_iterator begin() const
std::vector< TGraph > graphs
Definition: EcalPnGraphs.h:64
int ordinate[50]
Definition: EcalPnGraphs.h:62
const_iterator end() const
Log< level::Info, false > LogInfo
std::string digiProducer_
Definition: EcalPnGraphs.h:39
std::vector< int > listAllPns
Definition: EcalPnGraphs.h:57
~EcalPnGraphs() override
std::vector< std::string > ebs_
Definition: EcalPnGraphs.h:42
std::vector< int > listPns
Definition: EcalPnGraphs.h:56
Log< level::Warning, false > LogWarning
TFile * root_file
Definition: EcalPnGraphs.h:66
void endJob() override
EcalPnGraphs(const edm::ParameterSet &ps)
Definition: EcalPnGraphs.cc:29