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