Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:56:29

0001 #ifndef MILLEPEDEMONITOR_H
0002 #define MILLEPEDEMONITOR_H
0003 
0004 /// \class MillePedeMonitor
0005 ///
0006 /// monitoring of MillePedeAlignmentAlgorithm and its input tracks
0007 ///
0008 ///  \author    : Gero Flucke
0009 ///  date       : October 2006
0010 ///  $Revision: 1.14 $
0011 ///  $Date: 2010/10/26 20:52:23 $
0012 ///  (last update by $Author: flucke $)
0013 
0014 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h"
0015 
0016 #include "Alignment/ReferenceTrajectories/interface/ReferenceTrajectoryBase.h"
0017 
0018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0019 
0020 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
0021 
0022 #include <vector>
0023 #include <array>
0024 
0025 #include <TString.h>
0026 #include <TH1.h>
0027 #include <TH2.h>
0028 
0029 class TH1;
0030 class TH2;
0031 class TDirectory;
0032 class Trajectory;
0033 
0034 class TrackerTopology;
0035 
0036 namespace reco {
0037   class Track;
0038 }
0039 
0040 class Alignable;
0041 class AlignableDetOrUnitPtr;
0042 class TrajectoryStateOnSurface;
0043 
0044 /***************************************
0045 ****************************************/
0046 class MillePedeMonitor {
0047 public:
0048   // book histograms in constructor
0049   MillePedeMonitor(const TrackerTopology *tTopo, const char *rootFile = "trackMonitor.root");
0050   MillePedeMonitor(TDirectory *rootDir, const TrackerTopology *tTopo);
0051   // writes histograms in destructor
0052   ~MillePedeMonitor();  // non-virtual destructor: not intended to be parent class
0053 
0054   void fillTrack(const reco::Track *track);                                              //, const Trajectory *traj);
0055   void fillUsedTrack(const reco::Track *track, unsigned int nHitX, unsigned int nHitY);  //, const Trajectory *traj);
0056   void fillRefTrajectory(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr);
0057   void fillDerivatives(const TransientTrackingRecHit::ConstRecHitPointer &recHit,
0058                        const float *localDerivs,
0059                        unsigned int nLocal,
0060                        const float *globalDerivs,
0061                        unsigned int nGlobal,
0062                        const int *labels);
0063   void fillResiduals(const TransientTrackingRecHit::ConstRecHitPointer &recHit,
0064                      const TrajectoryStateOnSurface &tsos,
0065                      unsigned int nHit,
0066                      float residuum,
0067                      float sigma,
0068                      bool isY);
0069   void fillFrameToFrame(const AlignableDetOrUnitPtr &aliDet, const Alignable *ali);
0070 
0071   void fillCorrelations2D(float corr, const TransientTrackingRecHit::ConstRecHitPointer &hit);
0072 
0073   void fillPxbSurveyHistsChi2(const float &chi2);
0074   void fillPxbSurveyHistsLocalPars(const float &a0, const float &a1, const float &S, const float &phi);
0075 
0076 private:
0077   bool init(TDirectory *directory);
0078   bool equidistLogBins(double *bins, int nBins, double first, double last) const;
0079   void fillResidualHists(const std::vector<TH1 *> &hists, float phiSensToNorm, float residuum, float sigma);
0080   void fillResidualHitHists(
0081       const std::vector<TH1 *> &hists, float angle, float residuum, float sigma, unsigned int nHit);
0082   void fillTrack(const reco::Track *track, std::vector<TH1 *> &trackHists1D, std::vector<TH2 *> &trackHists2D);
0083 
0084   template <class OBJECT_TYPE>
0085   int GetIndex(const std::vector<OBJECT_TYPE *> &vec, const TString &name);
0086   template <class OBJECT_TYPE>
0087   std::vector<OBJECT_TYPE *> cloneHists(const std::vector<OBJECT_TYPE *> &orgs,
0088                                         const TString &namAd,
0089                                         const TString &titAd) const;
0090   template <class OBJECT_TYPE>
0091   void addToDirectory(const std::vector<OBJECT_TYPE *> &objs, TDirectory *dir) const;
0092   template <typename T, size_t SIZE>
0093   std::array<int, SIZE> indexArray1D(const std::vector<T *> &hists, const char *title);
0094   template <typename T, size_t SIZE>
0095   std::array<std::array<int, SIZE>, SIZE> indexArray2D(const std::vector<T *> &hists, const char *title);
0096 
0097   TDirectory *myRootDir;
0098   bool myDeleteDir;
0099 
0100   std::vector<TH1 *> myTrackHists1D;  // all input tracks
0101   std::vector<TH2 *> myTrackHists2D;
0102   std::vector<TH1 *> myUsedTrackHists1D;  // tracks used, i.e. tranferred to pede
0103   std::vector<TH2 *> myUsedTrackHists2D;
0104   std::vector<TH1 *> myTrajectoryHists1D;
0105   std::vector<TH2 *> myTrajectoryHists2D;
0106   std::vector<TH2 *> myDerivHists2D;
0107   std::vector<TH2 *> myResidHists2D;
0108   std::vector<std::vector<TH1 *> > myResidHistsVec1DX;  ///[0]=all [1]=TPB [2]=TPE [3]=TIB [4]=TID [5]=TOB [6]=TEC
0109   std::vector<std::vector<TH1 *> > myResidHistsVec1DY;  ///[0]=all [1]=TPB [2]=TPE [3]=TIB [4]=TID [5]=TOB [6]=TEC
0110   std::vector<TH1 *> myResidHitHists1DX;
0111   std::vector<TH1 *> myResidHitHists1DY;
0112   std::vector<TH2 *> myFrame2FrameHists2D;
0113   std::vector<TH1 *> myCorrHists;       // correlations
0114   std::vector<TH1 *> myPxbSurveyHists;  // correlations
0115 
0116   const TrackerTopology *trackerTopology;
0117 };
0118 
0119 template <class OBJECT_TYPE>
0120 int MillePedeMonitor::GetIndex(const std::vector<OBJECT_TYPE *> &vec, const TString &name) {
0121   int result = 0;
0122   for (typename std::vector<OBJECT_TYPE *>::const_iterator iter = vec.begin(), iterEnd = vec.end(); iter != iterEnd;
0123        ++iter, ++result) {
0124     if (*iter && (*iter)->GetName() == name)
0125       return result;
0126   }
0127   edm::LogError("Alignment") << "@SUB=MillePedeMonitor::GetIndex"
0128                              << " could not find " << name;
0129   return -1;
0130 }
0131 
0132 template <class OBJECT_TYPE>
0133 std::vector<OBJECT_TYPE *> MillePedeMonitor::cloneHists(const std::vector<OBJECT_TYPE *> &orgs,
0134                                                         const TString &namAd,
0135                                                         const TString &titAd) const {
0136   // OBJECT_TYPE required to have methods Clone(const char*), GetName(), SetTitle(const char*) and GetTitle()
0137   std::vector<OBJECT_TYPE *> result;
0138   for (typename std::vector<OBJECT_TYPE *>::const_iterator iter = orgs.begin(), iterEnd = orgs.end(); iter != iterEnd;
0139        ++iter) {
0140     if (!(*iter))
0141       continue;
0142     result.push_back(static_cast<OBJECT_TYPE *>((*iter)->Clone(namAd + (*iter)->GetName())));
0143     if (result.back())
0144       result.back()->SetTitle((*iter)->GetTitle() + titAd);
0145     else
0146       edm::LogError("Alignment") << "@SUB=MillePedeMonitor::cloneHists"
0147                                  << "out of memory?";
0148   }
0149 
0150   return result;
0151 }
0152 
0153 template <class OBJECT_TYPE>
0154 void MillePedeMonitor::addToDirectory(const std::vector<OBJECT_TYPE *> &obs, TDirectory *dir) const {
0155   // OBJECT_TYPE is required to have method SetDirectory(TDirectory *dir)
0156   for (typename std::vector<OBJECT_TYPE *>::const_iterator iter = obs.begin(), iterEnd = obs.end(); iter != iterEnd;
0157        ++iter) {
0158     if (*iter)
0159       (*iter)->SetDirectory(dir);
0160   }
0161 }
0162 
0163 template <typename T, size_t SIZE>
0164 std::array<int, SIZE> MillePedeMonitor::indexArray1D(const std::vector<T *> &hists, const char *title) {
0165   std::array<int, SIZE> result{};
0166   for (size_t i = 0; i < SIZE; ++i) {
0167     result[i] = this->GetIndex(hists, Form(title, i));
0168   }
0169   return result;
0170 }
0171 
0172 template <typename T, size_t SIZE>
0173 std::array<std::array<int, SIZE>, SIZE> MillePedeMonitor::indexArray2D(const std::vector<T *> &hists,
0174                                                                        const char *title) {
0175   std::array<std::array<int, SIZE>, SIZE> result{};
0176   for (size_t i = 0; i < SIZE; ++i) {
0177     for (size_t j = 0; j < SIZE; ++j) {
0178       result[i][j] = this->GetIndex(hists, Form(title, i, j));
0179     }
0180   }
0181   return result;
0182 }
0183 
0184 #endif