Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

edm::Range interface to TrackCandidate recHits #30361

Merged
merged 2 commits into from
Jun 26, 2020
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion DQM/TrackingMonitor/src/TrackBuildingAnalyzer.cc
Original file line number Diff line number Diff line change
Expand Up @@ -601,7 +601,7 @@ void TrackBuildingAnalyzer::analyze(const edm::Event& iEvent,
//double qoverp = tsAtClosestApproachTrackCand.trackStateAtPCA().charge()/p.mag();
//double theta = p.theta();
//double lambda = M_PI/2-p.theta();
double numberOfHits = candidate.recHits().second - candidate.recHits().first;
double numberOfHits = candidate.nRecHits();
double dxy = (-v.x() * sin(p.phi()) + v.y() * cos(p.phi()));

double dz = v.z() - (v.x() * p.x() + v.y() * p.y()) / p.perp() * p.z() / p.perp();
Expand Down
6 changes: 3 additions & 3 deletions DataFormats/TrackCandidate/interface/TrackCandidate.h
Original file line number Diff line number Diff line change
Expand Up @@ -6,6 +6,7 @@
#include "DataFormats/TrackReco/interface/TrajectoryStopReasons.h"
#include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
#include "DataFormats/TrajectorySeed/interface/TrajectorySeed.h"
#include "FWCore/Utilities/interface/Range.h"

#include <utility>

Expand All @@ -23,8 +24,6 @@ only the second is compulsory,the other three can be empty / not present
class TrackCandidate {
public:
typedef edm::OwnVector<TrackingRecHit> RecHitContainer;
typedef RecHitContainer::const_iterator const_iterator;
typedef std::pair<const_iterator, const_iterator> range;

TrackCandidate()
: rh_(), seed_(), state_(), seedRef_(), nLoops_(0), stopReason_((uint8_t)StopReason::UNINITIALIZED) {}
Expand Down Expand Up @@ -55,7 +54,8 @@ class TrackCandidate {

PTrajectoryStateOnDet const& trajectoryStateOnDet() const { return state_; }

range recHits() const { return std::make_pair(rh_.begin(), rh_.end()); }
edm::Range<RecHitContainer::const_iterator> recHits() const { return {rh_.begin(), rh_.end()}; }
auto nRecHits() const { return rh_.size(); }

TrajectorySeed const& seed() const { return seed_; }

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -131,7 +131,6 @@ void TrackProducerWithSCAssociation::produce(edm::Event& theEvent, const edm::Ev
for (TrackCandidateCollection::const_iterator i = theTCCollection->begin(); i != theTCCollection->end(); i++) {
const TrackCandidate* theTC = &(*i);
PTrajectoryStateOnDet state = theTC->trajectoryStateOnDet();
const TrackCandidate::range& recHitVec = theTC->recHits();
const TrajectorySeed& seed = theTC->seed();

//convert PTrajectoryStateOnDet to TrajectoryStateOnSurface
Expand All @@ -148,8 +147,8 @@ void TrackProducerWithSCAssociation::produce(edm::Event& theEvent, const edm::Ev

float ndof = 0;

for (edm::OwnVector<TrackingRecHit>::const_iterator i = recHitVec.first; i != recHitVec.second; i++) {
hits.push_back(theBuilder.product()->build(&(*i)));
for (auto const& recHit : theTC->recHits()) {
hits.push_back(theBuilder.product()->build(&recHit));
}

//build Track
Expand Down
5 changes: 2 additions & 3 deletions RecoTracker/DebugTools/plugins/TestHits.cc
Original file line number Diff line number Diff line change
Expand Up @@ -204,7 +204,6 @@ void TestHits::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)

const TrackCandidate* theTC = &(*i);
PTrajectoryStateOnDet state = theTC->trajectoryStateOnDet();
const TrackCandidate::range& recHitVec = theTC->recHits();

//convert PTrajectoryStateOnDet to TrajectoryStateOnSurface

Expand All @@ -218,8 +217,8 @@ void TestHits::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
//convert the TrackingRecHit vector to a TransientTrackingRecHit vector
TransientTrackingRecHit::RecHitContainer hits;

for (edm::OwnVector<TrackingRecHit>::const_iterator i = recHitVec.first; i != recHitVec.second; i++) {
hits.push_back(theBuilder->build(&(*i)));
for (auto const& recHit : theTC->recHits()) {
hits.push_back(theBuilder->build(&recHit));
}

//call the fitter
Expand Down
5 changes: 2 additions & 3 deletions RecoTracker/DebugTools/plugins/TestSmoothHits.cc
Original file line number Diff line number Diff line change
Expand Up @@ -205,7 +205,6 @@ void TestSmoothHits::analyze(const edm::Event& iEvent, const edm::EventSetup& iS

const TrackCandidate* theTC = &(*i);
PTrajectoryStateOnDet state = theTC->trajectoryStateOnDet();
const TrackCandidate::range& recHitVec = theTC->recHits();

//convert PTrajectoryStateOnDet to TrajectoryStateOnSurface

Expand All @@ -219,8 +218,8 @@ void TestSmoothHits::analyze(const edm::Event& iEvent, const edm::EventSetup& iS
//convert the TrackingRecHit vector to a TransientTrackingRecHit vector
TransientTrackingRecHit::RecHitContainer hits;

for (edm::OwnVector<TrackingRecHit>::const_iterator i = recHitVec.first; i != recHitVec.second; i++) {
hits.push_back(theBuilder->build(&(*i)));
for (auto const& recHit : theTC->recHits()) {
hits.push_back(theBuilder->build(&recHit));
}

//call the fitter
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -411,10 +411,9 @@ namespace reco {

#ifdef EDM_ML_DEBUG
LogDebug("CosmicTrackSplitter") << " dumping the hits now: ";
for (TrackCandidate::range hitR = cand.recHits(); hitR.first != hitR.second; ++hitR.first) {
auto const &tmp = *hitR.first;
LogTrace("CosmicTrackSplitter") << " hit detid = " << hitR.first->geographicalId().rawId()
<< ", type = " << typeid(tmp).name();
for (auto const &hit : cand.recHits()) {
LogTrace("CosmicTrackSplitter") << " hit detid = " << hit.geographicalId().rawId()
<< ", type = " << typeid(hit).name();
}
#endif

Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -159,7 +159,7 @@ namespace {
continue; // at least "loose" ( FIXME: take cut value from CutSelector)

// if( ChiSquaredProbability(matchedTrack.chi2(),matchedTrack.ndof()) < minTrkProbCut_)continue;
int dHits = (cand.recHits().second - cand.recHits().first) - matchedTrack.recHitsSize();
int dHits = cand.nRecHits() - matchedTrack.recHitsSize();
if (dHits > diffHitsCut_)
continue;
matches.push_back(std::array<int, 3>{{i, candidateComponents[cInd].first, candidateComponents[cInd].second}});
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -76,18 +76,17 @@ void TrackCandidateTopBottomHitFilter::produce(edm::Event& iEvent, const edm::Ev
auto pOut = std::make_unique<TrackCandidateCollection>();
for (TrackCandidateCollection::const_iterator it = pIn->begin(); it != pIn->end(); ++it) {
PTrajectoryStateOnDet state = it->trajectoryStateOnDet();
TrackCandidate::range oldhits = it->recHits();
TrajectorySeed seed = it->seed();
TrackCandidate::RecHitContainer hits;
for (TrackCandidate::RecHitContainer::const_iterator hit = oldhits.first; hit != oldhits.second; ++hit) {
if (hit->isValid()) {
double hitY = theBuilder->build(&*hit)->globalPosition().y();
for (auto const& hit : it->recHits()) {
if (hit.isValid()) {
double hitY = theBuilder->build(&hit)->globalPosition().y();
if (seedY * hitY > 0)
hits.push_back(hit->clone());
hits.push_back(hit.clone());
else
break;
} else
hits.push_back(hit->clone());
hits.push_back(hit.clone());
}
if (hits.size() >= 3) {
TrackCandidate newTC(hits, seed, state);
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -45,7 +45,6 @@ void TrackProducerAlgorithm<T>::runWithCandidate(const TrackingGeometry* theG,
int ntc = 0;
for (auto const& theTC : theTCCollection) {
PTrajectoryStateOnDet const& state = theTC.trajectoryStateOnDet();
const TrackCandidate::range& recHitVec = theTC.recHits();
const TrajectorySeed& seed = theTC.seed();

//convert PTrajectoryStateOnDet to TrajectoryStateOnSurface
Expand All @@ -61,8 +60,8 @@ void TrackProducerAlgorithm<T>::runWithCandidate(const TrackingGeometry* theG,

float ndof = 0;
// we assume are transient...
for (auto i = recHitVec.first; i != recHitVec.second; ++i) {
hits.push_back((*i).cloneSH()); // major waste, will be soon fitted and recloned...
for (auto const& recHit : theTC.recHits()) {
hits.push_back(recHit.cloneSH()); // major waste, will be soon fitted and recloned...
}

stopReason_ = theTC.stopReason();
Expand Down
12 changes: 1 addition & 11 deletions RecoTracker/TrackProducer/plugins/FakeTrackProducers.cc
Original file line number Diff line number Diff line change
Expand Up @@ -14,7 +14,6 @@
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "FWCore/Utilities/interface/Range.h"

#include "DataFormats/TrackReco/interface/Track.h"
#include "DataFormats/TrackReco/interface/TrackExtra.h"
Expand All @@ -34,15 +33,6 @@
#include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
#include "Geometry/Records/interface/TrackerTopologyRcd.h"
#include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
#include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"

namespace {
TrajectorySeed::RecHitRange getHits(const TrajectorySeed &seed) { return seed.recHits(); }
TrajectorySeed::RecHitRange getHits(const TrackCandidate &seed) {
auto range = seed.recHits();
return edm::Range{range.first, range.second};
}
} // namespace

template <class T>
class FakeTrackProducer : public edm::stream::EDProducer<> {
Expand Down Expand Up @@ -113,7 +103,7 @@ void FakeTrackProducer<T>::produce(edm::Event &iEvent, const edm::EventSetup &iS
reco::Track::Vector p(gp.x(), gp.y(), gp.z());
int charge = state.localParameters().charge();
out->push_back(reco::Track(1.0, 1.0, x, p, charge, reco::Track::CovarianceMatrix()));
auto hits = getHits(mu);
auto hits = mu.recHits();
out->back().appendHits(hits.begin(), hits.end(), ttopo);
// Now Track Extra
const TrackingRecHit *hit0 = &*hits.begin();
Expand Down