CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlignmentMonitorGeneric.cc
Go to the documentation of this file.
3 
6 #include "TObject.h"
7 
8 #include <TString.h>
9 
11  AlignmentMonitorBase(cfg, "AlignmentMonitorGeneric")
12 {
13 }
14 
16 {
17  std::vector<std::string> residNames; // names of residual histograms
18 
19  residNames.push_back("x hit residuals pos track");
20  residNames.push_back("x hit residuals neg track");
21  residNames.push_back("y hit residuals pos track");
22  residNames.push_back("y hit residuals neg track");
23 
24  const std::vector<Alignable*>& alignables = pStore()->alignables();
25 
26  unsigned int nAlignable = alignables.size();
27  unsigned int nResidName = residNames.size();
28 
29  for (unsigned int i = 0; i < nAlignable; ++i)
30  {
31  const Alignable* ali = alignables[i];
32 
33  Hist1Ds& hists = m_resHists[ali];
34 
35  hists.resize(nResidName, 0);
36 
37  align::ID id = ali->id();
39 
40  for (unsigned int n = 0; n < nResidName; ++n)
41  {
42  const std::string& name = residNames[n];
43 
44  TString histName(name.c_str());
45  histName += Form("_%s_%d", AlignableObjectId::idToString(type), id);
46  histName.ReplaceAll(" ", "");
47 
48  TString histTitle(name.c_str());
49  histTitle += Form(" for %s with ID %d (subdet %d)",
51  id, DetId(id).subdetId());
52 
53  hists[n] = book1D(std::string("/iterN/") + std::string(name) + std::string("/"), std::string(histName), std::string(histTitle), nBin_, -5., 5.);
54  }
55  }
56 
57  m_trkHists.resize(6, 0);
58 
59  m_trkHists[0] = book1D("/iterN/", "pt" , "track p_{t} (GeV)" , nBin_, 0.0,100.0);
60  m_trkHists[1] = book1D("/iterN/", "eta" , "track #eta" , nBin_, - 3.0, 3.0);
61  m_trkHists[2] = book1D("/iterN/", "phi" , "track #phi" , nBin_, -M_PI, M_PI);
62  m_trkHists[3] = book1D("/iterN/", "d0" , "track d0 (cm)" , nBin_, -0.02, 0.02);
63  m_trkHists[4] = book1D("/iterN/", "dz" , "track dz (cm)" , nBin_, -20.0, 20.0);
64  m_trkHists[5] = book1D("/iterN/", "chi2", "track #chi^{2}/dof", nBin_, 0.0, 20.0);
65 
66 }
67 
69  const edm::EventSetup&,
71 {
72  TrajectoryStateCombiner tsoscomb;
73 
74  for (unsigned int t = 0; t < tracks.size(); ++t)
75  {
76  const reco::Track* track = tracks[t].second;
77 
78  float charge = tracks[t].second->charge();
79 
80  const std::vector<TrajectoryMeasurement>& meass
81  = tracks[t].first->measurements();
82 
83  for (unsigned int m = 0; m < meass.size(); ++m)
84  {
85  const TrajectoryMeasurement& meas = meass[m];
86  const TransientTrackingRecHit& hit = *meas.recHit();
87 
88  if ( hit.isValid() )
89  {
90  const Alignable* ali = pNavigator()->alignableFromDetId( hit.geographicalId() );
91 
92  while (ali) {
93  std::map<const Alignable*, Hist1Ds>::iterator h = m_resHists.find(ali);
94  if ( h != m_resHists.end() )
95  {
97 
98  align::LocalVector res = tsos.localPosition() - hit.localPosition();
99  LocalError err1 = tsos.localError().positionError();
100  LocalError err2 = hit.localPositionError();
101 
102  float errX = std::sqrt( err1.xx() + err2.xx() );
103  float errY = std::sqrt( err1.yy() + err2.yy() );
104 
105  h->second[charge > 0 ? 0 : 1]->Fill(res.x() / errX);
106  h->second[charge > 0 ? 2 : 3]->Fill(res.y() / errY);
107  }
108  ali = ali->mother();
109  }
110  }
111  }
112 
113  m_trkHists[0]->Fill( track->pt() );
114  m_trkHists[1]->Fill( track->eta() );
115  m_trkHists[2]->Fill( track->phi() );
116  m_trkHists[3]->Fill( track->d0() );
117  m_trkHists[4]->Fill( track->dz() );
118  m_trkHists[5]->Fill( track->normalizedChi2() );
119  }
120 }
121 
124 
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:180
int i
Definition: DBlmapReader.cc:9
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:121
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:109
uint32_t ID
Definition: Definitions.h:26
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:137
double charge(const std::vector< uint8_t > &Ampls)
LocalError positionError() const
std::vector< TH1F * > Hist1Ds
int iEvent
Definition: GenABIO.cc:243
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:139
float yy() const
Definition: LocalError.h:26
T sqrt(T t)
Definition: SSEVec.h:48
double pt() const
track transverse momentum
Definition: TrackBase.h:129
virtual StructureType alignableObjectId() const =0
Return the alignable type identifier.
AlignmentMonitorGeneric(const edm::ParameterSet &)
const LocalTrajectoryError & localError() const
std::vector< ConstTrajTrackPair > ConstTrajTrackPairCollection
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)
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:125
virtual void event(const edm::Event &, const edm::EventSetup &, const ConstTrajTrackPairCollection &)
Called for each event (by &quot;run()&quot;): may be reimplemented.
Definition: DetId.h:18
#define M_PI
Definition: BFit3D.cc:3
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
tuple tracks
Definition: testEve_cfg.py:39
virtual LocalError localPositionError() const =0
bool isValid() const
virtual void book()
Book or retrieve histograms; MUST be reimplemented.
#define DEFINE_EDM_PLUGIN(factory, type, name)
DetId geographicalId() const
static const char * idToString(align::StructureType type)
static const unsigned int nBin_
virtual LocalPoint localPosition() const =0
AlignableDetOrUnitPtr alignableFromDetId(const DetId &detid)
Returns AlignableDetOrUnitPtr corresponding to given DetId.
Alignable * mother() const
Return pointer to container alignable (if any)
Definition: Alignable.h:85
const align::Alignables & alignables(void) const
get all alignables
TrajectoryStateOnSurface const & backwardPredictedState() const
Access to backward predicted state (from smoother)