CMS 3D CMS Logo

CandOneToOneDeltaRMatcher Class Reference

Inheritance diagram for CandOneToOneDeltaRMatcher:

edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 CandOneToOneDeltaRMatcher (const edm::ParameterSet &)
 ~CandOneToOneDeltaRMatcher ()

Private Member Functions

std::vector< intAlgoBruteForce (int, int)
std::vector< intAlgoSwitchMethod (int, int)
double lenght (std::vector< int >)
void produce (edm::Event &, const edm::EventSetup &)

Private Attributes

std::string algoMethod_
std::vector< std::vector< float > > AllDist
edm::InputTag matched_
edm::InputTag source_


Detailed Description

Definition at line 17 of file CandOneToOneDeltaRMatcher.cc.


Constructor & Destructor Documentation

CandOneToOneDeltaRMatcher::CandOneToOneDeltaRMatcher ( const edm::ParameterSet cfg  ) 

Definition at line 59 of file CandOneToOneDeltaRMatcher.cc.

00059                                                                                :
00060   source_( cfg.getParameter<InputTag>( "src" ) ),
00061   matched_( cfg.getParameter<InputTag>( "matched" ) ),
00062   algoMethod_( cfg.getParameter<string>( "algoMethod" ) ) {
00063   produces<CandViewMatchMap>("src2mtc");
00064   produces<CandViewMatchMap>("mtc2src");
00065 }

CandOneToOneDeltaRMatcher::~CandOneToOneDeltaRMatcher (  ) 

Definition at line 67 of file CandOneToOneDeltaRMatcher.cc.

00067                                                       {
00068 }


Member Function Documentation

vector< int > CandOneToOneDeltaRMatcher::AlgoBruteForce ( int  nMin,
int  nMax 
) [private]

Definition at line 219 of file CandOneToOneDeltaRMatcher.cc.

References i1, lenght(), and stdcomb::next_combination().

Referenced by produce().

00219                                                                           {
00220 
00221   vector<int> ca;
00222   vector<int> cb;
00223   vector<int> bestCB;
00224   float totalDeltaR=0;
00225   float BestTotalDeltaR=1000;
00226 
00227   for(int i1=0; i1<nMax; i1++) ca.push_back(i1);
00228   for(int i1=0; i1<nMin; i1++) cb.push_back(i1);
00229 
00230   do
00231     {
00232       //do your processing on the new combination here
00233       for(int cnt=0;cnt<TMath::Factorial(nMin); cnt++)
00234         {
00235           totalDeltaR = lenght(cb);
00236           if ( totalDeltaR < BestTotalDeltaR ) {
00237             BestTotalDeltaR = totalDeltaR;
00238             bestCB=cb;
00239           }
00240           next_permutation( cb.begin() , cb.end() );
00241         }
00242     }
00243   while(next_combination( ca.begin() , ca.end() , cb.begin() , cb.end() ));
00244   
00245   return bestCB;
00246 }

vector< int > CandOneToOneDeltaRMatcher::AlgoSwitchMethod ( int  nMin,
int  nMax 
) [private]

Definition at line 259 of file CandOneToOneDeltaRMatcher.cc.

References AllDist, i1, i2, and i3.

Referenced by produce().

00259                                                                             {
00260 
00261   vector<int> bestCB;
00262   for(int i1=0; i1<nMin; i1++) {
00263     int minInd=0;
00264     for(int i2=1; i2<nMax; i2++) if( AllDist[i1][i2] < AllDist[i1][minInd] ) minInd = i2; 
00265     bestCB.push_back(minInd);
00266   }
00267 
00268   bool inside = true;
00269   while( inside ) {
00270     inside = false;
00271     for(int i1=0;i1<nMin;i1++){
00272       for(int i2=i1+1;i2<nMin;i2++){
00273         if ( bestCB[i1] == bestCB[i2] ) {
00274           inside = true;
00275           if ( AllDist[i1][(bestCB[i1])] <= AllDist[i2][(bestCB[i2])]) {
00276             AllDist[i2][(bestCB[i2])]= 1000;
00277             int minInd=0;
00278             for(int i3=1; i3<nMax; i3++) if( AllDist[i2][i3] < AllDist[i2][minInd] ) minInd = i3; 
00279             bestCB[i2]= minInd;
00280           }  else {
00281             AllDist[i1][(bestCB[i1])]= 1000;
00282             int minInd=0;
00283             for(int i3=1; i3<nMax; i3++) if( AllDist[i1][i3] < AllDist[i1][minInd] ) minInd = i3; 
00284             bestCB[i1]= minInd;
00285           }
00286         } // End if
00287       } 
00288     }
00289   } // End while
00290 
00291   return bestCB;
00292 
00293 }

