CMS 3D CMS Logo

AlignmentMuonHIPTrajectorySelector.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: AlignmentMuonHIPTrajectorySelector
4 // Class: AlignmentMuonHIPTrajectorySelector
5 //
13 //
14 // Original Author: Jim Pivarski
15 // Created: Wed Feb 20 10:56:46 CST 2008
16 // $Id: AlignmentMuonHIPTrajectorySelector.cc,v 1.4 2010/01/06 15:26:10 mussgill Exp $
17 //
18 //
19 
20 // system include files
21 #include <memory>
22 #include <map>
23 #include <fstream>
24 
25 // user include files
35 
49 #include "TH1F.h"
50 
51 //
52 // class decleration
53 //
54 
56 public:
58  ~AlignmentMuonHIPTrajectorySelector() override = default;
59 
60  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
61 
62 private:
63  void produce(edm::Event&, const edm::EventSetup&) override;
64 
65  // ---------- member data --------------------------------
67  const double m_minPt;
69  const int m_minTrackerDOF;
70  const double m_maxMuonResidual;
71 
72  const bool m_hists;
76 };
77 
78 //
79 // constants, enums and typedefs
80 //
81 
82 //
83 // static data member definitions
84 //
85 
86 //
87 // constructors and destructor
88 //
90  : m_input(iConfig.getParameter<edm::InputTag>("input")),
91  m_minPt(iConfig.getParameter<double>("minPt")),
92  m_maxTrackerForwardRedChi2(iConfig.getParameter<double>("maxTrackerForwardRedChi2")),
93  m_minTrackerDOF(iConfig.getParameter<int>("minTrackerDOF")),
94  m_maxMuonResidual(iConfig.getParameter<double>("maxMuonResidual")),
95  m_hists(iConfig.getParameter<bool>("hists")),
96  mapToken_(consumes<TrajTrackAssociationCollection>(m_input)),
97  m_pt(nullptr),
98  m_tracker_forwardredchi2(nullptr),
99  m_tracker_dof(nullptr) {
100  produces<TrajTrackAssociationCollection>();
101 
102  if (m_hists) {
104  m_pt = fs->make<TH1F>("pt", "Transverse momentum (GeV)", 100, 0., 100.);
106  fs->make<TH1F>("trackerForwardRedChi2", "forward-biased reduced chi2 in tracker", 100, 0., 5.);
107  m_tracker_dof = fs->make<TH1F>("trackerDOF", "DOF in tracker", 61, -0.5, 60.5);
108  m_resid_before = fs->make<TH1F>("residBefore", "muon residuals before cut (cm)", 100, -20, 20);
109  m_resid_after = fs->make<TH1F>("residAfter", "muon residuals after cut (cm)", 100, -20, 20);
110  }
111 }
112 
113 //
114 // member functions
115 //
118  desc.setUnknown();
119  descriptions.addDefault(desc);
120 }
121 
122 // ------------ method called to produce the data ------------
124  // input
125  const edm::Handle<TrajTrackAssociationCollection>& originalTrajTrackMap = iEvent.getHandle(mapToken_);
126 
127  // output
128  auto newTrajTrackMap = std::make_unique<TrajTrackAssociationCollection>();
129 
130  TrajectoryStateCombiner tsoscomb;
131 
132  for (TrajTrackAssociationCollection::const_iterator iPair = originalTrajTrackMap->begin();
133  iPair != originalTrajTrackMap->end();
134  ++iPair) {
135  if (m_hists) {
136  m_pt->Fill((*(*iPair).val).pt());
137  }
138 
139  if ((*(*iPair).val).pt() > m_minPt) {
140  std::vector<TrajectoryMeasurement> measurements = (*(*iPair).key).measurements();
141 
142  bool has_bad_residual = false;
143 
144  double tracker_forwardchi2 = 0.;
145  double tracker_dof = 0.;
146  for (std::vector<TrajectoryMeasurement>::const_iterator im = measurements.begin(); im != measurements.end();
147  ++im) {
148  const TrajectoryMeasurement& meas = *im;
149  auto hit = &(*meas.recHit());
150  const DetId id = hit->geographicalId();
151 
152  if (hit->isValid() && id.det() == DetId::Tracker) {
153  if (hit->dimension() == 1) {
154  double residual = meas.forwardPredictedState().localPosition().x() - hit->localPosition().x();
155  double error2 =
156  meas.forwardPredictedState().localError().positionError().xx() + hit->localPositionError().xx();
157 
158  tracker_forwardchi2 += residual * residual / error2;
159  tracker_dof += 1.;
160  } else if (hit->dimension() == 2) {
161  double residualx = meas.forwardPredictedState().localPosition().x() - hit->localPosition().x();
162  double residualy = meas.forwardPredictedState().localPosition().y() - hit->localPosition().y();
163  double errorxx2 =
164  meas.forwardPredictedState().localError().positionError().xx() + hit->localPositionError().xx();
165  double errorxy2 =
166  meas.forwardPredictedState().localError().positionError().xy() + hit->localPositionError().xy();
167  double erroryy2 =
168  meas.forwardPredictedState().localError().positionError().yy() + hit->localPositionError().yy();
169 
170  tracker_forwardchi2 +=
171  (residualx * residualx + residualy * residualy) / (errorxx2 + 2. * errorxy2 + erroryy2);
172  tracker_dof += 2.;
173  }
174  } // end if a tracker hit
175 
176  if (hit->isValid() && id.det() == DetId::Muon &&
177  (id.subdetId() == MuonSubdetId::DT || id.subdetId() == MuonSubdetId::CSC)) {
179  tsoscomb.combine(meas.forwardPredictedState(), meas.backwardPredictedState());
180  double residual = tsosc.localPosition().x() - hit->localPosition().x();
181  m_resid_before->Fill(residual);
182  if (fabs(residual) > m_maxMuonResidual) {
183  has_bad_residual = true;
184  }
185  } // end if a muon hit
186  }
187  tracker_dof -= 5.;
188  double tracker_forwardredchi2 = tracker_forwardchi2 / tracker_dof;
189 
190  if (m_hists) {
191  m_tracker_forwardredchi2->Fill(tracker_forwardredchi2);
192  m_tracker_dof->Fill(tracker_dof);
193 
194  for (std::vector<TrajectoryMeasurement>::const_iterator im = measurements.begin(); im != measurements.end();
195  ++im) {
196  const TrajectoryMeasurement& meas = *im;
197  auto hit = &(*meas.recHit());
198  const DetId id = hit->geographicalId();
199 
200  if (!has_bad_residual) {
201  if (hit->isValid() && id.det() == DetId::Muon &&
202  (id.subdetId() == MuonSubdetId::DT || id.subdetId() == MuonSubdetId::CSC)) {
204  tsoscomb.combine(meas.forwardPredictedState(), meas.backwardPredictedState());
205  double residual = tsosc.localPosition().x() - hit->localPosition().x();
206  m_resid_after->Fill(residual);
207  }
208  } // end if residuals pass cut
209  } // end second loop over hits
210  } // end if filling histograms
211 
212  if (tracker_forwardredchi2 < m_maxTrackerForwardRedChi2 && tracker_dof >= m_minTrackerDOF && !has_bad_residual) {
213  newTrajTrackMap->insert((*iPair).key, (*iPair).val);
214  } // end if passes tracker cuts
215  } // end if passes pT cut
216  } // end loop over original trajTrackMap
217 
218  // put it in the Event
219  iEvent.put(std::move(newTrajTrackMap));
220 }
221 
222 //define this as a plug-in
const LocalTrajectoryError & localError() const
TrajectoryStateOnSurface const & forwardPredictedState() const
Access to forward predicted state (from fitter or builder)
const edm::EDGetTokenT< TrajTrackAssociationCollection > mapToken_
LocalError positionError() const
const_iterator end() const
last iterator over the map (read only)
float yy() const
Definition: LocalError.h:24
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
int iEvent
Definition: GenABIO.cc:224
void addDefault(ParameterSetDescription const &psetDescription)
TrajectoryStateOnSurface const & backwardPredictedState() const
Access to backward predicted state (from smoother)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
float xy() const
Definition: LocalError.h:23
Definition: DetId.h:17
const_iterator begin() const
first iterator over the map (read only)
void produce(edm::Event &, const edm::EventSetup &) override
HLT enums.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
~AlignmentMuonHIPTrajectorySelector() override=default
static constexpr int DT
Definition: MuonSubdetId.h:11
static constexpr int CSC
Definition: MuonSubdetId.h:12
float xx() const
Definition: LocalError.h:22
def move(src, dest)
Definition: eostools.py:511
ConstRecHitPointer const & recHit() const