CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Chi2ChargeMeasurementEstimatorESProducer.cc
Go to the documentation of this file.
1 
12 
14 
15 #include<limits>
16 
18 
20 
21 
22 namespace {
23 
25 public:
26 
32  template<typename... Args>
34  float minGoodPixelCharge, float minGoodStripCharge,
36  Args && ...args) :
38  minGoodPixelCharge_(minGoodPixelCharge),
39  minGoodStripCharge_(minGoodStripCharge),
40  pTChargeCutThreshold2_( pTChargeCutThreshold>=0.? pTChargeCutThreshold*pTChargeCutThreshold :std::numeric_limits<float>::max())
41  {}
42 
43 
44  bool preFilter(const TrajectoryStateOnSurface& ts,
45  const MeasurementEstimator::OpaquePayload & opay) const override;
46 
47 
48  virtual Chi2ChargeMeasurementEstimator* clone() const {
49  return new Chi2ChargeMeasurementEstimator(*this);
50  }
51 private:
52 
53  const float minGoodPixelCharge_;
54  const float minGoodStripCharge_;
55  const float pTChargeCutThreshold2_;
56 
57  bool checkClusterCharge(DetId id, SiStripCluster const & cluster, const TrajectoryStateOnSurface& ts) const {
58  return siStripClusterTools::chargePerCM(id, cluster, ts.localParameters() ) > minGoodStripCharge_;
59 
60  }
61 
62 
63 };
64 
66  const MeasurementEstimator::OpaquePayload & opay) const {
67 
68  // what we got?
69  if (opay.tag != ClusterFilterPayload::myTag) return true; // not mine...
70 
71  auto const & clf = reinterpret_cast<ClusterFilterPayload const &>(opay);
72 
73  if (ts.globalMomentum().perp2()>pTChargeCutThreshold2_) return true;
74 
75  DetId detid = clf.detId;
76  uint32_t subdet = detid.subdetId();
77 
78  if (subdet>2) {
79  return checkClusterCharge(detid, *clf.cluster[0],ts) && ( nullptr==clf.cluster[1] || checkClusterCharge(detid, *clf.cluster[1],ts) ) ;
80 
81  }
82 
83  /* pixel charge not implemented as not used...
84  auto const & thit = static_cast<const SiPixelRecHit &>(hit);
85  thit.cluster()->charge() ...
86 
87  */
88 
89  return true;
90 }
91 
92 }
93 
94 
95 
99 #include <boost/shared_ptr.hpp>
100 
102 
103 
104 namespace {
105 
106 
107 class Chi2ChargeMeasurementEstimatorESProducer: public edm::ESProducer{
108  public:
109  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
110  Chi2ChargeMeasurementEstimatorESProducer(const edm::ParameterSet & p);
111  virtual ~Chi2ChargeMeasurementEstimatorESProducer();
112  boost::shared_ptr<Chi2MeasurementEstimatorBase> produce(const TrackingComponentsRecord &);
113 
114  private:
115  boost::shared_ptr<Chi2MeasurementEstimatorBase> m_estimator;
116  const edm::ParameterSet m_pset;
117 };
118 
119 void
120 Chi2ChargeMeasurementEstimatorESProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
121 
123  desc.add<std::string>("ComponentName","Chi2Charge");
124  desc.add<double>("pTChargeCutThreshold",-1.);
126  desc.add<edm::ParameterSetDescription>("clusterChargeCut", descCCC);
127  descriptions.add("Chi2ChargeMeasurementEstimatorDefault", desc);
128 }
129 
130 
131 
132 Chi2ChargeMeasurementEstimatorESProducer::Chi2ChargeMeasurementEstimatorESProducer(const edm::ParameterSet & pset) :
133  m_pset(pset)
134 {
135  std::string const & myname = pset.getParameter<std::string>("ComponentName");
136  setWhatProduced(this,myname);
137 }
138 
139 Chi2ChargeMeasurementEstimatorESProducer::~Chi2ChargeMeasurementEstimatorESProducer() {}
140 
141 boost::shared_ptr<Chi2MeasurementEstimatorBase>
142 Chi2ChargeMeasurementEstimatorESProducer::produce(const TrackingComponentsRecord & iRecord){
143 
144  auto maxChi2 = m_pset.getParameter<double>("MaxChi2");
145  auto nSigma = m_pset.getParameter<double>("nSigma");
146  auto maxDis = m_pset.getParameter<double>("MaxDisplacement");
147  auto maxSag = m_pset.getParameter<double>("MaxSagitta");
148  auto minTol = m_pset.getParameter<double>("MinimalTolerance");
149  auto minpt = m_pset.getParameter<double>("MinPtForHitRecoveryInGluedDet");
150  auto minGoodPixelCharge = 0;
151  auto minGoodStripCharge = clusterChargeCut(m_pset);
152  auto pTChargeCutThreshold= m_pset.getParameter<double>("pTChargeCutThreshold");
153 
154  m_estimator = boost::shared_ptr<Chi2MeasurementEstimatorBase>(
156  maxChi2,nSigma, maxDis, maxSag, minTol,minpt) );
157 
158  return m_estimator;
159 }
160 
161 }
162 
163 
165 DEFINE_FWK_EVENTSETUP_MODULE(Chi2ChargeMeasurementEstimatorESProducer);
166 
T getParameter(std::string const &) const
virtual Chi2MeasurementEstimator * clone() const
virtual bool preFilter(const TrajectoryStateOnSurface &, OpaquePayload const &) const
const LocalTrajectoryParameters & localParameters() const
float chargePerCM(DetId detid, Iter a, Iter b)
float clusterChargeCut(const edm::ParameterSet &conf, const char *name="clusterChargeCut")
T perp2() const
Definition: PV3DBase.h:71
edm::ParameterSetDescription getFilledConfigurationDescription4CCC()
tuple preFilter
Definition: gather_cfg.py:50
edm::ParameterSetDescription getFilledConfigurationDescription()
ParameterDescriptionBase * add(U const &iLabel, T const &value)
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
Definition: DetId.h:18
#define DEFINE_FWK_EVENTSETUP_MODULE(type)
Definition: ModuleFactory.h:60
GlobalVector globalMomentum() const
static constexpr int myTag