Skip to content

Commit

Permalink
Merge pull request #22401 from swagata87/scouting_addTrkColl
Browse files Browse the repository at this point in the history
Scouting - add track collection
  • Loading branch information
cmsbuild committed Mar 6, 2018
2 parents 880944b + 0783d7f commit b25c0ef
Show file tree
Hide file tree
Showing 4 changed files with 174 additions and 0 deletions.
63 changes: 63 additions & 0 deletions DataFormats/Scouting/interface/ScoutingTrack.h
@@ -0,0 +1,63 @@
#ifndef DataFormats_ScoutingTrack_h
#define DataFormats_ScoutingTrack_h

#include <vector>

//class for holding track information, for use in data scouting
class ScoutingTrack
{
public:
//constructor with values for all data fields
ScoutingTrack(float tk_pt, float tk_eta, float tk_phi, float tk_chi2, float tk_ndof, int tk_charge, float tk_dxy, float tk_dz, int tk_nValidPixelHits, int tk_nTrackerLayersWithMeasurement, int tk_nValidStripHits, float tk_qoverp, float tk_lambda, float tk_dxy_Error, float tk_dz_Error, float tk_qoverp_Error, float tk_lambda_Error, float tk_phi_Error, float tk_dsz, float tk_dsz_Error):
tk_pt_(tk_pt), tk_eta_(tk_eta), tk_phi_(tk_phi), tk_chi2_(tk_chi2), tk_ndof_(tk_ndof), tk_charge_(tk_charge), tk_dxy_(tk_dxy), tk_dz_(tk_dz), tk_nValidPixelHits_(tk_nValidPixelHits), tk_nTrackerLayersWithMeasurement_(tk_nTrackerLayersWithMeasurement), tk_nValidStripHits_(tk_nValidStripHits), tk_qoverp_(tk_qoverp), tk_lambda_(tk_lambda), tk_dxy_Error_(tk_dxy_Error), tk_dz_Error_(tk_dz_Error), tk_qoverp_Error_(tk_qoverp_Error), tk_lambda_Error_(tk_lambda_Error), tk_phi_Error_(tk_phi_Error), tk_dsz_(tk_dsz), tk_dsz_Error_(tk_dsz_Error) {}
//default constructor
ScoutingTrack():tk_pt_(0), tk_eta_(0), tk_phi_(0), tk_chi2_(0), tk_ndof_(0), tk_charge_(0), tk_dxy_(0), tk_dz_(0), tk_nValidPixelHits_(0), tk_nTrackerLayersWithMeasurement_(0), tk_nValidStripHits_(0), tk_qoverp_(0), tk_lambda_(0), tk_dxy_Error_(0), tk_dz_Error_(0), tk_qoverp_Error_(0), tk_lambda_Error_(0), tk_phi_Error_(0), tk_dsz_(0), tk_dsz_Error_(0) {}

//accessor functions
float tk_pt() const { return tk_pt_; }
float tk_eta() const { return tk_eta_; }
float tk_phi() const { return tk_phi_; }
float tk_chi2() const { return tk_chi2_; }
float tk_ndof() const { return tk_ndof_; }
int tk_charge() const { return tk_charge_; }
float tk_dxy() const { return tk_dxy_; }
float tk_dz() const { return tk_dz_; }
int tk_nValidPixelHits() const { return tk_nValidPixelHits_; }
int tk_nTrackerLayersWithMeasurement() const { return tk_nTrackerLayersWithMeasurement_; }
int tk_nValidStripHits() const { return tk_nValidStripHits_; }
float tk_qoverp() const { return tk_qoverp_; }
float tk_lambda() const { return tk_lambda_; }
float tk_dxy_Error() const { return tk_dxy_Error_; }
float tk_dz_Error() const { return tk_dz_Error_; }
float tk_qoverp_Error() const { return tk_qoverp_Error_; }
float tk_lambda_Error() const { return tk_lambda_Error_; }
float tk_phi_Error() const { return tk_phi_Error_; }
float tk_dsz() const { return tk_dsz_; }
float tk_dsz_Error() const { return tk_dsz_Error_; }

private:
float tk_pt_;
float tk_eta_;
float tk_phi_;
float tk_chi2_;
float tk_ndof_;
int tk_charge_;
float tk_dxy_;
float tk_dz_;
int tk_nValidPixelHits_;
int tk_nTrackerLayersWithMeasurement_;
int tk_nValidStripHits_;
float tk_qoverp_;
float tk_lambda_;
float tk_dxy_Error_;
float tk_dz_Error_;
float tk_qoverp_Error_;
float tk_lambda_Error_;
float tk_phi_Error_;
float tk_dsz_;
float tk_dsz_Error_;
};

typedef std::vector<ScoutingTrack> ScoutingTrackCollection;

