CMS 3D CMS Logo

Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes

PF_PU_FirstVertexTracks Class Reference

#include <PF_PU_FirstVertexTracks.h>

Inheritance diagram for PF_PU_FirstVertexTracks:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 PF_PU_FirstVertexTracks (const edm::ParameterSet &)
 ~PF_PU_FirstVertexTracks ()

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)

Private Member Functions

virtual void produce (edm::Event &, const edm::EventSetup &)
virtual bool TrackMatch (reco::Track, reco::Track)

Private Attributes

edm::InputTag input_AssociationMap_
edm::InputTag input_AssociationType_
edm::InputTag input_generalTracksCollection_
int input_MinQuality_
edm::InputTag input_VertexCollection_

Detailed Description

Definition at line 40 of file PF_PU_FirstVertexTracks.h.


Constructor & Destructor Documentation

PF_PU_FirstVertexTracks::PF_PU_FirstVertexTracks ( const edm::ParameterSet iConfig) [explicit]

Definition at line 61 of file PF_PU_FirstVertexTracks.cc.

References gather_cfg::cout, and edm::ParameterSet::getParameter().

{
   //now do what ever other initialization is needed

        input_AssociationType_ = iConfig.getParameter<edm::InputTag>("AssociationType");

        input_AssociationMap_ = iConfig.getParameter<InputTag>("AssociationMap");

        input_generalTracksCollection_ = iConfig.getParameter<InputTag>("TrackCollection");

        input_VertexCollection_ = iConfig.getParameter<InputTag>("VertexCollection");

        input_MinQuality_ = iConfig.getParameter<int>("MinQuality");

   //register your products

        if ( input_AssociationType_.label() == "TracksToVertex" ) {
          produces<TrackCollection>("T2V");
        } else {
          if ( input_AssociationType_.label() == "VertexToTracks" ) {
            produces<TrackCollection>("V2T");
          } else {
            if ( input_AssociationType_.label() == "Both" ) {
              produces<TrackCollection>("T2V");
              produces<TrackCollection>("V2T");
            } else {
              std::cout << "No correct InputTag for AssociationType!" << std::endl;
              std::cout << "Won't produce any TrackCollection!" << std::endl;
            }
          }
        }

}
PF_PU_FirstVertexTracks::~PF_PU_FirstVertexTracks ( )

Definition at line 96 of file PF_PU_FirstVertexTracks.cc.

{
 
   // do anything here that needs to be done at desctruction time
   // (e.g. close files, deallocate resources etc.)

}

Member Function Documentation

void PF_PU_FirstVertexTracks::fillDescriptions ( edm::ConfigurationDescriptions descriptions) [static]

Reimplemented from edm::EDProducer.

Definition at line 240 of file PF_PU_FirstVertexTracks.cc.

References edm::ConfigurationDescriptions::addDefault(), and edm::ParameterSetDescription::setUnknown().

                                                                                    {
  //The following says we do not know what parameters are allowed so do no validation
  // Please change this to state exactly what you do use, even if it is no parameters
  edm::ParameterSetDescription desc;
  desc.setUnknown();
  descriptions.addDefault(desc);
}
void PF_PU_FirstVertexTracks::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [private, virtual]

Implements edm::EDProducer.

Definition at line 111 of file PF_PU_FirstVertexTracks.cc.

References gather_cfg::cout, edm::first(), edm::Event::getByLabel(), edm::Event::put(), edm::second(), and findQualityFiles::size.

