CMS 3D CMS Logo

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 <array>
24 
25 #include <TString.h>
26 #include <TH1.h>
27 #include <TH2.h>
28 
29 class TH1;
30 class TH2;
31 class TDirectory;
32 class Trajectory;
33 
34 class TrackerTopology;
35 
36 namespace reco {
37  class Track;
38 }
39 
40 class Alignable;
43 
44 /***************************************
45 ****************************************/
47 {
48  public:
49  // book histograms in constructor
50  MillePedeMonitor(const TrackerTopology* tTopo,const char *rootFile = "trackMonitor.root");
51  MillePedeMonitor(TDirectory *rootDir, const TrackerTopology* tTopo);
52  // writes histograms in destructor
53  ~MillePedeMonitor(); // non-virtual destructor: not intended to be parent class
54 
55  void fillTrack(const reco::Track *track);//, const Trajectory *traj);
56  void fillUsedTrack(const reco::Track *track, unsigned int nHitX, unsigned int nHitY);//, const Trajectory *traj);
57  void fillRefTrajectory(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr);
58  void fillDerivatives(const TransientTrackingRecHit::ConstRecHitPointer &recHit,
59  const float *localDerivs, unsigned int nLocal,
60  const float *globalDerivs, unsigned int nGlobal, const int *labels);
61  void fillResiduals(const TransientTrackingRecHit::ConstRecHitPointer &recHit,
62  const TrajectoryStateOnSurface &tsos, unsigned int nHit,
63  float residuum, float sigma, bool isY);
64  void fillFrameToFrame(const AlignableDetOrUnitPtr &aliDet, const Alignable *ali);
65 
66  void fillCorrelations2D(float corr, const TransientTrackingRecHit::ConstRecHitPointer &hit);
67 
68  void fillPxbSurveyHistsChi2(const float &chi2);
69  void fillPxbSurveyHistsLocalPars(const float &a0, const float &a1, const float &S, const float &phi);
70 
71  private:
72  bool init(TDirectory *directory);
73  bool equidistLogBins(double* bins, int nBins, double first, double last) const;
74  void fillResidualHists(const std::vector<TH1*> &hists, float phiSensToNorm,
75  float residuum, float sigma);
76  void fillResidualHitHists(const std::vector<TH1*> &hists, float angle,
77  float residuum, float sigma, unsigned int nHit);
78  void fillTrack(const reco::Track *track, std::vector<TH1*> &trackHists1D,
79  std::vector<TH2*> &trackHists2D);
80 
81  template <class OBJECT_TYPE>
82  int GetIndex(const std::vector<OBJECT_TYPE*> &vec, const TString &name);
83  template <class OBJECT_TYPE>
84  std::vector<OBJECT_TYPE*> cloneHists(const std::vector<OBJECT_TYPE*> &orgs,
85  const TString &namAd, const TString &titAd) const;
86  template <class OBJECT_TYPE>
87  void addToDirectory(const std::vector<OBJECT_TYPE*> &objs, TDirectory *dir) const;
88  template <typename T, size_t SIZE>
89  std::array<int,SIZE> indexArray1D(const std::vector<T*>& hists, const char* title);
90  template <typename T, size_t SIZE>
91  std::array<std::array<int,SIZE>,SIZE> indexArray2D(const std::vector<T*>& hists,
92  const char* title);
93 
94  TDirectory *myRootDir;
95  bool myDeleteDir;
96 
97  std::vector<TH1*> myTrackHists1D; // all input tracks
98  std::vector<TH2*> myTrackHists2D;
99  std::vector<TH1*> myUsedTrackHists1D; // tracks used, i.e. tranferred to pede
100  std::vector<TH2*> myUsedTrackHists2D;
101  std::vector<TH1*> myTrajectoryHists1D;
102  std::vector<TH2*> myTrajectoryHists2D;
103  std::vector<TH2*> myDerivHists2D;
104  std::vector<TH2*> myResidHists2D;
105  std::vector<std::vector<TH1*> > myResidHistsVec1DX;
106  std::vector<std::vector<TH1*> > myResidHistsVec1DY;
107  std::vector<TH1*> myResidHitHists1DX;
108  std::vector<TH1*> myResidHitHists1DY;
109  std::vector<TH2*> myFrame2FrameHists2D;
110  std::vector<TH1*> myCorrHists; // correlations
111  std::vector<TH1*> myPxbSurveyHists; // correlations
112 
114 };
115 
116 template <class OBJECT_TYPE>
117 int MillePedeMonitor::GetIndex(const std::vector<OBJECT_TYPE*> &vec, const TString &name)
118 {
119  int result = 0;
120  for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = vec.begin(), iterEnd = vec.end();
121  iter != iterEnd; ++iter, ++result) {
122  if (*iter && (*iter)->GetName() == name) return result;
123  }
124  edm::LogError("Alignment") << "@SUB=MillePedeMonitor::GetIndex" << " could not find " << name;
125  return -1;
126 }
127 
128 template <class OBJECT_TYPE>
129 std::vector<OBJECT_TYPE*> MillePedeMonitor::cloneHists(const std::vector<OBJECT_TYPE*> &orgs,
130  const TString &namAd,
131  const TString &titAd) const
132 {
133  // OBJECT_TYPE required to have methods Clone(const char*), GetName(), SetTitle(const char*) and GetTitle()
134  std::vector<OBJECT_TYPE*> result;
135  for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = orgs.begin(), iterEnd = orgs.end();
136  iter != iterEnd; ++iter) {
137  if (!(*iter)) continue;
138  result.push_back(static_cast<OBJECT_TYPE*>((*iter)->Clone(namAd + (*iter)->GetName())));
139  if (result.back()) result.back()->SetTitle((*iter)->GetTitle() + titAd);
140  else edm::LogError("Alignment") <<"@SUB=MillePedeMonitor::cloneHists"
141  << "out of memory?";
142  }
143 
144  return result;
145 }
146 
147 
148 template <class OBJECT_TYPE>
149 void MillePedeMonitor::addToDirectory(const std::vector<OBJECT_TYPE*> &obs,
150  TDirectory *dir) const
151 {
152  // OBJECT_TYPE is required to have method SetDirectory(TDirectory *dir)
153  for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = obs.begin(), iterEnd = obs.end();
154  iter != iterEnd; ++iter) {
155  if (*iter) (*iter)->SetDirectory(dir);
156  }
157 }
158 
159 template <typename T, size_t SIZE>
160 std::array<int,SIZE>
161 MillePedeMonitor::indexArray1D(const std::vector<T*>& hists, const char* title) {
162  std::array<int,SIZE> result{};
163  for (size_t i = 0; i < SIZE; ++i) {
164  result[i] = this->GetIndex(hists, Form(title, i));
165  }
166  return result;
167 }
168 
169 template <typename T, size_t SIZE>
170 std::array<std::array<int,SIZE>,SIZE>
171 MillePedeMonitor::indexArray2D(const std::vector<T*>& hists, const char* title) {
172  std::array<std::array<int,SIZE>,SIZE> result{};
173  for (size_t i = 0; i < SIZE; ++i) {
174  for (size_t j = 0; j < SIZE; ++j) {
175  result[i][j] = this->GetIndex(hists, Form(title, i, j));
176  }
177  }
178  return result;
179 }
180 
181 #endif
std::array< std::array< int, SIZE >, SIZE > indexArray2D(const std::vector< T * > &hists, const char *title)
std::vector< TH1 * > myTrackHists1D
int GetIndex(const std::vector< OBJECT_TYPE * > &vec, const TString &name)
std::vector< TH2 * > myDerivHists2D
std::vector< std::vector< TH1 * > > myResidHistsVec1DX
int init
Definition: HydjetWrapper.h:67
std::vector< TH1 * > myTrajectoryHists1D
const TrackerTopology * trackerTopology
std::vector< TH1 * > myResidHitHists1DY
std::vector< OBJECT_TYPE * > cloneHists(const std::vector< OBJECT_TYPE * > &orgs, const TString &namAd, const TString &titAd) const
std::array< int, SIZE > indexArray1D(const std::vector< T * > &hists, const char *title)
void addToDirectory(const std::vector< OBJECT_TYPE * > &objs, TDirectory *dir) const
std::vector< TH2 * > myTrackHists2D
std::vector< TH2 * > myResidHists2D
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
std::vector< TH1 * > myCorrHists
std::vector< TH1 * > myPxbSurveyHists
TDirectory * myRootDir
std::vector< TH1 * > myUsedTrackHists1D
JetCorrectorParameters corr
Definition: classes.h:5
std::vector< TH2 * > myTrajectoryHists2D
std::vector< std::vector< TH1 * > > myResidHistsVec1DY
[0]=all [1]=TPB [2]=TPE [3]=TIB [4]=TID [5]=TOB [6]=TEC
double S(const TLorentzVector &, const TLorentzVector &)
Definition: Particle.cc:99
std::vector< TH1 * > myResidHitHists1DX
[0]=all [1]=TPB [2]=TPE [3]=TIB [4]=TID [5]=TOB [6]=TEC
fixed size matrix
std::vector< TH2 * > myUsedTrackHists2D
dbl *** dir
Definition: mlp_gen.cc:35
std::vector< TH2 * > myFrame2FrameHists2D
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11