CMS 3D CMS Logo

PFConversionAlgo Class Reference

#include <RecoParticleFlow/PFAlgo/interface/PFConversionAlgo.h>

List of all members.

Public Member Functions

std::vector< reco::PFCandidateconversionCandidates ()
bool isConversionValidCandidate (const reco::PFBlockRef &blockRef, std::vector< bool > &active)
 PFConversionAlgo ()
 ~PFConversionAlgo ()

Private Types

typedef std::multimap
< unsigned, std::vector
< unsigned > > 
AssMap

Private Member Functions

void runPFConversion (const reco::PFBlockRef &blockRef, std::vector< bool > &active)
void setActive (const reco::PFBlockRef &blockRef, AssMap &assToConv, std::vector< bool > &active)
void setCandidates (const reco::PFBlockRef &blockref, AssMap &assToConv)
bool setLinks (const reco::PFBlockRef &blockRef, AssMap &assToConv, std::vector< bool > &active)

Private Attributes

std::vector< reco::PFCandidateconversionCandidate_
bool isvalid_


Detailed Description

Definition at line 11 of file PFConversionAlgo.h.


Member Typedef Documentation

typedef std::multimap<unsigned, std::vector<unsigned> > PFConversionAlgo::AssMap [private]

Definition at line 32 of file PFConversionAlgo.h.


Constructor & Destructor Documentation

PFConversionAlgo::PFConversionAlgo (  ) 

Definition at line 20 of file PFConversionAlgo.cc.

00020                                     {
00021 
00022 }

PFConversionAlgo::~PFConversionAlgo (  )  [inline]

Definition at line 20 of file PFConversionAlgo.h.

00020 {;};


Member Function Documentation

std::vector<reco::PFCandidate> PFConversionAlgo::conversionCandidates (  )  [inline]

Definition at line 32 of file PFConversionAlgo.h.

References conversionCandidate_.

00032 {return conversionCandidate_;};

bool PFConversionAlgo::isConversionValidCandidate ( const reco::PFBlockRef blockRef,
std::vector< bool > &  active 
) [inline]

Definition at line 23 of file PFConversionAlgo.h.

References isvalid_, and runPFConversion().

00025   {
00026     isvalid_=false;
00027     runPFConversion(blockRef,active);
00028     return isvalid_;
00029   };

void PFConversionAlgo::runPFConversion ( const reco::PFBlockRef blockRef,
std::vector< bool > &  active 
) [private]

Definition at line 24 of file PFConversionAlgo.cc.

References conversionCandidate_, reco::PFCandidate::gamma, reco::PFCandidate::GAMMA_TO_GAMMACONV, isvalid_, it, setActive(), setCandidates(), and setLinks().

Referenced by isConversionValidCandidate().

00025                                                                 {
00026 
00027   //  std::cout << " PFConversionAlgo::RunPFConversion " << std::endl;
00028 
00029   AssMap elemAssociatedToConv;
00030   
00031   bool blockHasConversion = setLinks(blockRef,elemAssociatedToConv, active );
00032   
00033   if (  blockHasConversion ) {
00034     conversionCandidate_.clear();
00035     setCandidates(blockRef,elemAssociatedToConv);
00036     if (conversionCandidate_.size() > 0 ){
00037       isvalid_ = true;
00038       //cout << " There is a candidate " << endl;
00039       // if there is at least a candidate the active vector is modified
00040       // setting = false all the elements used to build the candidate
00041       setActive(blockRef,elemAssociatedToConv, active);
00042       // this is just debug to check that all is going fine. Will take it out later      
00043       std::vector<reco::PFCandidate>::iterator it;
00044       for ( it = conversionCandidate_.begin(); it != conversionCandidate_.end(); ++it )  {
00045         reco::PFCandidate::ParticleType type = (*it).particleId();
00046         bool isConverted =  (*it).flag( reco::PFCandidate::GAMMA_TO_GAMMACONV );
00047         if ( type == reco::PFCandidate::gamma && isConverted ) {
00048           //      cout<<"Conversion PFCandidate!"<<  *it <<  endl;
00049         }
00050       }
00051     }
00052     
00053   } // conversion was found in the block
00054 
00055 }

void PFConversionAlgo::setActive ( const reco::PFBlockRef blockRef,
AssMap assToConv,
std::vector< bool > &  active 
) [private]

