CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/DQM/TrackingMonitor/src/TrackAnalyzer.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2012/03/28 22:59:43 $
00005  *  $Revision: 1.25 $
00006  *  \author Suchandra Dutta , Giorgia Mila
00007  */
00008 
00009 #include "DataFormats/TrackReco/interface/Track.h"
00010 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00011 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00012 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00013 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00014 #include "MagneticField/Engine/interface/MagneticField.h"
00015 
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 #include "FWCore/Utilities/interface/InputTag.h"
00018 #include "DQMServices/Core/interface/DQMStore.h"
00019 #include "DQM/TrackingMonitor/interface/TrackAnalyzer.h"
00020 #include <string>
00021 #include "TMath.h"
00022 
00023 TrackAnalyzer::TrackAnalyzer(const edm::ParameterSet& iConfig) 
00024     : conf_( iConfig )
00025     , doTrackerSpecific_                   ( conf_.getParameter<bool>("doTrackerSpecific") )
00026     , doAllPlots_                          ( conf_.getParameter<bool>("doAllPlots") )
00027     , doBSPlots_                           ( conf_.getParameter<bool>("doBeamSpotPlots") )
00028     , doGoodTrackPlots_                    ( conf_.getParameter<bool>("doGoodTrackPlots") )
00029     , doDCAPlots_                          ( conf_.getParameter<bool>("doDCAPlots") )
00030     , doGeneralPropertiesPlots_            ( conf_.getParameter<bool>("doGeneralPropertiesPlots") )
00031     , doMeasurementStatePlots_             ( conf_.getParameter<bool>("doMeasurementStatePlots") )
00032     , doHitPropertiesPlots_                ( conf_.getParameter<bool>("doHitPropertiesPlots") )
00033     , doRecHitVsPhiVsEtaPerTrack_          ( conf_.getParameter<bool>("doRecHitVsPhiVsEtaPerTrack") )
00034     , doLayersVsPhiVsEtaPerTrack_          ( conf_.getParameter<bool>("doLayersVsPhiVsEtaPerTrack") )
00035     , doGoodTrackRecHitVsPhiVsEtaPerTrack_ ( conf_.getParameter<bool>("doGoodTrackRecHitVsPhiVsEtaPerTrack") )
00036     , doGoodTrackLayersVsPhiVsEtaPerTrack_ ( conf_.getParameter<bool>("doGoodTrackLayersVsPhiVsEtaPerTrack") )
00037     , doGoodTrack2DChi2Plots_              ( conf_.getParameter<bool>("doGoodTrack2DChi2Plots") )
00038     , doThetaPlots_                        ( conf_.getParameter<bool>("doThetaPlots") )
00039     , doTrackPxPyPlots_                    ( conf_.getParameter<bool>("doTrackPxPyPlots") )
00040     , doDCAwrt000Plots_                    ( conf_.getParameter<bool>("doDCAwrt000Plots") )
00041     , doLumiAnalysis_                      ( conf_.getParameter<bool>("doLumiAnalysis") )
00042     , doTestPlots_                         ( conf_.getParameter<bool>("doTestPlots") )
00043     , NumberOfRecHitsPerTrack(NULL)
00044     , NumberOfRecHitsFoundPerTrack(NULL)
00045     , NumberOfRecHitsLostPerTrack(NULL)
00046     , NumberOfLayersPerTrack(NULL)
00047     , NumberOfRecHitVsPhiVsEtaPerTrack(NULL)
00048     , NumberOfLayersVsPhiVsEtaPerTrack(NULL)
00049     , Chi2(NULL)
00050     , Chi2Prob(NULL)
00051     , Chi2oNDF(NULL)
00052     , DistanceOfClosestApproach(NULL)
00053     , DistanceOfClosestApproachToBS(NULL)
00054     , DistanceOfClosestApproachVsTheta(NULL)
00055     , DistanceOfClosestApproachVsPhi(NULL)
00056     , DistanceOfClosestApproachToBSVsPhi(NULL)
00057     , DistanceOfClosestApproachVsEta(NULL)
00058     , xPointOfClosestApproach(NULL)
00059     , xPointOfClosestApproachVsZ0wrt000(NULL)
00060     , xPointOfClosestApproachVsZ0wrtBS(NULL)
00061     , yPointOfClosestApproach(NULL)
00062     , yPointOfClosestApproachVsZ0wrt000(NULL)
00063     , yPointOfClosestApproachVsZ0wrtBS(NULL)
00064     , zPointOfClosestApproach(NULL)
00065     , zPointOfClosestApproachVsPhi(NULL)
00066     , algorithm(NULL)
00067      // TESTING MEs
00068     , TESTDistanceOfClosestApproachToBS(NULL)
00069     , TESTDistanceOfClosestApproachToBSVsPhi(NULL)
00070                             // add by mia in order to deal w/ LS transitions
00071     , Chi2oNDF_lumiFlag(NULL)
00072     , NumberOfRecHitsPerTrack_lumiFlag(NULL)
00073     , GoodTrackChi2oNDF_lumiFlag(NULL)
00074     , GoodTrackNumberOfRecHitsPerTrack_lumiFlag(NULL)
00075 
00076     , NumberOfTOBRecHitsPerTrack(NULL)
00077     , NumberOfTOBRecHitsPerTrackVsPhiProfile(NULL)
00078     , NumberOfTOBRecHitsPerTrackVsEtaProfile(NULL)
00079     , NumberOfTOBLayersPerTrack(NULL)
00080     , NumberOfTOBLayersPerTrackVsPhiProfile(NULL)
00081     , NumberOfTOBLayersPerTrackVsEtaProfile(NULL)
00082 
00083     , NumberOfTIBRecHitsPerTrack(NULL)
00084     , NumberOfTIBRecHitsPerTrackVsPhiProfile(NULL)
00085     , NumberOfTIBRecHitsPerTrackVsEtaProfile(NULL)
00086     , NumberOfTIBLayersPerTrack(NULL)
00087     , NumberOfTIBLayersPerTrackVsPhiProfile(NULL)
00088     , NumberOfTIBLayersPerTrackVsEtaProfile(NULL)
00089 
00090     , NumberOfTIDRecHitsPerTrack(NULL)
00091     , NumberOfTIDRecHitsPerTrackVsPhiProfile(NULL)
00092     , NumberOfTIDRecHitsPerTrackVsEtaProfile(NULL)
00093     , NumberOfTIDLayersPerTrack(NULL)
00094     , NumberOfTIDLayersPerTrackVsPhiProfile(NULL)
00095     , NumberOfTIDLayersPerTrackVsEtaProfile(NULL)
00096 
00097     , NumberOfTECRecHitsPerTrack(NULL)
00098     , NumberOfTECRecHitsPerTrackVsPhiProfile(NULL)
00099     , NumberOfTECRecHitsPerTrackVsEtaProfile(NULL)
00100     , NumberOfTECLayersPerTrack(NULL)
00101     , NumberOfTECLayersPerTrackVsPhiProfile(NULL)
00102     , NumberOfTECLayersPerTrackVsEtaProfile(NULL)
00103 
00104     , NumberOfPixBarrelRecHitsPerTrack(NULL)
00105     , NumberOfPixBarrelRecHitsPerTrackVsPhiProfile(NULL)
00106     , NumberOfPixBarrelRecHitsPerTrackVsEtaProfile(NULL)
00107     , NumberOfPixBarrelLayersPerTrack(NULL)
00108     , NumberOfPixBarrelLayersPerTrackVsPhiProfile(NULL)
00109     , NumberOfPixBarrelLayersPerTrackVsEtaProfile(NULL)
00110 
00111     , NumberOfPixEndcapRecHitsPerTrack(NULL)
00112     , NumberOfPixEndcapRecHitsPerTrackVsPhiProfile(NULL)
00113     , NumberOfPixEndcapRecHitsPerTrackVsEtaProfile(NULL)
00114     , NumberOfPixEndcapLayersPerTrack(NULL)
00115     , NumberOfPixEndcapLayersPerTrackVsPhiProfile(NULL)
00116     , NumberOfPixEndcapLayersPerTrackVsEtaProfile(NULL)
00117 
00118     , GoodTrackNumberOfRecHitVsPhiVsEtaPerTrack(NULL)
00119     , GoodTrackNumberOfLayersVsPhiVsEtaPerTrack(NULL)
00120     , GoodTrackNumberOfRecHitsPerTrackVsPhiProfile(NULL)
00121     , GoodTrackNumberOfRecHitsPerTrackVsEtaProfile(NULL)
00122     , GoodTrackNumberOfFoundRecHitsPerTrackVsPhiProfile(NULL)
00123     , GoodTrackNumberOfFoundRecHitsPerTrackVsEtaProfile(NULL)
00124     , GoodTrackChi2oNDF(NULL)
00125     , GoodTrackChi2Prob(NULL)
00126     , GoodTrackChi2oNDFVsPhi(NULL)
00127     , GoodTrackChi2ProbVsPhi(NULL)
00128     , GoodTrackChi2oNDFVsEta(NULL)
00129     , GoodTrackChi2ProbVsEta(NULL)
00130     , GoodTrackNumberOfRecHitsPerTrack(NULL)
00131     , GoodTrackNumberOfFoundRecHitsPerTrack(NULL)
00132     , GoodTrackAlgorithm(NULL)
00133 {
00134 
00135   //  std::cout << "TrackAnalyzer::TrackAnalyzer() - doGoodTrackPlots_ = "  << doGoodTrackPlots_ << std::endl;
00136 
00137 }
00138 
00139 TrackAnalyzer::~TrackAnalyzer() 
00140 { 
00141 }
00142 
00143 void TrackAnalyzer::beginJob(DQMStore * dqmStore_) 
00144 {
00145 
00146     // parameters from the configuration
00147     std::string QualName       = conf_.getParameter<std::string>("Quality");
00148     std::string AlgoName       = conf_.getParameter<std::string>("AlgoName");
00149     std::string MEFolderName   = conf_.getParameter<std::string>("FolderName"); 
00150     std::string MEBSFolderName = conf_.getParameter<std::string>("BSFolderName"); 
00151 
00152     // use the AlgoName and Quality Name 
00153     std::string CatagoryName = QualName != "" ? AlgoName + "_" + QualName : AlgoName;
00154 
00155     // get binning from the configuration
00156     double RecHitMin   = conf_.getParameter<double>("RecHitMin");
00157     double RecHitMax   = conf_.getParameter<double>("RecHitMax");
00158 
00159     int    TKHitBin     = conf_.getParameter<int>(   "RecHitBin");
00160     double TKHitMin     = conf_.getParameter<double>("RecHitMin");
00161     double TKHitMax     = conf_.getParameter<double>("RecHitMax");
00162 
00163     int    TKLostBin    = conf_.getParameter<int>(   "RecLostBin");
00164     double TKLostMin    = conf_.getParameter<double>("RecLostMin");
00165     double TKLostMax    = conf_.getParameter<double>("RecLostMax");
00166 
00167     int    TKLayBin     = conf_.getParameter<int>(   "RecLayBin");
00168     double TKLayMin     = conf_.getParameter<double>("RecLayMin");
00169     double TKLayMax     = conf_.getParameter<double>("RecLayMax");
00170 
00171     int    Chi2Bin      = conf_.getParameter<int>(   "Chi2Bin");
00172     double Chi2Min      = conf_.getParameter<double>("Chi2Min");
00173     double Chi2Max      = conf_.getParameter<double>("Chi2Max");
00174 
00175     int    Chi2NDFBin   = conf_.getParameter<int>(   "Chi2NDFBin");
00176     double Chi2NDFMin   = conf_.getParameter<double>("Chi2NDFMin");
00177     double Chi2NDFMax   = conf_.getParameter<double>("Chi2NDFMax");
00178 
00179     int    Chi2ProbBin  = conf_.getParameter<int>(   "Chi2ProbBin");
00180     double Chi2ProbMin  = conf_.getParameter<double>("Chi2ProbMin");
00181     double Chi2ProbMax  = conf_.getParameter<double>("Chi2ProbMax");
00182 
00183     int    PhiBin       = conf_.getParameter<int>(   "PhiBin");
00184     double PhiMin       = conf_.getParameter<double>("PhiMin");
00185     double PhiMax       = conf_.getParameter<double>("PhiMax");
00186 
00187     int    EtaBin       = conf_.getParameter<int>(   "EtaBin");
00188     double EtaMin       = conf_.getParameter<double>("EtaMin");
00189     double EtaMax       = conf_.getParameter<double>("EtaMax");
00190 
00191     int    ThetaBin     = conf_.getParameter<int>(   "ThetaBin");
00192     double ThetaMin     = conf_.getParameter<double>("ThetaMin");
00193     double ThetaMax     = conf_.getParameter<double>("ThetaMax");
00194 
00195     int    DxyBin       = conf_.getParameter<int>(   "DxyBin");
00196     double DxyMin       = conf_.getParameter<double>("DxyMin");
00197     double DxyMax       = conf_.getParameter<double>("DxyMax");
00198 
00199     int    VXBin        = conf_.getParameter<int>(   "VXBin");
00200     double VXMin        = conf_.getParameter<double>("VXMin");
00201     double VXMax        = conf_.getParameter<double>("VXMax");
00202 
00203     int    VYBin        = conf_.getParameter<int>(   "VYBin");
00204     double VYMin        = conf_.getParameter<double>("VYMin");
00205     double VYMax        = conf_.getParameter<double>("VYMax");
00206 
00207     int    VZBin        = conf_.getParameter<int>(   "VZBin");
00208     double VZMin        = conf_.getParameter<double>("VZMin");
00209     double VZMax        = conf_.getParameter<double>("VZMax");
00210 
00211     int    X0Bin        = conf_.getParameter<int>(   "X0Bin");
00212     double X0Min        = conf_.getParameter<double>("X0Min");
00213     double X0Max        = conf_.getParameter<double>("X0Max");
00214 
00215     int    Y0Bin        = conf_.getParameter<int>(   "Y0Bin");
00216     double Y0Min        = conf_.getParameter<double>("Y0Min");
00217     double Y0Max        = conf_.getParameter<double>("Y0Max");
00218 
00219     int    Z0Bin        = conf_.getParameter<int>(   "Z0Bin");
00220     double Z0Min        = conf_.getParameter<double>("Z0Min");
00221     double Z0Max        = conf_.getParameter<double>("Z0Max");
00222 
00223     dqmStore_->setCurrentFolder(MEFolderName);
00224 
00225     // book the Hit Property histograms
00226     // ---------------------------------------------------------------------------------//
00227 
00228     if ( doHitPropertiesPlots_ || doAllPlots_ ){
00229 
00230       dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties");
00231       
00232       histname = "NumberOfRecHitsPerTrack_";
00233       NumberOfRecHitsPerTrack = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, TKHitBin, TKHitMin, TKHitMax);
00234       NumberOfRecHitsPerTrack->setAxisTitle("Number of all RecHits of each Track");
00235       NumberOfRecHitsPerTrack->setAxisTitle("Number of Tracks", 2);
00236 
00237       histname = "NumberOfRecHitsFoundPerTrack_";
00238       NumberOfRecHitsFoundPerTrack = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, TKHitBin, TKHitMin, TKHitMax);
00239       NumberOfRecHitsFoundPerTrack->setAxisTitle("Number of RecHits found for each Track");
00240       NumberOfRecHitsFoundPerTrack->setAxisTitle("Number of Tracks", 2);
00241 
00242       histname = "NumberOfRecHitsLostPerTrack_";
00243       NumberOfRecHitsLostPerTrack = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, TKLostBin, TKLostMin, TKLostMax);
00244       NumberOfRecHitsLostPerTrack->setAxisTitle("Number of RecHits lost for each Track");
00245       NumberOfRecHitsLostPerTrack->setAxisTitle("Number of Tracks", 2);
00246 
00247       histname = "NumberOfLayersPerTrack_";
00248       NumberOfLayersPerTrack = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, TKLayBin, TKLayMin, TKLayMax);
00249       NumberOfLayersPerTrack->setAxisTitle("Number of Layers of each Track", 1);
00250       NumberOfLayersPerTrack->setAxisTitle("Number of Tracks", 2);
00251       
00252       if ( doRecHitVsPhiVsEtaPerTrack_ || doAllPlots_ ){
00253         
00254         histname = "NumberOfRecHitVsPhiVsEtaPerTrack_";
00255         NumberOfRecHitVsPhiVsEtaPerTrack = dqmStore_->bookProfile2D(histname+CatagoryName, histname+CatagoryName, 
00256                                                                     EtaBin, EtaMin, EtaMax, PhiBin, PhiMin, PhiMax, 0, 40., "");
00257         NumberOfRecHitVsPhiVsEtaPerTrack->setAxisTitle("Track #eta ", 1);
00258         NumberOfRecHitVsPhiVsEtaPerTrack->setAxisTitle("Track #phi ", 2);
00259       }
00260       if ( doLayersVsPhiVsEtaPerTrack_ || doAllPlots_ ){
00261         
00262         histname = "NumberOfLayersVsPhiVsEtaPerTrack_";
00263         NumberOfLayersVsPhiVsEtaPerTrack = dqmStore_->bookProfile2D(histname+CatagoryName, histname+CatagoryName, 
00264                                                                     EtaBin, EtaMin, EtaMax, PhiBin, PhiMin, PhiMax, 0, 40., "");
00265         NumberOfLayersVsPhiVsEtaPerTrack->setAxisTitle("Track #eta ", 1);
00266         NumberOfLayersVsPhiVsEtaPerTrack->setAxisTitle("Track #phi ", 2);
00267       }
00268     }
00269 
00270 
00271     // book the General Property histograms
00272     // ---------------------------------------------------------------------------------//
00273 
00274     if (doGeneralPropertiesPlots_ || doAllPlots_){
00275 
00276       dqmStore_->setCurrentFolder(MEFolderName+"/GeneralProperties");
00277 
00278       histname = "Chi2_";
00279       Chi2 = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, Chi2Bin, Chi2Min, Chi2Max);
00280       Chi2->setAxisTitle("Track #chi^{2}"  ,1);
00281       Chi2->setAxisTitle("Number of Tracks",2);
00282 
00283       histname = "Chi2Prob_";
00284       Chi2Prob = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, Chi2ProbBin, Chi2ProbMin, Chi2ProbMax);
00285       Chi2Prob->setAxisTitle("Track #chi^{2} probability",1);
00286       Chi2Prob->setAxisTitle("Number of Tracks"        ,2);
00287       
00288       histname = "Chi2oNDF_";
00289       Chi2oNDF = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, Chi2NDFBin, Chi2NDFMin, Chi2NDFMax);
00290       Chi2oNDF->setAxisTitle("Track #chi^{2}/ndf",1);
00291       Chi2oNDF->setAxisTitle("Number of Tracks"  ,2);
00292       
00293       histname = "xPointOfClosestApproach_";
00294       xPointOfClosestApproach = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, VXBin, VXMin, VXMax);
00295       xPointOfClosestApproach->setAxisTitle("x component of Track PCA to beam line (cm)",1);
00296       xPointOfClosestApproach->setAxisTitle("Number of Tracks",2);
00297       
00298       histname = "yPointOfClosestApproach_";
00299       yPointOfClosestApproach = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, VYBin, VYMin, VYMax);
00300       yPointOfClosestApproach->setAxisTitle("y component of Track PCA to beam line (cm)",1);
00301       yPointOfClosestApproach->setAxisTitle("Number of Tracks",2);
00302 
00303       histname = "zPointOfClosestApproach_";
00304       zPointOfClosestApproach = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, VZBin, VZMin, VZMax);
00305       zPointOfClosestApproach->setAxisTitle("z component of Track PCA to beam line (cm)",1);
00306       zPointOfClosestApproach->setAxisTitle("Number of Tracks",2);
00307       
00308       // See DataFormats/TrackReco/interface/TrackBase.h for track algorithm enum definition
00309       // http://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/DataFormats/TrackReco/interface/TrackBase.h?view=log
00310       histname = "algorithm_";
00311       algorithm = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, 32, 0., 32.);
00312       algorithm->setAxisTitle("Tracking algorithm",1);
00313       algorithm->setAxisTitle("Number of Tracks",2);
00314       
00315     }
00316     
00317     // book LS analysis related histograms
00318     // -----------------------------------
00319     if ( doLumiAnalysis_ ) {
00320       // add by Mia in order to deal w/ LS transitions  
00321       dqmStore_->setCurrentFolder(MEFolderName+"/LSanalysis");
00322 
00323       histname = "NumberOfRecHitsPerTrack_lumiFlag_";
00324       NumberOfRecHitsPerTrack_lumiFlag = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, TKHitBin, TKHitMin, TKHitMax);
00325       NumberOfRecHitsPerTrack_lumiFlag->setAxisTitle("Number of all RecHits of each Track");
00326       NumberOfRecHitsPerTrack_lumiFlag->setAxisTitle("Number of Tracks", 2);
00327 
00328       histname = "Chi2oNDF_lumiFlag_";
00329       Chi2oNDF_lumiFlag = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, Chi2NDFBin, Chi2NDFMin, Chi2NDFMax);
00330       Chi2oNDF_lumiFlag->setAxisTitle("Track #chi^{2}/ndf",1);
00331       Chi2oNDF_lumiFlag->setAxisTitle("Number of Tracks"  ,2);
00332 
00333       histname = "GoodTrackNumberOfRecHitsPerTrack_lumiFlag_";
00334       GoodTrackNumberOfRecHitsPerTrack_lumiFlag = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, TKHitBin, TKHitMin, TKHitMax);
00335       GoodTrackNumberOfRecHitsPerTrack_lumiFlag->setAxisTitle("Number of all RecHits of each Good Track");
00336       GoodTrackNumberOfRecHitsPerTrack_lumiFlag->setAxisTitle("Number of Good Tracks", 2);
00337       
00338       histname = "GoodTrackChi2oNDF_lumiFlag_";
00339       GoodTrackChi2oNDF_lumiFlag = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, Chi2NDFBin, Chi2NDFMin, Chi2NDFMax);
00340       GoodTrackChi2oNDF_lumiFlag->setAxisTitle("Good Track #chi^{2}/ndf",1);
00341       GoodTrackChi2oNDF_lumiFlag->setAxisTitle("Number of Good Tracks"  ,2);
00342     }
00343 
00344 
00345     // book the Beam Spot related histograms
00346     // ---------------------------------------------------------------------------------//
00347     
00348     if(doBSPlots_ || doAllPlots_)
00349       {
00350         //        dqmStore_->setCurrentFolder(MEBSFolderName);
00351         dqmStore_->setCurrentFolder(MEFolderName+"/GeneralProperties");
00352 
00353         histname = "DistanceOfClosestApproachToBS_";
00354         DistanceOfClosestApproachToBS = dqmStore_->book1D(histname+CatagoryName,histname+CatagoryName,DxyBin,DxyMin,DxyMax);
00355         DistanceOfClosestApproachToBS->setAxisTitle("Track d_{xy} wrt beam spot (cm)",1);
00356         DistanceOfClosestApproachToBS->setAxisTitle("Number of Tracks",2);
00357         
00358         histname = "DistanceOfClosestApproachToBSVsPhi_";
00359         DistanceOfClosestApproachToBSVsPhi = dqmStore_->bookProfile(histname+CatagoryName,histname+CatagoryName, PhiBin, PhiMin, PhiMax, DxyBin, DxyMin, DxyMax,"");
00360         DistanceOfClosestApproachToBSVsPhi->getTH1()->SetBit(TH1::kCanRebin);
00361         DistanceOfClosestApproachToBSVsPhi->setAxisTitle("Track #phi",1);
00362         DistanceOfClosestApproachToBSVsPhi->setAxisTitle("Track d_{xy} wrt beam spot (cm)",2);
00363         
00364         histname = "xPointOfClosestApproachVsZ0wrt000_";
00365         xPointOfClosestApproachVsZ0wrt000 = dqmStore_->bookProfile(histname+CatagoryName, histname+CatagoryName, Z0Bin, Z0Min, Z0Max, X0Bin, X0Min, X0Max,"");
00366         xPointOfClosestApproachVsZ0wrt000->setAxisTitle("d_{z} (cm)",1);
00367         xPointOfClosestApproachVsZ0wrt000->setAxisTitle("x component of Track PCA to beam line (cm)",2);
00368         
00369         histname = "yPointOfClosestApproachVsZ0wrt000_";
00370         yPointOfClosestApproachVsZ0wrt000 = dqmStore_->bookProfile(histname+CatagoryName, histname+CatagoryName, Z0Bin, Z0Min, Z0Max, Y0Bin, Y0Min, Y0Max,"");
00371         yPointOfClosestApproachVsZ0wrt000->setAxisTitle("d_{z} (cm)",1);
00372         yPointOfClosestApproachVsZ0wrt000->setAxisTitle("y component of Track PCA to beam line (cm)",2);
00373 
00374         histname = "xPointOfClosestApproachVsZ0wrtBS_";
00375         xPointOfClosestApproachVsZ0wrtBS = dqmStore_->bookProfile(histname+CatagoryName, histname+CatagoryName, Z0Bin, Z0Min, Z0Max, X0Bin, X0Min, X0Max,"");
00376         xPointOfClosestApproachVsZ0wrtBS->setAxisTitle("d_{z} w.r.t. Beam Spot  (cm)",1);
00377         xPointOfClosestApproachVsZ0wrtBS->setAxisTitle("x component of Track PCA to BS (cm)",2);
00378         
00379         histname = "yPointOfClosestApproachVsZ0wrtBS_";
00380         yPointOfClosestApproachVsZ0wrtBS = dqmStore_->bookProfile(histname+CatagoryName, histname+CatagoryName, Z0Bin, Z0Min, Z0Max, Y0Bin, Y0Min, Y0Max,"");
00381         yPointOfClosestApproachVsZ0wrtBS->setAxisTitle("d_{z} w.r.t. Beam Spot (cm)",1);
00382         yPointOfClosestApproachVsZ0wrtBS->setAxisTitle("y component of Track PCA to BS (cm)",2);
00383 
00384         histname = "zPointOfClosestApproachVsPhi_";
00385         zPointOfClosestApproachVsPhi = dqmStore_->bookProfile(histname+CatagoryName, histname+CatagoryName, PhiBin, PhiMin, PhiMax, VZBin, VZMin, VZMax, "");
00386         zPointOfClosestApproachVsPhi->setAxisTitle("Track #phi",1);
00387         zPointOfClosestApproachVsPhi->setAxisTitle("y component of Track PCA to beam line (cm)",2);
00388 
00389         if (doTestPlots_) {
00390           
00391           histname = "TESTDistanceOfClosestApproachToBS_";
00392           TESTDistanceOfClosestApproachToBS = dqmStore_->book1D(histname+CatagoryName,histname+CatagoryName,DxyBin,DxyMin,DxyMax);
00393           TESTDistanceOfClosestApproachToBS->setAxisTitle("Track d_{xy} wrt beam spot (cm)",1);
00394           TESTDistanceOfClosestApproachToBS->setAxisTitle("Number of Tracks",2);
00395           
00396           histname = "TESTDistanceOfClosestApproachToBSVsPhi_";
00397           TESTDistanceOfClosestApproachToBSVsPhi = dqmStore_->bookProfile(histname+CatagoryName,histname+CatagoryName, PhiBin, PhiMin, PhiMax, DxyBin, DxyMin, DxyMax,"");
00398           TESTDistanceOfClosestApproachToBSVsPhi->getTH1()->SetBit(TH1::kCanRebin);
00399           TESTDistanceOfClosestApproachToBSVsPhi->setAxisTitle("Track #phi",1);
00400           TESTDistanceOfClosestApproachToBSVsPhi->setAxisTitle("Track d_{xy} wrt beam spot (cm)",2);
00401 
00402         }
00403 
00404       }
00405     
00406     // book the Profile plots for DCA related histograms
00407     // ---------------------------------------------------------------------------------//
00408     if(doDCAPlots_ || doAllPlots_)
00409       {
00410         if (doDCAwrt000Plots_) {
00411           if (doThetaPlots_) {
00412             dqmStore_->setCurrentFolder(MEFolderName+"/GeneralProperties");
00413             histname = "DistanceOfClosestApproachVsTheta_";
00414             DistanceOfClosestApproachVsTheta = dqmStore_->bookProfile(histname+CatagoryName,histname+CatagoryName, ThetaBin, ThetaMin, ThetaMax, DxyMin,DxyMax,"");
00415             DistanceOfClosestApproachVsTheta->setAxisTitle("Track #theta",1);
00416             DistanceOfClosestApproachVsTheta->setAxisTitle("Track d_{xy} wrt (0,0,0) (cm)",2);
00417           }
00418           
00419           histname = "DistanceOfClosestApproachVsEta_";
00420           DistanceOfClosestApproachVsEta = dqmStore_->bookProfile(histname+CatagoryName,histname+CatagoryName, EtaBin, EtaMin, EtaMax, DxyMin, DxyMax,"");
00421           DistanceOfClosestApproachVsEta->setAxisTitle("Track #eta",1);
00422           DistanceOfClosestApproachVsEta->setAxisTitle("Track d_{xy} wrt (0,0,0) (cm)",2);
00423           // temporary patch in order to put back those MEs in Muon Workspace
00424 
00425           histname = "DistanceOfClosestApproach_";
00426           DistanceOfClosestApproach = dqmStore_->book1D(histname+CatagoryName,histname+CatagoryName,DxyBin,DxyMin,DxyMax);
00427           DistanceOfClosestApproach->setAxisTitle("Track d_{xy} wrt (0,0,0) (cm)",1);
00428           DistanceOfClosestApproach->setAxisTitle("Number of Tracks",2);
00429           
00430           histname = "DistanceOfClosestApproachVsPhi_";
00431           DistanceOfClosestApproachVsPhi = dqmStore_->bookProfile(histname+CatagoryName,histname+CatagoryName, PhiBin, PhiMin, PhiMax, DxyMin,DxyMax,"");
00432           DistanceOfClosestApproachVsPhi->getTH1()->SetBit(TH1::kCanRebin);
00433           DistanceOfClosestApproachVsPhi->setAxisTitle("Track #phi",1);
00434           DistanceOfClosestApproachVsPhi->setAxisTitle("Track d_{xy} wrt (0,0,0) (cm)",2);
00435         }
00436       }
00437 
00438 
00439     // book tracker specific related histograms
00440     // ---------------------------------------------------------------------------------//
00441     if(doTrackerSpecific_ || doAllPlots_) 
00442       {
00443         doTrackerSpecificInitialization(dqmStore_);
00444       }
00445     
00446     // book state related histograms
00447     // ---------------------------------------------------------------------------------//
00448     if (doMeasurementStatePlots_ || doAllPlots_){
00449 
00450       std::string StateName = conf_.getParameter<std::string>("MeasurementState");
00451       
00452       if (StateName == "All") 
00453         {
00454           bookHistosForState("OuterSurface", dqmStore_);
00455           bookHistosForState("InnerSurface", dqmStore_);
00456           bookHistosForState("ImpactPoint" , dqmStore_);
00457         } 
00458       else if 
00459         (   
00460          StateName != "OuterSurface" && 
00461          StateName != "InnerSurface" && 
00462          StateName != "ImpactPoint" &&
00463          StateName != "default" 
00464          ) 
00465         {
00466           bookHistosForState("default", dqmStore_);
00467         }
00468       else
00469         {
00470           bookHistosForState(StateName, dqmStore_);
00471         }
00472     }
00473     
00474     // book histos for good tracks (HP + Pt>1GeV)
00475     // ---------------------------------------------------------------------------------//
00476 
00477     if ( doGoodTrackPlots_ || doAllPlots_ ) {
00478 
00479       dqmStore_->setCurrentFolder(MEFolderName+"/GeneralProperties/GoodTracks");
00480 
00481       histname = "GoodTrackChi2oNDF_";
00482       GoodTrackChi2oNDF = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, Chi2NDFBin, Chi2NDFMin, Chi2NDFMax);
00483       GoodTrackChi2oNDF->setAxisTitle("Good Track #chi^{2}/ndf",1);
00484       GoodTrackChi2oNDF->setAxisTitle("Number of Good Tracks"  ,2);
00485 
00486       histname = "GoodTrackChi2Prob_";
00487       GoodTrackChi2Prob = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, Chi2ProbBin, Chi2ProbMin, Chi2ProbMax);
00488       GoodTrackChi2Prob->setAxisTitle("Good Track #chi^{2} probability",1);
00489       GoodTrackChi2Prob->setAxisTitle("Number of Good Tracks"  ,2);
00490 
00491       histname = "GoodTrackChi2oNDFVsPhi_";
00492       GoodTrackChi2oNDFVsPhi = dqmStore_->bookProfile(histname+CatagoryName, histname+CatagoryName, PhiBin, PhiMin, PhiMax, Chi2NDFMin, Chi2NDFMax);
00493       GoodTrackChi2oNDFVsPhi->setAxisTitle("Good Tracks #phi"  ,1);
00494       GoodTrackChi2oNDFVsPhi->setAxisTitle("Good Track #chi^{2}/ndf",2);
00495 
00496       histname = "GoodTrackChi2ProbVsPhi_";
00497       GoodTrackChi2ProbVsPhi = dqmStore_->bookProfile(histname+CatagoryName, histname+CatagoryName, PhiBin, PhiMin, PhiMax, Chi2ProbMin, Chi2ProbMax);
00498       GoodTrackChi2ProbVsPhi->setAxisTitle("Good Tracks #phi"  ,1);
00499       GoodTrackChi2ProbVsPhi->setAxisTitle("Good Track #chi^{2} probability",2);
00500 
00501       histname = "GoodTrackChi2oNDFVsEta_";
00502       GoodTrackChi2oNDFVsEta = dqmStore_->bookProfile(histname+CatagoryName, histname+CatagoryName, EtaBin, EtaMin, EtaMax, Chi2NDFMin, Chi2NDFMax);
00503       GoodTrackChi2oNDFVsEta->setAxisTitle("Good Tracks #eta"  ,1);
00504       GoodTrackChi2oNDFVsEta->setAxisTitle("Good Track #chi^{2}/ndf",2);
00505 
00506       histname = "GoodTrackChi2ProbVsEta_";
00507       GoodTrackChi2ProbVsEta = dqmStore_->bookProfile(histname+CatagoryName, histname+CatagoryName, EtaBin, EtaMin, EtaMax, Chi2ProbMin, Chi2ProbMax);
00508       GoodTrackChi2ProbVsEta->setAxisTitle("Good Tracks #eta"  ,1);
00509       GoodTrackChi2ProbVsEta->setAxisTitle("Good Track #chi^{2} probability",2);
00510 
00511       histname = "GoodTrackAlgorithm_";
00512       GoodTrackAlgorithm = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, 32, 0., 32.);
00513       GoodTrackAlgorithm->setAxisTitle("Good Track tracking algorithm",1);
00514       GoodTrackAlgorithm->setAxisTitle("Number of Good Tracks",2);
00515 
00516 
00517       dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties/GoodTracks");
00518 
00519       histname = "GoodTrackNumberOfRecHitsPerTrack_";
00520       GoodTrackNumberOfRecHitsPerTrack = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, TKHitBin, TKHitMin, TKHitMax);
00521       GoodTrackNumberOfRecHitsPerTrack->setAxisTitle("Number of all RecHits of each Good Track");
00522       GoodTrackNumberOfRecHitsPerTrack->setAxisTitle("Number of Good Tracks", 2);
00523       
00524       histname = "GoodTrackNumberOfRecHitsFoundPerTrack_";
00525       GoodTrackNumberOfFoundRecHitsPerTrack = dqmStore_->book1D(histname+CatagoryName, histname+CatagoryName, TKHitBin, TKHitMin, TKHitMax);
00526       GoodTrackNumberOfFoundRecHitsPerTrack->setAxisTitle("Number of found RecHits of each Good Track");
00527       GoodTrackNumberOfFoundRecHitsPerTrack->setAxisTitle("Number of Good Tracks", 2);
00528       
00529       if ( doGoodTrackRecHitVsPhiVsEtaPerTrack_ || doAllPlots_ ){
00530         
00531         histname = "GoodTrackNumberOfRecHitVsPhiVsEtaPerTrack_";
00532         GoodTrackNumberOfRecHitVsPhiVsEtaPerTrack = dqmStore_->bookProfile2D(histname+CatagoryName, histname+CatagoryName, 
00533                                                                              EtaBin, EtaMin, EtaMax, PhiBin, PhiMin, PhiMax, 0, 40., "");
00534         GoodTrackNumberOfRecHitVsPhiVsEtaPerTrack->setAxisTitle("Good Track #eta ", 1);
00535         GoodTrackNumberOfRecHitVsPhiVsEtaPerTrack->setAxisTitle("Good Track #phi ", 2);
00536       }
00537       
00538       if ( doGoodTrackLayersVsPhiVsEtaPerTrack_ || doAllPlots_ ){
00539         
00540         histname = "GoodTrackNumberOfLayersVsPhiVsEtaPerTrack_";
00541         GoodTrackNumberOfLayersVsPhiVsEtaPerTrack = dqmStore_->bookProfile2D(histname+CatagoryName, histname+CatagoryName, 
00542                                                                              EtaBin, EtaMin, EtaMax, PhiBin, PhiMin, PhiMax, 0, 40., "");
00543         GoodTrackNumberOfLayersVsPhiVsEtaPerTrack->setAxisTitle("Good Track #eta ", 1);
00544         GoodTrackNumberOfLayersVsPhiVsEtaPerTrack->setAxisTitle("Good Track #phi ", 2);
00545       }
00546       // rechits
00547       histname = "GoodTrackNumberOfRecHitsPerTrackVsPhiProfile_";
00548       GoodTrackNumberOfRecHitsPerTrackVsPhiProfile = dqmStore_->bookProfile(histname+CatagoryName, histname+CatagoryName, PhiBin, PhiMin, PhiMax, RecHitMin, RecHitMax,"");
00549       GoodTrackNumberOfRecHitsPerTrackVsPhiProfile->setAxisTitle("Good Track #phi",1);
00550       GoodTrackNumberOfRecHitsPerTrackVsPhiProfile->setAxisTitle("Number of  RecHits of each Track",2);
00551 
00552       histname = "GoodTrackNumberOfRecHitsPerTrackVsEtaProfile_";
00553       GoodTrackNumberOfRecHitsPerTrackVsEtaProfile = dqmStore_->bookProfile(histname+CatagoryName, histname+CatagoryName, EtaBin, EtaMin, EtaMax, RecHitMin, RecHitMax,"");
00554       GoodTrackNumberOfRecHitsPerTrackVsEtaProfile->setAxisTitle("Good Track #eta",1);
00555       GoodTrackNumberOfRecHitsPerTrackVsEtaProfile->setAxisTitle("Number of  RecHits of each Track",2);
00556 
00557       histname = "GoodTrackNumberOfFoundRecHitsPerTrackVsPhiProfile_";
00558       GoodTrackNumberOfFoundRecHitsPerTrackVsPhiProfile = dqmStore_->bookProfile(histname+CatagoryName, histname+CatagoryName, PhiBin, PhiMin, PhiMax, RecHitMin, RecHitMax,"");
00559       GoodTrackNumberOfFoundRecHitsPerTrackVsPhiProfile->setAxisTitle("Good Track #phi",1);
00560       GoodTrackNumberOfFoundRecHitsPerTrackVsPhiProfile->setAxisTitle("Number of found RecHits of each Track",2);
00561 
00562       histname = "GoodTrackNumberOfFoundRecHitsPerTrackVsEtaProfile_";
00563       GoodTrackNumberOfFoundRecHitsPerTrackVsEtaProfile = dqmStore_->bookProfile(histname+CatagoryName, histname+CatagoryName, EtaBin, EtaMin, EtaMax, RecHitMin, RecHitMax,"");
00564       GoodTrackNumberOfFoundRecHitsPerTrackVsEtaProfile->setAxisTitle("Good Track #eta",1);
00565       GoodTrackNumberOfFoundRecHitsPerTrackVsEtaProfile->setAxisTitle("Number of found RecHits of each Track",2);
00566 
00567     }
00568    
00569 }
00570 
00571 // -- Analyse
00572 // ---------------------------------------------------------------------------------//
00573 void TrackAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup, const reco::Track& track)
00574 {
00575 
00576   if ( doHitPropertiesPlots_ || doAllPlots_ ){
00577     // rec hits
00578     // originally it was track.recHitsSize(), but it is exploiting extra
00579     // therefore it is moved to track.hitPattern().numberOfHits()
00580     NumberOfRecHitsPerTrack->Fill(track.hitPattern().numberOfHits());
00581     NumberOfRecHitsFoundPerTrack->Fill(track.numberOfValidHits());
00582     NumberOfRecHitsLostPerTrack->Fill(track.numberOfLostHits());
00583 
00584     // 2D plots    
00585     if ( doRecHitVsPhiVsEtaPerTrack_ || doAllPlots_ )
00586       NumberOfRecHitVsPhiVsEtaPerTrack->Fill(track.eta(),track.phi(),track.hitPattern().numberOfHits());    
00587 
00588     // layers
00589     NumberOfLayersPerTrack->Fill(track.hitPattern().trackerLayersWithMeasurement());
00590     // 2D plots    
00591     if ( doLayersVsPhiVsEtaPerTrack_ || doAllPlots_ )
00592       NumberOfLayersVsPhiVsEtaPerTrack->Fill(track.eta(),track.phi(),track.hitPattern().trackerLayersWithMeasurement());
00593 
00594 
00595   }
00596 
00597   if (doGeneralPropertiesPlots_ || doAllPlots_){
00598     // fitting
00599     Chi2->Fill(track.chi2());
00600     Chi2Prob->Fill(TMath::Prob(track.chi2(),(int)track.ndof()));
00601     Chi2oNDF->Fill(track.normalizedChi2());
00602 
00603     // DCA
00604     // temporary patch in order to put back those MEs in Muon Workspace 
00605     if (doDCAwrt000Plots_) {
00606       DistanceOfClosestApproach->Fill(track.dxy());
00607       DistanceOfClosestApproachVsPhi->Fill(track.phi(), track.dxy());
00608     }
00609     /*
00610       http://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/DataFormats/TrackReco/interface/TrackBase.h?view=markup
00611       vertex() is DEPRECATED !!
00612     xPointOfClosestApproach->Fill(track.vertex().x());
00613     yPointOfClosestApproach->Fill(track.vertex().y());
00614     zPointOfClosestApproach->Fill(track.vertex().z());
00615     */
00616 
00617     // PCA
00618     xPointOfClosestApproach->Fill(track.referencePoint().x());
00619     yPointOfClosestApproach->Fill(track.referencePoint().y());
00620     zPointOfClosestApproach->Fill(track.referencePoint().z());
00621 
00622     // algorithm
00623     algorithm->Fill(static_cast<double>(track.algo()));
00624 
00625   }
00626 
00627   if ( doLumiAnalysis_ ) {
00628     NumberOfRecHitsPerTrack_lumiFlag -> Fill(track.hitPattern().numberOfHits());    
00629     Chi2oNDF_lumiFlag->Fill(track.normalizedChi2());
00630   }
00631 
00632   if(doBSPlots_ || doAllPlots_)
00633     {
00634       edm::InputTag bsSrc = conf_.getParameter< edm::InputTag >("beamSpot");
00635 
00636       edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00637         iEvent.getByLabel(bsSrc,recoBeamSpotHandle);
00638         reco::BeamSpot bs = *recoBeamSpotHandle;      
00639 
00640         DistanceOfClosestApproachToBS->Fill(track.dxy(bs.position()));
00641         DistanceOfClosestApproachToBSVsPhi->Fill(track.phi(), track.dxy(bs.position()));
00642         zPointOfClosestApproachVsPhi->Fill(track.phi(),track.vz());
00643         xPointOfClosestApproachVsZ0wrt000->Fill(track.dz(),track.vx());
00644         yPointOfClosestApproachVsZ0wrt000->Fill(track.dz(),track.vy());
00645         xPointOfClosestApproachVsZ0wrtBS->Fill(track.dz(bs.position()),(track.vx()-bs.position(track.vz()).x()));
00646         yPointOfClosestApproachVsZ0wrtBS->Fill(track.dz(bs.position()),(track.vy()-bs.position(track.vz()).y()));
00647         if (doTestPlots_) {
00648           TESTDistanceOfClosestApproachToBS->Fill(track.dxy(bs.position(track.vz())));
00649           TESTDistanceOfClosestApproachToBSVsPhi->Fill(track.phi(), track.dxy(bs.position(track.vz())));
00650         }
00651     }
00652 
00653     if(doDCAPlots_ || doAllPlots_)
00654       {
00655         if (doDCAwrt000Plots_) {
00656           if (doThetaPlots_) {
00657             DistanceOfClosestApproachVsTheta->Fill(track.theta(), track.d0());
00658           }
00659           DistanceOfClosestApproachVsEta->Fill(track.eta(), track.d0());
00660         }  
00661       }
00662     //Tracker Specific Histograms
00663     if(doTrackerSpecific_ || doAllPlots_) 
00664       {
00665         doTrackerSpecificFillHists(track);
00666       }
00667 
00668     if (doMeasurementStatePlots_ || doAllPlots_){
00669       std::string StateName = conf_.getParameter<std::string>("MeasurementState");
00670       if (StateName == "All") 
00671         {
00672           fillHistosForState(iSetup, track, std::string("OuterSurface"));
00673           fillHistosForState(iSetup, track, std::string("InnerSurface"));
00674           fillHistosForState(iSetup, track, std::string("ImpactPoint"));
00675         } 
00676       else if 
00677         (   
00678          StateName != "OuterSurface" && 
00679          StateName != "InnerSurface" && 
00680          StateName != "ImpactPoint" &&
00681          StateName != "default" 
00682          ) 
00683         {
00684           fillHistosForState(iSetup, track, std::string("default"));
00685         }
00686       else
00687         {
00688           fillHistosForState(iSetup, track, StateName);
00689         }
00690     }
00691 
00692     // Good Tracks plots
00693     if ( track.quality(reco::TrackBase::highPurity) && track.pt() > 1. ) {
00694       if ( doGoodTrackPlots_ || doAllPlots_ ) {
00695         GoodTrackChi2oNDF->Fill(track.normalizedChi2());
00696         GoodTrackChi2Prob->Fill(TMath::Prob(track.chi2(),(int)track.ndof()));
00697         GoodTrackChi2oNDFVsPhi->Fill(track.phi(),track.normalizedChi2());
00698         GoodTrackChi2ProbVsPhi->Fill(track.phi(),TMath::Prob(track.chi2(),(int)track.ndof()));
00699         GoodTrackChi2oNDFVsEta->Fill(track.eta(),track.normalizedChi2());
00700         GoodTrackChi2ProbVsEta->Fill(track.eta(),TMath::Prob(track.chi2(),(int)track.ndof()));
00701 
00702         // originally it was track.recHitsSize(), but it is exploiting extra
00703         // therefore it is moved to track.hitPattern().numberOfHits()
00704         GoodTrackNumberOfRecHitsPerTrack->Fill(track.hitPattern().numberOfHits());
00705         GoodTrackNumberOfFoundRecHitsPerTrack->Fill(track.numberOfValidHits());
00706         // algorithm
00707         GoodTrackAlgorithm->Fill(static_cast<double>(track.algo()));
00708         GoodTrackNumberOfFoundRecHitsPerTrackVsPhiProfile->Fill(track.phi(),track.numberOfValidHits());
00709         GoodTrackNumberOfFoundRecHitsPerTrackVsEtaProfile->Fill(track.eta(),track.numberOfValidHits());
00710         GoodTrackNumberOfRecHitsPerTrackVsPhiProfile->Fill(track.phi(),track.hitPattern().numberOfHits());
00711         GoodTrackNumberOfRecHitsPerTrackVsEtaProfile->Fill(track.eta(),track.hitPattern().numberOfHits());
00712 
00713         if ( doGoodTrackRecHitVsPhiVsEtaPerTrack_ || doAllPlots_ ) 
00714           GoodTrackNumberOfRecHitVsPhiVsEtaPerTrack->Fill(track.eta(),track.phi(),track.hitPattern().numberOfHits());
00715 
00716         if ( doGoodTrackLayersVsPhiVsEtaPerTrack_ || doAllPlots_ ) 
00717           GoodTrackNumberOfLayersVsPhiVsEtaPerTrack->Fill(track.eta(),track.phi(),track.hitPattern().trackerLayersWithMeasurement());
00718       }
00719       if ( doLumiAnalysis_ ) {
00720         GoodTrackChi2oNDF_lumiFlag                -> Fill(track.normalizedChi2()); 
00721         GoodTrackNumberOfRecHitsPerTrack_lumiFlag -> Fill(track.hitPattern().numberOfHits());
00722       }
00723     }
00724 
00725 }
00726 
00727 
00728 // book histograms at differnt measurement points
00729 // ---------------------------------------------------------------------------------//
00730 void TrackAnalyzer::bookHistosForState(std::string sname, DQMStore * dqmStore_) 
00731 {
00732 
00733 
00734     // parameters from the configuration
00735     std::string QualName       = conf_.getParameter<std::string>("Quality");
00736     std::string AlgoName       = conf_.getParameter<std::string>("AlgoName");
00737     std::string MEFolderName   = conf_.getParameter<std::string>("FolderName"); 
00738 
00739     // use the AlgoName and Quality Name 
00740     std::string CatagoryName = QualName != "" ? AlgoName + "_" + QualName : AlgoName;
00741 
00742     // get binning from the configuration
00743     double Chi2NDFMin = conf_.getParameter<double>("Chi2NDFMin");
00744     double Chi2NDFMax = conf_.getParameter<double>("Chi2NDFMax");
00745 
00746     int    RecHitBin   = conf_.getParameter<int>(   "RecHitBin");
00747     double RecHitMin   = conf_.getParameter<double>("RecHitMin");
00748     double RecHitMax   = conf_.getParameter<double>("RecHitMax");
00749 
00750     int    RecLayBin   = conf_.getParameter<int>(   "RecHitBin");
00751     double RecLayMin   = conf_.getParameter<double>("RecHitMin");
00752     double RecLayMax   = conf_.getParameter<double>("RecHitMax");
00753 
00754 
00755     int    PhiBin     = conf_.getParameter<int>(   "PhiBin");
00756     double PhiMin     = conf_.getParameter<double>("PhiMin");
00757     double PhiMax     = conf_.getParameter<double>("PhiMax");
00758 
00759     int    EtaBin     = conf_.getParameter<int>(   "EtaBin");
00760     double EtaMin     = conf_.getParameter<double>("EtaMin");
00761     double EtaMax     = conf_.getParameter<double>("EtaMax");
00762 
00763     int    ThetaBin   = conf_.getParameter<int>(   "ThetaBin");
00764     double ThetaMin   = conf_.getParameter<double>("ThetaMin");
00765     double ThetaMax   = conf_.getParameter<double>("ThetaMax");
00766 
00767     int    TrackQBin  = conf_.getParameter<int>(   "TrackQBin");
00768     double TrackQMin  = conf_.getParameter<double>("TrackQMin");
00769     double TrackQMax  = conf_.getParameter<double>("TrackQMax");
00770 
00771     int    TrackPtBin = conf_.getParameter<int>(   "TrackPtBin");
00772     double TrackPtMin = conf_.getParameter<double>("TrackPtMin");
00773     double TrackPtMax = conf_.getParameter<double>("TrackPtMax");
00774 
00775     int    TrackPBin  = conf_.getParameter<int>(   "TrackPBin");
00776     double TrackPMin  = conf_.getParameter<double>("TrackPMin");
00777     double TrackPMax  = conf_.getParameter<double>("TrackPMax");
00778 
00779     int    TrackPxBin = conf_.getParameter<int>(   "TrackPxBin");
00780     double TrackPxMin = conf_.getParameter<double>("TrackPxMin");
00781     double TrackPxMax = conf_.getParameter<double>("TrackPxMax");
00782 
00783     int    TrackPyBin = conf_.getParameter<int>(   "TrackPyBin");
00784     double TrackPyMin = conf_.getParameter<double>("TrackPyMin");
00785     double TrackPyMax = conf_.getParameter<double>("TrackPyMax");
00786 
00787     int    TrackPzBin = conf_.getParameter<int>(   "TrackPzBin");
00788     double TrackPzMin = conf_.getParameter<double>("TrackPzMin");
00789     double TrackPzMax = conf_.getParameter<double>("TrackPzMax");
00790 
00791     int    ptErrBin   = conf_.getParameter<int>(   "ptErrBin");
00792     double ptErrMin   = conf_.getParameter<double>("ptErrMin");
00793     double ptErrMax   = conf_.getParameter<double>("ptErrMax");
00794 
00795     int    pxErrBin   = conf_.getParameter<int>(   "pxErrBin");
00796     double pxErrMin   = conf_.getParameter<double>("pxErrMin");
00797     double pxErrMax   = conf_.getParameter<double>("pxErrMax");
00798 
00799     int    pyErrBin   = conf_.getParameter<int>(   "pyErrBin");
00800     double pyErrMin   = conf_.getParameter<double>("pyErrMin");
00801     double pyErrMax   = conf_.getParameter<double>("pyErrMax");
00802 
00803     int    pzErrBin   = conf_.getParameter<int>(   "pzErrBin");
00804     double pzErrMin   = conf_.getParameter<double>("pzErrMin");
00805     double pzErrMax   = conf_.getParameter<double>("pzErrMax");
00806 
00807     int    pErrBin    = conf_.getParameter<int>(   "pErrBin");
00808     double pErrMin    = conf_.getParameter<double>("pErrMin");
00809     double pErrMax    = conf_.getParameter<double>("pErrMax");
00810 
00811     int    phiErrBin  = conf_.getParameter<int>(   "phiErrBin");
00812     double phiErrMin  = conf_.getParameter<double>("phiErrMin");
00813     double phiErrMax  = conf_.getParameter<double>("phiErrMax");
00814 
00815     int    etaErrBin  = conf_.getParameter<int>(   "etaErrBin");
00816     double etaErrMin  = conf_.getParameter<double>("etaErrMin");
00817     double etaErrMax  = conf_.getParameter<double>("etaErrMax");
00818 
00819 
00820     dqmStore_->setCurrentFolder(MEFolderName);
00821 
00822     TkParameterMEs tkmes;
00823 
00824     std::string histTag = (sname == "default") ? CatagoryName : sname + "_" + CatagoryName;
00825 
00826     if(doAllPlots_)
00827     {
00828 
00829         //        COMMENTED BEACUSE THERE IS ALREADY THE PROFILE !!! (blablaProfile)
00830         /*
00831         // hit properties
00832         dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties");
00833         // rechits
00834         histname = "NumberOfRecHitsPerTrackVsPhi_" + histTag;
00835         tkmes.NumberOfRecHitsPerTrackVsPhi = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, RecHitMin, RecHitMax,"");
00836         tkmes.NumberOfRecHitsPerTrackVsPhi->setAxisTitle("Track #phi",1);
00837         tkmes.NumberOfRecHitsPerTrackVsPhi->setAxisTitle("Number of found RecHits of each Track",2);
00838 
00839         if (doThetaPlots_) {
00840           histname = "NumberOfRecHitsPerTrackVsTheta_" + histTag;
00841           tkmes.NumberOfRecHitsPerTrackVsTheta = dqmStore_->bookProfile(histname, histname, ThetaBin, ThetaMin, ThetaMax, RecHitMin, RecHitMax,"");
00842           tkmes.NumberOfRecHitsPerTrackVsTheta->setAxisTitle("Track #theta",1);
00843           tkmes.NumberOfRecHitsPerTrackVsTheta->setAxisTitle("Number of found RecHits of each Track",2);
00844         }
00845         histname = "NumberOfRecHitsPerTrackVsEta_" + histTag;
00846         tkmes.NumberOfRecHitsPerTrackVsEta = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, RecHitMin, RecHitMax,"");
00847         tkmes.NumberOfRecHitsPerTrackVsEta->setAxisTitle("Track #eta",1);
00848         tkmes.NumberOfRecHitsPerTrackVsEta->setAxisTitle("Number of found RecHits of each Track",2);
00849 
00850         histname = "NumberOfLayersPerTrackVsPhi_" + histTag;
00851         tkmes.NumberOfLayersPerTrackVsPhi = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, RecLayMin, RecLayMax,"");
00852         tkmes.NumberOfLayersPerTrackVsPhi->setAxisTitle("Track #phi",1);
00853         tkmes.NumberOfLayersPerTrackVsPhi->setAxisTitle("Number of Layers of each Track",2);
00854 
00855         if (doThetaPlots_) {
00856           histname = "NumberOfLayersPerTrackVsTheta_" + histTag;
00857           tkmes.NumberOfLayersPerTrackVsTheta = dqmStore_->bookProfile(histname, histname, ThetaBin, ThetaMin, ThetaMax, RecLayMin, RecLayMax,"");
00858           tkmes.NumberOfLayersPerTrackVsTheta->setAxisTitle("Track #theta",1);
00859           tkmes.NumberOfLayersPerTrackVsTheta->setAxisTitle("Number of Layers of each Track",2);
00860         }
00861         histname = "NumberOfLayersPerTrackVsEta_" + histTag;
00862         tkmes.NumberOfLayersPerTrackVsEta = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, RecLayMin, RecLayMax,"");
00863         tkmes.NumberOfLayersPerTrackVsEta->setAxisTitle("Track #eta",1);
00864         tkmes.NumberOfLayersPerTrackVsEta->setAxisTitle("Number of Layers of each Track",2);
00865         */
00866 
00867         // general properties
00868         dqmStore_->setCurrentFolder(MEFolderName+"/GeneralProperties");
00869 
00870         if (doThetaPlots_) {
00871           histname = "Chi2oNDFVsTheta_" + histTag;
00872           tkmes.Chi2oNDFVsTheta = dqmStore_->bookProfile(histname, histname, ThetaBin, ThetaMin, ThetaMax, Chi2NDFMin, Chi2NDFMax,"");
00873           tkmes.Chi2oNDFVsTheta->setAxisTitle("Track #theta",1);
00874           tkmes.Chi2oNDFVsTheta->setAxisTitle("Track #chi^{2}/ndf",2);
00875         }
00876         histname = "Chi2oNDFVsPhi_" + histTag;
00877         tkmes.Chi2oNDFVsPhi   = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, Chi2NDFMin, Chi2NDFMax,"");
00878         tkmes.Chi2oNDFVsPhi->setAxisTitle("Track #phi",1);
00879         tkmes.Chi2oNDFVsPhi->setAxisTitle("Track #chi^{2}/ndf",2);
00880 
00881         histname = "Chi2oNDFVsEta_" + histTag;
00882         tkmes.Chi2oNDFVsEta   = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, Chi2NDFMin, Chi2NDFMax,"");
00883         tkmes.Chi2oNDFVsEta->setAxisTitle("Track #eta",1);
00884         tkmes.Chi2oNDFVsEta->setAxisTitle("Track #chi^{2}/ndf",2);
00885     }
00886 
00887     // general properties
00888     dqmStore_->setCurrentFolder(MEFolderName+"/GeneralProperties");
00889 
00890     histname = "TrackP_" + histTag;
00891 
00892     tkmes.TrackP = dqmStore_->book1D(histname, histname, TrackPBin, TrackPMin, TrackPMax);
00893     tkmes.TrackP->setAxisTitle("Track |p| (GeV/c)", 1);
00894     tkmes.TrackP->setAxisTitle("Number of Tracks",2);
00895 
00896     histname = "TrackPt_" + histTag;
00897     tkmes.TrackPt = dqmStore_->book1D(histname, histname, TrackPtBin, TrackPtMin, TrackPtMax);
00898     tkmes.TrackPt->setAxisTitle("Track p_{T} (GeV/c)", 1);
00899     tkmes.TrackPt->setAxisTitle("Number of Tracks",2);
00900 
00901     if (doTrackPxPyPlots_) {
00902       histname = "TrackPx_" + histTag;
00903       tkmes.TrackPx = dqmStore_->book1D(histname, histname, TrackPxBin, TrackPxMin, TrackPxMax);
00904       tkmes.TrackPx->setAxisTitle("Track p_{x} (GeV/c)", 1);
00905       tkmes.TrackPx->setAxisTitle("Number of Tracks",2);
00906 
00907       histname = "TrackPy_" + histTag;
00908       tkmes.TrackPy = dqmStore_->book1D(histname, histname, TrackPyBin, TrackPyMin, TrackPyMax);
00909       tkmes.TrackPy->setAxisTitle("Track p_{y} (GeV/c)", 1);
00910       tkmes.TrackPy->setAxisTitle("Number of Tracks",2);
00911     }
00912     histname = "TrackPz_" + histTag;
00913     tkmes.TrackPz = dqmStore_->book1D(histname, histname, TrackPzBin, TrackPzMin, TrackPzMax);
00914     tkmes.TrackPz->setAxisTitle("Track p_{z} (GeV/c)", 1);
00915     tkmes.TrackPz->setAxisTitle("Number of Tracks",2);
00916 
00917     histname = "TrackPhi_" + histTag;
00918     tkmes.TrackPhi = dqmStore_->book1D(histname, histname, PhiBin, PhiMin, PhiMax);
00919     tkmes.TrackPhi->setAxisTitle("Track #phi", 1);
00920     tkmes.TrackPhi->setAxisTitle("Number of Tracks",2);
00921 
00922     histname = "TrackEta_" + histTag;
00923     tkmes.TrackEta = dqmStore_->book1D(histname, histname, EtaBin, EtaMin, EtaMax);
00924     tkmes.TrackEta->setAxisTitle("Track #eta", 1);
00925     tkmes.TrackEta->setAxisTitle("Number of Tracks",2);
00926 
00927     if (doThetaPlots_) {  
00928       histname = "TrackTheta_" + histTag;
00929       tkmes.TrackTheta = dqmStore_->book1D(histname, histname, ThetaBin, ThetaMin, ThetaMax);
00930       tkmes.TrackTheta->setAxisTitle("Track #theta", 1);
00931       tkmes.TrackTheta->setAxisTitle("Number of Tracks",2);
00932     }
00933     histname = "TrackQ_" + histTag;
00934     tkmes.TrackQ = dqmStore_->book1D(histname, histname, TrackQBin, TrackQMin, TrackQMax);
00935     tkmes.TrackQ->setAxisTitle("Track Charge", 1);
00936     tkmes.TrackQ->setAxisTitle("Number of Tracks",2);
00937 
00938     histname = "TrackPErrOverP_" + histTag;
00939     tkmes.TrackPErr = dqmStore_->book1D(histname, histname, pErrBin, pErrMin, pErrMax);
00940     tkmes.TrackPErr->setAxisTitle("error(p)/p", 1);
00941     tkmes.TrackPErr->setAxisTitle("Number of Tracks",2);
00942 
00943     histname = "TrackPtErrOverPt_" + histTag;
00944     tkmes.TrackPtErr = dqmStore_->book1D(histname, histname, ptErrBin, ptErrMin, ptErrMax);
00945     tkmes.TrackPtErr->setAxisTitle("error(p_{T})/p_{T}", 1);
00946     tkmes.TrackPtErr->setAxisTitle("Number of Tracks",2);
00947 
00948     histname = "TrackPtErrOverPtVsEta_" + histTag;
00949     tkmes.TrackPtErrVsEta = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, ptErrMin, ptErrMax);
00950     tkmes.TrackPtErrVsEta->setAxisTitle("Track #eta",1);
00951     tkmes.TrackPtErrVsEta->setAxisTitle("error(p_{T})/p_{T}", 2);
00952 
00953     if (doTrackPxPyPlots_) {
00954       histname = "TrackPxErrOverPx_" + histTag;
00955       tkmes.TrackPxErr = dqmStore_->book1D(histname, histname, pxErrBin, pxErrMin, pxErrMax);
00956       tkmes.TrackPxErr->setAxisTitle("error(p_{x})/p_{x}", 1);
00957       tkmes.TrackPxErr->setAxisTitle("Number of Tracks",2);
00958       
00959       histname = "TrackPyErrOverPy_" + histTag;
00960       tkmes.TrackPyErr = dqmStore_->book1D(histname, histname, pyErrBin, pyErrMin, pyErrMax);
00961       tkmes.TrackPyErr->setAxisTitle("error(p_{y})/p_{y}", 1);
00962       tkmes.TrackPyErr->setAxisTitle("Number of Tracks",2);
00963     }
00964     histname = "TrackPzErrOverPz_" + histTag;
00965     tkmes.TrackPzErr = dqmStore_->book1D(histname, histname, pzErrBin, pzErrMin, pzErrMax);
00966     tkmes.TrackPzErr->setAxisTitle("error(p_{z})/p_{z}", 1);
00967     tkmes.TrackPzErr->setAxisTitle("Number of Tracks",2);
00968 
00969     histname = "TrackPhiErr_" + histTag;
00970     tkmes.TrackPhiErr = dqmStore_->book1D(histname, histname, phiErrBin, phiErrMin, phiErrMax);
00971     tkmes.TrackPhiErr->setAxisTitle("error(#phi)");
00972     tkmes.TrackPhiErr->setAxisTitle("Number of Tracks",2);
00973 
00974     histname = "TrackEtaErr_" + histTag;
00975     tkmes.TrackEtaErr = dqmStore_->book1D(histname, histname, etaErrBin, etaErrMin, etaErrMax);
00976     tkmes.TrackEtaErr->setAxisTitle("error(#eta)");
00977     tkmes.TrackEtaErr->setAxisTitle("Number of Tracks",2);
00978 
00979     // rec hit profiles
00980     histname = "NumberOfRecHitsPerTrackVsPhiProfile_" + histTag;
00981     tkmes.NumberOfRecHitsPerTrackVsPhiProfile = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, RecHitBin, RecHitMin, RecHitMax,"");
00982     tkmes.NumberOfRecHitsPerTrackVsPhiProfile->setAxisTitle("Track #phi",1);
00983     tkmes.NumberOfRecHitsPerTrackVsPhiProfile->setAxisTitle("Number of found RecHits of each Track",2);
00984 
00985     if (doThetaPlots_) {
00986       histname = "NumberOfRecHitsPerTrackVsThetaProfile_" + histTag;
00987       tkmes.NumberOfRecHitsPerTrackVsThetaProfile = dqmStore_->bookProfile(histname, histname, ThetaBin, ThetaMin, ThetaMax, RecHitBin, RecHitMin, RecHitMax,"");
00988       tkmes.NumberOfRecHitsPerTrackVsThetaProfile->setAxisTitle("Track #phi",1);
00989       tkmes.NumberOfRecHitsPerTrackVsThetaProfile->setAxisTitle("Number of found RecHits of each Track",2);
00990     }
00991     histname = "NumberOfRecHitsPerTrackVsEtaProfile_" + histTag;
00992     tkmes.NumberOfRecHitsPerTrackVsEtaProfile = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, RecHitBin, RecHitMin, RecHitMax,"");
00993     tkmes.NumberOfRecHitsPerTrackVsEtaProfile->setAxisTitle("Track #eta",1);
00994     tkmes.NumberOfRecHitsPerTrackVsEtaProfile->setAxisTitle("Number of found RecHits of each Track",2);
00995 
00996     histname = "NumberOfLayersPerTrackVsPhiProfile_" + histTag;
00997     tkmes.NumberOfLayersPerTrackVsPhiProfile = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, RecLayBin, RecLayMin, RecLayMax,"");
00998     tkmes.NumberOfLayersPerTrackVsPhiProfile->setAxisTitle("Track #phi",1);
00999     tkmes.NumberOfLayersPerTrackVsPhiProfile->setAxisTitle("Number of Layers of each Track",2);
01000 
01001     if (doThetaPlots_) {
01002       histname = "NumberOfLayersPerTrackVsThetaProfile_" + histTag;
01003       tkmes.NumberOfLayersPerTrackVsThetaProfile = dqmStore_->bookProfile(histname, histname, ThetaBin, ThetaMin, ThetaMax, RecLayBin, RecLayMin, RecLayMax,"");
01004       tkmes.NumberOfLayersPerTrackVsThetaProfile->setAxisTitle("Track #phi",1);
01005       tkmes.NumberOfLayersPerTrackVsThetaProfile->setAxisTitle("Number of Layers of each Track",2);
01006     }
01007     histname = "NumberOfLayersPerTrackVsEtaProfile_" + histTag;
01008     tkmes.NumberOfLayersPerTrackVsEtaProfile = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, RecLayBin, RecLayMin, RecLayMax,"");
01009     tkmes.NumberOfLayersPerTrackVsEtaProfile->setAxisTitle("Track #eta",1);
01010     tkmes.NumberOfLayersPerTrackVsEtaProfile->setAxisTitle("Number of Layers of each Track",2);
01011 
01012     if ( doGoodTrackPlots_ || doAllPlots_ ) {
01013 
01014       dqmStore_->setCurrentFolder(MEFolderName+"/GeneralProperties/GoodTracks");
01015 
01016       histname = "GoodTrackPt_" + histTag;
01017       tkmes.GoodTrackPt = dqmStore_->book1D(histname, histname, TrackPtBin, TrackPtMin, TrackPtMax);
01018       tkmes.GoodTrackPt->setAxisTitle("Good Track p_{T} (GeV/c)", 1);
01019       tkmes.GoodTrackPt->setAxisTitle("Number of Tracks",2);
01020       
01021       histname = "GoodTrackEta_" + histTag;
01022       tkmes.GoodTrackEta = dqmStore_->book1D(histname, histname, EtaBin, EtaMin, EtaMax);
01023       tkmes.GoodTrackEta->setAxisTitle("Good Track #eta", 1);
01024       tkmes.GoodTrackEta->setAxisTitle("Number of Tracks",2);
01025       
01026       histname = "GoodTrackPhi_" + histTag;
01027       tkmes.GoodTrackPhi = dqmStore_->book1D(histname, histname, PhiBin, PhiMin, PhiMax);
01028       tkmes.GoodTrackPhi->setAxisTitle("Good Track #phi", 1);
01029       tkmes.GoodTrackPhi->setAxisTitle("Number of Tracks",2);
01030 
01031     }
01032 
01033     // now put the MEs in the map
01034     TkParameterMEMap.insert( std::make_pair(sname, tkmes) );
01035 }
01036 
01037 
01038 // fill histograms at differnt measurement points
01039 // ---------------------------------------------------------------------------------//
01040 void TrackAnalyzer::fillHistosForState(const edm::EventSetup& iSetup, const reco::Track & track, std::string sname) 
01041 {
01042     //get the kinematic parameters
01043     double p, px, py, pz, pt, theta, phi, eta, q;
01044     double pxerror, pyerror, pzerror, pterror, perror, phierror, etaerror;
01045 
01046     bool isHighPurity = track.quality(reco::TrackBase::highPurity);
01047 
01048     if (sname == "default") 
01049     {
01050         p     = track.p();
01051         px    = track.px();
01052         py    = track.py();
01053         pz    = track.pz();
01054         pt    = track.pt();
01055         phi   = track.phi();
01056         theta = track.theta();
01057         eta   = track.eta();
01058         q     = track.charge();
01059 
01060         //        pterror  = (pt) ? track.ptError()/pt : 0.0;
01061         pterror  = (pt) ? track.ptError()/(pt*pt) : 0.0;
01062         pxerror  = -1.0;
01063         pyerror  = -1.0;
01064         pzerror  = -1.0;
01065         perror   = -1.0;
01066         phierror = track.phiError();
01067         etaerror = track.etaError();
01068 
01069     } 
01070     else 
01071     {
01072         edm::ESHandle<TransientTrackBuilder> theB;
01073         iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
01074         reco::TransientTrack TransTrack = theB->build(track);
01075 
01076         TrajectoryStateOnSurface TSOS;
01077 
01078         if      (sname == "OuterSurface")  TSOS = TransTrack.outermostMeasurementState();
01079         else if (sname == "InnerSurface")  TSOS = TransTrack.innermostMeasurementState();
01080         else if (sname == "ImpactPoint")   TSOS = TransTrack.impactPointState();
01081 
01082         p     = TSOS.globalMomentum().mag();
01083         px    = TSOS.globalMomentum().x();
01084         py    = TSOS.globalMomentum().y();
01085         pz    = TSOS.globalMomentum().z();
01086         pt    = TSOS.globalMomentum().perp();
01087         phi   = TSOS.globalMomentum().phi();
01088         theta = TSOS.globalMomentum().theta();
01089         eta   = TSOS.globalMomentum().eta();
01090         q     = TSOS.charge();
01091 
01092         //get the error of the kinimatic parameters
01093         AlgebraicSymMatrix66 errors = TSOS.cartesianError().matrix();
01094         double partialPterror = errors(3,3)*pow(TSOS.globalMomentum().x(),2) + errors(4,4)*pow(TSOS.globalMomentum().y(),2);
01095         pterror  = sqrt(partialPterror)/TSOS.globalMomentum().perp();
01096         pxerror  = sqrt(errors(3,3))/TSOS.globalMomentum().x();
01097         pyerror  = sqrt(errors(4,4))/TSOS.globalMomentum().y();
01098         pzerror  = sqrt(errors(5,5))/TSOS.globalMomentum().z();
01099         perror   = sqrt(partialPterror+errors(5,5)*pow(TSOS.globalMomentum().z(),2))/TSOS.globalMomentum().mag();
01100         phierror = sqrt(TSOS.curvilinearError().matrix()(2,2));
01101         etaerror = sqrt(TSOS.curvilinearError().matrix()(1,1))*fabs(sin(TSOS.globalMomentum().theta()));
01102 
01103     }
01104 
01105     std::map<std::string, TkParameterMEs>::iterator iPos = TkParameterMEMap.find(sname); 
01106     if (iPos != TkParameterMEMap.end()) 
01107     {
01108         TkParameterMEs tkmes = iPos->second;
01109 
01110         // momentum
01111         tkmes.TrackP->Fill(p);
01112         if (doTrackPxPyPlots_) {
01113           tkmes.TrackPx->Fill(px);
01114           tkmes.TrackPy->Fill(py);
01115         }
01116         tkmes.TrackPz->Fill(pz);
01117         tkmes.TrackPt->Fill(pt);
01118 
01119         // angles
01120         tkmes.TrackPhi->Fill(phi);
01121         tkmes.TrackEta->Fill(eta);
01122         if (doThetaPlots_) {
01123           tkmes.TrackTheta->Fill(theta);
01124         }
01125         tkmes.TrackQ->Fill(q);
01126 
01127         // errors
01128         tkmes.TrackPtErr->Fill(pterror);
01129         tkmes.TrackPtErrVsEta->Fill(eta,pterror);
01130         if (doTrackPxPyPlots_) {
01131           tkmes.TrackPxErr->Fill(pxerror);
01132           tkmes.TrackPyErr->Fill(pyerror);
01133         }
01134         tkmes.TrackPzErr->Fill(pzerror);
01135         tkmes.TrackPErr->Fill(perror);
01136         tkmes.TrackPhiErr->Fill(phierror);
01137         tkmes.TrackEtaErr->Fill(etaerror);
01138 
01139         // rec hits 
01140         tkmes.NumberOfRecHitsPerTrackVsPhiProfile->Fill(phi,    track.hitPattern().numberOfValidHits());
01141         if (doThetaPlots_) {
01142           tkmes.NumberOfRecHitsPerTrackVsThetaProfile->Fill(theta,track.hitPattern().numberOfValidHits());
01143         }
01144         tkmes.NumberOfRecHitsPerTrackVsEtaProfile->Fill(eta,    track.hitPattern().numberOfValidHits());
01145 
01146 
01147         // rec layers 
01148         tkmes.NumberOfLayersPerTrackVsPhiProfile->Fill(phi,     track.hitPattern().stripLayersWithMeasurement());
01149         if (doThetaPlots_) {
01150           tkmes.NumberOfLayersPerTrackVsThetaProfile->Fill(theta, track.hitPattern().stripLayersWithMeasurement());
01151         }
01152         tkmes.NumberOfLayersPerTrackVsEtaProfile->Fill(eta,     track.hitPattern().stripLayersWithMeasurement());
01153 
01154         if(doAllPlots_)
01155         {
01156           /*
01157             // COMMENTED because the same MEs are already fill below !!!!
01158             // hit related
01159             tkmes.NumberOfRecHitsPerTrackVsPhi->Fill(phi, track.found());
01160             if (doThetaPlots_) {
01161               tkmes.NumberOfRecHitsPerTrackVsTheta->Fill(theta, track.found());
01162             }
01163             tkmes.NumberOfRecHitsPerTrackVsEta->Fill(eta, track.found());
01164           */
01165 
01166             // general properties
01167             if (doThetaPlots_) {
01168               tkmes.Chi2oNDFVsTheta->Fill(theta, track.normalizedChi2());
01169             }
01170             tkmes.Chi2oNDFVsPhi->Fill(phi, track.normalizedChi2());
01171             tkmes.Chi2oNDFVsEta->Fill(eta, track.normalizedChi2());
01172 
01173             // COMMENTED because there are already those quantity in blablaProfile !!
01174             /*
01175             // rec hits 
01176             tkmes.NumberOfRecHitsPerTrackVsPhi->Fill(phi,        track.hitPattern().numberOfValidHits());
01177             if (doThetaPlots_) {
01178               tkmes.NumberOfRecHitsPerTrackVsTheta->Fill(theta,    track.hitPattern().numberOfValidHits());
01179             }
01180             tkmes.NumberOfRecHitsPerTrackVsEta->Fill(eta,        track.hitPattern().numberOfValidHits());
01181             */      
01182 
01183         }
01184 
01185         if ( doGoodTrackPlots_ || doAllPlots_ ) {
01186           if ( isHighPurity && pt > 1. ) {
01187             tkmes.GoodTrackPt->Fill(pt);
01188             tkmes.GoodTrackEta->Fill(eta);
01189             tkmes.GoodTrackPhi->Fill(phi);
01190           }
01191         }
01192 
01193     }
01194 }
01195 
01196 
01197 void TrackAnalyzer::doTrackerSpecificInitialization(DQMStore * dqmStore_) 
01198 {
01199 
01200     // parameters from the configuration
01201     std::string QualName     = conf_.getParameter<std::string>("Quality");
01202     std::string AlgoName     = conf_.getParameter<std::string>("AlgoName");
01203     std::string MEFolderName = conf_.getParameter<std::string>("FolderName"); 
01204 
01205     // use the AlgoName and Quality Name 
01206     std::string CatagoryName = QualName != "" ? AlgoName + "_" + QualName : AlgoName;
01207 
01208     // get binning from the configuration
01209     int    TOBHitBin = conf_.getParameter<int>(   "TOBHitBin");
01210     double TOBHitMin = conf_.getParameter<double>("TOBHitMin");
01211     double TOBHitMax = conf_.getParameter<double>("TOBHitMax");
01212 
01213     int    TIBHitBin = conf_.getParameter<int>(   "TIBHitBin");
01214     double TIBHitMin = conf_.getParameter<double>("TIBHitMin");
01215     double TIBHitMax = conf_.getParameter<double>("TIBHitMax");
01216 
01217     int    TIDHitBin = conf_.getParameter<int>(   "TIDHitBin");
01218     double TIDHitMin = conf_.getParameter<double>("TIDHitMin");
01219     double TIDHitMax = conf_.getParameter<double>("TIDHitMax");
01220 
01221     int    TECHitBin = conf_.getParameter<int>(   "TECHitBin");
01222     double TECHitMin = conf_.getParameter<double>("TECHitMin");
01223     double TECHitMax = conf_.getParameter<double>("TECHitMax");
01224 
01225     int    PXBHitBin = conf_.getParameter<int>(   "PXBHitBin");
01226     double PXBHitMin = conf_.getParameter<double>("PXBHitMin");
01227     double PXBHitMax = conf_.getParameter<double>("PXBHitMax");
01228 
01229     int    PXFHitBin = conf_.getParameter<int>(   "PXFHitBin");
01230     double PXFHitMin = conf_.getParameter<double>("PXFHitMin");
01231     double PXFHitMax = conf_.getParameter<double>("PXFHitMax");
01232 
01233     int    TOBLayBin = conf_.getParameter<int>(   "TOBLayBin");
01234     double TOBLayMin = conf_.getParameter<double>("TOBLayMin");
01235     double TOBLayMax = conf_.getParameter<double>("TOBLayMax");
01236 
01237     int    TIBLayBin = conf_.getParameter<int>(   "TIBLayBin");
01238     double TIBLayMin = conf_.getParameter<double>("TIBLayMin");
01239     double TIBLayMax = conf_.getParameter<double>("TIBLayMax");
01240 
01241     int    TIDLayBin = conf_.getParameter<int>(   "TIDLayBin");
01242     double TIDLayMin = conf_.getParameter<double>("TIDLayMin");
01243     double TIDLayMax = conf_.getParameter<double>("TIDLayMax");
01244 
01245     int    TECLayBin = conf_.getParameter<int>(   "TECLayBin");
01246     double TECLayMin = conf_.getParameter<double>("TECLayMin");
01247     double TECLayMax = conf_.getParameter<double>("TECLayMax");
01248 
01249     int    PXBLayBin = conf_.getParameter<int>(   "PXBLayBin");
01250     double PXBLayMin = conf_.getParameter<double>("PXBLayMin");
01251     double PXBLayMax = conf_.getParameter<double>("PXBLayMax");
01252 
01253     int    PXFLayBin = conf_.getParameter<int>(   "PXFLayBin");
01254     double PXFLayMin = conf_.getParameter<double>("PXFLayMin");
01255     double PXFLayMax = conf_.getParameter<double>("PXFLayMax");
01256 
01257 
01258     int    PhiBin     = conf_.getParameter<int>(   "PhiBin");
01259     double PhiMin     = conf_.getParameter<double>("PhiMin");
01260     double PhiMax     = conf_.getParameter<double>("PhiMax");
01261 
01262     int    EtaBin     = conf_.getParameter<int>(   "EtaBin");
01263     double EtaMin     = conf_.getParameter<double>("EtaMin");
01264     double EtaMax     = conf_.getParameter<double>("EtaMax");
01265 
01266     // book hit property histograms
01267     // ---------------------------------------------------------------------------------//
01268     dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties");
01269 
01270 
01271     TkParameterMEs tkmes;
01272 
01273     // TOB hits properties
01274     dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties/TOB");
01275 
01276     histname = "NumberOfTOBRecHitsPerTrack_" + CatagoryName;
01277     NumberOfTOBRecHitsPerTrack = dqmStore_->book1D(histname, histname, TOBHitBin, TOBHitMin, TOBHitMax);
01278     NumberOfTOBRecHitsPerTrack->setAxisTitle("Number of TOB found RecHits of each Track",1);
01279     NumberOfTOBRecHitsPerTrack->setAxisTitle("Number of Tracks", 2);
01280 
01281     histname = "NumberOfTOBRecHitsPerTrackVsPhiProfile_" + CatagoryName;
01282     NumberOfTOBRecHitsPerTrackVsPhiProfile = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, TOBHitBin, TOBHitMin, TOBHitMax,"");
01283     NumberOfTOBRecHitsPerTrackVsPhiProfile->setAxisTitle("Track #phi",1);
01284     NumberOfTOBRecHitsPerTrackVsPhiProfile->setAxisTitle("Number of TOB found RecHits of each Track",2);
01285 
01286     histname = "NumberOfTOBRecHitsPerTrackVsEtaProfile_" + CatagoryName;
01287     NumberOfTOBRecHitsPerTrackVsEtaProfile = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, TOBHitBin, TOBHitMin, TOBHitMax,"");
01288     NumberOfTOBRecHitsPerTrackVsEtaProfile->setAxisTitle("Track #eta",1);
01289     NumberOfTOBRecHitsPerTrackVsEtaProfile->setAxisTitle("Number of TOB found RecHits of each Track",2);
01290 
01291     histname = "NumberOfTOBLayersPerTrack_" + CatagoryName;
01292     NumberOfTOBLayersPerTrack = dqmStore_->book1D(histname, histname, TOBLayBin, TOBLayMin, TOBLayMax);
01293     NumberOfTOBLayersPerTrack->setAxisTitle("Number of TOB Layers of each Track",1);
01294     NumberOfTOBLayersPerTrack->setAxisTitle("Number of Tracks", 2);
01295 
01296     histname = "NumberOfTOBLayersPerTrackVsPhiProfile_" + CatagoryName;
01297     NumberOfTOBLayersPerTrackVsPhiProfile = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, TOBLayBin, TOBLayMin, TOBLayMax,"");
01298     NumberOfTOBLayersPerTrackVsPhiProfile->setAxisTitle("Track #phi",1);
01299     NumberOfTOBLayersPerTrackVsPhiProfile->setAxisTitle("Number of TOB Layers of each Track",2);
01300 
01301     histname = "NumberOfTOBLayersPerTrackVsEtaProfile_" + CatagoryName;
01302     NumberOfTOBLayersPerTrackVsEtaProfile = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, TOBLayBin, TOBLayMin, TOBLayMax,"");
01303     NumberOfTOBLayersPerTrackVsEtaProfile->setAxisTitle("Track #eta",1);
01304     NumberOfTOBLayersPerTrackVsEtaProfile->setAxisTitle("Number of TOB Layers of each Track",2);
01305 
01306 
01307     // TIB hits properties
01308     dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties/TIB");
01309 
01310     histname = "NumberOfTIBRecHitsPerTrack_" + CatagoryName;
01311     NumberOfTIBRecHitsPerTrack = dqmStore_->book1D(histname, histname, TIBHitBin, TIBHitMin, TIBHitMax);
01312     NumberOfTIBRecHitsPerTrack->setAxisTitle("Number of TIB found RecHits of each Track",1);
01313     NumberOfTIBRecHitsPerTrack->setAxisTitle("Number of Tracks", 2);
01314 
01315     histname = "NumberOfTIBRecHitsPerTrackVsPhiProfile_" + CatagoryName;
01316     NumberOfTIBRecHitsPerTrackVsPhiProfile = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, TIBHitBin, TIBHitMin, TIBHitMax,"");
01317     NumberOfTIBRecHitsPerTrackVsPhiProfile->setAxisTitle("Track #phi",1);
01318     NumberOfTIBRecHitsPerTrackVsPhiProfile->setAxisTitle("Number of TIB found RecHits of each Track",2);
01319 
01320     histname = "NumberOfTIBRecHitsPerTrackVsEtaProfile_" + CatagoryName;
01321     NumberOfTIBRecHitsPerTrackVsEtaProfile = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, TIBHitBin, TIBHitMin, TIBHitMax,"");
01322     NumberOfTIBRecHitsPerTrackVsEtaProfile->setAxisTitle("Track #eta",1);
01323     NumberOfTIBRecHitsPerTrackVsEtaProfile->setAxisTitle("Number of TIB found RecHits of each Track",2);
01324 
01325     histname = "NumberOfTIBLayersPerTrack_" + CatagoryName;
01326     NumberOfTIBLayersPerTrack = dqmStore_->book1D(histname, histname, TIBLayBin, TIBLayMin, TIBLayMax);
01327     NumberOfTIBLayersPerTrack->setAxisTitle("Number of TIB Layers of each Track",1);
01328     NumberOfTIBLayersPerTrack->setAxisTitle("Number of Tracks", 2);
01329 
01330     histname = "NumberOfTIBLayersPerTrackVsPhiProfile_" + CatagoryName;
01331     NumberOfTIBLayersPerTrackVsPhiProfile = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, TIBLayBin, TIBLayMin, TIBLayMax,"");
01332     NumberOfTIBLayersPerTrackVsPhiProfile->setAxisTitle("Track #phi",1);
01333     NumberOfTIBLayersPerTrackVsPhiProfile->setAxisTitle("Number of TIB Layers of each Track",2);
01334 
01335     histname = "NumberOfTIBLayersPerTrackVsEtaProfile_" + CatagoryName;
01336     NumberOfTIBLayersPerTrackVsEtaProfile = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, TIBLayBin, TIBLayMin, TIBLayMax,"");
01337     NumberOfTIBLayersPerTrackVsEtaProfile->setAxisTitle("Track #eta",1);
01338     NumberOfTIBLayersPerTrackVsEtaProfile->setAxisTitle("Number of TIB Layers of each Track",2);
01339 
01340 
01341     // TID hit properties
01342     dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties/TID");
01343 
01344     histname = "NumberOfTIDRecHitsPerTrack_" + CatagoryName;
01345     NumberOfTIDRecHitsPerTrack = dqmStore_->book1D(histname, histname, TIDHitBin, TIDHitMin, TIDHitMax);
01346     NumberOfTIDRecHitsPerTrack->setAxisTitle("Number of TID found RecHits of each Track",1);
01347     NumberOfTIDRecHitsPerTrack->setAxisTitle("Number of Tracks", 2);
01348 
01349     histname = "NumberOfTIDRecHitsPerTrackVsPhiProfile_" + CatagoryName;
01350     NumberOfTIDRecHitsPerTrackVsPhiProfile = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, TIDHitBin, TIDHitMin, TIDHitMax,"");
01351     NumberOfTIDRecHitsPerTrackVsPhiProfile->setAxisTitle("Track #phi",1);
01352     NumberOfTIDRecHitsPerTrackVsPhiProfile->setAxisTitle("Number of TID found RecHits of each Track",2);
01353 
01354     histname = "NumberOfTIDRecHitsPerTrackVsEtaProfile_" + CatagoryName;
01355     NumberOfTIDRecHitsPerTrackVsEtaProfile = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, TIDHitBin, TIDHitMin, TIDHitMax,"");
01356     NumberOfTIDRecHitsPerTrackVsEtaProfile->setAxisTitle("Track #eta",1);
01357     NumberOfTIDRecHitsPerTrackVsEtaProfile->setAxisTitle("Number of TID found RecHits of each Track",2);
01358 
01359     histname = "NumberOfTIDLayersPerTrack_" + CatagoryName;
01360     NumberOfTIDLayersPerTrack = dqmStore_->book1D(histname, histname, TIDLayBin, TIDLayMin, TIDLayMax);
01361     NumberOfTIDLayersPerTrack->setAxisTitle("Number of TID Layers of each Track",1);
01362     NumberOfTIDLayersPerTrack->setAxisTitle("Number of Tracks", 2);
01363 
01364     histname = "NumberOfTIDLayersPerTrackVsPhiProfile_" + CatagoryName;
01365     NumberOfTIDLayersPerTrackVsPhiProfile = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, TIDLayBin, TIDLayMin, TIDLayMax,"");
01366     NumberOfTIDLayersPerTrackVsPhiProfile->setAxisTitle("Track #phi",1);
01367     NumberOfTIDLayersPerTrackVsPhiProfile->setAxisTitle("Number of TID Layers of each Track",2);
01368 
01369     histname = "NumberOfTIDLayersPerTrackVsEtaProfile_" + CatagoryName;
01370     NumberOfTIDLayersPerTrackVsEtaProfile = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, TIDLayBin, TIDLayMin, TIDLayMax,"");
01371     NumberOfTIDLayersPerTrackVsEtaProfile->setAxisTitle("Track #eta",1);
01372     NumberOfTIDLayersPerTrackVsEtaProfile->setAxisTitle("Number of TID Layers of each Track",2);
01373 
01374 
01375     // TEC hits properties    
01376     dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties/TEC");
01377 
01378     histname = "NumberOfTECRecHitsPerTrack_"+ CatagoryName;
01379     NumberOfTECRecHitsPerTrack = dqmStore_->book1D(histname, histname, TECHitBin, TECHitMin, TECHitMax);
01380     NumberOfTECRecHitsPerTrack->setAxisTitle("Number of TEC found RecHits of each Track",1);
01381     NumberOfTECRecHitsPerTrack->setAxisTitle("Number of Tracks", 2);
01382 
01383     histname = "NumberOfTECRecHitsPerTrackVsPhiProfile_" + CatagoryName;
01384     NumberOfTECRecHitsPerTrackVsPhiProfile = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, TECHitBin, TECHitMin, TECHitMax,"");
01385     NumberOfTECRecHitsPerTrackVsPhiProfile->setAxisTitle("Track #phi",1);
01386     NumberOfTECRecHitsPerTrackVsPhiProfile->setAxisTitle("Number of TEC found RecHits of each Track",2);
01387 
01388     histname = "NumberOfTECRecHitsPerTrackVsEtaProfile_" + CatagoryName;
01389     NumberOfTECRecHitsPerTrackVsEtaProfile = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, TECHitBin, TECHitMin, TECHitMax,"");
01390     NumberOfTECRecHitsPerTrackVsEtaProfile->setAxisTitle("Track #eta",1);
01391     NumberOfTECRecHitsPerTrackVsEtaProfile->setAxisTitle("Number of TEC found RecHits of each Track",2);
01392 
01393     histname = "NumberOfTECLayersPerTrack_"+ CatagoryName;
01394     NumberOfTECLayersPerTrack = dqmStore_->book1D(histname, histname, TECLayBin, TECLayMin, TECLayMax);
01395     NumberOfTECLayersPerTrack->setAxisTitle("Number of TEC Layers of each Track",1);
01396     NumberOfTECLayersPerTrack->setAxisTitle("Number of Tracks", 2);
01397 
01398     histname = "NumberOfTECLayersPerTrackVsPhiProfile_" + CatagoryName;
01399     NumberOfTECLayersPerTrackVsPhiProfile = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, TECLayBin, TECLayMin, TECLayMax,"");
01400     NumberOfTECLayersPerTrackVsPhiProfile->setAxisTitle("Track #phi",1);
01401     NumberOfTECLayersPerTrackVsPhiProfile->setAxisTitle("Number of TEC Layers of each Track",2);
01402 
01403     histname = "NumberOfTECLayersPerTrackVsEtaProfile_" + CatagoryName;
01404     NumberOfTECLayersPerTrackVsEtaProfile = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, TECLayBin, TECLayMin, TECLayMax,"");
01405     NumberOfTECLayersPerTrackVsEtaProfile->setAxisTitle("Track #eta",1);
01406     NumberOfTECLayersPerTrackVsEtaProfile->setAxisTitle("Number of TEC Layers of each Track",2);
01407 
01408 
01409     // PixBarrel hits properties
01410     dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties/PixBarrel");
01411 
01412     histname = "NumberOfPixBarrelRecHitsPerTrack_" + CatagoryName;
01413     NumberOfPixBarrelRecHitsPerTrack = dqmStore_->book1D(histname, histname, PXBHitBin, PXBHitMin, PXBHitMax);
01414     NumberOfPixBarrelRecHitsPerTrack->setAxisTitle("Number of Pixel Barrel found RecHits of each Track",1);
01415     NumberOfPixBarrelRecHitsPerTrack->setAxisTitle("Number of Tracks", 2);
01416 
01417     histname = "NumberOfPixBarrelRecHitsPerTrackVsPhiProfile_" + CatagoryName;
01418     NumberOfPixBarrelRecHitsPerTrackVsPhiProfile = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, PXBHitBin, PXBHitMin, PXBHitMax,"");
01419     NumberOfPixBarrelRecHitsPerTrackVsPhiProfile->setAxisTitle("Track #phi",1);
01420     NumberOfPixBarrelRecHitsPerTrackVsPhiProfile->setAxisTitle("Number of Pixel Barrel found RecHits of each Track",2);
01421 
01422     histname = "NumberOfPixBarrelRecHitsPerTrackVsEtaProfile_" + CatagoryName;
01423     NumberOfPixBarrelRecHitsPerTrackVsEtaProfile = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, PXBHitBin, PXBHitMin, PXBHitMax,"");
01424     NumberOfPixBarrelRecHitsPerTrackVsEtaProfile->setAxisTitle("Track #eta",1);
01425     NumberOfPixBarrelRecHitsPerTrackVsEtaProfile->setAxisTitle("Number of Pixel Barrel found RecHits of each Track",2);
01426 
01427     histname = "NumberOfPixBarrelLayersPerTrack_" + CatagoryName;
01428     NumberOfPixBarrelLayersPerTrack = dqmStore_->book1D(histname, histname, PXBLayBin, PXBLayMin, PXBLayMax);
01429     NumberOfPixBarrelLayersPerTrack->setAxisTitle("Number of Pixel Barrel Layers of each Track",1);
01430     NumberOfPixBarrelLayersPerTrack->setAxisTitle("Number of Tracks", 2);
01431 
01432     histname = "NumberOfPixBarrelLayersPerTrackVsPhiProfile_" + CatagoryName;
01433     NumberOfPixBarrelLayersPerTrackVsPhiProfile = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, PXBLayBin, PXBLayMin, PXBLayMax,"");
01434     NumberOfPixBarrelLayersPerTrackVsPhiProfile->setAxisTitle("Track #phi",1);
01435     NumberOfPixBarrelLayersPerTrackVsPhiProfile->setAxisTitle("Number of Pixel Barrel Layers of each Track",2);
01436 
01437     histname = "NumberOfPixBarrelLayersPerTrackVsEtaProfile_" + CatagoryName;
01438     NumberOfPixBarrelLayersPerTrackVsEtaProfile = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, PXBLayBin, PXBLayMin, PXBLayMax,"");
01439     NumberOfPixBarrelLayersPerTrackVsEtaProfile->setAxisTitle("Track #eta",1);
01440     NumberOfPixBarrelLayersPerTrackVsEtaProfile->setAxisTitle("Number of Pixel Barrel Layers of each Track",2);
01441 
01442 
01443     // PixEndcap hits profiles
01444     dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties/PixEndcap");
01445 
01446     histname = "NumberOfPixEndcapRecHitsPerTrack_" + CatagoryName;
01447     NumberOfPixEndcapRecHitsPerTrack = dqmStore_->book1D(histname, histname, PXFHitBin, PXFHitMin, PXFHitMax);
01448     NumberOfPixEndcapRecHitsPerTrack->setAxisTitle("Number of Pixel Endcap found RecHits of each Track",1);
01449     NumberOfPixEndcapRecHitsPerTrack->setAxisTitle("Number of Tracks", 2);
01450 
01451     histname = "NumberOfPixEndcapRecHitsPerTrackVsPhiProfile_" + CatagoryName;
01452     NumberOfPixEndcapRecHitsPerTrackVsPhiProfile = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, PXBHitBin, PXBHitMin, PXBHitMax,"");
01453     NumberOfPixEndcapRecHitsPerTrackVsPhiProfile->setAxisTitle("Track #phi",1);
01454     NumberOfPixEndcapRecHitsPerTrackVsPhiProfile->setAxisTitle("Number of Pixel Endcap found RecHits of each Track",2);
01455 
01456     histname = "NumberOfPixEndcapRecHitsPerTrackVsEtaProfile_" + CatagoryName;
01457     NumberOfPixEndcapRecHitsPerTrackVsEtaProfile = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, PXBHitBin, PXBHitMin, PXBHitMax,"");
01458     NumberOfPixEndcapRecHitsPerTrackVsEtaProfile->setAxisTitle("Track #eta",1);
01459     NumberOfPixEndcapRecHitsPerTrackVsEtaProfile->setAxisTitle("Number of Pixel Endcap found RecHits of each Track",2);
01460 
01461     histname = "NumberOfPixEndcapLayersPerTrack_" + CatagoryName;
01462     NumberOfPixEndcapLayersPerTrack = dqmStore_->book1D(histname, histname, PXFLayBin, PXFLayMin, PXFLayMax);
01463     NumberOfPixEndcapLayersPerTrack->setAxisTitle("Number of Pixel Endcap Layers of each Track",1);
01464     NumberOfPixEndcapLayersPerTrack->setAxisTitle("Number of Tracks", 2);
01465 
01466     histname = "NumberOfPixEndcapLayersPerTrackVsPhiProfile_" + CatagoryName;
01467     NumberOfPixEndcapLayersPerTrackVsPhiProfile = dqmStore_->bookProfile(histname, histname, PhiBin, PhiMin, PhiMax, PXBLayBin, PXBLayMin, PXBLayMax,"");
01468     NumberOfPixEndcapLayersPerTrackVsPhiProfile->setAxisTitle("Track #phi",1);
01469     NumberOfPixEndcapLayersPerTrackVsPhiProfile->setAxisTitle("Number of Pixel Endcap Layers of each Track",2);
01470 
01471     histname = "NumberOfPixEndcapLayersPerTrackVsEtaProfile_" + CatagoryName;
01472     NumberOfPixEndcapLayersPerTrackVsEtaProfile = dqmStore_->bookProfile(histname, histname, EtaBin, EtaMin, EtaMax, PXBLayBin, PXBLayMin, PXBLayMax,"");
01473     NumberOfPixEndcapLayersPerTrackVsEtaProfile->setAxisTitle("Track #eta",1);
01474     NumberOfPixEndcapLayersPerTrackVsEtaProfile->setAxisTitle("Number of Pixel Endcap Layers of each Track",2);
01475 
01476 }
01477 
01478 
01479 void TrackAnalyzer::doTrackerSpecificFillHists(const reco::Track & track) 
01480 {
01481   double phi   = track.phi();
01482   double eta   = track.eta();
01483 
01484   //Fill TIB Layers and RecHits
01485   NumberOfTIBRecHitsPerTrack->Fill(track.hitPattern().numberOfValidStripTIBHits()); 
01486   NumberOfTIBRecHitsPerTrackVsPhiProfile->Fill(phi,    track.hitPattern().numberOfValidStripTIBHits());
01487   NumberOfTIBRecHitsPerTrackVsEtaProfile->Fill(eta,    track.hitPattern().numberOfValidStripTIBHits());
01488 
01489   NumberOfTIBLayersPerTrack->Fill(track.hitPattern().stripTIBLayersWithMeasurement());
01490   NumberOfTIBLayersPerTrackVsPhiProfile->Fill(phi,     track.hitPattern().stripTIBLayersWithMeasurement());
01491   NumberOfTIBLayersPerTrackVsEtaProfile->Fill(eta,     track.hitPattern().stripTIBLayersWithMeasurement());
01492 
01493   //Fill TOB Layers and RecHits
01494   NumberOfTOBRecHitsPerTrack->Fill(track.hitPattern().numberOfValidStripTOBHits());
01495   NumberOfTOBRecHitsPerTrackVsPhiProfile->Fill(phi,    track.hitPattern().numberOfValidStripTOBHits());
01496   NumberOfTOBRecHitsPerTrackVsEtaProfile->Fill(eta,    track.hitPattern().numberOfValidStripTOBHits());
01497 
01498   NumberOfTOBLayersPerTrack->Fill(track.hitPattern().stripTOBLayersWithMeasurement());
01499   NumberOfTOBLayersPerTrackVsPhiProfile->Fill(phi,     track.hitPattern().stripTOBLayersWithMeasurement());
01500   NumberOfTOBLayersPerTrackVsEtaProfile->Fill(eta,     track.hitPattern().stripTOBLayersWithMeasurement());
01501   
01502 
01503   //Fill TID Layers and RecHits
01504   NumberOfTIDRecHitsPerTrack->Fill(track.hitPattern().numberOfValidStripTIDHits());
01505   NumberOfTIDRecHitsPerTrackVsPhiProfile->Fill(phi,    track.hitPattern().numberOfValidStripTIDHits());
01506   NumberOfTIDRecHitsPerTrackVsEtaProfile->Fill(eta,    track.hitPattern().numberOfValidStripTIDHits());
01507 
01508   NumberOfTIDLayersPerTrack->Fill(track.hitPattern().stripTIDLayersWithMeasurement());
01509   NumberOfTIDLayersPerTrackVsPhiProfile->Fill(phi,     track.hitPattern().stripTIDLayersWithMeasurement());
01510   NumberOfTIDLayersPerTrackVsEtaProfile->Fill(eta,     track.hitPattern().stripTIDLayersWithMeasurement());
01511 
01512   
01513   //Fill TEC Layers and RecHits
01514   NumberOfTECRecHitsPerTrack->Fill(track.hitPattern().numberOfValidStripTECHits());
01515   NumberOfTECRecHitsPerTrackVsPhiProfile->Fill(phi,     track.hitPattern().numberOfValidStripTECHits());
01516   NumberOfTECRecHitsPerTrackVsEtaProfile->Fill(eta,     track.hitPattern().numberOfValidStripTECHits());
01517 
01518   NumberOfTECLayersPerTrack->Fill(track.hitPattern().stripTECLayersWithMeasurement());
01519   NumberOfTECLayersPerTrackVsPhiProfile->Fill(phi,     track.hitPattern().stripTECLayersWithMeasurement());
01520   NumberOfTECLayersPerTrackVsEtaProfile->Fill(eta,     track.hitPattern().stripTECLayersWithMeasurement());
01521 
01522   //Fill PixBarrel Layers and RecHits
01523   NumberOfPixBarrelRecHitsPerTrack->Fill(track.hitPattern().numberOfValidPixelBarrelHits());
01524   NumberOfPixBarrelRecHitsPerTrackVsPhiProfile->Fill(phi,    track.hitPattern().numberOfValidPixelBarrelHits());
01525   NumberOfPixBarrelRecHitsPerTrackVsEtaProfile->Fill(eta,    track.hitPattern().numberOfValidPixelBarrelHits());
01526 
01527   NumberOfPixBarrelLayersPerTrack->Fill(track.hitPattern().pixelBarrelLayersWithMeasurement());
01528   NumberOfPixBarrelLayersPerTrackVsPhiProfile->Fill(phi,     track.hitPattern().pixelBarrelLayersWithMeasurement());
01529   NumberOfPixBarrelLayersPerTrackVsEtaProfile->Fill(eta,     track.hitPattern().pixelBarrelLayersWithMeasurement());
01530 
01531   //Fill PixEndcap Layers and RecHits
01532   NumberOfPixEndcapRecHitsPerTrack->Fill(track.hitPattern().numberOfValidPixelEndcapHits());
01533   NumberOfPixEndcapRecHitsPerTrackVsPhiProfile->Fill(phi,    track.hitPattern().numberOfValidPixelEndcapHits());
01534   NumberOfPixEndcapRecHitsPerTrackVsEtaProfile->Fill(eta,    track.hitPattern().numberOfValidPixelEndcapHits());
01535 
01536   NumberOfPixEndcapLayersPerTrack->Fill(track.hitPattern().pixelEndcapLayersWithMeasurement());
01537   NumberOfPixEndcapLayersPerTrackVsPhiProfile->Fill(phi,     track.hitPattern().pixelEndcapLayersWithMeasurement());
01538   NumberOfPixEndcapLayersPerTrackVsEtaProfile->Fill(eta,     track.hitPattern().pixelEndcapLayersWithMeasurement());
01539 
01540 }
01541 //
01542 // -- Set Lumi Flag
01543 //
01544 void TrackAnalyzer::setLumiFlag() { 
01545   /*
01546   if (Chi2oNDF) Chi2oNDF->setLumiFlag();
01547   if (NumberOfRecHitsPerTrack) NumberOfRecHitsPerTrack->setLumiFlag();
01548   if (GoodTrackChi2oNDF) GoodTrackChi2oNDF->setLumiFlag();
01549   if (GoodTrackNumberOfRecHitsPerTrack) GoodTrackNumberOfRecHitsPerTrack->setLumiFlag();
01550   */
01551   if ( Chi2oNDF_lumiFlag                         ) Chi2oNDF_lumiFlag                         -> setLumiFlag();
01552   if ( NumberOfRecHitsPerTrack_lumiFlag          ) NumberOfRecHitsPerTrack_lumiFlag          -> setLumiFlag();
01553   if ( GoodTrackChi2oNDF_lumiFlag                ) GoodTrackChi2oNDF_lumiFlag                -> setLumiFlag();
01554   if ( GoodTrackNumberOfRecHitsPerTrack_lumiFlag ) GoodTrackNumberOfRecHitsPerTrack_lumiFlag -> setLumiFlag();
01555 }
01556 //
01557 // -- Apply SoftReset 
01558 //
01559 void TrackAnalyzer::doSoftReset(DQMStore * dqmStore_) {
01560   dqmStore_->softReset(Chi2oNDF);
01561   dqmStore_->softReset(NumberOfRecHitsPerTrack);
01562   dqmStore_->softReset(GoodTrackChi2oNDF);
01563   dqmStore_->softReset(GoodTrackNumberOfRecHitsPerTrack);
01564 }
01565 //
01566 // -- Apply Reset 
01567 //
01568 void TrackAnalyzer::doReset(DQMStore * dqmStore_) {
01569   if ( Chi2oNDF_lumiFlag                         ) Chi2oNDF_lumiFlag                         -> Reset();
01570   if ( NumberOfRecHitsPerTrack_lumiFlag          ) NumberOfRecHitsPerTrack_lumiFlag          -> Reset();
01571   if ( GoodTrackChi2oNDF_lumiFlag                ) GoodTrackChi2oNDF_lumiFlag                -> Reset();
01572   if ( GoodTrackNumberOfRecHitsPerTrack_lumiFlag ) GoodTrackNumberOfRecHitsPerTrack_lumiFlag -> Reset();
01573 }
01574 //
01575 // -- Remove SoftReset
01576 //
01577 void TrackAnalyzer::undoSoftReset(DQMStore * dqmStore_) {
01578   dqmStore_->disableSoftReset(Chi2oNDF);
01579   dqmStore_->disableSoftReset(NumberOfRecHitsPerTrack);
01580   dqmStore_->disableSoftReset(GoodTrackChi2oNDF);
01581   dqmStore_->disableSoftReset(GoodTrackNumberOfRecHitsPerTrack);
01582 }
01583 
01584