{

        auto_ptr<TrackCollection> t2v_firstvertextracks(new TrackCollection() );
        auto_ptr<TrackCollection> v2t_firstvertextracks(new TrackCollection() );

        bool t2vassmap = false;
        bool v2tassmap = false;
  
        //get the input vertex<->general track association map
        Handle<TrackToVertexAssMap> t2vAM;
        Handle<VertexToTrackAssMap> v2tAM;
        
        string asstype = input_AssociationType_.label();

        if ( ( asstype == "TracksToVertex" ) || ( asstype == "Both" ) ) {
          if ( iEvent.getByLabel(input_AssociationMap_, t2vAM ) ) {
            t2vassmap = true;
          }
        }

        if ( ( asstype == "VertexToTracks" ) || ( asstype == "Both" ) ) {
          if ( iEvent.getByLabel(input_AssociationMap_, v2tAM ) ) {
            v2tassmap = true;
          }
        }

        if ( !t2vassmap && !v2tassmap ) {
          cout << "No input collection could be found" << endl;
          return;
        }
 
        //get the input track collection
        Handle<TrackCollection> input_trckcollH;
        iEvent.getByLabel(input_generalTracksCollection_,input_trckcollH);

        if ( t2vassmap ){

          const TrackQualityPairVector trckcoll = t2vAM->begin()->val;

          //get the tracks associated to the first vertex and store them in a track collection
          for (unsigned int trckcoll_ite = 0; trckcoll_ite < trckcoll.size(); trckcoll_ite++){

            float quality = trckcoll[trckcoll_ite].second;

            if ( quality>=input_MinQuality_ ) { 

              TrackRef AMtrkref = trckcoll[trckcoll_ite].first;

              for(unsigned int index_input_trck=0; index_input_trck<input_trckcollH->size(); index_input_trck++){

                TrackRef input_trackref = TrackRef(input_trckcollH,index_input_trck);

                if( TrackMatch(*AMtrkref,*input_trackref) ){

                  t2v_firstvertextracks->push_back(*AMtrkref);
                  break;
                              
                } 

              }

            }

          }

          iEvent.put( t2v_firstvertextracks, "T2V" );

        } 

        if ( v2tassmap ) {
 
          //get the input vertex collection
          Handle<VertexCollection> input_vtxcollH;
          iEvent.getByLabel(input_VertexCollection_,input_vtxcollH);

          VertexRef firstVertexRef(input_vtxcollH,0);

          VertexToTrackAssMap::const_iterator v2t_ite;

          for(v2t_ite=v2tAM->begin(); v2t_ite!=v2tAM->end(); v2t_ite++){

            TrackRef AMtrkref = v2t_ite->key;

            for(unsigned int index_input_trck=0; index_input_trck<input_trckcollH->size(); index_input_trck++){

              TrackRef input_trackref = TrackRef(input_trckcollH,index_input_trck);

              if(TrackMatch(*AMtrkref,*input_trackref)){
    
                for(unsigned v_ite = 0; v_ite<(v2t_ite->val).size(); v_ite++){

                  VertexRef vtxref = (v2t_ite->val)[v_ite].first;
                  float quality = (v2t_ite->val)[v_ite].second;

                  if ( (vtxref==firstVertexRef) && (quality>=input_MinQuality_) ){
                    v2t_firstvertextracks->push_back(*AMtrkref);
                  }

                }

              }

            }

          }

          iEvent.put( v2t_firstvertextracks, "V2T" );

        }

}
bool PF_PU_FirstVertexTracks::TrackMatch ( reco::Track  track1,
reco::Track  track2 
) [private, virtual]

Definition at line 225 of file PF_PU_FirstVertexTracks.cc.

References eta(), AlCaHLTBitMon_ParallelJobs::p, and phi.

{

        return (
          (track1).eta()  == (track2).eta() &&
          (track1).phi()  == (track2).phi() &&
          (track1).chi2() == (track2).chi2() &&
          (track1).ndof() == (track2).ndof() &&
          (track1).p()    == (track2).p()
        );

}

Member Data Documentation

Definition at line 55 of file PF_PU_FirstVertexTracks.h.

Definition at line 53 of file PF_PU_FirstVertexTracks.h.

Definition at line 56 of file PF_PU_FirstVertexTracks.h.

Definition at line 59 of file PF_PU_FirstVertexTracks.h.

Definition at line 57 of file PF_PU_FirstVertexTracks.h.