#endif
2 changes: 2 additions & 0 deletions DataFormats/Scouting/src/classes.h
@@ -1,6 +1,7 @@
#include "DataFormats/Scouting/interface/ScoutingCaloJet.h"
#include "DataFormats/Scouting/interface/ScoutingPFJet.h"
#include "DataFormats/Scouting/interface/ScoutingParticle.h"
#include "DataFormats/Scouting/interface/ScoutingTrack.h"
#include "DataFormats/Scouting/interface/ScoutingVertex.h"
#include "DataFormats/Scouting/interface/ScoutingElectron.h"
#include "DataFormats/Scouting/interface/ScoutingMuon.h"
Expand All @@ -18,5 +19,6 @@ namespace DataFormats_Scouting {
edm::Wrapper<ScoutingElectronCollection> sc5;
edm::Wrapper<ScoutingMuonCollection> sc6;
edm::Wrapper<ScoutingPhotonCollection> sc7;
edm::Wrapper<ScoutingTrackCollection> sc8;
};
}
5 changes: 5 additions & 0 deletions DataFormats/Scouting/src/classes_def.xml
Expand Up @@ -24,18 +24,23 @@
<class name="ScoutingPhoton" ClassVersion="2">
<version ClassVersion="2" checksum="2797410539"/>
</class>
<class name="ScoutingTrack" ClassVersion="2">
<version ClassVersion="2" checksum="1342928770"/>
</class>
<class name="std::vector<ScoutingCaloJet>"/>
<class name="std::vector<ScoutingPFJet>"/>
<class name="std::vector<ScoutingParticle>"/>
<class name="std::vector<ScoutingVertex>"/>
<class name="std::vector<ScoutingElectron>"/>
<class name="std::vector<ScoutingMuon>"/>
<class name="std::vector<ScoutingPhoton>"/>
<class name="std::vector<ScoutingTrack>"/>
<class name="edm::Wrapper<std::vector<ScoutingCaloJet> >"/>
<class name="edm::Wrapper<std::vector<ScoutingPFJet> >"/>
<class name="edm::Wrapper<std::vector<ScoutingParticle> >"/>
<class name="edm::Wrapper<std::vector<ScoutingVertex> >"/>
<class name="edm::Wrapper<std::vector<ScoutingElectron> >"/>
<class name="edm::Wrapper<std::vector<ScoutingMuon> >"/>
<class name="edm::Wrapper<std::vector<ScoutingPhoton> >"/>
<class name="edm::Wrapper<std::vector<ScoutingTrack> >"/>
</lcgdict>
104 changes: 104 additions & 0 deletions HLTrigger/Muon/plugins/HLTScoutingTrackProducer.cc
@@ -0,0 +1,104 @@
// -*- C++ -*-
//
// Package: HLTrigger/Muon
// Class: HLTScoutingTrackProducer
//
/**\class HLTScoutingTrackProducer HLTScoutingTrackProducer.cc HLTScoutingTrackProducer.cc
Description: Producer for Scouting Tracks
*/
//

#include <memory>
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/global/EDProducer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "DataFormats/Common/interface/AssociationMap.h"
#include "DataFormats/Common/interface/getRef.h"
#include "DataFormats/Common/interface/ValueMap.h"
#include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
#include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
#include "DataFormats/TrackReco/interface/HitPattern.h"
#include "DataFormats/TrackReco/interface/Track.h"
#include "DataFormats/TrackReco/interface/TrackFwd.h"
#include "DataFormats/VertexReco/interface/Vertex.h"
#include "DataFormats/VertexReco/interface/VertexFwd.h"
#include "DataFormats/Scouting/interface/ScoutingMuon.h"
#include "DataFormats/Scouting/interface/ScoutingTrack.h"
#include "DataFormats/Scouting/interface/ScoutingVertex.h"
#include "DataFormats/MuonReco/interface/MuonTrackLinks.h"
#include "DataFormats/MuonReco/interface/MuonFwd.h"

class HLTScoutingTrackProducer : public edm::global::EDProducer<> {

public:
explicit HLTScoutingTrackProducer(const edm::ParameterSet&);
~HLTScoutingTrackProducer();

static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

private:
void produce(edm::StreamID sid, edm::Event & iEvent, edm::EventSetup const & setup)
const final;

const edm::EDGetTokenT<reco::TrackCollection> otherTrackCollection_;


};

//
// constructors and destructor
//
HLTScoutingTrackProducer::HLTScoutingTrackProducer(const edm::ParameterSet& iConfig):
otherTrackCollection_(consumes<reco::TrackCollection>
(iConfig.getParameter<edm::InputTag>("OtherTracks")))
{
//register products
produces<ScoutingTrackCollection>();
}

HLTScoutingTrackProducer::~HLTScoutingTrackProducer() = default;

// ------------ method called to produce the data ------------
void HLTScoutingTrackProducer::produce(edm::StreamID sid, edm::Event & iEvent,
edm::EventSetup const & setup) const
{
using namespace edm;

std::unique_ptr<ScoutingTrackCollection> outTrack(new ScoutingTrackCollection());

Handle<reco::TrackCollection> otherTrackCollection;

if(iEvent.getByToken(otherTrackCollection_, otherTrackCollection)){
// Produce tracks in event
for (auto &trk : *otherTrackCollection) {
outTrack->emplace_back(trk.pt(), trk.eta(), trk.phi(),trk.chi2(), trk.ndof(),
trk.charge(), trk.dxy(), trk.dz(), trk.hitPattern().numberOfValidPixelHits(),
trk.hitPattern().trackerLayersWithMeasurement(), trk.hitPattern().numberOfValidStripHits(),
trk.qoverp(), trk.lambda(), trk.dxyError(), trk.dzError(),
trk.qoverpError(),
trk.lambdaError(),
trk.phiError(),
trk.dsz(),
trk.dszError()
);

}
}

iEvent.put(std::move(outTrack));
}

// ------------ method fills 'descriptions' with the allowed parameters for the module ------------
void HLTScoutingTrackProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
edm::ParameterSetDescription desc;
desc.add<edm::InputTag>("OtherTracks", edm::InputTag("hltPixelTracksL3MuonNoVtx"));
descriptions.add("hltScoutingTrackProducer", desc);
}

// declare this class as a framework plugin
#include "FWCore/Framework/interface/MakerMacros.h"
DEFINE_FWK_MODULE(HLTScoutingTrackProducer);

0 comments on commit b25c0ef

Please sign in to comment.