CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/Alignment/OfflineValidation/plugins/TrackerOfflineValidation.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    TrackerOfflineValidation
00004 // Class:      TrackerOfflineValidation
00005 // 
00013 //
00014 // Original Author:  Erik Butz
00015 //         Created:  Tue Dec 11 14:03:05 CET 2007
00016 // $Id: TrackerOfflineValidation.cc,v 1.44 2010/11/19 12:45:53 flucke Exp $
00017 //
00018 //
00019 
00020 // system include files
00021 #include <memory>
00022 #include <map>
00023 #include <sstream>
00024 #include <math.h>
00025 #include <utility>
00026 #include <vector>
00027 #include <iostream>
00028 
00029 // ROOT includes
00030 #include "TH1.h"
00031 #include "TH2.h"
00032 #include "TProfile.h"
00033 #include "TFile.h"
00034 #include "TTree.h"
00035 #include "TF1.h"
00036 #include "TMath.h"
00037 
00038 // user include files
00039 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00040 #include "FWCore/Framework/interface/Frameworkfwd.h"
00041 #include "FWCore/Framework/interface/EDAnalyzer.h"
00042 #include "FWCore/Framework/interface/EventSetup.h"
00043 #include "FWCore/Framework/interface/Event.h"
00044 #include "FWCore/Framework/interface/MakerMacros.h"
00045 #include "FWCore/Framework/interface/ESHandle.h"
00046 #include "FWCore/ServiceRegistry/interface/Service.h"
00047 
00048 #include "DataFormats/DetId/interface/DetId.h"
00049 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00050 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00051 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00052 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00053 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00054 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00055 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00056 #include "DataFormats/Math/interface/deltaPhi.h"
00057 
00058 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00059 #include "CommonTools/Utils/interface/TFileDirectory.h"
00060 
00061 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00062 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00063 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00064 
00065 #include "DQMServices/Core/interface/DQMStore.h"
00066 #include "DQMServices/Core/interface/MonitorElement.h"
00067 
00068 #include "Alignment/CommonAlignment/interface/Alignable.h"
00069 #include "Alignment/CommonAlignment/interface/Utilities.h"
00070 #include "Alignment/CommonAlignment/interface/AlignableObjectId.h"
00071 #include "Alignment/TrackerAlignment/interface/AlignableTracker.h"
00072 #include "Alignment/TrackerAlignment/interface/TrackerAlignableId.h"
00073 
00074 #include "Alignment/OfflineValidation/interface/TrackerValidationVariables.h"
00075 #include "Alignment/OfflineValidation/interface/TkOffTreeVariables.h"
00076 
00077 //
00078 // class declaration
00079 //
00080 class TrackerOfflineValidation : public edm::EDAnalyzer {
00081 public:
00082   explicit TrackerOfflineValidation(const edm::ParameterSet&);
00083   ~TrackerOfflineValidation();
00084   
00085   enum HistogrammType { XResidual, NormXResidual, 
00086                         YResidual, /*NormYResidual, */
00087                         XprimeResidual, NormXprimeResidual, 
00088                         YprimeResidual, NormYprimeResidual,
00089                         XResidualProfile, YResidualProfile };
00090   
00091 private:
00092 
00093   struct ModuleHistos{
00094     ModuleHistos() :  ResHisto(), NormResHisto(), ResYHisto(), /*NormResYHisto(),*/
00095                       ResXprimeHisto(), NormResXprimeHisto(), 
00096                       ResYprimeHisto(), NormResYprimeHisto() {} 
00097     TH1* ResHisto;
00098     TH1* NormResHisto;
00099     TH1* ResYHisto;
00100     /* TH1* NormResYHisto; */
00101     TH1* ResXprimeHisto;
00102     TH1* NormResXprimeHisto;
00103     TH1* ResYprimeHisto;
00104     TH1* NormResYprimeHisto;
00105 
00106     TH1* LocalX;
00107     TH1* LocalY;
00108     TProfile* ResXvsXProfile;
00109     TProfile* ResXvsYProfile;
00110     TProfile* ResYvsXProfile;
00111     TProfile* ResYvsYProfile;
00112   };
00113   
00114   // container struct to organize collection of histogramms during endJob
00115   struct SummaryContainer{
00116     SummaryContainer() : sumXResiduals_(), summaryXResiduals_(), 
00117                          sumNormXResiduals_(), summaryNormXResiduals_(),
00118                          sumYResiduals_(), summaryYResiduals_() ,
00119                          sumNormYResiduals_(), summaryNormYResiduals_() {}
00120     
00121     TH1* sumXResiduals_;
00122     TH1* summaryXResiduals_;
00123     TH1* sumNormXResiduals_;
00124     TH1* summaryNormXResiduals_;
00125     TH1* sumYResiduals_;
00126     TH1* summaryYResiduals_;
00127     TH1* sumNormYResiduals_;
00128     TH1* summaryNormYResiduals_;
00129   };
00130   
00131   
00132   struct DirectoryWrapper{
00133     DirectoryWrapper(const DirectoryWrapper& upDir,const std::string& newDir,
00134                      const std::string& basedir,bool useDqmMode)
00135       : tfd(0),
00136         dqmMode(useDqmMode),
00137         theDbe(0) {
00138       if (newDir.length()!=0){
00139         if(upDir.directoryString.length()!=0)directoryString=upDir.directoryString+"/"+newDir;
00140         else directoryString = newDir;
00141       }
00142       else
00143         directoryString=upDir.directoryString;
00144 
00145       if (!dqmMode){
00146         if (newDir.length()==0) tfd.reset(&(*upDir.tfd));
00147         else
00148           tfd.reset(new TFileDirectory(upDir.tfd->mkdir(newDir)));
00149       }
00150       else {
00151         theDbe=edm::Service<DQMStore>().operator->();
00152       }
00153     }
00154     
00155     DirectoryWrapper(const std::string& newDir,const std::string& basedir,bool useDqmMode)
00156       : tfd(0),
00157         dqmMode(useDqmMode),
00158         theDbe(0) {
00159       if (!dqmMode){
00160         edm::Service<TFileService> fs;
00161         if (newDir.length()==0){
00162           tfd.reset(new TFileDirectory(static_cast<TFileDirectory&>(*fs)));
00163         }
00164         else {
00165           tfd.reset(new TFileDirectory(fs->mkdir(newDir)));
00166           directoryString=newDir;
00167         }
00168       }
00169       else {
00170         if (newDir.length()!=0){
00171           if(basedir.length()!=0)directoryString=basedir+"/"+newDir;
00172           else directoryString = newDir;
00173         }
00174         else directoryString=basedir;
00175         theDbe=edm::Service<DQMStore>().operator->();
00176       }
00177     }
00178     // Generalization of Histogram Booking; allows switch between TFileService and DQMStore
00179     template <typename T> TH1* make(const char* name,const char* title,int nBinX,double minBinX,double maxBinX);
00180     template <typename T> TH1* make(const char* name,const char* title,int nBinX,double *xBins);//variable bin size in x for profile histo 
00181     template <typename T> TH1* make(const char* name,const char* title,int nBinX,double minBinX,double maxBinX,int nBinY,double minBinY,double maxBinY);
00182     template <typename T> TH1* make(const char* name,const char* title,int nBinX,double minBinX,double maxBinX,double minBinY,double maxBinY);  // at present not used
00183     
00184     std::auto_ptr<TFileDirectory> tfd;
00185     std::string directoryString;
00186     const bool dqmMode;
00187     DQMStore* theDbe;
00188   };
00189   
00190   
00191   // 
00192   // ------------- private member function -------------
00193   // 
00194   virtual void analyze(const edm::Event&, const edm::EventSetup&);
00195   virtual void endJob();
00196   
00197   virtual void checkBookHists(const edm::EventSetup& setup);
00198 
00199   void bookGlobalHists(DirectoryWrapper& tfd);
00200   void bookDirHists(DirectoryWrapper& tfd, const Alignable& ali, const AlignableObjectId& aliobjid);
00201   void bookHists(DirectoryWrapper& tfd, const Alignable& ali, align::StructureType type, int i, 
00202                  const AlignableObjectId& aliobjid);
00203  
00204   void collateSummaryHists( DirectoryWrapper& tfd, const Alignable& ali, int i, 
00205                             const AlignableObjectId& aliobjid, 
00206                             std::vector<TrackerOfflineValidation::SummaryContainer>& vLevelProfiles);
00207   
00208   void fillTree(TTree& tree, const std::map<int, TrackerOfflineValidation::ModuleHistos>& moduleHist_, 
00209                 TkOffTreeVariables& treeMem, const TrackerGeometry& tkgeom);
00210   
00211   TrackerOfflineValidation::SummaryContainer bookSummaryHists(DirectoryWrapper& tfd, 
00212                                                               const Alignable& ali, 
00213                                                               align::StructureType type, int i, 
00214                                                               const AlignableObjectId& aliobjid); 
00215 
00216   ModuleHistos& getHistStructFromMap(const DetId& detid); 
00217 
00218   bool isBarrel(uint32_t subDetId);
00219   bool isEndCap(uint32_t subDetId);
00220   bool isPixel(uint32_t subDetId);
00221   bool isDetOrDetUnit(align::StructureType type);
00222 
00223   TH1* bookTH1F(bool isTransient, DirectoryWrapper& tfd, const char* histName, const char* histTitle, 
00224                 int nBinsX, double lowX, double highX);
00225 
00226   TProfile* bookTProfile(bool isTransient, DirectoryWrapper& tfd, const char* histName, const char* histTitle, 
00227                          int nBinsX, double lowX, double highX);
00228 
00229   TProfile* bookTProfile(bool isTransient, DirectoryWrapper& tfd, const char* histName, const char* histTitle, 
00230                          int nBinsX, double lowX, double highX, double lowY, double highY);
00231 
00232   void getBinning(uint32_t subDetId, TrackerOfflineValidation::HistogrammType residualtype, 
00233                   int& nBinsX, double& lowerBoundX, double& upperBoundX);
00234 
00235   void summarizeBinInContainer(int bin, SummaryContainer& targetContainer, 
00236                                SummaryContainer& sourceContainer);
00237 
00238   void summarizeBinInContainer(int bin, uint32_t subDetId, SummaryContainer& targetContainer, 
00239                                ModuleHistos& sourceContainer);
00240 
00241   void setSummaryBin(int bin, TH1* targetHist, TH1* sourceHist);
00242     
00243   float Fwhm(const TH1* hist) const;
00244   std::pair<float,float> fitResiduals(TH1* hist) const; //, float meantmp, float rmstmp);
00245   float getMedian( const TH1* hist) const; 
00246   
00247   // From MillePedeAlignmentMonitor: Get Index for Arbitary vector<class> by name
00248   template <class OBJECT_TYPE> int GetIndex(const std::vector<OBJECT_TYPE*>& vec, const TString& name);
00249   
00250   
00251   // ---------- member data ---------------------------
00252 
00253   const edm::ParameterSet parSet_;
00254   edm::ESHandle<TrackerGeometry> tkGeom_;
00255   const TrackerGeometry *bareTkGeomPtr_; // ugly hack to book hists only once, but check 
00256 
00257   // parameters from cfg to steer
00258   const bool lCoorHistOn_;
00259   const bool moduleLevelHistsTransient_;
00260   const bool moduleLevelProfiles_;
00261   const bool overlappOn_;
00262   const bool stripYResiduals_;
00263   const bool useFwhm_;
00264   const bool useFit_;
00265   const bool useOverflowForRMS_;
00266   const bool dqmMode_;
00267   const std::string moduleDirectory_;
00268   
00269   std::map< std::pair<uint32_t, uint32_t>, TH1*> hOverlappResidual;
00270   
00271   // a vector to keep track which pointers should be deleted at the very end
00272   std::vector<TH1*> vDeleteObjects_;
00273 
00274   std::vector<TH1*> vTrackHistos_;
00275   std::vector<TH1*> vTrackProfiles_;
00276   std::vector<TH1*> vTrack2DHistos_;
00277   
00278   std::map<int,TrackerOfflineValidation::ModuleHistos> mPxbResiduals_;
00279   std::map<int,TrackerOfflineValidation::ModuleHistos> mPxeResiduals_;
00280   std::map<int,TrackerOfflineValidation::ModuleHistos> mTibResiduals_;
00281   std::map<int,TrackerOfflineValidation::ModuleHistos> mTidResiduals_;
00282   std::map<int,TrackerOfflineValidation::ModuleHistos> mTobResiduals_;
00283   std::map<int,TrackerOfflineValidation::ModuleHistos> mTecResiduals_;
00284 };
00285 
00286 
00287 //
00288 // constants, enums and typedefs
00289 //
00290 
00291 //
00292 // static data member definitions
00293 //
00294 
00295 template <class OBJECT_TYPE>  
00296 int TrackerOfflineValidation::GetIndex(const std::vector<OBJECT_TYPE*> &vec, const TString &name)
00297 {
00298   int result = 0;
00299   for (typename std::vector<OBJECT_TYPE*>::const_iterator iter = vec.begin(), iterEnd = vec.end();
00300        iter != iterEnd; ++iter, ++result) {
00301     if (*iter && (*iter)->GetName() == name) return result;
00302   }
00303   edm::LogError("Alignment") << "@SUB=TrackerOfflineValidation::GetIndex" << " could not find " << name;
00304   return -1;
00305 }
00306 
00307 
00308 template <> TH1* TrackerOfflineValidation::DirectoryWrapper::make<TH1F>(const char* name,const char* title,int nBinX,double minBinX,double maxBinX){
00309   if(dqmMode){theDbe->setCurrentFolder(directoryString); return theDbe->book1D(name,title,nBinX,minBinX,maxBinX)->getTH1();}
00310   else{return tfd->make<TH1F>(name,title,nBinX,minBinX,maxBinX);}
00311 }
00312 
00313 template <> TH1* TrackerOfflineValidation::DirectoryWrapper::make<TProfile>(const char* name,const char* title,int nBinX,double *xBins){
00314   if(dqmMode){
00315     theDbe->setCurrentFolder(directoryString);
00316     //DQM profile requires y-bins for construction... using TProfile creator by hand...
00317     TProfile *tmpProfile=new TProfile(name,title,nBinX,xBins);
00318     tmpProfile->SetDirectory(0);
00319     return theDbe->bookProfile(name,tmpProfile)->getTH1();
00320   }
00321   else{return tfd->make<TProfile>(name,title,nBinX,xBins);}
00322 }
00323 
00324 template <> TH1* TrackerOfflineValidation::DirectoryWrapper::make<TProfile>(const char* name,const char* title,int nBinX,double minBinX,double maxBinX){
00325   if(dqmMode){
00326     theDbe->setCurrentFolder(directoryString);
00327     //DQM profile requires y-bins for construction... using TProfile creator by hand...
00328     TProfile *tmpProfile=new TProfile(name,title,nBinX,minBinX,maxBinX);
00329     tmpProfile->SetDirectory(0);
00330     return theDbe->bookProfile(name,tmpProfile)->getTH1();
00331   }
00332   else{return tfd->make<TProfile>(name,title,nBinX,minBinX,maxBinX);}
00333 }
00334 
00335 template <> TH1* TrackerOfflineValidation::DirectoryWrapper::make<TProfile>(const char* name ,const char* title,int nbinX,double minX ,double maxX,double minY,double maxY){
00336   if(dqmMode){
00337     theDbe->setCurrentFolder(directoryString);
00338     int dummy(0); // DQMProfile wants Y channels... does not use them!
00339     return (theDbe->bookProfile(name,title,nbinX,minX,maxX,dummy,minY,maxY)->getTH1());
00340   }
00341   else{
00342     return tfd->make<TProfile>(name,title,nbinX,minX,maxX,minY,maxY);
00343   }
00344 }
00345 
00346 template <> TH1* TrackerOfflineValidation::DirectoryWrapper::make<TH2F>(const char* name,const char* title,int nBinX,double minBinX,double maxBinX,int nBinY,double minBinY,double maxBinY){
00347   if(dqmMode){theDbe->setCurrentFolder(directoryString); return theDbe->book2D(name,title,nBinX,minBinX,maxBinX,nBinY,minBinY,maxBinY)->getTH1();}
00348   else{return tfd->make<TH2F>(name,title,nBinX,minBinX,maxBinX,nBinY,minBinY,maxBinY);}
00349 }
00350 
00351 
00352 //
00353 // constructors and destructor
00354 //
00355 TrackerOfflineValidation::TrackerOfflineValidation(const edm::ParameterSet& iConfig)
00356   : parSet_(iConfig), bareTkGeomPtr_(0), lCoorHistOn_(parSet_.getParameter<bool>("localCoorHistosOn")),
00357     moduleLevelHistsTransient_(parSet_.getParameter<bool>("moduleLevelHistsTransient")),
00358     moduleLevelProfiles_(parSet_.getParameter<bool>("moduleLevelProfiles")),
00359     overlappOn_(parSet_.getParameter<bool>("overlappOn")), 
00360     stripYResiduals_(parSet_.getParameter<bool>("stripYResiduals")), 
00361     useFwhm_(parSet_.getParameter<bool>("useFwhm")),
00362     useFit_(parSet_.getParameter<bool>("useFit")),
00363     useOverflowForRMS_(parSet_.getParameter<bool>("useOverflowForRMS")),
00364     dqmMode_(parSet_.getParameter<bool>("useInDqmMode")),
00365     moduleDirectory_(parSet_.getParameter<std::string>("moduleDirectoryInOutput"))
00366 {
00367 }
00368 
00369 
00370 TrackerOfflineValidation::~TrackerOfflineValidation()
00371 {
00372    // do anything here that needs to be done at desctruction time
00373    // (e.g. close files, deallocate resources etc.)
00374   for( std::vector<TH1*>::const_iterator it = vDeleteObjects_.begin(), itEnd = vDeleteObjects_.end(); 
00375        it != itEnd;
00376        ++it) delete *it;
00377 }
00378 
00379 
00380 //
00381 // member functions
00382 //
00383 
00384 
00385 // ------------ method called once each job just before starting event loop  ------------
00386 void
00387 TrackerOfflineValidation::checkBookHists(const edm::EventSetup& es)
00388 {
00389   es.get<TrackerDigiGeometryRecord>().get( tkGeom_ );
00390   const TrackerGeometry *newBareTkGeomPtr = &(*tkGeom_);
00391   if (newBareTkGeomPtr == bareTkGeomPtr_) return; // already booked hists, nothing changed
00392 
00393   if (!bareTkGeomPtr_) { // pointer not yet set: called the first time => book hists
00394     AlignableObjectId aliobjid;
00395     
00396     // construct alignable tracker to get access to alignable hierarchy 
00397     AlignableTracker aliTracker(&(*tkGeom_));
00398     
00399     edm::LogInfo("TrackerOfflineValidation") << "There are " << newBareTkGeomPtr->detIds().size()
00400                                              << " dets in the Geometry record.\n"
00401                                              << "Out of these "<<newBareTkGeomPtr->detUnitIds().size()
00402                                              <<" are detUnits";
00403     
00404     // Book Histogramms for global track quantities
00405     std::string globDir("GlobalTrackVariables");
00406     DirectoryWrapper trackglobal(globDir,moduleDirectory_,dqmMode_);
00407     this->bookGlobalHists(trackglobal);
00408     
00409     // recursively book histogramms on lowest level
00410     DirectoryWrapper tfdw("",moduleDirectory_,dqmMode_);
00411     this->bookDirHists(tfdw, aliTracker, aliobjid);
00412   }
00413   else { // histograms booked, but changed TrackerGeometry?
00414     edm::LogWarning("GeometryChange") << "@SUB=checkBookHists"
00415                                       << "TrackerGeometry changed, but will not re-book hists!";
00416   }
00417   bareTkGeomPtr_ = newBareTkGeomPtr;
00418 }
00419 
00420 
00421 void 
00422 TrackerOfflineValidation::bookGlobalHists(DirectoryWrapper& tfd )
00423 {
00424 
00425   vTrackHistos_.push_back(tfd.make<TH1F>("h_tracketa",
00426                                          "Track #eta;#eta_{Track};Number of Tracks",
00427                                          90,-3.,3.));
00428   vTrackHistos_.push_back(tfd.make<TH1F>("h_trackphi",
00429                                          "Track #phi;#phi_{Track};Number of Tracks",
00430                                          90,-3.15,3.15));
00431   vTrackHistos_.push_back(tfd.make<TH1F>("h_trackNumberOfValidHits",
00432                                          "Track # of valid hits;# of valid hits _{Track};Number of Tracks",
00433                                         40,0.,40.));
00434   vTrackHistos_.push_back(tfd.make<TH1F>("h_trackNumberOfLostHits",
00435                                          "Track # of lost hits;# of lost hits _{Track};Number of Tracks",
00436                                         10,0.,10.));
00437   vTrackHistos_.push_back(tfd.make<TH1F>("h_curvature",
00438                                          "Curvature #kappa;#kappa_{Track};Number of Tracks",
00439                                          100,-.05,.05));
00440   vTrackHistos_.push_back(tfd.make<TH1F>("h_curvature_pos",
00441                                          "Curvature |#kappa| Positive Tracks;|#kappa_{pos Track}|;Number of Tracks",
00442                                          100,.0,.05));
00443   vTrackHistos_.push_back(tfd.make<TH1F>("h_curvature_neg",
00444                                          "Curvature |#kappa| Negative Tracks;|#kappa_{neg Track}|;Number of Tracks",
00445                                          100,.0,.05));
00446   vTrackHistos_.push_back(tfd.make<TH1F>("h_diff_curvature",
00447                                          "Curvature |#kappa| Tracks Difference;|#kappa_{Track}|;# Pos Tracks - # Neg Tracks",
00448                                          100,.0,.05));
00449   vTrackHistos_.push_back(tfd.make<TH1F>("h_chi2",
00450                                          "#chi^{2};#chi^{2}_{Track};Number of Tracks",
00451                                          500,-0.01,500.));
00452   vTrackHistos_.push_back(tfd.make<TH1F>("h_chi2Prob",
00453                                          "#chi^{2} probability;#chi^{2}prob_{Track};Number of Tracks",
00454                                         100,0.0,1.));          
00455   vTrackHistos_.push_back(tfd.make<TH1F>("h_normchi2",
00456                                          "#chi^{2}/ndof;#chi^{2}/ndof;Number of Tracks",
00457                                          100,-0.01,10.));     
00458   vTrackHistos_.push_back(tfd.make<TH1F>("h_pt",
00459                                          "p_{T}^{track};p_{T}^{track} [GeV];Number of Tracks",
00460                                          250,0.,250));           
00461   vTrackHistos_.push_back(tfd.make<TH1F>("h_ptResolution",
00462                                          "#delta{p_{T}/p_{T}^{track}};#delta_{p_{T}/p_{T}^{track}};Number of Tracks",
00463                                          100,0.,0.5));           
00464 
00465   vTrackProfiles_.push_back(tfd.make<TProfile>("p_d0_vs_phi",
00466                                                "Transverse Impact Parameter vs. #phi;#phi_{Track};#LT d_{0} #GT [cm]",
00467                                                100,-3.15,3.15));
00468   vTrackProfiles_.push_back(tfd.make<TProfile>("p_dz_vs_phi",
00469                                                "Longitudinal Impact Parameter vs. #phi;#phi_{Track};#LT d_{z} #GT [cm]",
00470                                                100,-3.15,3.15));
00471   vTrackProfiles_.push_back(tfd.make<TProfile>("p_d0_vs_eta",
00472                                                "Transverse Impact Parameter vs. #eta;#eta_{Track};#LT d_{0} #GT [cm]",
00473                                                100,-3.15,3.15));
00474   vTrackProfiles_.push_back(tfd.make<TProfile>("p_dz_vs_eta",
00475                                                "Longitudinal Impact Parameter vs. #eta;#eta_{Track};#LT d_{z} #GT [cm]",
00476                                                100,-3.15,3.15));
00477   vTrackProfiles_.push_back(tfd.make<TProfile>("p_chi2_vs_phi",
00478                                                "#chi^{2} vs. #phi;#phi_{Track};#LT #chi^{2} #GT",
00479                                                100,-3.15,3.15));
00480   vTrackProfiles_.push_back(tfd.make<TProfile>("p_chi2Prob_vs_phi",
00481                                                "#chi^{2} probablility vs. #phi;#phi_{Track};#LT #chi^{2} probability#GT",
00482                                                100,-3.15,3.15));
00483   vTrackProfiles_.push_back(tfd.make<TProfile>("p_normchi2_vs_phi",
00484                                                "#chi^{2}/ndof vs. #phi;#phi_{Track};#LT #chi^{2}/ndof #GT",
00485                                                100,-3.15,3.15));
00486   vTrackProfiles_.push_back(tfd.make<TProfile>("p_chi2_vs_eta",
00487                                                "#chi^{2} vs. #eta;#eta_{Track};#LT #chi^{2} #GT",
00488                                                100,-3.15,3.15));
00489   //variable binning for chi2/ndof vs. pT
00490   double xBins[19]={0.,0.15,0.5,1.,1.5,2.,2.5,3.,3.5,4.,4.5,5.,7.,10.,15.,25.,40.,100.,200.};
00491   vTrackProfiles_.push_back(tfd.make<TProfile>("p_normchi2_vs_pt",
00492                                                "norm #chi^{2} vs. p_{T}_{Track}; p_{T}_{Track};#LT #chi^{2}/ndof #GT",
00493                                                18,xBins));
00494 
00495   vTrackProfiles_.push_back(tfd.make<TProfile>("p_normchi2_vs_p",
00496                                                "#chi^{2}/ndof vs. p_{Track};p_{Track};#LT #chi^{2}/ndof #GT",
00497                                                 18,xBins));
00498   vTrackProfiles_.push_back(tfd.make<TProfile>("p_chi2Prob_vs_eta",
00499                                                "#chi^{2} probability vs. #eta;#eta_{Track};#LT #chi^{2} probability #GT",
00500                                                100,-3.15,3.15));   
00501   vTrackProfiles_.push_back(tfd.make<TProfile>("p_normchi2_vs_eta",
00502                                                "#chi^{2}/ndof vs. #eta;#eta_{Track};#LT #chi^{2}/ndof #GT",
00503                                                100,-3.15,3.15));
00504   vTrackProfiles_.push_back(tfd.make<TProfile>("p_kappa_vs_phi",
00505                                                "#kappa vs. #phi;#phi_{Track};#kappa",
00506                                                100,-3.15,3.15));
00507   vTrackProfiles_.push_back(tfd.make<TProfile>("p_kappa_vs_eta",
00508                                                "#kappa vs. #eta;#eta_{Track};#kappa",
00509                                                100,-3.15,3.15));
00510   vTrackProfiles_.push_back(tfd.make<TProfile>("p_ptResolution_vs_phi",
00511                                                "#delta_{p_{T}}/p_{T}^{track};#phi^{track};#delta_{p_{T}}/p_{T}^{track}",
00512                                                100, -3.15,3.15));
00513   vTrackProfiles_.push_back(tfd.make<TProfile>("p_ptResolution_vs_eta",
00514                                                "#delta_{p_{T}}/p_{T}^{track};#eta^{track};#delta_{p_{T}}/p_{T}^{track}",
00515                                                100, -3.15,3.15));
00516 
00517   vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_d0_vs_phi",
00518                                            "Transverse Impact Parameter vs. #phi;#phi_{Track};d_{0} [cm]",
00519                                            100, -3.15, 3.15, 100,-1.,1.) );
00520   vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_dz_vs_phi",
00521                                            "Longitudinal Impact Parameter vs. #phi;#phi_{Track};d_{z} [cm]",
00522                                            100, -3.15, 3.15, 100,-100.,100.));
00523   vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_d0_vs_eta",
00524                                            "Transverse Impact Parameter vs. #eta;#eta_{Track};d_{0} [cm]",
00525                                            100, -3.15, 3.15, 100,-1.,1.));
00526   vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_dz_vs_eta",
00527                                            "Longitudinal Impact Parameter vs. #eta;#eta_{Track};d_{z} [cm]",
00528                                            100, -3.15, 3.15, 100,-100.,100.));
00529   vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_chi2_vs_phi",
00530                                            "#chi^{2} vs. #phi;#phi_{Track};#chi^{2}",
00531                                            100, -3.15, 3.15, 500, 0., 500.));
00532   vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_chi2Prob_vs_phi",
00533                                            "#chi^{2} probability vs. #phi;#phi_{Track};#chi^{2} probability",
00534                                           100, -3.15, 3.15, 100, 0., 1.));
00535   vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_normchi2_vs_phi",
00536                                            "#chi^{2}/ndof vs. #phi;#phi_{Track};#chi^{2}/ndof",
00537                                            100, -3.15, 3.15, 100, 0., 10.));
00538   vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_chi2_vs_eta",
00539                                            "#chi^{2} vs. #eta;#eta_{Track};#chi^{2}",
00540                                            100, -3.15, 3.15, 500, 0., 500.));  
00541   vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_chi2Prob_vs_eta",
00542                                            "#chi^{2} probaility vs. #eta;#eta_{Track};#chi^{2} probability",
00543                                            100, -3.15, 3.15, 100, 0., 1.));  
00544   vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_normchi2_vs_eta",
00545                                            "#chi^{2}/ndof vs. #eta;#eta_{Track};#chi^{2}/ndof",
00546                                            100,-3.15,3.15, 100, 0., 10.));
00547   vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_kappa_vs_phi",
00548                                            "#kappa vs. #phi;#phi_{Track};#kappa",
00549                                            100,-3.15,3.15, 100, .0,.05));
00550   vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_kappa_vs_eta",
00551                                            "#kappa vs. #eta;#eta_{Track};#kappa",
00552                                            100,-3.15,3.15, 100, .0,.05));
00553   vTrack2DHistos_.push_back(tfd.make<TH2F>("h2_normchi2_vs_kappa",
00554                                            "#kappa vs. #chi^{2}/ndof;#chi^{2}/ndof;#kappa",
00555                                            100,0.,10, 100,-.03,.03));
00556 
00557   /****************** Definition of 2-D Histos of ResX vs momenta ****************************/
00558   vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_pixB",
00559                                            "#momentum vs. #resX in pixB;#momentum;#resX",
00560                                            15,0.,15., 200, -0.1,0.1));                     
00561   vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_pixE",
00562                                            "#momentum vs. #resX in pixE;#momentum;#resX",
00563                                            15,0.,15., 200, -0.1,0.1)); 
00564   vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_TIB",
00565                                            "#momentum vs. #resX in TIB;#momentum;#resX",
00566                                            15,0.,15., 200, -0.1,0.1)); 
00567   vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_TID",
00568                                            "#momentum vs. #resX in TID;#momentum;#resX",
00569                                            15,0.,15., 200, -0.1,0.1)); 
00570   vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_TOB",
00571                                            "#momentum vs. #resX in TOB;#momentum;#resX",
00572                                            15,0.,15., 200, -0.1,0.1));
00573   vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resXprime_TEC",
00574                                            "#momentum vs. #resX in TEC;#momentum;#resX",
00575                                            15,0.,15., 200, -0.1,0.1));  
00576 
00577   /****************** Definition of 2-D Histos of ResY vs momenta ****************************/
00578   vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resYprime_pixB",
00579                                            "#momentum vs. #resY in pixB;#momentum;#resY",
00580                                            15,0.,15., 200, -0.1,0.1));                     
00581   vTrack2DHistos_.push_back(tfd.make<TH2F>("p_vs_resYprime_pixE",
00582                                            "#momentum vs. #resY in pixE;#momentum;#resY",
00583                                            15,0.,15., 200, -0.1,0.1)); 
00584 
00585 }
00586 
00587 
00588 void
00589 TrackerOfflineValidation::bookDirHists(DirectoryWrapper& tfd, const Alignable& ali, const AlignableObjectId& aliobjid)
00590 {
00591   std::vector<Alignable*> alivec(ali.components());
00592   for(int i=0, iEnd = ali.components().size();i < iEnd; ++i) {
00593     std::string structurename  = aliobjid.typeToName((alivec)[i]->alignableObjectId());
00594     LogDebug("TrackerOfflineValidation") << "StructureName = " << structurename;
00595     std::stringstream dirname;
00596     dirname << structurename;
00597     // add no suffix counter to Strip and Pixel, just aesthetics
00598     if (structurename != "Strip" && structurename != "Pixel") dirname << "_" << i+1;
00599 
00600     if (structurename.find("Endcap",0) != std::string::npos ) {
00601       DirectoryWrapper f(tfd,dirname.str(),moduleDirectory_,dqmMode_);
00602       bookHists(f, *(alivec)[i], ali.alignableObjectId() , i, aliobjid);
00603       bookDirHists( f, *(alivec)[i], aliobjid);
00604     } else if( !(this->isDetOrDetUnit( (alivec)[i]->alignableObjectId()) )
00605               || alivec[i]->components().size() > 1) {      
00606       DirectoryWrapper f(tfd,dirname.str(),moduleDirectory_,dqmMode_);
00607       bookHists(tfd, *(alivec)[i], ali.alignableObjectId() , i, aliobjid);
00608       bookDirHists( f, *(alivec)[i], aliobjid);
00609     } else {
00610       bookHists(tfd, *(alivec)[i], ali.alignableObjectId() , i, aliobjid);
00611     }
00612   }
00613 }
00614 
00615 
00616 void 
00617 TrackerOfflineValidation::bookHists(DirectoryWrapper& tfd, const Alignable& ali, align::StructureType type, int i, const AlignableObjectId& aliobjid)
00618 {
00619   TrackerAlignableId aliid;
00620   const DetId id = ali.id();
00621 
00622   // comparing subdetandlayer to subdetIds gives a warning at compile time
00623   // -> subdetandlayer could also be pair<uint,uint> but this has to be adapted
00624   // in AlignableObjId 
00625   std::pair<int,int> subdetandlayer = aliid.typeAndLayerFromDetId(id);
00626 
00627   align::StructureType subtype = align::invalid;
00628   
00629   // are we on or just above det, detunit level respectively?
00630   if (type == align::AlignableDetUnit )subtype = type;
00631   else if( this->isDetOrDetUnit(ali.alignableObjectId()) ) subtype = ali.alignableObjectId();
00632   
00633   // construct histogramm title and name
00634   std::stringstream histoname, histotitle, normhistoname, normhistotitle, 
00635     yhistoname, yhistotitle,
00636     xprimehistoname, xprimehistotitle, normxprimehistoname, normxprimehistotitle,
00637     yprimehistoname, yprimehistotitle, normyprimehistoname, normyprimehistotitle,
00638     localxname, localxtitle, localyname, localytitle,
00639     resxvsxprofilename, resxvsxprofiletitle, resyvsxprofilename, resyvsxprofiletitle,
00640     resxvsyprofilename, resxvsyprofiletitle, resyvsyprofilename, resyvsyprofiletitle; 
00641   
00642   std::string wheel_or_layer;
00643 
00644   if( this->isEndCap(static_cast<uint32_t>(subdetandlayer.first)) ) wheel_or_layer = "_wheel_";
00645   else if ( this->isBarrel(static_cast<uint32_t>(subdetandlayer.first)) ) wheel_or_layer = "_layer_";
00646   else edm::LogWarning("TrackerOfflineValidation") << "@SUB=TrackerOfflineValidation::bookHists" 
00647                                                    << "Unknown subdetid: " <<  subdetandlayer.first;     
00648   
00649   histoname << "h_residuals_subdet_" << subdetandlayer.first 
00650             << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00651   yhistoname << "h_y_residuals_subdet_" << subdetandlayer.first 
00652              << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00653   xprimehistoname << "h_xprime_residuals_subdet_" << subdetandlayer.first 
00654                   << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00655   yprimehistoname << "h_yprime_residuals_subdet_" << subdetandlayer.first 
00656                   << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00657   normhistoname << "h_normresiduals_subdet_" << subdetandlayer.first 
00658                 << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00659   normxprimehistoname << "h_normxprimeresiduals_subdet_" << subdetandlayer.first 
00660                       << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00661   normyprimehistoname << "h_normyprimeresiduals_subdet_" << subdetandlayer.first 
00662                       << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00663   histotitle << "X Residual for module " << id.rawId() << ";x_{tr} - x_{hit} [cm]";
00664   yhistotitle << "Y Residual for module " << id.rawId() << ";y_{tr} - y_{hit} [cm]";
00665   normhistotitle << "Normalized Residual for module " << id.rawId() << ";x_{tr} - x_{hit}/#sigma";
00666   xprimehistotitle << "X' Residual for module " << id.rawId() << ";(x_{tr} - x_{hit})' [cm]";
00667   normxprimehistotitle << "Normalized X' Residual for module " << id.rawId() << ";(x_{tr} - x_{hit})'/#sigma";
00668   yprimehistotitle << "Y' Residual for module " << id.rawId() << ";(y_{tr} - y_{hit})' [cm]";
00669   normyprimehistotitle << "Normalized Y' Residual for module " << id.rawId() << ";(y_{tr} - y_{hit})'/#sigma";
00670   
00671   if ( moduleLevelProfiles_ ) {
00672     localxname << "h_localx_subdet_" << subdetandlayer.first 
00673                << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00674     localyname << "h_localy_subdet_" << subdetandlayer.first 
00675                << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00676     localxtitle << "u local for module " << id.rawId() << "; u_{tr,r}";
00677     localytitle << "v local for module " << id.rawId() << "; v_{tr,r}";
00678 
00679     resxvsxprofilename << "p_residuals_x_vs_x_subdet_" << subdetandlayer.first 
00680                        << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00681     resyvsxprofilename << "p_residuals_y_vs_x_subdet_" << subdetandlayer.first 
00682                        << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00683     resxvsyprofilename << "p_residuals_x_vs_y_subdet_" << subdetandlayer.first 
00684                        << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00685     resyvsyprofilename << "p_residuals_y_vs_y_subdet_" << subdetandlayer.first 
00686                        << wheel_or_layer << subdetandlayer.second << "_module_" << id.rawId();
00687     resxvsxprofiletitle << "U Residual vs u for module " << id.rawId() << "; u_{tr,r} ;(u_{tr} - u_{hit})/tan#alpha [cm]";
00688     resyvsxprofiletitle << "V Residual vs u for module " << id.rawId() << "; u_{tr,r} ;(v_{tr} - v_{hit})/tan#beta  [cm]";
00689     resxvsyprofiletitle << "U Residual vs v for module " << id.rawId() << "; v_{tr,r} ;(u_{tr} - u_{hit})/tan#alpha [cm]";
00690     resyvsyprofiletitle << "V Residual vs v for module " << id.rawId() << "; v_{tr,r} ;(v_{tr} - v_{hit})/tan#beta  [cm]";
00691   }
00692   
00693   if( this->isDetOrDetUnit( subtype ) ) {
00694     ModuleHistos &histStruct = this->getHistStructFromMap(id);
00695     int nbins = 0;
00696     double xmin = 0., xmax = 0.;
00697     double ymin = -0.1, ymax = 0.1;
00698 
00699     // do not allow transient hists in DQM mode
00700     bool moduleLevelHistsTransient(moduleLevelHistsTransient_);
00701     if (dqmMode_) moduleLevelHistsTransient = false;
00702     
00703     // decide via cfg if hists in local coordinates should be booked 
00704     if(lCoorHistOn_) {
00705       this->getBinning(id.subdetId(), XResidual, nbins, xmin, xmax);
00706       histStruct.ResHisto = this->bookTH1F(moduleLevelHistsTransient, tfd, 
00707                                            histoname.str().c_str(),histotitle.str().c_str(),                 
00708                                            nbins, xmin, xmax);
00709       this->getBinning(id.subdetId(), NormXResidual, nbins, xmin, xmax);
00710       histStruct.NormResHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
00711                                                normhistoname.str().c_str(),normhistotitle.str().c_str(),
00712                                                nbins, xmin, xmax);
00713     } 
00714     this->getBinning(id.subdetId(), XprimeResidual, nbins, xmin, xmax);
00715     histStruct.ResXprimeHisto = this->bookTH1F(moduleLevelHistsTransient, tfd, 
00716                                                xprimehistoname.str().c_str(),xprimehistotitle.str().c_str(),
00717                                                nbins, xmin, xmax);
00718     this->getBinning(id.subdetId(), NormXprimeResidual, nbins, xmin, xmax);
00719     histStruct.NormResXprimeHisto = this->bookTH1F(moduleLevelHistsTransient, tfd, 
00720                                                    normxprimehistoname.str().c_str(),normxprimehistotitle.str().c_str(),
00721                                                    nbins, xmin, xmax);
00722 
00723     if ( moduleLevelProfiles_ ) {
00724       this->getBinning(id.subdetId(), XResidualProfile, nbins, xmin, xmax);      
00725 
00726       histStruct.LocalX = this->bookTH1F(moduleLevelHistsTransient, tfd,
00727                                          localxname.str().c_str(),localxtitle.str().c_str(),
00728                                          nbins, xmin, xmax);
00729       histStruct.LocalY = this->bookTH1F(moduleLevelHistsTransient, tfd,
00730                                          localyname.str().c_str(),localytitle.str().c_str(),
00731                                          nbins, xmin, xmax);
00732       histStruct.ResXvsXProfile = this->bookTProfile(moduleLevelHistsTransient, tfd,
00733                                                      resxvsxprofilename.str().c_str(),resxvsxprofiletitle.str().c_str(),
00734                                                      nbins, xmin, xmax, ymin, ymax);
00735       histStruct.ResXvsYProfile = this->bookTProfile(moduleLevelHistsTransient, tfd,
00736                                                      resxvsyprofilename.str().c_str(),resxvsyprofiletitle.str().c_str(),
00737                                                      nbins, xmin, xmax, ymin, ymax);
00738     }
00739 
00740     if( this->isPixel(subdetandlayer.first) || stripYResiduals_ ) {
00741       this->getBinning(id.subdetId(), YprimeResidual, nbins, xmin, xmax);
00742       histStruct.ResYprimeHisto = this->bookTH1F(moduleLevelHistsTransient, tfd,
00743                                                  yprimehistoname.str().c_str(),yprimehistotitle.str().c_str(),
00744                                                  nbins, xmin, xmax);
00745       if (lCoorHistOn_) { // un-primed y-residual
00746         this->getBinning(id.subdetId(), YResidual, nbins, xmin, xmax);
00747         histStruct.ResYHisto = this->bookTH1F(moduleLevelHistsTransient, tfd, 
00748                                               yhistoname.str().c_str(), yhistotitle.str().c_str(),
00749                                               nbins, xmin, xmax);
00750       }
00751       this->getBinning(id.subdetId(), NormYprimeResidual, nbins, xmin, xmax);
00752       histStruct.NormResYprimeHisto = this->bookTH1F(moduleLevelHistsTransient, tfd, 
00753                                                      normyprimehistoname.str().c_str(),normyprimehistotitle.str().c_str(),
00754                                                      nbins, xmin, xmax);
00755       // Here we could add un-primed normalised y-residuals if(lCoorHistOn_)...
00756       if ( moduleLevelProfiles_ ) {
00757         this->getBinning(id.subdetId(), YResidualProfile, nbins, xmin, xmax);      
00758         
00759         histStruct.ResYvsXProfile = this->bookTProfile(moduleLevelHistsTransient, tfd,
00760                                                        resyvsxprofilename.str().c_str(),resyvsxprofiletitle.str().c_str(),
00761                                                        nbins, xmin, xmax, ymin, ymax);
00762         histStruct.ResYvsYProfile = this->bookTProfile(moduleLevelHistsTransient, tfd,
00763                                                        resyvsyprofilename.str().c_str(),resyvsyprofiletitle.str().c_str(),
00764                                                        nbins, xmin, xmax, ymin, ymax);
00765       }
00766     }
00767   }
00768 }
00769 
00770 
00771 TH1* TrackerOfflineValidation::bookTH1F(bool isTransient, DirectoryWrapper& tfd, const char* histName, const char* histTitle, 
00772                                         int nBinsX, double lowX, double highX)
00773 {
00774   if (isTransient) {
00775     vDeleteObjects_.push_back(new TH1F(histName, histTitle, nBinsX, lowX, highX));
00776     return vDeleteObjects_.back(); // return last element of vector
00777   } 
00778   else
00779     return tfd.make<TH1F>(histName, histTitle, nBinsX, lowX, highX);
00780 }
00781 
00782 TProfile* TrackerOfflineValidation::bookTProfile(bool isTransient, DirectoryWrapper& tfd, const char* histName, const char* histTitle, 
00783                                                  int nBinsX, double lowX, double highX)
00784 {
00785   if (isTransient) {
00786     TProfile * profile = new TProfile(histName, histTitle, nBinsX, lowX, highX);
00787     vDeleteObjects_.push_back(profile);
00788     return profile;
00789   }
00790   else
00791     return (TProfile*)tfd.make<TProfile>(histName, histTitle, nBinsX, lowX, highX);
00792 }
00793 
00794 
00795 TProfile* TrackerOfflineValidation::bookTProfile(bool isTransient, DirectoryWrapper& tfd, const char* histName, const char* histTitle, 
00796                                                  int nBinsX, double lowX, double highX, double lowY, double highY)
00797 {
00798   if (isTransient) {
00799     TProfile * profile = new TProfile(histName, histTitle, nBinsX, lowX, highX, lowY, highY);
00800     vDeleteObjects_.push_back(profile);
00801     return profile;
00802   }
00803   else
00804     return (TProfile*)tfd.make<TProfile>(histName, histTitle, nBinsX, lowX, highX, lowY, highY);
00805 }
00806 
00807 bool TrackerOfflineValidation::isBarrel(uint32_t subDetId)
00808 {
00809   return (subDetId == StripSubdetector::TIB ||
00810           subDetId == StripSubdetector::TOB ||
00811           subDetId == PixelSubdetector::PixelBarrel );
00812 }
00813 
00814 
00815 bool TrackerOfflineValidation::isEndCap(uint32_t subDetId)
00816 {
00817   return ( subDetId == StripSubdetector::TID ||
00818            subDetId == StripSubdetector::TEC ||
00819            subDetId == PixelSubdetector::PixelEndcap );
00820 }
00821 
00822 
00823 bool TrackerOfflineValidation::isPixel(uint32_t subDetId)
00824 {
00825   return (subDetId == PixelSubdetector::PixelBarrel ||
00826           subDetId == PixelSubdetector::PixelEndcap );
00827 }
00828 
00829 
00830 bool TrackerOfflineValidation::isDetOrDetUnit(align::StructureType type)
00831 {
00832   return ( type == align::AlignableDet ||
00833            type == align::AlignableDetUnit );
00834 }
00835 
00836 
00837 void 
00838 TrackerOfflineValidation::getBinning(uint32_t subDetId, 
00839                                      TrackerOfflineValidation::HistogrammType residualType, 
00840                                      int& nBinsX, double& lowerBoundX, double& upperBoundX)
00841 {
00842   // determine if 
00843   const bool isPixel = this->isPixel(subDetId);
00844   
00845   edm::ParameterSet binningPSet;
00846   
00847   switch(residualType) 
00848     {
00849     case XResidual :
00850       if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1XResPixelModules");                
00851       else binningPSet        = parSet_.getParameter<edm::ParameterSet>("TH1XResStripModules");                
00852       break;
00853     case NormXResidual : 
00854       if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormXResPixelModules");             
00855       else binningPSet        = parSet_.getParameter<edm::ParameterSet>("TH1NormXResStripModules");                
00856       break;
00857     case XprimeResidual :
00858       if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1XprimeResPixelModules");                
00859       else binningPSet        = parSet_.getParameter<edm::ParameterSet>("TH1XprimeResStripModules");                
00860       break;
00861     case NormXprimeResidual :
00862       if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormXprimeResPixelModules");
00863       else binningPSet        = parSet_.getParameter<edm::ParameterSet>("TH1NormXprimeResStripModules");
00864       break;
00865     case YResidual : // borrow y-residual binning from yprime
00866     case YprimeResidual :
00867       if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1YResPixelModules");                
00868       else binningPSet        = parSet_.getParameter<edm::ParameterSet>("TH1YResStripModules");                
00869       break; 
00870       /* case NormYResidual :*/
00871     case NormYprimeResidual :
00872       if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TH1NormYResPixelModules");             
00873       else binningPSet        = parSet_.getParameter<edm::ParameterSet>("TH1NormYResStripModules");  
00874       break;
00875     case XResidualProfile :
00876       if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TProfileXResPixelModules");                
00877       else binningPSet        = parSet_.getParameter<edm::ParameterSet>("TProfileXResStripModules");                
00878       break;
00879     case YResidualProfile :
00880       if(isPixel) binningPSet = parSet_.getParameter<edm::ParameterSet>("TProfileYResPixelModules");                
00881       else binningPSet        = parSet_.getParameter<edm::ParameterSet>("TProfileYResStripModules");                
00882       break;
00883     }
00884   nBinsX      = binningPSet.getParameter<int32_t>("Nbinx");                    
00885   lowerBoundX = binningPSet.getParameter<double>("xmin");                      
00886   upperBoundX = binningPSet.getParameter<double>("xmax");     
00887 }
00888 
00889 
00890 void 
00891 TrackerOfflineValidation::setSummaryBin(int bin, TH1* targetHist, TH1* sourceHist)
00892 {
00893   if(targetHist && sourceHist) {
00894     targetHist->SetBinContent(bin, sourceHist->GetMean(1));
00895     if(useFwhm_) targetHist->SetBinError(bin, Fwhm(sourceHist)/2.);
00896     else targetHist->SetBinError(bin, sourceHist->GetRMS(1) );
00897   }
00898   else return;
00899 }
00900 
00901 
00902 void 
00903 TrackerOfflineValidation::summarizeBinInContainer( int bin, SummaryContainer& targetContainer, 
00904                                                    SummaryContainer& sourceContainer)
00905 {
00906   this->setSummaryBin(bin, targetContainer.summaryXResiduals_, sourceContainer.sumXResiduals_);
00907   this->setSummaryBin(bin, targetContainer.summaryNormXResiduals_, sourceContainer.sumNormXResiduals_);
00908   // If no y-residual hists, just returns:
00909   this->setSummaryBin(bin, targetContainer.summaryYResiduals_, sourceContainer.sumYResiduals_);
00910   this->setSummaryBin(bin, targetContainer.summaryNormYResiduals_, sourceContainer.sumNormYResiduals_);
00911 }
00912 
00913 
00914 void 
00915 TrackerOfflineValidation::summarizeBinInContainer( int bin, uint32_t subDetId, 
00916                                                    SummaryContainer& targetContainer, 
00917                                                    ModuleHistos& sourceContainer)
00918 {
00919   // takes two summary Containers and sets summaryBins for all histograms
00920   this->setSummaryBin(bin, targetContainer.summaryXResiduals_, sourceContainer.ResXprimeHisto);
00921   this->setSummaryBin(bin, targetContainer.summaryNormXResiduals_, sourceContainer.NormResXprimeHisto);
00922   if( this->isPixel(subDetId) || stripYResiduals_ ) {
00923     this->setSummaryBin(bin, targetContainer.summaryYResiduals_, sourceContainer.ResYprimeHisto);
00924     this->setSummaryBin(bin, targetContainer.summaryNormYResiduals_, sourceContainer.NormResYprimeHisto);
00925   }
00926 }
00927 
00928 
00929 TrackerOfflineValidation::ModuleHistos& 
00930 TrackerOfflineValidation::getHistStructFromMap(const DetId& detid)
00931 {
00932   // get a struct with histograms from the respective map
00933   // if no object exist, the reference is automatically created by the map
00934   // throw exception if non-tracker id is passed
00935   uint subdetid = detid.subdetId();
00936   if(subdetid == PixelSubdetector::PixelBarrel ) {
00937     return mPxbResiduals_[detid.rawId()];
00938   } else if (subdetid == PixelSubdetector::PixelEndcap) {
00939     return mPxeResiduals_[detid.rawId()];
00940   } else if(subdetid  == StripSubdetector::TIB) {
00941     return mTibResiduals_[detid.rawId()];
00942   } else if(subdetid  == StripSubdetector::TID) {
00943     return mTidResiduals_[detid.rawId()];
00944   } else if(subdetid  == StripSubdetector::TOB) {
00945     return mTobResiduals_[detid.rawId()];
00946   } else if(subdetid  == StripSubdetector::TEC) {
00947     return mTecResiduals_[detid.rawId()];
00948   } else {
00949     throw cms::Exception("Geometry Error")
00950       << "[TrackerOfflineValidation] Error, tried to get reference for non-tracker subdet " << subdetid 
00951       << " from detector " << detid.det();
00952     return mPxbResiduals_[0];
00953   }
00954 }
00955 
00956 
00957 // ------------ method called to for each event  ------------
00958 void
00959 TrackerOfflineValidation::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00960 {
00961   if (useOverflowForRMS_)TH1::StatOverflows(kTRUE);
00962   this->checkBookHists(iSetup); // check whether hists are booked and do so if not yet done
00963   
00964   TrackerValidationVariables avalidator_(iSetup,parSet_);
00965     
00966   std::vector<TrackerValidationVariables::AVTrackStruct> vTrackstruct;
00967   avalidator_.fillTrackQuantities(iEvent, vTrackstruct);
00968   
00969   for (std::vector<TrackerValidationVariables::AVTrackStruct>::const_iterator itT = vTrackstruct.begin();        
00970        itT != vTrackstruct.end();
00971        ++itT) {
00972     
00973     // Fill 1D track histos
00974     static const int etaindex = this->GetIndex(vTrackHistos_,"h_tracketa");
00975     vTrackHistos_[etaindex]->Fill(itT->eta);
00976     static const int phiindex = this->GetIndex(vTrackHistos_,"h_trackphi");
00977     vTrackHistos_[phiindex]->Fill(itT->phi);
00978     static const int numOfValidHitsindex = this->GetIndex(vTrackHistos_,"h_trackNumberOfValidHits");
00979     vTrackHistos_[numOfValidHitsindex]->Fill(itT->numberOfValidHits);
00980     static const int numOfLostHitsindex = this->GetIndex(vTrackHistos_,"h_trackNumberOfLostHits");
00981     vTrackHistos_[numOfLostHitsindex]->Fill(itT->numberOfLostHits);
00982     static const int kappaindex = this->GetIndex(vTrackHistos_,"h_curvature");
00983     vTrackHistos_[kappaindex]->Fill(itT->kappa);
00984     static const int kappaposindex = this->GetIndex(vTrackHistos_,"h_curvature_pos");
00985     if (itT->charge > 0)
00986       vTrackHistos_[kappaposindex]->Fill(fabs(itT->kappa));
00987     static const int kappanegindex = this->GetIndex(vTrackHistos_,"h_curvature_neg");
00988     if (itT->charge < 0)
00989       vTrackHistos_[kappanegindex]->Fill(fabs(itT->kappa));
00990     static const int normchi2index = this->GetIndex(vTrackHistos_,"h_normchi2");
00991     vTrackHistos_[normchi2index]->Fill(itT->normchi2);
00992     static const int chi2index = this->GetIndex(vTrackHistos_,"h_chi2");
00993     vTrackHistos_[chi2index]->Fill(itT->chi2);
00994     static const int chi2Probindex = this->GetIndex(vTrackHistos_,"h_chi2Prob");
00995     vTrackHistos_[chi2Probindex]->Fill(itT->chi2Prob);
00996     static const int ptindex = this->GetIndex(vTrackHistos_,"h_pt");
00997     vTrackHistos_[ptindex]->Fill(itT->pt);
00998     if (itT->ptError != 0.) {
00999       static const int ptResolutionindex = this->GetIndex(vTrackHistos_,"h_ptResolution");
01000       vTrackHistos_[ptResolutionindex]->Fill(itT->ptError/itT->pt);
01001     }
01002     // Fill track profiles
01003     static const int d0phiindex = this->GetIndex(vTrackProfiles_,"p_d0_vs_phi");
01004     vTrackProfiles_[d0phiindex]->Fill(itT->phi,itT->d0);
01005     static const int dzphiindex = this->GetIndex(vTrackProfiles_,"p_dz_vs_phi");
01006     vTrackProfiles_[dzphiindex]->Fill(itT->phi,itT->dz);
01007     static const int d0etaindex = this->GetIndex(vTrackProfiles_,"p_d0_vs_eta");
01008     vTrackProfiles_[d0etaindex]->Fill(itT->eta,itT->d0);
01009     static const int dzetaindex = this->GetIndex(vTrackProfiles_,"p_dz_vs_eta");
01010     vTrackProfiles_[dzetaindex]->Fill(itT->eta,itT->dz);
01011     static const int chiphiindex = this->GetIndex(vTrackProfiles_,"p_chi2_vs_phi");
01012     vTrackProfiles_[chiphiindex]->Fill(itT->phi,itT->chi2);
01013     static const int chiProbphiindex = this->GetIndex(vTrackProfiles_,"p_chi2Prob_vs_phi");
01014     vTrackProfiles_[chiProbphiindex]->Fill(itT->phi,itT->chi2Prob);
01015     static const int normchiphiindex = this->GetIndex(vTrackProfiles_,"p_normchi2_vs_phi");
01016     vTrackProfiles_[normchiphiindex]->Fill(itT->phi,itT->normchi2);
01017     static const int chietaindex = this->GetIndex(vTrackProfiles_,"p_chi2_vs_eta");
01018     vTrackProfiles_[chietaindex]->Fill(itT->eta,itT->chi2);
01019     static const int normchiptindex = this->GetIndex(vTrackProfiles_,"p_normchi2_vs_pt");
01020     vTrackProfiles_[normchiptindex]->Fill(itT->pt,itT->normchi2);
01021     static const int normchipindex = this->GetIndex(vTrackProfiles_,"p_normchi2_vs_p");
01022     vTrackProfiles_[normchipindex]->Fill(itT->p,itT->normchi2);
01023     static const int chiProbetaindex = this->GetIndex(vTrackProfiles_,"p_chi2Prob_vs_eta");
01024     vTrackProfiles_[chiProbetaindex]->Fill(itT->eta,itT->chi2Prob);
01025     static const int normchietaindex = this->GetIndex(vTrackProfiles_,"p_normchi2_vs_eta");
01026     vTrackProfiles_[normchietaindex]->Fill(itT->eta,itT->normchi2);
01027     static const int kappaphiindex = this->GetIndex(vTrackProfiles_,"p_kappa_vs_phi");
01028     vTrackProfiles_[kappaphiindex]->Fill(itT->phi,itT->kappa);
01029     static const int kappaetaindex = this->GetIndex(vTrackProfiles_,"p_kappa_vs_eta");
01030     vTrackProfiles_[kappaetaindex]->Fill(itT->eta,itT->kappa);
01031     static const int ptResphiindex = this->GetIndex(vTrackProfiles_,"p_ptResolution_vs_phi");
01032     vTrackProfiles_[ptResphiindex]->Fill(itT->phi,itT->ptError/itT->pt);
01033     static const int ptResetaindex = this->GetIndex(vTrackProfiles_,"p_ptResolution_vs_eta");
01034     vTrackProfiles_[ptResetaindex]->Fill(itT->eta,itT->ptError/itT->pt);
01035 
01036     // Fill 2D track histos
01037     static const int d0phiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_d0_vs_phi");
01038     vTrack2DHistos_[d0phiindex_2d]->Fill(itT->phi,itT->d0);
01039     static const int dzphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_dz_vs_phi");
01040     vTrack2DHistos_[dzphiindex_2d]->Fill(itT->phi,itT->dz);
01041     static const int d0etaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_d0_vs_eta");
01042     vTrack2DHistos_[d0etaindex_2d]->Fill(itT->eta,itT->d0);
01043     static const int dzetaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_dz_vs_eta");
01044     vTrack2DHistos_[dzetaindex_2d]->Fill(itT->eta,itT->dz);
01045     static const int chiphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_chi2_vs_phi");
01046     vTrack2DHistos_[chiphiindex_2d]->Fill(itT->phi,itT->chi2);
01047     static const int chiProbphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_chi2Prob_vs_phi");
01048     vTrack2DHistos_[chiProbphiindex_2d]->Fill(itT->phi,itT->chi2Prob);
01049     static const int normchiphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_normchi2_vs_phi");
01050     vTrack2DHistos_[normchiphiindex_2d]->Fill(itT->phi,itT->normchi2);
01051     static const int chietaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_chi2_vs_eta");
01052     vTrack2DHistos_[chietaindex_2d]->Fill(itT->eta,itT->chi2);
01053     static const int chiProbetaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_chi2Prob_vs_eta");
01054     vTrack2DHistos_[chiProbetaindex_2d]->Fill(itT->eta,itT->chi2Prob);
01055     static const int normchietaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_normchi2_vs_eta");
01056     vTrack2DHistos_[normchietaindex_2d]->Fill(itT->eta,itT->normchi2);
01057     static const int kappaphiindex_2d = this->GetIndex(vTrack2DHistos_,"h2_kappa_vs_phi");
01058     vTrack2DHistos_[kappaphiindex_2d]->Fill(itT->phi,itT->kappa);
01059     static const int kappaetaindex_2d = this->GetIndex(vTrack2DHistos_,"h2_kappa_vs_eta");
01060     vTrack2DHistos_[kappaetaindex_2d]->Fill(itT->eta,itT->kappa);
01061     static const int normchi2kappa_2d = this->GetIndex(vTrack2DHistos_,"h2_normchi2_vs_kappa");
01062     vTrack2DHistos_[normchi2kappa_2d]->Fill(itT->normchi2,itT->kappa);
01063 
01064     // hit quantities: residuals, normalized residuals
01065     for (std::vector<TrackerValidationVariables::AVHitStruct>::const_iterator itH = itT->hits.begin();
01066          itH != itT->hits.end();
01067          ++itH) {
01068       
01069       DetId detid(itH->rawDetId);
01070       ModuleHistos &histStruct = this->getHistStructFromMap(detid);
01071 
01072       // fill histos in local coordinates if set in cf
01073       if (lCoorHistOn_) {
01074         histStruct.ResHisto->Fill(itH->resX);
01075         if(itH->resErrX != 0) histStruct.NormResHisto->Fill(itH->resX/itH->resErrX);
01076         if (this->isPixel(detid.subdetId()) || stripYResiduals_ ) {
01077           histStruct.ResYHisto->Fill(itH->resY);
01078           // here add un-primed normalised y-residuals if wanted
01079         }
01080       }
01081       if (itH->resXprime != -999.) {
01082         histStruct.ResXprimeHisto->Fill(itH->resXprime);
01083 
01084         /******************************* Fill 2-D histo ResX vs momenta *****************************/
01085         if (detid.subdetId() == PixelSubdetector::PixelBarrel) { 
01086           static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_pixB");
01087           vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
01088         }       
01089         if (detid.subdetId() == PixelSubdetector::PixelEndcap) { 
01090           static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_pixE");
01091           vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
01092         }
01093         if (detid.subdetId() == StripSubdetector::TIB) { 
01094           static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_TIB");
01095           vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
01096         }
01097         if (detid.subdetId() == StripSubdetector::TID) { 
01098           static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_TID");
01099           vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
01100         }
01101         if (detid.subdetId() == StripSubdetector::TOB) { 
01102           static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_TOB");
01103           vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
01104         }
01105         if (detid.subdetId() == StripSubdetector::TEC) {   
01106           static const int resXvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resXprime_TEC");
01107           vTrack2DHistos_[resXvsPindex_2d]->Fill(itT->p,itH->resXprime);
01108         }
01109         /******************************************/
01110 
01111         if ( moduleLevelProfiles_ && itH->inside ) {
01112           float tgalpha = tan(itH->localAlpha);
01113           if ( fabs(tgalpha)!=0 ){
01114             histStruct.LocalX->Fill(itH->localXnorm, tgalpha*tgalpha); 
01115             histStruct.LocalY->Fill(itH->localYnorm, tgalpha*tgalpha); 
01116             histStruct.ResXvsXProfile->Fill(itH->localXnorm, itH->resX/tgalpha, tgalpha*tgalpha); 
01117             histStruct.ResXvsYProfile->Fill(itH->localYnorm, itH->resX/tgalpha, tgalpha*tgalpha); 
01118           }
01119         }
01120 
01121         if(itH->resXprimeErr != 0 && itH->resXprimeErr != -999 ) {      
01122           histStruct.NormResXprimeHisto->Fill(itH->resXprime/itH->resXprimeErr);
01123         }
01124       }     
01125       
01126       if (itH->resYprime != -999.) {
01127         if (this->isPixel(detid.subdetId()) || stripYResiduals_ ) {
01128           histStruct.ResYprimeHisto->Fill(itH->resYprime);
01129 
01130           /******************************* Fill 2-D histo ResY vs momenta *****************************/
01131           if (detid.subdetId() == PixelSubdetector::PixelBarrel) { 
01132             static const int resYvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resYprime_pixB");
01133             vTrack2DHistos_[resYvsPindex_2d]->Fill(itT->p,itH->resYprime);
01134           }     
01135           if (detid.subdetId() == PixelSubdetector::PixelEndcap) { 
01136             static const int resYvsPindex_2d = this->GetIndex(vTrack2DHistos_,"p_vs_resYprime_pixE");
01137             vTrack2DHistos_[resYvsPindex_2d]->Fill(itT->p,itH->resYprime);
01138           }
01139           /******************************************/
01140 
01141           if ( moduleLevelProfiles_ && itH->inside ) {
01142             float tgbeta = tan(itH->localBeta);
01143             if ( fabs(tgbeta)!=0 ){
01144               histStruct.ResYvsXProfile->Fill(itH->localXnorm, itH->resY/tgbeta, tgbeta*tgbeta); 
01145               histStruct.ResYvsYProfile->Fill(itH->localYnorm, itH->resY/tgbeta, tgbeta*tgbeta); 
01146             }
01147           }
01148 
01149           if (itH->resYprimeErr != 0 && itH->resYprimeErr != -999. ) {  
01150             histStruct.NormResYprimeHisto->Fill(itH->resYprime/itH->resYprimeErr);
01151           } 
01152         }
01153       }
01154       
01155       if (overlappOn_) {
01156         std::pair<uint32_t,uint32_t> tmp_pair(std::make_pair(itH->rawDetId, itH->overlapres.first));
01157         if (itH->overlapres.first != 0) {
01158           if (hOverlappResidual[tmp_pair]) {
01159             hOverlappResidual[tmp_pair]->Fill(itH->overlapres.second);
01160           } else if (hOverlappResidual[std::make_pair(itH->overlapres.first, itH->rawDetId)]) {
01161             hOverlappResidual[std::make_pair( itH->overlapres.first, itH->rawDetId)]->Fill(itH->overlapres.second);
01162           } else {
01163             std::string dirname("OverlappResiduals");
01164             DirectoryWrapper tfd(dirname.c_str(),moduleDirectory_,dqmMode_);
01165             hOverlappResidual[tmp_pair] = tfd.make<TH1F>(Form("hOverlappResidual_%d_%d",tmp_pair.first,tmp_pair.second),
01166                                                          "Overlapp Residuals",100,-50,50);
01167             hOverlappResidual[tmp_pair]->Fill(itH->overlapres.second);
01168           }
01169         }
01170       }
01171       
01172     } // finish loop over hit quantities
01173   } // finish loop over track quantities
01174 
01175   if (useOverflowForRMS_) TH1::StatOverflows(kFALSE);  
01176 }
01177 
01178 
01179 // ------------ method called once each job just after ending the event loop  ------------
01180 void 
01181 TrackerOfflineValidation::endJob()
01182 {
01183 
01184   if (!tkGeom_.product()) return;
01185 
01186   AlignableTracker aliTracker(&(*tkGeom_));
01187   
01188   AlignableObjectId aliobjid;
01189   
01190   static const int kappadiffindex = this->GetIndex(vTrackHistos_,"h_diff_curvature");
01191   vTrackHistos_[kappadiffindex]->Add(vTrackHistos_[this->GetIndex(vTrackHistos_,"h_curvature_neg")],
01192                                      vTrackHistos_[this->GetIndex(vTrackHistos_,"h_curvature_pos")],-1,1);
01193  
01194   // Collate Information for Subdetectors
01195   // create summary histogramms recursively
01196   std::vector<TrackerOfflineValidation::SummaryContainer> vTrackerprofiles;
01197   DirectoryWrapper f("",moduleDirectory_,dqmMode_);
01198   this->collateSummaryHists(f,(aliTracker), 0, aliobjid, vTrackerprofiles);
01199   
01200   if (dqmMode_) return;
01201   // Should be excluded in dqmMode, since TTree is not usable
01202   // In dqmMode tree operations are are sourced out to the additional module TrackerOfflineValidationSummary
01203   
01204   edm::Service<TFileService> fs;
01205   TTree *tree = fs->make<TTree>("TkOffVal","TkOffVal");
01206  
01207   TkOffTreeVariables *treeMemPtr = new TkOffTreeVariables;
01208   // We create branches for all members of 'TkOffTreeVariables' (even if not needed).
01209   // This works because we have a dictionary for 'TkOffTreeVariables'
01210   // (see src/classes_def.xml and src/classes.h):
01211   tree->Branch("TkOffTreeVariables", &treeMemPtr); // address of pointer!
01212  
01213   this->fillTree(*tree, mPxbResiduals_, *treeMemPtr, *tkGeom_);
01214   this->fillTree(*tree, mPxeResiduals_, *treeMemPtr, *tkGeom_);
01215   this->fillTree(*tree, mTibResiduals_, *treeMemPtr, *tkGeom_);
01216   this->fillTree(*tree, mTidResiduals_, *treeMemPtr, *tkGeom_);
01217   this->fillTree(*tree, mTobResiduals_, *treeMemPtr, *tkGeom_);
01218   this->fillTree(*tree, mTecResiduals_, *treeMemPtr, *tkGeom_);
01219 
01220   delete treeMemPtr; treeMemPtr = 0;
01221 }
01222 
01223 
01224 void
01225 TrackerOfflineValidation::collateSummaryHists( DirectoryWrapper& tfd, const Alignable& ali, int i, 
01226                                                const AlignableObjectId& aliobjid, 
01227                                                std::vector<TrackerOfflineValidation::SummaryContainer>& vLevelProfiles)
01228 {
01229   std::vector<Alignable*> alivec(ali.components());
01230   if( this->isDetOrDetUnit((alivec)[0]->alignableObjectId()) ) return;
01231   
01232   for(int iComp=0, iCompEnd = ali.components().size();iComp < iCompEnd; ++iComp) {
01233     std::vector< TrackerOfflineValidation::SummaryContainer > vProfiles;        
01234     std::string structurename  = aliobjid.typeToName((alivec)[iComp]->alignableObjectId());
01235     
01236     LogDebug("TrackerOfflineValidation") << "StructureName = " << structurename;
01237     std::stringstream dirname;
01238     dirname << structurename;
01239     
01240     // add no suffix counter to strip and pixel -> just aesthetics
01241     if (structurename != "Strip" && structurename != "Pixel") dirname << "_" << iComp+1;
01242     
01243     if(  !(this->isDetOrDetUnit( (alivec)[iComp]->alignableObjectId()) )
01244          || (alivec)[0]->components().size() > 1 ) {
01245       DirectoryWrapper f(tfd,dirname.str(),moduleDirectory_,dqmMode_);
01246       this->collateSummaryHists( f, *(alivec)[iComp], i, aliobjid, vProfiles);
01247       vLevelProfiles.push_back(this->bookSummaryHists(tfd, *(alivec[iComp]), ali.alignableObjectId(), iComp+1, aliobjid));
01248       TH1 *hY = vLevelProfiles[iComp].sumYResiduals_;
01249       TH1 *hNormY = vLevelProfiles[iComp].sumNormYResiduals_;
01250       for(uint n = 0; n < vProfiles.size(); ++n) {
01251         this->summarizeBinInContainer(n+1, vLevelProfiles[iComp], vProfiles[n] );
01252         vLevelProfiles[iComp].sumXResiduals_->Add(vProfiles[n].sumXResiduals_);
01253         vLevelProfiles[iComp].sumNormXResiduals_->Add(vProfiles[n].sumNormXResiduals_);
01254         if (hY)     hY->Add(vProfiles[n].sumYResiduals_);         // only if existing
01255         if (hNormY) hNormY->Add(vProfiles[n].sumNormYResiduals_); // dito (pxl, stripYResiduals_)
01256       }
01257       if(dqmMode_)continue;  // No fits in dqmMode
01258       //add fit values to stat box
01259       this->fitResiduals(vLevelProfiles[iComp].sumXResiduals_);
01260       this->fitResiduals(vLevelProfiles[iComp].sumNormXResiduals_);
01261       if (hY)     this->fitResiduals(hY);     // only if existing (pixel or stripYResiduals_)
01262       if (hNormY) this->fitResiduals(hNormY); // dito
01263     } else {
01264       // nothing to be done for det or detunits
01265       continue;
01266     }
01267   }
01268 }
01269 
01270 
01271 TrackerOfflineValidation::SummaryContainer 
01272 TrackerOfflineValidation::bookSummaryHists(DirectoryWrapper& tfd, const Alignable& ali, 
01273                                            align::StructureType type, int i, 
01274                                            const AlignableObjectId& aliobjid)
01275 {
01276   const uint aliSize = ali.components().size();
01277   const align::StructureType alitype = ali.alignableObjectId();
01278   const align::StructureType subtype = ali.components()[0]->alignableObjectId();
01279   const char *aliTypeName = aliobjid.typeToName(alitype).c_str(); // lifetime of char* OK
01280   const char *aliSubtypeName = aliobjid.typeToName(subtype).c_str();
01281   const char *typeName = aliobjid.typeToName(type).c_str();
01282 
01283   const DetId aliDetId = ali.id(); 
01284   // y residuals only if pixel or specially requested for strip:
01285   const bool bookResidY = this->isPixel(aliDetId.subdetId()) || stripYResiduals_;
01286 
01287   SummaryContainer sumContainer;
01288   
01289   // Book summary hists with one bin per component, 
01290   // but special case for Det with two DetUnit that we want to summarize one level up 
01291   // (e.g. in TOBRods with 12 bins for 6 stereo and 6 rphi DetUnit.)
01292   //    component of ali is not Det or Det with just one components
01293   const uint subcompSize = ali.components()[0]->components().size();
01294   if (subtype != align::AlignableDet || subcompSize == 1) { // Det with 1 comp. should not exist anymore...
01295     const TString title(Form("Summary for substructures in %s %d;%s;",aliTypeName,i,aliSubtypeName));
01296     
01297     sumContainer.summaryXResiduals_ = tfd.make<TH1F>(Form("h_summaryX%s_%d",aliTypeName,i), 
01298                                                      title + "#LT #Delta x' #GT",
01299                                                      aliSize, 0.5, aliSize+0.5);
01300     sumContainer.summaryNormXResiduals_ = tfd.make<TH1F>(Form("h_summaryNormX%s_%d",aliTypeName,i), 
01301                                                          title + "#LT #Delta x'/#sigma #GT",
01302                                                          aliSize,0.5,aliSize+0.5);
01303     
01304     if (bookResidY) {
01305       sumContainer.summaryYResiduals_ = tfd.make<TH1F>(Form("h_summaryY%s_%d",aliTypeName,i), 
01306                                                        title + "#LT #Delta y' #GT",
01307                                                        aliSize, 0.5, aliSize+0.5);
01308       sumContainer.summaryNormYResiduals_ = tfd.make<TH1F>(Form("h_summaryNormY%s_%d",aliTypeName,i), 
01309                                                            title + "#LT #Delta y'/#sigma #GT",
01310                                                            aliSize,0.5,aliSize+0.5);
01311     }
01312     
01313   } else if (subtype == align::AlignableDet && subcompSize > 1) { // fixed: was aliSize before
01314     if (subcompSize != 2) { // strange... expect only 2 DetUnits in DS layers
01315       // this 2 is hardcoded factor 2 in binning below and also assumed later on
01316       edm::LogError("Alignment") << "@SUB=bookSummaryHists"
01317                                  << "Det with " << subcompSize << " components";
01318     }
01319     // title contains x-title
01320     const TString title(Form("Summary for substructures in %s %d;%s;", aliTypeName, i,
01321                              aliobjid.typeToName(ali.components()[0]->components()[0]->alignableObjectId()).c_str()));
01322     
01323     sumContainer.summaryXResiduals_ 
01324       = tfd.make<TH1F>(Form("h_summaryX%s_%d", aliTypeName, i), 
01325                        title + "#LT #Delta x' #GT", (2*aliSize), 0.5, 2*aliSize+0.5);
01326     sumContainer.summaryNormXResiduals_ 
01327       = tfd.make<TH1F>(Form("h_summaryNormX%s_%d", aliTypeName, i), 
01328                        title + "#LT #Delta x'/#sigma #GT", (2*aliSize), 0.5, 2*aliSize+0.5);
01329 
01330     if (bookResidY) {
01331       sumContainer.summaryYResiduals_ 
01332         = tfd.make<TH1F>(Form("h_summaryY%s_%d", aliTypeName, i), 
01333                          title + "#LT #Delta y' #GT", (2*aliSize), 0.5, 2*aliSize+0.5);
01334       sumContainer.summaryNormYResiduals_ 
01335         = tfd.make<TH1F>(Form("h_summaryNormY%s_%d", aliTypeName, i), 
01336                          title + "#LT #Delta y'/#sigma #GT", (2*aliSize), 0.5, 2*aliSize+0.5);
01337     }
01338 
01339   } else {
01340     edm::LogError("TrackerOfflineValidation") << "@SUB=TrackerOfflineValidation::bookSummaryHists" 
01341                                               << "No summary histogramm for hierarchy level " 
01342                                               << aliTypeName << " in subdet " << aliDetId.subdetId();
01343   }
01344 
01345   // Now book hists that just sum up the residual histograms from lower levels.
01346   // Axis title is copied from lowest level module of structure.
01347   // Should be safe that y-hists are only touched if non-null pointers...
01348   int nbins = 0;
01349   double xmin = 0., xmax = 0.;
01350   const TString sumTitle(Form("Residual for %s %d in %s;", aliTypeName, i, typeName));
01351   const ModuleHistos &xTitHists = this->getHistStructFromMap(aliDetId); // for x-axis titles
01352   this->getBinning(aliDetId.subdetId(), XprimeResidual, nbins, xmin, xmax);
01353   
01354   sumContainer.sumXResiduals_ = tfd.make<TH1F>(Form("h_Xprime_%s_%d", aliTypeName, i),
01355                                                sumTitle + xTitHists.ResXprimeHisto->GetXaxis()->GetTitle(),
01356                                                nbins, xmin, xmax);
01357   
01358   this->getBinning(aliDetId.subdetId(), NormXprimeResidual, nbins, xmin, xmax);
01359   sumContainer.sumNormXResiduals_ = tfd.make<TH1F>(Form("h_NormXprime_%s_%d",aliTypeName,i), 
01360                                                    sumTitle + xTitHists.NormResXprimeHisto->GetXaxis()->GetTitle(),
01361                                                    nbins, xmin, xmax);
01362   if (bookResidY) {
01363     this->getBinning(aliDetId.subdetId(), YprimeResidual, nbins, xmin, xmax);
01364     sumContainer.sumYResiduals_ = tfd.make<TH1F>(Form("h_Yprime_%s_%d",aliTypeName,i), 
01365                                                  sumTitle + xTitHists.ResYprimeHisto->GetXaxis()->GetTitle(),
01366                                                  nbins, xmin, xmax);
01367     
01368     this->getBinning(aliDetId.subdetId(), NormYprimeResidual, nbins, xmin, xmax);
01369     sumContainer.sumNormYResiduals_ = tfd.make<TH1F>(Form("h_NormYprime_%s_%d",aliTypeName,i), 
01370                                                      sumTitle + xTitHists.NormResYprimeHisto->GetXaxis()->GetTitle(),
01371                                                      nbins, xmin, xmax);
01372   }
01373   
01374   // If we are at the lowest level, we already sum up and fill the summary.
01375 
01376   // special case I: For DetUnits and Dets with only one subcomponent start filling summary histos
01377   if( ( subtype == align::AlignableDet && subcompSize == 1) || subtype  == align::AlignableDetUnit ) {
01378     for(uint k = 0; k < aliSize; ++k) {
01379       DetId detid = ali.components()[k]->id();
01380       ModuleHistos &histStruct = this->getHistStructFromMap(detid);
01381       this->summarizeBinInContainer(k+1, detid.subdetId() ,sumContainer, histStruct );
01382       sumContainer.sumXResiduals_->Add(histStruct.ResXprimeHisto);
01383       sumContainer.sumNormXResiduals_->Add(histStruct.NormResXprimeHisto);
01384       if( this->isPixel(detid.subdetId()) || stripYResiduals_ ) {
01385         sumContainer.sumYResiduals_->Add(histStruct.ResYprimeHisto);
01386         sumContainer.sumNormYResiduals_->Add(histStruct.NormResYprimeHisto);
01387       }
01388     }
01389   } else if( subtype == align::AlignableDet && subcompSize > 1) { // fixed: was aliSize before
01390     // special case II: Fill summary histos for dets with two detunits 
01391     for(uint k = 0; k < aliSize; ++k) {
01392       for(uint j = 0; j < subcompSize; ++j) { // assumes all have same size (as binning does)
01393         DetId detid = ali.components()[k]->components()[j]->id();
01394         ModuleHistos &histStruct = this->getHistStructFromMap(detid);   
01395         this->summarizeBinInContainer(2*k+j+1, detid.subdetId() ,sumContainer, histStruct );
01396         sumContainer.sumXResiduals_->Add( histStruct.ResXprimeHisto);
01397         sumContainer.sumNormXResiduals_->Add( histStruct.NormResXprimeHisto);
01398         if( this->isPixel(detid.subdetId()) || stripYResiduals_ ) {
01399           sumContainer.sumYResiduals_->Add( histStruct.ResYprimeHisto);
01400           sumContainer.sumNormYResiduals_->Add( histStruct.NormResYprimeHisto);
01401         }
01402       }
01403     }
01404   }
01405   return sumContainer;
01406 }
01407 
01408 
01409 float 
01410 TrackerOfflineValidation::Fwhm (const TH1* hist) const
01411 {
01412   float max = hist->GetMaximum();
01413   int left = -1, right = -1;
01414   for(unsigned int i = 1, iEnd = hist->GetNbinsX(); i <= iEnd; ++i) {
01415     if(hist->GetBinContent(i) < max/2. && hist->GetBinContent(i+1) > max/2. && left == -1) {
01416       if(max/2. - hist->GetBinContent(i) < hist->GetBinContent(i+1) - max/2.) {
01417         left = i;
01418         ++i;
01419       } else {
01420         left = i+1;
01421         ++i;
01422       }
01423     }
01424     if(left != -1 && right == -1) {
01425       if(hist->GetBinContent(i) > max/2. && hist->GetBinContent(i+1) < max/2.) {
01426         if( hist->GetBinContent(i) - max/2. < max/2. - hist->GetBinContent(i+1)) {
01427           right = i;
01428         } else {
01429           right = i+1;
01430         }
01431         
01432       }
01433     }
01434   }
01435   return hist->GetXaxis()->GetBinCenter(right) - hist->GetXaxis()->GetBinCenter(left);
01436 }
01437 
01438 
01439 void 
01440 TrackerOfflineValidation::fillTree(TTree& tree,
01441                                    const std::map<int, TrackerOfflineValidation::ModuleHistos>& moduleHist_,
01442                                    TkOffTreeVariables &treeMem, const TrackerGeometry& tkgeom)
01443 {
01444  
01445   for(std::map<int, TrackerOfflineValidation::ModuleHistos>::const_iterator it = moduleHist_.begin(), 
01446         itEnd= moduleHist_.end(); it != itEnd;++it ) { 
01447     treeMem.clear(); // make empty/default
01448     
01449     //variables concerning the tracker components/hierarchy levels
01450     DetId detId_ = it->first;
01451     treeMem.moduleId = detId_;
01452     treeMem.subDetId = detId_.subdetId();
01453     treeMem.isDoubleSide =0;
01454 
01455     if(treeMem.subDetId == PixelSubdetector::PixelBarrel){
01456       PXBDetId pxbId(detId_);
01457       unsigned int whichHalfBarrel(1), rawId(detId_.rawId());  //DetId does not know about halfBarrels is PXB ...
01458       if( (rawId>=302056964 && rawId<302059300) || (rawId>=302123268 && rawId<302127140) || (rawId>=302189572 && rawId<302194980) )whichHalfBarrel=2;
01459       treeMem.layer = pxbId.layer(); 
01460       treeMem.half = whichHalfBarrel;
01461       treeMem.rod = pxbId.ladder();     // ... so, ladder is not per halfBarrel-Layer, but per barrel-layer!
01462       treeMem.module = pxbId.module();
01463     } else if(treeMem.subDetId == PixelSubdetector::PixelEndcap){
01464       PXFDetId pxfId(detId_); 
01465       unsigned int whichHalfCylinder(1), rawId(detId_.rawId());  //DetId does not kmow about halfCylinders in PXF
01466       if( (rawId>=352394500 && rawId<352406032) || (rawId>=352460036 && rawId<352471568) || (rawId>=344005892 && rawId<344017424) || (rawId>=344071428 && rawId<344082960) )whichHalfCylinder=2;
01467       treeMem.layer = pxfId.disk(); 
01468       treeMem.side = pxfId.side();
01469       treeMem.half = whichHalfCylinder;
01470       treeMem.blade = pxfId.blade(); 
01471       treeMem.panel = pxfId.panel();
01472       treeMem.module = pxfId.module();
01473     } else if(treeMem.subDetId == StripSubdetector::TIB){
01474       TIBDetId tibId(detId_); 
01475       unsigned int whichHalfShell(1), rawId(detId_.rawId());  //DetId does not kmow about halfShells in TIB
01476        if ( (rawId>=369120484 && rawId<369120688) || (rawId>=369121540 && rawId<369121776) || (rawId>=369136932 && rawId<369137200) || (rawId>=369137988 && rawId<369138288) ||
01477             (rawId>=369153396 && rawId<369153744) || (rawId>=369154436 && rawId<369154800) || (rawId>=369169844 && rawId<369170256) || (rawId>=369170900 && rawId<369171344) ||
01478             (rawId>=369124580 && rawId<369124784) || (rawId>=369125636 && rawId<369125872) || (rawId>=369141028 && rawId<369141296) || (rawId>=369142084 && rawId<369142384) ||
01479             (rawId>=369157492 && rawId<369157840) || (rawId>=369158532 && rawId<369158896) || (rawId>=369173940 && rawId<369174352) || (rawId>=369174996 && rawId<369175440) ) whichHalfShell=2;
01480       treeMem.layer = tibId.layer(); 
01481       treeMem.side = tibId.string()[0];
01482       treeMem.half = whichHalfShell;
01483       treeMem.rod = tibId.string()[2]; 
01484       treeMem.outerInner = tibId.string()[1]; 
01485       treeMem.module = tibId.module();
01486       treeMem.isStereo = tibId.stereo();
01487       treeMem.isDoubleSide = tibId.isDoubleSide();
01488     } else if(treeMem.subDetId == StripSubdetector::TID){
01489       TIDDetId tidId(detId_); 
01490       treeMem.layer = tidId.wheel(); 
01491       treeMem.side = tidId.side();
01492       treeMem.ring = tidId.ring(); 
01493       treeMem.outerInner = tidId.module()[0]; 
01494       treeMem.module = tidId.module()[1];
01495       treeMem.isStereo = tidId.stereo();
01496       treeMem.isDoubleSide = tidId.isDoubleSide();
01497     } else if(treeMem.subDetId == StripSubdetector::TOB){
01498       TOBDetId tobId(detId_); 
01499       treeMem.layer = tobId.layer(); 
01500       treeMem.side = tobId.rod()[0];
01501       treeMem.rod = tobId.rod()[1]; 
01502       treeMem.module = tobId.module();
01503       treeMem.isStereo = tobId.stereo();
01504       treeMem.isDoubleSide = tobId.isDoubleSide();
01505     } else if(treeMem.subDetId == StripSubdetector::TEC) {
01506       TECDetId tecId(detId_); 
01507       treeMem.layer = tecId.wheel(); 
01508       treeMem.side  = tecId.side();
01509       treeMem.ring  = tecId.ring(); 
01510       treeMem.petal = tecId.petal()[1]; 
01511       treeMem.outerInner = tecId.petal()[0];
01512       treeMem.module = tecId.module();
01513       treeMem.isStereo = tecId.stereo();
01514       treeMem.isDoubleSide = tecId.isDoubleSide(); 
01515     }
01516     
01517     //variables concerning the tracker geometry
01518     const Surface::PositionType &gPModule = tkgeom.idToDet(detId_)->position();
01519     treeMem.posPhi = gPModule.phi();
01520     treeMem.posEta = gPModule.eta();
01521     treeMem.posR   = gPModule.perp();
01522     treeMem.posX   = gPModule.x();
01523     treeMem.posY   = gPModule.y();
01524     treeMem.posZ   = gPModule.z();
01525  
01526     const Surface& surface =  tkgeom.idToDet(detId_)->surface();
01527     
01528     //global Orientation of local coordinate system of dets/detUnits   
01529     LocalPoint  lUDirection(1.,0.,0.), lVDirection(0.,1.,0.), lWDirection(0.,0.,1.);
01530     GlobalPoint gUDirection = surface.toGlobal(lUDirection),
01531                 gVDirection = surface.toGlobal(lVDirection),
01532                 gWDirection = surface.toGlobal(lWDirection);
01533     double dR(999.), dPhi(999.), dZ(999.);
01534     if(treeMem.subDetId==PixelSubdetector::PixelBarrel || treeMem.subDetId==StripSubdetector::TIB || treeMem.subDetId==StripSubdetector::TOB){
01535       dR = gWDirection.perp() - gPModule.perp();
01536       dPhi = deltaPhi(gUDirection.phi(),gPModule.phi());
01537       dZ = gVDirection.z() - gPModule.z();
01538       if(dZ>=0.)treeMem.rOrZDirection = 1; else treeMem.rOrZDirection = -1;
01539     }else if(treeMem.subDetId==PixelSubdetector::PixelEndcap){
01540       dR = gUDirection.perp() - gPModule.perp();
01541       dPhi = deltaPhi(gVDirection.phi(),gPModule.phi());
01542       dZ = gWDirection.z() - gPModule.z();
01543       if(dR>=0.)treeMem.rOrZDirection = 1; else treeMem.rOrZDirection = -1;
01544     }else if(treeMem.subDetId==StripSubdetector::TID || treeMem.subDetId==StripSubdetector::TEC){
01545       dR = gVDirection.perp() - gPModule.perp();
01546       dPhi = deltaPhi(gUDirection.phi(),gPModule.phi());
01547       dZ = gWDirection.z() - gPModule.z();
01548       if(dR>=0.)treeMem.rOrZDirection = 1; else treeMem.rOrZDirection = -1;
01549     }
01550     if(dR>=0.)treeMem.rDirection = 1; else treeMem.rDirection = -1;
01551     if(dPhi>=0.)treeMem.phiDirection = 1; else treeMem.phiDirection = -1;
01552     if(dZ>=0.)treeMem.zDirection = 1; else treeMem.zDirection = -1;
01553     
01554     //mean and RMS values (extracted from histograms(Xprime on module level)
01555     treeMem.entries = static_cast<UInt_t>(it->second.ResXprimeHisto->GetEntries());
01556     treeMem.meanX = it->second.ResXprimeHisto->GetMean();
01557     treeMem.rmsX  = it->second.ResXprimeHisto->GetRMS();
01558     //treeMem.sigmaX = Fwhm(it->second.ResXprimeHisto)/2.355;
01559     
01560     if (useFit_) {
01561       //call fit function which returns mean and sigma from the fit
01562       //for absolute residuals
01563       std::pair<float,float> fitResult1 = this->fitResiduals(it->second.ResXprimeHisto);
01564       treeMem.fitMeanX = fitResult1.first;
01565       treeMem.fitSigmaX = fitResult1.second;
01566       //for normalized residuals
01567       std::pair<float,float> fitResult2 = this->fitResiduals(it->second.NormResXprimeHisto);
01568       treeMem.fitMeanNormX = fitResult2.first;
01569       treeMem.fitSigmaNormX = fitResult2.second;
01570     }
01571     
01572     //get median for absolute residuals
01573     treeMem.medianX   = this->getMedian(it->second.ResXprimeHisto);
01574 
01575     int numberOfBins=it->second.ResXprimeHisto->GetNbinsX();
01576     treeMem.numberOfUnderflows = it->second.ResXprimeHisto->GetBinContent(0);
01577     treeMem.numberOfOverflows = it->second.ResXprimeHisto->GetBinContent(numberOfBins+1);
01578     treeMem.numberOfOutliers =  it->second.ResXprimeHisto->GetBinContent(0)+it->second.ResXprimeHisto->GetBinContent(numberOfBins+1);
01579     
01580     //mean and RMS values (extracted from histograms(normalized Xprime on module level)
01581     treeMem.meanNormX = it->second.NormResXprimeHisto->GetMean();
01582     treeMem.rmsNormX = it->second.NormResXprimeHisto->GetRMS();
01583 
01584     double stats[20];
01585     it->second.NormResXprimeHisto->GetStats(stats);
01586     // GF  treeMem.chi2PerDofX = stats[3]/(stats[0]-1);
01587     if (stats[0]) treeMem.chi2PerDofX = stats[3]/stats[0];
01588     
01589     treeMem.sigmaNormX = Fwhm(it->second.NormResXprimeHisto)/2.355;
01590     treeMem.histNameX = it->second.ResXprimeHisto->GetName();
01591     treeMem.histNameNormX = it->second.NormResXprimeHisto->GetName();
01592     
01593     // fill tree variables in local coordinates if set in cfg
01594     if(lCoorHistOn_) {
01595       treeMem.meanLocalX = it->second.ResHisto->GetMean();
01596       treeMem.rmsLocalX = it->second.ResHisto->GetRMS();
01597       treeMem.meanNormLocalX = it->second.NormResHisto->GetMean();
01598       treeMem.rmsNormLocalX = it->second.NormResHisto->GetRMS();
01599 
01600       treeMem.histNameLocalX = it->second.ResHisto->GetName();
01601       treeMem.histNameNormLocalX = it->second.NormResHisto->GetName();
01602       if (it->second.ResYHisto) treeMem.histNameLocalY = it->second.ResYHisto->GetName();
01603     }
01604 
01605     // mean and RMS values in local y (extracted from histograms(normalized Yprime on module level)
01606     // might exist in pixel only
01607     if (it->second.ResYprimeHisto) {//(stripYResiduals_){
01608       TH1 *h = it->second.ResYprimeHisto;
01609       treeMem.meanY = h->GetMean();
01610       treeMem.rmsY  = h->GetRMS();
01611       
01612       if (useFit_) { // fit function which returns mean and sigma from the fit
01613         std::pair<float,float> fitMeanSigma = this->fitResiduals(h);
01614         treeMem.fitMeanY  = fitMeanSigma.first;
01615         treeMem.fitSigmaY = fitMeanSigma.second;
01616       }
01617       
01618       //get median for absolute residuals
01619       treeMem.medianY   = this->getMedian(h);
01620 
01621       treeMem.histNameY = h->GetName();
01622     }
01623     if (it->second.NormResYprimeHisto) {
01624       TH1 *h = it->second.NormResYprimeHisto;
01625       treeMem.meanNormY = h->GetMean();
01626       treeMem.rmsNormY  = h->GetRMS();
01627       h->GetStats(stats); // stats buffer defined above
01628       if (stats[0]) treeMem.chi2PerDofY = stats[3]/stats[0];
01629 
01630       if (useFit_) { // fit function which returns mean and sigma from the fit
01631         std::pair<float,float> fitMeanSigma = this->fitResiduals(h);
01632         treeMem.fitMeanNormY  = fitMeanSigma.first;
01633         treeMem.fitSigmaNormY = fitMeanSigma.second;
01634       }
01635       treeMem.histNameNormY = h->GetName();
01636     }
01637     tree.Fill();
01638   }
01639 }
01640 
01641 
01642 std::pair<float,float> 
01643 TrackerOfflineValidation::fitResiduals(TH1* hist) const
01644 {
01645   std::pair<float,float> fitResult(9999., 9999.);
01646   if (!hist || hist->GetEntries() < 20) return fitResult;
01647 
01648   float mean  = hist->GetMean();
01649   float sigma = hist->GetRMS();
01650 
01651   try { // for < CMSSW_2_2_0 since ROOT warnings from fit are converted to exceptions
01652     // Remove the try/catch for more recent CMSSW!
01653     // first fit: two RMS around mean
01654     TF1 func("tmp", "gaus", mean - 2.*sigma, mean + 2.*sigma); 
01655     if (0 == hist->Fit(&func,"QNR")) { // N: do not blow up file by storing fit!
01656       mean  = func.GetParameter(1);
01657       sigma = func.GetParameter(2);
01658       // second fit: three sigma of first fit around mean of first fit
01659       func.SetRange(mean - 3.*sigma, mean + 3.*sigma);
01660       // I: integral gives more correct results if binning is too wide
01661       // L: Likelihood can treat empty bins correctly (if hist not weighted...)
01662       if (0 == hist->Fit(&func, "Q0LR")) {
01663         if (hist->GetFunction(func.GetName())) { // Take care that it is later on drawn:
01664           hist->GetFunction(func.GetName())->ResetBit(TF1::kNotDraw);
01665         }
01666         fitResult.first = func.GetParameter(1);
01667         fitResult.second = func.GetParameter(2);
01668       }
01669     }
01670   } catch (cms::Exception const & e) {
01671     edm::LogWarning("Alignment") << "@SUB=TrackerOfflineValidation::fitResiduals"
01672                                  << "Caught this exception during ROOT fit: "
01673                                  << e.what();
01674   }
01675   return fitResult;
01676 }
01677 
01678 
01679 float 
01680 TrackerOfflineValidation::getMedian(const TH1* histo) const
01681 {
01682   float median = 999;
01683   int nbins = histo->GetNbinsX();
01684 
01685   //extract median from histogram
01686   double *x = new double[nbins];
01687   double *y = new double[nbins];
01688   for (int j = 0; j < nbins; j++) {
01689     x[j] = histo->GetBinCenter(j+1);
01690     y[j] = histo->GetBinContent(j+1);
01691   }
01692   median = TMath::Median(nbins, x, y);
01693   
01694   delete[] x; x = 0;
01695   delete [] y; y = 0;  
01696 
01697   return median;
01698 
01699 }
01700 //define this as a plug-in
01701 DEFINE_FWK_MODULE(TrackerOfflineValidation);