CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/CommonTools/RecoUtils/plugins/PF_PU_FirstVertexTracks.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    PF_PU_AssoMap
00004 // Class:      PF_PU_FirstVertexTracks
00005 // 
00010 //
00011 // Original Author:  Matthias Geisler
00012 //         Created:  Wed Apr 18 14:48:37 CEST 2012
00013 // $Id: PF_PU_FirstVertexTracks.cc,v 1.2 2013/05/23 15:41:36 gartung Exp $
00014 //
00015 //
00016 #include "CommonTools/RecoUtils/interface/PF_PU_FirstVertexTracks.h"
00017 
00018 // system include files
00019 #include <vector>
00020 #include <string>
00021 
00022 // user include files
00023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00024 
00025 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00026 #include "DataFormats/Common/interface/AssociationMap.h"
00027 #include "DataFormats/Common/interface/Handle.h"
00028 #include "DataFormats/Common/interface/OneToManyWithQuality.h"
00029 #include "DataFormats/Common/interface/OneToManyWithQualityGeneric.h"
00030 #include "DataFormats/Common/interface/View.h"
00031 
00032 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00033 #include "DataFormats/TrackReco/interface/TrackBase.h"
00034 #include "DataFormats/VertexReco/interface/Vertex.h"
00035 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00036 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00037 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00038 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00039 
00040 //
00041 // constants, enums and typedefs
00042 //
00043    
00044 using namespace edm;
00045 using namespace std;
00046 using namespace reco;
00047 
00048 typedef AssociationMap<OneToManyWithQuality< VertexCollection, TrackCollection, int> > TrackToVertexAssMap;
00049 typedef AssociationMap<OneToManyWithQuality< TrackCollection, VertexCollection, int> > VertexToTrackAssMap;
00050 
00051 typedef vector<pair<TrackRef, int> > TrackQualityPairVector;
00052 
00053 
00054 //
00055 // static data member definitions
00056 //
00057 
00058 //
00059 // constructors and destructor
00060 //
00061 PF_PU_FirstVertexTracks::PF_PU_FirstVertexTracks(const edm::ParameterSet& iConfig)
00062 {
00063    //now do what ever other initialization is needed
00064 
00065         input_AssociationType_ = iConfig.getParameter<edm::InputTag>("AssociationType");
00066 
00067         input_AssociationMap_ = iConfig.getParameter<InputTag>("AssociationMap");
00068 
00069         input_generalTracksCollection_ = iConfig.getParameter<InputTag>("TrackCollection");
00070 
00071         input_VertexCollection_ = iConfig.getParameter<InputTag>("VertexCollection");
00072 
00073         input_MinQuality_ = iConfig.getParameter<int>("MinQuality");
00074 
00075    //register your products
00076 
00077         if ( input_AssociationType_.label() == "TracksToVertex" ) {
00078           produces<TrackCollection>("T2V");
00079         } else {
00080           if ( input_AssociationType_.label() == "VertexToTracks" ) {
00081             produces<TrackCollection>("V2T");
00082           } else {
00083             if ( input_AssociationType_.label() == "Both" ) {
00084               produces<TrackCollection>("T2V");
00085               produces<TrackCollection>("V2T");
00086             } else {
00087               std::cout << "No correct InputTag for AssociationType!" << std::endl;
00088               std::cout << "Won't produce any TrackCollection!" << std::endl;
00089             }
00090           }
00091         }
00092 
00093 }
00094 
00095 
00096 PF_PU_FirstVertexTracks::~PF_PU_FirstVertexTracks()
00097 {
00098  
00099    // do anything here that needs to be done at desctruction time
00100    // (e.g. close files, deallocate resources etc.)
00101 
00102 }
00103 
00104 
00105 //
00106 // member functions
00107 //
00108 
00109 // ------------ method called to produce the data  ------------
00110 void
00111 PF_PU_FirstVertexTracks::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00112 {
00113 
00114         auto_ptr<TrackCollection> t2v_firstvertextracks(new TrackCollection() );
00115         auto_ptr<TrackCollection> v2t_firstvertextracks(new TrackCollection() );
00116 
00117         bool t2vassmap = false;
00118         bool v2tassmap = false;
00119   
00120         //get the input vertex<->general track association map
00121         Handle<TrackToVertexAssMap> t2vAM;
00122         Handle<VertexToTrackAssMap> v2tAM;
00123         
00124         string asstype = input_AssociationType_.label();
00125 
00126         if ( ( asstype == "TracksToVertex" ) || ( asstype == "Both" ) ) {
00127           if ( iEvent.getByLabel(input_AssociationMap_, t2vAM ) ) {
00128             t2vassmap = true;
00129           }
00130         }
00131 
00132         if ( ( asstype == "VertexToTracks" ) || ( asstype == "Both" ) ) {
00133           if ( iEvent.getByLabel(input_AssociationMap_, v2tAM ) ) {
00134             v2tassmap = true;
00135           }
00136         }
00137 
00138         if ( !t2vassmap && !v2tassmap ) {
00139           cout << "No input collection could be found" << endl;
00140           return;
00141         }
00142  
00143         //get the input track collection
00144         Handle<TrackCollection> input_trckcollH;
00145         iEvent.getByLabel(input_generalTracksCollection_,input_trckcollH);
00146 
00147         if ( t2vassmap ){
00148 
00149           const TrackQualityPairVector trckcoll = t2vAM->begin()->val;
00150 
00151           //get the tracks associated to the first vertex and store them in a track collection
00152           for (unsigned int trckcoll_ite = 0; trckcoll_ite < trckcoll.size(); trckcoll_ite++){
00153 
00154             float quality = trckcoll[trckcoll_ite].second;
00155 
00156             if ( quality>=input_MinQuality_ ) { 
00157 
00158               TrackRef AMtrkref = trckcoll[trckcoll_ite].first;
00159 
00160               for(unsigned int index_input_trck=0; index_input_trck<input_trckcollH->size(); index_input_trck++){
00161 
00162                 TrackRef input_trackref = TrackRef(input_trckcollH,index_input_trck);
00163 
00164                 if( TrackMatch(*AMtrkref,*input_trackref) ){
00165 
00166                   t2v_firstvertextracks->push_back(*AMtrkref);
00167                   break;
00168                               
00169                 } 
00170 
00171               }
00172 
00173             }
00174 
00175           }
00176 
00177           iEvent.put( t2v_firstvertextracks, "T2V" );
00178 
00179         } 
00180 
00181         if ( v2tassmap ) {
00182  
00183           //get the input vertex collection
00184           Handle<VertexCollection> input_vtxcollH;
00185           iEvent.getByLabel(input_VertexCollection_,input_vtxcollH);
00186 
00187           VertexRef firstVertexRef(input_vtxcollH,0);
00188 
00189           VertexToTrackAssMap::const_iterator v2t_ite;
00190 
00191           for(v2t_ite=v2tAM->begin(); v2t_ite!=v2tAM->end(); v2t_ite++){
00192 
00193             TrackRef AMtrkref = v2t_ite->key;
00194 
00195             for(unsigned int index_input_trck=0; index_input_trck<input_trckcollH->size(); index_input_trck++){
00196 
00197               TrackRef input_trackref = TrackRef(input_trckcollH,index_input_trck);
00198 
00199               if(TrackMatch(*AMtrkref,*input_trackref)){
00200     
00201                 for(unsigned v_ite = 0; v_ite<(v2t_ite->val).size(); v_ite++){
00202 
00203                   VertexRef vtxref = (v2t_ite->val)[v_ite].first;
00204                   float quality = (v2t_ite->val)[v_ite].second;
00205 
00206                   if ( (vtxref==firstVertexRef) && (quality>=input_MinQuality_) ){
00207                     v2t_firstvertextracks->push_back(*AMtrkref);
00208                   }
00209 
00210                 }
00211 
00212               }
00213 
00214             }
00215 
00216           }
00217 
00218           iEvent.put( v2t_firstvertextracks, "V2T" );
00219 
00220         }
00221 
00222 }
00223 
00224 bool 
00225 PF_PU_FirstVertexTracks::TrackMatch(const Track& track1,const Track& track2)
00226 {
00227 
00228         return (
00229           (track1).eta()  == (track2).eta() &&
00230           (track1).phi()  == (track2).phi() &&
00231           (track1).chi2() == (track2).chi2() &&
00232           (track1).ndof() == (track2).ndof() &&
00233           (track1).p()    == (track2).p()
00234         );
00235 
00236 }
00237 
00238 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
00239 void
00240 PF_PU_FirstVertexTracks::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
00241   //The following says we do not know what parameters are allowed so do no validation
00242   // Please change this to state exactly what you do use, even if it is no parameters
00243   edm::ParameterSetDescription desc;
00244   desc.setUnknown();
00245   descriptions.addDefault(desc);
00246 }
00247 
00248 //define this as a plug-in
00249 DEFINE_FWK_MODULE(PF_PU_FirstVertexTracks);