CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CandOneToOneDeltaRMatcher.cc
Go to the documentation of this file.
1 /* \class CandOneToOneDeltaRMatcher
2  *
3  * Producer for simple match map
4  * to match two collections of candidate
5  * with one-to-One matching
6  * minimizing Sum(DeltaR)
7  *
8  */
9 
12 
14 
15 #include<vector>
16 #include<iostream>
17 
19  public:
22  private:
23  void produce( edm::Event&, const edm::EventSetup& ) override;
24  double lenght( const std::vector<int>& );
25  std::vector<int> AlgoBruteForce(int, int);
26  std::vector<int> AlgoSwitchMethod(int, int);
27 
30  std::vector < std::vector<float> > AllDist;
32 
33 };
34 
36 
44 
49 
50 
51 #include <Math/VectorUtil.h>
52 #include <TMath.h>
53 
54 using namespace edm;
55 using namespace std;
56 using namespace reco;
57 using namespace ROOT::Math::VectorUtil;
58 using namespace stdcomb;
59 
61  sourceToken_( consumes<CandidateView>( cfg.getParameter<InputTag>( "src" ) ) ),
62  matchedToken_( consumes<CandidateView>( cfg.getParameter<InputTag>( "matched" ) ) ),
63  algoMethod_( cfg.getParameter<string>( "algoMethod" ) ) {
64  produces<CandViewMatchMap>("src2mtc");
65  produces<CandViewMatchMap>("mtc2src");
66 }
67 
69 }
70 
72 
74  Handle<CandidateView> matched;
75  evt.getByToken( sourceToken_, source ) ;
76  evt.getByToken( matchedToken_, matched ) ;
77 
78  edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "======== Source Collection =======";
79  for( CandidateView::const_iterator c = source->begin(); c != source->end(); ++c ) {
80  edm::LogVerbatim("CandOneToOneDeltaRMatcher") << " pt source " << c->pt() << " " << c->eta() << " " << c->phi() << endl;
81  }
82  edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "======== Matched Collection =======";
83  for( CandidateView::const_iterator c = matched->begin(); c != matched->end(); ++c ) {
84  edm::LogVerbatim("CandOneToOneDeltaRMatcher") << " pt source " << c->pt() << " " << c->eta() << " " << c->phi() << endl;
85  }
86 
87  const int nSrc = source->size();
88  const int nMtc = matched->size();
89 
90  const int nMin = min( source->size() , matched->size() );
91  const int nMax = max( source->size() , matched->size() );
92  if( nMin < 1 ) return;
93 
94  if( nSrc <= nMtc ) {
95  for(CandidateView::const_iterator iSr = source->begin();
96  iSr != source->end();
97  iSr++) {
98  vector <float> tempAllDist;
99  for(CandidateView::const_iterator iMt = matched->begin();
100  iMt != matched->end();
101  iMt++) {
102  tempAllDist.push_back(DeltaR( iSr->p4() , iMt->p4() ) );
103  }
104  AllDist.push_back(tempAllDist);
105  tempAllDist.clear();
106  }
107  } else {
108  for(CandidateView::const_iterator iMt = matched->begin();
109  iMt != matched->end();
110  iMt++) {
111  vector <float> tempAllDist;
112  for(CandidateView::const_iterator iSr = source->begin();
113  iSr != source->end();
114  iSr++) {
115  tempAllDist.push_back(DeltaR( iSr->p4() , iMt->p4() ) );
116  }
117  AllDist.push_back(tempAllDist);
118  tempAllDist.clear();
119  }
120  }
121 
122  /*
123  edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "======== The DeltaR Matrix =======";
124  for(int m0=0; m0<nMin; m0++) {
125  // for(int m1=0; m1<nMax; m1++) {
126  edm::LogVerbatim("CandOneToOneDeltaRMatcher") << setprecision(2) << fixed << (m1 AllDist[m0][m1] ;
127  //}
128  edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "\n";
129  }
130  */
131 
132  // Loop size if Brute Force
133  int nLoopToDo = (int) ( TMath::Factorial(nMax) / TMath::Factorial(nMax - nMin) );
134  edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "nLoop:" << nLoopToDo << endl;
135  edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "Choosen Algo is:" << algoMethod_ ;
136  vector<int> bestCB;
137 
138  // Algo is Brute Force
139  if( algoMethod_ == "BruteForce") {
140 
141  bestCB = AlgoBruteForce(nMin,nMax);
142 
143  // Algo is Switch Method
144  } else if( algoMethod_ == "SwitchMode" ) {
145 
146  bestCB = AlgoSwitchMethod(nMin,nMax);
147 
148  // Algo is Brute Force if nLoop < 10000
149  } else if( algoMethod_ == "MixMode" ) {
150 
151  if( nLoopToDo < 10000 ) {
152  bestCB = AlgoBruteForce(nMin,nMax);
153  } else {
154  bestCB = AlgoSwitchMethod(nMin,nMax);
155  }
156 
157  } else {
158  throw cms::Exception("OneToOne Constructor") << "wrong matching method in ParameterSet";
159  }
160 
161  for(int i1=0; i1<nMin; i1++) edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "min: " << i1 << " " << bestCB[i1] << " " << AllDist[i1][bestCB[i1]];
162 
163 /*
164  auto_ptr<CandViewMatchMap> matchMapSrMt( new CandViewMatchMap( CandViewMatchMap::ref_type( CandidateRefProd( source ),
165  CandidateRefProd( matched ) ) ) );
166  auto_ptr<CandViewMatchMap> matchMapMtSr( new CandViewMatchMap( CandViewMatchMap::ref_type( CandidateRefProd( matched ),
167  CandidateRefProd( source ) ) ) );
168 */
169 
170  auto_ptr<CandViewMatchMap> matchMapSrMt( new CandViewMatchMap() );
171  auto_ptr<CandViewMatchMap> matchMapMtSr( new CandViewMatchMap() );
172 
173  for( int c = 0; c != nMin; c ++ ) {
174  if( source->size() <= matched->size() ) {
175  matchMapSrMt->insert( source ->refAt(c ), matched->refAt(bestCB[c] ) );
176  matchMapMtSr->insert( matched->refAt(bestCB[c] ), source ->refAt(c ) );
177  } else {
178  matchMapSrMt->insert( source ->refAt(bestCB[c] ), matched->refAt(c ) );
179  matchMapMtSr->insert( matched->refAt(c ), source ->refAt(bestCB[c] ) );
180  }
181  }
182 
183 /*
184  for( int c = 0; c != nMin; c ++ ) {
185  if( source->size() <= matched->size() ) {
186  matchMapSrMt->insert( CandidateRef( source, c ), CandidateRef( matched, bestCB[c] ) );
187  matchMapMtSr->insert( CandidateRef( matched, bestCB[c] ), CandidateRef( source, c ) );
188  } else {
189  matchMapSrMt->insert( CandidateRef( source, bestCB[c] ), CandidateRef( matched, c ) );
190  matchMapMtSr->insert( CandidateRef( matched, c ), CandidateRef( source, bestCB[c] ) );
191  }
192  }
193 */
194  evt.put( matchMapSrMt, "src2mtc" );
195  evt.put( matchMapMtSr, "mtc2src" );
196 
197  AllDist.clear();
198 }
199 
200 
201 double CandOneToOneDeltaRMatcher::lenght(const vector<int>& best) {
202  double myLenght=0;
203  int row=0;
204  for(vector<int>::const_iterator it=best.begin(); it!=best.end(); it++ ) {
205  myLenght+=AllDist[row][*it];
206  row++;
207  }
208  return myLenght;
209 }
210 
211 // this is the Brute Force Algorithm
212 // All the possible combination are checked
213 // The best one is always found
214 // Be carefull when you have high values for nMin and nMax --> the combinatorial could explode!
215 // Sum(DeltaR) is minimized -->
216 // 0.1 - 0.2 - 1.0 - 1.5 is lower than
217 // 0.1 - 0.2 - 0.3 - 3.0
218 // Which one do you prefer? --> BruteForce select always the first
219 
220 vector<int> CandOneToOneDeltaRMatcher::AlgoBruteForce( int nMin, int nMax ) {
221 
222  vector<int> ca;
223  vector<int> cb;
224  vector<int> bestCB;
225  float totalDeltaR=0;
226  float BestTotalDeltaR=1000;
227 
228  for(int i1=0; i1<nMax; i1++) ca.push_back(i1);
229  for(int i1=0; i1<nMin; i1++) cb.push_back(i1);
230 
231  do
232  {
233  //do your processing on the new combination here
234  for(int cnt=0;cnt<TMath::Factorial(nMin); cnt++)
235  {
236  totalDeltaR = lenght(cb);
237  if ( totalDeltaR < BestTotalDeltaR ) {
238  BestTotalDeltaR = totalDeltaR;
239  bestCB=cb;
240  }
241  next_permutation( cb.begin() , cb.end() );
242  }
243  }
244  while(next_combination( ca.begin() , ca.end() , cb.begin() , cb.end() ));
245 
246  return bestCB;
247 }
248 
249 // This method (Developed originally by Daniele Benedetti) check for the best combination
250 // choosing the minimum DeltaR for each line in AllDist matrix
251 // If no repeated row is found: ie (line,col)=(1,3) and (2,3) --> same as BruteForce
252 // If repetition --> set the higher DeltaR between the 2 repetition to 1000 and re-check best combination
253 // Iterate until no repetition
254 // No guaranted minimum for Sum(DeltaR)
255 // If you have:
256 // 0.1 - 0.2 - 1.0 - 1.5 is lower than
257 // 0.1 - 0.2 - 0.3 - 3.0
258 // SwitchMethod normally select the second solution
259 
260 vector<int> CandOneToOneDeltaRMatcher::AlgoSwitchMethod( int nMin, int nMax ) {
261 
262  vector<int> bestCB;
263  for(int i1=0; i1<nMin; i1++) {
264  int minInd=0;
265  for(int i2=1; i2<nMax; i2++) if( AllDist[i1][i2] < AllDist[i1][minInd] ) minInd = i2;
266  bestCB.push_back(minInd);
267  }
268 
269  bool inside = true;
270  while( inside ) {
271  inside = false;
272  for(int i1=0;i1<nMin;i1++){
273  for(int i2=i1+1;i2<nMin;i2++){
274  if ( bestCB[i1] == bestCB[i2] ) {
275  inside = true;
276  if ( AllDist[i1][(bestCB[i1])] <= AllDist[i2][(bestCB[i2])]) {
277  AllDist[i2][(bestCB[i2])]= 1000;
278  int minInd=0;
279  for(int i3=1; i3<nMax; i3++) if( AllDist[i2][i3] < AllDist[i2][minInd] ) minInd = i3;
280  bestCB[i2]= minInd;
281  } else {
282  AllDist[i1][(bestCB[i1])]= 1000;
283  int minInd=0;
284  for(int i3=1; i3<nMax; i3++) if( AllDist[i1][i3] < AllDist[i1][minInd] ) minInd = i3;
285  bestCB[i1]= minInd;
286  }
287  } // End if
288  }
289  }
290  } // End while
291 
292  return bestCB;
293 
294 }
295 
298 
edm::EDGetTokenT< reco::CandidateView > sourceToken_
tuple cfg
Definition: looper.py:259
std::vector< std::vector< float > > AllDist
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< int > AlgoSwitchMethod(int, int)
std::vector< int > AlgoBruteForce(int, int)
CandOneToOneDeltaRMatcher(const edm::ParameterSet &)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
T min(T a, T b)
Definition: MathUtil.h:58
void produce(edm::Event &, const edm::EventSetup &) override
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85
double lenght(const std::vector< int > &)
bool next_combination(BidIt n_begin, BidIt n_end, BidIt r_begin, BidIt r_end)
Definition: combination.h:22
static std::string const source
Definition: EdmProvDump.cc:42
edm::AssociationMap< edm::OneToOneGeneric< reco::CandidateView, reco::CandidateView > > CandViewMatchMap
Definition: CandMatchMap.h:25
edm::EDGetTokenT< reco::CandidateView > matchedToken_