Definition at line 356 of file PFConversionAlgo.cc.

References i, j, and size.

Referenced by runPFConversion().

00357                                                                                           {
00358 
00359   // Lock tracks and clusters belonging to the conversion
00360   for (  std::multimap<unsigned, std::vector<unsigned> >::iterator i=elemAssociatedToConv.begin(); 
00361          i!=  elemAssociatedToConv.end(); ++i) {
00362     uint iConvTrk =  i->first;
00363     active[iConvTrk]=false;
00364     //std::cout << "  PFConversionAlgo::setActive locking all elements linked to a conversion " << std::endl;
00365     for (unsigned int j=0; j< (i->second).size(); ++j ) {
00366       active[(i->second)[j]]=false;
00367       //std::cout << " Locking element " << (i->second)[j] << std::endl;
00368     }
00369   }
00370 
00371   
00372   return;
00373 }

void PFConversionAlgo::setCandidates ( const reco::PFBlockRef blockref,
AssMap assToConv 
) [private]

Get the momentum of the parent track pair

for a first try just simple cuts

Build candidate

Definition at line 247 of file PFConversionAlgo.cc.

References funct::abs(), reco::PFCandidate::addElementInBlock(), parseConfig::block, reco::PFBlockElementCluster::clusterRef(), conversionCandidate_, deltaPhi(), reco::PFBlockElement::ECAL, bookConverter::elements, reco::PFBlock::elements(), reco::PFCluster::energy(), reco::PFCandidate::gamma, reco::PFCandidate::GAMMA_TO_GAMMACONV, i, j, reco::PFBlock::linkData(), pi, reco::PFCandidate::setConversionRef(), reco::PFCandidate::setEcalEnergy(), reco::PFCandidate::setFlag(), reco::PFCandidate::setHcalEnergy(), reco::PFCandidate::setPs1Energy(), reco::PFCandidate::setPs2Energy(), size, reco::PFBlockElement::TRACK, and reco::PFBlockElementTrack::trackRef().

Referenced by runPFConversion().

