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 
30 //=============================================================================
32 //=============================================================================
33 
34  digiProducer_ = ps.getParameter<std::string>("digiProducer");
35  fileName = ps.getUntrackedParameter<std::string >("fileName", std::string("toto") );
36 
37  first_Pn = 0;
38 
39  listPns = ps.getUntrackedParameter<std::vector<int> >("listPns", std::vector<int>());
40  numPn = ps.getUntrackedParameter< int >("numPn");
41 
42  std::vector<int> listDefaults;
43  listDefaults.push_back(-1);
44  feds_ = ps.getUntrackedParameter<std::vector<int> > ("requestedFeds",listDefaults);
45  bool fedIsGiven = false;
46 
47  std::vector<std::string> ebDefaults;
48  ebDefaults.push_back("none");
49  ebs_ = ps.getUntrackedParameter<std::vector<std::string> >("requestedEbs",ebDefaults);
50 
51  //FEDs and EBs
52  if ( feds_[0] != -1 ) {
53  edm::LogInfo("EcalPnGraphs") << "FED id is given! Goining to beginJob! ";
54  fedIsGiven = true;
55  }else {
56  feds_.clear();
57  if ( ebs_[0] !="none" ) {
58  //EB id is given and convert to FED id
59  fedMap = new EcalFedMap();
60  for (std::vector<std::string>::const_iterator ebItr = ebs_.begin();
61  ebItr!= ebs_.end(); ++ebItr) {
62  feds_.push_back(fedMap->getFedFromSlice(*ebItr));
63  }
64  delete fedMap;
65  } else {
66  //Select all FEDs in the Event
67  for ( int i=601; i<655; ++i){
68  feds_.push_back(i);
69  }
70  }
71  }
72 
73  // consistency checks checks
74  inputIsOk = true;
75  //check with FEDs
76  if ( fedIsGiven ) {
77  std::vector<int>::iterator fedIter;
78  for (fedIter = feds_.begin(); fedIter!=feds_.end(); ++fedIter) {
79  if ( (*fedIter) < 601 || (*fedIter) > 654) {
80  std::cout << "[EcalPnGraphs] pn number : " << (*fedIter) << " found in listFeds. "
81  << " Valid range is [601-654]. Returning." << std::endl;
82  inputIsOk = false;
83  return;
84  }
85  }
86  }
87 
88  //Check with Pns
89  if ( listPns[0] != -1 ) {
90  std::vector<int>::iterator intIter;
91  for (intIter = listPns.begin(); intIter != listPns.end(); intIter++) {
92  if ( ((*intIter) < 1) || (10 < (*intIter)) ) {
93  std::cout << "[EcalPnGraphs] pn number : " << (*intIter) << " found in listPns. "
94  << " Valid range is 1-10. Returning." << std::endl;
95  inputIsOk = false;
96  return;
97  }
98  if (!first_Pn ) first_Pn = (*intIter);
99  }
100  } else {
101  listPns.clear();
102  listPns.push_back(5);
103  listPns.push_back(6);
104  }
105 
106  // setting the abcissa array once for all
107  for (int i=0; i<50; i++) abscissa[i] = i;
108 
109  // local event counter (in general different from LV1)
110  eventCounter =0;
111 }
112 
113 
114 //=============================================================================
116 //=============================================================================
117  //delete *;
118 }
119 
120 //=============================================================================
122 //=============================================================================
123  edm::LogInfo("EcalPhGraphs") << "entering beginJob! " ;
124 }
125 
126 //=============================================================================
128 //=============================================================================
129 
130  eventCounter++;
131  if (!inputIsOk) 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  {
147  int ipn = (*pn_it);
148  int hpn = numPn;
149 
150  for (int u = (-hpn) ; u<=hpn; u++){
151  int ipn_c = ipn + u;
152  if (ipn_c < 1 || ipn_c > 10) 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 << " 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()) continue;
180 
181  for ( int i=0; i< (*pnItr).size() && i<50; ++i ) {
182  ordinate[i] = (*pnItr).sample(i).adc();
183  }
184  //make grapn of ph digis
185  TGraph oneGraph(50, abscissa,ordinate);
187  title = "Graph_ev" + intToString( eventCounter )
188  + "_FED" + intToString( FEDid )
189  + "_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  //
200  // outputs the number into the string stream and then flushes
201  // the buffer (makes sure the output is put into the stream)
202  //
203  std::ostringstream myStream; //creates an ostringstream object
204  myStream << num << std::flush;
205 
206  return(myStream.str()); //returns the string form of the stringstream object
207 }
208 
209 //=============================================================================
211 //=============================================================================
212  fileName += ( std::string("_Pn") + intToString(first_Pn) );
213  fileName += ".graph.root";
214 
215  root_file = new TFile( fileName.c_str() , "RECREATE" );
216  std::vector<TGraph>::iterator gr_it;
217  for ( gr_it = graphs.begin(); gr_it != graphs.end(); gr_it++ ) (*gr_it).Write();
218  root_file->Close();
219 
220  edm::LogInfo("EcalPnGraphs") << "DONE!.... " ;
221 }
222 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int abscissa[50]
Definition: EcalPnGraphs.h:67
int getFedFromSlice(std::string)
Definition: EcalFedMap.cc:104
std::vector< EcalPnDiodeDigi >::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:46
std::string fileName
Definition: EcalPnGraphs.h:58
void analyze(const edm::Event &e, const edm::EventSetup &c) override
void beginJob() override
EcalFedMap * fedMap
Definition: EcalPnGraphs.h:38
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:464
std::vector< TGraph > graphs
Definition: EcalPnGraphs.h:70
int ordinate[50]
Definition: EcalPnGraphs.h:68
const_iterator end() const
std::string digiProducer_
Definition: EcalPnGraphs.h:44
std::vector< int > listAllPns
Definition: EcalPnGraphs.h:63
~EcalPnGraphs() override
std::vector< std::string > ebs_
Definition: EcalPnGraphs.h:47
std::vector< int > listPns
Definition: EcalPnGraphs.h:62
TFile * root_file
Definition: EcalPnGraphs.h:72
const_iterator begin() const
void endJob() override
EcalPnGraphs(const edm::ParameterSet &ps)
Definition: EcalPnGraphs.cc:31