CMS 3D CMS Logo

TrackerSeedValidator.cc

Go to the documentation of this file.
00001 #include "Validation/RecoTrack/interface/TrackerSeedValidator.h"
00002 #include "Validation/Tools/interface/FitSlicesYTool.h"
00003 
00004 #include "FWCore/Framework/interface/MakerMacros.h"
00005 
00006 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00007 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00008 #include "SimDataFormats/Vertex/interface/SimVertexContainer.h"
00009 #include "SimTracker/TrackAssociation/interface/TrackAssociatorByChi2.h"
00010 #include "SimTracker/TrackAssociation/interface/TrackAssociatorByHits.h"
00011 #include "SimTracker/TrackerHitAssociation/interface/TrackerHitAssociator.h"
00012 #include "SimTracker/Records/interface/TrackAssociatorRecord.h"
00013 
00014 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00015 
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00018 #include "TrackingTools/PatternTools/interface/TrajectoryStateClosestToBeamLineBuilder.h"
00019 #include "TrackingTools/TrajectoryState/interface/PerigeeConversions.h"
00020 #include <TF1.h>
00021 #include "TrackingTools/Records/interface/TransientRecHitRecord.h" 
00022 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00023 
00024 using namespace std;
00025 using namespace edm;
00026 
00027 void TrackerSeedValidator::beginRun(edm::Run const&, edm::EventSetup const& setup) {
00028   setup.get<IdealMagneticFieldRecord>().get(theMF);  
00029   setup.get<TransientRecHitRecord>().get(builderName,theTTRHBuilder);
00030 
00031   int j=0;
00032   for (unsigned int ww=0;ww<associators.size();ww++){
00033     for (unsigned int www=0;www<label.size();www++){
00034 
00035       dbe_->cd();
00036       InputTag algo = label[www];
00037       string dirName="RecoTrackV/Seed/";
00038       if (algo.process()!="")
00039         dirName+=algo.process()+"_";
00040       if(algo.label()!="")
00041         dirName+=algo.label()+"_";
00042       if(algo.instance()!="")
00043         dirName+=algo.instance()+"_";      
00044       if (dirName.find("Seeds")<dirName.length()){
00045         dirName.replace(dirName.find("Seeds"),6,"");
00046       }
00047       string assoc= associators[ww];
00048       if (assoc.find("Track")<assoc.length()){
00049         assoc.replace(assoc.find("Track"),5,"");
00050       }
00051       dirName+=assoc;
00052       dbe_->setCurrentFolder(dirName.c_str());
00053 
00054       setUpVectors();
00055 
00056       dbe_->goUp();
00057       string subDirName = dirName + "/simulation";
00058       dbe_->setCurrentFolder(subDirName.c_str());
00059       h_ptSIM.push_back( dbe_->book1D("ptSIM", "generated p_{t}", 5500, 0, 110 ) );
00060       h_etaSIM.push_back( dbe_->book1D("etaSIM", "generated pseudorapidity", 500, -2.5, 2.5 ) );
00061       h_tracksSIM.push_back( dbe_->book1D("tracksSIM","number of simluated tracks",100,-0.5,99.5) );
00062       h_vertposSIM.push_back( dbe_->book1D("vertposSIM","Transverse position of sim vertices",1000,-0.5,10000.5) );
00063       
00064       dbe_->cd();
00065       dbe_->setCurrentFolder(dirName.c_str());
00066       h_tracks.push_back( dbe_->book1D("seeds","number of reconstructed seeds",20,-0.5,19.5) );
00067       h_fakes.push_back( dbe_->book1D("fakes","number of fake reco seeds",20,-0.5,19.5) );
00068       h_charge.push_back( dbe_->book1D("charge","charge",3,-1.5,1.5) );
00069       h_hits.push_back( dbe_->book1D("hits", "number of hits per seed", 30, -0.5, 29.5 ) );
00070 
00071       h_effic.push_back( dbe_->book1D("effic","efficiency vs #eta",nint,min,max) );
00072       h_efficPt.push_back( dbe_->book1D("efficPt","efficiency vs pT",nintpT,minpT,maxpT) );
00073       h_fakerate.push_back( dbe_->book1D("fakerate","fake rate vs #eta",nint,min,max) );
00074       h_fakeratePt.push_back( dbe_->book1D("fakeratePt","fake rate vs pT",nintpT,minpT,maxpT) );
00075       h_effic_vs_hit.push_back( dbe_->book1D("effic_vs_hit","effic vs hit",nintHit,minHit,maxHit) );
00076       h_fake_vs_hit.push_back( dbe_->book1D("fakerate_vs_hit","fake rate vs hit",nintHit,minHit,maxHit) );
00077 
00078       h_recoeta.push_back( dbe_->book1D("num_reco_eta","N of reco seed vs eta",nint,min,max) );
00079       h_assoceta.push_back( dbe_->book1D("num_assoc(simToReco)_eta","N of associated seeds (simToReco) vs eta",nint,min,max) );
00080       h_assoc2eta.push_back( dbe_->book1D("num_assoc(recoToSim)_eta","N of associated (recoToSim) seeds vs eta",nint,min,max) );
00081       h_simuleta.push_back( dbe_->book1D("num_simul_eta","N of simulated tracks vs eta",nint,min,max) );
00082       h_recopT.push_back( dbe_->book1D("num_reco_pT","N of reco seed vs pT",nintpT,minpT,maxpT) );
00083       h_assocpT.push_back( dbe_->book1D("num_assoc(simToReco)_pT","N of associated seeds (simToReco) vs pT",nintpT,minpT,maxpT) );
00084       h_assoc2pT.push_back( dbe_->book1D("num_assoc(recoToSim)_pT","N of associated (recoToSim) seeds vs pT",nintpT,minpT,maxpT) );
00085       h_simulpT.push_back( dbe_->book1D("num_simul_pT","N of simulated tracks vs pT",nintpT,minpT,maxpT) );
00086       
00087       h_eta.push_back( dbe_->book1D("eta", "pseudorapidity residue", 1000, -0.1, 0.1 ) );
00088       h_pt.push_back( dbe_->book1D("pullPt", "pull of p_{t}", 100, -10, 10 ) );
00089       h_pullTheta.push_back( dbe_->book1D("pullTheta","pull of #theta parameter",250,-25,25) );
00090       h_pullPhi.push_back( dbe_->book1D("pullPhi","pull of #phi parameter",250,-25,25) );
00091       h_pullDxy.push_back( dbe_->book1D("pullDxy","pull of dxy parameter",250,-25,25) );
00092       h_pullDz.push_back( dbe_->book1D("pullDz","pull of dz parameter",250,-25,25) );
00093       h_pullQoverp.push_back( dbe_->book1D("pullQoverp","pull of qoverp parameter",250,-25,25) );
00094       
00095       if (associators[ww]=="TrackAssociatorByHits"){
00096         h_assocFraction.push_back( dbe_->book1D("assocFraction","fraction of shared hits",200,0,2) );
00097         h_assocSharedHit.push_back(dbe_->book1D("assocSharedHit","number of shared hits",20,0,20));
00098       }
00099 
00100       //chi2_vs_nhits.push_back( dbe_->book2D("chi2_vs_nhits","#chi^{2} vs nhits",25,0,25,100,0,10) );
00101       nrec_vs_nsim.push_back( dbe_->book2D("nrec_vs_nsim","nrec vs nsim",20,-0.5,19.5,20,-0.5,19.5) );
00102 
00103       nhits_vs_eta.push_back( dbe_->book2D("nhits_vs_eta","nhits vs eta",nint,min,max,25,0,25) );
00104       h_hits_eta.push_back( dbe_->book1D("hits_eta","mean #hits vs eta",nint,min,max) );
00105 
00106       j++;
00107     }
00108   }
00109   edm::ESHandle<TrackAssociatorBase> theAssociator;
00110   for (unsigned int w=0;w<associators.size();w++) {
00111     setup.get<TrackAssociatorRecord>().get(associators[w],theAssociator);
00112     associator.push_back( theAssociator.product() );
00113   }
00114 }
00115 
00116 void TrackerSeedValidator::analyze(const edm::Event& event, const edm::EventSetup& setup){
00117 
00118   edm::LogInfo("TrackValidator") << "\n====================================================" << "\n"
00119                                  << "Analyzing new event" << "\n"
00120                                  << "====================================================\n" << "\n";
00121   
00122   edm::Handle<TrackingParticleCollection>  TPCollectionHeff ;
00123   event.getByLabel(label_tp_effic,TPCollectionHeff);
00124   const TrackingParticleCollection tPCeff = *(TPCollectionHeff.product());
00125   
00126   edm::Handle<TrackingParticleCollection>  TPCollectionHfake ;
00127   event.getByLabel(label_tp_fake,TPCollectionHfake);
00128   const TrackingParticleCollection tPCfake = *(TPCollectionHfake.product());
00129 
00130   if (tPCeff.size()==0) {edm::LogInfo("TrackValidator") << "TP Collection for efficiency studies has size = 0! Skipping Event." ; return;}
00131   if (tPCfake.size()==0) {edm::LogInfo("TrackValidator") << "TP Collection for fake rate studies has size = 0! Skipping Event." ; return;}
00132 
00133   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00134   event.getByLabel(bsSrc,recoBeamSpotHandle);
00135   reco::BeamSpot bs = *recoBeamSpotHandle;      
00136   
00137   int w=0;
00138   for (unsigned int ww=0;ww<associators.size();ww++){
00139     for (unsigned int www=0;www<label.size();www++){
00140       edm::LogVerbatim("TrackValidator") << "Analyzing " 
00141                                          << label[www].process()<<":"
00142                                          << label[www].label()<<":"
00143                                          << label[www].instance()<<" with "
00144                                          << associators[ww].c_str() <<"\n";
00145       //
00146       //get collections from the event
00147       //
00148       edm::Handle<edm::View<TrajectorySeed> > seedCollection;
00149       event.getByLabel(label[www], seedCollection);
00150       if (seedCollection->size()==0) {
00151         edm::LogInfo("TrackValidator") << "SeedCollection size = 0!" ; 
00152         continue;
00153       }
00154  
00155       //associate seeds
00156       LogTrace("TrackValidator") << "Calling associateRecoToSim method" << "\n";
00157       reco::RecoToSimCollectionSeed recSimColl=associator[ww]->associateRecoToSim(seedCollection,
00158                                                                                   TPCollectionHfake,
00159                                                                                   &event);
00160       LogTrace("TrackValidator") << "Calling associateSimToReco method" << "\n";
00161       reco::SimToRecoCollectionSeed simRecColl=associator[ww]->associateSimToReco(seedCollection,
00162                                                                                   TPCollectionHeff, 
00163                                                                                   &event);
00164       
00165       //
00166       //fill simulation histograms
00167       //compute number of seeds per eta interval
00168       //
00169       edm::LogVerbatim("TrackValidator") << "\n# of TrackingParticles: " << tPCeff.size() << "\n";
00170       int ats = 0;
00171       int st=0;
00172       for (TrackingParticleCollection::size_type i=0; i<tPCeff.size(); i++){
00173         TrackingParticleRef tp(TPCollectionHeff, i);
00174         if (tp->charge()==0) continue;
00175         st++;
00176         h_ptSIM[w]->Fill(sqrt(tp->momentum().perp2()));
00177         h_etaSIM[w]->Fill(tp->momentum().eta());
00178         h_vertposSIM[w]->Fill(sqrt(tp->vertex().perp2()));
00179 
00180         std::vector<std::pair<edm::RefToBase<TrajectorySeed>, double> > rt;
00181         if(simRecColl.find(tp) != simRecColl.end()){
00182           rt = simRecColl[tp];
00183           if (rt.size()!=0) {
00184             ats++;
00185             edm::LogVerbatim("TrackValidator") << "TrackingParticle #" << st 
00186                                                << " with pt=" << sqrt(tp->momentum().perp2()) 
00187                                                << " associated with quality:" << rt.begin()->second <<"\n";
00188           }
00189         }else{
00190           edm::LogVerbatim("TrackValidator") << "TrackingParticle #" << st
00191                                              << " with pt=" << sqrt(tp->momentum().perp2())
00192                                              << " NOT associated to any TrajectorySeed" << "\n";
00193         }
00194 
00195         for (unsigned int f=0; f<etaintervals[w].size()-1; f++){
00196           if (getEta(tp->momentum().eta())>etaintervals[w][f]&&
00197               getEta(tp->momentum().eta())<etaintervals[w][f+1]) {
00198             totSIMeta[w][f]++;
00199             if (rt.size()!=0) {
00200               totASSeta[w][f]++;
00201             }
00202           }
00203         } // END for (unsigned int f=0; f<etaintervals[w].size()-1; f++){
00204         
00205         for (unsigned int f=0; f<pTintervals[w].size()-1; f++){
00206           if (getPt(sqrt(tp->momentum().perp2()))>pTintervals[w][f]&&
00207               getPt(sqrt(tp->momentum().perp2()))<pTintervals[w][f+1]) {
00208             totSIMpT[w][f]++;
00209             if (rt.size()!=0) {
00210               totASSpT[w][f]++;
00211             }
00212           }
00213         } // END for (unsigned int f=0; f<pTintervals[w].size()-1; f++){
00214 
00215         totSIM_hit[w][tp->matchedHit()]++;//FIXME
00216         if (rt.size()!=0) totASS_hit[w][tp->matchedHit()]++;
00217       }
00218       if (st!=0) h_tracksSIM[w]->Fill(st);
00219       
00220 
00221       //
00222       //fill reconstructed seed histograms
00223       // 
00224       edm::LogVerbatim("TrackValidator") << "\n# of TrajectorySeeds with "
00225                                          << label[www].process()<<":"
00226                                          << label[www].label()<<":"
00227                                          << label[www].instance()
00228                                          << ": " << seedCollection->size() << "\n";
00229       int at=0;
00230       int rT=0;
00231       TrajectoryStateTransform tsTransform;
00232       TrajectoryStateClosestToBeamLineBuilder tscblBuilder;
00233       PerigeeConversions tspConverter;
00234       for(TrajectorySeedCollection::size_type i=0; i<seedCollection->size(); ++i){
00235         edm::RefToBase<TrajectorySeed> seed(seedCollection, i);
00236         rT++;
00237 
00238         //get parameters and errors from the seed state
00239         TransientTrackingRecHit::RecHitPointer recHit = theTTRHBuilder->build(&*(seed->recHits().second-1));
00240         TrajectoryStateOnSurface state = tsTransform.transientState( seed->startingState(), recHit->surface(), theMF.product());
00241 
00242         TrajectoryStateClosestToBeamLine tsAtClosestApproachSeed = tscblBuilder(*state.freeState(),bs);//as in TrackProducerAlgorithm
00243         GlobalPoint vSeed1 = tsAtClosestApproachSeed.trackStateAtPCA().position();
00244         GlobalVector pSeed = tsAtClosestApproachSeed.trackStateAtPCA().momentum();
00245         GlobalPoint vSeed(vSeed1.x()-bs.x0(),vSeed1.y()-bs.y0(),vSeed1.z()-bs.z0());
00246 
00247         double etaSeed = state.globalMomentum().eta();
00248         double ptSeed = sqrt(state.globalMomentum().perp2());
00249         double pmSeed = sqrt(state.globalMomentum().mag2());
00250         double pzSeed = state.globalMomentum().z();
00251         double numberOfHitsSeed = seed->recHits().second-seed->recHits().first;
00252         double qoverpSeed = tsAtClosestApproachSeed.trackStateAtPCA().charge()/pSeed.mag();
00253         double thetaSeed  = pSeed.theta();
00254         double lambdaSeed = M_PI/2-thetaSeed;
00255         double phiSeed    = pSeed.phi();
00256         double dxySeed    = (-vSeed.x()*sin(pSeed.phi())+vSeed.y()*cos(pSeed.phi()));
00257         double dzSeed     = vSeed.z() - (vSeed.x()*pSeed.x()+vSeed.y()*pSeed.y())/pSeed.perp() * pSeed.z()/pSeed.perp();
00258 
00259         PerigeeTrajectoryError seedPerigeeErrors = tspConverter.ftsToPerigeeError(tsAtClosestApproachSeed.trackStateAtPCA());
00260         double qoverpErrorSeed = tsAtClosestApproachSeed.trackStateAtPCA().curvilinearError().matrix().At(0,0);
00261         double lambdaErrorSeed = seedPerigeeErrors.thetaError();//=theta error
00262         double phiErrorSeed = seedPerigeeErrors.phiError();
00263         double dxyErrorSeed = seedPerigeeErrors.transverseImpactParameterError();
00264         double dzErrorSeed = seedPerigeeErrors.longitudinalImpactParameterError();
00265         double ptErrorSeed = (state.charge()!=0) ?  sqrt(
00266                ptSeed*ptSeed*pmSeed*pmSeed/state.charge()/state.charge() * tsAtClosestApproachSeed.trackStateAtPCA().curvilinearError().matrix().At(0,0)
00267                + 2*ptSeed*pmSeed/state.charge()*pzSeed * tsAtClosestApproachSeed.trackStateAtPCA().curvilinearError().matrix().At(0,1)
00268                + pzSeed*pzSeed * tsAtClosestApproachSeed.trackStateAtPCA().curvilinearError().matrix().At(1,1) ) : 1.e6;
00269         
00270         std::vector<std::pair<TrackingParticleRef, double> > tp;
00271         if(recSimColl.find(seed) != recSimColl.end()) {
00272           tp = recSimColl[seed];
00273           if (tp.size()!=0) {
00274             at++;
00275             edm::LogVerbatim("TrackValidator") << "TrajectorySeed #" << rT << " associated with quality:" << tp.begin()->second <<"\n";
00276           }
00277         } else {
00278           edm::LogVerbatim("TrackValidator") << "TrajectorySeed #" << rT << " NOT associated to any TrackingParticle" << "\n";            
00279         }
00280         
00281         //Compute fake rate vs eta
00282         for (unsigned int f=0; f<etaintervals[w].size()-1; f++){
00283           if (getEta(etaSeed)>etaintervals[w][f]&&
00284               getEta(etaSeed)<etaintervals[w][f+1]) {
00285             totRECeta[w][f]++; 
00286             if (tp.size()!=0) {
00287               totASS2eta[w][f]++;
00288             }           
00289           }
00290         }
00291         
00292         for (unsigned int f=0; f<pTintervals[w].size()-1; f++){
00293           if (getPt(ptSeed)>pTintervals[w][f]&&
00294               getPt(ptSeed)<pTintervals[w][f+1]) {
00295             totRECpT[w][f]++; 
00296             if (tp.size()!=0) {
00297               totASS2pT[w][f]++;
00298             }         
00299           }
00300         }
00301         totREC_hit[w][seed->nHits()]++;
00302         if (tp.size()!=0) totASS2_hit[w][seed->nHits()]++;
00303         
00304         //Fill other histos
00305         try{
00306           if (tp.size()==0) continue;
00307         
00308           TrackingParticleRef tpr = tp.begin()->first;
00309           const SimTrack * assocTrack = &(*tpr->g4Track_begin());
00310         
00311           if (associators[ww]=="TrackAssociatorByHits"){
00312             double fraction = tp.begin()->second;
00313             h_assocFraction[www]->Fill(fraction);
00314             h_assocSharedHit[www]->Fill(fraction*numberOfHitsSeed);
00315           }
00316     
00317           h_hits[w]->Fill(numberOfHitsSeed);
00318           h_charge[w]->Fill( state.charge() );
00319         
00320           //compute tracking particle parameters at point of closest approach to the beamline
00321           edm::ESHandle<MagneticField> theMF;
00322           setup.get<IdealMagneticFieldRecord>().get(theMF);
00323           FreeTrajectoryState 
00324             ftsAtProduction(GlobalPoint(tpr->vertex().x(),tpr->vertex().y(),tpr->vertex().z()),
00325                             GlobalVector(assocTrack->momentum().x(),assocTrack->momentum().y(),assocTrack->momentum().z()),
00326                             TrackCharge(tpr->charge()),
00327                             theMF.product());
00328           TrajectoryStateClosestToBeamLine tsAtClosestApproach = tscblBuilder(ftsAtProduction,bs);//as in TrackProducerAlgorithm
00329           GlobalPoint v1 = tsAtClosestApproach.trackStateAtPCA().position();
00330           GlobalVector p = tsAtClosestApproach.trackStateAtPCA().momentum();
00331           GlobalPoint v(v1.x()-bs.x0(),v1.y()-bs.y0(),v1.z()-bs.z0());
00332 
00333           double qoverpSim = tsAtClosestApproach.trackStateAtPCA().charge()/p.mag();
00334           double lambdaSim = M_PI/2-p.theta();
00335           double phiSim    = p.phi();
00336           double dxySim    = (-v.x()*sin(p.phi())+v.y()*cos(p.phi()));
00337           double dzSim     = v.z() - (v.x()*p.x()+v.y()*p.y())/p.perp() * p.z()/p.perp();
00338 
00339           // eta residue; pt, k, theta, phi, dxy, dz pulls
00340           double qoverpPull=(qoverpSeed-qoverpSim)/qoverpErrorSeed;
00341           double thetaPull=(lambdaSeed-lambdaSim)/lambdaErrorSeed;
00342           double phiPull=(phiSeed-phiSim)/phiErrorSeed;
00343           double dxyPull=(dxySeed-dxySim)/dxyErrorSeed;
00344           double dzPull=(dzSeed-dzSim)/dzErrorSeed;
00345 
00346           double contrib_Qoverp = ((qoverpSeed-qoverpSim)/qoverpErrorSeed)*
00347             ((qoverpSeed-qoverpSim)/qoverpErrorSeed)/5;
00348           double contrib_dxy = ((dxySeed-dxySim)/dxyErrorSeed)*((dxySeed-dxySim)/dxyErrorSeed)/5;
00349           double contrib_dz = ((dzSeed-dzSim)/dzErrorSeed)*((dzSeed-dzSim)/dzErrorSeed)/5;
00350           double contrib_theta = ((lambdaSeed-lambdaSim)/lambdaErrorSeed)*
00351             ((lambdaSeed-lambdaSim)/lambdaErrorSeed)/5;
00352           double contrib_phi = ((phiSeed-phiSim)/phiErrorSeed)*
00353             ((phiSeed-phiSim)/phiErrorSeed)/5;
00354           LogTrace("TrackValidatorTEST") << "assocChi2=" << tp.begin()->second << "\n"
00355                                          << "" <<  "\n"
00356                                          << "ptREC=" << ptSeed << "\n"
00357                                          << "etaREC=" << etaSeed << "\n"
00358                                          << "qoverpREC=" << qoverpSeed << "\n"
00359                                          << "dxyREC=" << dxySeed << "\n"
00360                                          << "dzREC=" << dzSeed << "\n"
00361                                          << "thetaREC=" << thetaSeed << "\n"
00362                                          << "phiREC=" << phiSeed << "\n"
00363                                          << "" <<  "\n"
00364                                          << "qoverpError()=" << qoverpErrorSeed << "\n"
00365                                          << "dxyError()=" << dxyErrorSeed << "\n"
00366                                          << "dzError()=" << dzErrorSeed << "\n"
00367                                          << "thetaError()=" << lambdaErrorSeed << "\n"
00368                                          << "phiError()=" << phiErrorSeed << "\n"
00369                                          << "" <<  "\n"
00370                                          << "ptSIM=" << sqrt(assocTrack->momentum().perp2()) << "\n"
00371                                          << "etaSIM=" << assocTrack->momentum().Eta() << "\n"    
00372                                          << "qoverpSIM=" << qoverpSim << "\n"
00373                                          << "dxySIM=" << dxySim << "\n"
00374                                          << "dzSIM=" << dzSim << "\n"
00375                                          << "thetaSIM=" << M_PI/2-lambdaSim << "\n"
00376                                          << "phiSIM=" << phiSim << "\n"
00377                                          << "" << "\n"
00378                                          << "contrib_Qoverp=" << contrib_Qoverp << "\n"
00379                                          << "contrib_dxy=" << contrib_dxy << "\n"
00380                                          << "contrib_dz=" << contrib_dz << "\n"
00381                                          << "contrib_theta=" << contrib_theta << "\n"
00382                                          << "contrib_phi=" << contrib_phi << "\n"
00383                                          << "" << "\n"
00384                                          <<"chi2PULL="<<contrib_Qoverp+contrib_dxy+contrib_dz+contrib_theta+contrib_phi<<"\n";
00385           
00386           h_pullQoverp[w]->Fill(qoverpPull);
00387           h_pullTheta[w]->Fill(thetaPull);
00388           h_pullPhi[w]->Fill(phiPull);
00389           h_pullDxy[w]->Fill(dxyPull);
00390           h_pullDz[w]->Fill(dzPull);
00391 
00392           double ptres=ptSeed-sqrt(assocTrack->momentum().perp2()); 
00393           //double etares=etaSeed-assocTrack->momentum().pseudoRapidity();
00394           double etares=etaSeed-assocTrack->momentum().Eta();
00395 
00396           h_pt[w]->Fill(ptres/ptErrorSeed);
00397           h_eta[w]->Fill(etares);
00398 
00399           //#hit vs eta: fill 2D histos
00400           nhits_vs_eta[w]->Fill(getEta(etaSeed),numberOfHitsSeed);
00401 
00402         } catch (cms::Exception e){
00403           LogTrace("TrackValidator") << "exception found: " << e.what() << "\n";
00404         }
00405       }
00406       if (at!=0) h_tracks[w]->Fill(at);
00407       h_fakes[w]->Fill(rT-at);
00408       edm::LogVerbatim("TrackValidator") << "Total Simulated: " << st << "\n"
00409                                          << "Total Associated (simToReco): " << ats << "\n"
00410                                          << "Total Reconstructed: " << rT << "\n"
00411                                          << "Total Associated (recoToSim): " << at << "\n"
00412                                          << "Total Fakes: " << rT-at << "\n";
00413       nrec_vs_nsim[w]->Fill(rT,st);
00414       w++;
00415     }
00416   }
00417 }
00418 
00419 void TrackerSeedValidator::endRun(edm::Run const&, edm::EventSetup const&) {
00420   LogTrace("TrackValidator") << "TrackerSeedValidator::endRun()";
00421   int w=0;
00422   for (unsigned int ww=0;ww<associators.size();ww++){
00423     for (unsigned int www=0;www<label.size();www++){
00424 
00425       doProfileX(nhits_vs_eta[w],h_hits_eta[w]);    
00426 
00427       //effic&fake
00428       fillPlotFromVectors(h_effic[w],totASSeta[w],totSIMeta[w],"effic");
00429       fillPlotFromVectors(h_fakerate[w],totASS2eta[w],totRECeta[w],"fakerate");
00430       fillPlotFromVectors(h_efficPt[w],totASSpT[w],totSIMpT[w],"effic");
00431       fillPlotFromVectors(h_fakeratePt[w],totASS2pT[w],totRECpT[w],"fakerate");
00432       fillPlotFromVectors(h_effic_vs_hit[w],totASS_hit[w],totSIM_hit[w],"effic");
00433       fillPlotFromVectors(h_fake_vs_hit[w],totASS2_hit[w],totREC_hit[w],"fakerate");
00434 
00435       fillPlotFromVector(h_recoeta[w],totRECeta[w]);
00436       fillPlotFromVector(h_simuleta[w],totSIMeta[w]);
00437       fillPlotFromVector(h_assoceta[w],totASSeta[w]);
00438       fillPlotFromVector(h_assoc2eta[w],totASS2eta[w]);
00439 
00440       fillPlotFromVector(h_recopT[w],totRECpT[w]);
00441       fillPlotFromVector(h_simulpT[w],totSIMpT[w]);
00442       fillPlotFromVector(h_assocpT[w],totASSpT[w]);
00443       fillPlotFromVector(h_assoc2pT[w],totASS2pT[w]);
00444 
00445       w++;
00446     }
00447   }
00448   if ( out.size() != 0 && dbe_ ) dbe_->save(out);
00449 }
00450 
00451 
00452 
00453 

Generated on Tue Jun 9 17:49:39 2009 for CMSSW by  doxygen 1.5.4