CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MkFitSeedConverter.cc
Go to the documentation of this file.
2 
6 
12 
16 
20 
23 
26 
28 
29 // ROOT
30 #include "Math/SVector.h"
31 #include "Math/SMatrix.h"
32 
33 // mkFit includes
36 
38 public:
40  ~MkFitSeedConverter() override = default;
41 
42  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
43 
44 private:
45  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
46 
48  const TrackerTopology& ttopo,
49  const TransientTrackingRecHitBuilder& ttrhBuilder,
50  const MagneticField& mf,
51  const MkFitGeometry& mkFitGeom) const;
52 
53  using SVector3 = ROOT::Math::SVector<float, 3>;
54  using SMatrixSym33 = ROOT::Math::SMatrix<float, 3, 3, ROOT::Math::MatRepSym<float, 3>>;
55  using SMatrixSym66 = ROOT::Math::SMatrix<float, 6, 6, ROOT::Math::MatRepSym<float, 6>>;
56 
63 };
64 
66  : seedToken_{consumes<edm::View<TrajectorySeed>>(iConfig.getParameter<edm::InputTag>("seeds"))},
67  ttrhBuilderToken_{esConsumes<TransientTrackingRecHitBuilder, TransientRecHitRecord>(
68  iConfig.getParameter<edm::ESInputTag>("ttrhBuilder"))},
69  ttopoToken_{esConsumes<TrackerTopology, TrackerTopologyRcd>()},
70  mkFitGeomToken_{esConsumes<MkFitGeometry, TrackerRecoGeometryRecord>()},
71  mfToken_{esConsumes<MagneticField, IdealMagneticFieldRecord>()},
72  putToken_{produces<MkFitSeedWrapper>()} {}
73 
76 
77  desc.add("seeds", edm::InputTag{"initialStepSeeds"});
78  desc.add("ttrhBuilder", edm::ESInputTag{"", "WithTrackAngle"});
79 
80  descriptions.addWithDefaultLabel(desc);
81 }
82 
84  iEvent.emplace(putToken_,
85  convertSeeds(iEvent.get(seedToken_),
86  iSetup.getData(ttopoToken_),
87  iSetup.getData(ttrhBuilderToken_),
88  iSetup.getData(mfToken_),
89  iSetup.getData(mkFitGeomToken_)));
90 }
91 
93  const TrackerTopology& ttopo,
94  const TransientTrackingRecHitBuilder& ttrhBuilder,
95  const MagneticField& mf,
96  const MkFitGeometry& mkFitGeom) const {
98  ret.reserve(seeds.size());
99 
100  auto isPlusSide = [&ttopo](const DetId& detid) {
101  return ttopo.side(detid) == static_cast<unsigned>(TrackerDetSide::PosEndcap);
102  };
103 
104  int seed_index = 0;
105  for (const auto& seed : seeds) {
106  auto const& hitRange = seed.recHits();
107  const auto lastRecHit = ttrhBuilder.build(&*(hitRange.end() - 1));
108  const auto tsos = trajectoryStateTransform::transientState(seed.startingState(), lastRecHit->surface(), &mf);
109  const auto& stateGlobal = tsos.globalParameters();
110  const auto& gpos = stateGlobal.position();
111  const auto& gmom = stateGlobal.momentum();
112  SVector3 pos(gpos.x(), gpos.y(), gpos.z());
113  SVector3 mom(gmom.x(), gmom.y(), gmom.z());
114 
115  const auto& cov = tsos.curvilinearError().matrix();
116  SMatrixSym66 err; //fill a sub-matrix, mkfit::TrackState will convert internally
117  for (int i = 0; i < 5; ++i) {
118  for (int j = i; j < 5; ++j) {
119  err.At(i, j) = cov[i][j];
120  }
121  }
122 
123  mkfit::TrackState state(tsos.charge(), pos, mom, err);
124  state.convertFromGlbCurvilinearToCCS();
125  ret.emplace_back(state, 0, seed_index, 0, nullptr);
126  LogTrace("MkFitSeedConverter") << "Inserted seed with index " << seed_index;
127 
128  // Add hits
129  for (auto const& recHit : hitRange) {
130  if (not trackerHitRTTI::isFromDet(recHit)) {
131  throw cms::Exception("Assert") << "Encountered a seed with a hit which is not trackerHitRTTI::isFromDet()";
132  }
133  auto& baseTrkRecHit = static_cast<const BaseTrackerRecHit&>(recHit);
134  if (!baseTrkRecHit.isMatched()) {
135  const auto& clusterRef = baseTrkRecHit.firstClusterRef();
136  const auto detId = recHit.geographicalId();
137  const auto ilay = mkFitGeom.layerNumberConverter().convertLayerNumber(
138  detId.subdetId(), ttopo.layer(detId), false, ttopo.isStereo(detId), isPlusSide(detId));
139  LogTrace("MkFitSeedConverter") << " adding hit detid " << detId.rawId() << " index " << clusterRef.index()
140  << " ilay " << ilay;
141  ret.back().addHitIdx(clusterRef.index(), ilay, 0); // per-hit chi2 is not known
142  } else {
143  auto& matched2D = dynamic_cast<const SiStripMatchedRecHit2D&>(recHit);
144  const OmniClusterRef* const clRefs[2] = {&matched2D.monoClusterRef(), &matched2D.stereoClusterRef()};
145  const DetId detIds[2] = {matched2D.monoId(), matched2D.stereoId()};
146  for (int ii = 0; ii < 2; ++ii) {
147  const auto& detId = detIds[ii];
148  const auto ilay = mkFitGeom.layerNumberConverter().convertLayerNumber(
149  detId.subdetId(), ttopo.layer(detId), false, ttopo.isStereo(detId), isPlusSide(detId));
150  LogTrace("MkFitSeedConverter") << " adding matched hit detid " << detId.rawId() << " index "
151  << clRefs[ii]->index() << " ilay " << ilay;
152  ret.back().addHitIdx(clRefs[ii]->index(), ilay, 0); // per-hit chi2 is not known
153  }
154  }
155  }
156  ++seed_index;
157  }
158  return ret;
159 }
160 
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
tuple ret
prodAgent to be discontinued
mkfit::LayerNumberConverter const & layerNumberConverter() const
Definition: MkFitGeometry.h:32
bool isFromDet(TrackingRecHit const &hit)
const edm::ESGetToken< MkFitGeometry, TrackerRecoGeometryRecord > mkFitGeomToken_
const edm::EDGetTokenT< edm::View< TrajectorySeed > > seedToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
unsigned int side(const DetId &id) const
ROOT::Math::SMatrix< float, 3, 3, ROOT::Math::MatRepSym< float, 3 >> SMatrixSym33
bool isStereo(const DetId &id) const
size_type size() const
const edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > ttopoToken_
mkfit::TrackVec convertSeeds(const edm::View< TrajectorySeed > &seeds, const TrackerTopology &ttopo, const TransientTrackingRecHitBuilder &ttrhBuilder, const MagneticField &mf, const MkFitGeometry &mkFitGeom) const
int ii
Definition: cuy.py:589
#define LogTrace(id)
MkFitSeedConverter(edm::ParameterSet const &iConfig)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
int iEvent
Definition: GenABIO.cc:224
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
const edm::ESGetToken< TransientTrackingRecHitBuilder, TransientRecHitRecord > ttrhBuilderToken_
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const edm::EDPutTokenT< MkFitSeedWrapper > putToken_
Definition: DetId.h:17
std::vector< Track > TrackVec
TrajectoryStateOnSurface transientState(const PTrajectoryStateOnDet &ts, const Surface *surface, const MagneticField *field)
OrphanHandle< PROD > emplace(EDPutTokenT< PROD > token, Args &&...args)
puts a new product
Definition: Event.h:433
const GlobalTrajectoryParameters & globalParameters() const
int convertLayerNumber(int det, int lay, bool useMatched, int isStereo, bool posZ) const
virtual OmniClusterRef const & firstClusterRef() const =0
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
unsigned int layer(const DetId &id) const
unsigned int index() const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
ROOT::Math::SVector< float, 3 > SVector3
ROOT::Math::SMatrix< float, 6, 6, ROOT::Math::MatRepSym< float, 6 >> SMatrixSym66
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > mfToken_
~MkFitSeedConverter() override=default