CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/Validation/RecoTrack/src/MTVHistoProducerAlgoForTracker.cc

Go to the documentation of this file.
00001 #include "Validation/RecoTrack/interface/MTVHistoProducerAlgoForTracker.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 #include "DataFormats/TrackReco/interface/Track.h"
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00005 
00006 #include "DataFormats/TrackReco/interface/DeDxData.h"
00007 #include "DataFormats/Common/interface/ValueMap.h"
00008 #include "DataFormats/Common/interface/Ref.h"
00009 #include <DataFormats/TrackReco/interface/TrackFwd.h>
00010 
00011 #include "DQMServices/ClientConfig/interface/FitSlicesYTool.h"
00012 
00013 #include "TMath.h"
00014 #include <TF1.h>
00015 
00016 using namespace std;
00017 
00018 MTVHistoProducerAlgoForTracker::MTVHistoProducerAlgoForTracker(const edm::ParameterSet& pset): MTVHistoProducerAlgo(pset){
00019   //parameters for _vs_eta plots
00020   minEta  = pset.getParameter<double>("minEta");  
00021   maxEta  = pset.getParameter<double>("maxEta");
00022   nintEta = pset.getParameter<int>("nintEta");
00023   useFabsEta = pset.getParameter<bool>("useFabsEta");
00024 
00025   //parameters for _vs_pt plots
00026   minPt  = pset.getParameter<double>("minPt");
00027   maxPt  = pset.getParameter<double>("maxPt");
00028   nintPt = pset.getParameter<int>("nintPt");
00029   useInvPt = pset.getParameter<bool>("useInvPt");
00030   useLogPt = pset.getUntrackedParameter<bool>("useLogPt",false);
00031 
00032   //parameters for _vs_Hit plots
00033   minHit  = pset.getParameter<double>("minHit");
00034   maxHit  = pset.getParameter<double>("maxHit");
00035   nintHit = pset.getParameter<int>("nintHit");
00036   
00037   //parameters for _vs_phi plots
00038   minPhi  = pset.getParameter<double>("minPhi");
00039   maxPhi  = pset.getParameter<double>("maxPhi");
00040   nintPhi = pset.getParameter<int>("nintPhi");
00041 
00042   //parameters for _vs_Dxy plots
00043   minDxy  = pset.getParameter<double>("minDxy");
00044   maxDxy  = pset.getParameter<double>("maxDxy");
00045   nintDxy = pset.getParameter<int>("nintDxy");
00046     
00047   //parameters for _vs_Dz plots
00048   minDz  = pset.getParameter<double>("minDz");
00049   maxDz  = pset.getParameter<double>("maxDz");
00050   nintDz = pset.getParameter<int>("nintDz");
00051 
00052   //parameters for _vs_ProductionVertexTransvPosition plots
00053   minVertpos  = pset.getParameter<double>("minVertpos");
00054   maxVertpos  = pset.getParameter<double>("maxVertpos");
00055   nintVertpos = pset.getParameter<int>("nintVertpos");
00056 
00057   //parameters for _vs_ProductionVertexZPosition plots
00058   minZpos  = pset.getParameter<double>("minZpos");
00059   maxZpos  = pset.getParameter<double>("maxZpos");
00060   nintZpos = pset.getParameter<int>("nintZpos");
00061 
00062   //parameters for dE/dx plots
00063   minDeDx  = pset.getParameter<double>("minDeDx");
00064   maxDeDx  = pset.getParameter<double>("maxDeDx");
00065   nintDeDx = pset.getParameter<int>("nintDeDx");
00066   
00067   //parameters for Pileup plots
00068   minVertcount  = pset.getParameter<double>("minVertcount");
00069   maxVertcount  = pset.getParameter<double>("maxVertcount");
00070   nintVertcount = pset.getParameter<int>("nintVertcount");
00071 
00072   //parameters for resolution plots
00073   ptRes_rangeMin = pset.getParameter<double>("ptRes_rangeMin");
00074   ptRes_rangeMax = pset.getParameter<double>("ptRes_rangeMax");
00075   ptRes_nbin = pset.getParameter<int>("ptRes_nbin");
00076 
00077   phiRes_rangeMin = pset.getParameter<double>("phiRes_rangeMin");
00078   phiRes_rangeMax = pset.getParameter<double>("phiRes_rangeMax");
00079   phiRes_nbin = pset.getParameter<int>("phiRes_nbin");
00080 
00081   cotThetaRes_rangeMin = pset.getParameter<double>("cotThetaRes_rangeMin");
00082   cotThetaRes_rangeMax = pset.getParameter<double>("cotThetaRes_rangeMax");
00083   cotThetaRes_nbin = pset.getParameter<int>("cotThetaRes_nbin");
00084 
00085   dxyRes_rangeMin = pset.getParameter<double>("dxyRes_rangeMin");
00086   dxyRes_rangeMax = pset.getParameter<double>("dxyRes_rangeMax");
00087   dxyRes_nbin = pset.getParameter<int>("dxyRes_nbin");
00088 
00089   dzRes_rangeMin = pset.getParameter<double>("dzRes_rangeMin");
00090   dzRes_rangeMax = pset.getParameter<double>("dzRes_rangeMax");
00091   dzRes_nbin = pset.getParameter<int>("dzRes_nbin");
00092 
00093 
00094   //--- tracking particle selectors for efficiency measurements
00095   using namespace edm;
00096 
00097   ParameterSet generalTpSelectorPSet = pset.getParameter<ParameterSet>("generalTpSelector");
00098   ParameterSet TpSelectorForEfficiencyVsEtaPSet = pset.getParameter<ParameterSet>("TpSelectorForEfficiencyVsEta");
00099   ParameterSet TpSelectorForEfficiencyVsPhiPSet = pset.getParameter<ParameterSet>("TpSelectorForEfficiencyVsPhi");
00100   ParameterSet TpSelectorForEfficiencyVsPtPSet = pset.getParameter<ParameterSet>("TpSelectorForEfficiencyVsPt");
00101   ParameterSet TpSelectorForEfficiencyVsVTXRPSet = pset.getParameter<ParameterSet>("TpSelectorForEfficiencyVsVTXR");
00102   ParameterSet TpSelectorForEfficiencyVsVTXZPSet = pset.getParameter<ParameterSet>("TpSelectorForEfficiencyVsVTXZ");
00103 
00104   ParameterSet TpSelectorForEfficiencyVsConPSet = TpSelectorForEfficiencyVsEtaPSet;
00105   Entry name("signalOnly",false,true);
00106   TpSelectorForEfficiencyVsConPSet.insert(true,"signalOnly",name);
00107   
00108   using namespace reco::modules;
00109   generalTpSelector               = new TrackingParticleSelector(ParameterAdapter<TrackingParticleSelector>::make(generalTpSelectorPSet));
00110   TpSelectorForEfficiencyVsEta    = new TrackingParticleSelector(ParameterAdapter<TrackingParticleSelector>::make(TpSelectorForEfficiencyVsEtaPSet));
00111   TpSelectorForEfficiencyVsCon    = new TrackingParticleSelector(ParameterAdapter<TrackingParticleSelector>::make(TpSelectorForEfficiencyVsConPSet));
00112   TpSelectorForEfficiencyVsPhi    = new TrackingParticleSelector(ParameterAdapter<TrackingParticleSelector>::make(TpSelectorForEfficiencyVsPhiPSet));
00113   TpSelectorForEfficiencyVsPt     = new TrackingParticleSelector(ParameterAdapter<TrackingParticleSelector>::make(TpSelectorForEfficiencyVsPtPSet));
00114   TpSelectorForEfficiencyVsVTXR   = new TrackingParticleSelector(ParameterAdapter<TrackingParticleSelector>::make(TpSelectorForEfficiencyVsVTXRPSet));
00115   TpSelectorForEfficiencyVsVTXZ   = new TrackingParticleSelector(ParameterAdapter<TrackingParticleSelector>::make(TpSelectorForEfficiencyVsVTXZPSet));
00116 
00117   // fix for the LogScale by Ryan
00118   if(useLogPt){
00119     maxPt=log10(maxPt);
00120     if(minPt > 0){
00121       minPt=log10(minPt);
00122     }
00123     else{
00124       edm::LogWarning("MultiTrackValidator") 
00125         << "minPt = "
00126         << minPt << " <= 0 out of range while requesting log scale.  Using minPt = 0.1.";
00127       minPt=log10(0.1);
00128     }
00129   }
00130 
00131 }
00132 
00133 MTVHistoProducerAlgoForTracker::~MTVHistoProducerAlgoForTracker(){
00134   delete generalTpSelector;
00135   delete TpSelectorForEfficiencyVsEta;
00136   delete TpSelectorForEfficiencyVsPhi;
00137   delete TpSelectorForEfficiencyVsPt;
00138   delete TpSelectorForEfficiencyVsVTXR;
00139   delete TpSelectorForEfficiencyVsVTXZ;
00140 }
00141 
00142 
00143 void MTVHistoProducerAlgoForTracker::setUpVectors(){
00144   std::vector<double> etaintervalsv;
00145   std::vector<double> phiintervalsv;
00146   std::vector<double> pTintervalsv;
00147   std::vector<double> dxyintervalsv;
00148   std::vector<double> dzintervalsv;
00149   std::vector<double> vertposintervalsv;
00150   std::vector<double> zposintervalsv;
00151   std::vector<double> vertcountintervalsv;
00152 
00153   std::vector<int>    totSIMveta,totASSveta,totASS2veta,totloopveta,totmisidveta,totASS2vetaSig,totRECveta;
00154   std::vector<int>    totSIMvpT,totASSvpT,totASS2vpT,totRECvpT,totloopvpT,totmisidvpT;
00155   std::vector<int>    totSIMv_hit,totASSv_hit,totASS2v_hit,totRECv_hit,totloopv_hit,totmisidv_hit;
00156   std::vector<int>    totSIMv_phi,totASSv_phi,totASS2v_phi,totRECv_phi,totloopv_phi,totmisidv_phi;
00157   std::vector<int>    totSIMv_dxy,totASSv_dxy,totASS2v_dxy,totRECv_dxy,totloopv_dxy,totmisidv_dxy;
00158   std::vector<int>    totSIMv_dz,totASSv_dz,totASS2v_dz,totRECv_dz,totloopv_dz,totmisidv_dz;
00159 
00160   std::vector<int>    totSIMv_vertpos,totASSv_vertpos,totSIMv_zpos,totASSv_zpos; 
00161   std::vector<int>    totSIMv_vertcount,totASSv_vertcount,totRECv_vertcount,totASS2v_vertcount; 
00162   std::vector<int>    totRECv_algo;
00163 
00164   double step=(maxEta-minEta)/nintEta;
00165   //std::ostringstream title,name; ///BM, what is this?
00166   etaintervalsv.push_back(minEta);
00167   for (int k=1;k<nintEta+1;k++) {
00168     double d=minEta+k*step;
00169     etaintervalsv.push_back(d);
00170     totSIMveta.push_back(0);
00171     totASSveta.push_back(0);
00172     totASS2veta.push_back(0);
00173     totloopveta.push_back(0);
00174     totmisidveta.push_back(0);
00175     totASS2vetaSig.push_back(0);
00176     totRECveta.push_back(0);
00177   }   
00178   etaintervals.push_back(etaintervalsv);
00179   totSIMeta.push_back(totSIMveta);
00180   totCONeta.push_back(totASSveta);
00181   totASSeta.push_back(totASSveta);
00182   totASS2eta.push_back(totASS2veta);
00183   totloopeta.push_back(totloopveta);
00184   totmisideta.push_back(totmisidveta);
00185   totASS2etaSig.push_back(totASS2vetaSig);
00186   totRECeta.push_back(totRECveta);
00187   totFOMT_eta.push_back(totASSveta);
00188     
00189   totASS2_itpu_eta_entire.push_back(totASS2veta);
00190   totASS2_itpu_eta_entire_signal.push_back(totASS2vetaSig);
00191   totASS2_ootpu_eta_entire.push_back(totASS2veta);
00192 
00193   for (size_t i = 0; i < 15; i++) {
00194     totRECv_algo.push_back(0);
00195   }
00196   totREC_algo.push_back(totRECv_algo);
00197 
00198   double stepPt = (maxPt-minPt)/nintPt;
00199   pTintervalsv.push_back(minPt);
00200   for (int k=1;k<nintPt+1;k++) {
00201     double d=0;
00202     if(useLogPt)d=pow(10,minPt+k*stepPt);
00203     else d=minPt+k*stepPt;
00204     pTintervalsv.push_back(d);
00205     totSIMvpT.push_back(0);
00206     totASSvpT.push_back(0);
00207     totASS2vpT.push_back(0);
00208     totRECvpT.push_back(0);
00209     totloopvpT.push_back(0);
00210     totmisidvpT.push_back(0);
00211   }
00212   pTintervals.push_back(pTintervalsv);
00213   totSIMpT.push_back(totSIMvpT);
00214   totASSpT.push_back(totASSvpT);
00215   totASS2pT.push_back(totASS2vpT);
00216   totRECpT.push_back(totRECvpT);
00217   totlooppT.push_back(totloopvpT);
00218   totmisidpT.push_back(totmisidvpT);
00219   
00220   for (int k=1;k<nintHit+1;k++) {
00221     totSIMv_hit.push_back(0);
00222     totASSv_hit.push_back(0);
00223     totASS2v_hit.push_back(0);
00224     totRECv_hit.push_back(0);
00225     totloopv_hit.push_back(0);
00226     totmisidv_hit.push_back(0);
00227   }
00228   totSIM_hit.push_back(totSIMv_hit);
00229   totASS_hit.push_back(totASSv_hit);
00230   totASS2_hit.push_back(totASS2v_hit);
00231   totREC_hit.push_back(totRECv_hit);
00232   totloop_hit.push_back(totloopv_hit);
00233   totmisid_hit.push_back(totmisidv_hit);
00234   
00235   double stepPhi = (maxPhi-minPhi)/nintPhi;
00236   phiintervalsv.push_back(minPhi);
00237   for (int k=1;k<nintPhi+1;k++) {
00238     double d=minPhi+k*stepPhi;
00239     phiintervalsv.push_back(d);
00240     totSIMv_phi.push_back(0);
00241     totASSv_phi.push_back(0);
00242     totASS2v_phi.push_back(0);
00243     totRECv_phi.push_back(0);
00244     totloopv_phi.push_back(0);
00245     totmisidv_phi.push_back(0);
00246   }
00247   phiintervals.push_back(phiintervalsv);
00248   totSIM_phi.push_back(totSIMv_phi);
00249   totASS_phi.push_back(totASSv_phi);
00250   totASS2_phi.push_back(totASS2v_phi);
00251   totREC_phi.push_back(totRECv_phi);
00252   totloop_phi.push_back(totloopv_phi);
00253   totmisid_phi.push_back(totmisidv_phi);
00254   
00255   double stepDxy = (maxDxy-minDxy)/nintDxy;
00256   dxyintervalsv.push_back(minDxy);
00257   for (int k=1;k<nintDxy+1;k++) {
00258     double d=minDxy+k*stepDxy;
00259     dxyintervalsv.push_back(d);
00260     totSIMv_dxy.push_back(0);
00261     totASSv_dxy.push_back(0);
00262     totASS2v_dxy.push_back(0);
00263     totRECv_dxy.push_back(0);
00264     totloopv_dxy.push_back(0);
00265     totmisidv_dxy.push_back(0);
00266   }
00267   dxyintervals.push_back(dxyintervalsv);
00268   totSIM_dxy.push_back(totSIMv_dxy);
00269   totASS_dxy.push_back(totASSv_dxy);
00270   totASS2_dxy.push_back(totASS2v_dxy);
00271   totREC_dxy.push_back(totRECv_dxy);
00272   totloop_dxy.push_back(totloopv_dxy);
00273   totmisid_dxy.push_back(totmisidv_dxy);
00274   
00275   
00276   double stepDz = (maxDz-minDz)/nintDz;
00277   dzintervalsv.push_back(minDz);
00278   for (int k=1;k<nintDz+1;k++) {
00279     double d=minDz+k*stepDz;
00280     dzintervalsv.push_back(d);
00281     totSIMv_dz.push_back(0);
00282     totASSv_dz.push_back(0);
00283     totASS2v_dz.push_back(0);
00284     totRECv_dz.push_back(0);
00285     totloopv_dz.push_back(0);
00286     totmisidv_dz.push_back(0);
00287   }
00288   dzintervals.push_back(dzintervalsv);
00289   totSIM_dz.push_back(totSIMv_dz);
00290   totASS_dz.push_back(totASSv_dz);
00291   totASS2_dz.push_back(totASS2v_dz);
00292   totREC_dz.push_back(totRECv_dz);
00293   totloop_dz.push_back(totloopv_dz);
00294   totmisid_dz.push_back(totmisidv_dz);
00295   
00296   double stepVertpos = (maxVertpos-minVertpos)/nintVertpos;
00297   vertposintervalsv.push_back(minVertpos);
00298   for (int k=1;k<nintVertpos+1;k++) {
00299     double d=minVertpos+k*stepVertpos;
00300     vertposintervalsv.push_back(d);
00301     totSIMv_vertpos.push_back(0);
00302     totASSv_vertpos.push_back(0);
00303   }
00304   vertposintervals.push_back(vertposintervalsv);
00305   totSIM_vertpos.push_back(totSIMv_vertpos);
00306   totASS_vertpos.push_back(totASSv_vertpos);
00307     
00308   double stepZpos = (maxZpos-minZpos)/nintZpos;
00309   zposintervalsv.push_back(minZpos);
00310   for (int k=1;k<nintZpos+1;k++) {
00311     double d=minZpos+k*stepZpos;
00312     zposintervalsv.push_back(d);
00313     totSIMv_zpos.push_back(0);
00314     totASSv_zpos.push_back(0);
00315   }
00316   zposintervals.push_back(zposintervalsv);
00317   totSIM_zpos.push_back(totSIMv_zpos);
00318   totCONzpos.push_back(totSIMv_zpos);
00319   totASS_zpos.push_back(totASSv_zpos);
00320   totSIM_vertz_entire.push_back(totSIMv_zpos);
00321   totASS_vertz_entire.push_back(totASSv_zpos);
00322   totSIM_vertz_barrel.push_back(totSIMv_zpos);
00323   totASS_vertz_barrel.push_back(totASSv_zpos);
00324   totSIM_vertz_fwdpos.push_back(totSIMv_zpos);
00325   totASS_vertz_fwdpos.push_back(totASSv_zpos);
00326   totSIM_vertz_fwdneg.push_back(totSIMv_zpos);
00327   totASS_vertz_fwdneg.push_back(totASSv_zpos);
00328 
00329   double stepVertcount=(maxVertcount-minVertcount)/nintVertcount;
00330   vertcountintervalsv.push_back(minVertcount);
00331   for (int k=1;k<nintVertcount+1;k++) {
00332     double d=minVertcount+k*stepVertcount;
00333     vertcountintervalsv.push_back(d);
00334     totSIMv_vertcount.push_back(0);
00335     totASSv_vertcount.push_back(0);
00336     totASS2v_vertcount.push_back(0);
00337     totRECv_vertcount.push_back(0);
00338   }   
00339   vertcountintervals.push_back(vertcountintervalsv);
00340   totSIM_vertcount_entire.push_back(totSIMv_vertcount);
00341   totCONvertcount.push_back(totSIMv_vertcount);
00342   totASS_vertcount_entire.push_back(totASSv_vertcount);
00343   totASS2_vertcount_entire.push_back(totASS2v_vertcount);
00344   totASS2_vertcount_entire_signal.push_back(totASS2v_vertcount);
00345   totREC_vertcount_entire.push_back(totRECv_vertcount);
00346   totSIM_vertcount_barrel.push_back(totSIMv_vertcount);
00347   totASS_vertcount_barrel.push_back(totASSv_vertcount);
00348   totASS2_vertcount_barrel.push_back(totASS2v_vertcount);
00349   totREC_vertcount_barrel.push_back(totRECv_vertcount);
00350   totSIM_vertcount_fwdpos.push_back(totSIMv_vertcount);
00351   totASS_vertcount_fwdpos.push_back(totASSv_vertcount);
00352   totASS2_vertcount_fwdpos.push_back(totASS2v_vertcount);
00353   totREC_vertcount_fwdpos.push_back(totRECv_vertcount);
00354   totSIM_vertcount_fwdneg.push_back(totSIMv_vertcount);
00355   totASS_vertcount_fwdneg.push_back(totASSv_vertcount);
00356   totASS2_vertcount_fwdneg.push_back(totASS2v_vertcount);
00357   totREC_vertcount_fwdneg.push_back(totRECv_vertcount);
00358   totFOMT_vertcount.push_back(totASSv_vertcount);
00359     
00360   totASS2_itpu_vertcount_entire.push_back(totASS2v_vertcount);
00361   totASS2_itpu_vertcount_entire_signal.push_back(totASS2v_vertcount);
00362 
00363   totASS2_ootpu_entire.push_back(totASS2v_vertcount);
00364   totREC_ootpu_entire.push_back(totRECv_vertcount);
00365   totASS2_ootpu_barrel.push_back(totASS2v_vertcount);
00366   totREC_ootpu_barrel.push_back(totRECv_vertcount);
00367   totASS2_ootpu_fwdpos.push_back(totASS2v_vertcount);
00368   totREC_ootpu_fwdpos.push_back(totRECv_vertcount);
00369   totASS2_ootpu_fwdneg.push_back(totASS2v_vertcount);
00370   totREC_ootpu_fwdneg.push_back(totRECv_vertcount);
00371 
00372 }
00373 
00374 void MTVHistoProducerAlgoForTracker::bookSimHistos(){
00375   h_ptSIM.push_back( dbe_->book1D("ptSIM", "generated p_{t}", 5500, 0, 110 ) );
00376   h_etaSIM.push_back( dbe_->book1D("etaSIM", "generated pseudorapidity", 500, -2.5, 2.5 ) );
00377   h_tracksSIM.push_back( dbe_->book1D("tracksSIM","number of simulated tracks",200,-0.5,99.5) );
00378   h_vertposSIM.push_back( dbe_->book1D("vertposSIM","Transverse position of sim vertices",100,0.,120.) );  
00379   h_bunchxSIM.push_back( dbe_->book1D("bunchxSIM", "bunch crossing", 22, -5, 5 ) );
00380 }
00381 
00382 
00383 void MTVHistoProducerAlgoForTracker::bookRecoHistos(){
00384   h_tracks.push_back( dbe_->book1D("tracks","number of reconstructed tracks",401,-0.5,400.5) );
00385   h_fakes.push_back( dbe_->book1D("fakes","number of fake reco tracks",101,-0.5,100.5) );
00386   h_charge.push_back( dbe_->book1D("charge","charge",3,-1.5,1.5) );
00387   
00388   h_hits.push_back( dbe_->book1D("hits", "number of hits per track", nintHit,minHit,maxHit ) );
00389   h_losthits.push_back( dbe_->book1D("losthits", "number of lost hits per track", nintHit,minHit,maxHit) );
00390   h_nchi2.push_back( dbe_->book1D("chi2", "normalized #chi^{2}", 200, 0, 20 ) );
00391   h_nchi2_prob.push_back( dbe_->book1D("chi2_prob", "normalized #chi^{2} probability",100,0,1));
00392 
00393   h_algo.push_back( dbe_->book1D("h_algo","Tracks by algo",15,0.0,15.0) );
00394 
00396   h_recoeta.push_back( dbe_->book1D("num_reco_eta","N of reco track vs eta",nintEta,minEta,maxEta) );
00397   h_assoceta.push_back( dbe_->book1D("num_assoc(simToReco)_eta","N of associated tracks (simToReco) vs eta",nintEta,minEta,maxEta) );
00398   h_assoc2eta.push_back( dbe_->book1D("num_assoc(recoToSim)_eta","N of associated (recoToSim) tracks vs eta",nintEta,minEta,maxEta) );
00399   h_simuleta.push_back( dbe_->book1D("num_simul_eta","N of simulated tracks vs eta",nintEta,minEta,maxEta) );
00400   h_loopereta.push_back( dbe_->book1D("num_duplicate_eta","N of associated (recoToSim) duplicate tracks vs eta",nintEta,minEta,maxEta) );
00401   h_misideta.push_back( dbe_->book1D("num_chargemisid_eta","N of associated (recoToSim) charge misIDed tracks vs eta",nintEta,minEta,maxEta) );
00402   h_recopT.push_back( dbe_->book1D("num_reco_pT","N of reco track vs pT",nintPt,minPt,maxPt) );
00403   h_assocpT.push_back( dbe_->book1D("num_assoc(simToReco)_pT","N of associated tracks (simToReco) vs pT",nintPt,minPt,maxPt) );
00404   h_assoc2pT.push_back( dbe_->book1D("num_assoc(recoToSim)_pT","N of associated (recoToSim) tracks vs pT",nintPt,minPt,maxPt) );
00405   h_simulpT.push_back( dbe_->book1D("num_simul_pT","N of simulated tracks vs pT",nintPt,minPt,maxPt) );
00406   h_looperpT.push_back( dbe_->book1D("num_duplicate_pT","N of associated (recoToSim) duplicate tracks vs pT",nintPt,minPt,maxPt) );
00407   h_misidpT.push_back( dbe_->book1D("num_chargemisid_pT","N of associated (recoToSim) charge misIDed tracks vs pT",nintPt,minPt,maxPt) );
00408   //
00409   h_recohit.push_back( dbe_->book1D("num_reco_hit","N of reco track vs hit",nintHit,minHit,maxHit) );
00410   h_assochit.push_back( dbe_->book1D("num_assoc(simToReco)_hit","N of associated tracks (simToReco) vs hit",nintHit,minHit,maxHit) );
00411   h_assoc2hit.push_back( dbe_->book1D("num_assoc(recoToSim)_hit","N of associated (recoToSim) tracks vs hit",nintHit,minHit,maxHit) );
00412   h_simulhit.push_back( dbe_->book1D("num_simul_hit","N of simulated tracks vs hit",nintHit,minHit,maxHit) );
00413   h_looperhit.push_back( dbe_->book1D("num_duplicate_hit","N of associated (recoToSim) duplicate tracks vs hit",nintHit,minHit,maxHit) );
00414   h_misidhit.push_back( dbe_->book1D("num_chargemisid_hit","N of associated (recoToSim) charge misIDed tracks vs hit",nintHit,minHit,maxHit) );
00415   //
00416   h_recophi.push_back( dbe_->book1D("num_reco_phi","N of reco track vs phi",nintPhi,minPhi,maxPhi) );
00417   h_assocphi.push_back( dbe_->book1D("num_assoc(simToReco)_phi","N of associated tracks (simToReco) vs phi",nintPhi,minPhi,maxPhi) );
00418   h_assoc2phi.push_back( dbe_->book1D("num_assoc(recoToSim)_phi","N of associated (recoToSim) tracks vs phi",nintPhi,minPhi,maxPhi) );
00419   h_simulphi.push_back( dbe_->book1D("num_simul_phi","N of simulated tracks vs phi",nintPhi,minPhi,maxPhi) );
00420   h_looperphi.push_back( dbe_->book1D("num_duplicate_phi","N of associated (recoToSim) duplicate tracks vs phi",nintPhi,minPhi,maxPhi) );
00421   h_misidphi.push_back( dbe_->book1D("num_chargemisid_phi","N of associated (recoToSim) charge misIDed tracks vs phi",nintPhi,minPhi,maxPhi) );
00422   
00423   h_recodxy.push_back( dbe_->book1D("num_reco_dxy","N of reco track vs dxy",nintDxy,minDxy,maxDxy) );
00424   h_assocdxy.push_back( dbe_->book1D("num_assoc(simToReco)_dxy","N of associated tracks (simToReco) vs dxy",nintDxy,minDxy,maxDxy) );
00425   h_assoc2dxy.push_back( dbe_->book1D("num_assoc(recoToSim)_dxy","N of associated (recoToSim) tracks vs dxy",nintDxy,minDxy,maxDxy) );
00426   h_simuldxy.push_back( dbe_->book1D("num_simul_dxy","N of simulated tracks vs dxy",nintDxy,minDxy,maxDxy) );
00427   h_looperdxy.push_back( dbe_->book1D("num_duplicate_dxy","N of associated (recoToSim) looper tracks vs dxy",nintDxy,minDxy,maxDxy) );
00428   h_misiddxy.push_back( dbe_->book1D("num_chargemisid_dxy","N of associated (recoToSim) charge misIDed tracks vs dxy",nintDxy,minDxy,maxDxy) );
00429   
00430   h_recodz.push_back( dbe_->book1D("num_reco_dz","N of reco track vs dz",nintDz,minDz,maxDz) );
00431   h_assocdz.push_back( dbe_->book1D("num_assoc(simToReco)_dz","N of associated tracks (simToReco) vs dz",nintDz,minDz,maxDz) );
00432   h_assoc2dz.push_back( dbe_->book1D("num_assoc(recoToSim)_dz","N of associated (recoToSim) tracks vs dz",nintDz,minDz,maxDz) );
00433   h_simuldz.push_back( dbe_->book1D("num_simul_dz","N of simulated tracks vs dz",nintDz,minDz,maxDz) );
00434   h_looperdz.push_back( dbe_->book1D("num_duplicate_dz","N of associated (recoToSim) looper tracks vs dz",nintDz,minDz,maxDz) );
00435   h_misiddz.push_back( dbe_->book1D("num_chargemisid_versus_dz","N of associated (recoToSim) charge misIDed tracks vs dz",nintDz,minDz,maxDz) );
00436   
00437   h_assocvertpos.push_back( dbe_->book1D("num_assoc(simToReco)_vertpos",
00438                                          "N of associated tracks (simToReco) vs transverse vert position",       
00439                                          nintVertpos,minVertpos,maxVertpos) );
00440   h_simulvertpos.push_back( dbe_->book1D("num_simul_vertpos","N of simulated tracks vs transverse vert position",
00441                                          nintVertpos,minVertpos,maxVertpos) );
00442   
00443   h_assoczpos.push_back( dbe_->book1D("num_assoc(simToReco)_zpos","N of associated tracks (simToReco) vs z vert position",
00444                                       nintZpos,minZpos,maxZpos) );
00445   h_simulzpos.push_back( dbe_->book1D("num_simul_zpos","N of simulated tracks vs z vert position",nintZpos,minZpos,maxZpos) );
00446   
00447 
00448   h_reco_vertcount_entire.push_back( dbe_->book1D("num_reco_vertcount_entire","N of reco tracks vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00449   h_assoc_vertcount_entire.push_back( dbe_->book1D("num_assoc(simToReco)_vertcount_entire","N of associated tracks (simToReco) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00450   h_assoc2_vertcount_entire.push_back( dbe_->book1D("num_assoc(recoToSim)_vertcount_entire","N of associated (recoToSim) tracks vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00451   h_simul_vertcount_entire.push_back( dbe_->book1D("num_simul_vertcount_entire","N of simulated tracks vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00452 
00453   h_reco_vertcount_barrel.push_back( dbe_->book1D("num_reco_vertcount_barrel","N of reco tracks in barrel vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00454   h_assoc_vertcount_barrel.push_back( dbe_->book1D("num_assoc(simToReco)_vertcount_barrel","N of associated tracks (simToReco) in barrel vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00455   h_assoc2_vertcount_barrel.push_back( dbe_->book1D("num_assoc(recoToSim)_vertcount_barrel","N of associated (recoToSim) tracks in barrel vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00456   h_simul_vertcount_barrel.push_back( dbe_->book1D("num_simul_vertcount_barrel","N of simulated tracks in barrel vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00457 
00458   h_reco_vertcount_fwdpos.push_back( dbe_->book1D("num_reco_vertcount_fwdpos","N of reco tracks in endcap(+) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00459   h_assoc_vertcount_fwdpos.push_back( dbe_->book1D("num_assoc(simToReco)_vertcount_fwdpos","N of associated tracks (simToReco) in endcap(+) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00460   h_assoc2_vertcount_fwdpos.push_back( dbe_->book1D("num_assoc(recoToSim)_vertcount_fwdpos","N of associated (recoToSim) tracks in endcap(+) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00461   h_simul_vertcount_fwdpos.push_back( dbe_->book1D("num_simul_vertcount_fwdpos","N of simulated tracks in endcap(+) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00462 
00463   h_reco_vertcount_fwdneg.push_back( dbe_->book1D("num_reco_vertcount_fwdneg","N of reco tracks in endcap(-) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00464   h_assoc_vertcount_fwdneg.push_back( dbe_->book1D("num_assoc(simToReco)_vertcount_fwdneg","N of associated tracks (simToReco) in endcap(-) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00465   h_assoc2_vertcount_fwdneg.push_back( dbe_->book1D("num_assoc(recoToSim)_vertcount_fwdneg","N of associated (recoToSim) tracks in endcap(-) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00466   h_simul_vertcount_fwdneg.push_back( dbe_->book1D("num_simul_vertcount_fwdneg","N of simulated tracks in endcap(-) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00467 
00468   h_assoc_vertz_entire.push_back( dbe_->book1D("num_assoc(simToReco)_vertz_entire","N of associated tracks (simToReco) in entire vs z of primary intercation vertex",nintZpos,minZpos,maxZpos) );
00469   h_simul_vertz_entire.push_back( dbe_->book1D("num_simul_vertz_entire","N of simulated tracks in entire vs N of pileup vertices",nintZpos,minZpos,maxZpos) );
00470 
00471   h_assoc_vertz_barrel.push_back( dbe_->book1D("num_assoc(simToReco)_vertz_barrel","N of associated tracks (simToReco) in barrel vs z of primary intercation vertex",nintZpos,minZpos,maxZpos) );
00472   h_simul_vertz_barrel.push_back( dbe_->book1D("num_simul_vertz_barrel","N of simulated tracks in barrel vs N of pileup vertices",nintZpos,minZpos,maxZpos) );
00473 
00474   h_assoc_vertz_fwdpos.push_back( dbe_->book1D("num_assoc(simToReco)_vertz_fwdpos","N of associated tracks (simToReco) in endcap(+) vs z of primary interaction vertex",nintZpos,minZpos,maxZpos) );
00475   h_simul_vertz_fwdpos.push_back( dbe_->book1D("num_simul_vertz_fwdpos","N of simulated tracks in endcap(+) vs z of primary interaction vertex",nintZpos,minZpos,maxZpos) );
00476 
00477   h_assoc_vertz_fwdneg.push_back( dbe_->book1D("num_assoc(simToReco)_vertz_fwdneg","N of associated tracks (simToReco) in endcap(-) vs N of pileup vertices",nintZpos,minZpos,maxZpos) );
00478   h_simul_vertz_fwdneg.push_back( dbe_->book1D("num_simul_vertz_fwdneg","N of simulated tracks in endcap(-) vs z of primary interaction vertex",nintZpos,minZpos,maxZpos) );
00479 /*
00480   h_assoc2_itpu_eta.push_back( dbe_->book1D("num_assoc(recoToSim)_itpu_eta_entire","N of associated tracks (simToReco) from in time pileup vs eta",nintEta,minEta,maxEta) );
00481   h_assoc2_itpu_sig_eta.push_back( dbe_->book1D("num_assoc(recoToSim)_itpu_eta_entire_signal","N of associated tracks (simToReco) from in time pileup vs eta",nintEta,minEta,maxEta) );
00482   h_assoc2_itpu_vertcount.push_back( dbe_->book1D("num_assoc(recoToSim)_itpu_vertcount_entire","N of associated tracks (simToReco) from in time pileup vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00483   h_assoc2_itpu_sig_vertcount.push_back( dbe_->book1D("num_assoc(recoToSim)_itpu_vertcount_entire_signal","N of associated tracks (simToReco) from in time pileup vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );*/
00484 
00485   h_reco_ootpu_eta.push_back( dbe_->book1D("num_reco_ootpu_eta_entire","N of reco tracks vs eta",nintEta,minEta,maxEta) );
00486   h_reco_ootpu_vertcount.push_back( dbe_->book1D("num_reco_ootpu_vertcount_entire","N of reco tracks vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00487 
00488   h_assoc2_ootpu_entire.push_back( dbe_->book1D("num_assoc(recoToSim)_ootpu_eta_entire","N of associated tracks (simToReco) from out of time pileup vs eta",nintEta,minEta,maxEta) );
00489   h_assoc2_ootpu_vertcount.push_back( dbe_->book1D("num_assoc(recoToSim)_ootpu_vertcount_entire","N of associated tracks (simToReco) from out of time pileup vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00490 
00491   h_reco_ootpu_entire.push_back( dbe_->book1D("num_reco_ootpu_entire","N of reco tracks vs z of primary interaction vertex",nintVertcount,minVertcount,maxVertcount) );
00492   h_assoc2_ootpu_barrel.push_back( dbe_->book1D("num_assoc(recoToSim)_ootpu_barrel","N of associated tracks (simToReco) from out of time pileup in barrel vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00493   h_reco_ootpu_barrel.push_back( dbe_->book1D("num_reco_ootpu_barrel","N of reco tracks in barrel vs z of primary interaction vertex",nintVertcount,minVertcount,maxVertcount) );
00494   h_assoc2_ootpu_fwdpos.push_back( dbe_->book1D("num_assoc(recoToSim)_ootpu_fwdpos","N of associated tracks (simToReco) from out of time pileup in endcap(+) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00495   h_reco_ootpu_fwdpos.push_back( dbe_->book1D("num_reco_ootpu_fwdpos","N of reco tracks in endcap(+) vs z of primary interaction vertex",nintVertcount,minVertcount,maxVertcount) );
00496   h_assoc2_ootpu_fwdneg.push_back( dbe_->book1D("num_assoc(recoToSim)_ootpu_fwdneg","N of associated tracks (simToReco) from out of time pileup in endcap(-) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00497   h_reco_ootpu_fwdneg.push_back( dbe_->book1D("num_reco_ootpu_fwdneg","N of reco tracks in endcap(-) vs z of primary interaction vertex",nintVertcount,minVertcount,maxVertcount) );
00498 
00499 
00501   
00502   h_eta.push_back( dbe_->book1D("eta", "pseudorapidity residue", 1000, -0.1, 0.1 ) );
00503   h_pt.push_back( dbe_->book1D("pullPt", "pull of p_{t}", 100, -10, 10 ) );
00504   h_pullTheta.push_back( dbe_->book1D("pullTheta","pull of #theta parameter",250,-25,25) );
00505   h_pullPhi.push_back( dbe_->book1D("pullPhi","pull of #phi parameter",250,-25,25) );
00506   h_pullDxy.push_back( dbe_->book1D("pullDxy","pull of dxy parameter",250,-25,25) );
00507   h_pullDz.push_back( dbe_->book1D("pullDz","pull of dz parameter",250,-25,25) );
00508   h_pullQoverp.push_back( dbe_->book1D("pullQoverp","pull of qoverp parameter",250,-25,25) );
00509   
00510   /* TO BE FIXED -----------
00511   if (associators[ww]=="TrackAssociatorByChi2"){
00512     h_assochi2.push_back( dbe_->book1D("assocChi2","track association #chi^{2}",1000000,0,100000) );
00513     h_assochi2_prob.push_back(dbe_->book1D("assocChi2_prob","probability of association #chi^{2}",100,0,1));
00514   } else if (associators[ww]=="TrackAssociatorByHits"){
00515     h_assocFraction.push_back( dbe_->book1D("assocFraction","fraction of shared hits",200,0,2) );
00516     h_assocSharedHit.push_back(dbe_->book1D("assocSharedHit","number of shared hits",20,0,20));
00517   }
00518   */
00519   h_assocFraction.push_back( dbe_->book1D("assocFraction","fraction of shared hits",200,0,2) );
00520   h_assocSharedHit.push_back(dbe_->book1D("assocSharedHit","number of shared hits",41,-0.5,40.5));
00521   // ----------------------
00522 
00523   chi2_vs_nhits.push_back( dbe_->book2D("chi2_vs_nhits","#chi^{2} vs nhits",25,0,25,100,0,10) );
00524   
00525   etares_vs_eta.push_back( dbe_->book2D("etares_vs_eta","etaresidue vs eta",nintEta,minEta,maxEta,200,-0.1,0.1) );
00526   nrec_vs_nsim.push_back( dbe_->book2D("nrec_vs_nsim","nrec vs nsim",401,-0.5,400.5,401,-0.5,400.5) );
00527   
00528   chi2_vs_eta.push_back( dbe_->book2D("chi2_vs_eta","chi2_vs_eta",nintEta,minEta,maxEta, 200, 0, 20 ));
00529   chi2_vs_phi.push_back( dbe_->book2D("chi2_vs_phi","#chi^{2} vs #phi",nintPhi,minPhi,maxPhi, 200, 0, 20 ) );
00530   
00531   nhits_vs_eta.push_back( dbe_->book2D("nhits_vs_eta","nhits vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00532   nPXBhits_vs_eta.push_back( dbe_->book2D("nPXBhits_vs_eta","# PXB its vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00533   nPXFhits_vs_eta.push_back( dbe_->book2D("nPXFhits_vs_eta","# PXF hits vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00534   nTIBhits_vs_eta.push_back( dbe_->book2D("nTIBhits_vs_eta","# TIB hits vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00535   nTIDhits_vs_eta.push_back( dbe_->book2D("nTIDhits_vs_eta","# TID hits vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00536   nTOBhits_vs_eta.push_back( dbe_->book2D("nTOBhits_vs_eta","# TOB hits vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00537   nTEChits_vs_eta.push_back( dbe_->book2D("nTEChits_vs_eta","# TEC hits vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00538 
00539   nLayersWithMeas_vs_eta.push_back( dbe_->book2D("nLayersWithMeas_vs_eta","# Layers with measurement vs eta",
00540                                                 nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00541   nPXLlayersWithMeas_vs_eta.push_back( dbe_->book2D("nPXLlayersWithMeas_vs_eta","# PXL Layers with measurement vs eta",
00542                                                    nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00543   nSTRIPlayersWithMeas_vs_eta.push_back( dbe_->book2D("nSTRIPlayersWithMeas_vs_eta","# STRIP Layers with measurement vs eta",
00544                                                       nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00545   nSTRIPlayersWith1dMeas_vs_eta.push_back( dbe_->book2D("nSTRIPlayersWith1dMeas_vs_eta","# STRIP Layers with 1D measurement vs eta",
00546                                                         nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00547   nSTRIPlayersWith2dMeas_vs_eta.push_back( dbe_->book2D("nSTRIPlayersWith2dMeas_vs_eta","# STRIP Layers with 2D measurement vs eta",
00548                                                        nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00549   
00550   nhits_vs_phi.push_back( dbe_->book2D("nhits_vs_phi","#hits vs #phi",nintPhi,minPhi,maxPhi,nintHit,minHit,maxHit) );
00551   
00552   nlosthits_vs_eta.push_back( dbe_->book2D("nlosthits_vs_eta","nlosthits vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00553   
00554   //resolution of track parameters
00555   //                       dPt/Pt    cotTheta        Phi            TIP            LIP
00556   // log10(pt)<0.5        100,0.1    240,0.08     100,0.015      100,0.1000    150,0.3000
00557   // 0.5<log10(pt)<1.5    100,0.1    120,0.01     100,0.003      100,0.0100    150,0.0500
00558   // >1.5                 100,0.3    100,0.005    100,0.0008     100,0.0060    120,0.0300
00559   
00560   ptres_vs_eta.push_back(dbe_->book2D("ptres_vs_eta","ptres_vs_eta",
00561                                       nintEta,minEta,maxEta, ptRes_nbin, ptRes_rangeMin, ptRes_rangeMax));
00562   
00563   ptres_vs_phi.push_back( dbe_->book2D("ptres_vs_phi","p_{t} res vs #phi",
00564                                        nintPhi,minPhi,maxPhi, ptRes_nbin, ptRes_rangeMin, ptRes_rangeMax));
00565   
00566   ptres_vs_pt.push_back(dbe_->book2D("ptres_vs_pt","ptres_vs_pt",nintPt,minPt,maxPt, ptRes_nbin, ptRes_rangeMin, ptRes_rangeMax));
00567   
00568   cotThetares_vs_eta.push_back(dbe_->book2D("cotThetares_vs_eta","cotThetares_vs_eta",
00569                                             nintEta,minEta,maxEta,cotThetaRes_nbin, cotThetaRes_rangeMin, cotThetaRes_rangeMax));
00570 
00571   
00572   cotThetares_vs_pt.push_back(dbe_->book2D("cotThetares_vs_pt","cotThetares_vs_pt",
00573                                            nintPt,minPt,maxPt, cotThetaRes_nbin, cotThetaRes_rangeMin, cotThetaRes_rangeMax));      
00574 
00575 
00576   phires_vs_eta.push_back(dbe_->book2D("phires_vs_eta","phires_vs_eta",
00577                                        nintEta,minEta,maxEta, phiRes_nbin, phiRes_rangeMin, phiRes_rangeMax));
00578 
00579   phires_vs_pt.push_back(dbe_->book2D("phires_vs_pt","phires_vs_pt",
00580                                       nintPt,minPt,maxPt, phiRes_nbin, phiRes_rangeMin, phiRes_rangeMax));
00581 
00582   phires_vs_phi.push_back(dbe_->book2D("phires_vs_phi","#phi res vs #phi",
00583                                        nintPhi,minPhi,maxPhi,phiRes_nbin, phiRes_rangeMin, phiRes_rangeMax));
00584 
00585   dxyres_vs_eta.push_back(dbe_->book2D("dxyres_vs_eta","dxyres_vs_eta",
00586                                        nintEta,minEta,maxEta,dxyRes_nbin, dxyRes_rangeMin, dxyRes_rangeMax));
00587   
00588   dxyres_vs_pt.push_back( dbe_->book2D("dxyres_vs_pt","dxyres_vs_pt",
00589                                        nintPt,minPt,maxPt,dxyRes_nbin, dxyRes_rangeMin, dxyRes_rangeMax));
00590   
00591   dzres_vs_eta.push_back(dbe_->book2D("dzres_vs_eta","dzres_vs_eta",
00592                                       nintEta,minEta,maxEta,dzRes_nbin, dzRes_rangeMin, dzRes_rangeMax));
00593 
00594   dzres_vs_pt.push_back(dbe_->book2D("dzres_vs_pt","dzres_vs_pt",nintPt,minPt,maxPt,dzRes_nbin, dzRes_rangeMin, dzRes_rangeMax));
00595   
00596   ptmean_vs_eta_phi.push_back(dbe_->bookProfile2D("ptmean_vs_eta_phi","mean p_{t} vs #eta and #phi",
00597                                                   nintPhi,minPhi,maxPhi,nintEta,minEta,maxEta,1000,0,1000));
00598   phimean_vs_eta_phi.push_back(dbe_->bookProfile2D("phimean_vs_eta_phi","mean #phi vs #eta and #phi",
00599                                                    nintPhi,minPhi,maxPhi,nintEta,minEta,maxEta,nintPhi,minPhi,maxPhi));
00600   
00601   //pulls of track params vs eta: to be used with fitslicesytool
00602   dxypull_vs_eta.push_back(dbe_->book2D("dxypull_vs_eta","dxypull_vs_eta",nintEta,minEta,maxEta,100,-10,10));
00603   ptpull_vs_eta.push_back(dbe_->book2D("ptpull_vs_eta","ptpull_vs_eta",nintEta,minEta,maxEta,100,-10,10)); 
00604   dzpull_vs_eta.push_back(dbe_->book2D("dzpull_vs_eta","dzpull_vs_eta",nintEta,minEta,maxEta,100,-10,10)); 
00605   phipull_vs_eta.push_back(dbe_->book2D("phipull_vs_eta","phipull_vs_eta",nintEta,minEta,maxEta,100,-10,10)); 
00606   thetapull_vs_eta.push_back(dbe_->book2D("thetapull_vs_eta","thetapull_vs_eta",nintEta,minEta,maxEta,100,-10,10));
00607   
00608   //      h_ptshiftetamean.push_back( dbe_->book1D("h_ptshifteta_Mean","<#deltapT/pT>[%] vs #eta",nintEta,minEta,maxEta) ); 
00609   
00610 
00611   //pulls of track params vs phi
00612   ptpull_vs_phi.push_back(dbe_->book2D("ptpull_vs_phi","p_{t} pull vs #phi",nintPhi,minPhi,maxPhi,100,-10,10)); 
00613   phipull_vs_phi.push_back(dbe_->book2D("phipull_vs_phi","#phi pull vs #phi",nintPhi,minPhi,maxPhi,100,-10,10)); 
00614   thetapull_vs_phi.push_back(dbe_->book2D("thetapull_vs_phi","#theta pull vs #phi",nintPhi,minPhi,maxPhi,100,-10,10));
00615 
00616   
00617   nrecHit_vs_nsimHit_sim2rec.push_back( dbe_->book2D("nrecHit_vs_nsimHit_sim2rec","nrecHit vs nsimHit (Sim2RecAssoc)",
00618                                                      nintHit,minHit,maxHit, nintHit,minHit,maxHit ));
00619   nrecHit_vs_nsimHit_rec2sim.push_back( dbe_->book2D("nrecHit_vs_nsimHit_rec2sim","nrecHit vs nsimHit (Rec2simAssoc)",
00620                                                      nintHit,minHit,maxHit, nintHit,minHit,maxHit ));
00621 
00622   // dE/dx stuff
00623   // FIXME: it would be nice to have an array
00624   h_dedx_estim1.push_back( dbe_->book1D("h_dedx_estim1","dE/dx estimator 1",nintDeDx,minDeDx,maxDeDx) ); 
00625   h_dedx_estim2.push_back( dbe_->book1D("h_dedx_estim2","dE/dx estimator 2",nintDeDx,minDeDx,maxDeDx) ); 
00626   h_dedx_nom1.push_back( dbe_->book1D("h_dedx_nom1","dE/dx number of measurements",nintHit,minHit,maxHit) ); 
00627   h_dedx_nom2.push_back( dbe_->book1D("h_dedx_nom2","dE/dx number of measurements",nintHit,minHit,maxHit) ); 
00628   h_dedx_sat1.push_back( dbe_->book1D("h_dedx_sat1","dE/dx number of measurements with saturation",nintHit,minHit,maxHit) ); 
00629   h_dedx_sat2.push_back( dbe_->book1D("h_dedx_sat2","dE/dx number of measurements with saturation",nintHit,minHit,maxHit) ); 
00630 
00631   // PU special stuff
00632   h_con_eta.push_back( dbe_->book1D("num_con_eta","N of PU tracks vs eta",nintEta,minEta,maxEta) );
00633   h_con_vertcount.push_back( dbe_->book1D("num_con_vertcount","N of PU tracks vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00634   h_con_zpos.push_back( dbe_->book1D("num_con_zpos","N of PU tracks vs z of primary interaction vertex",nintZpos,minZpos,maxZpos) );
00635 
00636 
00637   if(useLogPt){
00638     BinLogX(dzres_vs_pt.back()->getTH2F());
00639     BinLogX(dxyres_vs_pt.back()->getTH2F());
00640     BinLogX(phires_vs_pt.back()->getTH2F());
00641     BinLogX(cotThetares_vs_pt.back()->getTH2F());
00642     BinLogX(ptres_vs_pt.back()->getTH2F());
00643     BinLogX(h_looperpT.back()->getTH1F());
00644     BinLogX(h_misidpT.back()->getTH1F());
00645     BinLogX(h_recopT.back()->getTH1F());
00646     BinLogX(h_assocpT.back()->getTH1F());
00647     BinLogX(h_assoc2pT.back()->getTH1F());
00648     BinLogX(h_simulpT.back()->getTH1F());
00649   }  
00650 }
00651 
00652 void MTVHistoProducerAlgoForTracker::bookRecoHistosForStandaloneRunning(){
00653   h_effic.push_back( dbe_->book1D("effic","efficiency vs #eta",nintEta,minEta,maxEta) );
00654   h_efficPt.push_back( dbe_->book1D("efficPt","efficiency vs pT",nintPt,minPt,maxPt) );
00655   h_effic_vs_hit.push_back( dbe_->book1D("effic_vs_hit","effic vs hit",nintHit,minHit,maxHit) );
00656   h_effic_vs_phi.push_back( dbe_->book1D("effic_vs_phi","effic vs phi",nintPhi,minPhi,maxPhi) );
00657   h_effic_vs_dxy.push_back( dbe_->book1D("effic_vs_dxy","effic vs dxy",nintDxy,minDxy,maxDxy) );
00658   h_effic_vs_dz.push_back( dbe_->book1D("effic_vs_dz","effic vs dz",nintDz,minDz,maxDz) );
00659   h_effic_vs_vertpos.push_back( dbe_->book1D("effic_vs_vertpos","effic vs vertpos",nintVertpos,minVertpos,maxVertpos) );
00660   h_effic_vs_zpos.push_back( dbe_->book1D("effic_vs_zpos","effic vs zpos",nintZpos,minZpos,maxZpos) );
00661   h_effic_vertcount_entire.push_back( dbe_->book1D("effic_vertcount_entire","efficiency vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00662   h_effic_vertcount_barrel.push_back( dbe_->book1D("effic_vertcount_barrel","efficiency in barrel vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00663   h_effic_vertcount_fwdpos.push_back( dbe_->book1D("effic_vertcount_fwdpos","efficiency in endcap(+) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00664   h_effic_vertcount_fwdneg.push_back( dbe_->book1D("effic_vertcount_fwdneg","efficiency in endcap(-) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00665   h_effic_vertz_entire.push_back( dbe_->book1D("effic_vertz_entire","efficiency vs z of primary interaction vertex",nintZpos,minZpos,maxZpos) );
00666   h_effic_vertz_barrel.push_back( dbe_->book1D("effic_vertz_barrel","efficiency in barrel vs z of primary interaction vertex",nintZpos,minZpos,maxZpos) );
00667   h_effic_vertz_fwdpos.push_back( dbe_->book1D("effic_vertz_fwdpos","efficiency in endcap(+) vs z of primary interaction vertex",nintZpos,minZpos,maxZpos) );
00668   h_effic_vertz_fwdneg.push_back( dbe_->book1D("effic_vertz_fwdneg","efficiency in endcap(-) vs z of primary interaction vertex",nintZpos,minZpos,maxZpos) );
00669 
00670   h_looprate.push_back( dbe_->book1D("duplicatesRate","loop rate vs #eta",nintEta,minEta,maxEta) );
00671   h_misidrate.push_back( dbe_->book1D("chargeMisIdRate","misid rate vs #eta",nintEta,minEta,maxEta) );
00672   h_fakerate.push_back( dbe_->book1D("fakerate","fake rate vs #eta",nintEta,minEta,maxEta) );
00673   h_loopratepT.push_back( dbe_->book1D("duplicatesRate_Pt","loop rate vs pT",nintPt,minPt,maxPt) );
00674   h_misidratepT.push_back( dbe_->book1D("chargeMisIdRate_Pt","misid rate vs pT",nintPt,minPt,maxPt) );
00675   h_fakeratePt.push_back( dbe_->book1D("fakeratePt","fake rate vs pT",nintPt,minPt,maxPt) );
00676   h_loopratehit.push_back( dbe_->book1D("duplicatesRate_hit","loop rate vs hit",nintHit,minHit,maxHit) );
00677   h_misidratehit.push_back( dbe_->book1D("chargeMisIdRate_hit","misid rate vs hit",nintHit,minHit,maxHit) );
00678   h_fake_vs_hit.push_back( dbe_->book1D("fakerate_vs_hit","fake rate vs hit",nintHit,minHit,maxHit) );
00679   h_loopratephi.push_back( dbe_->book1D("duplicatesRate_phi","loop rate vs #phi",nintPhi,minPhi,maxPhi) );
00680   h_misidratephi.push_back( dbe_->book1D("chargeMisIdRate_phi","misid rate vs #phi",nintPhi,minPhi,maxPhi) );
00681   h_fake_vs_phi.push_back( dbe_->book1D("fakerate_vs_phi","fake vs #phi",nintPhi,minPhi,maxPhi) );
00682   h_loopratedxy.push_back( dbe_->book1D("duplicatesRate_dxy","loop rate vs dxy",nintDxy,minDxy,maxDxy) );
00683   h_misidratedxy.push_back( dbe_->book1D("chargeMisIdRate_dxy","misid rate vs dxy",nintDxy,minDxy,maxDxy) );
00684   h_fake_vs_dxy.push_back( dbe_->book1D("fakerate_vs_dxy","fake rate vs dxy",nintDxy,minDxy,maxDxy) );
00685   h_loopratedz.push_back( dbe_->book1D("duplicatesRate_dz","loop rate vs dz",nintDz,minDz,maxDz) );
00686   h_misidratedz.push_back( dbe_->book1D("chargeMisIdRate_dz","misid rate vs dz",nintDz,minDz,maxDz) );
00687   h_fake_vs_dz.push_back( dbe_->book1D("fakerate_vs_dz","fake vs dz",nintDz,minDz,maxDz) );
00688   h_fakerate_vertcount_entire.push_back( dbe_->book1D("fakerate_vertcount_entire","fake rate vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00689   h_fakerate_vertcount_barrel.push_back( dbe_->book1D("fakerate_vertcount_barrel","fake rate in barrel vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00690   h_fakerate_vertcount_fwdpos.push_back( dbe_->book1D("fakerate_vertcount_fwdpos","fake rate in endcap(+) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00691   h_fakerate_vertcount_fwdneg.push_back( dbe_->book1D("fakerate_vertcount_fwdneg","fake rate in endcap(-) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00692   h_fakerate_ootpu_entire.push_back( dbe_->book1D("fakerate_ootpu_eta_entire","Out of time pileup fake rate vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00693   h_fakerate_ootpu_barrel.push_back( dbe_->book1D("fakerate_ootpu_barrel","Out of time pileup fake rate in barrel vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00694   h_fakerate_ootpu_fwdpos.push_back( dbe_->book1D("fakerate_ootpu_fwdpos","Out of time pileup fake rate in endcap(+) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00695   h_fakerate_ootpu_fwdneg.push_back( dbe_->book1D("fakerate_ootpu_fwdneg","Out of time pileup fake rate in endcap(-) vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00696 
00697   h_fomt_eta.push_back( dbe_->book1D("fomt","fraction of misreconstructed tracks vs #eta",nintEta,minEta,maxEta) );
00698   h_fomt_sig_eta.push_back( dbe_->book1D("fomt_signal","fraction of misreconstructed tracks vs #eta",nintEta,minEta,maxEta) );
00699   h_fomt_vertcount.push_back( dbe_->book1D("fomt_vertcount","fraction of misreconstructed tracks vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00700   h_fomt_sig_vertcount.push_back( dbe_->book1D("fomt_vertcount_signal","fraction of misreconstructed tracks vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00701   h_fomt_ootpu_eta.push_back( dbe_->book1D("fomt_ootpu_eta","Out of time pileup fraction of misreconstructed tracks vs #eta",nintEta,minEta,maxEta) );
00702   h_fomt_ootpu_vertcount.push_back( dbe_->book1D("fomt_ootpu_vertcount","Out of time pileup fraction of misreconstructed tracks vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00703   h_fomt_itpu_eta.push_back( dbe_->book1D("fomt_itpu_eta","In time pileup fraction of misreconstructed tracks vs eta",nintEta,minEta,maxEta) );
00704   h_fomt_itpu_vertcount.push_back( dbe_->book1D("fomt_itpu_vertcount","In time pileup fraction of misreconstructed tracks vs N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00705 
00706   h_effic_PU_eta.push_back( dbe_->book1D("effic_PU_eta","PU efficiency vs #eta",nintEta,minEta,maxEta) );
00707   h_effic_PU_vertcount.push_back( dbe_->book1D("effic_PU_vertcount","PU efficiency s N of pileup vertices",nintVertcount,minVertcount,maxVertcount) );
00708   h_effic_PU_zpos.push_back( dbe_->book1D("effic_PU_zpos","PU efficiency vs z of primary interaction vertex",nintZpos,minZpos,maxZpos) );
00709 
00710   h_chi2meanhitsh.push_back( dbe_->bookProfile("chi2mean_vs_nhits","mean #chi^{2} vs nhits",25,0,25,100,0,10) );
00711   h_chi2meanh.push_back( dbe_->bookProfile("chi2mean","mean #chi^{2} vs #eta",nintEta,minEta,maxEta, 200, 0, 20) );
00712   h_chi2mean_vs_phi.push_back( dbe_->bookProfile("chi2mean_vs_phi","mean of #chi^{2} vs #phi",nintPhi,minPhi,maxPhi, 200, 0, 20) );
00713 
00714   h_hits_eta.push_back( dbe_->bookProfile("hits_eta","mean #hits vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00715   h_PXBhits_eta.push_back( dbe_->bookProfile("PXBhits_eta","mean # PXB hits vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00716   h_PXFhits_eta.push_back( dbe_->bookProfile("PXFhits_eta","mean # PXF hits vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00717   h_TIBhits_eta.push_back( dbe_->bookProfile("TIBhits_eta","mean # TIB hits vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00718   h_TIDhits_eta.push_back( dbe_->bookProfile("TIDhits_eta","mean # TID hits vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00719   h_TOBhits_eta.push_back( dbe_->bookProfile("TOBhits_eta","mean # TOB hits vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00720   h_TEChits_eta.push_back( dbe_->bookProfile("TEChits_eta","mean # TEC hits vs eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00721 
00722   h_LayersWithMeas_eta.push_back(dbe_->bookProfile("LayersWithMeas_eta","mean # LayersWithMeas vs eta",
00723                            nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00724   h_PXLlayersWithMeas_eta.push_back(dbe_->bookProfile("PXLlayersWith2dMeas_eta","mean # PXLlayersWithMeas vs eta",
00725                               nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00726   h_STRIPlayersWithMeas_eta.push_back(dbe_->bookProfile("STRIPlayersWithMeas_eta","mean # STRIPlayersWithMeas vs eta",
00727                             nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00728   h_STRIPlayersWith1dMeas_eta.push_back(dbe_->bookProfile("STRIPlayersWith1dMeas_eta","mean # STRIPlayersWith1dMeas vs eta",
00729                               nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00730   h_STRIPlayersWith2dMeas_eta.push_back(dbe_->bookProfile("STRIPlayersWith2dMeas_eta","mean # STRIPlayersWith2dMeas vs eta",
00731                               nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00732   h_hits_phi.push_back( dbe_->bookProfile("hits_phi","mean #hits vs #phi",nintPhi,minPhi,maxPhi, nintHit,minHit,maxHit) );
00733   h_losthits_eta.push_back( dbe_->bookProfile("losthits_eta","losthits_eta",nintEta,minEta,maxEta,nintHit,minHit,maxHit) );
00734 
00735   h_ptrmsh.push_back( dbe_->book1D("ptres_vs_eta_Sigma","#sigma(#deltap_{t}/p_{t}) vs #eta",nintEta,minEta,maxEta) );
00736   h_ptmeanhPhi.push_back( dbe_->book1D("ptres_vs_phi_Mean","mean of p_{t} resolution vs #phi",nintPhi,minPhi,maxPhi));
00737   h_ptrmshPhi.push_back( dbe_->book1D("ptres_vs_phi_Sigma","#sigma(#deltap_{t}/p_{t}) vs #phi",nintPhi,minPhi,maxPhi) );
00738   h_ptmeanhPt.push_back( dbe_->book1D("ptres_vs_pt_Mean","mean of p_{t} resolution vs p_{t}",nintPt,minPt,maxPt));
00739   h_ptrmshPt.push_back( dbe_->book1D("ptres_vs_pt_Sigma","#sigma(#deltap_{t}/p_{t}) vs pT",nintPt,minPt,maxPt) );
00740   h_cotThetameanh.push_back( dbe_->book1D("cotThetares_vs_eta_Mean","#sigma(cot(#theta)) vs #eta Mean",nintEta,minEta,maxEta) );
00741   h_cotThetarmsh.push_back( dbe_->book1D("cotThetares_vs_eta_Sigma","#sigma(cot(#theta)) vs #eta Sigma",nintEta,minEta,maxEta) );
00742   h_cotThetameanhPt.push_back( dbe_->book1D("cotThetares_vs_pt_Mean","#sigma(cot(#theta)) vs pT Mean",nintPt,minPt,maxPt) );
00743   h_cotThetarmshPt.push_back( dbe_->book1D("cotThetares_vs_pt_Sigma","#sigma(cot(#theta)) vs pT Sigma",nintPt,minPt,maxPt) );
00744   h_phimeanh.push_back(dbe_->book1D("phires_vs_eta_Mean","mean of #phi res vs #eta",nintEta,minEta,maxEta));
00745   h_phirmsh.push_back( dbe_->book1D("phires_vs_eta_Sigma","#sigma(#delta#phi) vs #eta",nintEta,minEta,maxEta) );
00746   h_phimeanhPt.push_back(dbe_->book1D("phires_vs_pt_Mean","mean of #phi res vs pT",nintPt,minPt,maxPt));
00747   h_phirmshPt.push_back( dbe_->book1D("phires_vs_pt_Sigma","#sigma(#delta#phi) vs pT",nintPt,minPt,maxPt) );
00748   h_phimeanhPhi.push_back(dbe_->book1D("phires_vs_phi_Mean","mean of #phi res vs #phi",nintPhi,minPhi,maxPhi));
00749   h_phirmshPhi.push_back( dbe_->book1D("phires_vs_phi_Sigma","#sigma(#delta#phi) vs #phi",nintPhi,minPhi,maxPhi) );
00750   h_dxymeanh.push_back( dbe_->book1D("dxyres_vs_eta_Mean","mean of dxyres vs #eta",nintEta,minEta,maxEta) );
00751   h_dxyrmsh.push_back( dbe_->book1D("dxyres_vs_eta_Sigma","#sigma(#deltadxy) vs #eta",nintEta,minEta,maxEta) );
00752   h_dxymeanhPt.push_back( dbe_->book1D("dxyres_vs_pt_Mean","mean of dxyres vs pT",nintPt,minPt,maxPt) );
00753   h_dxyrmshPt.push_back( dbe_->book1D("dxyres_vs_pt_Sigma","#sigmadxy vs pT",nintPt,minPt,maxPt) );
00754   h_dzmeanh.push_back( dbe_->book1D("dzres_vs_eta_Mean","mean of dzres vs #eta",nintEta,minEta,maxEta) );
00755   h_dzrmsh.push_back( dbe_->book1D("dzres_vs_eta_Sigma","#sigma(#deltadz) vs #eta",nintEta,minEta,maxEta) );
00756   h_dzmeanhPt.push_back( dbe_->book1D("dzres_vs_pt_Mean","mean of dzres vs pT",nintPt,minPt,maxPt) );
00757   h_dzrmshPt.push_back( dbe_->book1D("dzres_vs_pt_Sigma","#sigma(#deltadz vs pT",nintPt,minPt,maxPt) );
00758   h_dxypulletamean.push_back( dbe_->book1D("h_dxypulleta_Mean","mean of dxy pull vs #eta",nintEta,minEta,maxEta) ); 
00759   h_ptpulletamean.push_back( dbe_->book1D("h_ptpulleta_Mean","mean of p_{t} pull vs #eta",nintEta,minEta,maxEta) ); 
00760   h_dzpulletamean.push_back( dbe_->book1D("h_dzpulleta_Mean","mean of dz pull vs #eta",nintEta,minEta,maxEta) ); 
00761   h_phipulletamean.push_back( dbe_->book1D("h_phipulleta_Mean","mean of #phi pull vs #eta",nintEta,minEta,maxEta) ); 
00762   h_thetapulletamean.push_back( dbe_->book1D("h_thetapulleta_Mean","mean of #theta pull vs #eta",nintEta,minEta,maxEta) );
00763   h_dxypulleta.push_back( dbe_->book1D("h_dxypulleta_Sigma","#sigma of dxy pull vs #eta",nintEta,minEta,maxEta) ); 
00764   h_ptpulleta.push_back( dbe_->book1D("h_ptpulleta_Sigma","#sigma of p_{t} pull vs #eta",nintEta,minEta,maxEta) ); 
00765   h_dzpulleta.push_back( dbe_->book1D("h_dzpulleta_Sigma","#sigma of dz pull vs #eta",nintEta,minEta,maxEta) ); 
00766   h_phipulleta.push_back( dbe_->book1D("h_phipulleta_Sigma","#sigma of #phi pull vs #eta",nintEta,minEta,maxEta) ); 
00767   h_thetapulleta.push_back( dbe_->book1D("h_thetapulleta_Sigma","#sigma of #theta pull vs #eta",nintEta,minEta,maxEta) );
00768   h_ptshifteta.push_back( dbe_->book1D("ptres_vs_eta_Mean","<#deltapT/pT>[%] vs #eta",nintEta,minEta,maxEta) ); 
00769   h_ptpullphimean.push_back( dbe_->book1D("h_ptpullphi_Mean","mean of p_{t} pull vs #phi",nintPhi,minPhi,maxPhi) ); 
00770   h_phipullphimean.push_back( dbe_->book1D("h_phipullphi_Mean","mean of #phi pull vs #phi",nintPhi,minPhi,maxPhi) );
00771   h_thetapullphimean.push_back( dbe_->book1D("h_thetapullphi_Mean","mean of #theta pull vs #phi",nintPhi,minPhi,maxPhi) );
00772   h_ptpullphi.push_back( dbe_->book1D("h_ptpullphi_Sigma","#sigma of p_{t} pull vs #phi",nintPhi,minPhi,maxPhi) ); 
00773   h_phipullphi.push_back( dbe_->book1D("h_phipullphi_Sigma","#sigma of #phi pull vs #phi",nintPhi,minPhi,maxPhi) );
00774   h_thetapullphi.push_back( dbe_->book1D("h_thetapullphi_Sigma","#sigma of #theta pull vs #phi",nintPhi,minPhi,maxPhi) );
00775   
00776   if(useLogPt){
00777     BinLogX(h_dzmeanhPt.back()->getTH1F());
00778     BinLogX(h_dzrmshPt.back()->getTH1F());
00779     BinLogX(h_dxymeanhPt.back()->getTH1F());
00780     BinLogX(h_dxyrmshPt.back()->getTH1F());
00781     BinLogX(h_phimeanhPt.back()->getTH1F());
00782     BinLogX(h_phirmshPt.back()->getTH1F());
00783     BinLogX(h_cotThetameanhPt.back()->getTH1F());
00784     BinLogX(h_cotThetarmshPt.back()->getTH1F());
00785     BinLogX(h_ptmeanhPt.back()->getTH1F());
00786     BinLogX(h_ptrmshPt.back()->getTH1F());
00787     BinLogX(h_efficPt.back()->getTH1F());
00788     BinLogX(h_fakeratePt.back()->getTH1F());
00789     BinLogX(h_loopratepT.back()->getTH1F());
00790     BinLogX(h_misidratepT.back()->getTH1F());
00791   }    
00792 }
00793 
00794 void MTVHistoProducerAlgoForTracker::fill_generic_simTrack_histos(int count,
00795                                                                   ParticleBase::Vector momentumTP,
00796                                                                   ParticleBase::Point vertexTP,
00797                                   int bx){
00798   h_ptSIM[count]->Fill(sqrt(momentumTP.perp2()));
00799   h_etaSIM[count]->Fill(momentumTP.eta());
00800   h_vertposSIM[count]->Fill(sqrt(vertexTP.perp2()));
00801   h_bunchxSIM[count]->Fill(bx);
00802 }
00803 
00804 
00805 
00806 // TO BE FIXED USING PLAIN HISTOGRAMS INSTEAD OF RE-IMPLEMENTATION OF HISTOGRAMS (i.d. vectors<int/double>)
00807 void MTVHistoProducerAlgoForTracker::fill_recoAssociated_simTrack_histos(int count,
00808                                                                          const TrackingParticle& tp,
00809                                                                          ParticleBase::Vector momentumTP,
00810                                                                          ParticleBase::Point vertexTP,
00811                                                                          double dxySim, double dzSim, int nSimHits,
00812                                                                          const reco::Track* track,
00813                                      int numVertices, double vertz){
00814   bool isMatched = track;
00815 
00816   if((*TpSelectorForEfficiencyVsEta)(tp)){
00817     //effic vs hits
00818     int nSimHitsInBounds = std::min((int)nSimHits,int(maxHit-1));
00819     totSIM_hit[count][nSimHitsInBounds]++;
00820     if(isMatched) {
00821       totASS_hit[count][nSimHitsInBounds]++;
00822       nrecHit_vs_nsimHit_sim2rec[count]->Fill( track->numberOfValidHits(),nSimHits);
00823     }
00824 
00825     //effic vs eta
00826     for (unsigned int f=0; f<etaintervals[count].size()-1; f++){
00827       if (getEta(momentumTP.eta())>etaintervals[count][f]&&
00828           getEta(momentumTP.eta())<etaintervals[count][f+1]) {
00829         totSIMeta[count][f]++;
00830         if (isMatched) {
00831           totASSeta[count][f]++;
00832         }
00833       }
00834 
00835     } // END for (unsigned int f=0; f<etaintervals[w].size()-1; f++){
00836 
00837     //effic vs num pileup vertices
00838     for (unsigned int f=0; f<vertcountintervals[count].size()-1; f++){
00839       if (numVertices == vertcountintervals[count][f]) {
00840         totSIM_vertcount_entire[count][f]++;
00841         if (isMatched) {
00842           totASS_vertcount_entire[count][f]++;
00843         }
00844       }
00845       if (numVertices == vertcountintervals[count][f] && momentumTP.eta() <= 0.9 && momentumTP.eta() >= -0.9) {
00846         totSIM_vertcount_barrel[count][f]++;
00847         if (isMatched) {
00848           totASS_vertcount_barrel[count][f]++;
00849         }
00850       }
00851       if (numVertices == vertcountintervals[count][f] && momentumTP.eta() > 0.9) {
00852         totSIM_vertcount_fwdpos[count][f]++;
00853         if (isMatched) {
00854           totASS_vertcount_fwdpos[count][f]++;
00855         }
00856       }
00857       if (numVertices == vertcountintervals[count][f] && momentumTP.eta() < -0.9) {
00858         totSIM_vertcount_fwdneg[count][f]++;
00859         if (isMatched) {
00860           totASS_vertcount_fwdneg[count][f]++;
00861         }
00862       }
00863     }
00864 
00865   }
00866 
00867   if((*TpSelectorForEfficiencyVsPhi)(tp)){
00868     for (unsigned int f=0; f<phiintervals[count].size()-1; f++){
00869       if (momentumTP.phi() > phiintervals[count][f]&&
00870           momentumTP.phi() <phiintervals[count][f+1]) {
00871         totSIM_phi[count][f]++;
00872         if (isMatched) {
00873           totASS_phi[count][f]++;
00874         }
00875       }
00876     } // END for (unsigned int f=0; f<phiintervals[count].size()-1; f++){
00877   }
00878         
00879   if((*TpSelectorForEfficiencyVsPt)(tp)){
00880     for (unsigned int f=0; f<pTintervals[count].size()-1; f++){
00881       if (getPt(sqrt(momentumTP.perp2()))>pTintervals[count][f]&&
00882           getPt(sqrt(momentumTP.perp2()))<pTintervals[count][f+1]) {
00883         totSIMpT[count][f]++;
00884         if (isMatched) {
00885           totASSpT[count][f]++;
00886         }
00887       }
00888     } // END for (unsigned int f=0; f<pTintervals[count].size()-1; f++){
00889   }     
00890 
00891   if((*TpSelectorForEfficiencyVsVTXR)(tp)){
00892     for (unsigned int f=0; f<dxyintervals[count].size()-1; f++){
00893       if (dxySim>dxyintervals[count][f]&&
00894           dxySim<dxyintervals[count][f+1]) {
00895         totSIM_dxy[count][f]++;
00896         if (isMatched) {
00897           totASS_dxy[count][f]++;
00898         }
00899       }
00900     } // END for (unsigned int f=0; f<dxyintervals[count].size()-1; f++){
00901 
00902     for (unsigned int f=0; f<vertposintervals[count].size()-1; f++){
00903       if (sqrt(vertexTP.perp2())>vertposintervals[count][f]&&
00904           sqrt(vertexTP.perp2())<vertposintervals[count][f+1]) {
00905         totSIM_vertpos[count][f]++;
00906         if (isMatched) {
00907           totASS_vertpos[count][f]++;
00908         }
00909       }
00910     } // END for (unsigned int f=0; f<vertposintervals[count].size()-1; f++){
00911   }
00912 
00913   if((*TpSelectorForEfficiencyVsVTXZ)(tp)){
00914     for (unsigned int f=0; f<dzintervals[count].size()-1; f++){
00915       if (dzSim>dzintervals[count][f]&&
00916           dzSim<dzintervals[count][f+1]) {
00917         totSIM_dz[count][f]++;
00918         if (isMatched) {
00919           totASS_dz[count][f]++;
00920         }
00921       }
00922     } // END for (unsigned int f=0; f<dzintervals[count].size()-1; f++){
00923 
00924   
00925     for (unsigned int f=0; f<zposintervals[count].size()-1; f++){
00926         if (vertexTP.z()>zposintervals[count][f]&&vertexTP.z()<zposintervals[count][f+1]) {
00927                 totSIM_zpos[count][f]++;
00928                 if (isMatched) totASS_zpos[count][f]++;
00929         }
00930         if (vertz>zposintervals[count][f]&&vertz<zposintervals[count][f+1]) {
00931                 totSIM_vertz_entire[count][f]++;
00932                 if (isMatched) totASS_vertz_entire[count][f]++;
00933         }
00934         if (vertz>zposintervals[count][f]&&vertz<zposintervals[count][f+1] && fabs(momentumTP.eta())<0.9) {
00935                 totSIM_vertz_barrel[count][f]++;
00936                 if (isMatched) totASS_vertz_barrel[count][f]++;
00937         }
00938         if (vertz>zposintervals[count][f]&&vertz<zposintervals[count][f+1] && momentumTP.eta()>0.9) {
00939                 totSIM_vertz_fwdpos[count][f]++;
00940                 if (isMatched) totASS_vertz_fwdpos[count][f]++;
00941         }
00942         if (vertz>zposintervals[count][f]&&vertz<zposintervals[count][f+1] && momentumTP.eta()<-0.9) {
00943                 totSIM_vertz_fwdneg[count][f]++;
00944                 if (isMatched) totASS_vertz_fwdneg[count][f]++;
00945         }
00946     } // END for (unsigned int f=0; f<zposintervals[count].size()-1; f++){
00947   }
00948 
00949   //Special investigations for PU  
00950   if(((*TpSelectorForEfficiencyVsCon)(tp)) && (!((*TpSelectorForEfficiencyVsEta)(tp)))){
00951 
00952    //efficPU vs eta
00953     for (unsigned int f=0; f<etaintervals[count].size()-1; f++){
00954       if (getEta(momentumTP.eta())>etaintervals[count][f]&&
00955           getEta(momentumTP.eta())<etaintervals[count][f+1]) {
00956         totCONeta[count][f]++;
00957       }
00958     } // END for (unsigned int f=0; f<etaintervals[w].size()-1; f++){
00959 
00960     //efficPU vs num pileup vertices
00961     for (unsigned int f=0; f<vertcountintervals[count].size()-1; f++){
00962       if (numVertices == vertcountintervals[count][f]) {
00963         totCONvertcount[count][f]++;
00964       }
00965     } // END for (unsigned int f=0; f<vertcountintervals[count].size()-1; f++){
00966   
00967     for (unsigned int f=0; f<zposintervals[count].size()-1; f++){
00968       if (vertexTP.z()>zposintervals[count][f]&&vertexTP.z()<zposintervals[count][f+1]) {
00969         totCONzpos[count][f]++;
00970       }
00971     } // END for (unsigned int f=0; f<zposintervals[count].size()-1; f++){
00972 
00973   }
00974 
00975 }
00976 
00977 // dE/dx
00978 void MTVHistoProducerAlgoForTracker::fill_dedx_recoTrack_histos(int count, edm::RefToBase<reco::Track>& trackref, std::vector< edm::ValueMap<reco::DeDxData> > v_dEdx) {
00979 //void MTVHistoProducerAlgoForTracker::fill_dedx_recoTrack_histos(reco::TrackRef trackref, std::vector< edm::ValueMap<reco::DeDxData> > v_dEdx) {
00980   double dedx;
00981   int nom;
00982   int sat;
00983   edm::ValueMap<reco::DeDxData> dEdxTrack;
00984   for (unsigned int i=0; i<v_dEdx.size(); i++) {
00985     dEdxTrack = v_dEdx.at(i);
00986     dedx = dEdxTrack[trackref].dEdx(); 
00987     nom  = dEdxTrack[trackref].numberOfMeasurements();
00988     sat  = dEdxTrack[trackref].numberOfSaturatedMeasurements();
00989     if (i==0) {
00990       h_dedx_estim1[count]->Fill(dedx);
00991       h_dedx_nom1[count]->Fill(nom);
00992       h_dedx_sat1[count]->Fill(sat);
00993     } else if (i==1) {
00994       h_dedx_estim2[count]->Fill(dedx);
00995       h_dedx_nom2[count]->Fill(nom);
00996       h_dedx_sat2[count]->Fill(sat);
00997     }
00998   }
00999 }
01000 
01001 
01002 // TO BE FIXED USING PLAIN HISTOGRAMS INSTEAD OF RE-IMPLEMENTATION OF HISTOGRAMS (i.d. vectors<int/double>)
01003 void MTVHistoProducerAlgoForTracker::fill_generic_recoTrack_histos(int count,
01004                                                                    const reco::Track& track,
01005                                                                    math::XYZPoint bsPosition,
01006                                                                    bool isMatched,
01007                                                                    bool isSigMatched,
01008                                                                    bool isChargeMatched,
01009                                    int numAssocRecoTracks,
01010                                    int numVertices,
01011                                    int tpbunchcrossing, int nSimHits,
01012                                    double sharedFraction){
01013 
01014   //Fill track algo histogram
01015 
01016   if (track.algo()>=4 && track.algo()<=14) totREC_algo[count][track.algo()-4]++;  
01017   int sharedHits = sharedFraction *  track.numberOfValidHits();
01018 
01019   //Compute fake rate vs eta
01020   for (unsigned int f=0; f<etaintervals[count].size()-1; f++){
01021     if (getEta(track.momentum().eta())>etaintervals[count][f]&&
01022         getEta(track.momentum().eta())<etaintervals[count][f+1]) {
01023       totRECeta[count][f]++;
01024       if (isMatched) {
01025         totASS2eta[count][f]++;
01026         if (!isChargeMatched) totmisideta[count][f]++;
01027         if (numAssocRecoTracks>1) totloopeta[count][f]++;
01028         if (tpbunchcrossing==0) totASS2_itpu_eta_entire[count][f]++;
01029         if (tpbunchcrossing!=0) totASS2_ootpu_eta_entire[count][f]++;
01030         nrecHit_vs_nsimHit_rec2sim[count]->Fill( track.numberOfValidHits(),nSimHits);
01031         h_assocFraction[count]->Fill( sharedFraction);
01032         h_assocSharedHit[count]->Fill( sharedHits);
01033       }
01034       if (isSigMatched) {
01035         totASS2etaSig[count][f]++;
01036         if (tpbunchcrossing==0) totASS2_itpu_eta_entire_signal[count][f]++;
01037       }
01038     }
01039   } // End for (unsigned int f=0; f<etaintervals[count].size()-1; f++){
01040 
01041   for (unsigned int f=0; f<phiintervals[count].size()-1; f++){
01042     if (track.momentum().phi()>phiintervals[count][f]&&
01043         track.momentum().phi()<phiintervals[count][f+1]) {
01044       totREC_phi[count][f]++; 
01045       if (isMatched) {
01046         totASS2_phi[count][f]++;
01047         if (!isChargeMatched) totmisid_phi[count][f]++;
01048         if (numAssocRecoTracks>1) totloop_phi[count][f]++;
01049       }         
01050     }
01051   } // End for (unsigned int f=0; f<phiintervals[count].size()-1; f++){
01052 
01053         
01054   for (unsigned int f=0; f<pTintervals[count].size()-1; f++){
01055     if (getPt(sqrt(track.momentum().perp2()))>pTintervals[count][f]&&
01056         getPt(sqrt(track.momentum().perp2()))<pTintervals[count][f+1]) {
01057       totRECpT[count][f]++; 
01058       if (isMatched) {
01059         totASS2pT[count][f]++;
01060         if (!isChargeMatched) totmisidpT[count][f]++;
01061         if (numAssocRecoTracks>1) totlooppT[count][f]++;
01062       }       
01063     }
01064   } // End for (unsigned int f=0; f<pTintervals[count].size()-1; f++){
01065   
01066   for (unsigned int f=0; f<dxyintervals[count].size()-1; f++){
01067     if (track.dxy(bsPosition)>dxyintervals[count][f]&&
01068         track.dxy(bsPosition)<dxyintervals[count][f+1]) {
01069       totREC_dxy[count][f]++; 
01070       if (isMatched) {
01071         totASS2_dxy[count][f]++;
01072         if (!isChargeMatched) totmisid_dxy[count][f]++;
01073         if (numAssocRecoTracks>1) totloop_dxy[count][f]++;
01074       }       
01075     }
01076   } // End for (unsigned int f=0; f<dxyintervals[count].size()-1; f++){
01077   
01078   for (unsigned int f=0; f<dzintervals[count].size()-1; f++){
01079     if (track.dz(bsPosition)>dzintervals[count][f]&&
01080         track.dz(bsPosition)<dzintervals[count][f+1]) {
01081       totREC_dz[count][f]++; 
01082       if (isMatched) {
01083         totASS2_dz[count][f]++;
01084         if (!isChargeMatched) totmisid_dz[count][f]++;
01085         if (numAssocRecoTracks>1) totloop_dz[count][f]++;
01086       }       
01087     }
01088   } // End for (unsigned int f=0; f<dzintervals[count].size()-1; f++){
01089 
01090   int tmp = std::min((int)track.found(),int(maxHit-1));
01091   totREC_hit[count][tmp]++;
01092   if (isMatched) {
01093     totASS2_hit[count][tmp]++;
01094     if (!isChargeMatched) totmisid_hit[count][tmp]++;
01095     if (numAssocRecoTracks>1) totloop_hit[count][tmp]++;
01096   }
01097 
01098   for (unsigned int f=0; f<vertcountintervals[count].size()-1; f++){
01099     if (numVertices ==  vertcountintervals[count][f]) {
01100       totREC_vertcount_entire[count][f]++;
01101       totREC_ootpu_entire[count][f]++;
01102       if (isMatched) {
01103         totASS2_vertcount_entire[count][f]++;
01104         if (tpbunchcrossing==0) totASS2_itpu_vertcount_entire[count][f]++;
01105         if (tpbunchcrossing!=0) totASS2_ootpu_entire[count][f]++;
01106       }
01107       if (isSigMatched) {
01108         totASS2_vertcount_entire_signal[count][f]++;
01109         if (tpbunchcrossing==0) totASS2_itpu_vertcount_entire_signal[count][f]++;
01110       }
01111     }
01112     if (numVertices ==  vertcountintervals[count][f] && track.eta() <= 0.9 && track.eta() >= -0.9) {
01113       totREC_vertcount_barrel[count][f]++;
01114       totREC_ootpu_barrel[count][f]++;
01115       if (isMatched) {
01116         totASS2_vertcount_barrel[count][f]++;
01117         if (isMatched && tpbunchcrossing!=0) totASS2_ootpu_barrel[count][f]++;
01118       }
01119     }
01120     if (numVertices ==  vertcountintervals[count][f] && track.eta() > 0.9) {
01121       totREC_vertcount_fwdpos[count][f]++;
01122       totREC_ootpu_fwdpos[count][f]++;
01123       if (isMatched) {
01124         totASS2_vertcount_fwdpos[count][f]++;
01125         if (isMatched && tpbunchcrossing!=0) totASS2_ootpu_fwdpos[count][f]++;
01126       }
01127     }
01128     if (numVertices ==  vertcountintervals[count][f] && track.eta() < -0.9) {
01129       totREC_vertcount_fwdneg[count][f]++;
01130       totREC_ootpu_fwdneg[count][f]++;
01131       if (isMatched) {
01132         totASS2_vertcount_fwdneg[count][f]++;
01133         if (isMatched && tpbunchcrossing!=0) totASS2_ootpu_fwdneg[count][f]++;
01134       }
01135     }
01136 
01137   }
01138 
01139 }
01140 
01141 
01142 void MTVHistoProducerAlgoForTracker::fill_simAssociated_recoTrack_histos(int count,
01143                                                                          const reco::Track& track){
01144     //nchi2 and hits global distributions
01145     h_nchi2[count]->Fill(track.normalizedChi2());
01146     h_nchi2_prob[count]->Fill(TMath::Prob(track.chi2(),(int)track.ndof()));
01147     h_hits[count]->Fill(track.numberOfValidHits());
01148     h_losthits[count]->Fill(track.numberOfLostHits());
01149     chi2_vs_nhits[count]->Fill(track.numberOfValidHits(),track.normalizedChi2());
01150     h_charge[count]->Fill( track.charge() );
01151     
01152     //chi2 and #hit vs eta: fill 2D histos
01153     chi2_vs_eta[count]->Fill(getEta(track.eta()),track.normalizedChi2());
01154     nhits_vs_eta[count]->Fill(getEta(track.eta()),track.numberOfValidHits());
01155     nPXBhits_vs_eta[count]->Fill(getEta(track.eta()),track.hitPattern().numberOfValidPixelBarrelHits());
01156     nPXFhits_vs_eta[count]->Fill(getEta(track.eta()),track.hitPattern().numberOfValidPixelEndcapHits());
01157     nTIBhits_vs_eta[count]->Fill(getEta(track.eta()),track.hitPattern().numberOfValidStripTIBHits());
01158     nTIDhits_vs_eta[count]->Fill(getEta(track.eta()),track.hitPattern().numberOfValidStripTIDHits());
01159     nTOBhits_vs_eta[count]->Fill(getEta(track.eta()),track.hitPattern().numberOfValidStripTOBHits());
01160     nTEChits_vs_eta[count]->Fill(getEta(track.eta()),track.hitPattern().numberOfValidStripTECHits());
01161     nLayersWithMeas_vs_eta[count]->Fill(getEta(track.eta()),track.hitPattern().trackerLayersWithMeasurement());
01162     nPXLlayersWithMeas_vs_eta[count]->Fill(getEta(track.eta()),track.hitPattern().pixelLayersWithMeasurement());
01163     int LayersAll = track.hitPattern().stripLayersWithMeasurement();
01164     int Layers2D = track.hitPattern().numberOfValidStripLayersWithMonoAndStereo(); 
01165     int Layers1D = LayersAll - Layers2D;        
01166     nSTRIPlayersWithMeas_vs_eta[count]->Fill(getEta(track.eta()),LayersAll);
01167     nSTRIPlayersWith1dMeas_vs_eta[count]->Fill(getEta(track.eta()),Layers1D);
01168     nSTRIPlayersWith2dMeas_vs_eta[count]->Fill(getEta(track.eta()),Layers2D);
01169         
01170     nlosthits_vs_eta[count]->Fill(getEta(track.eta()),track.numberOfLostHits());
01171 }
01172 
01173 
01174 void MTVHistoProducerAlgoForTracker::fill_trackBased_histos(int count, int assTracks, int numRecoTracks, int numSimTracks){
01175 
01176         h_tracks[count]->Fill(assTracks);
01177         h_fakes[count]->Fill(numRecoTracks-assTracks);
01178         nrec_vs_nsim[count]->Fill(numRecoTracks,numSimTracks);
01179 
01180 }
01181 
01182 
01183 
01184 void MTVHistoProducerAlgoForTracker::fill_ResoAndPull_recoTrack_histos(int count,
01185                                                                        ParticleBase::Vector momentumTP,
01186                                                                        ParticleBase::Point vertexTP,
01187                                                                        int chargeTP,
01188                                                                        const reco::Track& track,
01189                                                                        math::XYZPoint bsPosition){
01190 
01191   // evaluation of TP parameters
01192   double qoverpSim = chargeTP/sqrt(momentumTP.x()*momentumTP.x()+momentumTP.y()*momentumTP.y()+momentumTP.z()*momentumTP.z());
01193   double lambdaSim = M_PI/2-momentumTP.theta();
01194   double phiSim    = momentumTP.phi();
01195   double dxySim    = (-vertexTP.x()*sin(momentumTP.phi())+vertexTP.y()*cos(momentumTP.phi()));
01196   double dzSim     = vertexTP.z() - (vertexTP.x()*momentumTP.x()+vertexTP.y()*momentumTP.y())/sqrt(momentumTP.perp2()) 
01197     * momentumTP.z()/sqrt(momentumTP.perp2());
01198 
01199           
01200   //  reco::Track::ParameterVector rParameters = track.parameters(); // UNUSED
01201   
01202   double qoverpRec(0);
01203   double qoverpErrorRec(0); 
01204   double ptRec(0);
01205   double ptErrorRec(0);
01206   double lambdaRec(0); 
01207   double lambdaErrorRec(0);
01208   double phiRec(0);
01209   double phiErrorRec(0);
01210 
01211   /* TO BE FIXED LATER  -----------
01212   //loop to decide whether to take gsfTrack (utilisation of mode-function) or common track
01213   const GsfTrack* gsfTrack(0);
01214   if(useGsf){
01215     gsfTrack = dynamic_cast<const GsfTrack*>(&(*track));
01216     if (gsfTrack==0) edm::LogInfo("TrackValidator") << "Trying to access mode for a non-GsfTrack";
01217   }
01218   
01219   if (gsfTrack) {
01220     // get values from mode
01221     getRecoMomentum(*gsfTrack, ptRec, ptErrorRec, qoverpRec, qoverpErrorRec, 
01222                     lambdaRec,lambdaErrorRec, phiRec, phiErrorRec); 
01223   }
01224          
01225   else {
01226     // get values from track (without mode) 
01227     getRecoMomentum(*track, ptRec, ptErrorRec, qoverpRec, qoverpErrorRec, 
01228                     lambdaRec,lambdaErrorRec, phiRec, phiErrorRec); 
01229   }
01230   */
01231   getRecoMomentum(track, ptRec, ptErrorRec, qoverpRec, qoverpErrorRec, 
01232                   lambdaRec,lambdaErrorRec, phiRec, phiErrorRec); 
01233   // -------------
01234 
01235   double ptError = ptErrorRec;  
01236   double ptres=ptRec-sqrt(momentumTP.perp2()); 
01237   double etares=track.eta()-momentumTP.Eta();
01238 
01239 
01240   double dxyRec    = track.dxy(bsPosition);
01241   double dzRec     = track.dz(bsPosition);
01242 
01243   // eta residue; pt, k, theta, phi, dxy, dz pulls
01244   double qoverpPull=(qoverpRec-qoverpSim)/qoverpErrorRec;
01245   double thetaPull=(lambdaRec-lambdaSim)/lambdaErrorRec;
01246   double phiPull=(phiRec-phiSim)/phiErrorRec;
01247   double dxyPull=(dxyRec-dxySim)/track.dxyError();
01248   double dzPull=(dzRec-dzSim)/track.dzError();
01249 
01250   double contrib_Qoverp = ((qoverpRec-qoverpSim)/qoverpErrorRec)*
01251     ((qoverpRec-qoverpSim)/qoverpErrorRec)/5;
01252   double contrib_dxy = ((dxyRec-dxySim)/track.dxyError())*((dxyRec-dxySim)/track.dxyError())/5;
01253   double contrib_dz = ((dzRec-dzSim)/track.dzError())*((dzRec-dzSim)/track.dzError())/5;
01254   double contrib_theta = ((lambdaRec-lambdaSim)/lambdaErrorRec)*
01255     ((lambdaRec-lambdaSim)/lambdaErrorRec)/5;
01256   double contrib_phi = ((phiRec-phiSim)/phiErrorRec)*
01257     ((phiRec-phiSim)/phiErrorRec)/5;
01258 
01259   LogTrace("TrackValidatorTEST") 
01260     //<< "assocChi2=" << tp.begin()->second << "\n"
01261     << "" <<  "\n"
01262     << "ptREC=" << ptRec << "\n" << "etaREC=" << track.eta() << "\n" << "qoverpREC=" << qoverpRec << "\n"
01263     << "dxyREC=" << dxyRec << "\n" << "dzREC=" << dzRec << "\n"
01264     << "thetaREC=" << track.theta() << "\n" << "phiREC=" << phiRec << "\n"
01265     << "" <<  "\n"
01266     << "qoverpError()=" << qoverpErrorRec << "\n" << "dxyError()=" << track.dxyError() << "\n"<< "dzError()=" 
01267     << track.dzError() << "\n"
01268     << "thetaError()=" << lambdaErrorRec << "\n" << "phiError()=" << phiErrorRec << "\n"
01269     << "" <<  "\n"
01270     << "ptSIM=" << sqrt(momentumTP.perp2()) << "\n"<< "etaSIM=" << momentumTP.Eta() << "\n"<< "qoverpSIM=" << qoverpSim << "\n"
01271     << "dxySIM=" << dxySim << "\n"<< "dzSIM=" << dzSim << "\n" << "thetaSIM=" << M_PI/2-lambdaSim << "\n" 
01272     << "phiSIM=" << phiSim << "\n"
01273     << "" << "\n"
01274     << "contrib_Qoverp=" << contrib_Qoverp << "\n"<< "contrib_dxy=" << contrib_dxy << "\n"<< "contrib_dz=" << contrib_dz << "\n"
01275     << "contrib_theta=" << contrib_theta << "\n"<< "contrib_phi=" << contrib_phi << "\n"
01276     << "" << "\n"
01277     <<"chi2PULL="<<contrib_Qoverp+contrib_dxy+contrib_dz+contrib_theta+contrib_phi<<"\n";
01278 
01279   h_pullQoverp[count]->Fill(qoverpPull);
01280   h_pullTheta[count]->Fill(thetaPull);
01281   h_pullPhi[count]->Fill(phiPull);
01282   h_pullDxy[count]->Fill(dxyPull);
01283   h_pullDz[count]->Fill(dzPull);
01284 
01285 
01286   h_pt[count]->Fill(ptres/ptError);
01287   h_eta[count]->Fill(etares);
01288   etares_vs_eta[count]->Fill(getEta(track.eta()),etares);
01289   
01290 
01291         
01292   //resolution of track params: fill 2D histos
01293   dxyres_vs_eta[count]->Fill(getEta(track.eta()),dxyRec-dxySim);
01294   ptres_vs_eta[count]->Fill(getEta(track.eta()),(ptRec-sqrt(momentumTP.perp2()))/ptRec);
01295   dzres_vs_eta[count]->Fill(getEta(track.eta()),dzRec-dzSim);
01296   phires_vs_eta[count]->Fill(getEta(track.eta()),phiRec-phiSim);
01297   cotThetares_vs_eta[count]->Fill(getEta(track.eta()),1/tan(M_PI*0.5-lambdaRec)-1/tan(M_PI*0.5-lambdaSim));         
01298   
01299   //same as before but vs pT
01300   dxyres_vs_pt[count]->Fill(getPt(ptRec),dxyRec-dxySim);
01301   ptres_vs_pt[count]->Fill(getPt(ptRec),(ptRec-sqrt(momentumTP.perp2()))/ptRec);
01302   dzres_vs_pt[count]->Fill(getPt(ptRec),dzRec-dzSim);
01303   phires_vs_pt[count]->Fill(getPt(ptRec),phiRec-phiSim);
01304   cotThetares_vs_pt[count]->Fill(getPt(ptRec),1/tan(M_PI*0.5-lambdaRec)-1/tan(M_PI*0.5-lambdaSim));      
01305         
01306   //pulls of track params vs eta: fill 2D histos
01307   dxypull_vs_eta[count]->Fill(getEta(track.eta()),dxyPull);
01308   ptpull_vs_eta[count]->Fill(getEta(track.eta()),ptres/ptError);
01309   dzpull_vs_eta[count]->Fill(getEta(track.eta()),dzPull);
01310   phipull_vs_eta[count]->Fill(getEta(track.eta()),phiPull);
01311   thetapull_vs_eta[count]->Fill(getEta(track.eta()),thetaPull);
01312         
01313   //plots vs phi
01314   nhits_vs_phi[count]->Fill(phiRec,track.numberOfValidHits());
01315   chi2_vs_phi[count]->Fill(phiRec,track.normalizedChi2());
01316   ptmean_vs_eta_phi[count]->Fill(phiRec,getEta(track.eta()),ptRec);
01317   phimean_vs_eta_phi[count]->Fill(phiRec,getEta(track.eta()),phiRec);
01318   ptres_vs_phi[count]->Fill(phiRec,(ptRec-sqrt(momentumTP.perp2()))/ptRec);
01319   phires_vs_phi[count]->Fill(phiRec,phiRec-phiSim); 
01320   ptpull_vs_phi[count]->Fill(phiRec,ptres/ptError);
01321   phipull_vs_phi[count]->Fill(phiRec,phiPull); 
01322   thetapull_vs_phi[count]->Fill(phiRec,thetaPull); 
01323 
01324 
01325 }
01326 
01327 
01328 
01329 void 
01330 MTVHistoProducerAlgoForTracker::getRecoMomentum (const reco::Track& track, double& pt, double& ptError,
01331                                                  double& qoverp, double& qoverpError, double& lambda,double& lambdaError,
01332                                                  double& phi, double& phiError ) const {
01333   pt = track.pt();
01334   ptError = track.ptError();
01335   qoverp = track.qoverp();
01336   qoverpError = track.qoverpError();
01337   lambda = track.lambda();
01338   lambdaError = track.lambdaError(); 
01339   phi = track.phi(); 
01340   phiError = track.phiError();
01341   //   cout <<"test1" << endl; 
01342   
01343 
01344 
01345 }
01346 
01347 void 
01348 MTVHistoProducerAlgoForTracker::getRecoMomentum (const reco::GsfTrack& gsfTrack, double& pt, double& ptError,
01349                                                  double& qoverp, double& qoverpError, double& lambda,double& lambdaError,
01350                                                  double& phi, double& phiError  ) const {
01351 
01352   pt = gsfTrack.ptMode();
01353   ptError = gsfTrack.ptModeError();
01354   qoverp = gsfTrack.qoverpMode();
01355   qoverpError = gsfTrack.qoverpModeError();
01356   lambda = gsfTrack.lambdaMode();
01357   lambdaError = gsfTrack.lambdaModeError(); 
01358   phi = gsfTrack.phiMode(); 
01359   phiError = gsfTrack.phiModeError();
01360   //   cout <<"test2" << endl;
01361 
01362 }
01363 
01364 double 
01365 MTVHistoProducerAlgoForTracker::getEta(double eta) {
01366   if (useFabsEta) return fabs(eta);
01367   else return eta;
01368 }
01369   
01370 double 
01371 MTVHistoProducerAlgoForTracker::getPt(double pt) {
01372   if (useInvPt && pt!=0) return 1/pt;
01373   else return pt;
01374 }
01375 
01376 
01377 void MTVHistoProducerAlgoForTracker::finalHistoFits(int counter){
01378   //resolution of track params: get sigma from 2D histos
01379   FitSlicesYTool fsyt_dxy(dxyres_vs_eta[counter]);
01380   fsyt_dxy.getFittedSigmaWithError(h_dxyrmsh[counter]);
01381   fsyt_dxy.getFittedMeanWithError(h_dxymeanh[counter]);
01382   FitSlicesYTool fsyt_dxyPt(dxyres_vs_pt[counter]);
01383   fsyt_dxyPt.getFittedSigmaWithError(h_dxyrmshPt[counter]);
01384   fsyt_dxyPt.getFittedMeanWithError(h_dxymeanhPt[counter]);
01385   FitSlicesYTool fsyt_pt(ptres_vs_eta[counter]);
01386   fsyt_pt.getFittedSigmaWithError(h_ptrmsh[counter]);
01387   fsyt_pt.getFittedMeanWithError(h_ptshifteta[counter]);      
01388   FitSlicesYTool fsyt_ptPt(ptres_vs_pt[counter]);
01389   fsyt_ptPt.getFittedSigmaWithError(h_ptrmshPt[counter]);
01390   fsyt_ptPt.getFittedMeanWithError(h_ptmeanhPt[counter]);
01391   FitSlicesYTool fsyt_ptPhi(ptres_vs_phi[counter]); 
01392   fsyt_ptPhi.getFittedSigmaWithError(h_ptrmshPhi[counter]);
01393   fsyt_ptPhi.getFittedMeanWithError(h_ptmeanhPhi[counter]);
01394   FitSlicesYTool fsyt_dz(dzres_vs_eta[counter]);
01395   fsyt_dz.getFittedSigmaWithError(h_dzrmsh[counter]);
01396   fsyt_dz.getFittedMeanWithError(h_dzmeanh[counter]);
01397   FitSlicesYTool fsyt_dzPt(dzres_vs_pt[counter]);
01398   fsyt_dzPt.getFittedSigmaWithError(h_dzrmshPt[counter]);
01399   fsyt_dzPt.getFittedMeanWithError(h_dzmeanhPt[counter]);
01400   FitSlicesYTool fsyt_phi(phires_vs_eta[counter]);
01401   fsyt_phi.getFittedSigmaWithError(h_phirmsh[counter]);
01402   fsyt_phi.getFittedMeanWithError(h_phimeanh[counter]);
01403   FitSlicesYTool fsyt_phiPt(phires_vs_pt[counter]);
01404   fsyt_phiPt.getFittedSigmaWithError(h_phirmshPt[counter]);
01405   fsyt_phiPt.getFittedMeanWithError(h_phimeanhPt[counter]);
01406   FitSlicesYTool fsyt_phiPhi(phires_vs_phi[counter]); 
01407   fsyt_phiPhi.getFittedSigmaWithError(h_phirmshPhi[counter]); 
01408   fsyt_phiPhi.getFittedMeanWithError(h_phimeanhPhi[counter]); 
01409   FitSlicesYTool fsyt_cotTheta(cotThetares_vs_eta[counter]);
01410   fsyt_cotTheta.getFittedSigmaWithError(h_cotThetarmsh[counter]);
01411   fsyt_cotTheta.getFittedMeanWithError(h_cotThetameanh[counter]);
01412   FitSlicesYTool fsyt_cotThetaPt(cotThetares_vs_pt[counter]);
01413   fsyt_cotThetaPt.getFittedSigmaWithError(h_cotThetarmshPt[counter]);
01414   fsyt_cotThetaPt.getFittedMeanWithError(h_cotThetameanhPt[counter]);
01415 
01416   //pulls of track params vs eta: get sigma from 2D histos
01417   FitSlicesYTool fsyt_dxyp(dxypull_vs_eta[counter]);
01418   fsyt_dxyp.getFittedSigmaWithError(h_dxypulleta[counter]);
01419   fsyt_dxyp.getFittedMeanWithError(h_dxypulletamean[counter]);
01420   FitSlicesYTool fsyt_ptp(ptpull_vs_eta[counter]);
01421   fsyt_ptp.getFittedSigmaWithError(h_ptpulleta[counter]);
01422   fsyt_ptp.getFittedMeanWithError(h_ptpulletamean[counter]);
01423   FitSlicesYTool fsyt_dzp(dzpull_vs_eta[counter]);
01424   fsyt_dzp.getFittedSigmaWithError(h_dzpulleta[counter]);
01425   fsyt_dzp.getFittedMeanWithError(h_dzpulletamean[counter]);
01426   FitSlicesYTool fsyt_phip(phipull_vs_eta[counter]);
01427   fsyt_phip.getFittedSigmaWithError(h_phipulleta[counter]);
01428   fsyt_phip.getFittedMeanWithError(h_phipulletamean[counter]);
01429   FitSlicesYTool fsyt_thetap(thetapull_vs_eta[counter]);
01430   fsyt_thetap.getFittedSigmaWithError(h_thetapulleta[counter]);
01431   fsyt_thetap.getFittedMeanWithError(h_thetapulletamean[counter]);
01432   //vs phi
01433   FitSlicesYTool fsyt_ptpPhi(ptpull_vs_phi[counter]);
01434   fsyt_ptpPhi.getFittedSigmaWithError(h_ptpullphi[counter]);
01435   fsyt_ptpPhi.getFittedMeanWithError(h_ptpullphimean[counter]);
01436   FitSlicesYTool fsyt_phipPhi(phipull_vs_phi[counter]);
01437   fsyt_phipPhi.getFittedSigmaWithError(h_phipullphi[counter]);
01438   fsyt_phipPhi.getFittedMeanWithError(h_phipullphimean[counter]);
01439   FitSlicesYTool fsyt_thetapPhi(thetapull_vs_phi[counter]);
01440   fsyt_thetapPhi.getFittedSigmaWithError(h_thetapullphi[counter]);
01441   fsyt_thetapPhi.getFittedMeanWithError(h_thetapullphimean[counter]);
01442   
01443   //effic&fake;
01444   for (unsigned int ite = 0;ite<totASSeta[counter].size();ite++) totFOMT_eta[counter][ite]=totASS2eta[counter][ite]-totASS2etaSig[counter][ite];
01445   for (unsigned int ite = 0;ite<totASS2_vertcount_entire[counter].size();ite++)  
01446     totFOMT_vertcount[counter][ite]=totASS2_vertcount_entire[counter][ite]-totASS2_vertcount_entire_signal[counter][ite];
01447   for (unsigned int ite = 0;ite<totASS2_itpu_eta_entire[counter].size();ite++) totASS2_itpu_eta_entire[counter][ite]-=totASS2_itpu_eta_entire_signal[counter][ite];
01448   for (unsigned int ite = 0;ite<totASS2_itpu_vertcount_entire[counter].size();ite++) totASS2_itpu_vertcount_entire[counter][ite]-=totASS2_itpu_vertcount_entire_signal[counter][ite];
01449 
01450   fillPlotFromVectors(h_effic[counter],totASSeta[counter],totSIMeta[counter],"effic");
01451   fillPlotFromVectors(h_fakerate[counter],totASS2eta[counter],totRECeta[counter],"fakerate");
01452   fillPlotFromVectors(h_looprate[counter],totloopeta[counter],totRECeta[counter],"effic");
01453   fillPlotFromVectors(h_misidrate[counter],totmisideta[counter],totRECeta[counter],"effic");
01454   fillPlotFromVectors(h_efficPt[counter],totASSpT[counter],totSIMpT[counter],"effic");
01455   fillPlotFromVectors(h_fakeratePt[counter],totASS2pT[counter],totRECpT[counter],"fakerate");
01456   fillPlotFromVectors(h_loopratepT[counter],totlooppT[counter],totRECpT[counter],"effic");
01457   fillPlotFromVectors(h_misidratepT[counter],totmisidpT[counter],totRECpT[counter],"effic");
01458   fillPlotFromVectors(h_effic_vs_hit[counter],totASS_hit[counter],totSIM_hit[counter],"effic");
01459   fillPlotFromVectors(h_fake_vs_hit[counter],totASS2_hit[counter],totREC_hit[counter],"fakerate");
01460   fillPlotFromVectors(h_loopratehit[counter],totloop_hit[counter],totREC_hit[counter],"effic");
01461   fillPlotFromVectors(h_misidratehit[counter],totmisid_hit[counter],totREC_hit[counter],"effic");
01462   fillPlotFromVectors(h_effic_vs_phi[counter],totASS_phi[counter],totSIM_phi[counter],"effic");
01463   fillPlotFromVectors(h_fake_vs_phi[counter],totASS2_phi[counter],totREC_phi[counter],"fakerate");
01464   fillPlotFromVectors(h_loopratephi[counter],totloop_phi[counter],totREC_phi[counter],"effic");
01465   fillPlotFromVectors(h_misidratephi[counter],totmisid_phi[counter],totREC_phi[counter],"effic");
01466   fillPlotFromVectors(h_effic_vs_dxy[counter],totASS_dxy[counter],totSIM_dxy[counter],"effic");
01467   fillPlotFromVectors(h_fake_vs_dxy[counter],totASS2_dxy[counter],totREC_dxy[counter],"fakerate");
01468   fillPlotFromVectors(h_loopratedxy[counter],totloop_dxy[counter],totREC_dxy[counter],"effic");
01469   fillPlotFromVectors(h_misidratedxy[counter],totmisid_dxy[counter],totREC_dxy[counter],"effic");
01470   fillPlotFromVectors(h_effic_vs_dz[counter],totASS_dz[counter],totSIM_dz[counter],"effic");
01471   fillPlotFromVectors(h_fake_vs_dz[counter],totASS2_dz[counter],totREC_dz[counter],"fakerate");
01472   fillPlotFromVectors(h_loopratedz[counter],totloop_dz[counter],totREC_dz[counter],"effic");
01473   fillPlotFromVectors(h_misidratedz[counter],totmisid_dz[counter],totREC_dz[counter],"effic");
01474   fillPlotFromVectors(h_effic_vs_vertpos[counter],totASS_vertpos[counter],totSIM_vertpos[counter],"effic");
01475   fillPlotFromVectors(h_effic_vs_zpos[counter],totASS_zpos[counter],totSIM_zpos[counter],"effic");
01476   fillPlotFromVectors(h_effic_vertcount_entire[counter],totASS_vertcount_entire[counter],totSIM_vertcount_entire[counter],"effic");
01477   fillPlotFromVectors(h_effic_vertcount_barrel[counter],totASS_vertcount_barrel[counter],totSIM_vertcount_barrel[counter],"effic");
01478   fillPlotFromVectors(h_effic_vertcount_fwdpos[counter],totASS_vertcount_fwdpos[counter],totSIM_vertcount_fwdpos[counter],"effic");
01479   fillPlotFromVectors(h_effic_vertcount_fwdneg[counter],totASS_vertcount_fwdneg[counter],totSIM_vertcount_fwdneg[counter],"effic");
01480   fillPlotFromVectors(h_fakerate_vertcount_entire[counter],totASS2_vertcount_entire[counter],totREC_vertcount_entire[counter],"fakerate");
01481   fillPlotFromVectors(h_fakerate_vertcount_barrel[counter],totASS2_vertcount_barrel[counter],totREC_vertcount_barrel[counter],"fakerate");
01482   fillPlotFromVectors(h_fakerate_vertcount_fwdpos[counter],totASS2_vertcount_fwdpos[counter],totREC_vertcount_fwdpos[counter],"fakerate");
01483   fillPlotFromVectors(h_fakerate_vertcount_fwdneg[counter],totASS2_vertcount_fwdneg[counter],totREC_vertcount_fwdneg[counter],"fakerate");
01484   fillPlotFromVectors(h_effic_vertz_entire[counter],totASS_vertz_entire[counter],totSIM_vertz_entire[counter],"effic");
01485   fillPlotFromVectors(h_effic_vertz_barrel[counter],totASS_vertz_barrel[counter],totSIM_vertz_barrel[counter],"effic");
01486   fillPlotFromVectors(h_effic_vertz_fwdpos[counter],totASS_vertz_fwdpos[counter],totSIM_vertz_fwdpos[counter],"effic");
01487   fillPlotFromVectors(h_effic_vertz_fwdneg[counter],totASS_vertz_fwdneg[counter],totSIM_vertz_fwdneg[counter],"effic");
01488   fillPlotFromVectors(h_fakerate_ootpu_entire[counter],totASS2_ootpu_entire[counter],totREC_ootpu_entire[counter],"effic");
01489   fillPlotFromVectors(h_fakerate_ootpu_barrel[counter],totASS2_ootpu_barrel[counter],totREC_ootpu_barrel[counter],"effic");
01490   fillPlotFromVectors(h_fakerate_ootpu_fwdpos[counter],totASS2_ootpu_fwdpos[counter],totREC_ootpu_fwdpos[counter],"effic");
01491   fillPlotFromVectors(h_fakerate_ootpu_fwdneg[counter],totASS2_ootpu_fwdneg[counter],totREC_ootpu_fwdneg[counter],"effic");
01492 
01493   fillPlotFromVectors(h_fomt_eta[counter],totFOMT_eta[counter],totASS2etaSig[counter],"pileup");
01494   fillPlotFromVectors(h_fomt_vertcount[counter],totFOMT_vertcount[counter],totASS2_vertcount_entire_signal[counter],"pileup");
01495   fillPlotFromVectors(h_fomt_sig_eta[counter],totASS2etaSig[counter],totRECeta[counter],"fakerate");
01496   fillPlotFromVectors(h_fomt_sig_vertcount[counter],totASS2_vertcount_entire_signal[counter],totREC_vertcount_entire[counter],"fakerate");
01497   fillPlotFromVectors(h_fomt_itpu_eta[counter],totASS2_itpu_eta_entire[counter],totRECeta[counter],"effic");
01498   fillPlotFromVectors(h_fomt_itpu_vertcount[counter],totASS2_itpu_vertcount_entire[counter],totREC_vertcount_entire[counter],"effic");
01499   fillPlotFromVectors(h_fomt_ootpu_eta[counter],totASS2_ootpu_eta_entire[counter],totRECeta[counter],"effic");
01500   fillPlotFromVectors(h_fomt_ootpu_vertcount[counter],totASS2_ootpu_entire[counter],totREC_ootpu_entire[counter],"effic");
01501 
01502   fillPlotFromVectors(h_effic_PU_eta[counter],totASSeta[counter],totCONeta[counter],"pileup");
01503   fillPlotFromVectors(h_effic_PU_vertcount[counter],totASS_vertcount_entire[counter],totCONvertcount[counter],"pileup");
01504   fillPlotFromVectors(h_effic_PU_zpos[counter],totASS_zpos[counter],totCONzpos[counter],"pileup");
01505 
01506 }
01507 
01508 void MTVHistoProducerAlgoForTracker::fillProfileHistosFromVectors(int counter){
01509   //chi2 and #hit vs eta: get mean from 2D histos
01510   doProfileX(chi2_vs_eta[counter],h_chi2meanh[counter]);
01511   doProfileX(nhits_vs_eta[counter],h_hits_eta[counter]);
01512   doProfileX(nPXBhits_vs_eta[counter],h_PXBhits_eta[counter]);
01513   doProfileX(nPXFhits_vs_eta[counter],h_PXFhits_eta[counter]);
01514   doProfileX(nTIBhits_vs_eta[counter],h_TIBhits_eta[counter]);
01515   doProfileX(nTIDhits_vs_eta[counter],h_TIDhits_eta[counter]);
01516   doProfileX(nTOBhits_vs_eta[counter],h_TOBhits_eta[counter]);
01517   doProfileX(nTEChits_vs_eta[counter],h_TEChits_eta[counter]);
01518 
01519   doProfileX(nLayersWithMeas_vs_eta[counter],h_LayersWithMeas_eta[counter]);
01520   doProfileX(nPXLlayersWithMeas_vs_eta[counter],h_PXLlayersWithMeas_eta[counter]);    
01521   doProfileX(nSTRIPlayersWithMeas_vs_eta[counter],h_STRIPlayersWithMeas_eta[counter]);    
01522   doProfileX(nSTRIPlayersWith1dMeas_vs_eta[counter],h_STRIPlayersWith1dMeas_eta[counter]);
01523   doProfileX(nSTRIPlayersWith2dMeas_vs_eta[counter],h_STRIPlayersWith2dMeas_eta[counter]);
01524 
01525 
01526 
01527   doProfileX(nlosthits_vs_eta[counter],h_losthits_eta[counter]);
01528   //vs phi
01529   doProfileX(chi2_vs_nhits[counter],h_chi2meanhitsh[counter]); 
01530   //      doProfileX(ptres_vs_eta[counter],h_ptresmean_vs_eta[counter]);
01531   //      doProfileX(phires_vs_eta[counter],h_phiresmean_vs_eta[counter]);
01532   doProfileX(chi2_vs_phi[counter],h_chi2mean_vs_phi[counter]);
01533   doProfileX(nhits_vs_phi[counter],h_hits_phi[counter]);
01534   //       doProfileX(ptres_vs_phi[counter],h_ptresmean_vs_phi[counter]);
01535   //       doProfileX(phires_vs_phi[counter],h_phiresmean_vs_phi[counter]);
01536 }
01537 
01538 void MTVHistoProducerAlgoForTracker::fillHistosFromVectors(int counter){
01539   fillPlotFromVector(h_algo[counter],totREC_algo[counter]);
01540 
01541   fillPlotFromVector(h_recoeta[counter],totRECeta[counter]);
01542   fillPlotFromVector(h_simuleta[counter],totSIMeta[counter]);
01543   fillPlotFromVector(h_assoceta[counter],totASSeta[counter]);
01544   fillPlotFromVector(h_assoc2eta[counter],totASS2eta[counter]);
01545   fillPlotFromVector(h_loopereta[counter],totloopeta[counter]);
01546   fillPlotFromVector(h_misideta[counter],totmisideta[counter]);
01547   
01548   fillPlotFromVector(h_recopT[counter],totRECpT[counter]);
01549   fillPlotFromVector(h_simulpT[counter],totSIMpT[counter]);
01550   fillPlotFromVector(h_assocpT[counter],totASSpT[counter]);
01551   fillPlotFromVector(h_assoc2pT[counter],totASS2pT[counter]);
01552   fillPlotFromVector(h_looperpT[counter],totlooppT[counter]);
01553   fillPlotFromVector(h_misidpT[counter],totmisidpT[counter]);
01554   
01555   fillPlotFromVector(h_recohit[counter],totREC_hit[counter]);
01556   fillPlotFromVector(h_simulhit[counter],totSIM_hit[counter]);
01557   fillPlotFromVector(h_assochit[counter],totASS_hit[counter]);
01558   fillPlotFromVector(h_assoc2hit[counter],totASS2_hit[counter]);
01559   fillPlotFromVector(h_looperhit[counter],totloop_hit[counter]);
01560   fillPlotFromVector(h_misidhit[counter],totmisid_hit[counter]);
01561   
01562   fillPlotFromVector(h_recophi[counter],totREC_phi[counter]);
01563   fillPlotFromVector(h_simulphi[counter],totSIM_phi[counter]);
01564   fillPlotFromVector(h_assocphi[counter],totASS_phi[counter]);
01565   fillPlotFromVector(h_assoc2phi[counter],totASS2_phi[counter]);
01566   fillPlotFromVector(h_looperphi[counter],totloop_phi[counter]);
01567   fillPlotFromVector(h_misidphi[counter],totmisid_phi[counter]);
01568   
01569   fillPlotFromVector(h_recodxy[counter],totREC_dxy[counter]);
01570   fillPlotFromVector(h_simuldxy[counter],totSIM_dxy[counter]);
01571   fillPlotFromVector(h_assocdxy[counter],totASS_dxy[counter]);
01572   fillPlotFromVector(h_assoc2dxy[counter],totASS2_dxy[counter]);
01573   fillPlotFromVector(h_looperdxy[counter],totloop_dxy[counter]);
01574   fillPlotFromVector(h_misiddxy[counter],totmisid_dxy[counter]);
01575 
01576   fillPlotFromVector(h_recodz[counter],totREC_dz[counter]);
01577   fillPlotFromVector(h_simuldz[counter],totSIM_dz[counter]);
01578   fillPlotFromVector(h_assocdz[counter],totASS_dz[counter]);
01579   fillPlotFromVector(h_assoc2dz[counter],totASS2_dz[counter]);
01580   fillPlotFromVector(h_looperdz[counter],totloop_dz[counter]);
01581   fillPlotFromVector(h_misiddz[counter],totmisid_dz[counter]);
01582   
01583   fillPlotFromVector(h_simulvertpos[counter],totSIM_vertpos[counter]);
01584   fillPlotFromVector(h_assocvertpos[counter],totASS_vertpos[counter]);
01585   
01586   fillPlotFromVector(h_simulzpos[counter],totSIM_zpos[counter]);
01587   fillPlotFromVector(h_assoczpos[counter],totASS_zpos[counter]);  
01588 
01589   fillPlotFromVector(h_reco_vertcount_entire[counter],totREC_vertcount_entire[counter]);
01590   fillPlotFromVector(h_simul_vertcount_entire[counter],totSIM_vertcount_entire[counter]);
01591   fillPlotFromVector(h_assoc_vertcount_entire[counter],totASS_vertcount_entire[counter]);
01592   fillPlotFromVector(h_assoc2_vertcount_entire[counter],totASS2_vertcount_entire[counter]);
01593   
01594   fillPlotFromVector(h_reco_vertcount_barrel[counter],totREC_vertcount_barrel[counter]);
01595   fillPlotFromVector(h_simul_vertcount_barrel[counter],totSIM_vertcount_barrel[counter]);
01596   fillPlotFromVector(h_assoc_vertcount_barrel[counter],totASS_vertcount_barrel[counter]);
01597   fillPlotFromVector(h_assoc2_vertcount_barrel[counter],totASS2_vertcount_barrel[counter]);
01598   
01599   fillPlotFromVector(h_reco_vertcount_fwdpos[counter],totREC_vertcount_fwdpos[counter]);
01600   fillPlotFromVector(h_simul_vertcount_fwdpos[counter],totSIM_vertcount_fwdpos[counter]);
01601   fillPlotFromVector(h_assoc_vertcount_fwdpos[counter],totASS_vertcount_fwdpos[counter]);
01602   fillPlotFromVector(h_assoc2_vertcount_fwdpos[counter],totASS2_vertcount_fwdpos[counter]);
01603   
01604   fillPlotFromVector(h_reco_vertcount_fwdneg[counter],totREC_vertcount_fwdneg[counter]);
01605   fillPlotFromVector(h_simul_vertcount_fwdneg[counter],totSIM_vertcount_fwdneg[counter]);
01606   fillPlotFromVector(h_assoc_vertcount_fwdneg[counter],totASS_vertcount_fwdneg[counter]);
01607   fillPlotFromVector(h_assoc2_vertcount_fwdneg[counter],totASS2_vertcount_fwdneg[counter]);
01608   
01609   fillPlotFromVector(h_simul_vertz_entire[counter],totSIM_vertz_entire[counter]);
01610   fillPlotFromVector(h_assoc_vertz_entire[counter],totASS_vertz_entire[counter]);
01611   
01612   fillPlotFromVector(h_simul_vertz_barrel[counter],totSIM_vertz_barrel[counter]);
01613   fillPlotFromVector(h_assoc_vertz_barrel[counter],totASS_vertz_barrel[counter]);
01614   
01615   fillPlotFromVector(h_simul_vertz_fwdpos[counter],totSIM_vertz_fwdpos[counter]);
01616   fillPlotFromVector(h_assoc_vertz_fwdpos[counter],totASS_vertz_fwdpos[counter]);
01617   
01618   fillPlotFromVector(h_simul_vertz_fwdneg[counter],totSIM_vertz_fwdneg[counter]);
01619   fillPlotFromVector(h_assoc_vertz_fwdneg[counter],totASS_vertz_fwdneg[counter]);
01620   
01621   fillPlotFromVector(h_reco_ootpu_entire[counter],totREC_ootpu_entire[counter]);
01622   fillPlotFromVector(h_assoc2_ootpu_entire[counter],totASS2_ootpu_entire[counter]);
01623   
01624   fillPlotFromVector(h_reco_ootpu_barrel[counter],totREC_ootpu_barrel[counter]);
01625   fillPlotFromVector(h_assoc2_ootpu_barrel[counter],totASS2_ootpu_barrel[counter]);
01626   
01627   fillPlotFromVector(h_reco_ootpu_fwdpos[counter],totREC_ootpu_fwdpos[counter]);
01628   fillPlotFromVector(h_assoc2_ootpu_fwdpos[counter],totASS2_ootpu_fwdpos[counter]);
01629   
01630   fillPlotFromVector(h_reco_ootpu_fwdneg[counter],totREC_ootpu_fwdneg[counter]);
01631   fillPlotFromVector(h_assoc2_ootpu_fwdneg[counter],totASS2_ootpu_fwdneg[counter]);
01632 
01633   fillPlotFromVector(h_con_eta[counter],totCONeta[counter]);
01634   fillPlotFromVector(h_con_vertcount[counter],totCONvertcount[counter]);
01635   fillPlotFromVector(h_con_zpos[counter],totCONzpos[counter]);
01636   
01637 }