CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/Validation/RecoTrack/plugins/TrackerSeedValidator.cc

Go to the documentation of this file.
00001 #include "Validation/RecoTrack/interface/TrackerSeedValidator.h"
00002 #include "DQMServices/ClientConfig/interface/FitSlicesYTool.h"
00003 
00004 #include "FWCore/Framework/interface/MakerMacros.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 
00007 #include "DataFormats/TrackReco/interface/Track.h"
00008 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00009 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00010 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00011 #include "SimDataFormats/Vertex/interface/SimVertexContainer.h"
00012 #include "SimTracker/TrackAssociation/interface/TrackAssociatorByChi2.h"
00013 #include "SimTracker/TrackAssociation/interface/QuickTrackAssociatorByHits.h"
00014 #include "SimTracker/TrackerHitAssociation/interface/TrackerHitAssociator.h"
00015 #include "SimTracker/Records/interface/TrackAssociatorRecord.h"
00016 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00017 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00018 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
00019 #include "TrackingTools/TrajectoryState/interface/PerigeeConversions.h"
00020 #include "TrackingTools/Records/interface/TransientRecHitRecord.h" 
00021 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00022 #include "Validation/RecoTrack/interface/MTVHistoProducerAlgoFactory.h"
00023 #include "SimDataFormats/TrackingAnalysis/interface/TrackingVertex.h"
00024 #include "SimDataFormats/TrackingAnalysis/interface/TrackingVertexContainer.h"
00025 #include "SimDataFormats/PileupSummaryInfo/interface/PileupSummaryInfo.h"
00026 #include "SimDataFormats/EncodedEventId/interface/EncodedEventId.h"
00027 #include "SimTracker/TrackAssociation/plugins/ParametersDefinerForTPESProducer.h"
00028 
00029 #include <TF1.h>
00030 
00031 using namespace std;
00032 using namespace edm;
00033 
00034 typedef edm::Ref<edm::HepMCProduct, HepMC::GenParticle > GenParticleRef;
00035 
00036 TrackerSeedValidator::TrackerSeedValidator(const edm::ParameterSet& pset):MultiTrackValidatorBase(pset){
00037   //theExtractor = IsoDepositExtractorFactory::get()->create( extractorName, extractorPSet);
00038 
00039   ParameterSet psetForHistoProducerAlgo = pset.getParameter<ParameterSet>("histoProducerAlgoBlock");
00040   string histoProducerAlgoName = psetForHistoProducerAlgo.getParameter<string>("ComponentName");
00041   histoProducerAlgo_ = MTVHistoProducerAlgoFactory::get()->create(histoProducerAlgoName ,psetForHistoProducerAlgo);
00042   histoProducerAlgo_->setDQMStore(dbe_);
00043 
00044   dirName_ = pset.getParameter<std::string>("dirName");
00045 
00046   tpSelector = TrackingParticleSelector(pset.getParameter<double>("ptMinTP"),
00047                                         pset.getParameter<double>("minRapidityTP"),
00048                                         pset.getParameter<double>("maxRapidityTP"),
00049                                         pset.getParameter<double>("tipTP"),
00050                                         pset.getParameter<double>("lipTP"),
00051                                         pset.getParameter<int>("minHitTP"),
00052                                         pset.getParameter<bool>("signalOnlyTP"),
00053                                         pset.getParameter<bool>("chargedOnlyTP"),
00054                                         pset.getParameter<bool>("stableOnlyTP"),
00055                                         pset.getParameter<std::vector<int> >("pdgIdTP"));
00056 
00057   runStandalone = pset.getParameter<bool>("runStandalone");
00058 
00059   builderName = pset.getParameter<std::string>("TTRHBuilder");
00060 }
00061 
00062 TrackerSeedValidator::~TrackerSeedValidator(){delete histoProducerAlgo_;}
00063 
00064 void TrackerSeedValidator::beginRun(edm::Run const&, edm::EventSetup const& setup) {
00065   setup.get<IdealMagneticFieldRecord>().get(theMF);  
00066   setup.get<TransientRecHitRecord>().get(builderName,theTTRHBuilder);
00067 
00068   for (unsigned int ww=0;ww<associators.size();ww++){
00069     for (unsigned int www=0;www<label.size();www++){
00070 
00071       dbe_->cd();
00072       InputTag algo = label[www];
00073       string dirName=dirName_;
00074       if (algo.process()!="")
00075         dirName+=algo.process()+"_";
00076       if(algo.label()!="")
00077         dirName+=algo.label()+"_";
00078       if(algo.instance()!="")
00079         dirName+=algo.instance()+"_";      
00080       //      if (dirName.find("Seeds")<dirName.length()){
00081       //        dirName.replace(dirName.find("Seeds"),6,"");
00082       //      }
00083       string assoc= associators[ww];
00084       if (assoc.find("Track")<assoc.length()){
00085         assoc.replace(assoc.find("Track"),5,"");
00086       }
00087       dirName+=assoc;
00088       std::replace(dirName.begin(), dirName.end(), ':', '_');
00089 
00090       dbe_->setCurrentFolder(dirName.c_str());
00091 
00092       
00093       // vector of vector initialization
00094       histoProducerAlgo_->initialize(); //TO BE FIXED. I'D LIKE TO AVOID THIS CALL
00095 
00096       dbe_->goUp(); //Is this really necessary ???
00097       string subDirName = dirName + "/simulation";
00098       dbe_->setCurrentFolder(subDirName.c_str());
00099 
00100       //Booking histograms concerning with simulated tracks
00101       histoProducerAlgo_->bookSimHistos();
00102 
00103       dbe_->cd();
00104       dbe_->setCurrentFolder(dirName.c_str());
00105 
00106       //Booking histograms concerning with reconstructed tracks
00107       histoProducerAlgo_->bookRecoHistos();
00108       if (runStandalone) histoProducerAlgo_->bookRecoHistosForStandaloneRunning();
00109     }//end loop www
00110     edm::ESHandle<TrackAssociatorBase> theAssociator;
00111     for (unsigned int w=0;w<associators.size();w++) {
00112       setup.get<TrackAssociatorRecord>().get(associators[w],theAssociator);
00113       associator.push_back( theAssociator.product() );
00114     }//end loop w
00115   }// end loop ww
00116 }
00117 
00118 void TrackerSeedValidator::analyze(const edm::Event& event, const edm::EventSetup& setup){
00119 
00120   edm::LogInfo("TrackValidator") << "\n====================================================" << "\n"
00121                                  << "Analyzing new event" << "\n"
00122                                  << "====================================================\n" << "\n";
00123   
00124   edm::ESHandle<ParametersDefinerForTP> parametersDefinerTP; 
00125   setup.get<TrackAssociatorRecord>().get(parametersDefiner,parametersDefinerTP);    
00126 
00127   edm::Handle<TrackingParticleCollection>  TPCollectionHeff ;
00128   event.getByLabel(label_tp_effic,TPCollectionHeff);
00129   const TrackingParticleCollection tPCeff = *(TPCollectionHeff.product());
00130   
00131   edm::Handle<TrackingParticleCollection>  TPCollectionHfake ;
00132   event.getByLabel(label_tp_fake,TPCollectionHfake);
00133   const TrackingParticleCollection tPCfake = *(TPCollectionHfake.product());
00134 
00135   if (tPCeff.size()==0) {edm::LogInfo("TrackValidator") << "TP Collection for efficiency studies has size = 0! Skipping Event." ; return;}
00136   if (tPCfake.size()==0) {edm::LogInfo("TrackValidator") << "TP Collection for fake rate studies has size = 0! Skipping Event." ; return;}
00137 
00138   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00139   event.getByLabel(bsSrc,recoBeamSpotHandle);
00140   reco::BeamSpot bs = *recoBeamSpotHandle;      
00141     
00142   edm::Handle< vector<PileupSummaryInfo> > puinfoH;
00143   event.getByLabel(label_pileupinfo,puinfoH);
00144   PileupSummaryInfo puinfo;      
00145   
00146   for (unsigned int puinfo_ite=0;puinfo_ite<(*puinfoH).size();++puinfo_ite){ 
00147     if ((*puinfoH)[puinfo_ite].getBunchCrossing()==0){
00148       puinfo=(*puinfoH)[puinfo_ite];
00149       break;
00150     }
00151   }
00152 
00153   edm::Handle<TrackingVertexCollection> tvH;
00154   event.getByLabel(label_tv,tvH);
00155   TrackingVertexCollection tv = *tvH;      
00156 
00157   int w=0;
00158   for (unsigned int ww=0;ww<associators.size();ww++){
00159     for (unsigned int www=0;www<label.size();www++){
00160       edm::LogVerbatim("TrackValidator") << "Analyzing " 
00161                                          << label[www].process()<<":"
00162                                          << label[www].label()<<":"
00163                                          << label[www].instance()<<" with "
00164                                          << associators[ww].c_str() <<"\n";
00165       //
00166       //get collections from the event
00167       //
00168       edm::Handle<edm::View<TrajectorySeed> > seedCollection;
00169       event.getByLabel(label[www], seedCollection);
00170       if (seedCollection->size()==0) {
00171         edm::LogInfo("TrackValidator") << "SeedCollection size = 0!" ; 
00172         continue;
00173       }
00174  
00175       //associate seeds
00176       LogTrace("TrackValidator") << "Calling associateRecoToSim method" << "\n";
00177       reco::RecoToSimCollectionSeed recSimColl=associator[ww]->associateRecoToSim(seedCollection,
00178                                                                                   TPCollectionHfake,
00179                                                                                   &event);
00180       LogTrace("TrackValidator") << "Calling associateSimToReco method" << "\n";
00181       reco::SimToRecoCollectionSeed simRecColl=associator[ww]->associateSimToReco(seedCollection,
00182                                                                                   TPCollectionHeff, 
00183                                                                                   &event);
00184       
00185       //
00186       //fill simulation histograms
00187       //compute number of seeds per eta interval
00188       //
00189       edm::LogVerbatim("TrackValidator") << "\n# of TrackingParticles: " << tPCeff.size() << "\n";
00190       int ats(0);         //This counter counts the number of simTracks that are "associated" to recoTracks
00191       int st(0);          //This counter counts the number of simulated tracks passing the MTV selection (i.e. tpSelector(tp) )
00192       unsigned sts(0);   //This counter counts the number of simTracks surviving the bunchcrossing cut 
00193       unsigned asts(0);  //This counter counts the number of simTracks that are "associated" to recoTracks surviving the bunchcrossing cut
00194       for (TrackingParticleCollection::size_type i=0; i<tPCeff.size(); i++){
00195         TrackingParticleRef tp(TPCollectionHeff, i);
00196 
00197         if (tp->charge()==0) continue;
00198 
00199         if(! tpSelector(*tp)) continue;
00200 
00201         ParticleBase::Vector momentumTP = tp->momentum();
00202         ParticleBase::Point vertexTP = tp->vertex();
00203         //Calcualte the impact parameters w.r.t. PCA
00204         ParticleBase::Vector momentum = parametersDefinerTP->momentum(event,setup,*tp);
00205         ParticleBase::Point vertex = parametersDefinerTP->vertex(event,setup,*tp);
00206         double dxySim = (-vertex.x()*sin(momentum.phi())+vertex.y()*cos(momentum.phi()));
00207         double dzSim = vertex.z() - (vertex.x()*momentum.x()+vertex.y()*momentum.y())/sqrt(momentum.perp2()) 
00208           * momentum.z()/sqrt(momentum.perp2());
00209         
00210         st++;
00211 
00212         histoProducerAlgo_->fill_generic_simTrack_histos(w,momentumTP,vertexTP, tp->eventId().bunchCrossing());
00213 
00214         const TrajectorySeed* matchedSeedPointer=0;
00215         std::vector<std::pair<edm::RefToBase<TrajectorySeed>, double> > rt;
00216         if(simRecColl.find(tp) != simRecColl.end()){
00217           rt = simRecColl[tp];
00218           if (rt.size()!=0) {
00219             ats++;
00220             matchedSeedPointer = rt.begin()->first.get();
00221             edm::LogVerbatim("TrackValidator") << "TrackingParticle #" << st 
00222                                                << " with pt=" << sqrt(tp->momentum().perp2()) 
00223                                                << " associated with quality:" << rt.begin()->second <<"\n";
00224           }
00225         }else{
00226           edm::LogVerbatim("TrackValidator") << "TrackingParticle #" << st
00227                                              << " with pt=" << sqrt(tp->momentum().perp2())
00228                                              << " NOT associated to any TrajectorySeed" << "\n";
00229         }
00230 
00231         std::vector<PSimHit> simhits=tp->trackPSimHit(DetId::Tracker);
00232         int nSimHits = simhits.end()-simhits.begin();
00233 
00234         double vtx_z_PU = tp->vertex().z();
00235         for (size_t j = 0; j < tv.size(); j++) {
00236             if (tp->eventId().event() == tv[j].eventId().event()) {
00237                 vtx_z_PU = tv[j].position().z();
00238                 break;
00239             }
00240         }
00241 
00242 
00243         //fixme convert seed into track
00244         reco::Track* matchedTrackPointer = 0;
00245         if (matchedSeedPointer) {
00246           TSCBLBuilderNoMaterial tscblBuilder;
00247           TransientTrackingRecHit::RecHitPointer recHit = theTTRHBuilder->build(&*(matchedSeedPointer->recHits().second-1));
00248           TrajectoryStateOnSurface state = trajectoryStateTransform::transientState( matchedSeedPointer->startingState(), recHit->surface(), theMF.product());
00249           TrajectoryStateClosestToBeamLine tsAtClosestApproachSeed = tscblBuilder(*state.freeState(),bs);//as in TrackProducerAlgorithm
00250           if(!(tsAtClosestApproachSeed.isValid())){
00251             edm::LogVerbatim("SeedValidator")<<"TrajectoryStateClosestToBeamLine not valid";
00252             continue;
00253           }
00254           const reco::TrackBase::Point vSeed1(tsAtClosestApproachSeed.trackStateAtPCA().position().x(),
00255                                               tsAtClosestApproachSeed.trackStateAtPCA().position().y(),
00256                                               tsAtClosestApproachSeed.trackStateAtPCA().position().z());
00257           const reco::TrackBase::Vector pSeed(tsAtClosestApproachSeed.trackStateAtPCA().momentum().x(),
00258                                               tsAtClosestApproachSeed.trackStateAtPCA().momentum().y(),
00259                                               tsAtClosestApproachSeed.trackStateAtPCA().momentum().z());
00260           //GlobalPoint vSeed(vSeed1.x()-bs.x0(),vSeed1.y()-bs.y0(),vSeed1.z()-bs.z0());
00261           PerigeeTrajectoryError seedPerigeeErrors = PerigeeConversions::ftsToPerigeeError(tsAtClosestApproachSeed.trackStateAtPCA());
00262           matchedTrackPointer = new reco::Track(0.,0., vSeed1, pSeed, 1, seedPerigeeErrors.covarianceMatrix());
00263           matchedTrackPointer->setHitPattern(matchedSeedPointer->recHits().first,matchedSeedPointer->recHits().second);
00264         }
00265 
00266         histoProducerAlgo_->fill_recoAssociated_simTrack_histos(w,*tp,tp->momentum(),tp->vertex(),dxySim,dzSim,nSimHits,
00267                                                                 matchedTrackPointer,puinfo.getPU_NumInteractions(), vtx_z_PU);
00268 
00269         sts++;
00270         if (matchedTrackPointer) asts++;
00271 
00272       } // End  for (TrackingParticleCollection::size_type i=0; i<tPCeff.size(); i++){
00273 
00274       //
00275       //fill reconstructed seed histograms
00276       // 
00277       edm::LogVerbatim("TrackValidator") << "\n# of TrajectorySeeds with "
00278                                          << label[www].process()<<":"
00279                                          << label[www].label()<<":"
00280                                          << label[www].instance()
00281                                          << ": " << seedCollection->size() << "\n";
00282       int sat(0); //This counter counts the number of recoTracks that are associated to SimTracks from Signal only
00283       int at(0); //This counter counts the number of recoTracks that are associated to SimTracks
00284       int rT(0); //This counter counts the number of recoTracks in general
00285       
00286       TSCBLBuilderNoMaterial tscblBuilder;
00287       for(TrajectorySeedCollection::size_type i=0; i<seedCollection->size(); ++i){
00288         edm::RefToBase<TrajectorySeed> seed(seedCollection, i);
00289         rT++;
00290 
00291         //get parameters and errors from the seed state
00292         TransientTrackingRecHit::RecHitPointer recHit = theTTRHBuilder->build(&*(seed->recHits().second-1));
00293         TrajectoryStateOnSurface state = trajectoryStateTransform::transientState( seed->startingState(), recHit->surface(), theMF.product());
00294         TrajectoryStateClosestToBeamLine tsAtClosestApproachSeed = tscblBuilder(*state.freeState(),bs);//as in TrackProducerAlgorithm
00295         if(!(tsAtClosestApproachSeed.isValid())){
00296           edm::LogVerbatim("SeedValidator")<<"TrajectoryStateClosestToBeamLine not valid";
00297           continue;
00298             }
00299         const reco::TrackBase::Point vSeed1(tsAtClosestApproachSeed.trackStateAtPCA().position().x(),
00300                                             tsAtClosestApproachSeed.trackStateAtPCA().position().y(),
00301                                             tsAtClosestApproachSeed.trackStateAtPCA().position().z());
00302         const reco::TrackBase::Vector pSeed(tsAtClosestApproachSeed.trackStateAtPCA().momentum().x(),
00303                                             tsAtClosestApproachSeed.trackStateAtPCA().momentum().y(),
00304                                             tsAtClosestApproachSeed.trackStateAtPCA().momentum().z());
00305         //GlobalPoint vSeed(vSeed1.x()-bs.x0(),vSeed1.y()-bs.y0(),vSeed1.z()-bs.z0());
00306         PerigeeTrajectoryError seedPerigeeErrors = PerigeeConversions::ftsToPerigeeError(tsAtClosestApproachSeed.trackStateAtPCA());
00307 
00308         //fixme
00309         reco::Track* trackFromSeed = new reco::Track(0.,0., vSeed1, pSeed, 1, seedPerigeeErrors.covarianceMatrix());
00310         trackFromSeed->setHitPattern(seed->recHits().first,seed->recHits().second);
00311 
00312         bool isSigSimMatched(false);
00313         bool isSimMatched(false);
00314         bool isChargeMatched(true);
00315         int numAssocSeeds = 0;
00316         int tpbx = 0;
00317         int nSimHits = 0;
00318         double sharedFraction = 0.;
00319         std::vector<std::pair<TrackingParticleRef, double> > tp;
00320         if(recSimColl.find(seed) != recSimColl.end()) {
00321           tp = recSimColl[seed];
00322           if (tp.size()!=0) {
00323 
00324             std::vector<PSimHit> simhits=tp[0].first->trackPSimHit(DetId::Tracker);
00325             nSimHits = simhits.end()-simhits.begin();
00326             sharedFraction = tp[0].second;
00327             isSimMatched = true;
00328             if (tp[0].first->charge() != seed->startingState().parameters().charge()) isChargeMatched = false;
00329             if(simRecColl.find(tp[0].first) != simRecColl.end()) numAssocSeeds = simRecColl[tp[0].first].size();
00330             //std::cout << numAssocRecoTracks << std::endl;
00331             tpbx = tp[0].first->eventId().bunchCrossing();
00332 
00333             at++;
00334 
00335             for (unsigned int tp_ite=0;tp_ite<tp.size();++tp_ite){ 
00336               TrackingParticle trackpart = *(tp[tp_ite].first);
00337               if ((trackpart.eventId().event() == 0) && (trackpart.eventId().bunchCrossing() == 0)){
00338                 isSigSimMatched = true;
00339                 sat++;
00340                 break;
00341               }
00342             }
00343 
00344 
00345             edm::LogVerbatim("SeedValidator") << "TrajectorySeed #" << rT << " associated with quality:" << tp.begin()->second <<"\n";
00346           }
00347         } else {
00348           edm::LogVerbatim("SeedValidator") << "TrajectorySeed #" << rT << " NOT associated to any TrackingParticle" << "\n";             
00349         }
00350 
00351         histoProducerAlgo_->fill_generic_recoTrack_histos(w,*trackFromSeed,bs.position(),isSimMatched,isSigSimMatched, 
00352                                                           isChargeMatched, numAssocSeeds, puinfo.getPU_NumInteractions(), 
00353                                                           tpbx, nSimHits, sharedFraction);
00354         
00355         //Fill other histos
00356         try{
00357           if (tp.size()==0) continue;
00358 
00359           histoProducerAlgo_->fill_simAssociated_recoTrack_histos(w,*trackFromSeed);
00360 
00361           TrackingParticleRef tpr = tp.begin()->first;
00362         
00363           //compute tracking particle parameters at point of closest approach to the beamline
00364           ParticleBase::Vector momentumTP = parametersDefinerTP->momentum(event,setup,*(tpr.get()));
00365           ParticleBase::Point vertexTP = parametersDefinerTP->vertex(event,setup,*(tpr.get()));                  
00366 
00367           //      LogTrace("SeedValidatorTEST") << "assocChi2=" << tp.begin()->second << "\n"
00368           //                                     << "" <<  "\n"
00369           //                                     << "ptREC=" << ptSeed << "\n"
00370           //                                     << "etaREC=" << etaSeed << "\n"
00371           //                                     << "qoverpREC=" << qoverpSeed << "\n"
00372           //                                     << "dxyREC=" << dxySeed << "\n"
00373           //                                     << "dzREC=" << dzSeed << "\n"
00374           //                                     << "thetaREC=" << thetaSeed << "\n"
00375           //                                     << "phiREC=" << phiSeed << "\n"
00376           //                                     << "" <<  "\n"
00377           //                                     << "qoverpError()=" << qoverpErrorSeed << "\n"
00378           //                                     << "dxyError()=" << dxyErrorSeed << "\n"
00379           //                                     << "dzError()=" << dzErrorSeed << "\n"
00380           //                                     << "thetaError()=" << lambdaErrorSeed << "\n"
00381           //                                     << "phiError()=" << phiErrorSeed << "\n"
00382           //                                     << "" <<  "\n"
00383           //                                     << "ptSIM=" << sqrt(assocTrack->momentum().perp2()) << "\n"
00384           //                                     << "etaSIM=" << assocTrack->momentum().Eta() << "\n"    
00385           //                                     << "qoverpSIM=" << qoverpSim << "\n"
00386           //                                     << "dxySIM=" << dxySim << "\n"
00387           //                                     << "dzSIM=" << dzSim << "\n"
00388           //                                     << "thetaSIM=" << M_PI/2-lambdaSim << "\n"
00389           //                                     << "phiSIM=" << phiSim << "\n"
00390           //                                     << "" << "\n"
00391           //                                     << "contrib_Qoverp=" << contrib_Qoverp << "\n"
00392           //                                     << "contrib_dxy=" << contrib_dxy << "\n"
00393           //                                     << "contrib_dz=" << contrib_dz << "\n"
00394           //                                     << "contrib_theta=" << contrib_theta << "\n"
00395           //                                     << "contrib_phi=" << contrib_phi << "\n"
00396           //                                     << "" << "\n"
00397           //                                     <<"chi2PULL="<<contrib_Qoverp+contrib_dxy+contrib_dz+contrib_theta+contrib_phi<<"\n";
00398 
00399           histoProducerAlgo_->fill_ResoAndPull_recoTrack_histos(w,momentumTP,vertexTP,tpr->charge(),
00400                                                                 *trackFromSeed,bs.position());
00401 
00402 
00403         } catch (cms::Exception e){
00404           LogTrace("SeedValidator") << "exception found: " << e.what() << "\n";
00405         }
00406       }// End of for(TrajectorySeedCollection::size_type i=0; i<seedCollection->size(); ++i)
00407 
00408       histoProducerAlgo_->fill_trackBased_histos(w,at,rT,st);
00409 
00410       edm::LogVerbatim("SeedValidator") << "Total Simulated: " << st << "\n"
00411                                          << "Total Associated (simToReco): " << ats << "\n"
00412                                          << "Total Reconstructed: " << rT << "\n"
00413                                          << "Total Associated (recoToSim): " << at << "\n"
00414                                          << "Total Fakes: " << rT-at << "\n";
00415       w++;
00416     }
00417   }
00418 }
00419 
00420 void TrackerSeedValidator::endRun(edm::Run const&, edm::EventSetup const&) {
00421   LogTrace("SeedValidator") << "TrackerSeedValidator::endRun()";
00422   int w=0;
00423   for (unsigned int ww=0;ww<associators.size();ww++){
00424     for (unsigned int www=0;www<label.size();www++){
00425 
00426       if (runStandalone) histoProducerAlgo_->finalHistoFits(w);
00427       if (runStandalone) histoProducerAlgo_->fillProfileHistosFromVectors(w);
00428       histoProducerAlgo_->fillHistosFromVectors(w);
00429 
00430       w++;
00431     }
00432   }
00433   if ( out.size() != 0 && dbe_ ) dbe_->save(out);
00434 }
00435 
00436 
00437