CMS 3D CMS Logo

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