CMS 3D CMS Logo

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