double CandOneToOneDeltaRMatcher::lenght ( std::vector< int  )  [private]

Referenced by AlgoBruteForce().

void CandOneToOneDeltaRMatcher::produce ( edm::Event evt,
const edm::EventSetup es 
) [private, virtual]

Implements edm::EDProducer.

Definition at line 70 of file CandOneToOneDeltaRMatcher.cc.

References AlgoBruteForce(), algoMethod_, AlgoSwitchMethod(), AllDist, c, lat::endl(), Exception, edm::Event::getByLabel(), i1, int, matched_, max, min, edm::Event::put(), source, and source_.

00070                                                                           {
00071   
00072   Handle<CandidateView> source;  
00073   Handle<CandidateView> matched;  
00074   evt.getByLabel( source_, source ) ;
00075   evt.getByLabel( matched_, matched ) ;
00076  
00077   edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "======== Source Collection =======";
00078   for( CandidateView::const_iterator c = source->begin(); c != source->end(); ++c ) {
00079     edm::LogVerbatim("CandOneToOneDeltaRMatcher") << " pt source  " << c->pt() << " " << c->eta() << " " << c->phi()  << endl;
00080   }    
00081   edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "======== Matched Collection =======";
00082   for( CandidateView::const_iterator c = matched->begin(); c != matched->end(); ++c ) {
00083     edm::LogVerbatim("CandOneToOneDeltaRMatcher") << " pt source  " << c->pt() << " " << c->eta() << " " << c->phi()  << endl;
00084   } 
00085 
00086   const int nSrc = source->size();
00087   const int nMtc = matched->size();
00088 
00089   const int nMin = min( source->size() , matched->size() );
00090   const int nMax = max( source->size() , matched->size() );
00091   if( nMin < 1 ) return;
00092 
00093   if( nSrc <= nMtc ) {
00094     for(CandidateView::const_iterator iSr  = source->begin();
00095         iSr != source->end();
00096         iSr++) {
00097       vector <float> tempAllDist;
00098       for(CandidateView::const_iterator iMt  = matched->begin();
00099           iMt != matched->end();
00100           iMt++) { 
00101         tempAllDist.push_back(DeltaR( iSr->p4() , iMt->p4() ) );
00102       }
00103       AllDist.push_back(tempAllDist);
00104       tempAllDist.clear();
00105     } 
00106   } else {
00107     for(CandidateView::const_iterator iMt  = matched->begin();
00108         iMt != matched->end();
00109         iMt++) {
00110       vector <float> tempAllDist;
00111       for(CandidateView::const_iterator iSr  = source->begin();
00112           iSr != source->end();
00113           iSr++) { 
00114         tempAllDist.push_back(DeltaR( iSr->p4() , iMt->p4() ) );
00115       }
00116       AllDist.push_back(tempAllDist);
00117       tempAllDist.clear();
00118     } 
00119   }
00120   
00121   /*
00122   edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "======== The DeltaR Matrix =======";
00123   for(int m0=0; m0<nMin; m0++) {
00124     //    for(int m1=0; m1<nMax; m1++) {
00125       edm::LogVerbatim("CandOneToOneDeltaRMatcher") << setprecision(2) << fixed << (m1 AllDist[m0][m1] ;
00126     //}
00127     edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "\n"; 
00128   }
00129   */
00130   
00131   // Loop size if Brute Force
00132   int nLoopToDo = (int) ( TMath::Factorial(nMax) / TMath::Factorial(nMax - nMin) );
00133   edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "nLoop:" << nLoopToDo << endl;
00134   edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "Choosen Algo is:" << algoMethod_ ;
00135   vector<int> bestCB;
00136 
00137   // Algo is Brute Force
00138   if( algoMethod_ == "BruteForce") {
00139 
00140     bestCB = AlgoBruteForce(nMin,nMax);
00141 
00142   // Algo is Switch Method
00143   } else if( algoMethod_ == "SwitchMode" ) {
00144 
00145     bestCB = AlgoSwitchMethod(nMin,nMax);
00146 
00147   // Algo is Brute Force if nLoop < 10000
00148   } else if( algoMethod_ == "MixMode" ) {
00149 
00150     if( nLoopToDo < 10000 ) {
00151       bestCB = AlgoBruteForce(nMin,nMax);
00152     } else { 
00153       bestCB = AlgoSwitchMethod(nMin,nMax);
00154     } 
00155 
00156   } else {
00157     throw cms::Exception("OneToOne Constructor") << "wrong matching method in ParameterSet";
00158   }
00159 
00160   for(int i1=0; i1<nMin; i1++) edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "min: " << i1 << " " << bestCB[i1] << " " << AllDist[i1][bestCB[i1]];
00161 
00162 /*
00163   auto_ptr<CandViewMatchMap> matchMapSrMt( new CandViewMatchMap( CandViewMatchMap::ref_type( CandidateRefProd( source  ),
00164                                                                                              CandidateRefProd( matched )  ) ) );
00165   auto_ptr<CandViewMatchMap> matchMapMtSr( new CandViewMatchMap( CandViewMatchMap::ref_type( CandidateRefProd( matched ),
00166                                                                                              CandidateRefProd( source  )  ) ) );
00167 */
00168 
00169   auto_ptr<CandViewMatchMap> matchMapSrMt( new CandViewMatchMap() );
00170   auto_ptr<CandViewMatchMap> matchMapMtSr( new CandViewMatchMap() );
00171 
00172   for( int c = 0; c != nMin; c ++ ) {
00173     if( source->size() <= matched->size() ) {
00174       matchMapSrMt->insert( source ->refAt(c         ), matched->refAt(bestCB[c] ) );
00175       matchMapMtSr->insert( matched->refAt(bestCB[c] ), source ->refAt(c         ) );
00176     } else {
00177       matchMapSrMt->insert( source ->refAt(bestCB[c] ), matched->refAt(c         ) );
00178       matchMapMtSr->insert( matched->refAt(c         ), source ->refAt(bestCB[c] ) );
00179     }
00180   }
00181 
00182 /*
00183   for( int c = 0; c != nMin; c ++ ) {
00184     if( source->size() <= matched->size() ) { 
00185       matchMapSrMt->insert( CandidateRef( source,  c         ), CandidateRef( matched, bestCB[c] ) ); 
00186       matchMapMtSr->insert( CandidateRef( matched, bestCB[c] ), CandidateRef( source, c          ) );
00187     } else {
00188       matchMapSrMt->insert( CandidateRef( source,  bestCB[c] ), CandidateRef( matched, c         ) );
00189       matchMapMtSr->insert( CandidateRef( matched, c         ), CandidateRef( source,  bestCB[c] ) );
00190     }
00191   }
00192 */
00193   evt.put( matchMapSrMt, "src2mtc" );
00194   evt.put( matchMapMtSr, "mtc2src" );
00195 
00196   AllDist.clear();
00197 }


Member Data Documentation

std::string CandOneToOneDeltaRMatcher::algoMethod_ [private]

Definition at line 30 of file CandOneToOneDeltaRMatcher.cc.

Referenced by produce().

std::vector< std::vector<float> > CandOneToOneDeltaRMatcher::AllDist [private]

Definition at line 29 of file CandOneToOneDeltaRMatcher.cc.

Referenced by AlgoSwitchMethod(), and produce().

edm::InputTag CandOneToOneDeltaRMatcher::matched_ [private]

Definition at line 28 of file CandOneToOneDeltaRMatcher.cc.

Referenced by produce().

edm::InputTag CandOneToOneDeltaRMatcher::source_ [private]

Definition at line 27 of file CandOneToOneDeltaRMatcher.cc.

Referenced by produce().


The documentation for this class was generated from the following file:
Generated on Tue Jun 9 18:16:08 2009 for CMSSW by  doxygen 1.5.4