Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:26:35

0001 // -*- C++ -*-
0002 //
0003 // Package:    TrackAssociator
0004 // Class:      TrivialExample
0005 // 
0006 /*
0007 
0008  Description: Trivial example to use get energy for a collection of ctfWithMaterialTracks
0009 
0010 */
0011 //
0012 // Original Author:  Dmytro Kovalskyi
0013 
0014 #include "FWCore/Framework/interface/EDAnalyzer.h"
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "FWCore/Framework/interface/EventSetup.h"
0017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0018 #include "FWCore/Framework/interface/MakerMacros.h"
0019 
0020 #include "DataFormats/TrackReco/interface/Track.h"
0021 #include "DataFormats/Common/interface/Handle.h"
0022 
0023 #include "TrackingTools/TrackAssociator/interface/TrackDetectorAssociator.h"
0024 #include "TrackingTools/TrackAssociator/interface/TrackAssociatorParameters.h"
0025 
0026 class TrivialExample : public edm::EDAnalyzer {
0027  public:
0028    explicit TrivialExample(const edm::ParameterSet&);
0029    virtual ~TrivialExample(){}
0030    virtual void analyze (const edm::Event&, const edm::EventSetup&);
0031 
0032  private:
0033    TrackDetectorAssociator trackAssociator_;
0034    TrackAssociatorParameters parameters_;
0035 };
0036 
0037 TrivialExample::TrivialExample(const edm::ParameterSet& iConfig)
0038 {
0039    // TrackAssociator parameters
0040    edm::ParameterSet parameters = iConfig.getParameter<edm::ParameterSet>("TrackAssociatorParameters");
0041    parameters_.loadParameters( parameters );
0042    trackAssociator_.useDefaultPropagator();
0043 }
0044 
0045 void TrivialExample::analyze( const edm::Event& iEvent, const edm::EventSetup& iSetup)
0046 {
0047    // get reco tracks 
0048    edm::Handle<reco::TrackCollection> recoTracks;
0049    iEvent.getByLabel("ctfWithMaterialTracks", recoTracks);
0050    if (! recoTracks.isValid() ) throw cms::Exception("FatalError") << "No reco tracks were found\n";
0051 
0052    for(reco::TrackCollection::const_iterator recoTrack = recoTracks->begin();
0053        recoTrack != recoTracks->end(); ++recoTrack){
0054        
0055       if (recoTrack->pt() < 2) continue; // skip low Pt tracks
0056 
0057       
0058       TrackDetMatchInfo info = trackAssociator_.associate(iEvent, iSetup, *recoTrack, parameters_);
0059       
0060       edm::LogVerbatim("TrackAssociator") << "\n-------------------------------------------------------\n Track (pt,eta,phi): " << 
0061     recoTrack->pt() << " , " << recoTrack->eta() << " , " << recoTrack->phi() ;
0062       edm::LogVerbatim("TrackAssociator") << "Ecal energy in crossed crystals based on RecHits: " << 
0063     info.crossedEnergy(TrackDetMatchInfo::EcalRecHits);
0064       edm::LogVerbatim("TrackAssociator") << "Ecal energy in 3x3 crystals based on RecHits: " << 
0065     info.nXnEnergy(TrackDetMatchInfo::EcalRecHits, 1);
0066       edm::LogVerbatim("TrackAssociator") << "Hcal energy in crossed towers based on RecHits: " << 
0067     info.crossedEnergy(TrackDetMatchInfo::HcalRecHits);
0068       edm::LogVerbatim("TrackAssociator") << "Hcal energy in 3x3 towers based on RecHits: " << 
0069     info.nXnEnergy(TrackDetMatchInfo::HcalRecHits, 1);
0070       edm::LogVerbatim("TrackAssociator") << "Number of muon segment matches: " << info.numberOfSegments();
0071 
0072    }
0073 }
0074 
0075 //define this as a plug-in
0076 DEFINE_FWK_MODULE(TrivialExample);