CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/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/04/24 17:42:42 $
00005  *  $Revision: 1.38 $
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       histname = "NumberOfGoodTracks_lumiFlag_" + CategoryName;
00247       NumberOfGoodTracks_lumiFlag = dqmStore_->book1D(histname, histname, TKNoBin, TKNoMin, TKNoMax);
00248       NumberOfGoodTracks_lumiFlag->setAxisTitle("Number of Good Tracks per Event", 1);
00249       NumberOfGoodTracks_lumiFlag->setAxisTitle("Number of Events", 2);
00250       
00251     }
00252 
00253     // book profile plots vs LS :  
00254     //---------------------------  
00255    
00256  
00257     if ( doProfilesVsLS_ || doAllPlots) {
00258 
00259       dqmStore_->setCurrentFolder(MEFolderName+"/GeneralProperties");
00260 
00261       histname = "NumberOfTracksVsLS_"+ CategoryName;
00262       NumberOfTracksVsLS = dqmStore_->bookProfile(histname,histname, LSBin,LSMin,LSMax, TKNoMin, (TKNoMax+0.5)*3.-0.5,"");
00263       NumberOfTracksVsLS->getTH1()->SetBit(TH1::kCanRebin);
00264       NumberOfTracksVsLS->setAxisTitle("#Lumi section",1);
00265       NumberOfTracksVsLS->setAxisTitle("Number of  Tracks",2);
00266 
00267       if (doGoodTrackPlots_ || doAllPlots ){
00268         dqmStore_->setCurrentFolder(MEFolderName+"/GeneralProperties/GoodTracks");
00269 
00270         histname = "NumberOfGoodTracksVsLS_"+ CategoryName;
00271         NumberOfGoodTracksVsLS = dqmStore_->bookProfile(histname,histname, LSBin,LSMin,LSMax, TKNoMin, TKNoMax,"");
00272         NumberOfGoodTracksVsLS->getTH1()->SetBit(TH1::kCanRebin);
00273         NumberOfGoodTracksVsLS->setAxisTitle("#Lumi section",1);
00274         NumberOfGoodTracksVsLS->setAxisTitle("Number of Good Tracks",2);
00275         
00276         histname = "GoodTracksFractionVsLS_"+ CategoryName;
00277         GoodTracksFractionVsLS = dqmStore_->bookProfile(histname,histname, LSBin,LSMin,LSMax,0,1.1,"");
00278         GoodTracksFractionVsLS->getTH1()->SetBit(TH1::kCanRebin);
00279         GoodTracksFractionVsLS->setAxisTitle("#Lumi section",1);
00280         GoodTracksFractionVsLS->setAxisTitle("Fraction of Good Tracks",2);
00281         
00282         histname = "GoodTracksNumberOfRecHitsPerTrackVsLS_" + CategoryName;
00283         GoodTracksNumberOfRecHitsPerTrackVsLS = dqmStore_->bookProfile(histname,histname, LSBin,LSMin,LSMax,0.,40.,"");
00284         GoodTracksNumberOfRecHitsPerTrackVsLS->getTH1()->SetBit(TH1::kCanRebin);
00285         GoodTracksNumberOfRecHitsPerTrackVsLS->setAxisTitle("#Lumi section",1);
00286         GoodTracksNumberOfRecHitsPerTrackVsLS->setAxisTitle("Mean number of RecHits per Good track",2);
00287       }
00288     }
00289 
00290     // book PU monitoring plots :  
00291     //---------------------------  
00292    
00293     if ( doPUmonitoring_ ) {
00294 
00295       for (size_t i=0; i<theVertexMonitor.size(); i++)
00296         theVertexMonitor[i]->beginJob(dqmStore_);
00297 
00298       dqmStore_->setCurrentFolder(MEFolderName+"/PUmonitoring");
00299       
00300       if ( doPlotsVsGoodPVtx_ ) {
00301         // get binning from the configuration
00302         int    GoodPVtxBin   = conf_.getParameter<int>("GoodPVtxBin");
00303         double GoodPVtxMin   = conf_.getParameter<double>("GoodPVtxMin");
00304         double GoodPVtxMax   = conf_.getParameter<double>("GoodPVtxMax");
00305 
00306         histname = "NumberOfTracksVsGoodPVtx";
00307         NumberOfTracksVsGoodPVtx = dqmStore_->bookProfile(histname,histname,GoodPVtxBin,GoodPVtxMin,GoodPVtxMax,TKNoMin, (TKNoMax+0.5)*3.-0.5,"");
00308         NumberOfTracksVsGoodPVtx->getTH1()->SetBit(TH1::kCanRebin);
00309         NumberOfTracksVsGoodPVtx->setAxisTitle("Number of PV",1);
00310         NumberOfTracksVsGoodPVtx->setAxisTitle("Mean number of Tracks per Event",2);
00311       
00312         histname = "NumberOfGoodTracksVsGoodPVtx";
00313         NumberOfGoodTracksVsGoodPVtx = dqmStore_->bookProfile(histname,histname,GoodPVtxBin,GoodPVtxMin,GoodPVtxMax, TKNoMin, TKNoMax,"");
00314         NumberOfGoodTracksVsGoodPVtx->getTH1()->SetBit(TH1::kCanRebin);
00315         NumberOfGoodTracksVsGoodPVtx->setAxisTitle("Number of PV",1);
00316         NumberOfGoodTracksVsGoodPVtx->setAxisTitle("Mean number of Good Tracks per Event",2);
00317 
00318         histname = "FractionOfGoodTracksVsGoodPVtx";
00319         FractionOfGoodTracksVsGoodPVtx = dqmStore_->bookProfile(histname,histname,GoodPVtxBin,GoodPVtxMin,GoodPVtxMax, TKNoMin, TKNoMax,"");
00320         FractionOfGoodTracksVsGoodPVtx->getTH1()->SetBit(TH1::kCanRebin);
00321         FractionOfGoodTracksVsGoodPVtx->setAxisTitle("Number of PV",1);
00322         FractionOfGoodTracksVsGoodPVtx->setAxisTitle("Mean fraction of Good Tracks per Event",2);
00323       }
00324       
00325       if ( doPlotsVsBXlumi_ ) {
00326         // get binning from the configuration
00327         edm::ParameterSet BXlumiParameters = conf_.getParameter<edm::ParameterSet>("BXlumiSetup");
00328         int    BXlumiBin   = BXlumiParameters.getParameter<int>("BXlumiBin");
00329         double BXlumiMin   = BXlumiParameters.getParameter<double>("BXlumiMin");
00330         double BXlumiMax   = BXlumiParameters.getParameter<double>("BXlumiMax");
00331       
00332         histname = "NumberOfTracksVsBXlumi_"+ CategoryName;
00333         NumberOfTracksVsBXlumi = dqmStore_->bookProfile(histname,histname, BXlumiBin,BXlumiMin,BXlumiMax, TKNoMin, TKNoMax,"");
00334         NumberOfTracksVsBXlumi->getTH1()->SetBit(TH1::kCanRebin);
00335         NumberOfTracksVsBXlumi->setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00336         NumberOfTracksVsBXlumi->setAxisTitle("Mean number of Good Tracks",2);
00337 
00338         histname = "NumberOfGoodTracksVsBXlumi_"+ CategoryName;
00339         NumberOfGoodTracksVsBXlumi = dqmStore_->bookProfile(histname,histname, BXlumiBin,BXlumiMin,BXlumiMax, TKNoMin, TKNoMax,"");
00340         NumberOfGoodTracksVsBXlumi->getTH1()->SetBit(TH1::kCanRebin);
00341         NumberOfGoodTracksVsBXlumi->setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00342         NumberOfGoodTracksVsBXlumi->setAxisTitle("Mean number of Good Tracks",2);
00343       
00344         histname = "FractionOfGoodTracksVsBXlumi_"+ CategoryName;
00345         FractionOfGoodTracksVsBXlumi = dqmStore_->bookProfile(histname,histname, BXlumiBin,BXlumiMin,BXlumiMax, TKNoMin, TKNoMax,"");
00346         FractionOfGoodTracksVsBXlumi->getTH1()->SetBit(TH1::kCanRebin);
00347         FractionOfGoodTracksVsBXlumi->setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00348         FractionOfGoodTracksVsBXlumi->setAxisTitle("Mean fraction of Good Tracks",2);
00349       
00350       }
00351     }
00352 
00353     theTrackAnalyzer->beginJob(dqmStore_);
00354 
00355     // book the Seed Property histograms
00356     // ---------------------------------------------------------------------------------//
00357 
00358     dqmStore_->setCurrentFolder(MEFolderName+"/TrackBuilding");
00359 
00360     doAllSeedPlots=conf_.getParameter<bool>("doSeedParameterHistos");
00361     doSeedNumberPlot=conf_.getParameter<bool>("doSeedNumberHisto");
00362     doSeedVsClusterPlot=conf_.getParameter<bool>("doSeedVsClusterHisto");
00363     //    if (doAllPlots) doAllSeedPlots=true;
00364 
00365     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"));
00366 
00367     edm::InputTag seedProducer   = conf_.getParameter<edm::InputTag>("SeedProducer");
00368 
00369     if (doAllSeedPlots || doSeedNumberPlot){
00370       histname = "NumberOfSeeds_"+ seedProducer.label() + "_"+ CategoryName;
00371       NumberOfSeeds = dqmStore_->book1D(histname, histname, TKNoSeedBin, TKNoSeedMin, TKNoSeedMax);
00372       NumberOfSeeds->setAxisTitle("Number of Seeds per Event", 1);
00373       NumberOfSeeds->setAxisTitle("Number of Events", 2);
00374     }
00375 
00376     if (doAllSeedPlots || doSeedVsClusterPlot){
00377 
00378       ClusterLabels=  conf_.getParameter<std::vector<std::string> >("ClusterLabels");
00379 
00380       std::vector<double> histoMin,histoMax;
00381       std::vector<int> histoBin; //these vectors are for max min and nbins in histograms 
00382 
00383       int    NClusPxBin  = conf_.getParameter<int>(   "NClusPxBin");
00384       double NClusPxMin  = conf_.getParameter<double>("NClusPxMin");
00385       double NClusPxMax = conf_.getParameter<double>("NClusPxMax");
00386       
00387       int    NClusStrBin = conf_.getParameter<int>(   "NClusStrBin");
00388       double NClusStrMin = conf_.getParameter<double>("NClusStrMin");
00389       double NClusStrMax = conf_.getParameter<double>("NClusStrMax");
00390 
00391       setMaxMinBin(histoMin,histoMax,histoBin,NClusStrMin,NClusStrMax,NClusStrBin,NClusPxMin,NClusPxMax,NClusPxBin);
00392      
00393       for (uint i=0; i<ClusterLabels.size(); i++){
00394         histname = "SeedsVsClusters_" + seedProducer.label() + "_Vs_" + ClusterLabels[i] + "_" + CategoryName;
00395         SeedsVsClusters.push_back(dynamic_cast<MonitorElement*>(dqmStore_->book2D(histname, histname, histoBin[i], histoMin[i], histoMax[i],
00396                                                                                   TKNoSeedBin, TKNoSeedMin, TKNoSeedMax)));
00397         SeedsVsClusters[i]->setAxisTitle("Number of Clusters", 1);
00398         SeedsVsClusters[i]->setAxisTitle("Number of Seeds", 2);
00399       }
00400     }
00401     
00402 
00403     doTkCandPlots=conf_.getParameter<bool>("doTrackCandHistos");
00404     //    if (doAllPlots) doTkCandPlots=true;
00405 
00406     if (doTkCandPlots){
00407 
00408       edm::InputTag tcProducer     = conf_.getParameter<edm::InputTag>("TCProducer");
00409 
00410       histname = "NumberOfTrackCandidates_"+ tcProducer.label() + "_"+ CategoryName;
00411       NumberOfTrackCandidates = dqmStore_->book1D(histname, histname, TCNoBin, TCNoMin, TCNoMax);
00412       NumberOfTrackCandidates->setAxisTitle("Number of Track Candidates per Event", 1);
00413       NumberOfTrackCandidates->setAxisTitle("Number of Event", 2);
00414     }
00415 
00416    
00417     theTrackBuildingAnalyzer->beginJob(dqmStore_);
00418     
00419 
00420     if (doLumiAnalysis) {
00421 //      if (NumberOfTracks) NumberOfTracks->setLumiFlag();
00422 //      if (NumberOfGoodTracks) NumberOfGoodTracks->setLumiFlag();
00423 //      if (FractionOfGoodTracks) FractionOfGoodTracks->setLumiFlag();
00424       if ( NumberOfTracks_lumiFlag       ) NumberOfTracks_lumiFlag       -> setLumiFlag();
00425       if ( NumberOfGoodTracks_lumiFlag   ) NumberOfGoodTracks_lumiFlag   -> setLumiFlag();
00426       theTrackAnalyzer->setLumiFlag();    
00427     }
00428 
00429     if (doTrackerSpecific_ || doAllPlots) {
00430 
00431       ClusterLabels=  conf_.getParameter<std::vector<std::string> >("ClusterLabels");
00432 
00433       std::vector<double> histoMin,histoMax;
00434       std::vector<int> histoBin; //these vectors are for max min and nbins in histograms 
00435 
00436       /*
00437       int    NClusPxBin  = conf_.getParameter<int>(   "NClusPxBin");
00438       double NClusPxMin  = conf_.getParameter<double>("NClusPxMin");
00439       double NClusPxMax = conf_.getParameter<double>("NClusPxMax");
00440       
00441       int    NClusStrBin = conf_.getParameter<int>(   "NClusStrBin");
00442       double NClusStrMin = conf_.getParameter<double>("NClusStrMin");
00443       double NClusStrMax = conf_.getParameter<double>("NClusStrMax");
00444       */
00445 
00446       /*
00447       int    NClus2DTotBin = conf_.getParameter<int>(   "NClus2DTotBin");
00448       double NClus2DTotMin = conf_.getParameter<double>("NClus2DTotMin");
00449       double NClus2DTotMax = conf_.getParameter<double>("NClus2DTotMax");
00450       */
00451 
00452       int    NClusStrBin = conf_.getParameter<int>(   "NClusStrBin");
00453       double NClusStrMin = conf_.getParameter<double>("NClusStrMin");
00454       double NClusStrMax = conf_.getParameter<double>("NClusStrMax");
00455 
00456       int    NClusPxBin = conf_.getParameter<int>(   "NClusPxBin");
00457       double NClusPxMin = conf_.getParameter<double>("NClusPxMin");
00458       double NClusPxMax = conf_.getParameter<double>("NClusPxMax");
00459 
00460       int    NTrk2DBin     = conf_.getParameter<int>(   "NTrk2DBin");
00461       double NTrk2DMin     = conf_.getParameter<double>("NTrk2DMin");
00462       double NTrk2DMax     = conf_.getParameter<double>("NTrk2DMax");
00463 
00464       //      setMaxMinBin(histoMin,histoMax,histoBin,NClusStrMin,NClusStrMax,NClusStrBin,NClusPxMin,NClusPxMax,NClusPxBin);
00465       setMaxMinBin(histoMin,histoMax,histoBin,
00466                    NClusStrMin,NClusStrMax,NClusStrBin,
00467                    NClusPxMin,  NClusPxMax,  NClusPxBin);
00468      
00469       /*
00470       int    NClusPxBin  = conf_.getParameter<int>(   "NClusPxBin");
00471       double NClusPxMin  = conf_.getParameter<double>("NClusPxMin");
00472       double NClusPxMax = conf_.getParameter<double>("NClusPxMax");
00473 
00474       int    NClusStrBin = conf_.getParameter<int>(   "NClusStrBin");
00475       double NClusStrMin = conf_.getParameter<double>("NClusStrMin");
00476       double NClusStrMax = conf_.getParameter<double>("NClusStrMax");
00477       */
00478 
00479       dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties");
00480 
00481       for (uint i=0; i<ClusterLabels.size(); i++){
00482 
00483         dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties");
00484         histname = "TracksVs" + ClusterLabels[i] + "Cluster_" + CategoryName;
00485         NumberOfTrkVsClusters.push_back(dynamic_cast<MonitorElement*>(dqmStore_->book2D(histname, histname,
00486                                                                                         histoBin[i], histoMin[i], histoMax[i],
00487                                                                                         NTrk2DBin,NTrk2DMin,NTrk2DMax
00488                                                                                         )));
00489         std::string title = "Number of " + ClusterLabels[i] + " Clusters";
00490         NumberOfTrkVsClusters[i]->setAxisTitle(title, 1);
00491         NumberOfTrkVsClusters[i]->setAxisTitle("Number of Seeds", 2);
00492 
00493         if (doGoodTrackPlots_ || doAllPlots ) {
00494 
00495           dqmStore_->setCurrentFolder(MEFolderName+"/HitProperties/GoodTracks");
00496 
00497           if(ClusterLabels[i].compare("Tot")==0){
00498             histname = "GoodTracksVs" + ClusterLabels[i] + "Cluster_" + CategoryName; 
00499             NumberOfGoodTrkVsClus = dqmStore_->book2D(histname,histname,
00500                                                       histoBin[i], histoMin[i], histoMax[i],
00501                                                       TKNoBin,TKNoMin,TKNoMax
00502                                                       );
00503             NumberOfGoodTrkVsClus->setAxisTitle("# of Clusters in (Pixel+Strip) Detectors", 1);
00504             NumberOfGoodTrkVsClus->setAxisTitle("Number of Good Tracks", 2);
00505           }
00506         }
00507       }
00508 
00509       /*
00510       histname = "TracksVsClusters_" + CategoryName; 
00511       NumberOfTrkVsClus = dqmStore_->book2D(histname,histname,
00512                                             NClus2DTotBin,NClus2DTotMin,NClus2DTotMax,
00513                                             NTrk2DBin,NTrk2DMin,NTrk2DMax
00514                                             );
00515       NumberOfTrkVsClus->setAxisTitle("# of Clusters in (Pixel+Strip) Detectors", 1);
00516       NumberOfTrkVsClus->setAxisTitle("Number of Tracks", 2);
00517       */
00518     }
00519 
00520     
00521 }
00522 
00523 // -- BeginRun
00524 //---------------------------------------------------------------------------------//
00525 void TrackingMonitor::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup)
00526 {
00527  
00528   // Initialize the GenericTriggerEventFlag
00529   if ( genTriggerEventFlag_->on() ) genTriggerEventFlag_->initRun( iRun, iSetup );
00530 }
00531 
00532 // - BeginLumi
00533 // ---------------------------------------------------------------------------------//
00534 void TrackingMonitor::beginLuminosityBlock(const edm::LuminosityBlock& lumi, const edm::EventSetup&  eSetup) {
00535 
00536   if (doLumiAnalysis) {
00537 //    dqmStore_->softReset(NumberOfTracks);
00538 //    dqmStore_->softReset(NumberOfGoodTracks);
00539 //    dqmStore_->softReset(FractionOfGoodTracks);
00540 //    theTrackAnalyzer->doSoftReset(dqmStore_);    
00541     if ( NumberOfTracks_lumiFlag       ) NumberOfTracks_lumiFlag       -> Reset();
00542     if ( NumberOfGoodTracks_lumiFlag   ) NumberOfGoodTracks_lumiFlag   -> Reset();
00543     theTrackAnalyzer->doReset(dqmStore_);    
00544   }
00545 }
00546 
00547 // -- Analyse
00548 // ---------------------------------------------------------------------------------//
00549 void TrackingMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) 
00550 {
00551     // Filter out events if Trigger Filtering is requested
00552     if (genTriggerEventFlag_->on()&& ! genTriggerEventFlag_->accept( iEvent, iSetup) ) return;
00553 
00554     // input tags for collections from the configuration
00555     edm::InputTag trackProducer  = conf_.getParameter<edm::InputTag>("TrackProducer");
00556     edm::InputTag seedProducer   = conf_.getParameter<edm::InputTag>("SeedProducer");
00557     edm::InputTag tcProducer     = conf_.getParameter<edm::InputTag>("TCProducer");
00558     edm::InputTag bsSrc          = conf_.getParameter<edm::InputTag>("beamSpot");
00559     edm::InputTag pvSrc_         = conf_.getParameter<edm::InputTag>("primaryVertex");
00560     std::string Quality = conf_.getParameter<std::string>("Quality");
00561     std::string Algo    = conf_.getParameter<std::string>("AlgoName");
00562 
00563     //  Analyse the tracks
00564     //  if the collection is empty, do not fill anything
00565     // ---------------------------------------------------------------------------------//
00566 
00567     // get the track collection
00568     edm::Handle<reco::TrackCollection> trackHandle;
00569     iEvent.getByLabel(trackProducer, trackHandle);
00570 
00571     if (trackHandle.isValid()) 
00572     {
00573 
00574        reco::TrackCollection trackCollection = *trackHandle;
00575         // calculate the mean # rechits and layers
00576         int totalNumTracks = 0, totalRecHits = 0, totalLayers = 0;
00577         int totalNumHPTracks = 0, totalNumPt1Tracks = 0, totalNumHPPt1Tracks = 0;
00578         double frac = 0.;
00579         for (reco::TrackCollection::const_iterator track = trackCollection.begin(); track!=trackCollection.end(); ++track) 
00580         {
00581 
00582             if( track->quality(reco::TrackBase::highPurity) ) {
00583               ++totalNumHPTracks;
00584               if ( track->pt() >= 1. ) {
00585                 ++totalNumHPPt1Tracks;
00586                 if ( doProfilesVsLS_ || doAllPlots)
00587                   if ( doGoodTrackPlots_ || doAllPlots ) {
00588                     GoodTracksNumberOfRecHitsPerTrackVsLS->Fill(static_cast<double>(iEvent.id().luminosityBlock()),track->recHitsSize());
00589                   }
00590               }
00591             }
00592             
00593             if ( track->pt() >= 1. ) ++totalNumPt1Tracks;
00594         
00595 
00596             if( Quality == "highPurity") 
00597             {
00598                 if( !track->quality(reco::TrackBase::highPurity) ) continue;
00599             }
00600             else if( Quality == "tight") 
00601             {
00602                 if( !track->quality(reco::TrackBase::tight) ) continue;
00603             }
00604             else if( Quality == "loose") 
00605             {
00606                 if( !track->quality(reco::TrackBase::loose) ) continue;
00607             }
00608             
00609             totalNumTracks++;
00610             totalRecHits    += track->numberOfValidHits();
00611             totalLayers     += track->hitPattern().trackerLayersWithMeasurement();
00612 
00613             // do analysis per track
00614             theTrackAnalyzer->analyze(iEvent, iSetup, *track);
00615         }
00616 
00617         if (totalNumPt1Tracks > 0) frac = static_cast<double>(totalNumHPPt1Tracks)/static_cast<double>(totalNumPt1Tracks);
00618 
00619         if (doGeneralPropertiesPlots_ || doAllPlots){
00620           NumberOfTracks       -> Fill(totalNumTracks);
00621           NumberOfGoodTracks   -> Fill(totalNumHPPt1Tracks);
00622           FractionOfGoodTracks -> Fill(frac);
00623         }
00624 
00625         if ( doProfilesVsLS_ || doAllPlots) {
00626           NumberOfTracksVsLS->Fill(static_cast<double>(iEvent.id().luminosityBlock()),totalNumTracks);
00627           if ( doGoodTrackPlots_ || doAllPlots ) {
00628             NumberOfGoodTracksVsLS->Fill(static_cast<double>(iEvent.id().luminosityBlock()),totalNumHPPt1Tracks);
00629             GoodTracksFractionVsLS->Fill(static_cast<double>(iEvent.id().luminosityBlock()),frac);
00630           }
00631         }
00632 
00633         if ( doLumiAnalysis ) {
00634           NumberOfTracks_lumiFlag       -> Fill(totalNumTracks);
00635           NumberOfGoodTracks_lumiFlag   -> Fill(totalNumHPPt1Tracks);
00636         }
00637 
00638         if (doGeneralPropertiesPlots_ || doAllPlots){
00639           if( totalNumTracks > 0 )
00640               {
00641                 double meanRecHits = static_cast<double>(totalRecHits) / static_cast<double>(totalNumTracks);
00642                 double meanLayers  = static_cast<double>(totalLayers)  / static_cast<double>(totalNumTracks);
00643                 NumberOfMeanRecHitsPerTrack -> Fill(meanRecHits);
00644                 NumberOfMeanLayersPerTrack  -> Fill(meanLayers);
00645               }
00646           }
00647     
00648         //  Analyse the Track Building variables 
00649         //  if the collection is empty, do not fill anything
00650         // ---------------------------------------------------------------------------------//
00651         
00652                             
00653             // fill the TrackCandidate info
00654             if (doTkCandPlots) 
00655               {
00656         
00657                 // magnetic field
00658                 edm::ESHandle<MagneticField> theMF;
00659                 iSetup.get<IdealMagneticFieldRecord>().get(theMF);  
00660                 
00661                 // get the beam spot
00662                 edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00663                 iEvent.getByLabel(bsSrc,recoBeamSpotHandle);
00664                 const reco::BeamSpot& bs = *recoBeamSpotHandle;      
00665                 
00666                 // get the candidate collection
00667                 edm::Handle<TrackCandidateCollection> theTCHandle;
00668                 iEvent.getByLabel(tcProducer, theTCHandle ); 
00669                 const TrackCandidateCollection& theTCCollection = *theTCHandle;
00670 
00671                 if (theTCHandle.isValid())
00672                   {
00673                     NumberOfTrackCandidates->Fill(theTCCollection.size());
00674                     iSetup.get<TransientRecHitRecord>().get(builderName,theTTRHBuilder);
00675                     for( TrackCandidateCollection::const_iterator cand = theTCCollection.begin(); cand != theTCCollection.end(); ++cand)
00676                       {
00677                         theTrackBuildingAnalyzer->analyze(iEvent, iSetup, *cand, bs, theMF, theTTRHBuilder);
00678                       }
00679                   }
00680                 else
00681                   {
00682                     edm::LogWarning("TrackingMonitor") << "No Track Candidates in the event.  Not filling associated histograms";
00683                   }
00684               }
00685     
00686             //plots for trajectory seeds
00687 
00688             if (doAllSeedPlots || doSeedNumberPlot || doSeedVsClusterPlot || runTrackBuildingAnalyzerForSeed){
00689             
00690             // get the seed collection
00691             edm::Handle<edm::View<TrajectorySeed> > seedHandle;
00692             iEvent.getByLabel(seedProducer, seedHandle);
00693             const edm::View<TrajectorySeed>& seedCollection = *seedHandle;
00694             
00695             // fill the seed info
00696             if (seedHandle.isValid()) 
00697               {
00698                 if(doAllSeedPlots || doSeedNumberPlot) NumberOfSeeds->Fill(seedCollection.size());
00699 
00700                 if(doAllSeedPlots || doSeedVsClusterPlot){
00701 
00702                   std::vector<int> NClus;
00703                   setNclus(iEvent,NClus);
00704                   for (uint  i=0; i< ClusterLabels.size(); i++){
00705                     SeedsVsClusters[i]->Fill(NClus[i],seedCollection.size());
00706                   }
00707                 }
00708 
00709                 if (doAllSeedPlots || runTrackBuildingAnalyzerForSeed){
00710 
00711                   //here duplication of mag field and be informations is needed to allow seed and track cand histos to be independent
00712                   // magnetic field
00713                   edm::ESHandle<MagneticField> theMF;
00714                   iSetup.get<IdealMagneticFieldRecord>().get(theMF);  
00715                   
00716                   // get the beam spot
00717                   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00718                   iEvent.getByLabel(bsSrc,recoBeamSpotHandle);
00719                   const reco::BeamSpot& bs = *recoBeamSpotHandle;      
00720                   
00721                   iSetup.get<TransientRecHitRecord>().get(builderName,theTTRHBuilder);
00722                   for(size_t i=0; i < seedHandle->size(); ++i)
00723                     {
00724                       edm::RefToBase<TrajectorySeed> seed(seedHandle, i);
00725                       theTrackBuildingAnalyzer->analyze(iEvent, iSetup, *seed, bs, theMF, theTTRHBuilder);
00726                     }
00727                 }
00728                 
00729               }
00730             else
00731               {
00732                 edm::LogWarning("TrackingMonitor") << "No Trajectory seeds in the event.  Not filling associated histograms";
00733               }
00734             }
00735 
00736 
00737          if (doTrackerSpecific_ || doAllPlots) 
00738           {
00739             std::vector<int> NClus;
00740             setNclus(iEvent,NClus);
00741             for (uint  i=0; i< ClusterLabels.size(); i++){
00742               NumberOfTrkVsClusters[i]->Fill(NClus[i],totalNumTracks);
00743             }
00744             if ( doGoodTrackPlots_ || doAllPlots ) {
00745               for (uint  i=0; i< ClusterLabels.size(); i++){
00746                 if(ClusterLabels[i].compare("Tot")==0){
00747                   NumberOfGoodTrkVsClus->Fill( NClus[i],totalNumHPPt1Tracks);
00748                 }
00749               }
00750             }
00751 
00752             /*
00753             edm::Handle< edmNew::DetSetVector<SiStripCluster> > strip_clusters;
00754             iEvent.getByLabel("siStripClusters", strip_clusters);
00755             edm::Handle< edmNew::DetSetVector<SiPixelCluster> > pixel_clusters;
00756             iEvent.getByLabel("siPixelClusters", pixel_clusters);
00757             if (strip_clusters.isValid() && pixel_clusters.isValid()) 
00758               {
00759                 unsigned int ncluster_pix   = (*pixel_clusters).dataSize(); 
00760                 unsigned int ncluster_strip = (*strip_clusters).dataSize(); 
00761                 double ratio = 0.0;
00762                 if ( ncluster_pix > 0) ratio = atan(ncluster_pix*1.0/ncluster_strip);
00763 
00764                 NumberOfStripClus->Fill(ncluster_strip);
00765                 NumberOfPixelClus->Fill(ncluster_pix);
00766                 RatioOfPixelAndStripClus->Fill(ratio);
00767 
00768                 NumberOfTrkVsClus->Fill( ncluster_strip+ncluster_pix,totalNumTracks);
00769             
00770                 if (doGoodTrackPlots_) {
00771                   NumberOfGoodTrkVsClus->Fill( ncluster_strip+ncluster_pix,totalNumHPPt1Tracks);
00772                 }
00773               }
00774             */
00775           }
00776 
00777          if ( doPUmonitoring_ ) {
00778            
00779            // do vertex monitoring
00780            for (size_t i=0; i<theVertexMonitor.size(); i++)
00781              theVertexMonitor[i]->analyze(iEvent, iSetup);
00782            
00783            if ( doPlotsVsGoodPVtx_ ) {
00784              
00785              edm::InputTag primaryVertexInputTag = conf_.getParameter<edm::InputTag>("primaryVertex");
00786              
00787              size_t totalNumGoodPV = 0;
00788              edm::Handle< reco::VertexCollection > pvHandle;
00789              iEvent.getByLabel(primaryVertexInputTag, pvHandle );
00790              if (pvHandle.isValid())
00791                {
00792                  
00793                  for (reco::VertexCollection::const_iterator pv = pvHandle->begin();
00794                       pv != pvHandle->end(); ++pv) {
00795                    
00796                    //--- pv fake (the pv collection should have size==1 and the pv==beam spot) 
00797                    if (pv->isFake() || pv->tracksSize()==0) continue;
00798                    
00799                    // definition of goodOfflinePrimaryVertex
00800                    if (pv->ndof() < 4. || pv->z() > 24.)  continue;
00801                    totalNumGoodPV++;
00802                  }
00803                  
00804                  NumberOfTracksVsGoodPVtx       -> Fill( totalNumGoodPV, totalNumTracks      );
00805                  NumberOfGoodTracksVsGoodPVtx   -> Fill( totalNumGoodPV, totalNumHPPt1Tracks );
00806                  FractionOfGoodTracksVsGoodPVtx -> Fill( totalNumGoodPV, frac                );
00807                }
00808            }
00809            
00810            if ( doPlotsVsBXlumi_ ) {
00811              
00812              double bxlumi = theLumiDetails_->getValue(iEvent);
00813              
00814              NumberOfTracksVsBXlumi       -> Fill( bxlumi, totalNumTracks      );
00815              NumberOfGoodTracksVsBXlumi   -> Fill( bxlumi, totalNumHPPt1Tracks );
00816              FractionOfGoodTracksVsBXlumi -> Fill( bxlumi, frac                );
00817            }
00818            
00819          }
00820          
00821     }
00822 
00823 }
00824 
00825 
00826 void TrackingMonitor::endRun(const edm::Run&, const edm::EventSetup&) 
00827 {
00828   if (doLumiAnalysis) {
00829     /*
00830     dqmStore_->disableSoftReset(NumberOfTracks);
00831     dqmStore_->disableSoftReset(NumberOfGoodTracks);
00832     dqmStore_->disableSoftReset(FractionOfGoodTracks);
00833     theTrackAnalyzer->undoSoftReset(dqmStore_);    
00834     */
00835   }
00836   
00837 }
00838 
00839 void TrackingMonitor::endJob(void) 
00840 {
00841     bool outputMEsInRootFile   = conf_.getParameter<bool>("OutputMEsInRootFile");
00842     std::string outputFileName = conf_.getParameter<std::string>("OutputFileName");
00843     if(outputMEsInRootFile)
00844     {
00845         dqmStore_->showDirStructure();
00846         dqmStore_->save(outputFileName);
00847     }
00848 }
00849 
00850 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) 
00851 {
00852   arrayMin.resize(ClusterLabels.size());
00853   arrayMax.resize(ClusterLabels.size());
00854   arrayBin.resize(ClusterLabels.size());
00855 
00856   for (uint i=0; i<ClusterLabels.size(); ++i){
00857 
00858     if (ClusterLabels[i].compare("Pix")==0) {arrayMin[i]=pmin; arrayMax[i]=pmax; arrayBin[i]=pbin;}
00859     else if(ClusterLabels[i].compare("Strip")==0){arrayMin[i]=smin; arrayMax[i]=smax; arrayBin[i]=sbin;}
00860     else if(ClusterLabels[i].compare("Tot")==0){arrayMin[i]=smin; arrayMax[i]=smax+pmax; arrayBin[i]=sbin;}
00861     else {edm::LogWarning("TrackingMonitor")  << "Cluster Label " << ClusterLabels[i] << " not defined, using strip parameters "; 
00862       arrayMin[i]=smin; arrayMax[i]=smax; arrayBin[i]=sbin;}
00863 
00864   }
00865     
00866 }
00867 
00868 void TrackingMonitor::setNclus(const edm::Event& iEvent,std::vector<int> &arrayNclus) 
00869 {
00870 
00871   int ncluster_pix=-1;
00872   int ncluster_strip=-1;
00873 
00874   edm::Handle< edmNew::DetSetVector<SiStripCluster> > strip_clusters;
00875   iEvent.getByLabel("siStripClusters", strip_clusters);
00876   edm::Handle< edmNew::DetSetVector<SiPixelCluster> > pixel_clusters;
00877   iEvent.getByLabel("siPixelClusters", pixel_clusters);
00878 
00879   if (strip_clusters.isValid() && pixel_clusters.isValid()) 
00880     {
00881                 ncluster_pix   = (*pixel_clusters).dataSize(); 
00882                 ncluster_strip = (*strip_clusters).dataSize(); 
00883     }
00884 
00885   arrayNclus.resize(ClusterLabels.size());
00886   for (uint i=0; i<ClusterLabels.size(); ++i){
00887 
00888     if (ClusterLabels[i].compare("Pix")==0) arrayNclus[i]=ncluster_pix ;
00889     else if(ClusterLabels[i].compare("Strip")==0) arrayNclus[i]=ncluster_strip;
00890     else if(ClusterLabels[i].compare("Tot")==0)arrayNclus[i]=ncluster_pix+ncluster_strip;
00891     else {edm::LogWarning("TrackingMonitor") << "Cluster Label " << ClusterLabels[i] << " not defined using stri parametrs ";
00892       arrayNclus[i]=ncluster_strip ;}
00893   }
00894     
00895 }
00896 DEFINE_FWK_MODULE(TrackingMonitor);