CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
MkFitEventOfHitsProducer.cc
Go to the documentation of this file.
2 
6 
9 
12 
17 
20 
26 
27 // mkFit includes
28 #include "mkFit/HitStructures.h"
29 #include "mkFit/MkStdSeqs.h"
30 #include "LayerNumberConverter.h"
31 
33 public:
35  ~MkFitEventOfHitsProducer() override = default;
36 
37  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
38 
39 private:
40  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
41 
42  void fill(const std::vector<const TrackingRecHit*>& hits,
43  mkfit::EventOfHits& eventOfHits,
44  const MkFitGeometry& mkFitGeom) const;
45 
56  const bool usePixelQualityDB_;
58 };
59 
61  : beamSpotToken_{consumes(iConfig.getParameter<edm::InputTag>("beamSpot"))},
62  pixelHitsToken_{consumes(iConfig.getParameter<edm::InputTag>("pixelHits"))},
63  stripHitsToken_{consumes(iConfig.getParameter<edm::InputTag>("stripHits"))},
64  pixelClusterIndexToHitToken_{consumes(iConfig.getParameter<edm::InputTag>("pixelHits"))},
65  stripClusterIndexToHitToken_{consumes(iConfig.getParameter<edm::InputTag>("stripHits"))},
66  mkFitGeomToken_{esConsumes()},
67  putToken_{produces<MkFitEventOfHits>()},
68  usePixelQualityDB_{iConfig.getParameter<bool>("usePixelQualityDB")},
69  useStripStripQualityDB_{iConfig.getParameter<bool>("useStripStripQualityDB")} {
70  if (useStripStripQualityDB_ || usePixelQualityDB_)
71  geomToken_ = esConsumes();
72  if (usePixelQualityDB_) {
74  }
75  if (useStripStripQualityDB_) {
76  stripQualityToken_ = esConsumes();
77  }
78 }
79 
82 
83  desc.add("beamSpot", edm::InputTag{"offlineBeamSpot"});
84  desc.add("pixelHits", edm::InputTag{"mkFitSiPixelHits"});
85  desc.add("stripHits", edm::InputTag{"mkFitSiStripHits"});
86  desc.add("usePixelQualityDB", true)->setComment("Use SiPixelQuality DB information");
87  desc.add("useStripStripQualityDB", true)->setComment("Use SiStrip quality DB information");
88 
89  descriptions.addWithDefaultLabel(desc);
90 }
91 
93  const auto& pixelHits = iEvent.get(pixelHitsToken_);
94  const auto& stripHits = iEvent.get(stripHitsToken_);
95  const auto& mkFitGeom = iSetup.getData(mkFitGeomToken_);
96 
97  auto eventOfHits = std::make_unique<mkfit::EventOfHits>(mkFitGeom.trackerInfo());
98  mkfit::StdSeq::Cmssw_LoadHits_Begin(*eventOfHits, {&pixelHits.hits(), &stripHits.hits()});
99 
101  std::vector<mkfit::DeadVec> deadvectors(mkFitGeom.layerNumberConverter().nLayers());
102  const auto& trackerGeom = iSetup.getData(geomToken_);
103 
104  if (usePixelQualityDB_) {
105  const auto& pixelQuality = iSetup.getData(pixelQualityToken_);
106  const auto& badPixels = pixelQuality.getBadComponentList();
107  for (const auto& bp : badPixels) {
108  const DetId detid(bp.DetID);
109  const auto& surf = trackerGeom.idToDet(detid)->surface();
110  bool isBarrel = (mkFitGeom.topology()->side(detid) == static_cast<unsigned>(TrackerDetSide::Barrel));
111  const auto ilay = mkFitGeom.mkFitLayerNumber(detid);
112  const auto q1 = isBarrel ? surf.zSpan().first : surf.rSpan().first;
113  const auto q2 = isBarrel ? surf.zSpan().second : surf.rSpan().second;
114  if (bp.errorType == 0)
115  deadvectors[ilay].push_back({surf.phiSpan().first, surf.phiSpan().second, q1, q2});
116  }
117  }
118 
120  const auto& siStripQuality = iSetup.getData(stripQualityToken_);
121  const auto& badStrips = siStripQuality.getBadComponentList();
122  for (const auto& bs : badStrips) {
123  const DetId detid(bs.detid);
124  const auto& surf = trackerGeom.idToDet(detid)->surface();
125  bool isBarrel = (mkFitGeom.topology()->side(detid) == static_cast<unsigned>(TrackerDetSide::Barrel));
126  const auto ilay = mkFitGeom.mkFitLayerNumber(detid);
127  const auto q1 = isBarrel ? surf.zSpan().first : surf.rSpan().first;
128  const auto q2 = isBarrel ? surf.zSpan().second : surf.rSpan().second;
129  if (bs.BadModule)
130  deadvectors[ilay].push_back({surf.phiSpan().first, surf.phiSpan().second, q1, q2});
131  else { //assume that BadApvs are filled in sync with BadFibers
132  auto const& topo = dynamic_cast<const StripTopology&>(trackerGeom.idToDet(detid)->topology());
133  int firstApv = -1;
134  int lastApv = -1;
135 
136  auto addRangeAPV = [&topo, &surf, &q1, &q2](int first, int last, mkfit::DeadVec& dv) {
137  auto const firstPoint = surf.toGlobal(topo.localPosition(first * sistrip::STRIPS_PER_APV));
138  auto const lastPoint = surf.toGlobal(topo.localPosition((last + 1) * sistrip::STRIPS_PER_APV));
139  float phi1 = firstPoint.phi();
140  float phi2 = lastPoint.phi();
141  if (reco::deltaPhi(phi1, phi2) > 0)
142  std::swap(phi1, phi2);
143  LogTrace("SiStripBadComponents")
144  << "insert bad range " << first << " to " << last << " " << phi1 << " " << phi2;
145  dv.push_back({phi1, phi2, q1, q2});
146  };
147 
148  const int nApvs = topo.nstrips() / sistrip::STRIPS_PER_APV;
149  for (int apv = 0; apv < nApvs; ++apv) {
150  const bool isBad = bs.BadApvs & (1 << apv);
151  if (isBad) {
152  LogTrace("SiStripBadComponents") << "bad apv " << apv << " on " << bs.detid;
153  if (lastApv == -1) {
154  firstApv = apv;
155  lastApv = apv;
156  } else if (lastApv + 1 == apv)
157  lastApv++;
158 
159  if (apv + 1 == nApvs)
160  addRangeAPV(firstApv, lastApv, deadvectors[ilay]);
161  } else if (firstApv != -1) {
162  addRangeAPV(firstApv, lastApv, deadvectors[ilay]);
163  //and reset
164  firstApv = -1;
165  lastApv = -1;
166  }
167  }
168  }
169  }
170  }
171  mkfit::StdSeq::LoadDeads(*eventOfHits, deadvectors);
172  }
173 
174  fill(iEvent.get(pixelClusterIndexToHitToken_).hits(), *eventOfHits, mkFitGeom);
175  fill(iEvent.get(stripClusterIndexToHitToken_).hits(), *eventOfHits, mkFitGeom);
176 
177  mkfit::StdSeq::Cmssw_LoadHits_End(*eventOfHits);
178 
179  auto const bs = iEvent.get(beamSpotToken_);
180  eventOfHits->SetBeamSpot(
181  mkfit::BeamSpot(bs.x0(), bs.y0(), bs.z0(), bs.sigmaZ(), bs.BeamWidthX(), bs.BeamWidthY(), bs.dxdz(), bs.dydz()));
182 
183  iEvent.emplace(putToken_, std::move(eventOfHits));
184 }
185 
186 void MkFitEventOfHitsProducer::fill(const std::vector<const TrackingRecHit*>& hits,
187  mkfit::EventOfHits& eventOfHits,
188  const MkFitGeometry& mkFitGeom) const {
189  for (int i = 0, end = hits.size(); i < end; ++i) {
190  const auto* hit = hits[i];
191  if (hit != nullptr) {
192  const auto ilay = mkFitGeom.mkFitLayerNumber(hit->geographicalId());
193  eventOfHits[ilay].RegisterHit(i);
194  }
195  }
196 }
197 
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
void setComment(std::string const &value)
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
MkFitEventOfHitsProducer(edm::ParameterSet const &iConfig)
~MkFitEventOfHitsProducer() override=default
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
#define LogTrace(id)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
bool getData(T &iHolder) const
Definition: EventSetup.h:128
const edm::EDPutTokenT< MkFitEventOfHits > putToken_
int iEvent
Definition: GenABIO.cc:224
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::ESGetToken< SiPixelQuality, SiPixelQualityRcd > pixelQualityToken_
def move
Definition: eostools.py:511
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
const edm::EDGetTokenT< MkFitHitWrapper > stripHitsToken_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: DetId.h:17
pixelQualityToken_(iC.esConsumes())
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > geomToken_
OrphanHandle< PROD > emplace(EDPutTokenT< PROD > token, Args &&...args)
puts a new product
Definition: Event.h:433
const edm::EDGetTokenT< MkFitClusterIndexToHit > pixelClusterIndexToHitToken_
int mkFitLayerNumber(DetId detId) const
const edm::ESGetToken< MkFitGeometry, TrackerRecoGeometryRecord > mkFitGeomToken_
const edm::EDGetTokenT< MkFitClusterIndexToHit > stripClusterIndexToHitToken_
Constants and enumerated types for FED/FEC systems.
const edm::EDGetTokenT< MkFitHitWrapper > pixelHitsToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::ESGetToken< SiStripQuality, SiStripQualityRcd > stripQualityToken_
static const uint16_t STRIPS_PER_APV
void fill(const std::vector< const TrackingRecHit * > &hits, mkfit::EventOfHits &eventOfHits, const MkFitGeometry &mkFitGeom) const
string end
Definition: dataset.py:937
tuple last
Definition: dqmdumpme.py:56
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283