CMS 3D CMS Logo

compareCands.h
Go to the documentation of this file.
1 #ifndef compareCands_h
2 #define compareCands_h
3 
13 #include "TH2.h"
14 #include "TH1.h"
16 
17 //first declare stuff about the class template
18 //notice that we pass our own defined struct of data with necessary mbx information (we can modify this in future without changing much)
19 //we also avoid messing about with class hierarchy...
20 template <class T>
21 class compareCands {
22  public:
23  compareCands(const T &data, const T &emu, const GctErrorAnalyzerMBxInfo &mbxparams);
24  ~compareCands();
25 
26  bool doCompare(TH1I *errorFlag_hist_, TH1I *mismatchD_Rank, TH2I *mismatchD_EtEtaPhi, TH1I *mismatchE_Rank, TH2I *mismatchE_EtEtaPhi);
27 
28  private:
31 
32 };
33 
34 //now the implementation
35 template<class T>
36 compareCands<T>::compareCands(const T &data, const T &emu, const GctErrorAnalyzerMBxInfo &mbxparams) :
37  data_(data),
38  emu_(emu),
39  mbxparams_(mbxparams)
40 {
41  //std::cout << "initialising..." << std::endl;
42 }
43 
44 template<class T>
46  //anything need destructing?
47 }
48 
49 template<class T>
50 bool compareCands<T>::doCompare(TH1I *errorFlag_hist_, TH1I *mismatchD_Rank, TH2I *mismatchD_EtEtaPhi, TH1I *mismatchE_Rank, TH2I *mismatchE_EtEtaPhi) {
51 
52  //this code has now been patched to be multiple_bx compliant. However, this still means that only 1 comparison will happen per event, and this has to be
53  //matched such that the RCTTrigBx(=0) data is run over the emulator and the EmuTrigBx analysis (Bx=0) corresponds to the GCTTrigBx (Bx=0) analysis
54  //These TrigBx parameters are set in the configuration to make things more flexible if things change later
55 
56  //define some temporary local variables
57  bool errorFlag=false;
58  unsigned int i=0, j=0;
59  std::vector<bool> matched(GCT_OBJECT_QUANTA);
60  //this makes a vector of GCT_OBJECT_QUANTA=4 bools, all set to false
61  //remember that pushing back will make the vector larger!
62 
63  for(i=0; i < data_->size(); i++) {
64 
65  //The first thing to check is that the BX of the data corresponds to the trig Bx (we expect these to be contiguous i.e. data sorted in Bx)
66  if(data_->at(i).bx() != mbxparams_.GCTTrigBx) continue;
67 
68  //If the data candidate has zero rank, move to the next data candidate
69  //since all the candidates are ranked in order of rank, this implies all the remaining data candidates also have rank = 0
70  if(data_->at(i).rank() == 0 ) continue;
71 
72  for(j=0; j < emu_->size(); j++) {
73 
74  //Again, the first thing to check in this loop is that the BX of the emulator data corresponds to the trig Bx
75  if(emu_->at(j).bx() != mbxparams_.EmuTrigBx) continue;
76 
77  if( data_->at(i).rank() == emu_->at(j).rank()
78  && data_->at(i).regionId().ieta() == emu_->at(j).regionId().ieta()
79  && data_->at(i).regionId().iphi() == emu_->at(j).regionId().iphi()
80  && matched.at((j % GCT_OBJECT_QUANTA)) == 0 ) {
81  //this means that the ith data candidate matches the jth emulator candidate
82  errorFlag_hist_->Fill(0); //fill the errorflag histo in the matched bin
83  matched.at((j % GCT_OBJECT_QUANTA)) = true; //set emulator candidate to matched so it doesn't get re-used
84  break; //matched the current data candidate, now move to the next
85  }
86 
87  if( (j % GCT_OBJECT_QUANTA) + 1 == GCT_OBJECT_QUANTA) {
88  errorFlag_hist_->Fill(1); //fill the errorflag histo in the unmatched data candidate bin
89  mismatchD_Rank->Fill(data_->at(i).rank()); //fill the rank histogram of mismatched data candidates
90  mismatchD_EtEtaPhi->Fill(data_->at(i).regionId().ieta(),data_->at(i).regionId().iphi(),data_->at(i).rank()); //fill the EtEtaPhi dist for mismatched candidates
91  errorFlag=true; //set the errorFlag to true
92  }
93 
94  }
95  }
96 
97 //loop over the matched boolean vector and see if there are any rank>0 unmatched emu candidates - if there are populate the histogram in the emulator mismatched bin
98  for(i=0; i<matched.size(); i++) {
99  //the first thing to check is that the matched flag for object i out of 0,1,2,3 (0 -> GCT_OBJECT_QUANTA-1) is not set - then we can check that the corresponding
100  //emulator candidates either have rank = 0 (which is good) or rank > 0 (which is bad)
101  if(matched.at(i)) continue;
102 
103  //now loop over the emulator candidates
104  for(j=0; j< emu_->size(); j++) {
105  //check that the bx of the emulator candidates is the trigbx
106  if(emu_->at(j).bx() != mbxparams_.EmuTrigBx) continue;
107 
108  //now check that the j%GCT_OBJECT_QUANTA is the same as the index of the false entry in the bool_matched vector so that we are looking at the right candidate
109  if((j%GCT_OBJECT_QUANTA == i) && (emu_->at(j).rank()>0) ) {
110  errorFlag_hist_->Fill(2); //increment emulator mismatched bin
111  mismatchE_Rank->Fill(emu_->at(j).rank()); //fill the rank histogram for unmatched emulator
112  mismatchE_EtEtaPhi->Fill(emu_->at(j).regionId().ieta(),emu_->at(j).regionId().iphi(),emu_->at(j).rank()); //fill EtEtaPhi for unmatched emu cands
113  errorFlag=true; //set the errorFlag (if it's not already)
114  }
115  }
116  }
117 
118  return errorFlag;
119 }
120 
121 
122 #endif
bool doCompare(TH1I *errorFlag_hist_, TH1I *mismatchD_Rank, TH2I *mismatchD_EtEtaPhi, TH1I *mismatchE_Rank, TH2I *mismatchE_EtEtaPhi)
Definition: compareCands.h:50
compareCands(const T &data, const T &emu, const GctErrorAnalyzerMBxInfo &mbxparams)
Definition: compareCands.h:36
const unsigned int GCT_OBJECT_QUANTA
GctErrorAnalyzerMBxInfo mbxparams_
Definition: compareCands.h:30
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
long double T