CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

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< int > AlgoBruteForce (int, int)
std::vector< int > AlgoSwitchMethod (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.

                                                                               :
  source_( cfg.getParameter<InputTag>( "src" ) ),
  matched_( cfg.getParameter<InputTag>( "matched" ) ),
  algoMethod_( cfg.getParameter<string>( "algoMethod" ) ) {
  produces<CandViewMatchMap>("src2mtc");
  produces<CandViewMatchMap>("mtc2src");
}
CandOneToOneDeltaRMatcher::~CandOneToOneDeltaRMatcher ( )

Definition at line 67 of file CandOneToOneDeltaRMatcher.cc.

                                                      {
}

Member Function Documentation

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

Definition at line 219 of file CandOneToOneDeltaRMatcher.cc.

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

Referenced by produce().

                                                                          {

  vector<int> ca;
  vector<int> cb;
  vector<int> bestCB;
  float totalDeltaR=0;
  float BestTotalDeltaR=1000;

  for(int i1=0; i1<nMax; i1++) ca.push_back(i1);
  for(int i1=0; i1<nMin; i1++) cb.push_back(i1);

  do
    {
      //do your processing on the new combination here
      for(int cnt=0;cnt<TMath::Factorial(nMin); cnt++)
        {
          totalDeltaR = lenght(cb);
          if ( totalDeltaR < BestTotalDeltaR ) {
            BestTotalDeltaR = totalDeltaR;
            bestCB=cb;
          }
          next_permutation( cb.begin() , cb.end() );
        }
    }
  while(next_combination( ca.begin() , ca.end() , cb.begin() , cb.end() ));
  
  return bestCB;
}
vector< int > CandOneToOneDeltaRMatcher::AlgoSwitchMethod ( int  nMin,
int  nMax 
) [private]

Definition at line 259 of file CandOneToOneDeltaRMatcher.cc.

References AllDist.

Referenced by produce().

                                                                            {

  vector<int> bestCB;
  for(int i1=0; i1<nMin; i1++) {
    int minInd=0;
    for(int i2=1; i2<nMax; i2++) if( AllDist[i1][i2] < AllDist[i1][minInd] ) minInd = i2; 
    bestCB.push_back(minInd);
  }

  bool inside = true;
  while( inside ) {
    inside = false;
    for(int i1=0;i1<nMin;i1++){
      for(int i2=i1+1;i2<nMin;i2++){
        if ( bestCB[i1] == bestCB[i2] ) {
          inside = true;
          if ( AllDist[i1][(bestCB[i1])] <= AllDist[i2][(bestCB[i2])]) {
            AllDist[i2][(bestCB[i2])]= 1000;
            int minInd=0;
            for(int i3=1; i3<nMax; i3++) if( AllDist[i2][i3] < AllDist[i2][minInd] ) minInd = i3; 
            bestCB[i2]= minInd;
          }  else {
            AllDist[i1][(bestCB[i1])]= 1000;
            int minInd=0;
            for(int i3=1; i3<nMax; i3++) if( AllDist[i1][i3] < AllDist[i1][minInd] ) minInd = i3; 
            bestCB[i1]= minInd;
          }
        } // End if
      } 
    }
  } // End while

  return bestCB;

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

Definition at line 200 of file CandOneToOneDeltaRMatcher.cc.

References AllDist.

Referenced by AlgoBruteForce().

                                                         {
  double myLenght=0;
  int row=0;
  for(vector<int>::iterator it=best.begin(); it!=best.end(); it++ ) {           
    myLenght+=AllDist[row][*it];
    row++;
  }
  return myLenght;
}
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, trackerHits::c, Exception, edm::Event::getByLabel(), matched_, max(), min, edm::Event::put(), LaserTracksInput_cfi::source, and source_.

                                                                          {
  
  Handle<CandidateView> source;  
  Handle<CandidateView> matched;  
  evt.getByLabel( source_, source ) ;
  evt.getByLabel( matched_, matched ) ;
 
  edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "======== Source Collection =======";
  for( CandidateView::const_iterator c = source->begin(); c != source->end(); ++c ) {
    edm::LogVerbatim("CandOneToOneDeltaRMatcher") << " pt source  " << c->pt() << " " << c->eta() << " " << c->phi()  << endl;
  }    
  edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "======== Matched Collection =======";
  for( CandidateView::const_iterator c = matched->begin(); c != matched->end(); ++c ) {
    edm::LogVerbatim("CandOneToOneDeltaRMatcher") << " pt source  " << c->pt() << " " << c->eta() << " " << c->phi()  << endl;
  } 

  const int nSrc = source->size();
  const int nMtc = matched->size();

  const int nMin = min( source->size() , matched->size() );
  const int nMax = max( source->size() , matched->size() );
  if( nMin < 1 ) return;

  if( nSrc <= nMtc ) {
    for(CandidateView::const_iterator iSr  = source->begin();
        iSr != source->end();
        iSr++) {
      vector <float> tempAllDist;
      for(CandidateView::const_iterator iMt  = matched->begin();
          iMt != matched->end();
          iMt++) { 
        tempAllDist.push_back(DeltaR( iSr->p4() , iMt->p4() ) );
      }
      AllDist.push_back(tempAllDist);
      tempAllDist.clear();
    } 
  } else {
    for(CandidateView::const_iterator iMt  = matched->begin();
        iMt != matched->end();
        iMt++) {
      vector <float> tempAllDist;
      for(CandidateView::const_iterator iSr  = source->begin();
          iSr != source->end();
          iSr++) { 
        tempAllDist.push_back(DeltaR( iSr->p4() , iMt->p4() ) );
      }
      AllDist.push_back(tempAllDist);
      tempAllDist.clear();
    } 
  }
  
  /*
  edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "======== The DeltaR Matrix =======";
  for(int m0=0; m0<nMin; m0++) {
    //    for(int m1=0; m1<nMax; m1++) {
      edm::LogVerbatim("CandOneToOneDeltaRMatcher") << setprecision(2) << fixed << (m1 AllDist[m0][m1] ;
    //}
    edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "\n"; 
  }
  */
  
  // Loop size if Brute Force
  int nLoopToDo = (int) ( TMath::Factorial(nMax) / TMath::Factorial(nMax - nMin) );
  edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "nLoop:" << nLoopToDo << endl;
  edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "Choosen Algo is:" << algoMethod_ ;
  vector<int> bestCB;

  // Algo is Brute Force
  if( algoMethod_ == "BruteForce") {

    bestCB = AlgoBruteForce(nMin,nMax);

  // Algo is Switch Method
  } else if( algoMethod_ == "SwitchMode" ) {

    bestCB = AlgoSwitchMethod(nMin,nMax);

  // Algo is Brute Force if nLoop < 10000
  } else if( algoMethod_ == "MixMode" ) {

    if( nLoopToDo < 10000 ) {
      bestCB = AlgoBruteForce(nMin,nMax);
    } else { 
      bestCB = AlgoSwitchMethod(nMin,nMax);
    } 

  } else {
    throw cms::Exception("OneToOne Constructor") << "wrong matching method in ParameterSet";
  }

  for(int i1=0; i1<nMin; i1++) edm::LogVerbatim("CandOneToOneDeltaRMatcher") << "min: " << i1 << " " << bestCB[i1] << " " << AllDist[i1][bestCB[i1]];

/*
  auto_ptr<CandViewMatchMap> matchMapSrMt( new CandViewMatchMap( CandViewMatchMap::ref_type( CandidateRefProd( source  ),
                                                                                             CandidateRefProd( matched )  ) ) );
  auto_ptr<CandViewMatchMap> matchMapMtSr( new CandViewMatchMap( CandViewMatchMap::ref_type( CandidateRefProd( matched ),
                                                                                             CandidateRefProd( source  )  ) ) );
*/

  auto_ptr<CandViewMatchMap> matchMapSrMt( new CandViewMatchMap() );
  auto_ptr<CandViewMatchMap> matchMapMtSr( new CandViewMatchMap() );

  for( int c = 0; c != nMin; c ++ ) {
    if( source->size() <= matched->size() ) {
      matchMapSrMt->insert( source ->refAt(c         ), matched->refAt(bestCB[c] ) );
      matchMapMtSr->insert( matched->refAt(bestCB[c] ), source ->refAt(c         ) );
    } else {
      matchMapSrMt->insert( source ->refAt(bestCB[c] ), matched->refAt(c         ) );
      matchMapMtSr->insert( matched->refAt(c         ), source ->refAt(bestCB[c] ) );
    }
  }

/*
  for( int c = 0; c != nMin; c ++ ) {
    if( source->size() <= matched->size() ) { 
      matchMapSrMt->insert( CandidateRef( source,  c         ), CandidateRef( matched, bestCB[c] ) ); 
      matchMapMtSr->insert( CandidateRef( matched, bestCB[c] ), CandidateRef( source, c          ) );
    } else {
      matchMapSrMt->insert( CandidateRef( source,  bestCB[c] ), CandidateRef( matched, c         ) );
      matchMapMtSr->insert( CandidateRef( matched, c         ), CandidateRef( source,  bestCB[c] ) );
    }
  }
*/
  evt.put( matchMapSrMt, "src2mtc" );
  evt.put( matchMapMtSr, "mtc2src" );

  AllDist.clear();
}

Member Data Documentation

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(), lenght(), and produce().

Definition at line 28 of file CandOneToOneDeltaRMatcher.cc.

Referenced by produce().

Definition at line 27 of file CandOneToOneDeltaRMatcher.cc.

Referenced by produce().