CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MillePedeMonitor.h
Go to the documentation of this file.
1 #ifndef MILLEPEDEMONITOR_H
2 #define MILLEPEDEMONITOR_H
3 
13 
15 
17 
19 
21 
22 #include <vector>
23 #include <TString.h>
24 #include <TH1.h>
25 #include <TH2.h>
26 
27 class TH1;
28 class TH2;
29 class TDirectory;
30 class Trajectory;
31 
32 namespace reco {
33  class Track;
34 }
35 
36 class Alignable;
39 
40 /***************************************
41 ****************************************/
43 {
44  public:
45  // book histograms in constructor
46  MillePedeMonitor(const char *rootFile = "trackMonitor.root");
47  MillePedeMonitor(TDirectory *rootDir);
48  // writes histograms in destructor
49  ~MillePedeMonitor(); // non-virtual destructor: not intended to be parent class
50 
51  void fillTrack(const reco::Track *track);//, const Trajectory *traj);
52  void fillUsedTrack(const reco::Track *track, unsigned int nHitX, unsigned int nHitY);//, const Trajectory *traj);
55  const float *localDerivs, unsigned int nLocal,
56  const float *globalDerivs, unsigned int nGlobal, const int *labels);
58  const TrajectoryStateOnSurface &tsos, unsigned int nHit,
59  float residuum, float sigma, bool isY);
60  void fillFrameToFrame(const AlignableDetOrUnitPtr &aliDet, const Alignable *ali);
61 
63 
64  void fillPxbSurveyHistsChi2(const float &chi2);
65  void fillPxbSurveyHistsLocalPars(const float &a0, const float &a1, const float &S, const float &phi);
66 
67  private:
68  bool init(TDirectory *directory);
69  bool equidistLogBins(double* bins, int nBins, double first, double last) const;
70  void fillResidualHists(const std::vector<TH1*> &hists, float phiSensToNorm,
71  float residuum, float sigma);
72  void fillResidualHitHists(const std::vector<TH1*> &hists, float angle,
73  float residuum, float sigma, unsigned int nHit);
74  void fillTrack(const reco::Track *track, std::vector<TH1*> &trackHists1D,
75  std::vector<TH2*> &trackHists2D);
76 
77  template <class OBJECT_TYPE>
78  int GetIndex(const std::vector<OBJECT_TYPE*> &vec, const TString &name);
79  template <class OBJECT_TYPE>
80  std::vector<OBJECT_TYPE*> cloneHists(const std::vector<OBJECT_TYPE*> &orgs,
81  const TString &namAd, const TString &titAd) const;
82  template <class OBJECT_TYPE>
83  void addToDirectory(const std::vector<OBJECT_TYPE*> &objs, TDirectory *dir) const;
84 
85  TDirectory *myRootDir;
86  bool myDeleteDir;
87 
88  std::vector<TH1*> myTrackHists1D; // all input tracks
89  std::vector<TH2*> myTrackHists2D;
90  std::vector<TH1*> myUsedTrackHists1D; // tracks used, i.e. tranferred to pede
91  std::vector<TH2*> myUsedTrackHists2D;
92  std::vector<TH1*> myTrajectoryHists1D;
93  std::vector<TH2*> myTrajectoryHists2D;
94  std::vector<TH2*> myDerivHists2D;
95  std::vector<TH2*> myResidHists2D;
96  std::vector<std::vector<TH1*> > myResidHistsVec1DX;
97  std::vector<std::vector<TH1*> > myResidHistsVec1DY;
98  std::vector<TH1*> myResidHitHists1DX;
99  std::vector<TH1*> myResidHitHists1DY;
100  std::vector<TH2*> myFrame2FrameHists2D;
101  std::vector<TH1*> myCorrHists; // correlations
102  std::vector<TH1*> myPxbSurveyHists; // correlations
103 
104 };
105 
106 template <class OBJECT_TYPE>
107 int MillePedeMonitor::GetIndex(const std::vector<OBJECT_TYPE*> &vec, const TString &name)
108 {
109  int result = 0;
110  for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = vec.begin(), iterEnd = vec.end();
111  iter != iterEnd; ++iter, ++result) {
112  if (*iter && (*iter)->GetName() == name) return result;
113  }
114  edm::LogError("Alignment") << "@SUB=MillePedeMonitor::GetIndex" << " could not find " << name;
115  return -1;
116 }
117 
118 template <class OBJECT_TYPE>
119 std::vector<OBJECT_TYPE*> MillePedeMonitor::cloneHists(const std::vector<OBJECT_TYPE*> &orgs,
120  const TString &namAd,
121  const TString &titAd) const
122 {
123  // OBJECT_TYPE required to have methods Clone(const char*), GetName(), SetTitle(const char*) and GetTitle()
124  std::vector<OBJECT_TYPE*> result;
125  for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = orgs.begin(), iterEnd = orgs.end();
126  iter != iterEnd; ++iter) {
127  if (!(*iter)) continue;
128  result.push_back(static_cast<OBJECT_TYPE*>((*iter)->Clone(namAd + (*iter)->GetName())));
129  if (result.back()) result.back()->SetTitle((*iter)->GetTitle() + titAd);
130  else edm::LogError("Alignment") <<"@SUB=MillePedeMonitor::cloneHists"
131  << "out of memory?";
132  }
133 
134  return result;
135 }
136 
137 
138 template <class OBJECT_TYPE>
139 void MillePedeMonitor::addToDirectory(const std::vector<OBJECT_TYPE*> &obs,
140  TDirectory *dir) const
141 {
142  // OBJECT_TYPE is required to have method SetDirectory(TDirectory *dir)
143  for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = obs.begin(), iterEnd = obs.end();
144  iter != iterEnd; ++iter) {
145  if (*iter) (*iter)->SetDirectory(dir);
146  }
147 }
148 
149 #endif
std::vector< TH1 * > myTrackHists1D
int GetIndex(const std::vector< OBJECT_TYPE * > &vec, const TString &name)
void fillCorrelations2D(float corr, const TransientTrackingRecHit::ConstRecHitPointer &hit)
std::vector< TH2 * > myDerivHists2D
std::vector< std::vector< TH1 * > > myResidHistsVec1DX
void fillFrameToFrame(const AlignableDetOrUnitPtr &aliDet, const Alignable *ali)
void fillPxbSurveyHistsLocalPars(const float &a0, const float &a1, const float &S, const float &phi)
std::vector< TH1 * > myTrajectoryHists1D
std::vector< TH1 * > myResidHitHists1DY
std::vector< OBJECT_TYPE * > cloneHists(const std::vector< OBJECT_TYPE * > &orgs, const TString &namAd, const TString &titAd) const
void fillResiduals(const TransientTrackingRecHit::ConstRecHitPointer &recHit, const TrajectoryStateOnSurface &tsos, unsigned int nHit, float residuum, float sigma, bool isY)
void fillPxbSurveyHistsChi2(const float &chi2)
void addToDirectory(const std::vector< OBJECT_TYPE * > &objs, TDirectory *dir) const
std::vector< TH2 * > myTrackHists2D
MillePedeMonitor(const char *rootFile="trackMonitor.root")
std::vector< TH2 * > myResidHists2D
void fillUsedTrack(const reco::Track *track, unsigned int nHitX, unsigned int nHitY)
std::vector< TH1 * > myCorrHists
std::vector< TH1 * > myPxbSurveyHists
tuple result
Definition: query.py:137
TDirectory * myRootDir
std::vector< TH1 * > myUsedTrackHists1D
bool init(TDirectory *directory)
bool first
Definition: L1TdeRCT.cc:79
JetCorrectorParameters corr
Definition: classes.h:9
std::vector< TH2 * > myTrajectoryHists2D
tuple labels
Definition: L1TDQM_cfg.py:62
std::vector< std::vector< TH1 * > > myResidHistsVec1DY
[0]=all [1]=TPB [2]=TPE [3]=TIB [4]=TID [5]=TOB [6]=TEC
bool equidistLogBins(double *bins, int nBins, double first, double last) const
void fillResidualHists(const std::vector< TH1 * > &hists, float phiSensToNorm, float residuum, float sigma)
void fillTrack(const reco::Track *track)
std::vector< TH1 * > myResidHitHists1DX
[0]=all [1]=TPB [2]=TPE [3]=TIB [4]=TID [5]=TOB [6]=TEC
std::vector< TH2 * > myUsedTrackHists2D
void fillRefTrajectory(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr)
void fillDerivatives(const TransientTrackingRecHit::ConstRecHitPointer &recHit, const float *localDerivs, unsigned int nLocal, const float *globalDerivs, unsigned int nGlobal, const int *labels)
dbl *** dir
Definition: mlp_gen.cc:35
std::vector< TH2 * > myFrame2FrameHists2D
void fillResidualHitHists(const std::vector< TH1 * > &hists, float angle, float residuum, float sigma, unsigned int nHit)
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11
Definition: DDAxes.h:10