00248                                                                      {
00249 
00250   
00251   vector<uint> elementsToAdd(0);
00252   const reco::PFBlock& block = *blockRef;
00253   PFBlock::LinkData linkData =  block.linkData();     
00254   const edm::OwnVector< reco::PFBlockElement >&  elements = block.elements();
00255 
00258   float EcalEne=0;
00259   float pairPx=0;
00260   float pairPy=0;
00261   float pairPz=0;
00262   const reco::PFBlockElementTrack*  elTrack=0;
00263   reco::TrackRef convTrackRef;
00264   for (  std::multimap<unsigned, std::vector<unsigned> >::iterator i=elemAssociatedToConv.begin(); 
00265          i!=  elemAssociatedToConv.end(); ++i) {
00266 
00267     uint iTrack =  i->first;
00268     elementsToAdd.push_back(iTrack);
00269     //    std::cout << " setCandidates adding track " << iTrack << " to block in PFCandiate " << std::endl;
00270     elTrack = dynamic_cast<const reco::PFBlockElementTrack*>((&elements[iTrack]));
00271     convTrackRef= elTrack->trackRef();     
00272     pairPx+=convTrackRef->innerMomentum().x();
00273     pairPy+=convTrackRef->innerMomentum().y();
00274     pairPz+=convTrackRef->innerMomentum().z();
00275     
00276     ConversionRef origConv = elements[iTrack].convRef();
00277     //    std::cout << " Ref to original conversions: track size " << origConv->tracks().size() << " SC energy " << origConv->caloCluster()[0]->energy() << std::endl;
00278     // std::cout  << " SC Et " <<  origConv->caloCluster()[0]->energy()/cosh(origConv->caloCluster()[0]->eta()) <<" eta " << origConv->caloCluster()[0]->eta() << " phi " << origConv->caloCluster()[0]->phi() <<  std::endl;
00279 
00280     unsigned int nEl=  (i->second).size();
00281     //    std::cout << " Number of elements connected " << nEl << std::endl;
00282     for (unsigned int j=0; j< nEl; ++j ) {
00283       unsigned int iEl = (i->second)[j];
00284       //std::cout  << " Adding element " << iEl << std::endl;
00285       PFBlockElement::Type typeassCalo = elements[iEl].type();
00286 
00288       if ( typeassCalo == reco::PFBlockElement::TRACK) {
00289         const reco::PFBlockElementTrack * elTrack = dynamic_cast<const reco::PFBlockElementTrack*>((&elements[iEl]));   
00290         elTrack = dynamic_cast<const reco::PFBlockElementTrack*>((&elements[iEl]));
00291         convTrackRef= elTrack->trackRef();     
00292         pairPx+=convTrackRef->innerMomentum().x();
00293         pairPy+=convTrackRef->innerMomentum().y();
00294         pairPz+=convTrackRef->innerMomentum().z();
00295       }
00296 
00297       if ( typeassCalo == reco::PFBlockElement::ECAL) {
00298         const reco::PFBlockElementCluster * clu = dynamic_cast<const reco::PFBlockElementCluster*>((&elements[iEl]));
00299         reco::PFCluster cl=*clu->clusterRef();
00300         EcalEne+= cl.energy();
00301       }
00302 
00303       elementsToAdd.push_back(iEl);
00304     }
00305 
00306     //    std::cout << " setCandidates EcalEne " << EcalEne << std::endl;
00307     reco::PFCandidate::ParticleType particleType = reco::PFCandidate::gamma;
00308 
00309     // define energy and momentum of the conversion. Momentum from the track pairs and Energy from 
00310     // the sum of the ecal cluster(s)
00311     math::XYZTLorentzVector  momentum;
00312     momentum.SetPxPyPzE(pairPx , pairPy, pairPz, EcalEne);
00313 
00315 
00316     float deltaCotTheta=origConv->pairCotThetaSeparation();
00317     float phiTk1=  origConv->tracks()[0]->innerMomentum().phi();
00318     float phiTk2=  origConv->tracks()[1]->innerMomentum().phi();
00319     float deltaPhi = phiTk1-phiTk2;
00320     if(deltaPhi >  pi) {deltaPhi = deltaPhi - twopi;}
00321     if(deltaPhi < -pi) {deltaPhi = deltaPhi + twopi;}
00322 
00324     if (  fabs(deltaCotTheta) < 0.05 && abs(deltaPhi<0.1) )  {
00325       
00326       
00328       reco::PFCandidate aCandidate = PFCandidate(0, momentum,particleType);
00329       for (uint elad=0; elad<elementsToAdd.size();elad++){
00330         aCandidate.addElementInBlock(blockRef,elementsToAdd[elad]);
00331       }
00332       
00333       
00334       aCandidate.setFlag( reco::PFCandidate::GAMMA_TO_GAMMACONV, true);
00335       aCandidate.setConversionRef(origConv);
00336       aCandidate.setEcalEnergy(EcalEne);
00337       aCandidate.setHcalEnergy(0.); 
00338       aCandidate.setPs1Energy(0.); 
00339       aCandidate.setPs2Energy(0.); 
00340       
00341       
00342       conversionCandidate_.push_back( aCandidate);       
00343     }  
00344     
00345 
00346   }
00347 
00348 
00349 
00350   return;
00351 }

bool PFConversionAlgo::setLinks ( const reco::PFBlockRef blockRef,
AssMap assToConv,
std::vector< bool > &  active 
) [private]

Definition at line 58 of file PFConversionAlgo.cc.

References reco::PFBlock::associatedElements(), parseConfig::block, GenMuonPlsPt100GeV_cfg::cout, reco::PFBlockElement::ECAL, bookConverter::elements, reco::PFBlock::elements(), lat::endl(), find(), reco::PFBlockElement::HCAL, i, index, j, reco::PFBlock::linkData(), reco::PFBlock::LINKTEST_ALL, reco::PFBlock::LINKTEST_CHI2, edm::OwnVector< T, P >::push_back(), edm::OwnVector< T, P >::size(), size, and reco::PFBlockElement::TRACK.

Referenced by runPFConversion().

