File indexing completed on 2024-04-06 12:28:12
0001 #ifndef TkStripMeasurementDet_H
0002 #define TkStripMeasurementDet_H
0003
0004 #include "TrackingTools/MeasurementDet/interface/MeasurementDet.h"
0005 #include "RecoTracker/MeasurementDet/interface/MeasurementTrackerEvent.h"
0006 #include "RecoTracker/MeasurementDet/src/TkMeasurementDetSet.h"
0007
0008 #include "RecoLocalTracker/ClusterParameterEstimator/interface/StripClusterParameterEstimator.h"
0009 #include "RecoLocalTracker/SiStripRecHitConverter/interface/StripCPE.h"
0010
0011 #include "DataFormats/TrackerRecHit2D/interface/OmniClusterRef.h"
0012 #include "DataFormats/SiStripCluster/interface/SiStripClusterCollection.h"
0013 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
0014 #include "DataFormats/Common/interface/DetSetVector.h"
0015 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0016 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
0017 #include "DataFormats/Common/interface/Handle.h"
0018 #include "CondFormats/SiStripObjects/interface/SiStripNoises.h"
0019 #include "CondFormats/SiStripObjects/interface/SiStripBadStrip.h"
0020 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
0021 #include "RecoTracker/TransientTrackingRecHit/interface/TSiStripRecHit2DLocalPos.h"
0022
0023 #include "TrackingTools/DetLayers/interface/MeasurementEstimator.h"
0024 #include "RecoTracker/MeasurementDet/interface/ClusterFilterPayload.h"
0025
0026 #include <tuple>
0027
0028 class TrackingRecHit;
0029
0030 class TkStripMeasurementDet;
0031
0032 struct dso_hidden TkStripRecHitIter {
0033 using detset = edmNew::DetSet<SiStripCluster>;
0034 using new_const_iterator = detset::const_iterator;
0035
0036 TkStripRecHitIter() {}
0037 TkStripRecHitIter(const TkStripMeasurementDet& imdet,
0038 const TrajectoryStateOnSurface& itsos,
0039 const MeasurementTrackerEvent& idata)
0040 : mdet(&imdet), tsos(&itsos), data(&idata) {}
0041
0042 TkStripRecHitIter(new_const_iterator ci,
0043 new_const_iterator ce,
0044 const TkStripMeasurementDet& imdet,
0045 const TrajectoryStateOnSurface& itsos,
0046 const MeasurementTrackerEvent& idata)
0047 : mdet(&imdet), tsos(&itsos), data(&idata), clusterI(ci), clusterE(ce) {}
0048
0049 const TkStripMeasurementDet* mdet = nullptr;
0050 const TrajectoryStateOnSurface* tsos = nullptr;
0051 const MeasurementTrackerEvent* data = nullptr;
0052
0053 new_const_iterator clusterI;
0054 new_const_iterator clusterE;
0055
0056 inline SiStripRecHit2D buildHit() const;
0057 inline void advance();
0058
0059 public:
0060 bool empty() const { return clusterI == clusterE; }
0061
0062 bool operator==(TkStripRecHitIter const& rh) { return clusterI == rh.clusterI; }
0063 bool operator!=(TkStripRecHitIter const& rh) { return clusterI != rh.clusterI; }
0064 bool operator<(TkStripRecHitIter const& rh) { return clusterI < rh.clusterI; }
0065
0066 TkStripRecHitIter& operator++() {
0067 advance();
0068 return *this;
0069 }
0070
0071 SiStripRecHit2D operator*() const { return buildHit(); }
0072 };
0073
0074 class dso_hidden TkStripMeasurementDet final : public MeasurementDet {
0075 public:
0076 typedef StripClusterParameterEstimator::LocalValues LocalValues;
0077 typedef StripClusterParameterEstimator::VLocalValues VLocalValues;
0078
0079 typedef SiStripRecHit2D::ClusterRef SiStripClusterRef;
0080
0081 typedef edmNew::DetSet<SiStripCluster> detset;
0082 typedef detset::const_iterator new_const_iterator;
0083
0084 typedef std::vector<SiStripCluster>::const_iterator const_iterator;
0085
0086 ~TkStripMeasurementDet() override {}
0087
0088 TkStripMeasurementDet(const GeomDet* gdet, StMeasurementConditionSet& conditionSet);
0089
0090 void setIndex(int i) { index_ = i; }
0091
0092 void setEmpty(StMeasurementDetSet& theDets) const { theDets.setEmpty(index()); }
0093
0094 bool isEmpty(const StMeasurementDetSet& theDets) const { return theDets.empty(index()); }
0095
0096 int index() const { return index_; }
0097
0098 unsigned int rawId() const { return conditionSet().id(index()); }
0099 unsigned char subId() const { return conditionSet().subId(index()); }
0100
0101 const detset& theSet(const StMeasurementDetSet& theDets) const { return theDets.detSet(index()); }
0102 const detset& detSet(const StMeasurementDetSet& theDets) const { return theDets.detSet(index()); }
0103
0104
0105 bool isActive(const MeasurementTrackerEvent& data) const override { return data.stripData().isActive(index()); }
0106
0107
0108 bool hasBadComponents(const TrajectoryStateOnSurface& tsos, const MeasurementTrackerEvent& data) const override {
0109 return false;
0110 }
0111
0112 std::tuple<TkStripRecHitIter, TkStripRecHitIter> hitRange(const TrajectoryStateOnSurface&,
0113 const MeasurementTrackerEvent& data) const;
0114 void advance(TkStripRecHitIter& hi) const;
0115 SiStripRecHit2D hit(TkStripRecHitIter const& hi) const;
0116
0117 RecHitContainer recHits(const TrajectoryStateOnSurface&, const MeasurementTrackerEvent& data) const override;
0118
0119 bool empty(const MeasurementTrackerEvent& data) const;
0120
0121 void simpleRecHits(const TrajectoryStateOnSurface& ts,
0122 const MeasurementTrackerEvent& data,
0123 std::vector<SiStripRecHit2D>& result) const;
0124 bool simpleRecHits(const TrajectoryStateOnSurface& ts,
0125 const MeasurementEstimator& est,
0126 const MeasurementTrackerEvent& data,
0127 std::vector<SiStripRecHit2D>& result) const;
0128
0129
0130 bool recHits(SimpleHitContainer& result,
0131 const TrajectoryStateOnSurface& stateOnThisDet,
0132 const MeasurementEstimator&,
0133 const MeasurementTrackerEvent& data) const override;
0134
0135
0136 bool recHits(const TrajectoryStateOnSurface& stateOnThisDet,
0137 const MeasurementEstimator& est,
0138 const MeasurementTrackerEvent& data,
0139 RecHitContainer& result,
0140 std::vector<float>& diffs) const override;
0141
0142 bool measurements(const TrajectoryStateOnSurface& stateOnThisDet,
0143 const MeasurementEstimator& est,
0144 const MeasurementTrackerEvent& data,
0145 TempMeasurements& result) const override;
0146
0147 const StripGeomDetUnit& specificGeomDet() const { return static_cast<StripGeomDetUnit const&>(fastGeomDet()); }
0148
0149 template <class ClusterRefT>
0150 TrackingRecHit::RecHitPointer buildRecHit(const ClusterRefT& cluster, const TrajectoryStateOnSurface& ltp) const {
0151 const GeomDetUnit& gdu(specificGeomDet());
0152 LocalValues lv = cpe()->localParameters(*cluster, gdu, ltp);
0153 return std::make_shared<SiStripRecHit2D>(lv.first, lv.second, fastGeomDet(), cluster);
0154 }
0155
0156 template <class ClusterRefT>
0157 void buildRecHits(const ClusterRefT& cluster,
0158 const TrajectoryStateOnSurface& ltp,
0159 const RecHitContainer& _res) const {
0160 RecHitContainer res = _res;
0161 const GeomDetUnit& gdu(specificGeomDet());
0162 VLocalValues vlv = cpe()->localParametersV(*cluster, gdu, ltp);
0163 for (VLocalValues::const_iterator it = vlv.begin(); it != vlv.end(); ++it)
0164 res.push_back(std::make_shared<SiStripRecHit2D>(it->first, it->second, fastGeomDet(), cluster));
0165 }
0166
0167 template <class ClusterRefT>
0168 bool filteredRecHits(const ClusterRefT& cluster,
0169 StripCPE::AlgoParam const& cpepar,
0170 const TrajectoryStateOnSurface& ltp,
0171 const MeasurementEstimator& est,
0172 const std::vector<bool>& skipClusters,
0173 RecHitContainer& result,
0174 std::vector<float>& diffs) const {
0175 if (isMasked(*cluster))
0176 return true;
0177 if (!accept(cluster, skipClusters))
0178 return true;
0179 if (!est.preFilter(ltp, ClusterFilterPayload(rawId(), &*cluster)))
0180 return true;
0181 auto const& vl = cpe()->localParameters(*cluster, cpepar);
0182 SiStripRecHit2D recHit(vl.first,
0183 vl.second,
0184 fastGeomDet(),
0185 cluster);
0186 std::pair<bool, double> diffEst = est.estimate(ltp, recHit);
0187 LogDebug("TkStripMeasurementDet") << " chi2=" << diffEst.second;
0188 if (diffEst.first) {
0189 result.push_back(std::make_shared<SiStripRecHit2D>(recHit));
0190 diffs.push_back(diffEst.second);
0191 }
0192 return diffEst.first;
0193 }
0194
0195 template <class ClusterRefT>
0196 bool filteredRecHits(const ClusterRefT& cluster,
0197 StripCPE::AlgoParam const& cpepar,
0198 const TrajectoryStateOnSurface& ltp,
0199 const MeasurementEstimator& est,
0200 const std::vector<bool>& skipClusters,
0201 std::vector<SiStripRecHit2D>& result) const {
0202 if (isMasked(*cluster))
0203 return true;
0204 if (!accept(cluster, skipClusters))
0205 return true;
0206 if (!est.preFilter(ltp, ClusterFilterPayload(rawId(), &*cluster)))
0207 return true;
0208 auto const& vl = cpe()->localParameters(*cluster, cpepar);
0209 result.emplace_back(vl.first, vl.second, fastGeomDet(), cluster);
0210 std::pair<bool, double> diffEst = est.estimate(ltp, result.back());
0211 LogDebug("TkStripMeasurementDet") << " chi2=" << diffEst.second;
0212 if (!diffEst.first)
0213 result.pop_back();
0214 return diffEst.first;
0215 }
0216
0217
0218 void setActiveThisPeriod(StMeasurementDetSet& theDets, bool active) { conditionSet().setActive(index(), active); }
0219
0220
0221
0222 void setActiveThisEvent(StMeasurementDetSet& theDets, bool active) const {
0223 theDets.setActiveThisEvent(index(), active);
0224 }
0225
0226
0227 bool hasAllGoodChannels() const { return (!hasAny128StripBad()) && badStripBlocks().empty(); }
0228
0229
0230 void set128StripStatus(bool good, int idx = -1) { conditionSet().set128StripStatus(index(), good, idx); }
0231
0232 typedef StMeasurementConditionSet::BadStripCuts BadStripCuts;
0233
0234
0235 bool testStrips(float utraj, float uerr) const;
0236
0237 typedef StMeasurementConditionSet::BadStripBlock BadStripBlock;
0238
0239 std::vector<BadStripBlock>& getBadStripBlocks() { return conditionSet().getBadStripBlocks(index()); }
0240 std::vector<BadStripBlock> const& badStripBlocks() const { return conditionSet().badStripBlocks(index()); }
0241
0242 bool maskBad128StripBlocks() const { return conditionSet().maskBad128StripBlocks(); }
0243
0244 private:
0245 using AClusters = StripClusterParameterEstimator::AClusters;
0246 using ALocalValues = StripClusterParameterEstimator::ALocalValues;
0247
0248 int index_;
0249 StMeasurementConditionSet* theDetConditions;
0250 StMeasurementConditionSet& conditionSet() { return *theDetConditions; }
0251 const StMeasurementConditionSet& conditionSet() const { return *theDetConditions; }
0252
0253 const StripCPE* cpe() const { return static_cast<const StripCPE*>(conditionSet().stripCPE()); }
0254
0255
0256 int totalStrips() const { return conditionSet().totalStrips(index()); }
0257 BadStripCuts const& badStripCuts() const { return conditionSet().badStripCuts(index()); }
0258
0259 bool hasAny128StripBad() const { return conditionSet().hasAny128StripBad(index()); }
0260
0261 inline bool isMasked(const SiStripCluster& cluster) const { return conditionSet().isMasked(index(), cluster); }
0262
0263 void buildSimpleRecHits(AClusters const& clusters,
0264 const MeasurementTrackerEvent& data,
0265 const detset& detSet,
0266 const TrajectoryStateOnSurface& ltp,
0267 std::vector<SiStripRecHit2D>& res) const {
0268 const GeomDetUnit& gdu(specificGeomDet());
0269 declareDynArray(LocalValues, clusters.size(), alv);
0270 cpe()->localParameters(clusters, alv, gdu, ltp.localParameters());
0271 res.reserve(alv.size());
0272 for (unsigned int i = 0; i < clusters.size(); ++i)
0273 res.emplace_back(alv[i].first, alv[i].second, gdu, detSet.makeRefTo(data.stripData().handle(), clusters[i]));
0274 }
0275
0276 public:
0277 inline bool accept(SiStripClusterRef const& r, const std::vector<bool>& skipClusters) const {
0278 return accept(r.key(), skipClusters);
0279 }
0280
0281 inline bool accept(unsigned int key, const std::vector<bool>& skipClusters) const {
0282 if (skipClusters.empty())
0283 return true;
0284 if (key >= skipClusters.size()) {
0285 LogDebug("TkStripMeasurementDet")
0286 << key << " is larger than: " << skipClusters.size()
0287 << "\n This must be a new cluster, and therefore should not be skiped most likely.";
0288 return true;
0289 }
0290 return (not(skipClusters[key]));
0291 }
0292 };
0293
0294 inline SiStripRecHit2D TkStripRecHitIter::buildHit() const { return mdet->hit(*this); }
0295 inline void TkStripRecHitIter::advance() { mdet->advance(*this); }
0296
0297 #endif