CMS 3D CMS Logo

SiStripClusters2ApproxClusters.cc
Go to the documentation of this file.
32 
33 #include <vector>
34 #include <memory>
35 
37 public:
39  void produce(edm::Event&, const edm::EventSetup&) override;
40 
41  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
42 
43 private:
46 
47  unsigned int maxNSat;
48  static constexpr double subclusterWindow_ = .7;
49  static constexpr double seedCutMIPs_ = .35;
50  static constexpr double seedCutSN_ = 7.;
51  static constexpr double subclusterCutMIPs_ = .45;
52  static constexpr double subclusterCutSN_ = 12.;
53 
56 
58 
61 
64 
67 };
68 
70  inputClusters = conf.getParameter<edm::InputTag>("inputClusters");
71  maxNSat = conf.getParameter<unsigned int>("maxSaturatedStrips");
72 
73  clusterToken = consumes<edmNew::DetSetVector<SiStripCluster> >(inputClusters);
74 
75  beamSpot_ = conf.getParameter<edm::InputTag>("beamSpot");
76  beamSpotToken_ = consumes<reco::BeamSpot>(beamSpot_);
77 
79 
82 
83  csfLabel_ = conf.getParameter<std::string>("clusterShapeHitFilterLabel");
85 
87  produces<SiStripApproximateClusterCollection>();
88 }
89 
91  const auto& clusterCollection = event.get(clusterToken);
92  auto result = std::make_unique<SiStripApproximateClusterCollection>();
93  result->reserve(clusterCollection.size(), clusterCollection.dataSize());
94 
95  auto const beamSpotHandle = event.getHandle(beamSpotToken_);
96  auto const& bs = beamSpotHandle.isValid() ? *beamSpotHandle : reco::BeamSpot();
97  if (not beamSpotHandle.isValid()) {
98  edm::LogError("SiStripClusters2ApproxClusters")
99  << "didn't find a valid beamspot with label \"" << beamSpot_.encode() << "\" -> using (0,0,0)";
100  }
101 
102  const auto& tkGeom = &iSetup.getData(tkGeomToken_);
103  const auto& theFilter = &iSetup.getData(csfToken_);
104  const auto& theNoise_ = &iSetup.getData(stripNoiseToken_);
105 
106  for (const auto& detClusters : clusterCollection) {
107  auto ff = result->beginDet(detClusters.id());
108 
109  unsigned int detId = detClusters.id();
110  const GeomDet* det = tkGeom->idToDet(detId);
111  double nApvs = detInfo_.getNumberOfApvsAndStripLength(detId).first;
112  double stripLength = detInfo_.getNumberOfApvsAndStripLength(detId).second;
113  double barycenter_ypos = 0.5 * stripLength;
114 
115  const StripGeomDetUnit* stripDet = dynamic_cast<const StripGeomDetUnit*>(det);
116  float mip = 3.9 / (sistrip::MeVperADCStrip / stripDet->surface().bounds().thickness());
117 
118  for (const auto& cluster : detClusters) {
119  const LocalPoint& lp = LocalPoint(((cluster.barycenter() * 10 / (sistrip::STRIPS_PER_APV * nApvs)) -
120  ((stripDet->surface().bounds().width()) * 0.5f)),
121  barycenter_ypos - (0.5f * stripLength),
122  0.);
123  const GlobalPoint& gpos = det->surface().toGlobal(lp);
124  GlobalPoint beamspot(bs.position().x(), bs.position().y(), bs.position().z());
125  const GlobalVector& gdir = gpos - beamspot;
126  const LocalVector& ldir = det->toLocal(gdir);
127 
128  int hitStrips;
129  float hitPredPos;
130  bool usable = theFilter->getSizes(detId, cluster, lp, ldir, hitStrips, hitPredPos);
131  // (almost) same logic as in StripSubClusterShapeTrajectoryFilter
132  bool isTrivial = (std::abs(hitPredPos) < 2.f && hitStrips <= 2);
133 
134  if (!usable || isTrivial) {
135  ff.push_back(SiStripApproximateCluster(cluster, maxNSat, hitPredPos, true));
136  } else {
137  bool peakFilter = false;
138  SlidingPeakFinder pf(std::max<int>(2, std::ceil(std::abs(hitPredPos) + subclusterWindow_)));
139  float mipnorm = mip / std::abs(ldir.z());
140  PeakFinderTest test(mipnorm,
141  detId,
142  cluster.firstStrip(),
143  theNoise_,
144  seedCutMIPs_,
145  seedCutSN_,
148  peakFilter = pf.apply(cluster.amplitudes(), test);
149 
150  ff.push_back(SiStripApproximateCluster(cluster, maxNSat, hitPredPos, peakFilter));
151  }
152  }
153  }
154 
155  event.put(std::move(result));
156 }
157 
160  desc.add<edm::InputTag>("inputClusters", edm::InputTag("siStripClusters"));
161  desc.add<unsigned int>("maxSaturatedStrips", 3);
162  desc.add<std::string>("clusterShapeHitFilterLabel", "ClusterShapeHitFilter"); // add CSF label
163  desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot")); // add BeamSpot tag
164  descriptions.add("SiStripClusters2ApproxClusters", desc);
165 }
166 
constexpr int32_t ceil(float num)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::string fullPath() const
Definition: FileInPath.cc:161
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
T z() const
Definition: PV3DBase.h:61
std::string encode() const
Definition: InputTag.cc:159
SiStripClusters2ApproxClusters(const edm::ParameterSet &conf)
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
edm::ESHandle< SiStripNoises > theNoise_
Log< level::Error, false > LogError
edm::ESGetToken< ClusterShapeHitFilter, CkfComponentsRecord > csfToken_
virtual float thickness() const =0
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
SiStripDetInfo read(std::string filePath)
void produce(edm::Event &, const edm::EventSetup &) override
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
const std::pair< unsigned short, double > getNumberOfApvsAndStripLength(uint32_t detId) const
Constants and enumerated types for FED/FEC systems.
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::ESGetToken< SiStripNoises, SiStripNoisesRcd > stripNoiseToken_
static const uint16_t STRIPS_PER_APV
static constexpr float MeVperADCStrip
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
static constexpr char const *const kDefaultFile
virtual float width() const =0
edm::EDGetTokenT< edmNew::DetSetVector< SiStripCluster > > clusterToken
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1
const Bounds & bounds() const
Definition: Surface.h:87