CMS 3D CMS Logo

MkFitInputConverter.cc
Go to the documentation of this file.
2 
6 
12 
16 
20 
23 
25 
26 // ROOT
27 #include "Math/SVector.h"
28 #include "Math/SMatrix.h"
29 
30 // mkFit includes
31 #include "Hit.h"
32 #include "Track.h"
33 #include "LayerNumberConverter.h"
34 
36 public:
37  explicit MkFitInputConverter(edm::ParameterSet const& iConfig);
38  ~MkFitInputConverter() override = default;
39 
40  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
41 
42 private:
43  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
44 
45  template <typename HitCollection>
46  void convertHits(const HitCollection& hits,
47  std::vector<mkfit::HitVec>& mkFitHits,
48  MkFitHitIndexMap& hitIndexMap,
49  int& totalHits,
50  const TrackerTopology& ttopo,
51  const TransientTrackingRecHitBuilder& ttrhBuilder,
52  const mkfit::LayerNumberConverter& lnc) const;
53 
54  bool passCCC(const SiStripRecHit2D& hit, const DetId hitId) const;
55  bool passCCC(const SiPixelRecHit& hit, const DetId hitId) const;
56 
58  const MkFitHitIndexMap& hitIndexMap,
59  const TransientTrackingRecHitBuilder& ttrhBuilder,
60  const MagneticField& mf) const;
61 
62  using SVector3 = ROOT::Math::SVector<float, 3>;
63  using SMatrixSym33 = ROOT::Math::SMatrix<float, 3, 3, ROOT::Math::MatRepSym<float, 3>>;
64  using SMatrixSym66 = ROOT::Math::SMatrix<float, 6, 6, ROOT::Math::MatRepSym<float, 6>>;
65 
74  const float minGoodStripCharge_;
75 };
76 
78  : pixelRecHitToken_{consumes<SiPixelRecHitCollection>(iConfig.getParameter<edm::InputTag>("pixelRecHits"))},
79  stripRphiRecHitToken_{
80  consumes<SiStripRecHit2DCollection>(iConfig.getParameter<edm::InputTag>("stripRphiRecHits"))},
81  stripStereoRecHitToken_{
82  consumes<SiStripRecHit2DCollection>(iConfig.getParameter<edm::InputTag>("stripStereoRecHits"))},
83  seedToken_{consumes<edm::View<TrajectorySeed>>(iConfig.getParameter<edm::InputTag>("seeds"))},
84  ttrhBuilderToken_{esConsumes<TransientTrackingRecHitBuilder, TransientRecHitRecord>(
85  iConfig.getParameter<edm::ESInputTag>("ttrhBuilder"))},
86  ttopoToken_{esConsumes<TrackerTopology, TrackerTopologyRcd>()},
87  mfToken_{esConsumes<MagneticField, IdealMagneticFieldRecord>()},
88  putToken_{produces<MkFitInputWrapper>()},
89  minGoodStripCharge_{static_cast<float>(
90  iConfig.getParameter<edm::ParameterSet>("minGoodStripCharge").getParameter<double>("value"))} {}
91 
94 
95  desc.add("pixelRecHits", edm::InputTag{"siPixelRecHits"});
96  desc.add("stripRphiRecHits", edm::InputTag{"siStripMatchedRecHits", "rphiRecHit"});
97  desc.add("stripStereoRecHits", edm::InputTag{"siStripMatchedRecHits", "stereoRecHit"});
98  desc.add("seeds", edm::InputTag{"initialStepSeeds"});
99  desc.add("ttrhBuilder", edm::ESInputTag{"", "WithTrackAngle"});
100 
102  descCCC.add<double>("value");
103  desc.add("minGoodStripCharge", descCCC);
104 
105  descriptions.add("mkFitInputConverterDefault", desc);
106 }
107 
109  mkfit::LayerNumberConverter lnc{mkfit::TkLayout::phase1};
110 
111  // Then import hits
112  const auto& ttrhBuilder = iSetup.getData(ttrhBuilderToken_);
113  const auto& ttopo = iSetup.getData(ttopoToken_);
114 
115  std::vector<mkfit::HitVec> mkFitHits(lnc.nLayers());
116  MkFitHitIndexMap hitIndexMap;
117  int totalHits = 0; // I need to have a global hit index in order to have the hit remapping working?
118  // Process strips first for better memory allocation pattern
119  convertHits(iEvent.get(stripRphiRecHitToken_), mkFitHits, hitIndexMap, totalHits, ttopo, ttrhBuilder, lnc);
120  convertHits(iEvent.get(stripStereoRecHitToken_), mkFitHits, hitIndexMap, totalHits, ttopo, ttrhBuilder, lnc);
121  convertHits(iEvent.get(pixelRecHitToken_), mkFitHits, hitIndexMap, totalHits, ttopo, ttrhBuilder, lnc);
122 
123  // Then import seeds
124  auto mkFitSeeds = convertSeeds(iEvent.get(seedToken_), hitIndexMap, ttrhBuilder, iSetup.getData(mfToken_));
125 
126  iEvent.emplace(putToken_, std::move(hitIndexMap), std::move(mkFitHits), std::move(mkFitSeeds), std::move(lnc));
127 }
128 
129 bool MkFitInputConverter::passCCC(const SiStripRecHit2D& hit, const DetId hitId) const {
130  return (siStripClusterTools::chargePerCM(hitId, hit.firstClusterRef().stripCluster()) > minGoodStripCharge_);
131 }
132 
133 bool MkFitInputConverter::passCCC(const SiPixelRecHit& hit, const DetId hitId) const { return true; }
134 
135 template <typename HitCollection>
137  std::vector<mkfit::HitVec>& mkFitHits,
138  MkFitHitIndexMap& hitIndexMap,
139  int& totalHits,
140  const TrackerTopology& ttopo,
141  const TransientTrackingRecHitBuilder& ttrhBuilder,
142  const mkfit::LayerNumberConverter& lnc) const {
143  if (hits.empty())
144  return;
145  auto isPlusSide = [&ttopo](const DetId& detid) {
146  return ttopo.side(detid) == static_cast<unsigned>(TrackerDetSide::PosEndcap);
147  };
148 
149  {
150  const DetId detid{hits.ids().back()};
151  const auto ilay =
152  lnc.convertLayerNumber(detid.subdetId(), ttopo.layer(detid), false, ttopo.isStereo(detid), isPlusSide(detid));
153  // Do initial reserves to minimize further memory allocations
154  const auto& lastClusterRef = hits.data().back().firstClusterRef();
155  hitIndexMap.resizeByClusterIndex(lastClusterRef.id(), lastClusterRef.index());
156  hitIndexMap.increaseLayerSize(ilay, hits.detsetSize(hits.ids().size() - 1));
157  }
158 
159  for (const auto& detset : hits) {
160  const DetId detid = detset.detId();
161  const auto subdet = detid.subdetId();
162  const auto layer = ttopo.layer(detid);
163  const auto isStereo = ttopo.isStereo(detid);
164  const auto ilay = lnc.convertLayerNumber(subdet, layer, false, isStereo, isPlusSide(detid));
165  hitIndexMap.increaseLayerSize(ilay, detset.size()); // to minimize memory allocations
166 
167  for (const auto& hit : detset) {
168  if (!passCCC(hit, detid))
169  continue;
170 
171  const auto& gpos = hit.globalPosition();
172  SVector3 pos(gpos.x(), gpos.y(), gpos.z());
173  const auto& gerr = hit.globalPositionError();
175  err.At(0, 0) = gerr.cxx();
176  err.At(1, 1) = gerr.cyy();
177  err.At(2, 2) = gerr.czz();
178  err.At(0, 1) = gerr.cyx();
179  err.At(0, 2) = gerr.czx();
180  err.At(1, 2) = gerr.czy();
181 
182  LogTrace("MkFitInputConverter") << "Adding hit detid " << detid.rawId() << " subdet " << subdet << " layer "
183  << layer << " isStereo " << isStereo << " zplus " << isPlusSide(detid) << " ilay "
184  << ilay;
185 
186  hitIndexMap.insert(hit.firstClusterRef().id(),
187  hit.firstClusterRef().index(),
188  MkFitHitIndexMap::MkFitHit{static_cast<int>(mkFitHits[ilay].size()), ilay},
189  &hit);
190  mkFitHits[ilay].emplace_back(pos, err, totalHits);
191  ++totalHits;
192  }
193  }
194 }
195 
197  const MkFitHitIndexMap& hitIndexMap,
198  const TransientTrackingRecHitBuilder& ttrhBuilder,
199  const MagneticField& mf) const {
201  ret.reserve(seeds.size());
202  int index = 0;
203  for (const auto& seed : seeds) {
204  const auto hitRange = seed.recHits();
205  const auto lastRecHit = ttrhBuilder.build(&*(hitRange.second - 1));
206  const auto tsos = trajectoryStateTransform::transientState(seed.startingState(), lastRecHit->surface(), &mf);
207  const auto& stateGlobal = tsos.globalParameters();
208  const auto& gpos = stateGlobal.position();
209  const auto& gmom = stateGlobal.momentum();
210  SVector3 pos(gpos.x(), gpos.y(), gpos.z());
211  SVector3 mom(gmom.x(), gmom.y(), gmom.z());
212 
213  const auto cartError = tsos.cartesianError(); // returns a temporary, so can't chain with the following line
214  const auto& cov = cartError.matrix();
216  for (int i = 0; i < 6; ++i) {
217  for (int j = i; j < 6; ++j) {
218  err.At(i, j) = cov[i][j];
219  }
220  }
221 
222  mkfit::TrackState state(tsos.charge(), pos, mom, err);
223  state.convertFromCartesianToCCS();
224  ret.emplace_back(state, 0, index, 0, nullptr);
225 
226  // Add hits
227  for (auto iHit = hitRange.first; iHit != hitRange.second; ++iHit) {
228  if (not trackerHitRTTI::isFromDet(*iHit)) {
229  throw cms::Exception("Assert") << "Encountered a seed with a hit which is not trackerHitRTTI::isFromDet()";
230  }
231  const auto& clusterRef = static_cast<const BaseTrackerRecHit&>(*iHit).firstClusterRef();
232  const auto& mkFitHit = hitIndexMap.mkFitHit(clusterRef.id(), clusterRef.index());
233  ret.back().addHitIdx(mkFitHit.index(), mkFitHit.layer(), 0); // per-hit chi2 is not known
234  }
235  ++index;
236  }
237  return ret;
238 }
239 
runTheMatrix.ret
ret
prodAgent to be discontinued
Definition: runTheMatrix.py:355
MkFitInputConverter::minGoodStripCharge_
const float minGoodStripCharge_
Definition: MkFitInputConverter.cc:74
edm::StreamID
Definition: StreamID.h:30
MkFitInputConverter::stripStereoRecHitToken_
edm::EDGetTokenT< SiStripRecHit2DCollection > stripStereoRecHitToken_
Definition: MkFitInputConverter.cc:68
TrackerTopology::side
unsigned int side(const DetId &id) const
Definition: TrackerTopology.cc:28
ZDCDigiStudy_cfi.HitCollection
HitCollection
Definition: ZDCDigiStudy_cfi.py:7
mps_fire.i
i
Definition: mps_fire.py:355
edm::ESInputTag
Definition: ESInputTag.h:87
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
hit::id
unsigned int id
Definition: SiStripHitEffFromCalibTree.cc:92
hfClusterShapes_cfi.hits
hits
Definition: hfClusterShapes_cfi.py:5
GlobalTrajectoryParameters::position
GlobalPoint position() const
Definition: GlobalTrajectoryParameters.h:60
MkFitInputConverter::SVector3
ROOT::Math::SVector< float, 3 > SVector3
Definition: MkFitInputConverter.cc:62
MkFitInputConverter::pixelRecHitToken_
edm::EDGetTokenT< SiPixelRecHitCollection > pixelRecHitToken_
Definition: MkFitInputConverter.cc:66
MkFitInputConverter::ttrhBuilderToken_
edm::ESGetToken< TransientTrackingRecHitBuilder, TransientRecHitRecord > ttrhBuilderToken_
Definition: MkFitInputConverter.cc:70
edm::EDGetTokenT
Definition: EDGetToken.h:33
MkFitInputConverter::SMatrixSym33
ROOT::Math::SMatrix< float, 3, 3, ROOT::Math::MatRepSym< float, 3 > > SMatrixSym33
Definition: MkFitInputConverter.cc:63
edm::EDPutTokenT< MkFitInputWrapper >
TrackerTopology
Definition: TrackerTopology.h:16
TransientRecHitRecord.h
pos
Definition: PixelAliasList.h:18
siStripClusterTools::chargePerCM
float chargePerCM(DetId detid, Iter a, Iter b)
Definition: SiStripClusterTools.h:29
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
TrackerTopology::layer
unsigned int layer(const DetId &id) const
Definition: TrackerTopology.cc:47
SiStripRecHit2D
Definition: SiStripRecHit2D.h:7
MkFitHitIndexMap::resizeByClusterIndex
void resizeByClusterIndex(edm::ProductID id, size_t clusterIndex)
Definition: MkFitHitIndexMap.cc:22
TransientTrackingRecHitBuilder::build
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
MkFitInputConverter::ttopoToken_
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > ttopoToken_
Definition: MkFitInputConverter.cc:71
MkFitInputConverter::MkFitInputConverter
MkFitInputConverter(edm::ParameterSet const &iConfig)
Definition: MkFitInputConverter.cc:77
MkFitInputConverter
Definition: MkFitInputConverter.cc:35
SiPixelRecHit
Our base class.
Definition: SiPixelRecHit.h:23
MkFitHitIndexMap::increaseLayerSize
void increaseLayerSize(int layer, size_t additionalSize)
Definition: MkFitHitIndexMap.cc:26
MkFitInputConverter::passCCC
bool passCCC(const SiStripRecHit2D &hit, const DetId hitId) const
Definition: MkFitInputConverter.cc:129
DetId
Definition: DetId.h:17
MakerMacros.h
TrackerTopology.h
TrackerTopologyRcd.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
TrackerTopology::isStereo
bool isStereo(const DetId &id) const
Definition: TrackerTopology.cc:158
IdealMagneticFieldRecord.h
MkFitHitIndexMap::mkFitHit
const MkFitHit & mkFitHit(edm::ProductID id, size_t clusterIndex) const
Get mkFit hit index and layer.
Definition: MkFitHitIndexMap.cc:51
TrajectorySeed.h
SiStripClusterTools.h
DetId::subdetId
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector's numbering enum)
Definition: DetId.h:48
edm::global::EDProducer
Definition: EDProducer.h:32
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
InitialStep_cff.seeds
seeds
Definition: InitialStep_cff.py:232
edm::View
Definition: CaloClusterFwd.h:14
mkfit::TrackVec
std::vector< Track > TrackVec
Definition: MkFitInputWrapper.h:14
SiPixelRecHitCollection.h
MkFitInputConverter::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: MkFitInputConverter.cc:92
edm::ParameterSet
Definition: ParameterSet.h:36
MkFitHitIndexMap::insert
void insert(edm::ProductID id, size_t clusterIndex, MkFitHit hit, const TrackingRecHit *hitPtr)
Definition: MkFitHitIndexMap.cc:33
Event.h
MkFitInputConverter::mfToken_
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > mfToken_
Definition: MkFitInputConverter.cc:72
runTheMatrix.err
err
Definition: runTheMatrix.py:288
MkFitInputConverter::~MkFitInputConverter
~MkFitInputConverter() override=default
MkFitInputConverter::seedToken_
edm::EDGetTokenT< edm::View< TrajectorySeed > > seedToken_
Definition: MkFitInputConverter.cc:69
iEvent
int iEvent
Definition: GenABIO.cc:224
trajectoryStateTransform::transientState
TrajectoryStateOnSurface transientState(const PTrajectoryStateOnDet &ts, const Surface *surface, const MagneticField *field)
Definition: TrajectoryStateTransform.cc:35
trackerHitRTTI.h
MkFitInputWrapper.h
SiStripRecHit2DCollection.h
MkFitInputConverter::convertSeeds
mkfit::TrackVec convertSeeds(const edm::View< TrajectorySeed > &seeds, const MkFitHitIndexMap &hitIndexMap, const TransientTrackingRecHitBuilder &ttrhBuilder, const MagneticField &mf) const
Definition: MkFitInputConverter.cc:196
trackerHitRTTI::isFromDet
bool isFromDet(TrackingRecHit const &hit)
Definition: trackerHitRTTI.h:35
MagneticField.h
TrackerDetSide.h
edm::EventSetup
Definition: EventSetup.h:57
edm::ESGetToken< TransientTrackingRecHitBuilder, TransientRecHitRecord >
edm::EventSetup::getData
bool getData(T &iHolder) const
Definition: EventSetup.h:113
TrackerDetSide::PosEndcap
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
eostools.move
def move(src, dest)
Definition: eostools.py:511
DetId::rawId
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
MkFitHitIndexMap
Definition: MkFitHitIndexMap.h:15
MkFitInputConverter::SMatrixSym66
ROOT::Math::SMatrix< float, 6, 6, ROOT::Math::MatRepSym< float, 6 > > SMatrixSym66
Definition: MkFitInputConverter.cc:64
MkFitInputConverter::stripRphiRecHitToken_
edm::EDGetTokenT< SiStripRecHit2DCollection > stripRphiRecHitToken_
Definition: MkFitInputConverter.cc:67
Exception
Definition: hltDiff.cc:246
TrajectoryStateTransform.h
MkFitInputConverter::produce
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
Definition: MkFitInputConverter.cc:108
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
TransientTrackingRecHitBuilder.h
LogTrace
#define LogTrace(id)
Definition: MessageLogger.h:671
TrajectoryStateOnSurface::globalParameters
const GlobalTrajectoryParameters & globalParameters() const
Definition: TrajectoryStateOnSurface.h:64
ParameterSet.h
MkFitInputConverter::convertHits
void convertHits(const HitCollection &hits, std::vector< mkfit::HitVec > &mkFitHits, MkFitHitIndexMap &hitIndexMap, int &totalHits, const TrackerTopology &ttopo, const TransientTrackingRecHitBuilder &ttrhBuilder, const mkfit::LayerNumberConverter &lnc) const
Definition: MkFitInputConverter.cc:136
MkFitHitIndexMap::MkFitHit
Definition: MkFitHitIndexMap.h:19
EDProducer.h
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
edm::Event
Definition: Event.h:73
MagneticField
Definition: MagneticField.h:19
TransientTrackingRecHitBuilder
Definition: TransientTrackingRecHitBuilder.h:6
edm::InputTag
Definition: InputTag.h:15
SurveyInfoScenario_cff.seed
seed
Definition: SurveyInfoScenario_cff.py:295
hit
Definition: SiStripHitEffFromCalibTree.cc:88
findQualityFiles.size
size
Write out results.
Definition: findQualityFiles.py:443
MkFitInputConverter::putToken_
edm::EDPutTokenT< MkFitInputWrapper > putToken_
Definition: MkFitInputConverter.cc:73