CMS 3D CMS Logo

AlignmentMonitorGeneric.cc
Go to the documentation of this file.
3 
7 #include "TObject.h"
8 
9 #include <TString.h>
10 
12  AlignmentMonitorBase(cfg, "AlignmentMonitorGeneric")
13 {
14 }
15 
17 {
18  std::vector<std::string> residNames; // names of residual histograms
19 
20  residNames.push_back("x hit residuals pos track");
21  residNames.push_back("x hit residuals neg track");
22  residNames.push_back("y hit residuals pos track");
23  residNames.push_back("y hit residuals neg track");
24 
25 
26  auto alignableObjectId =
28 
29  const std::vector<Alignable*>& alignables = pStore()->alignables();
30 
31  unsigned int nAlignable = alignables.size();
32  unsigned int nResidName = residNames.size();
33 
34  for (unsigned int i = 0; i < nAlignable; ++i)
35  {
36  const Alignable* ali = alignables[i];
37 
38  Hist1Ds& hists = m_resHists[ali];
39 
40  hists.resize(nResidName, nullptr);
41 
42  align::ID id = ali->id();
44 
45  for (unsigned int n = 0; n < nResidName; ++n)
46  {
47  const std::string& name = residNames[n];
48 
49  TString histName(name.c_str());
50  histName += Form("_%s_%d", alignableObjectId.idToString(type), id);
51  histName.ReplaceAll(" ", "");
52 
53  TString histTitle(name.c_str());
54  histTitle += Form(" for %s with ID %d (subdet %d)",
55  alignableObjectId.idToString(type),
56  id, DetId(id).subdetId());
57 
58  hists[n] = book1D(std::string("/iterN/") + std::string(name) + std::string("/"), std::string(histName.Data()), std::string(histTitle.Data()), nBin_, -5., 5.);
59  }
60  }
61 
62  m_trkHists.resize(6, nullptr);
63 
64  m_trkHists[0] = book1D("/iterN/", "pt" , "track p_{t} (GeV)" , nBin_, 0.0,100.0);
65  m_trkHists[1] = book1D("/iterN/", "eta" , "track #eta" , nBin_, - 3.0, 3.0);
66  m_trkHists[2] = book1D("/iterN/", "phi" , "track #phi" , nBin_, -M_PI, M_PI);
67  m_trkHists[3] = book1D("/iterN/", "d0" , "track d0 (cm)" , nBin_, -0.02, 0.02);
68  m_trkHists[4] = book1D("/iterN/", "dz" , "track dz (cm)" , nBin_, -20.0, 20.0);
69  m_trkHists[5] = book1D("/iterN/", "chi2", "track #chi^{2}/dof", nBin_, 0.0, 20.0);
70 
71 }
72 
74  const edm::EventSetup&,
76 {
77  TrajectoryStateCombiner tsoscomb;
78 
79  for (unsigned int t = 0; t < tracks.size(); ++t)
80  {
81  const reco::Track* track = tracks[t].second;
82 
83  float charge = tracks[t].second->charge();
84 
85  const std::vector<TrajectoryMeasurement>& meass
86  = tracks[t].first->measurements();
87 
88  for (unsigned int m = 0; m < meass.size(); ++m)
89  {
90  const TrajectoryMeasurement& meas = meass[m];
91  const TransientTrackingRecHit& hit = *meas.recHit();
92 
93  if ( hit.isValid() )
94  {
95  const Alignable* ali = pNavigator()->alignableFromDetId( hit.geographicalId() );
96 
97  while (ali) {
98  std::map<const Alignable*, Hist1Ds>::iterator h = m_resHists.find(ali);
99  if ( h != m_resHists.end() )
100  {
102 
104  LocalError err1 = tsos.localError().positionError();
105  LocalError err2 = hit.localPositionError(); // CPE+APE
106 
107  // subtract APEs from err2 (if existing) from covariance matrix
108  auto det = static_cast<const TrackerGeomDet*>(hit.det());
109  const auto localAPE = det->localAlignmentError();
110  if (localAPE.valid()) {
111  err2 = LocalError(err2.xx() - localAPE.xx(),
112  err2.xy() - localAPE.xy(),
113  err2.yy() - localAPE.yy());
114  }
115 
116  float errX = std::sqrt( err1.xx() + err2.xx() );
117  float errY = std::sqrt( err1.yy() + err2.yy() );
118 
119  h->second[charge > 0 ? 0 : 1]->Fill(res.x() / errX);
120  h->second[charge > 0 ? 2 : 3]->Fill(res.y() / errY);
121  }
122  ali = ali->mother();
123  }
124  }
125  }
126 
127  m_trkHists[0]->Fill( track->pt() );
128  m_trkHists[1]->Fill( track->eta() );
129  m_trkHists[2]->Fill( track->phi() );
130  m_trkHists[3]->Fill( track->d0() );
131  m_trkHists[4]->Fill( track->dz() );
132  m_trkHists[5]->Fill( track->normalizedChi2() );
133  }
134 }
135 
138 
AlignableMuon * pMuon()
type
Definition: HCALResponse.h:21
align::ID id() const
Return the ID of Alignable, i.e. DetId of &#39;first&#39; component GeomDet(Unit).
Definition: Alignable.h:189
float xx() const
Definition: LocalError.h:24
AlignableNavigator * pNavigator()
AlignmentParameterStore * pStore()
ConstRecHitPointer const & recHit() const
double d0() const
dxy parameter in perigee convention (d0 = -dxy)
Definition: TrackBase.h:597
std::map< const Alignable *, Hist1Ds > m_resHists
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:561
static AlignableObjectId commonObjectIdProvider(const AlignableObjectId &, const AlignableObjectId &)
uint32_t ID
Definition: Definitions.h:26
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:645
AlignableTracker * pTracker()
void event(const edm::Event &, const edm::EventSetup &, const ConstTrajTrackPairCollection &) override
Called for each event (by "run()"): may be reimplemented.
LocalError positionError() const
Definition: Electron.h:4
void book() override
Book or retrieve histograms; MUST be reimplemented.
std::vector< TH1F * > Hist1Ds
float xy() const
Definition: LocalError.h:25
int iEvent
Definition: GenABIO.cc:230
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:651
virtual StructureType alignableObjectId() const =0
Return the alignable type identifier.
float yy() const
Definition: LocalError.h:26
const GeomDet * det() const
T sqrt(T t)
Definition: SSEVec.h:18
double pt() const
track transverse momentum
Definition: TrackBase.h:621
AlignmentMonitorGeneric(const edm::ParameterSet &)
const LocalTrajectoryError & localError() const
std::vector< ConstTrajTrackPair > ConstTrajTrackPairCollection
virtual LocalPoint localPosition() const =0
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
TH1F * book1D(std::string dir, std::string name, std::string title, int nchX, double lowX, double highX)
TrajectoryStateOnSurface const & forwardPredictedState() const
Access to forward predicted state (from fitter or builder)
#define M_PI
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:609
Definition: DetId.h:18
bool isValid() const
#define DEFINE_EDM_PLUGIN(factory, type, name)
DetId geographicalId() const
virtual LocalError localPositionError() const =0
static const unsigned int nBin_
AlignableDetOrUnitPtr alignableFromDetId(const DetId &detid)
Returns AlignableDetOrUnitPtr corresponding to given DetId.
Alignable * mother() const
Return pointer to container alignable (if any)
Definition: Alignable.h:94
LocalError const & localAlignmentError() const
Return local alligment error.
const align::Alignables & alignables(void) const
get all alignables
TrajectoryStateOnSurface const & backwardPredictedState() const
Access to backward predicted state (from smoother)