CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoParticleFlow/PFProducer/interface/PFBlockAlgo.h

Go to the documentation of this file.
00001 #ifndef RecoParticleFlow_PFProducer_PFBlockAlgo_h
00002 #define RecoParticleFlow_PFProducer_PFBlockAlgo_h 
00003 
00004 #include <set>
00005 #include <vector>
00006 #include <iostream>
00007 
00008 // #include "FWCore/Framework/interface/Handle.h"
00009 #include "DataFormats/Common/interface/Handle.h"
00010 // #include "FWCore/Framework/interface/OrphanHandle.h"
00011 #include "DataFormats/Common/interface/OrphanHandle.h"
00012 
00013 
00014 #include "DataFormats/ParticleFlowReco/interface/PFClusterFwd.h"
00015 #include "DataFormats/ParticleFlowReco/interface/PFRecTrackFwd.h"
00016 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
00017 #include "DataFormats/ParticleFlowReco/interface/PFRecTrack.h"
00018 #include "DataFormats/ParticleFlowReco/interface/PFDisplacedTrackerVertex.h" // gouzevitch
00019 #include "DataFormats/ParticleFlowReco/interface/PFConversionFwd.h"
00020 #include "DataFormats/ParticleFlowReco/interface/PFConversion.h"
00021 #include "DataFormats/ParticleFlowReco/interface/PFV0Fwd.h"
00022 #include "DataFormats/ParticleFlowReco/interface/PFV0.h"
00023 #include "DataFormats/ParticleFlowReco/interface/GsfPFRecTrack.h"
00024 #include "DataFormats/ParticleFlowReco/interface/GsfPFRecTrackFwd.h"
00025 #include "DataFormats/ParticleFlowReco/interface/PFBrem.h"
00026 #include "DataFormats/ParticleFlowReco/interface/PFTrajectoryPoint.h"
00027 
00028 #include "DataFormats/ParticleFlowReco/interface/PFBlockElement.h"
00029 #include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
00030 #include "DataFormats/ParticleFlowReco/interface/PFBlockFwd.h"
00031 
00032 // #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00033 
00034 #include "RecoParticleFlow/PFProducer/interface/PFMuonAlgo.h"
00035 #include "RecoParticleFlow/PFProducer/interface/PhotonSelectorAlgo.h"
00036 #include "RecoParticleFlow/PFProducer/interface/PFBlockElementSCEqual.h"
00037 #include "RecoParticleFlow/PFClusterTools/interface/PFResolutionMap.h"
00038 
00039 #include "DataFormats/MuonReco/interface/MuonFwd.h"
00040 #include "DataFormats/MuonReco/interface/Muon.h"
00041 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
00042 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00043 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
00044 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00045 #include "DataFormats/EgammaReco/interface/ElectronSeedFwd.h"
00046 #include "DataFormats/EgammaReco/interface/ElectronSeed.h"
00047 #include "DataFormats/EgammaCandidates/interface/PhotonFwd.h"
00048 #include "DataFormats/EgammaCandidates/interface/Photon.h"
00049 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementSuperCluster.h"
00050 
00051 #include "RecoParticleFlow/PFClusterTools/interface/ClusterClusterMapping.h"
00052 
00053 #include "RecoParticleFlow/PFProducer/interface/PFBlockLink.h"
00054 
00055 #include <map>
00056 
00057 
00059 
00064 class PFBlockAlgo {
00065 
00066  public:
00067 
00068   PFBlockAlgo();
00069 
00070   ~PFBlockAlgo();
00071   
00072 
00073   void setParameters( std::vector<double>& DPtovPtCut, 
00074                       std::vector<unsigned>& NHitCut,
00075                       bool useConvBremPFRecTracks,
00076                       bool useIterTracking,
00077                       int nuclearInteractionsPurity,
00078                       bool useEGPhotons,
00079                       std::vector<double> & photonSelectionCuts
00080                       );
00081   
00082   typedef std::vector<bool> Mask;
00083 
00085   template< template<typename> class T>
00086     void  setInput(const T<reco::PFRecTrackCollection>&    trackh,
00087                    const T<reco::GsfPFRecTrackCollection>&    gsftrackh,
00088                    const T<reco::GsfPFRecTrackCollection>&    convbremgsftrackh,
00089                    const T<reco::MuonCollection>&    muonh,
00090                    const T<reco::PFDisplacedTrackerVertexCollection>&  nuclearh,
00091                    const T<reco::PFRecTrackCollection>&    nucleartrackh,
00092                    const T<reco::PFConversionCollection>&  conv,
00093                    const T<reco::PFV0Collection>&  v0,
00094                    const T<reco::PFClusterCollection>&  ecalh,
00095                    const T<reco::PFClusterCollection>&  hcalh,
00096                    const T<reco::PFClusterCollection>&  hfemh,
00097                    const T<reco::PFClusterCollection>&  hfhadh,
00098                    const T<reco::PFClusterCollection>&  psh,
00099                    const T<reco::PhotonCollection>&  egphh,
00100                    const Mask& trackMask = dummyMask_,
00101                    const Mask& gsftrackMask = dummyMask_,
00102                    const Mask& ecalMask = dummyMask_,
00103                    const Mask& hcalMask = dummyMask_,
00104                    const Mask& hfemMask = dummyMask_,              
00105                    const Mask& hfhadMask = dummyMask_,
00106                    const Mask& psMask = dummyMask_,
00107                    const Mask& phMask = dummyMask_); 
00108   
00110   template< template<typename> class T >
00111     void setInput(const T<reco::PFRecTrackCollection>&    trackh,
00112                   const T<reco::MuonCollection>&    muonh,
00113                   const T<reco::PFClusterCollection>&  ecalh,
00114                   const T<reco::PFClusterCollection>&  hcalh,
00115                   const T<reco::PFClusterCollection>&  hfemh,
00116                   const T<reco::PFClusterCollection>&  hfhadh,
00117                   const T<reco::PFClusterCollection>&  psh,
00118                   const Mask& trackMask = dummyMask_,
00119                   const Mask& ecalMask = dummyMask_,
00120                   const Mask& hcalMask = dummyMask_,
00121                   const Mask& psMask = dummyMask_ ) {
00122     T<reco::GsfPFRecTrackCollection> gsftrackh;
00123     T<reco::GsfPFRecTrackCollection> convbremgsftrackh;
00124     //T<reco::MuonCollection> muonh;
00125     T<reco::PFDisplacedTrackerVertexCollection> nuclearh;
00126     T<reco::PFRecTrackCollection>    nucleartrackh;
00127     T<reco::PFConversionCollection> convh;
00128     T<reco::PFV0Collection> v0;
00129     T<reco::PhotonCollection> phh;
00130     setInput<T>( trackh, gsftrackh, convbremgsftrackh, muonh, nuclearh, nucleartrackh, convh, v0, 
00131                  ecalh, hcalh, hfemh, hfhadh, psh, phh,
00132                  trackMask, ecalMask, hcalMask, psMask); 
00133   }
00134   
00136   template< template<typename> class T >
00137     void setInput(const T<reco::PFRecTrackCollection>&    trackh,
00138                   const T<reco::GsfPFRecTrackCollection>&    gsftrackh,
00139                   const T<reco::PFClusterCollection>&  ecalh,
00140                   const T<reco::PFClusterCollection>&  hcalh,
00141                   const T<reco::PFClusterCollection>&  psh,
00142                   const Mask& trackMask = dummyMask_,
00143                   const Mask& gsftrackMask = dummyMask_,
00144                   const Mask& ecalMask = dummyMask_,
00145                   const Mask& hcalMask = dummyMask_,
00146                   const Mask& psMask = dummyMask_ ) {
00147     T<reco::GsfPFRecTrackCollection> convbremgsftrackh;
00148     T<reco::MuonCollection> muonh;
00149     T<reco::PFDisplacedTrackerVertexCollection>  nuclearh;
00150     T<reco::PFRecTrackCollection>    nucleartrackh;
00151     T<reco::PFConversionCollection> convh;
00152     T<reco::PFV0Collection> v0;
00153     T<reco::PhotonCollection> egphh;
00154     setInput<T>( trackh, gsftrackh, convbremgsftrackh, muonh, nuclearh, nucleartrackh, convh, v0, ecalh, hcalh, psh, egphh,
00155                  trackMask, gsftrackMask,ecalMask, hcalMask, psMask); 
00156   }
00157   
00158   
00160   void setDebug( bool debug ) {debug_ = debug;}
00161   
00163   void findBlocks();
00164   
00165 
00167   /*   const  reco::PFBlockCollection& blocks() const {return *blocks_;} */
00168   const std::auto_ptr< reco::PFBlockCollection >& blocks() const 
00169     {return blocks_;}
00170   
00172   std::auto_ptr< reco::PFBlockCollection > transferBlocks() {return blocks_;}
00173 
00175   typedef std::list< reco::PFBlockElement* >::iterator IE;
00176   typedef std::list< reco::PFBlockElement* >::const_iterator IEC;  
00177   typedef reco::PFBlockCollection::const_iterator IBC;
00178   
00179  private:
00180   
00187   IE associate(IE next, IE last, std::vector<PFBlockLink>& links);
00188 
00191   void packLinks(reco::PFBlock& block, 
00192                  const std::vector<PFBlockLink>& links) const; 
00193 
00195   void checkDisplacedVertexLinks( reco::PFBlock& block ) const;
00196   
00200   void buildGraph(); 
00201 
00203   void link( const reco::PFBlockElement* el1, 
00204              const reco::PFBlockElement* el2, 
00205              PFBlockLink::Type& linktype, 
00206              reco::PFBlock::LinkTest& linktest,
00207              double& dist) const;
00208                  
00211   double testTrackAndPS(const reco::PFRecTrack& track,
00212                         const reco::PFCluster& ps) const;
00213   
00216   double testECALAndHCAL(const reco::PFCluster& ecal, 
00217                          const reco::PFCluster& hcal) const;
00218                          
00221   double testPS1AndPS2(const reco::PFCluster& ps1,
00222                        const reco::PFCluster& ps2) const;
00223 
00225   double testLinkBySuperCluster(const reco::PFClusterRef & elt1,
00226                                 const reco::PFClusterRef & elt2) const;   
00227 
00229   double testSuperClusterPFCluster(const reco::SuperClusterRef & sct1,
00230                                    const reco::PFClusterRef & elt2) const;
00231 
00235   void checkMaskSize( const reco::PFRecTrackCollection& tracks,
00236                       const reco::GsfPFRecTrackCollection& gsftracks,
00237                       const reco::PFClusterCollection&  ecals,
00238                       const reco::PFClusterCollection&  hcals,
00239                       const reco::PFClusterCollection&  hfems,
00240                       const reco::PFClusterCollection&  hfhads,
00241                       const reco::PFClusterCollection&  pss, 
00242                       const reco::PhotonCollection&  egphh, 
00243                       const Mask& trackMask,
00244                       const Mask& gsftrackMask, 
00245                       const Mask& ecalMask, 
00246                       const Mask& hcalMask,
00247                       const Mask& hfemMask,
00248                       const Mask& hfhadMask,                  
00249                       const Mask& psMask,
00250                       const Mask& phMask) const;
00251 
00253   // PFResolutionMap* openResolutionMap(const char* resMapName);
00254 
00256   bool goodPtResolution( const reco::TrackRef& trackref);
00257 
00258   double testLinkByVertex(const reco::PFBlockElement* elt1,
00259                           const reco::PFBlockElement* elt2) const;
00260 
00261   int muAssocToTrack( const reco::TrackRef& trackref,
00262                       const edm::Handle<reco::MuonCollection>& muonh) const;
00263   int muAssocToTrack( const reco::TrackRef& trackref,
00264                       const edm::OrphanHandle<reco::MuonCollection>& muonh) const;
00265 
00266   void fillFromPhoton(const reco::Photon & photon, reco::PFBlockElementSuperCluster * pfbe);
00267 
00268   std::auto_ptr< reco::PFBlockCollection >    blocks_;
00269   
00271   // std::vector< reco::PFCandidate >   particles_;
00272 
00273   // the test elements will be transferred to the blocks
00274   std::list< reco::PFBlockElement* >     elements_;
00275 
00276   static const Mask                      dummyMask_;
00277 
00279   std::vector<double> DPtovPtCut_;
00280   
00282   std::vector<unsigned> NHitCut_;
00283   
00285   bool useIterTracking_;
00286 
00288   bool useEGPhotons_;
00289 
00290   // This parameters defines the level of purity of
00291   // nuclear interactions choosen.
00292   // Level 1 is only high Purity sample labeled as isNucl
00293   // Level 2 isNucl + isNucl_Loose (2 secondary tracks vertices)
00294   // Level 3 isNucl + isNucl_Loose + isNucl_Kink
00295   //         (low purity sample made of 1 primary and 1 secondary track)
00296   // By default the level 1 is teh safest one.
00297   int nuclearInteractionsPurity_;
00298 
00300   bool  useConvBremPFRecTracks_;
00301 
00303   const PhotonSelectorAlgo * photonSelector_;
00305   std::vector<reco::SuperClusterRef > superClusters_;
00306 
00308   //  std::map<reco::PFClusterRef,int>  pfcRefSCMap_;
00309   std::vector<int> pfcSCVec_;
00310 
00312   std::vector<std::vector<reco::PFClusterRef> > scpfcRefs_;
00314   bool   debug_;
00315   
00316   friend std::ostream& operator<<(std::ostream&, const PFBlockAlgo&);
00317 
00318 };
00319 
00320 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementGsfTrack.h"
00321 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementBrem.h"
00322 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementTrack.h"
00323 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementCluster.h"
00324 #include "DataFormats/ParticleFlowReco/interface/PFLayer.h"
00325 
00326 template< template<typename> class T >
00327 void
00328 PFBlockAlgo::setInput(const T<reco::PFRecTrackCollection>&    trackh,
00329                       const T<reco::GsfPFRecTrackCollection>&    gsftrackh, 
00330                       const T<reco::GsfPFRecTrackCollection>&    convbremgsftrackh,
00331                       const T<reco::MuonCollection>&    muonh,
00332                       const T<reco::PFDisplacedTrackerVertexCollection>&  nuclearh,
00333                       const T<reco::PFRecTrackCollection>&    nucleartrackh,
00334                       const T<reco::PFConversionCollection>&  convh,
00335                       const T<reco::PFV0Collection>&  v0,
00336                       const T<reco::PFClusterCollection>&  ecalh,
00337                       const T<reco::PFClusterCollection>&  hcalh,
00338                       const T<reco::PFClusterCollection>&  hfemh,
00339                       const T<reco::PFClusterCollection>&  hfhadh,
00340                       const T<reco::PFClusterCollection>&  psh,
00341                       const T<reco::PhotonCollection>& egphh,
00342                       const Mask& trackMask,
00343                       const Mask& gsftrackMask,
00344                       const Mask& ecalMask,
00345                       const Mask& hcalMask,
00346                       const Mask& hfemMask,
00347                       const Mask& hfhadMask,
00348                       const Mask& psMask,
00349                       const Mask& phMask) {
00350 
00351 
00352   checkMaskSize( *trackh,
00353                  *gsftrackh,
00354                  *ecalh,
00355                  *hcalh,
00356                  *hfemh,
00357                  *hfhadh,
00358                  *psh,
00359                  *egphh,
00360                  trackMask,
00361                  gsftrackMask,
00362                  ecalMask,
00363                  hcalMask,
00364                  hfemMask,
00365                  hfhadMask,
00366                  psMask,
00367                  phMask);
00368 
00369   /*
00370   if (nucleartrackh.isValid()){
00371     for(unsigned i=0;i<nucleartrackh->size(); i++) {
00372       reco::PFRecTrackRef trackRef(nucleartrackh,i);
00373       std::cout << *trackRef << std::endl;
00374     }
00375   }
00376   */
00377 
00379   std::vector<reco::PFRecTrackRef> convBremPFRecTracks;
00380   convBremPFRecTracks.clear();
00381   // Super cluster mapping
00382   superClusters_.clear();
00383   scpfcRefs_.clear();
00384   pfcSCVec_.clear();
00385 
00386   if(gsftrackh.isValid() ) {
00387     const  reco::GsfPFRecTrackCollection PFGsfProd = *(gsftrackh.product());
00388     for(unsigned i=0;i<gsftrackh->size(); i++) {
00389       if( !gsftrackMask.empty() &&
00390           !gsftrackMask[i] ) continue;
00391       reco::GsfPFRecTrackRef refgsf(gsftrackh,i );   
00392    
00393       if((refgsf).isNull()) continue;
00394       reco::GsfTrackRef gsf=refgsf->gsfTrackRef();
00395 
00396       // retrieve and save the SC if ECAL-driven - Florian
00397       if(gsf->extra().isAvailable() && gsf->extra()->seedRef().isAvailable()) {
00398         reco::ElectronSeedRef seedRef=  gsf->extra()->seedRef().castTo<reco::ElectronSeedRef>();
00399         // check that the seed is valid
00400         if(seedRef.isAvailable() && seedRef->isEcalDriven()) {
00401           reco::SuperClusterRef scRef = seedRef->caloCluster().castTo<reco::SuperClusterRef>();
00402           if(scRef.isNonnull())   {          
00403             std::vector<reco::SuperClusterRef>::iterator itcheck=find(superClusters_.begin(),superClusters_.end(),scRef);
00404             // not present, add it
00405             if(itcheck==superClusters_.end())
00406               {
00407                 superClusters_.push_back(scRef);
00408                 reco::PFBlockElementSuperCluster * sce =
00409                   new reco::PFBlockElementSuperCluster(scRef);
00410                 sce->setFromGsfElectron(true);
00411                 elements_.push_back(sce);
00412               }
00413             else // it is already present, update the PFBE
00414               {
00415                 PFBlockElementSCEqual myEqual(scRef);
00416                 std::list<reco::PFBlockElement*>::iterator itcheck=find_if(elements_.begin(),elements_.end(),myEqual);
00417                 if(itcheck!=elements_.end())
00418                   {
00419                     reco::PFBlockElementSuperCluster* thePFBE=dynamic_cast<reco::PFBlockElementSuperCluster*>(*itcheck);
00420                     thePFBE->setFromGsfElectron(true);
00421                     //              std::cout << " Updating element to add electron information" << std::endl;
00422                   }
00423 //              else
00424 //                {
00425 //                  std::cout << " Missing element " << std::endl;
00426 //                }
00427               }
00428           }
00429         }
00430       }
00431 
00432       reco::PFBlockElement* gsfEl;
00433       
00434       const  std::vector<reco::PFTrajectoryPoint> 
00435         PfGsfPoint =  PFGsfProd[i].trajectoryPoints();
00436   
00437       unsigned int c_gsf=0;
00438       bool PassTracker = false;
00439       bool GetPout = false;
00440       unsigned int IndexPout = 0;
00441       
00442       typedef std::vector<reco::PFTrajectoryPoint>::const_iterator IP;
00443       for(IP itPfGsfPoint =  PfGsfPoint.begin();  
00444           itPfGsfPoint!= PfGsfPoint.end();itPfGsfPoint++) {
00445         
00446         if (itPfGsfPoint->isValid()){
00447           int layGsfP = itPfGsfPoint->layer();
00448           if (layGsfP == -1) PassTracker = true;
00449           if (PassTracker && layGsfP > 0 && GetPout == false) {
00450             IndexPout = c_gsf-1;
00451             GetPout = true;
00452           }
00453           //const math::XYZTLorentzVector GsfMoment = itPfGsfPoint->momentum();
00454           c_gsf++;
00455         }
00456       }
00457       math::XYZTLorentzVector pin = PfGsfPoint[0].momentum();      
00458       math::XYZTLorentzVector pout = PfGsfPoint[IndexPout].momentum();
00459 
00461       if(useConvBremPFRecTracks_) {
00462         const std::vector<reco::PFRecTrackRef>& temp_convBremPFRecTracks(refgsf->convBremPFRecTrackRef());
00463         if(temp_convBremPFRecTracks.size() > 0) {
00464           for(unsigned int iconv = 0; iconv <temp_convBremPFRecTracks.size(); iconv++) {
00465             convBremPFRecTracks.push_back(temp_convBremPFRecTracks[iconv]);
00466           }
00467         }
00468       }
00469       
00470       gsfEl = new reco::PFBlockElementGsfTrack(refgsf, pin, pout);
00471       
00472       elements_.push_back( gsfEl);
00473 
00474       std::vector<reco::PFBrem> pfbrem = refgsf->PFRecBrem();
00475       
00476       for (unsigned i2=0;i2<pfbrem.size(); i2++) {
00477         const double DP = pfbrem[i2].DeltaP();
00478         const double SigmaDP =  pfbrem[i2].SigmaDeltaP(); 
00479         const unsigned int TrajP = pfbrem[i2].indTrajPoint();
00480         if(TrajP == 99) continue;
00481 
00482         reco::PFBlockElement* bremEl;
00483         bremEl = new reco::PFBlockElementBrem(refgsf,DP,SigmaDP,TrajP);
00484         elements_.push_back(bremEl);
00485         
00486       }
00487     }
00488 
00489   }
00490 
00492   if(useEGPhotons_ && egphh.isValid()) {
00493     unsigned size=egphh->size();
00494     for(unsigned isc=0; isc<size; ++isc) {
00495       if(!phMask.empty() && !(phMask)[isc] ) continue;
00496       if(!photonSelector_->passPhotonSelection((*egphh)[isc])) continue;
00497       //      std::cout << " Selected a supercluster" << std::endl;
00498       // Add only the super clusters not already included 
00499       reco::SuperClusterRef scRef((*egphh)[isc].superCluster());
00500       std::vector<reco::SuperClusterRef>::iterator itcheck=find(superClusters_.begin(),superClusters_.end(),(*egphh)[isc].superCluster());
00501       if(itcheck==superClusters_.end())
00502         {
00503           superClusters_.push_back(scRef);
00504           reco::PFBlockElementSuperCluster * sce =
00505             new reco::PFBlockElementSuperCluster((*egphh)[isc].superCluster());
00506           fillFromPhoton((*egphh)[isc],sce);
00507           elements_.push_back(sce);
00508         }
00509       else
00510         {
00511           PFBlockElementSCEqual myEqual(scRef);
00512           std::list<reco::PFBlockElement*>::iterator itcheck=find_if(elements_.begin(),elements_.end(),myEqual);
00513           if(itcheck!=elements_.end())
00514             {
00515               reco::PFBlockElementSuperCluster* thePFBE=dynamic_cast<reco::PFBlockElementSuperCluster*>(*itcheck);
00516               fillFromPhoton((*egphh)[isc],thePFBE);
00517               thePFBE->setFromPhoton(true);
00518               //              std::cout << " Updating element to add Photon information " << photonSelector_->passPhotonSelection((*egphh)[isc]) << std::endl;
00519 
00520             }
00521 //        else
00522 //          {
00523 //            std::cout << " Missing element " << std::endl;
00524 //          }
00525         }
00526     }
00527   }
00528   
00529   // set the vector to the right size so to allow random access
00530   scpfcRefs_.resize(superClusters_.size());
00531 
00533 
00535 
00536   if(convh.isValid() ) {
00537     reco::PFBlockElement* trkFromConversionElement;
00538     for(unsigned i=0;i<convh->size(); i++) {
00539       reco::PFConversionRef convRef(convh,i);
00540 
00541       unsigned int trackSize=(convRef->pfTracks()).size();
00542       if ( convRef->pfTracks().size() < 2) continue;
00543       for(unsigned iTk=0;iTk<trackSize; iTk++) {
00544         
00545         reco::PFRecTrackRef compPFTkRef = convRef->pfTracks()[iTk];     
00546         trkFromConversionElement = new reco::PFBlockElementTrack(convRef->pfTracks()[iTk]);
00547         trkFromConversionElement->setConversionRef( convRef->originalConversion(), reco::PFBlockElement::T_FROM_GAMMACONV);
00548 
00549         elements_.push_back( trkFromConversionElement );
00550 
00551 
00552         if (debug_){
00553           std::cout << "PF Block Element from Conversion electron " << 
00554             (*trkFromConversionElement).trackRef().key() << std::endl;
00555           std::cout << *trkFromConversionElement << std::endl;
00556         }
00557        
00558       }     
00559     }  
00560   }
00561   
00562   
00564   
00566   
00567   if(v0.isValid() ) {
00568     reco::PFBlockElement* trkFromV0Element = 0;
00569     for(unsigned i=0;i<v0->size(); i++) {
00570       reco::PFV0Ref v0Ref( v0, i );
00571       unsigned int trackSize=(v0Ref->pfTracks()).size();
00572       for(unsigned iTk=0;iTk<trackSize; iTk++) {
00573 
00574         reco::PFRecTrackRef newPFRecTrackRef = (v0Ref->pfTracks())[iTk]; 
00575         reco::TrackBaseRef newTrackBaseRef(newPFRecTrackRef->trackRef());
00576         bool bNew = true;
00577         
00580         for(IE iel = elements_.begin(); iel != elements_.end(); iel++){
00581           reco::TrackBaseRef elemTrackBaseRef((*iel)->trackRef());
00582           if (newTrackBaseRef == elemTrackBaseRef){         
00583             trkFromV0Element = *iel;
00584             bNew = false;
00585             continue;
00586           }
00587         } 
00588 
00590         if (bNew) {
00591           trkFromV0Element = new reco::PFBlockElementTrack(v0Ref->pfTracks()[iTk]);
00592           elements_.push_back( trkFromV0Element );
00593         }
00594 
00595         trkFromV0Element->setV0Ref( v0Ref->originalV0(),
00596                                     reco::PFBlockElement::T_FROM_V0 );
00597           
00598         if (debug_){
00599           std::cout << "PF Block Element from V0 track New = " << bNew 
00600                     << (*trkFromV0Element).trackRef().key() << std::endl;
00601           std::cout << *trkFromV0Element << std::endl;
00602         }
00603 
00604         
00605       }
00606     }
00607   }
00608   
00610 
00612 
00613   if(nuclearh.isValid()) {
00614     reco::PFBlockElement* trkFromDisplacedVertexElement = 0;
00615     for(unsigned i=0;i<nuclearh->size(); i++) {
00616 
00617       const reco::PFDisplacedTrackerVertexRef dispacedVertexRef( nuclearh, i );
00618 
00619       //      std::cout << "Nuclear Interactions Purity " <<  nuclearInteractionsPurity_ << std::endl;
00620       //     dispacedVertexRef->displacedVertexRef()->Dump();
00621       //bool bIncludeVertices = true;
00622 
00623       
00624       bool bIncludeVertices = false; 
00625       bool bNucl = dispacedVertexRef->displacedVertexRef()->isNucl();
00626       bool bNucl_Loose = dispacedVertexRef->displacedVertexRef()->isNucl_Loose();
00627       bool bNucl_Kink = dispacedVertexRef->displacedVertexRef()->isNucl_Kink();
00628 
00629       if (nuclearInteractionsPurity_ >= 1) bIncludeVertices = bNucl;
00630       if (nuclearInteractionsPurity_ >= 2) bIncludeVertices = bIncludeVertices || bNucl_Loose;
00631       if (nuclearInteractionsPurity_ >= 3) bIncludeVertices = bIncludeVertices || bNucl_Kink;
00632 
00633       if (bIncludeVertices){
00634 
00635         unsigned int trackSize= dispacedVertexRef->pfRecTracks().size();
00636         if (debug_){
00637           std::cout << "" << std::endl;
00638           std::cout << "Displaced Vertex " << i << std::endl;
00639           dispacedVertexRef->displacedVertexRef()->Dump();
00640         }
00641         for(unsigned iTk=0;iTk < trackSize; iTk++) {
00642 
00643 
00644           // This peace of code looks weired at first but it seems to be necessary to let 
00645           // PFRooTEvent work properly. Since the track position called REPPoint is transient 
00646           // it has to be calculated and the PFRecTrack collection associted to Displaced Vertex is not
00647           // anymore the original collection. So here we match both collections if possible
00648           reco::PFRecTrackRef newPFRecTrackRef = dispacedVertexRef->pfRecTracks()[iTk]; 
00649           reco::TrackBaseRef constTrackBaseRef(newPFRecTrackRef->trackRef());
00650 
00651           
00652           if (nucleartrackh.isValid()){
00653             for(unsigned i=0;i<nucleartrackh->size(); i++) {
00654               reco::PFRecTrackRef transientPFRecTrackRef(nucleartrackh,i);
00655               reco::TrackBaseRef transientTrackBaseRef(transientPFRecTrackRef->trackRef());
00656               if (constTrackBaseRef==transientTrackBaseRef){
00657                 newPFRecTrackRef = transientPFRecTrackRef;
00658                 break;
00659               }
00660             }
00661           }
00662           reco::TrackBaseRef newTrackBaseRef(newPFRecTrackRef->trackRef());
00663           
00664 
00665 
00666           bool bNew = true;
00667           reco::PFBlockElement::TrackType blockType;
00668 
00671           for(IE iel = elements_.begin(); iel != elements_.end(); iel++){
00672             reco::TrackBaseRef elemTrackBaseRef((*iel)->trackRef());
00673             if (newTrackBaseRef == elemTrackBaseRef){
00674               trkFromDisplacedVertexElement = *iel;
00675               bNew = false;
00676               continue;
00677             }
00678           }
00679 
00680 
00682           if (bNew) { 
00683             
00684 
00685 
00686             trkFromDisplacedVertexElement = new reco::PFBlockElementTrack(newPFRecTrackRef);
00687             elements_.push_back( trkFromDisplacedVertexElement );
00688           }
00689 
00690           if (dispacedVertexRef->isIncomingTrack(newPFRecTrackRef)) 
00691             blockType = reco::PFBlockElement::T_TO_DISP;
00692           else if (dispacedVertexRef->isOutgoingTrack(newPFRecTrackRef)) 
00693             blockType = reco::PFBlockElement::T_FROM_DISP;
00694           else 
00695             blockType = reco::PFBlockElement::DEFAULT;
00696 
00698           trkFromDisplacedVertexElement->setDisplacedVertexRef( dispacedVertexRef, blockType );
00699 
00700 
00701           if (debug_){
00702             std::cout << "PF Block Element from DisplacedTrackingVertex track New = " << bNew
00703                       << (*trkFromDisplacedVertexElement).trackRef().key() << std::endl;
00704             std::cout << *trkFromDisplacedVertexElement << std::endl;
00705           }
00706         
00707         
00708         }
00709       }
00710     }  
00711 
00712     if (debug_) std::cout << "" << std::endl;
00713 
00714   }
00715 
00717 
00721 
00722   if(trackh.isValid() ) {
00723 
00724     if (debug_) std::cout << "Tracks already in from Displaced Vertices " << std::endl;
00725 
00726     Mask trackMaskVertex;
00727 
00728     for(unsigned i=0;i<trackh->size(); i++) {
00729       reco::PFRecTrackRef pfRefTrack( trackh,i );
00730       reco::TrackRef trackRef = pfRefTrack->trackRef();
00731 
00732       bool bMask = true;
00733       for(IE iel = elements_.begin(); iel != elements_.end(); iel++){
00734         reco::TrackRef elemTrackRef = (*iel)->trackRef();
00735         if( trackRef == elemTrackRef ) {
00736           if (debug_) std::cout << " " << trackRef.key();
00737           bMask = false; continue;
00738         }
00739       }
00740     
00741       trackMaskVertex.push_back(bMask);
00742     }
00743 
00744     if (debug_) std::cout << "" << std::endl;
00745 
00746     if (debug_) std::cout << "Additionnal tracks from main collection " << std::endl;
00747 
00748     for(unsigned i=0;i<trackh->size(); i++) {
00749 
00750 
00751       // this track has been disabled
00752       if( !trackMask.empty() && !trackMask[i] ) continue;
00753       
00754       reco::PFRecTrackRef ref( trackh,i );
00755 
00756       if (debug_) std::cout << " " << ref->trackRef().key();
00757 
00758       // Get the eventual muon associated to this track
00759       int muId_ = muAssocToTrack( ref->trackRef(), muonh );
00760       bool thisIsAPotentialMuon = false;
00761       if( muId_ != -1 ) {
00762         reco::MuonRef muonref( muonh, muId_ );
00763         thisIsAPotentialMuon = 
00764           PFMuonAlgo::isLooseMuon(muonref) || 
00765           PFMuonAlgo::isMuon(muonref);
00766       }
00767       // Reject bad tracks (except if identified as muon
00768       if( !thisIsAPotentialMuon && !goodPtResolution( ref->trackRef() ) ) continue;
00769 
00770       if (thisIsAPotentialMuon && debug_) std::cout << "Potential Muon P " <<  ref->trackRef()->p() 
00771                                                     << " pt " << ref->trackRef()->p() << std::endl; 
00772 
00773 
00774 
00775       reco::PFBlockElement* primaryElement = new reco::PFBlockElementTrack( ref );
00776 
00777       if( muId_ != -1 ) {
00778         // if a muon has been found
00779         reco::MuonRef muonref( muonh, muId_ );
00780 
00781         // If this track was already added to the collection, we just need to find the associated element and 
00782         // attach to it the reference
00783         if (!trackMaskVertex.empty() && !trackMaskVertex[i]){
00784           reco::TrackRef primaryTrackRef = ref->trackRef();
00785           for(IE iel = elements_.begin(); iel != elements_.end(); iel++){
00786             reco::TrackRef elemTrackRef = (*iel)->trackRef();
00787             if( primaryTrackRef == elemTrackRef ) {
00788               (*iel)->setMuonRef( muonref );
00789               if (debug_) std::cout << "One of the tracks identified in displaced vertices collections was spotted as muon" <<std:: endl;
00790             }
00791           }
00792         } else primaryElement->setMuonRef( muonref );
00793 
00794       }
00795 
00796       if (!trackMaskVertex.empty() && !trackMaskVertex[i]) continue;
00797 
00798       
00799       // set track type T_FROM_GAMMA for pfrectracks associated to conv brems
00800       if(useConvBremPFRecTracks_) {
00801         if(convBremPFRecTracks.size() > 0.) {
00802           for(unsigned int iconv = 0; iconv < convBremPFRecTracks.size(); iconv++) {
00803             if((*ref).trackRef() == (*convBremPFRecTracks[iconv]).trackRef()) {
00804               bool value = true;
00805               primaryElement->setTrackType(reco::PFBlockElement::T_FROM_GAMMACONV, value);
00806             }
00807           }
00808         }
00809       }
00810       elements_.push_back( primaryElement );
00811 
00812     }
00813 
00814     if (debug_) std::cout << " " << std::endl;
00815 
00816   }
00817 
00818  
00819   // -------------- GSF tracks and brems for Conversion Recovery ----------
00820    
00821   if(convbremgsftrackh.isValid() ) {
00822     
00823  
00824     const  reco::GsfPFRecTrackCollection ConvPFGsfProd = *(convbremgsftrackh.product());
00825     for(unsigned i=0;i<convbremgsftrackh->size(); i++) {
00826 
00827       reco::GsfPFRecTrackRef refgsf(convbremgsftrackh,i );   
00828       
00829       if((refgsf).isNull()) continue;
00830       
00831       reco::PFBlockElement* gsfEl;
00832       
00833       const  std::vector<reco::PFTrajectoryPoint> 
00834         PfGsfPoint =  ConvPFGsfProd[i].trajectoryPoints();
00835       
00836       unsigned int c_gsf=0;
00837       bool PassTracker = false;
00838       bool GetPout = false;
00839       unsigned int IndexPout = -1;
00840       
00841       typedef std::vector<reco::PFTrajectoryPoint>::const_iterator IP;
00842       for(IP itPfGsfPoint =  PfGsfPoint.begin();  
00843           itPfGsfPoint!= PfGsfPoint.end();itPfGsfPoint++) {
00844         
00845         if (itPfGsfPoint->isValid()){
00846           int layGsfP = itPfGsfPoint->layer();
00847           if (layGsfP == -1) PassTracker = true;
00848           if (PassTracker && layGsfP > 0 && GetPout == false) {
00849             IndexPout = c_gsf-1;
00850             GetPout = true;
00851           }
00852           //const math::XYZTLorentzVector GsfMoment = itPfGsfPoint->momentum();
00853           c_gsf++;
00854         }
00855       }
00856       math::XYZTLorentzVector pin = PfGsfPoint[0].momentum();      
00857       math::XYZTLorentzVector pout = PfGsfPoint[IndexPout].momentum();
00858       
00859       
00860     
00861       gsfEl = new reco::PFBlockElementGsfTrack(refgsf, pin, pout);
00862       
00863       bool valuegsf = true;
00864       // IMPORTANT SET T_FROM_GAMMACONV trackType() FOR CONVERSIONS
00865       gsfEl->setTrackType(reco::PFBlockElement::T_FROM_GAMMACONV, valuegsf);
00866 
00867       
00868 
00869       elements_.push_back( gsfEl);
00870       std::vector<reco::PFBrem> pfbrem = refgsf->PFRecBrem();
00871       
00872       for (unsigned i2=0;i2<pfbrem.size(); i2++) {
00873         const double DP = pfbrem[i2].DeltaP();
00874         const double SigmaDP =  pfbrem[i2].SigmaDeltaP(); 
00875         const unsigned int TrajP = pfbrem[i2].indTrajPoint();
00876         if(TrajP == 99) continue;
00877 
00878         reco::PFBlockElement* bremEl;
00879         bremEl = new reco::PFBlockElementBrem(refgsf,DP,SigmaDP,TrajP);
00880         elements_.push_back(bremEl);
00881         
00882       }
00883     }
00884   }
00885 
00886   
00887   // -------------- ECAL clusters ---------------------
00888 
00889 
00890   if(ecalh.isValid() ) {
00891     pfcSCVec_.resize(ecalh->size(),-1);
00892     for(unsigned i=0;i<ecalh->size(); i++)  {
00893 
00894       // this ecal cluster has been disabled
00895       if( !ecalMask.empty() &&
00896           !ecalMask[i] ) continue;
00897 
00898       reco::PFClusterRef ref( ecalh,i );
00899       reco::PFBlockElement* te
00900         = new reco::PFBlockElementCluster( ref,
00901                                            reco::PFBlockElement::ECAL);
00902       elements_.push_back( te );
00903       // Now mapping with Superclusters
00904       int scindex= ClusterClusterMapping::checkOverlap(*ref,superClusters_);
00905 
00906       if(scindex>=0)    {
00907           pfcSCVec_[ref.key()]=scindex;
00908           scpfcRefs_[scindex].push_back(ref);
00909         }
00910     }
00911   }
00912 
00913   // -------------- HCAL clusters ---------------------
00914 
00915   if(hcalh.isValid() ) {
00916     
00917     for(unsigned i=0;i<hcalh->size(); i++)  {
00918       
00919       // this hcal cluster has been disabled
00920       if( !hcalMask.empty() &&
00921           !hcalMask[i] ) continue;
00922       
00923       reco::PFClusterRef ref( hcalh,i );
00924       reco::PFBlockElement* th
00925         = new reco::PFBlockElementCluster( ref,
00926                                            reco::PFBlockElement::HCAL );
00927       elements_.push_back( th );
00928     }
00929   }
00930 
00931 
00932   // -------------- HFEM clusters ---------------------
00933 
00934   if(hfemh.isValid() ) {
00935     
00936     for(unsigned i=0;i<hfemh->size(); i++)  {
00937       
00938       // this hfem cluster has been disabled
00939       if( !hfemMask.empty() &&
00940           !hfemMask[i] ) continue;
00941       
00942       reco::PFClusterRef ref( hfemh,i );
00943       reco::PFBlockElement* th
00944         = new reco::PFBlockElementCluster( ref,
00945                                            reco::PFBlockElement::HFEM );
00946       elements_.push_back( th );
00947     }
00948   }
00949 
00950 
00951   // -------------- HFHAD clusters ---------------------
00952 
00953   if(hfhadh.isValid() ) {
00954     
00955     for(unsigned i=0;i<hfhadh->size(); i++)  {
00956       
00957       // this hfhad cluster has been disabled
00958       if( !hfhadMask.empty() &&
00959           !hfhadMask[i] ) continue;
00960       
00961       reco::PFClusterRef ref( hfhadh,i );
00962       reco::PFBlockElement* th
00963         = new reco::PFBlockElementCluster( ref,
00964                                            reco::PFBlockElement::HFHAD );
00965       elements_.push_back( th );
00966     }
00967   }
00968 
00969 
00970 
00971 
00972   // -------------- PS clusters ---------------------
00973 
00974   if(psh.isValid() ) {
00975     for(unsigned i=0;i<psh->size(); i++)  {
00976 
00977       // this ps cluster has been disabled
00978       if( !psMask.empty() &&
00979           !psMask[i] ) continue;
00980       reco::PFBlockElement::Type type = reco::PFBlockElement::NONE;
00981       reco::PFClusterRef ref( psh,i );
00982       // two types of elements:  PS1 (V) and PS2 (H) 
00983       // depending on layer:  PS1 or PS2
00984       switch(ref->layer()){
00985       case PFLayer::PS1:
00986         type = reco::PFBlockElement::PS1;
00987         break;
00988       case PFLayer::PS2:
00989         type = reco::PFBlockElement::PS2;
00990         break;
00991       default:
00992         break;
00993       }
00994       reco::PFBlockElement* tp
00995         = new reco::PFBlockElementCluster( ref,
00996                                            type );
00997       elements_.push_back( tp );
00998       
00999     }
01000   }
01001 }
01002 
01003 
01004 
01005 #endif
01006 
01007