CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/DQM/TrackingMonitor/plugins/TrackingMonitor.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2012/08/05 12:52:42 $
00005  *  $Revision: 1.39 $
00006  *  \author Suchandra Dutta , Giorgia Mila
00007  */
00008 
00009 #include "DataFormats/TrackReco/interface/Track.h"
00010 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00011 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00012 #include "DataFormats/TrackCandidate/interface/TrackCandidateCollection.h" 
00013 #include "DataFormats/TrackCandidate/interface/TrackCandidate.h" 
00014 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00015 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00016 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00017 #include "MagneticField/Engine/interface/MagneticField.h"
00018 
00019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00020 #include "FWCore/Utilities/interface/InputTag.h"
00021 #include "DQMServices/Core/interface/DQMStore.h"
00022 #include "DQM/TrackingMonitor/interface/TrackBuildingAnalyzer.h"
00023 #include "DQM/TrackingMonitor/interface/TrackAnalyzer.h"
00024 #include "DQM/TrackingMonitor/interface/VertexMonitor.h"
00025 #include "DQM/TrackingMonitor/plugins/TrackingMonitor.h"
00026 
00027 #include "FWCore/Framework/interface/ESHandle.h"
00028 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h" 
00029 
00030 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00031 #include "TrackingTools/Records/interface/TransientRecHitRecord.h" 
00032 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
00033 #include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
00034 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00035 #include "CommonTools/TriggerUtils/interface/GenericTriggerEventFlag.h"
00036 #include "DataFormats/Common/interface/DetSetVectorNew.h"
00037 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
00038 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
00039 
00040 #include "DataFormats/VertexReco/interface/Vertex.h"
00041 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00042 
00043 #include "DQM/TrackingMonitor/interface/GetLumi.h"
00044 
00045 #include <string>
00046 
00047 // TrackingMonitor 
00048 // ----------------------------------------------------------------------------------//
00049 
00050 TrackingMonitor::TrackingMonitor(const edm::ParameterSet& iConfig) 
00051     : dqmStore_( edm::Service<DQMStore>().operator->() )
00052     , conf_ ( iConfig )
00053     , theTrackAnalyzer( new TrackAnalyzer(conf_) )
00054     , theTrackBuildingAnalyzer( new TrackBuildingAnalyzer(conf_) )
00055     , NumberOfTracks(NULL)
00056     , NumberOfMeanRecHitsPerTrack(NULL)
00057     , NumberOfMeanLayersPerTrack(NULL)
00058     , NumberOfGoodTracks(NULL)
00059     , FractionOfGoodTracks(NULL)
00060     , NumberOfSeeds(NULL)
00061     , NumberOfTrackCandidates(NULL)
00062                                 // ADD by Mia
00063                                 /*
00064     , NumberOfPixelClus(NULL)
00065     , NumberOfStripClus(NULL)
00066     , RatioOfPixelAndStripClus(NULL)
00067     , NumberOfTrkVsClus(NULL)
00068     , NumberOfTrkVsStripClus(NULL)
00069     , NumberOfTrkVsPixelClus(NULL)
00070                                 */
00071     , NumberOfGoodTrkVsClus(NULL)
00072     , NumberOfTracksVsLS(NULL)
00073     , NumberOfGoodTracksVsLS(NULL)
00074     , GoodTracksFractionVsLS(NULL)
00075     , GoodTracksNumberOfRecHitsPerTrackVsLS(NULL)
00076                                 // ADD by Mia for PU monitoring
00077                                 // vertex plots to be moved in ad hoc class
00078     , NumberOfTracksVsGoodPVtx(NULL)
00079     , NumberOfTracksVsBXlumi(NULL)
00080     , NumberOfGoodTracksVsGoodPVtx(NULL)
00081     , NumberOfGoodTracksVsBXlumi(NULL)
00082     , FractionOfGoodTracksVsGoodPVtx(NULL)
00083     , FractionOfGoodTracksVsBXlumi(NULL)
00084                                 // ADD by Mia in order to deal with LS transitions
00085     , NumberOfTracks_lumiFlag(NULL)
00086     , NumberOfGoodTracks_lumiFlag(NULL)
00087 
00088     , builderName              ( conf_.getParameter<std::string>("TTRHBuilder"))
00089     , doTrackerSpecific_       ( conf_.getParameter<bool>("doTrackerSpecific") )
00090     , doLumiAnalysis           ( conf_.getParameter<bool>("doLumiAnalysis"))
00091     , doProfilesVsLS_          ( conf_.getParameter<bool>("doProfilesVsLS"))
00092     , doAllPlots               ( conf_.getParameter<bool>("doAllPlots"))
00093     , doGeneralPropertiesPlots_( conf_.getParameter<bool>("doGeneralPropertiesPlots"))
00094     , doHitPropertiesPlots_    ( conf_.getParameter<bool>("doHitPropertiesPlots"))
00095     , doGoodTrackPlots_        ( conf_.getParameter<bool>("doGoodTrackPlots") )
00096     , doPUmonitoring_          ( conf_.getParameter<bool>("doPUmonitoring") )
00097     , genTriggerEventFlag_(new GenericTriggerEventFlag(iConfig))
00098 {
00099 
00100   if ( doPUmonitoring_ ) {
00101 
00102     // get flag from the configuration
00103     doPlotsVsBXlumi_   = conf_.getParameter<bool>("doPlotsVsBXlumi");
00104     doPlotsVsGoodPVtx_ = conf_.getParameter<bool>("doPlotsVsGoodPVtx");
00105 
00106     if ( doPlotsVsBXlumi_ )
00107       theLumiDetails_ = new GetLumi( iConfig.getParameter<edm::ParameterSet>("BXlumiSetup") );
00108 
00109     std::vector<edm::InputTag> primaryVertexInputTags    = conf_.getParameter<std::vector<edm::InputTag> >("primaryVertexInputTags");
00110     std::vector<edm::InputTag> selPrimaryVertexInputTags = conf_.getParameter<std::vector<edm::InputTag> >("selPrimaryVertexInputTags");
00111     std::vector<std::string>   pvLabels                  = conf_.getParameter<std::vector<std::string> >  ("pvLabels");
00112 
00113     if (primaryVertexInputTags.size()==pvLabels.size() and primaryVertexInputTags.size()==selPrimaryVertexInputTags.size()) {
00114       for (size_t i=0; i<primaryVertexInputTags.size(); i++) {
00115         edm::InputTag iPVinputTag    = primaryVertexInputTags[i];
00116         edm::InputTag iSelPVinputTag = selPrimaryVertexInputTags[i];
00117         std::string   iPVlabel       = pvLabels[i];
00118         
00119         theVertexMonitor.push_back(new VertexMonitor(conf_,iPVinputTag,iSelPVinputTag,iPVlabel));
00120       }
00121     }
00122   }
00123 
00124 }
00125 
00126 
00127 TrackingMonitor::~TrackingMonitor() 
00128 {
00129   if (theTrackAnalyzer)          delete theTrackAnalyzer;
00130   if (theTrackBuildingAnalyzer)  delete theTrackBuildingAnalyzer;
00131   if ( doPUmonitoring_ )
00132     for (size_t i=0; i<theVertexMonitor.size(); i++)
00133       if (theVertexMonitor[i]) delete theVertexMonitor[i];
00134   if (genTriggerEventFlag_)      delete genTriggerEventFlag_;
00135 }
00136 
00137 
00138 void TrackingMonitor::beginJob(void) 
00139 {
00140 
00141     // parameters from the configuration
00142     std::string Quality      = conf_.getParameter<std::string>("Quality");
00143     std::string AlgoName     = conf_.getParameter<std::string>("AlgoName");
00144     std::string MEFolderName = conf_.getParameter<std::string>("FolderName"); 
00145 
00146     // test for the Quality veriable validity
00147     if( Quality != "")
00148     {
00149         if( Quality != "highPurity" && Quality != "tight" && Quality != "loose") 
00150         {
00151             edm::LogWarning("TrackingMonitor")  << "Qualty Name is invalid, using no quality criterea by default";
00152             Quality = "";
00153         }
00154     }
00155 
00156     // use the AlgoName and Quality Name
00157     std::string CategoryName = Quality != "" ? AlgoName + "_" + Quality : AlgoName;
00158 
00159     // get binning from the configuration
00160     int    TKNoBin     = conf_.getParameter<int>(   "TkSizeBin");
00161     double TKNoMin     = conf_.getParameter<double>("TkSizeMin");
00162     double TKNoMax     = conf_.getParameter<double>("TkSizeMax");
00163 
00164     int    TCNoBin     = conf_.getParameter<int>(   "TCSizeBin");
00165     double TCNoMin     = conf_.getParameter<double>("TCSizeMin");
00166     double TCNoMax     = conf_.getParameter<double>("TCSizeMax");
00167 
00168     int    TKNoSeedBin = conf_.getParameter<int>(   "TkSeedSizeBin");
00169     double TKNoSeedMin = conf_.getParameter<double>("TkSeedSizeMin");
00170     double TKNoSeedMax = conf_.getParameter<double>("TkSeedSizeMax");
00171 
00172     int    MeanHitBin  = conf_.getParameter<int>(   "MeanHitBin");
00173     double MeanHitMin  = conf_.getParameter<double>("MeanHitMin");
00174     double MeanHitMax  = conf_.getParameter<double>("MeanHitMax");
00175 
00176     int    MeanLayBin  = conf_.getParameter<int>(   "MeanLayBin");
00177     double MeanLayMin  = conf_.getParameter<double>("MeanLayMin");
00178     double MeanLayMax  = conf_.getParameter<double>("MeanLayMax");
00179 
00180     int LSBin = conf_.getParameter<int>(   "LSBin");
00181     int LSMin = conf_.getParameter<double>("LSMin");
00182     int LSMax = conf_.getParameter<double>("LSMax");
00183 
00184     std::string StateName = conf_.getParameter<std::string>("MeasurementState");
00185     if
00186     (
00187         StateName != "OuterSurface" &&
00188         StateName != "InnerSurface" &&
00189         StateName != "ImpactPoint"  &&
00190         StateName != "default"      &&
00191         StateName != "All"
00192     )
00193     {
00194         // print warning
00195         edm::LogWarning("TrackingMonitor")  << "State Name is invalid, using 'ImpactPoint' by default";
00196     }
00197 
00198     dqmStore_->setCurrentFolder(MEFolderName);
00199 
00200     // book the General Property histograms
00201     // ---------------------------------------------------------------------------------//
00202 
00203     if (doGeneralPropertiesPlots_ || doAllPlots){
00204      
00205       dqmStore_->setCurrentFolder(MEFolderName+"/GeneralProperties");
00206 
00207       histname = "NumberOfTracks_" + CategoryName;
00208       // MODIFY by Mia in order to cope w/ high multiplicity
00209       //      NumberOfTracks = dqmStore_->book1D(histname, histname, TKNoBin, TKNoMin, TKNoMax);
00210       NumberOfTracks = dqmStore_->book1D(histname, histname, 3*TKNoBin, TKNoMin, (TKNoMax+0.5)*3.-0.5);
00211       NumberOfTracks->setAxisTitle("Number of Tracks per Event", 1);
00212       NumberOfTracks->setAxisTitle("Number of Events", 2);
00213       
00214       histname = "NumberOfMeanRecHitsPerTrack_" + CategoryName;
00215       NumberOfMeanRecHitsPerTrack = dqmStore_->book1D(histname, histname, MeanHitBin, MeanHitMin, MeanHitMax);
00216       NumberOfMeanRecHitsPerTrack->setAxisTitle("Mean number of found RecHits per Track", 1);
00217       NumberOfMeanRecHitsPerTrack->setAxisTitle("Entries", 2);
00218       
00219       histname = "NumberOfMeanLayersPerTrack_" + CategoryName;
00220       NumberOfMeanLayersPerTrack = dqmStore_->book1D(histname, histname, MeanLayBin, MeanLayMin, MeanLayMax);
00221       NumberOfMeanLayersPerTrack->setAxisTitle("Mean number of Layers per Track", 1);
00222       NumberOfMeanLayersPerTrack->setAxisTitle("Entries", 2);
00223       
00224       dqmStore_->setCurrentFolder(MEFolderName+"/GeneralProperties/GoodTracks");
00225 
00226       histname = "NumberOfGoodTracks_" + CategoryName;
00227       NumberOfGoodTracks = dqmStore_->book1D(histname, histname, TKNoBin, TKNoMin, TKNoMax);
00228       NumberOfGoodTracks->setAxisTitle("Number of Good Tracks per Event", 1);
00229       NumberOfGoodTracks->setAxisTitle("Number of Events", 2);
00230       
00231       histname = "FractionOfGoodTracks_" + CategoryName;
00232       FractionOfGoodTracks = dqmStore_->book1D(histname, histname, 101, -0.005, 1.005);
00233       FractionOfGoodTracks->setAxisTitle("Fraction of High Purity Tracks (Tracks with Pt>1GeV)", 1);
00234       FractionOfGoodTracks->setAxisTitle("Entries", 2);
00235     }
00236 
00237     if ( doLumiAnalysis ) {
00238       // add by Mia in order to deal with LS transitions
00239       dqmStore_->setCurrentFolder(MEFolderName+"/LSanalysis");
00240 
00241       histname = "NumberOfTracks_lumiFlag_" + CategoryName;
00242       NumberOfTracks_lumiFlag = dqmStore_->book1D(histname, histname, 3*TKNoBin, TKNoMin, (TKNoMax+0.5)*3.-0.5);
00243       NumberOfTracks_lumiFlag->setAxisTitle("Number of Tracks per Event", 1);
00244       NumberOfTracks_lumiFlag->setAxisTitle("Number of Events", 2);
00245 
00246       if ( doGoodTrackPlots_ ) {
00247         histname = "NumberOfGoodTracks_lumiFlag_" + CategoryName;
00248         NumberOfGoodTracks_lumiFlag = dqmStore_->book1D(histname, histname, TKNoBin, TKNoMin, TKNoMax);
00249         NumberOfGoodTracks_lumiFlag->setAxisTitle("Number of Good Tracks per Event", 1);
00250         NumberOfGoodTracks_lumiFlag->setAxisTitle("Number of Events", 2);
00251       }
00252       
00253     }
00254 
00255     // book profile plots vs LS :  
00256     //---------------------------  
00257    
00258  
00259     if ( doProfilesVsLS_ || doAllPlots) {
00260 
00261       dqmStore_->setCurrentFolder(MEFolderName+"/GeneralProperties");
00262 
00263       histname = "NumberOfTracksVsLS_"+ CategoryName;
00264       NumberOfTracksVsLS = dqmStore_->bookProfile(histname,histname, LSBin,LSMin,LSMax, TKNoMin, (TKNoMax+0.5)*3.-0.5,"");
00265       NumberOfTracksVsLS->getTH1()->SetBit(TH1::kCanRebin);
00266       NumberOfTracksVsLS->setAxisTitle("#Lumi section",1);
00267       NumberOfTracksVsLS->setAxisTitle("Number of  Tracks",2);
00268 
00269       if (doGoodTrackPlots_ || doAllPlots ){
00270         dqmStore_->setCurrentFolder(MEFolderName+"/GeneralProperties/GoodTracks");
00271 
00272         histname = "NumberOfGoodTracksVsLS_"+ CategoryName;
00273         NumberOfGoodTracksVsLS = dqmStore_->bookProfile(histname,histname, LSBin,LSMin,LSMax, TKNoMin, TKNoMax,"");
00274         NumberOfGoodTracksVsLS->getTH1()->SetBit(TH1::kCanRebin);
00275         NumberOfGoodTracksVsLS->setAxisTitle("#Lumi section",1);
00276         NumberOfGoodTracksVsLS->setAxisTitle("Number of Good Tracks",2);
00277         
00278         histname = "GoodTracksFractionVsLS_"+ CategoryName;
00279         GoodTracksFractionVsLS = dqmStore_->bookProfile(histname,histname, LSBin,LSMin,LSMax,0,1.1,"");
00280         GoodTracksFractionVsLS->getTH1()->SetBit(TH1::kCanRebin);
00281         GoodTracksFractionVsLS->setAxisTitle("#Lumi section",1);
00282         GoodTracksFractionVsLS->setAxisTitle("Fraction of Good Tracks",2);
00283         
00284         histname = "GoodTracksNumberOfRecHitsPerTrackVsLS_" + CategoryName;
00285         GoodTracksNumberOfRecHitsPerTrackVsLS = dqmStore_->bookProfile(histname,histname, LSBin,LSMin,LSMax,0.,40.,"");
00286         GoodTracksNumberOfRecHitsPerTrackVsLS->getTH1()->SetBit(TH1::kCanRebin);
00287         GoodTracksNumberOfRecHitsPerTrackVsLS->setAxisTitle("#Lumi section",1);
00288         GoodTracksNumberOfRecHitsPerTrackVsLS->setAxisTitle("Mean number of RecHits per Good track",2);
00289       }
00290     }
00291 
00292     // book PU monitoring plots :  
00293     //---------------------------  
00294    
00295     if ( doPUmonitoring_ ) {
00296 
00297       for (size_t i=0; i<theVertexMonitor.size(); i++)
00298         theVertexMonitor[i]->beginJob(dqmStore_);
00299 
00300       dqmStore_->setCurrentFolder(MEFolderName+"/PUmonitoring");
00301       
00302       if ( doPlotsVsGoodPVtx_ ) {
00303         // get binning from the configuration
00304         int    GoodPVtxBin   = conf_.getParameter<int>("GoodPVtxBin");
00305         double GoodPVtxMin   = conf_.getParameter<double>("GoodPVtxMin");
00306         double GoodPVtxMax   = conf_.getParameter<double>("GoodPVtxMax");
00307 
00308         histname = "NumberOfTracksVsGoodPVtx";
00309         NumberOfTracksVsGoodPVtx = dqmStore_->bookProfile(histname,histname,GoodPVtxBin,GoodPVtxMin,GoodPVtxMax,TKNoMin, (TKNoMax+0.5)*3.-0.5,"");
00310         NumberOfTracksVsGoodPVtx->getTH1()->SetBit(TH1::kCanRebin);
00311         NumberOfTracksVsGoodPVtx->setAxisTitle("Number of PV",1);
00312         NumberOfTracksVsGoodPVtx->setAxisTitle("Mean number of Tracks per Event",2);
00313       
00314         histname = "NumberOfGoodTracksVsGoodPVtx";
00315         NumberOfGoodTracksVsGoodPVtx = dqmStore_->bookProfile(histname,histname,GoodPVtxBin,GoodPVtxMin,GoodPVtxMax, TKNoMin, TKNoMax,"");
00316         NumberOfGoodTracksVsGoodPVtx->getTH1()->SetBit(TH1::kCanRebin);
00317         NumberOfGoodTracksVsGoodPVtx->setAxisTitle("Number of PV",1);
00318         NumberOfGoodTracksVsGoodPVtx->setAxisTitle("Mean number of Good Tracks per Event",2);
00319 
00320         histname = "FractionOfGoodTracksVsGoodPVtx";
00321         FractionOfGoodTracksVsGoodPVtx = dqmStore_->bookProfile(histname,histname,GoodPVtxBin,GoodPVtxMin,GoodPVtxMax, TKNoMin, TKNoMax,"");
00322         FractionOfGoodTracksVsGoodPVtx->getTH1()->SetBit(TH1::kCanRebin);
00323         FractionOfGoodTracksVsGoodPVtx->setAxisTitle("Number of PV",1);
00324         FractionOfGoodTracksVsGoodPVtx->setAxisTitle("Mean fraction of Good Tracks per Event",2);
00325       }
00326       
00327       if ( doPlotsVsBXlumi_ ) {
00328         // get binning from the configuration
00329         edm::ParameterSet BXlumiParameters = conf_.getParameter<edm::ParameterSet>("BXlumiSetup");
00330         int    BXlumiBin   = BXlumiParameters.getParameter<int>("BXlumiBin");
00331         double BXlumiMin   = BXlumiParameters.getParameter<double>("BXlumiMin");
00332         double BXlumiMax   = BXlumiParameters.getParameter<double>("BXlumiMax");
00333       
00334         histname = "NumberOfTracksVsBXlumi_"+ CategoryName;
00335         NumberOfTracksVsBXlumi = dqmStore_->bookProfile(histname,histname, BXlumiBin,BXlumiMin,BXlumiMax, TKNoMin, TKNoMax,"");
00336         NumberOfTracksVsBXlumi->getTH1()->SetBit(TH1::kCanRebin);
00337         NumberOfTracksVsBXlumi->setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00338         NumberOfTracksVsBXlumi->setAxisTitle("Mean number of Good Tracks",2);
00339 
00340         histname = "NumberOfGoodTracksVsBXlumi_"+ CategoryName;
00341         NumberOfGoodTracksVsBXlumi = dqmStore_->bookProfile(histname,histname, BXlumiBin,BXlumiMin,BXlumiMax, TKNoMin, TKNoMax,"");
00342         NumberOfGoodTracksVsBXlumi->getTH1()->SetBit(TH1::kCanRebin);
00343         NumberOfGoodTracksVsBXlumi->setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00344         NumberOfGoodTracksVsBXlumi->setAxisTitle("Mean number of Good Tracks",2);
00345       
00346         histname = "FractionOfGoodTracksVsBXlumi_"+ CategoryName;
00347         FractionOfGoodTracksVsBXlumi = dqmStore_->bookProfile(histname,histname, BXlumiBin,BXlumiMin,BXlumiMax, TKNoMin, TKNoMax,"");
00348         FractionOfGoodTracksVsBXlumi->getTH1()->SetBit(TH1::kCanRebin);
00349         FractionOfGoodTracksVsBXlumi->setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00350         FractionOfGoodTracksVsBXlumi->setAxisTitle("Mean fraction of Good Tracks",2);
00351       
00352       }
00353     }
00354 
00355     theTrackAnalyzer->beginJob(dqmStore_);
00356 
00357     // book the Seed Property histograms
00358     // ---------------------------------------------------------------------------------//
00359 
00360     dqmStore_->setCurrentFolder(MEFolderName+"/TrackBuilding");
00361 
00362     doAllSeedPlots=conf_.getParameter<bool>("doSeedParameterHistos");
00363     doSeedNumberPlot=conf_.getParameter<bool>("doSeedNumberHisto");
00364     doSeedVsClusterPlot=conf_.getParameter<bool>("doSeedVsClusterHisto");
00365     //    if (doAllPlots) doAllSeedPlots=true;
00366 
00367     runTrackBuildingAnalyzerForSeed=(doAllSeedPlots || conf_.getParameter<bool>("doSeedPTHisto") ||conf_.getParameter<bool>("doSeedETAHisto") || conf_.getParameter<bool>("doSeedPHIHisto") || conf_.getParameter<bool>("doSeedPHIVsETAHisto") || conf_.getParameter<bool>("doSeedThetaHisto") || conf_.getParameter<bool>("doSeedQHisto") || conf_.getParameter<bool>("doSeedDxyHisto") || conf_.getParameter<bool>("doSeedDzHisto") || conf_.getParameter<bool>("doSeedNRecHitsHisto") || conf_.getParameter<bool>("doSeedNVsPhiProf")|| conf_.getParameter<bool>("doSeedNVsEtaProf"));
00368 
00369     edm::InputTag seedProducer   = conf_.getParameter<edm::InputTag>("SeedProducer");
00370 
00371     if (doAllSeedPlots || doSeedNumberPlot){
00372       histname = "NumberOfSeeds_"+ seedProducer.label() + "_"+ CategoryName;
00373       NumberOfSeeds = dqmStore_->book1D(histname, histname, TKNoSeedBin, TKNoSeedMin, TKNoSeedMax);
00374       NumberOfSeeds->setAxisTitle("Number of Seeds per Event", 1);
00375       NumberOfSeeds->setAxisTitle("Number of Events", 2);
00376     }
00377 
00378     if (doAllSeedPlots || doSeedVsClusterPlot){
00379 
00380       ClusterLabels=  conf_.getParameter<std::vector<std::string> >("ClusterLabels");
00381 
00382       std::vector<double> histoMin,histoMax;
00383       std::vector<int> histoBin; //these vectors are for max min and nbins in histograms 
00384 
00385       int    NClusPxBin  = conf_.getParameter<int>(   "NClusPxBin");
00386       double NClusPxMin  = conf_.getParameter<double>("NClusPxMin");
00387       double NClusPxMax = conf_.getParameter<double>("NClusPxMax");
00388       
00389       int    NClusStrBin = conf_.getParameter<int>(   "NClusStrBin");
00390       double NClusStrMin = conf_.getParameter<double>("NClusStrMin");
00391       double NClusStrMax = conf_.getParameter<double>("NClusStrMax");
00392 
00393       setMaxMinBin(histoMin,histoMax,histoBin,NClusStrMin,NClusStrMax,NClusStrBin,NClusPxMin,NClusPxMax,NClusPxBin);
00394      
00395       for (uint i=0; i<ClusterLabels.size(); i++){
00396         histname = "SeedsVsClusters_" + seedProducer.label() + "_Vs_" + ClusterLabels[i] + "_" + CategoryName;
00397         SeedsVsClusters.push_back(dynamic_cast<MonitorElement*>(dqmStore_->book2D(histname, histname, histoBin[i], histoMin[i], histoMax[i],
00398                                                                                   TKNoSeedBin, TKNoSeedMin, TKNoSeedMax)));
00399         SeedsVsClusters[i]->setAxisTitle("Number of Clusters", 1);
00400         SeedsVsClusters[i]->setAxisTitle("Number of Seeds", 2);
00401       }
00402     }
00403     
00404 
00405     doTkCandPlots=conf_.getParameter<bool>("doTrackCandHistos");
00406     //    if (doAllPlots) doTkCandPlots=true;
00407 
00408     if (doTkCandPlots){
00409 
00410       edm::InputTag tcProducer     = conf_.getParameter<edm::InputTag>("TCProducer");
00411 
00412       histname = "NumberOfTrackCandidates_"+ tcProducer.label() + "_"+ CategoryName;
00413       NumberOfTrackCandidates = dqmStore_->book1D(histname, histname, TCNoBin, TCNoMin, TCNoMax);
00414       NumberOfTrackCandidates->setAxisTitle("Number of Track Candidates per Event", 1);
00415       NumberOfTrackCandidates->setAxisTitle("Number of Event", 2);
00416     }
00417 
00418    
00419     theTrackBuildingAnalyzer->beginJob(dqmStore_);
00420     
00421 
00422     if (doLumiAnalysis) {
00423 //      if (NumberOfTracks) NumberOfTracks->setLumiFlag();
00424 //      if (NumberOfGoodTracks) NumberOfGoodTracks->setLumiFlag();
00425 //      if (FractionOfGoodTracks) FractionOfGoodTracks->setLumiFlag();
00426       if ( NumberOfTracks_lumiFlag       ) NumberOfTracks_lumiFlag       -> setLumiFlag();
00427       if ( NumberOfGoodTracks_lumiFlag   ) NumberOfGoodTracks_lumiFlag   -> setLumiFlag();
00428       theTrackAnalyzer->setLumiFlag();    
00429     }
00430 
00431     if (doTrackerSpecific_ || doAllPlots) {
00432 
00433       ClusterLabels=  conf_.getParameter<std::vector<std::string> >("ClusterLabels");
00434 
00435       std::vector<double> histoMin,histoMax;
00436       std::vector<int> histoBin; //these vectors are for max min and nbins in histograms 
00437 
00438       /*
00439       int    NClusPxBin  = conf_.getParameter<int>(   "NClusPxBin");
00440       double NClusPxMin  = conf_.getParameter<double>("NClusPxMin");
00441       double NClusPxMax = conf_.getParameter<double>("NClusPxMax");
00442       
00443       int    NClusStrBin = conf_.getParameter<int>(   "NClusStrBin");
00444       double NClusStrMin = conf_.getParameter<double>("NClusStrMin");
00445       double NClusStrMax = conf_.getParameter<double>("NClusStrMax");
00446       */
00447 
00448       /*
00449       int    NClus2DTotBin = conf_.getParameter<int>(   "NClus2DTotBin");
00450       double NClus2DTotMin = conf_.getParameter<double>("NClus2DTotMin");
00451       double NClus2DTotMax = conf_.getParameter<double>("NClus2DTotMax");
00452       */
00453 
00454       int    NClusStrBin = conf_.getParameter<int>(   "NClusStrBin");
00455       double NClusStrMin = conf_.getParameter<double>("NClusStrMin");
00456       double NClusStrMax = conf_.getParameter<double>("NClusStrMax");
00457 
00458       int    NClusPxBin = conf_.getParameter<int>(   "NClusPxBin");
00459       double NClusPxMin = conf_.getParameter<double>("NClusPxMin");
00460       double NClusPxMax = conf_.getParameter<double>("NClusPxMax");
00461 
00462       int    NTrk2DBin     = conf_.getParameter<int>(   "NTrk2DBin");
00463       double NTrk2DMin     = conf_.getParameter<double>("NTrk2DMin");
00464       double NTrk2DMax     = conf_.getParameter<double>("NTrk2DMax");
00465 
00466       //      setMaxMinBin(histoMin,histoMax,histoBin,NClusStrMin,NClusStrMax,NClusStrBin,NClusPxMin,NClusPxMax,NClusPxBin);
00467       setMaxMinBin(histoMin,histoMax,histoBin,
00468                    NClusStrMin,NClusStrMax,NClusStrBin,
00469                    NClusPxMin,  NClusPxMax,  NClusPxBin);
00470      
00471       /*
00472       int    NClusPxBin  = conf_.getParameter<int>(   "NClusPxBin");
00473       double NClusPxMin  = conf_.getParameter<double>("NClusPxMin");
00474       double NClusPxMax = conf_.getParameter<double>("NClusPxMax");
00475 
00476       int    NClusStrBin = conf_.getParameter<int>(   "NClusStrBin");
00477       double NClusStrMin = conf_.getParameter<double>("NClusStrMin");
00478       double NClusStrMax = conf_.getParameter<double>("NClusStrMax");
00479       */
00480 
00481       dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties");
00482 
00483       for (uint i=0; i<ClusterLabels.size(); i++){
00484 
00485         dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties");
00486         histname = "TracksVs" + ClusterLabels[i] + "Cluster_" + CategoryName;
00487         NumberOfTrkVsClusters.push_back(dynamic_cast<MonitorElement*>(dqmStore_->book2D(histname, histname,
00488                                                                                         histoBin[i], histoMin[i], histoMax[i],
00489                                                                                         NTrk2DBin,NTrk2DMin,NTrk2DMax
00490                                                                                         )));
00491         std::string title = "Number of " + ClusterLabels[i] + " Clusters";
00492         NumberOfTrkVsClusters[i]->setAxisTitle(title, 1);
00493         NumberOfTrkVsClusters[i]->setAxisTitle("Number of Seeds", 2);
00494 
00495         if (doGoodTrackPlots_ || doAllPlots ) {
00496 
00497           dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties/GoodTracks");
00498 
00499           if(ClusterLabels[i].compare("Tot")==0){
00500             histname = "GoodTracksVs" + ClusterLabels[i] + "Cluster_" + CategoryName; 
00501             NumberOfGoodTrkVsClus = dqmStore_->book2D(histname,histname,
00502                                                       histoBin[i], histoMin[i], histoMax[i],
00503                                                       TKNoBin,TKNoMin,TKNoMax
00504                                                       );
00505             NumberOfGoodTrkVsClus->setAxisTitle("# of Clusters in (Pixel+Strip) Detectors", 1);
00506             NumberOfGoodTrkVsClus->setAxisTitle("Number of Good Tracks", 2);
00507           }
00508         }
00509       }
00510 
00511       /*
00512       histname = "TracksVsClusters_" + CategoryName; 
00513       NumberOfTrkVsClus = dqmStore_->book2D(histname,histname,
00514                                             NClus2DTotBin,NClus2DTotMin,NClus2DTotMax,
00515                                             NTrk2DBin,NTrk2DMin,NTrk2DMax
00516                                             );
00517       NumberOfTrkVsClus->setAxisTitle("# of Clusters in (Pixel+Strip) Detectors", 1);
00518       NumberOfTrkVsClus->setAxisTitle("Number of Tracks", 2);
00519       */
00520     }
00521 
00522     
00523 }
00524 
00525 // -- BeginRun
00526 //---------------------------------------------------------------------------------//
00527 void TrackingMonitor::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup)
00528 {
00529  
00530   // Initialize the GenericTriggerEventFlag
00531   if ( genTriggerEventFlag_->on() ) genTriggerEventFlag_->initRun( iRun, iSetup );
00532 }
00533 
00534 // - BeginLumi
00535 // ---------------------------------------------------------------------------------//
00536 void TrackingMonitor::beginLuminosityBlock(const edm::LuminosityBlock& lumi, const edm::EventSetup&  eSetup) {
00537 
00538   if (doLumiAnalysis) {
00539 //    dqmStore_->softReset(NumberOfTracks);
00540 //    dqmStore_->softReset(NumberOfGoodTracks);
00541 //    dqmStore_->softReset(FractionOfGoodTracks);
00542 //    theTrackAnalyzer->doSoftReset(dqmStore_);    
00543     if ( NumberOfTracks_lumiFlag       ) NumberOfTracks_lumiFlag       -> Reset();
00544     if ( doGoodTrackPlots_ )
00545       if ( NumberOfGoodTracks_lumiFlag   ) NumberOfGoodTracks_lumiFlag   -> Reset();
00546     theTrackAnalyzer->doReset(dqmStore_);    
00547   }
00548 }
00549 
00550 // -- Analyse
00551 // ---------------------------------------------------------------------------------//
00552 void TrackingMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) 
00553 {
00554     // Filter out events if Trigger Filtering is requested
00555     if (genTriggerEventFlag_->on()&& ! genTriggerEventFlag_->accept( iEvent, iSetup) ) return;
00556 
00557     // input tags for collections from the configuration
00558     edm::InputTag trackProducer  = conf_.getParameter<edm::InputTag>("TrackProducer");
00559     edm::InputTag seedProducer   = conf_.getParameter<edm::InputTag>("SeedProducer");
00560     edm::InputTag tcProducer     = conf_.getParameter<edm::InputTag>("TCProducer");
00561     edm::InputTag bsSrc          = conf_.getParameter<edm::InputTag>("beamSpot");
00562     edm::InputTag pvSrc_         = conf_.getParameter<edm::InputTag>("primaryVertex");
00563     std::string Quality = conf_.getParameter<std::string>("Quality");
00564     std::string Algo    = conf_.getParameter<std::string>("AlgoName");
00565 
00566     //  Analyse the tracks
00567     //  if the collection is empty, do not fill anything
00568     // ---------------------------------------------------------------------------------//
00569 
00570     // get the track collection
00571     edm::Handle<reco::TrackCollection> trackHandle;
00572     iEvent.getByLabel(trackProducer, trackHandle);
00573 
00574     if (trackHandle.isValid()) 
00575     {
00576 
00577        reco::TrackCollection trackCollection = *trackHandle;
00578         // calculate the mean # rechits and layers
00579         int totalNumTracks = 0, totalRecHits = 0, totalLayers = 0;
00580         int totalNumHPTracks = 0, totalNumPt1Tracks = 0, totalNumHPPt1Tracks = 0;
00581         double frac = 0.;
00582         for (reco::TrackCollection::const_iterator track = trackCollection.begin(); track!=trackCollection.end(); ++track) 
00583         {
00584 
00585             if( track->quality(reco::TrackBase::highPurity) ) {
00586               ++totalNumHPTracks;
00587               if ( track->pt() >= 1. ) {
00588                 ++totalNumHPPt1Tracks;
00589                 if ( doProfilesVsLS_ || doAllPlots)
00590                   if ( doGoodTrackPlots_ || doAllPlots ) {
00591                     GoodTracksNumberOfRecHitsPerTrackVsLS->Fill(static_cast<double>(iEvent.id().luminosityBlock()),track->recHitsSize());
00592                   }
00593               }
00594             }
00595             
00596             if ( track->pt() >= 1. ) ++totalNumPt1Tracks;
00597         
00598 
00599             if( Quality == "highPurity") 
00600             {
00601                 if( !track->quality(reco::TrackBase::highPurity) ) continue;
00602             }
00603             else if( Quality == "tight") 
00604             {
00605                 if( !track->quality(reco::TrackBase::tight) ) continue;
00606             }
00607             else if( Quality == "loose") 
00608             {
00609                 if( !track->quality(reco::TrackBase::loose) ) continue;
00610             }
00611             
00612             totalNumTracks++;
00613             totalRecHits    += track->numberOfValidHits();
00614             totalLayers     += track->hitPattern().trackerLayersWithMeasurement();
00615 
00616             // do analysis per track
00617             theTrackAnalyzer->analyze(iEvent, iSetup, *track);
00618         }
00619 
00620         if (totalNumPt1Tracks > 0) frac = static_cast<double>(totalNumHPPt1Tracks)/static_cast<double>(totalNumPt1Tracks);
00621 
00622         if (doGeneralPropertiesPlots_ || doAllPlots){
00623           NumberOfTracks       -> Fill(totalNumTracks);
00624           NumberOfGoodTracks   -> Fill(totalNumHPPt1Tracks);
00625           FractionOfGoodTracks -> Fill(frac);
00626         }
00627 
00628         if ( doProfilesVsLS_ || doAllPlots) {
00629           NumberOfTracksVsLS->Fill(static_cast<double>(iEvent.id().luminosityBlock()),totalNumTracks);
00630           if ( doGoodTrackPlots_ || doAllPlots ) {
00631             NumberOfGoodTracksVsLS->Fill(static_cast<double>(iEvent.id().luminosityBlock()),totalNumHPPt1Tracks);
00632             GoodTracksFractionVsLS->Fill(static_cast<double>(iEvent.id().luminosityBlock()),frac);
00633           }
00634         }
00635 
00636         if ( doLumiAnalysis ) {
00637           NumberOfTracks_lumiFlag       -> Fill(totalNumTracks);
00638           if ( doGoodTrackPlots_ )
00639             NumberOfGoodTracks_lumiFlag   -> Fill(totalNumHPPt1Tracks);
00640         }
00641 
00642         if (doGeneralPropertiesPlots_ || doAllPlots){
00643           if( totalNumTracks > 0 )
00644               {
00645                 double meanRecHits = static_cast<double>(totalRecHits) / static_cast<double>(totalNumTracks);
00646                 double meanLayers  = static_cast<double>(totalLayers)  / static_cast<double>(totalNumTracks);
00647                 NumberOfMeanRecHitsPerTrack -> Fill(meanRecHits);
00648                 NumberOfMeanLayersPerTrack  -> Fill(meanLayers);
00649               }
00650           }
00651     
00652         //  Analyse the Track Building variables 
00653         //  if the collection is empty, do not fill anything
00654         // ---------------------------------------------------------------------------------//
00655         
00656                             
00657             // fill the TrackCandidate info
00658             if (doTkCandPlots) 
00659               {
00660         
00661                 // magnetic field
00662                 edm::ESHandle<MagneticField> theMF;
00663                 iSetup.get<IdealMagneticFieldRecord>().get(theMF);  
00664                 
00665                 // get the beam spot
00666                 edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00667                 iEvent.getByLabel(bsSrc,recoBeamSpotHandle);
00668                 const reco::BeamSpot& bs = *recoBeamSpotHandle;      
00669                 
00670                 // get the candidate collection
00671                 edm::Handle<TrackCandidateCollection> theTCHandle;
00672                 iEvent.getByLabel(tcProducer, theTCHandle ); 
00673                 const TrackCandidateCollection& theTCCollection = *theTCHandle;
00674 
00675                 if (theTCHandle.isValid())
00676                   {
00677                     NumberOfTrackCandidates->Fill(theTCCollection.size());
00678                     iSetup.get<TransientRecHitRecord>().get(builderName,theTTRHBuilder);
00679                     for( TrackCandidateCollection::const_iterator cand = theTCCollection.begin(); cand != theTCCollection.end(); ++cand)
00680                       {
00681                         theTrackBuildingAnalyzer->analyze(iEvent, iSetup, *cand, bs, theMF, theTTRHBuilder);
00682                       }
00683                   }
00684                 else
00685                   {
00686                     edm::LogWarning("TrackingMonitor") << "No Track Candidates in the event.  Not filling associated histograms";
00687                   }
00688               }
00689     
00690             //plots for trajectory seeds
00691 
00692             if (doAllSeedPlots || doSeedNumberPlot || doSeedVsClusterPlot || runTrackBuildingAnalyzerForSeed){
00693             
00694             // get the seed collection
00695             edm::Handle<edm::View<TrajectorySeed> > seedHandle;
00696             iEvent.getByLabel(seedProducer, seedHandle);
00697             const edm::View<TrajectorySeed>& seedCollection = *seedHandle;
00698             
00699             // fill the seed info
00700             if (seedHandle.isValid()) 
00701               {
00702                 if(doAllSeedPlots || doSeedNumberPlot) NumberOfSeeds->Fill(seedCollection.size());
00703 
00704                 if(doAllSeedPlots || doSeedVsClusterPlot){
00705 
00706                   std::vector<int> NClus;
00707                   setNclus(iEvent,NClus);
00708                   for (uint  i=0; i< ClusterLabels.size(); i++){
00709                     SeedsVsClusters[i]->Fill(NClus[i],seedCollection.size());
00710                   }
00711                 }
00712 
00713                 if (doAllSeedPlots || runTrackBuildingAnalyzerForSeed){
00714 
00715                   //here duplication of mag field and be informations is needed to allow seed and track cand histos to be independent
00716                   // magnetic field
00717                   edm::ESHandle<MagneticField> theMF;
00718                   iSetup.get<IdealMagneticFieldRecord>().get(theMF);  
00719                   
00720                   // get the beam spot
00721                   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00722                   iEvent.getByLabel(bsSrc,recoBeamSpotHandle);
00723                   const reco::BeamSpot& bs = *recoBeamSpotHandle;      
00724                   
00725                   iSetup.get<TransientRecHitRecord>().get(builderName,theTTRHBuilder);
00726                   for(size_t i=0; i < seedHandle->size(); ++i)
00727                     {
00728                       edm::RefToBase<TrajectorySeed> seed(seedHandle, i);
00729                       theTrackBuildingAnalyzer->analyze(iEvent, iSetup, *seed, bs, theMF, theTTRHBuilder);
00730                     }
00731                 }
00732                 
00733               }
00734             else
00735               {
00736                 edm::LogWarning("TrackingMonitor") << "No Trajectory seeds in the event.  Not filling associated histograms";
00737               }
00738             }
00739 
00740 
00741          if (doTrackerSpecific_ || doAllPlots) 
00742           {
00743             std::vector<int> NClus;
00744             setNclus(iEvent,NClus);
00745             for (uint  i=0; i< ClusterLabels.size(); i++){
00746               NumberOfTrkVsClusters[i]->Fill(NClus[i],totalNumTracks);
00747             }
00748             if ( doGoodTrackPlots_ || doAllPlots ) {
00749               for (uint  i=0; i< ClusterLabels.size(); i++){
00750                 if(ClusterLabels[i].compare("Tot")==0){
00751                   NumberOfGoodTrkVsClus->Fill( NClus[i],totalNumHPPt1Tracks);
00752                 }
00753               }
00754             }
00755 
00756             /*
00757             edm::Handle< edmNew::DetSetVector<SiStripCluster> > strip_clusters;
00758             iEvent.getByLabel("siStripClusters", strip_clusters);
00759             edm::Handle< edmNew::DetSetVector<SiPixelCluster> > pixel_clusters;
00760             iEvent.getByLabel("siPixelClusters", pixel_clusters);
00761             if (strip_clusters.isValid() && pixel_clusters.isValid()) 
00762               {
00763                 unsigned int ncluster_pix   = (*pixel_clusters).dataSize(); 
00764                 unsigned int ncluster_strip = (*strip_clusters).dataSize(); 
00765                 double ratio = 0.0;
00766                 if ( ncluster_pix > 0) ratio = atan(ncluster_pix*1.0/ncluster_strip);
00767 
00768                 NumberOfStripClus->Fill(ncluster_strip);
00769                 NumberOfPixelClus->Fill(ncluster_pix);
00770                 RatioOfPixelAndStripClus->Fill(ratio);
00771 
00772                 NumberOfTrkVsClus->Fill( ncluster_strip+ncluster_pix,totalNumTracks);
00773             
00774                 if (doGoodTrackPlots_) {
00775                   NumberOfGoodTrkVsClus->Fill( ncluster_strip+ncluster_pix,totalNumHPPt1Tracks);
00776                 }
00777               }
00778             */
00779           }
00780 
00781          if ( doPUmonitoring_ ) {
00782            
00783            // do vertex monitoring
00784            for (size_t i=0; i<theVertexMonitor.size(); i++)
00785              theVertexMonitor[i]->analyze(iEvent, iSetup);
00786            
00787            if ( doPlotsVsGoodPVtx_ ) {
00788              
00789              edm::InputTag primaryVertexInputTag = conf_.getParameter<edm::InputTag>("primaryVertex");
00790              
00791              size_t totalNumGoodPV = 0;
00792              edm::Handle< reco::VertexCollection > pvHandle;
00793              iEvent.getByLabel(primaryVertexInputTag, pvHandle );
00794              if (pvHandle.isValid())
00795                {
00796                  
00797                  for (reco::VertexCollection::const_iterator pv = pvHandle->begin();
00798                       pv != pvHandle->end(); ++pv) {
00799                    
00800                    //--- pv fake (the pv collection should have size==1 and the pv==beam spot) 
00801                    if (pv->isFake() || pv->tracksSize()==0) continue;
00802                    
00803                    // definition of goodOfflinePrimaryVertex
00804                    if (pv->ndof() < 4. || pv->z() > 24.)  continue;
00805                    totalNumGoodPV++;
00806                  }
00807                  
00808                  NumberOfTracksVsGoodPVtx       -> Fill( totalNumGoodPV, totalNumTracks      );
00809                  NumberOfGoodTracksVsGoodPVtx   -> Fill( totalNumGoodPV, totalNumHPPt1Tracks );
00810                  FractionOfGoodTracksVsGoodPVtx -> Fill( totalNumGoodPV, frac                );
00811                }
00812            }
00813            
00814            if ( doPlotsVsBXlumi_ ) {
00815              
00816              double bxlumi = theLumiDetails_->getValue(iEvent);
00817              
00818              NumberOfTracksVsBXlumi       -> Fill( bxlumi, totalNumTracks      );
00819              NumberOfGoodTracksVsBXlumi   -> Fill( bxlumi, totalNumHPPt1Tracks );
00820              FractionOfGoodTracksVsBXlumi -> Fill( bxlumi, frac                );
00821            }
00822            
00823          }
00824          
00825     }
00826 
00827 }
00828 
00829 
00830 void TrackingMonitor::endRun(const edm::Run&, const edm::EventSetup&) 
00831 {
00832   if (doLumiAnalysis) {
00833     /*
00834     dqmStore_->disableSoftReset(NumberOfTracks);
00835     dqmStore_->disableSoftReset(NumberOfGoodTracks);
00836     dqmStore_->disableSoftReset(FractionOfGoodTracks);
00837     theTrackAnalyzer->undoSoftReset(dqmStore_);    
00838     */
00839   }
00840   
00841 }
00842 
00843 void TrackingMonitor::endJob(void) 
00844 {
00845     bool outputMEsInRootFile   = conf_.getParameter<bool>("OutputMEsInRootFile");
00846     std::string outputFileName = conf_.getParameter<std::string>("OutputFileName");
00847     if(outputMEsInRootFile)
00848     {
00849         dqmStore_->showDirStructure();
00850         dqmStore_->save(outputFileName);
00851     }
00852 }
00853 
00854 void TrackingMonitor::setMaxMinBin(std::vector<double> &arrayMin,  std::vector<double> &arrayMax, std::vector<int> &arrayBin, double smin, double smax, int sbin, double pmin, double pmax, int pbin) 
00855 {
00856   arrayMin.resize(ClusterLabels.size());
00857   arrayMax.resize(ClusterLabels.size());
00858   arrayBin.resize(ClusterLabels.size());
00859 
00860   for (uint i=0; i<ClusterLabels.size(); ++i){
00861 
00862     if (ClusterLabels[i].compare("Pix")==0) {arrayMin[i]=pmin; arrayMax[i]=pmax; arrayBin[i]=pbin;}
00863     else if(ClusterLabels[i].compare("Strip")==0){arrayMin[i]=smin; arrayMax[i]=smax; arrayBin[i]=sbin;}
00864     else if(ClusterLabels[i].compare("Tot")==0){arrayMin[i]=smin; arrayMax[i]=smax+pmax; arrayBin[i]=sbin;}
00865     else {edm::LogWarning("TrackingMonitor")  << "Cluster Label " << ClusterLabels[i] << " not defined, using strip parameters "; 
00866       arrayMin[i]=smin; arrayMax[i]=smax; arrayBin[i]=sbin;}
00867 
00868   }
00869     
00870 }
00871 
00872 void TrackingMonitor::setNclus(const edm::Event& iEvent,std::vector<int> &arrayNclus) 
00873 {
00874 
00875   int ncluster_pix=-1;
00876   int ncluster_strip=-1;
00877 
00878   edm::Handle< edmNew::DetSetVector<SiStripCluster> > strip_clusters;
00879   iEvent.getByLabel("siStripClusters", strip_clusters);
00880   edm::Handle< edmNew::DetSetVector<SiPixelCluster> > pixel_clusters;
00881   iEvent.getByLabel("siPixelClusters", pixel_clusters);
00882 
00883   if (strip_clusters.isValid() && pixel_clusters.isValid()) 
00884     {
00885                 ncluster_pix   = (*pixel_clusters).dataSize(); 
00886                 ncluster_strip = (*strip_clusters).dataSize(); 
00887     }
00888 
00889   arrayNclus.resize(ClusterLabels.size());
00890   for (uint i=0; i<ClusterLabels.size(); ++i){
00891 
00892     if (ClusterLabels[i].compare("Pix")==0) arrayNclus[i]=ncluster_pix ;
00893     else if(ClusterLabels[i].compare("Strip")==0) arrayNclus[i]=ncluster_strip;
00894     else if(ClusterLabels[i].compare("Tot")==0)arrayNclus[i]=ncluster_pix+ncluster_strip;
00895     else {edm::LogWarning("TrackingMonitor") << "Cluster Label " << ClusterLabels[i] << " not defined using stri parametrs ";
00896       arrayNclus[i]=ncluster_strip ;}
00897   }
00898     
00899 }
00900 DEFINE_FWK_MODULE(TrackingMonitor);