00060                                                           {
00061 
00062   bool conversionFound = false;
00063   typedef std::multimap<double, unsigned>::iterator IE;
00064   const reco::PFBlock& block = *blockRef;
00065   //  std::cout << " PFConversionAlgo::setLinks block " << block << std::endl;
00066   const edm::OwnVector< reco::PFBlockElement >&  elements = block.elements();
00067   PFBlock::LinkData linkData =  block.linkData();    
00068 
00069   // this method looks in all elements in a block, idenitifies those 
00070   // blonging to a conversions and stores them in a local map 
00071   // so that in the further code the search needs not be done again
00072 
00073   unsigned convTrack1Ind=100;
00074   unsigned convTrack2Ind=100;
00075   unsigned convEcal1Ind=200;
00076   unsigned convEcal2Ind=200;
00077   unsigned convHcal1Ind=200;
00078   unsigned convHcal2Ind=200;
00079 
00080   for(unsigned iElem=0; iElem<elements.size(); iElem++) {
00081     bool trackFromConv = elements[iElem].trackType( PFBlockElement::T_FROM_GAMMACONV);
00082 
00083     if (!active[iElem]) continue;  
00084     if ( !trackFromConv ) continue;
00085     conversionFound = true;    
00086 
00087     //    std::cout << " Track " <<  iElem << " is from conversion" << std::endl;
00088     convTrack1Ind= iElem;
00089 
00090     std::multimap<unsigned, std::vector<unsigned> >::iterator found = elemAssociatedToConv.find(iElem);
00091     if ( found!= elemAssociatedToConv.end()) {
00092       //      std::cout << " Track " <<  iElem << " has already been included " << std::endl;
00093       continue;
00094     }
00095     
00096     
00097     bool alreadyTaken=false;
00098     for (  std::multimap<unsigned, std::vector<unsigned> >::iterator i=elemAssociatedToConv.begin(); 
00099            i!=  elemAssociatedToConv.end(); ++i) {
00100       for (unsigned int j=0; j< (i->second).size(); ++j ) {
00101         if ( iElem == (i->second)[j] )  alreadyTaken=true;
00102       }
00103     }
00104     if ( alreadyTaken ) {
00105       //      std::cout << " iElem " << iElem << " already taken" <<  std::endl;
00106       continue;
00107     }
00108 
00109  
00110     vector<unsigned> assElements(0);     
00111     vector<unsigned>::iterator iVec;
00112 
00113     std::multimap<double, unsigned> ecalElems;
00114     block.associatedElements(iElem ,  linkData,
00115                               ecalElems ,
00116                               reco::PFBlockElement::ECAL );
00117 
00118     std::multimap<double, unsigned> hcalElems;
00119     block.associatedElements( iElem ,  linkData,
00120                               hcalElems,
00121                               reco::PFBlockElement::HCAL,
00122                               reco::PFBlock::LINKTEST_ALL );
00123     
00124     std::multimap<double, unsigned> trackElems;
00125     block.associatedElements( iElem,  linkData,
00126                               trackElems ,
00127                               reco::PFBlockElement::TRACK,
00128                               reco::PFBlock::LINKTEST_CHI2);
00129 
00130 
00131     if(trackElems.empty() ) {
00132       //      std::cout<<"PFConversionAlgo::setLinks no track element connected to track "<<iElem<<std::endl;
00133     }
00134     
00135     if(ecalElems.empty() ) {
00136       //  std::cout<<"PFConversionAlgo::setLinks no ecal element connected to track "<<iElem<<std::endl;
00137     }
00138     
00139     if(hcalElems.empty() ) {
00140       //  std::cout<<"PFConversionAlgo::setLinks no hcal element connected to track "<<iElem<<std::endl;
00141     }
00142     
00143     //std::cout<<"PFConversionAlgo::setLinks now looping on elements associated to the track"<<std::endl;
00144 
00145 
00146 
00147     //    std::cout<<"  look at linked hcal clusters"<<std::endl;
00148     for(IE iTk = hcalElems.begin(); iTk != hcalElems.end(); ++iTk ) {
00149       double   chi2  = iTk->first;
00150       unsigned index = iTk->second;
00151       PFBlockElement::Type type = elements[index].type();
00152       if ( type ==  reco::PFBlockElement::HCAL) {
00153         // link track-ecal is found 
00154         convHcal1Ind=index;
00155         //      std::cout << " Hcal-Track link found with " << convHcal1Ind << std::endl;
00156         //      if ( index< 100)  assElements.push_back(index);
00157       }
00158     }
00159 
00160 
00161 
00162     //    std::cout<<"  look at linked ecal clusters"<<std::endl;
00163     for(IE iTk = ecalElems.begin(); iTk != ecalElems.end(); ++iTk ) {
00164       double   chi2  = iTk->first;
00165       unsigned index = iTk->second;
00166       PFBlockElement::Type type = elements[index].type();
00167       if ( type ==  reco::PFBlockElement::ECAL) {
00168         // link track-ecal is found 
00169         convEcal1Ind=index;
00170         //std::cout << " Ecal-Track link found with " << convEcal1Ind << std::endl;
00171         iVec = find ( assElements.begin(), assElements.end(), index) ;        
00172         if ( index< 100 && iVec == assElements.end() )  assElements.push_back(index);
00173       }
00174     }
00175 
00176 
00177     //    std::cout<<"PFConversionAlgo::setLinks  look at linked tracks"<<std::endl;    
00178     for(IE iTk = trackElems.begin(); iTk != trackElems.end(); ++iTk ) {
00179       double   chi2  = iTk->first;
00180       unsigned index = iTk->second;
00181       PFBlockElement::Type type = elements[index].type();
00182         // link track-track is found 
00183         convTrack2Ind=index;
00184         if ( index< 100)  assElements.push_back(index);
00185         //std::cout << " Track-Track link found with " << convTrack2Ind << std::endl;
00186         std::multimap<double, unsigned> ecalElems2;
00187         block.associatedElements(convTrack2Ind ,  linkData,
00188                                  ecalElems2 ,
00189                                  reco::PFBlockElement::ECAL );
00190 
00191 
00192         for(IE iTk = ecalElems2.begin(); iTk != ecalElems2.end(); ++iTk ) {
00193           double   chi2  = iTk->first;
00194           unsigned index = iTk->second;
00195           PFBlockElement::Type type = elements[index].type();
00196           if ( type ==  reco::PFBlockElement::ECAL) {
00197             convEcal2Ind=index;
00198             //      std::cout << " 2nd ecal track link found betwtenn track " << convTrack2Ind  << " and Ecal " << convEcal2Ind << std::endl;
00199             iVec = find ( assElements.begin(), assElements.end(), index) ;        
00200             if ( index< 100 && iVec== assElements.end() )  assElements.push_back(index);
00201           
00202           }
00203         }
00204 
00205         std::multimap<double, unsigned> hcalElems2;
00206         block.associatedElements(convTrack2Ind ,  linkData,
00207                                  hcalElems2 ,
00208                                  reco::PFBlockElement::HCAL );
00209         for(IE iTk = hcalElems.begin(); iTk != hcalElems.end(); ++iTk ) {
00210           double   chi2  = iTk->first;
00211           unsigned index = iTk->second;
00212           PFBlockElement::Type type = elements[index].type();
00213           if ( type ==  reco::PFBlockElement::HCAL) {
00214             // link track-ecal is found 
00215             convHcal2Ind=index;
00216             std::cout << " Hcal-Track link found with " << convHcal2Ind << std::endl;
00217             //   if ( index< 100)  assElements.push_back(index);
00218           }
00219         }
00220 
00221     }
00222 
00223 
00224     elemAssociatedToConv.insert(make_pair(convTrack1Ind, assElements));
00225 
00226     //std::cout << " PFConversionAlgo::setLink map size " << elemAssociatedToConv.size() << std::endl;
00227     for (  std::multimap<unsigned, std::vector<unsigned> >::iterator i=elemAssociatedToConv.begin(); 
00228            i!=  elemAssociatedToConv.end(); ++i) {
00229       //std::cout << " links found for " << i->first << std::endl;
00230       //std::cout << " elements " << (i->second).size() << std::endl;
00231       for (unsigned int j=0; j< (i->second).size(); ++j ) {
00232         unsigned int iEl = (i->second)[j];
00233         //std::cout  << " ass element " << iEl << std::endl;
00234       }
00235     }
00236 
00237 
00238 
00239   } // end loop over elements of the block looking for conversion track
00240 
00241 
00242   return conversionFound;
00243 
00244 }


Member Data Documentation

std::vector<reco::PFCandidate> PFConversionAlgo::conversionCandidate_ [private]

Definition at line 49 of file PFConversionAlgo.h.

Referenced by conversionCandidates(), runPFConversion(), and setCandidates().

bool PFConversionAlgo::isvalid_ [private]

Definition at line 50 of file PFConversionAlgo.h.

Referenced by isConversionValidCandidate(), and runPFConversion().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:29:42 2009 for CMSSW by  doxygen 1.5.4