CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/Alignment/MillePedeAlignmentAlgorithm/interface/MillePedeMonitor.h

Go to the documentation of this file.
00001 #ifndef MILLEPEDEMONITOR_H
00002 #define MILLEPEDEMONITOR_H
00003 
00013 
00014 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h"
00015 
00016 #include "Alignment/ReferenceTrajectories/interface/ReferenceTrajectoryBase.h"
00017 
00018 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00019 
00020 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
00021 
00022 #include <vector>
00023 #include <TString.h>
00024 #include <TH1.h>
00025 #include <TH2.h>
00026 
00027 class TH1;
00028 class TH2;
00029 class TDirectory;
00030 class Trajectory;
00031 
00032 namespace reco {
00033   class Track;
00034 }
00035 
00036 class Alignable;
00037 class AlignableDetOrUnitPtr;
00038 class TrajectoryStateOnSurface;
00039 
00040 /***************************************
00041 ****************************************/
00042 class MillePedeMonitor
00043 {
00044  public:
00045   // book histograms in constructor
00046   MillePedeMonitor(const char *rootFile = "trackMonitor.root");
00047   MillePedeMonitor(TDirectory *rootDir);
00048   // writes histograms in destructor
00049   ~MillePedeMonitor(); // non-virtual destructor: not intended to be parent class
00050 
00051   void fillTrack(const reco::Track *track);//, const Trajectory *traj);
00052   void fillUsedTrack(const reco::Track *track, unsigned int nHitX, unsigned int nHitY);//, const Trajectory *traj);
00053   void fillRefTrajectory(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr);
00054   void fillDerivatives(const TransientTrackingRecHit::ConstRecHitPointer &recHit,
00055                        const float *localDerivs, unsigned int nLocal,
00056                        const float *globalDerivs, unsigned int nGlobal, const int *labels);
00057   void fillResiduals(const TransientTrackingRecHit::ConstRecHitPointer &recHit,
00058                      const TrajectoryStateOnSurface &tsos, unsigned int nHit,
00059                      float residuum, float sigma, bool isY);
00060   void fillFrameToFrame(const AlignableDetOrUnitPtr &aliDet, const Alignable *ali);
00061 
00062   void fillCorrelations2D(float corr, const TransientTrackingRecHit::ConstRecHitPointer &hit);
00063   
00064   void fillPxbSurveyHistsChi2(const float &chi2);
00065   void fillPxbSurveyHistsLocalPars(const float &a0, const float &a1, const float &S, const float &phi);
00066 
00067  private:
00068   bool init(TDirectory *directory);
00069   bool equidistLogBins(double* bins, int nBins, double first, double last) const;
00070   void fillResidualHists(const std::vector<TH1*> &hists, float phiSensToNorm,
00071                          float residuum, float sigma);
00072   void fillResidualHitHists(const std::vector<TH1*> &hists, float angle,
00073                             float residuum, float sigma, unsigned int nHit);
00074   void fillTrack(const reco::Track *track, std::vector<TH1*> &trackHists1D,
00075                  std::vector<TH2*> &trackHists2D);
00076 
00077   template <class OBJECT_TYPE>  
00078     int GetIndex(const std::vector<OBJECT_TYPE*> &vec, const TString &name);
00079   template <class OBJECT_TYPE>  
00080   std::vector<OBJECT_TYPE*> cloneHists(const std::vector<OBJECT_TYPE*> &orgs,
00081                                        const TString &namAd, const TString &titAd) const;
00082   template <class OBJECT_TYPE>  
00083   void addToDirectory(const std::vector<OBJECT_TYPE*> &objs, TDirectory *dir) const;
00084 
00085   TDirectory *myRootDir;
00086   bool        myDeleteDir; 
00087 
00088   std::vector<TH1*> myTrackHists1D; // all input tracks 
00089   std::vector<TH2*> myTrackHists2D;
00090   std::vector<TH1*> myUsedTrackHists1D; // tracks used, i.e. tranferred to pede
00091   std::vector<TH2*> myUsedTrackHists2D;
00092   std::vector<TH1*> myTrajectoryHists1D;
00093   std::vector<TH2*> myTrajectoryHists2D;
00094   std::vector<TH2*> myDerivHists2D;
00095   std::vector<TH2*> myResidHists2D;
00096   std::vector<std::vector<TH1*> > myResidHistsVec1DX;
00097   std::vector<std::vector<TH1*> > myResidHistsVec1DY;
00098   std::vector<TH1*> myResidHitHists1DX;
00099   std::vector<TH1*> myResidHitHists1DY;
00100   std::vector<TH2*> myFrame2FrameHists2D;
00101   std::vector<TH1*> myCorrHists; // correlations
00102   std::vector<TH1*> myPxbSurveyHists; // correlations
00103 
00104 };
00105 
00106 template <class OBJECT_TYPE>  
00107 int MillePedeMonitor::GetIndex(const std::vector<OBJECT_TYPE*> &vec, const TString &name)
00108 {
00109   int result = 0;
00110   for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = vec.begin(), iterEnd = vec.end();
00111        iter != iterEnd; ++iter, ++result) {
00112     if (*iter && (*iter)->GetName() == name) return result;
00113   }
00114   edm::LogError("Alignment") << "@SUB=MillePedeMonitor::GetIndex" << " could not find " << name;
00115   return -1;
00116 }
00117 
00118 template <class OBJECT_TYPE>  
00119 std::vector<OBJECT_TYPE*> MillePedeMonitor::cloneHists(const std::vector<OBJECT_TYPE*> &orgs,
00120                                                        const TString &namAd,
00121                                                        const TString &titAd) const
00122 {
00123   // OBJECT_TYPE required to have methods Clone(const char*), GetName(), SetTitle(const char*) and GetTitle()
00124   std::vector<OBJECT_TYPE*> result;
00125   for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = orgs.begin(), iterEnd = orgs.end();
00126        iter != iterEnd; ++iter) {
00127     if (!(*iter)) continue;
00128     result.push_back(static_cast<OBJECT_TYPE*>((*iter)->Clone(namAd + (*iter)->GetName())));
00129     if (result.back()) result.back()->SetTitle((*iter)->GetTitle() + titAd);
00130     else edm::LogError("Alignment") <<"@SUB=MillePedeMonitor::cloneHists"
00131                                     << "out of memory?";
00132   }
00133   
00134   return result;
00135 }
00136 
00137 
00138 template <class OBJECT_TYPE>  
00139 void MillePedeMonitor::addToDirectory(const std::vector<OBJECT_TYPE*> &obs,
00140                                            TDirectory *dir) const
00141 {
00142   // OBJECT_TYPE is required to have method SetDirectory(TDirectory *dir)
00143   for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = obs.begin(), iterEnd = obs.end();
00144        iter != iterEnd; ++iter) {
00145     if (*iter) (*iter)->SetDirectory(dir);
00146   }
00147 }
00148 
00149 #endif