CMS 3D CMS Logo

Digi2Raw2Digi.cc
Go to the documentation of this file.
5 
12 
15 
22 #include <vector>
23 #include <utility>
24 #include <ostream>
25 #include <string>
26 #include <algorithm>
27 #include <cmath>
28 
29 
30 // Qunatities of interest in :
31 // DataFormats/HcalDigi/test/HcalDigiDump.cc - example of Digi dumping...
32 // DataFormats/HcalDigi/interface/HcalQIESample.h - adc, capid etc.
33 // DataFormats/HcalDigi/interface/HBHEDataFrame.h -
34 // zsMarkAndPass, zsUnsuppressed etc.
35 
36 template<class Digi>
37 
38 
40 
41  typename edm::Handle<edm::SortedCollection<Digi> > digiCollection1;
43  typename edm::Handle<edm::SortedCollection<Digi> > digiCollection2;
45 
46  if(unsuppressed) { // ZDC
47  iEvent.getByToken (tok1, digiCollection1);
48  }
49  else iEvent.getByToken (tok1, digiCollection1);
50 
51  iEvent.getByToken (tok2, digiCollection2);
52 
53  int size1 = 0;
54  int size2 = 0;
55 
56  for (digiItr1=digiCollection1->begin();digiItr1!=digiCollection1->end();digiItr1++) {
57  size1++;
58  }
59 
60  for (digiItr2=digiCollection2->begin();digiItr2!=digiCollection2->end();digiItr2++) {
61  size2++;
62  }
63 
64  //std::cout << "Digi collections size1 = "<< size1
65  // << " size2 = " << size2 << std::endl;
66 
67 
68  // CYCLE over first DIGI collection ======================================
69 
70  for (digiItr1=digiCollection1->begin();digiItr1!=digiCollection1->end();digiItr1++) {
71  HcalGenericDetId HcalGenDetId(digiItr1->id());
72  int tsize = (*digiItr1).size();
73  int match = 0;
74 
75  if(HcalGenDetId.isHcalZDCDetId()){
76  //for zdc
77  HcalZDCDetId element(digiItr1->id());
78  int zside = element.zside();
79  int section = element.section();
80  int channel = element.channel();
81  int gsub = HcalGenDetId.genericSubdet();
82 
83  if(section==3){// lumi section not reconstructed
84  size2++;
85  match = 1;
86  goto lumi;
87  }
88 
89 
90  //std::cout<< " Zdc genSubdet="<< gsub << " zside=" <<zside
91  // << " section= "<< section << " channel " <<channel
92  // <<std::endl;
93 
94  for (digiItr2=digiCollection2->begin();digiItr2!=digiCollection2->end();digiItr2++) {
95  HcalZDCDetId element2(digiItr2->id());
96 
97  //int zside2 = element2.zside();
98  //int section2 = element2.section();
99  //int channel2 = element2.channel();
100  //int gsub2 = HcalGenDetId.genericSubdet();
101 
102  //std::cout<< " Zdc genSubdet="<<gsub2
103  // << " zside=" <<zside2
104  // << " section= "<<section2
105  // <<" channel "<<channel2
106  // <<std::endl;
107 
108  if(element == element2) {
109  match = 1;
110  int identical = 1;
111  for (int i=0; i<tsize; i++) {
112  double adc = (*digiItr1)[i].adc();
113  int capid = (*digiItr1)[i].capid();
114  // std::cout << std::endl << " capid1=" << capid
115  // << " adc1=" << adc
116  // << std::endl;
117  double adc2 = (*digiItr2)[i].adc();
118  int capid2 = (*digiItr2)[i].capid();
119  // std::cout << " capid2=" << capid2 << " adc2=" << adc2
120  // << std::endl;
121  if( capid != capid2 || adc != adc2) {
122  std::cout << "===> PROBLEM !!! gebsubdet=" << gsub
123  << " zside=" <<zside
124  << " section= "<< section << " channel " <<channel
125  << std::endl;
126  std::cout << " capid1["<< i << "]=" << capid
127  << " adc1["<< i << "]=" << adc
128  << " capid2["<< i << "]=" << capid2
129  << " adc2["<< i << "]=" << adc2
130  << std::endl;
131  identical = 0;
132  meStatus->Fill(1.);
133  break;
134  }
135  } // end of DataFrames array
136  if(identical) meStatus->Fill(0.);
137  break; // matched HcalZDCID is processed,
138  // go to next (primary collection) cell
139  }
140  } // end of cycle over 2d DIGI collection
141  lumi:
142  if (!match) {
143  meStatus->Fill(2.);
144  std::cout << "===> PROBLEM !!! gsubdet=" << gsub
145  << " zside=" <<zside
146  << " section= "<< section << " channel " <<channel
147  << " HcalZDCId match is not found !!!"
148  << std::endl;
149  }
150 
151  }
152  else{
153  //for Hcal subdetectors
154  HcalDetId cell(digiItr1->id());
155  int depth = cell.depth();
156  int iphi = cell.iphi()-1;
157  int ieta = cell.ieta();
158  int sub = cell.subdet();
159  // if(ieta > 0) ieta--;
160  // std::cout << " Cell subdet=" << sub << " ieta=" << ieta
161  // << " inphi=" << iphi << " depth=" << depth << std::endl;
162 
163  // CYCLE over second DIGI collection ======================================
164  for (digiItr2=digiCollection2->begin();digiItr2!=digiCollection2->end();digiItr2++) {
165 
166  HcalDetId cell2(digiItr2->id());
167 
168  if( cell == cell2) {
169  match = 1;
170  int identical = 1;
171  for (int i=0; i<tsize; i++) {
172  double adc = (*digiItr1)[i].adc();
173  int capid = (*digiItr1)[i].capid();
174  // std::cout << std::endl << " capid1=" << capid
175  // << " adc1=" << adc
176  // << std::endl;
177  double adc2 = (*digiItr2)[i].adc();
178  int capid2 = (*digiItr2)[i].capid();
179  // std::cout << " capid2=" << capid2 << " adc2=" << adc2
180  // << std::endl;
181  if( capid != capid2 || adc != adc2) {
182  std::cout << "===> PROBLEM !!! subdet=" << sub << " ieta="
183  << ieta << " inphi=" << iphi << " depth=" << depth
184  << std::endl;
185  std::cout << " capid1["<< i << "]=" << capid
186  << " adc1["<< i << "]=" << adc
187  << " capid2["<< i << "]=" << capid2
188  << " adc2["<< i << "]=" << adc2
189  << std::endl;
190  identical = 0;
191  meStatus->Fill(1.);
192  break;
193  }
194  } // end of DataFrames array
195  if(identical) meStatus->Fill(0.);
196  break; // matched HcalID is processed,
197  // go to next (primary collection) cell
198  }
199  } // end of cycle over 2d DIGI collection
200  if (!match) {
201  meStatus->Fill(2.);
202  std::cout << "===> PROBLEM !!! subdet=" << sub << " ieta="
203  << ieta << " inphi=" << iphi << " depth=" << depth
204  << " HcalID match is not found !!!"
205  << std::endl;
206  }
207  }
208  } // end of cycle over 1st DIGI collection
209 
210  if (size1 != size2) {
211  meStatus->Fill(3.);
212  std::cout << "===> PROBLEM !!! Different size of Digi collections : "
213  << size1 << " and " << size2
214  << std::endl;
215  }
216 }
217 
218 
220  : inputTag1_(iConfig.getParameter<edm::InputTag>("digiLabel1")),
221  inputTag2_(iConfig.getParameter<edm::InputTag>("digiLabel2")),
222  outputFile_(iConfig.getUntrackedParameter<std::string>("outputFile"))
223 {
224 
225  // register for data access
226  tok_hbhe1_ = consumes<edm::SortedCollection<HBHEDataFrame> >(inputTag1_);
227  tok_hbhe2_ = consumes<edm::SortedCollection<HBHEDataFrame> >(inputTag2_);
228  tok_ho1_ = consumes<edm::SortedCollection<HODataFrame> >(inputTag1_);
229  tok_ho2_ = consumes<edm::SortedCollection<HODataFrame> >(inputTag2_);
230  tok_hf1_ = consumes<edm::SortedCollection<HFDataFrame> >(inputTag1_);
231  tok_hf2_ = consumes<edm::SortedCollection<HFDataFrame> >(inputTag2_);
232  tok_zdc1_ = consumes<edm::SortedCollection<ZDCDataFrame> >(edm::InputTag("simHcalUnsuppressedDigis"));
233  tok_zdc2_ = consumes<edm::SortedCollection<ZDCDataFrame> >(inputTag2_);
234 
235 
236  // DQM ROOT output
237  if ( outputFile_.size() != 0 ) {
238  edm::LogInfo("OutputInfo")
239  << " Hcal RecHit Task histograms will be saved to '"
240  << outputFile_.c_str() << "'";
241  } else {
242  edm::LogInfo("OutputInfo")
243  << " Hcal RecHit Task histograms will NOT be saved";
244  }
245 
246 
247 }
248 
250 {
251 
252  ibooker.setCurrentFolder("Digi2Raw2DigiV/Digi2Raw2DigiTask");
253 
254 
255  // const char * sub = hcalselector_.c_str();
256  char histo[100];
257 
258  sprintf (histo, "Digi2Raw2Digi_status") ;
259  // bins: 1)full match 2)ID match, not content 3) no match
260  // 4) number of events with diff number of Digis
261  meStatus = ibooker.book1D(histo, histo, 5, 0., 5.);
262 
263 }
264 
266 
267 
269 {
270  unsuppressed = 0;
271 
272  // std::cout << "=== HBHE ==================" << std::endl;
273  compare<HBHEDataFrame>(iEvent,iSetup,tok_hbhe1_,tok_hbhe2_);
274 
275  // std::cout << "=== HO ====================" << std::endl;
276  compare<HODataFrame>(iEvent,iSetup,tok_ho1_,tok_ho2_);
277 
278  // std::cout << "=== HF ====================" << std::endl;
279  compare<HFDataFrame>(iEvent,iSetup,tok_hf1_,tok_hf2_);
280 
281 
282  // std::cout << "=== ZDC ===================" << std::endl;
283  unsuppressed = 1;
284  compare<ZDCDataFrame>(iEvent,iSetup,tok_zdc1_,tok_zdc2_);
285 
286 
287  // std::cout << "=== CASTOR ================" << std::endl;
288  // compare<CastorDataFrame>(iEvent,iSetup);
289 
290 }
291 
292 
int adc(sample_type sample)
get the ADC sample (12 bits)
Digi2Raw2Digi(const edm::ParameterSet &)
std::string outputFile_
Definition: Digi2Raw2Digi.h:41
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< edm::SortedCollection< HBHEDataFrame > > tok_hbhe2_
Definition: Digi2Raw2Digi.h:33
std::vector< T >::const_iterator const_iterator
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< edm::SortedCollection< ZDCDataFrame > > tok_zdc2_
Definition: Digi2Raw2Digi.h:39
virtual void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &)
edm::EDGetTokenT< edm::SortedCollection< HFDataFrame > > tok_hf2_
Definition: Digi2Raw2Digi.h:37
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
edm::EDGetTokenT< edm::SortedCollection< ZDCDataFrame > > tok_zdc1_
Definition: Digi2Raw2Digi.h:38
void compare(const edm::Event &, const edm::EventSetup &, const edm::EDGetTokenT< edm::SortedCollection< Digi > > &tok1, const edm::EDGetTokenT< edm::SortedCollection< Digi > > &tok2)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:277
edm::EDGetTokenT< edm::SortedCollection< HBHEDataFrame > > tok_hbhe1_
Definition: Digi2Raw2Digi.h:32
edm::EDGetTokenT< edm::SortedCollection< HODataFrame > > tok_ho1_
Definition: Digi2Raw2Digi.h:34
edm::EDGetTokenT< edm::SortedCollection< HFDataFrame > > tok_hf1_
Definition: Digi2Raw2Digi.h:36
HLT enums.
edm::EDGetTokenT< edm::SortedCollection< HODataFrame > > tok_ho2_
Definition: Digi2Raw2Digi.h:35
edm::InputTag inputTag2_
Definition: Digi2Raw2Digi.h:30
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
edm::InputTag inputTag1_
Definition: Digi2Raw2Digi.h:29
virtual void analyze(const edm::Event &, const edm::EventSetup &)
Definition: Run.h:42
MonitorElement * meStatus
Definition: Digi2Raw2Digi.h:43