CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/SimTracker/TrackHistory/plugins/SVTagInfoValidationAnalyzer.cc

Go to the documentation of this file.
00001 #include <map>
00002 #include <string>
00003 #include <cmath>
00004 #include <Math/Functions.h>
00005 #include <Math/SVector.h>
00006 #include <Math/SMatrix.h>
00007 #include "TH1F.h"
00008 #include "TH1D.h"
00009 #include "TH2D.h"
00010 #include "TMath.h"
00011 #include "DataFormats/BTauReco/interface/SecondaryVertexTagInfo.h"
00012 
00013 #include "FWCore/Framework/interface/Frameworkfwd.h"
00014 #include "FWCore/Framework/interface/EDAnalyzer.h"
00015 #include "FWCore/Framework/interface/MakerMacros.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 #include "FWCore/ServiceRegistry/interface/Service.h"
00018 
00019 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00020 
00021 #include "SimTracker/TrackHistory/interface/VertexClassifierByProxy.h"
00022 #include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
00023 #include "DataFormats/Math/interface/Vector.h"
00024 #include "TROOT.h"
00025 #include "Math/VectorUtil.h"
00026 #include <TVector3.h>
00027 #include <Math/GenVector/PxPyPzE4D.h>
00028 #include <Math/GenVector/PxPyPzM4D.h>
00029 #include "DataFormats/Math/interface/LorentzVector.h"
00030 //
00031 // class decleration
00032 //
00033 using namespace reco;
00034 using namespace std;
00035 using namespace edm;
00036 
00037 class SVTagInfoValidationAnalyzer : public edm::EDAnalyzer
00038 {
00039 
00040 public:
00041 
00042     explicit SVTagInfoValidationAnalyzer(const edm::ParameterSet&);
00043 
00044 private:
00045 
00046     virtual void analyze(const edm::Event&, const edm::EventSetup&);
00047     virtual void endJob() ;
00048     // Member data
00049 
00050     VertexClassifierByProxy<reco::SecondaryVertexTagInfoCollection> classifier_;
00051 
00052     Int_t numberVertexClassifier_;
00053     edm::InputTag trackingTruth_;
00054     edm::InputTag svTagInfoProducer_;
00055 
00056     // Get the file service
00057     edm::Service<TFileService> fs_;
00058 
00059     Int_t n_event;
00060     Int_t rs_total_nall ; 
00061     Int_t rs_total_nsv ;
00062     Int_t rs_total_nbv ;
00063     Int_t rs_total_nbsv ;
00064     Int_t rs_total_ncv ;
00065     Int_t rs_total_nlv ;
00066     Int_t total_nfake ;
00067 
00068     Int_t sr_total_nall ; 
00069     Int_t sr_total_nsv ;
00070     Int_t sr_total_nbv ;
00071     Int_t sr_total_nbsv ;
00072     Int_t sr_total_ncv ;
00073     Int_t sr_total_nlv ;
00074     Int_t total_nmiss ;
00075         
00076 
00077   // Bookeeping of all the histograms per category
00078   void bookRecoToSim(std::string const &);
00079     void bookSimToReco(std::string const &);
00080 
00081     // Fill all histogram per category
00082     void fillRecoToSim(std::string const &, reco::Vertex const &, TrackingVertexRef const &);
00083     void fillSimToReco(std::string const &, reco::VertexBaseRef const &, TrackingVertexRef const &);
00084 
00085     // Histogram handlers
00086     std::map<std::string, TH1D *> TH1Index_;
00087 
00088 };
00089 
00090 
00091 SVTagInfoValidationAnalyzer::SVTagInfoValidationAnalyzer(const edm::ParameterSet& config) : classifier_(config)
00092 {
00093   //Initialize counters
00094     n_event = 0;
00095     rs_total_nall = 0; 
00096     rs_total_nsv = 0;
00097     rs_total_nbv = 0;
00098     rs_total_nbsv = 0;
00099     rs_total_ncv = 0;
00100     rs_total_nlv = 0;
00101     total_nfake = 0;
00102 
00103     sr_total_nall = 0; 
00104     sr_total_nsv = 0;
00105     sr_total_nbv = 0;
00106     sr_total_nbsv = 0;
00107     sr_total_ncv = 0;
00108     sr_total_nlv = 0;
00109     total_nmiss = 0;
00110 
00111     // Get the track collection
00112     svTagInfoProducer_ = config.getUntrackedParameter<edm::InputTag> ( "svTagInfoProducer" );
00113 
00114     // Name of the traking pariticle collection
00115     trackingTruth_ = config.getUntrackedParameter<edm::InputTag> ( "trackingTruth" );
00116 
00117     // Number of track categories
00118     numberVertexClassifier_ = VertexCategories::Unknown+1;
00119 
00120     // Define histogram for counting categories
00121     TH1Index_["VertexClassifier"] = fs_->make<TH1D>(
00122                                         "VertexClassifier",
00123                                         "Frequency for the different track categories",
00124                                         numberVertexClassifier_,
00125                                         -0.5,
00126                                         numberVertexClassifier_ - 0.5
00127                                     );
00128     
00129     //--- RecoToSim
00130     TH1Index_["rs_All_MatchQuality"]= fs_->make<TH1D>( "rs_All_MatchQuality", "Quality of Match", 51, -0.01, 1.01 );
00131     TH1Index_["rs_All_FlightDistance2d"]= fs_->make<TH1D>( "rs_All_FlightDistance2d", "Transverse flight distance [cm]", 100, 0, 5 );
00132     TH1Index_["rs_SecondaryVertex_FlightDistance2d"]= fs_->make<TH1D>( "rs_SecondaryVertex_FlightDistance2d", "Transverse flight distance [cm]", 100, 0, 5 );
00133     TH1Index_["rs_BSV_FlightDistance2d"]= fs_->make<TH1D>( "rs_BSV_FlightDistance2d", "Transverse flight distance [cm]", 100, 0, 5 );
00134     TH1Index_["rs_BWeakDecay_FlightDistance2d"]= fs_->make<TH1D>( "rs_BWeakDecay_FlightDistance2d", "Transverse flight distance [cm]", 100, 0, 5 );
00135     TH1Index_["rs_CWeakDecay_FlightDistance2d"]= fs_->make<TH1D>( "rs_CWeakDecay_FlightDistance2d", "Transverse flight distance [cm]", 100, 0, 5 );
00136     TH1Index_["rs_Light_FlightDistance2d"]= fs_->make<TH1D>( "rs_Light_FlightDistance2d", "Transverse flight distance [cm]", 100, 0, 5 );
00137 
00138     TH1Index_["rs_All_nRecVtx"]= fs_->make<TH1D>( "rs_All_nRecVtx", "Number of Vertices per event", 11, -0.5, 10.5 );
00139     TH1Index_["rs_SecondaryVertex_nRecVtx"]= fs_->make<TH1D>( "rs_SecondaryVertex_nRecVtx", "Number of Vertices per event", 11, -0.5, 10.5 );
00140     TH1Index_["rs_BSV_nRecVtx"]= fs_->make<TH1D>( "rs_BSV_nRecVtx", "Number of Vertices per event", 11, -0.5, 10.5 );
00141     TH1Index_["rs_BWeakDecay_nRecVtx"]= fs_->make<TH1D>( "rs_BWeakDecay_nRecVtx", "Number of Vertices per event", 11, -0.5, 10.5 );
00142     TH1Index_["rs_CWeakDecay_nRecVtx"]= fs_->make<TH1D>( "rs_CWeakDecay_nRecVtx", "Number of Vertices per event", 11, -0.5, 10.5 );
00143     TH1Index_["rs_Light_nRecVtx"]= fs_->make<TH1D>( "rs_Light_nRecVtx", "Number of Vertices per event", 11, -0.5, 10.5 );
00144 
00145 
00146     //--- SimToReco
00147     TH1Index_["sr_All_MatchQuality"]= fs_->make<TH1D>( "sr_All_MatchQuality", "Quality of Match", 51, -0.01, 1.01);
00148     TH1Index_["sr_All_nRecVtx"]= fs_->make<TH1D>( "sr_All_nRecVtx", "Number of Vertices per event", 11, -0.5, 10.5 );
00149     TH1Index_["sr_SecondaryVertex_nRecVtx"]= fs_->make<TH1D>( "sr_SecondaryVertex_nRecVtx", "Number of Vertices per event", 11, -0.5, 10.5 );
00150     TH1Index_["sr_BSV_nRecVtx"]= fs_->make<TH1D>( "sr_BSV_nRecVtx", "Number of Vertices per event", 11, -0.5, 10.5 );
00151     TH1Index_["sr_BWeakDecay_nRecVtx"]= fs_->make<TH1D>( "sr_BWeakDecay_nRecVtx", "Number of Vertices per event", 11, -0.5, 10.5 );
00152     TH1Index_["sr_CWeakDecay_nRecVtx"]= fs_->make<TH1D>( "sr_CWeakDecay_nRecVtx", "Number of Vertices per event", 11, -0.5, 10.5 );
00153     TH1Index_["sr_Light_nRecVtx"]= fs_->make<TH1D>( "sr_Light_nRecVtx", "Number of Vertices per event", 11, -0.5, 10.5 );  
00154 
00155 
00156     // Set the proper categories names
00157     for (Int_t i = 0; i < numberVertexClassifier_; ++i)
00158         TH1Index_["VertexClassifier"]->GetXaxis()->SetBinLabel(i+1, VertexCategories::Names[i]);
00159 
00160     // book histograms
00161     bookRecoToSim("rs_All");
00162     bookRecoToSim("rs_SecondaryVertex");
00163     bookRecoToSim("rs_BSV");
00164     bookRecoToSim("rs_BWeakDecay");
00165     bookRecoToSim("rs_CWeakDecay");
00166     bookRecoToSim("rs_Light");
00167 
00168     bookSimToReco("sr_All"); 
00169     bookSimToReco("sr_SecondaryVertex");
00170     bookSimToReco("sr_BSV");  
00171     bookSimToReco("sr_BWeakDecay");  
00172     bookSimToReco("sr_CWeakDecay");  
00173     bookSimToReco("sr_Light");  
00174 }   
00175 
00176 
00177 void SVTagInfoValidationAnalyzer::analyze(const edm::Event& event, const edm::EventSetup& setup)
00178 {
00179   ++n_event;
00180 
00181   std::cout << "*** Analyzing " << event.id() << " n_event = " << n_event << std::endl << std::endl;
00182 
00183   // Set the classifier for a new event
00184   classifier_.newEvent(event, setup);
00185   
00186 
00187   // Vertex collection
00188   edm::Handle<reco::SecondaryVertexTagInfoCollection> svTagInfoCollection;
00189   event.getByLabel(svTagInfoProducer_, svTagInfoCollection);
00190 
00191   // Get a constant reference to the track history associated to the classifier
00192   VertexHistory const & tracer = classifier_.history();
00193 
00194   cout << "* Event " << n_event << " ; svTagInfoCollection->size() = " << svTagInfoCollection->size() << endl; 
00195   
00196   int rs_nall = 0; 
00197   int rs_nsv = 0;
00198   int rs_nbv = 0;
00199   int rs_nbsv = 0;
00200   int rs_ncv = 0;
00201   int rs_nlv = 0;
00202   int nfake = 0;
00203 
00204   int sr_nall = 0; 
00205   int sr_nsv = 0;
00206   int sr_nbv = 0;
00207   int sr_nbsv = 0;
00208   int sr_ncv = 0;
00209   int sr_nlv = 0;
00210   int nmiss = 0;
00211 
00212   // Loop over the svTagInfo collection.
00213   for (std::size_t index = 0; index < svTagInfoCollection->size(); ++index){
00214 
00215     reco::SecondaryVertexTagInfoRef svTagInfo(svTagInfoCollection, index);
00216 
00217     // Loop over the vertexes in svTagInfo
00218     for ( std::size_t vindex = 0; vindex < svTagInfo->nVertices(); ++vindex ){
00219 
00220  
00221       // Classify the vertices
00222       classifier_.evaluate(svTagInfo, vindex);
00223 
00224       //quality of the match
00225       double rs_quality = tracer.quality();
00226  
00227 
00228       // Fill the histogram with the categories
00229       for (Int_t i = 0; i != numberVertexClassifier_; ++i) {
00230 
00231         if ( classifier_.is( (VertexCategories::Category) i ) ) {
00232 
00233           TH1Index_["VertexClassifier"]->Fill(i);
00234         }
00235       }
00236       if ( !classifier_.is(VertexCategories::Fake) ) {
00237 
00238         cout << "R2S: MatchQuality = " << rs_quality << endl;
00239 
00240         TH1Index_["rs_All_MatchQuality"]->Fill( rs_quality );
00241         fillRecoToSim("rs_All", svTagInfo->secondaryVertex(vindex), tracer.simVertex());
00242         TH1Index_["rs_All_FlightDistance2d"]->Fill( svTagInfo->flightDistance( vindex, true ).value() );
00243         rs_nall++;
00244 
00245         if ( classifier_.is(VertexCategories::SecondaryVertex) ) {
00246 
00247           //cout << "    R2S VertexCategories::SecondaryVertex" << endl;
00248           fillRecoToSim("rs_SecondaryVertex", svTagInfo->secondaryVertex(vindex), tracer.simVertex());
00249           TH1Index_["rs_SecondaryVertex_FlightDistance2d"]->Fill( svTagInfo->flightDistance( vindex, true ).value() );
00250           rs_nsv++;
00251         }
00252         
00253         if ( classifier_.is(VertexCategories::BWeakDecay) ) {
00254 
00255           //cout << "    R2S VertexCategories::BWeakDecay" << endl;
00256           fillRecoToSim("rs_BWeakDecay", svTagInfo->secondaryVertex(vindex), tracer.simVertex());
00257           TH1Index_["rs_BWeakDecay_FlightDistance2d"]->Fill( svTagInfo->flightDistance( vindex, true ).value() );
00258           rs_nbv++;
00259 
00260           if ( classifier_.is(VertexCategories::SecondaryVertex) ) {
00261             //cout << "    R2S VertexCategories::BWeakDecay SecondaryVertex" << endl;
00262             fillRecoToSim("rs_BSV", svTagInfo->secondaryVertex(vindex), tracer.simVertex());
00263             TH1Index_["rs_BSV_FlightDistance2d"]->Fill( svTagInfo->flightDistance( vindex, true ).value() );
00264             rs_nbsv++;   
00265 
00266           }          
00267         }//BWeakDecay
00268 
00269         else if ( classifier_.is(VertexCategories::CWeakDecay) ) {
00270 
00271           //cout << "    R2S VertexCategories::CWeakDecay" << endl;
00272           fillRecoToSim("rs_CWeakDecay", svTagInfo->secondaryVertex(vindex), tracer.simVertex());
00273           TH1Index_["rs_CWeakDecay_FlightDistance2d"]->Fill( svTagInfo->flightDistance( vindex, true ).value() );
00274           rs_ncv++;
00275 
00276         }                
00277         else {
00278           //cout << "    R2S Light (rest of categories)" << endl;
00279           fillRecoToSim("rs_Light", svTagInfo->secondaryVertex(vindex), tracer.simVertex());
00280           TH1Index_["rs_Light_FlightDistance2d"]->Fill( svTagInfo->flightDistance( vindex, true ).value() );
00281           rs_nlv++;
00282         }
00283       }//end if classifier
00284 
00285       else {
00286         cout << "    VertexCategories::Fake!!" << endl;
00287         nfake++;
00288       }
00289 
00290    }//end loop over vertices in svTagInfo
00291 
00292   }//loop over svTagInfo
00293 
00294   TH1Index_["rs_All_nRecVtx"]->Fill( rs_nall );
00295   TH1Index_["rs_SecondaryVertex_nRecVtx"]->Fill( rs_nsv );
00296   TH1Index_["rs_BWeakDecay_nRecVtx"]->Fill( rs_nbv );
00297   TH1Index_["rs_BSV_nRecVtx"]->Fill( rs_nbsv );
00298   TH1Index_["rs_CWeakDecay_nRecVtx"]->Fill( rs_ncv );  
00299   TH1Index_["rs_Light_nRecVtx"]->Fill( rs_nlv );  
00300   cout << endl;  
00301 
00302   //----------------------------------------------------------------
00303   // SIM TO RECO!
00304 
00305   // Vertex collection
00306   edm::Handle<TrackingVertexCollection>  TVCollection;
00307   event.getByLabel(trackingTruth_, TVCollection);
00308     
00309   // Loop over the TV collection.
00310   for (std::size_t index = 0; index < TVCollection->size(); ++index){
00311         
00312     TrackingVertexRef trackingVertex(TVCollection, index);
00313 
00314     classifier_.evaluate(trackingVertex);
00315         
00316     double sr_quality = tracer.quality();
00317  
00318     if ( classifier_.is(VertexCategories::Reconstructed) ) {
00319 
00320       cout << "S2R: MatchQuality = " << sr_quality << endl;
00321 
00322       //cout << "    S2R VertexCategories::Reconstructed" << endl;
00323       TH1Index_["sr_All_MatchQuality"]->Fill( sr_quality );      
00324       fillSimToReco("sr_All", tracer.recoVertex(), trackingVertex);
00325       sr_nall++;
00326 
00327       if ( classifier_.is(VertexCategories::SecondaryVertex) ) {
00328 
00329         //cout << "    S2R VertexCategories::Reconstructed::SecondaryVertex" << endl;
00330         fillSimToReco("sr_SecondaryVertex", tracer.recoVertex(), trackingVertex);
00331         sr_nsv++;
00332       }
00333 
00334       if ( classifier_.is(VertexCategories::BWeakDecay) ) {
00335 
00336         //cout << "    S2R VertexCategories::Reconstructed::BWeakDecay" << endl;
00337         fillSimToReco("sr_BWeakDecay", tracer.recoVertex(), trackingVertex);
00338         sr_nbv++;
00339 
00340         if ( classifier_.is(VertexCategories::SecondaryVertex) ) {
00341 
00342           //cout << "    S2R VertexCategories::Reconstructed::BWeakDecay SecondaryVertex" << endl;
00343          fillSimToReco("sr_BSV", tracer.recoVertex(), trackingVertex);
00344          sr_nbsv++;
00345         }
00346     
00347       }//BWeakDecay
00348 
00349       else if ( classifier_.is(VertexCategories::CWeakDecay) ) {
00350 
00351         //cout << "    S2R VertexCategories::CWeakDecay" << endl;
00352         fillSimToReco("sr_CWeakDecay", tracer.recoVertex(), trackingVertex);
00353         sr_ncv++;
00354       }
00355  
00356       else {
00357 
00358         //cout << "    S2R Light (rest of categories)" << endl;
00359         fillSimToReco("sr_Light", tracer.recoVertex(), trackingVertex);
00360         sr_nlv++;
00361       }
00362 
00363     }//Reconstructed
00364     else {
00365       //cout << "##### Not reconstructed!" << endl;
00366       nmiss++;
00367     }
00368 
00369   }//TVCollection.size()
00370 
00371   TH1Index_["sr_All_nRecVtx"]->Fill( sr_nall );
00372   TH1Index_["sr_SecondaryVertex_nRecVtx"]->Fill( sr_nsv );
00373   TH1Index_["sr_BWeakDecay_nRecVtx"]->Fill( sr_nbv );
00374   TH1Index_["sr_BSV_nRecVtx"]->Fill( sr_nbsv );
00375   TH1Index_["sr_CWeakDecay_nRecVtx"]->Fill( sr_ncv );  
00376   TH1Index_["sr_Light_nRecVtx"]->Fill( rs_nlv );  
00377 
00378   rs_total_nall += rs_nall;  
00379   rs_total_nsv += rs_nsv;
00380   rs_total_nbv += rs_nbv;
00381   rs_total_nbsv += rs_nbsv;
00382   rs_total_ncv += rs_ncv ;
00383   rs_total_nlv += rs_nlv;
00384   total_nfake += nfake;
00385 
00386   sr_total_nall += sr_nall; 
00387   sr_total_nsv += sr_nsv;
00388   sr_total_nbv += sr_nbv;
00389   sr_total_nbsv += sr_nbsv;
00390   sr_total_ncv += sr_ncv;
00391   sr_total_nlv += sr_nlv;
00392   total_nmiss += nmiss;
00393 
00394 
00395 }
00396 
00397 
00398 void SVTagInfoValidationAnalyzer::bookRecoToSim(std::string const & prefix){
00399   // Book pull histograms
00400 
00401   std::string name = prefix + "_Pullx";
00402   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, -10., 10.);
00403   name = prefix + "_Pully";
00404   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, -10., 10.);
00405   name = prefix + "_Pullz";
00406   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, -10., 10.);
00407 
00408    name = prefix + "_Resx";
00409   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, -0.05, 0.05);
00410   name = prefix + "_Resy";
00411   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, -0.05, 0.05); 
00412   name = prefix + "_Resz";
00413   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, -0.05, 0.05); 
00414 
00415   name = prefix + "_Chi2Norm";
00416   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, 0, 10.); 
00417   name = prefix + "_Chi2Prob";
00418   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, 0., 1.); 
00419  
00420   name = prefix + "_nRecTrks";
00421   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 501, -0.5, 500.5);
00422 
00423   name = prefix + "_AverageTrackWeight";
00424   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, -0.1, 1.1);
00425 
00426   name = prefix + "_Mass";
00427   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 65, 0., 6.5);
00428 
00429   name = prefix + "_RecPt";
00430   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 2000, 0., 1000.);
00431 
00432   name = prefix + "_RecEta";
00433   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 200, -3., 3.);
00434 
00435   name = prefix + "_RecCharge";
00436   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 21, -0.5, 20.5);
00437 
00438   name = prefix + "_RecTrackPt";
00439   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 2000, 0., 1000.);
00440 
00441   name = prefix + "_RecTrackEta";
00442   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 200, -3., 3.);
00443 
00444   name = prefix + "_nSimTrks";
00445   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 501, -0.5, 500.5); 
00446 
00447   name = prefix + "_SimPt";
00448   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 2000, 0., 1000.);
00449 
00450   name = prefix + "_SimEta";
00451   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 200, -3., 3.);
00452 
00453   name = prefix + "_SimCharge";
00454   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 21, -0.5, 20.5);
00455   
00456   name = prefix + "_SimTrackPt";
00457   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 500, 0., 500.);
00458 
00459   name = prefix + "_SimTrackEta";
00460   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 200, -3., 3.);
00461 
00462 
00463 }
00464 
00465 
00466 void SVTagInfoValidationAnalyzer::bookSimToReco(std::string const & prefix){
00467   // Book pull histograms
00468 
00469   std::string name = prefix + "_Pullx";
00470   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, -10., 10.);
00471   name = prefix + "_Pully";
00472   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, -10., 10.);
00473   name = prefix + "_Pullz";
00474   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, -10., 10.);
00475 
00476    name = prefix + "_Resx";
00477   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, -0.05, 0.05);
00478   name = prefix + "_Resy";
00479   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, -0.05, 0.05); 
00480   name = prefix + "_Resz";
00481   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, -0.05, 0.05); 
00482 
00483   name = prefix + "_Chi2Norm";
00484   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, 0, 10.); 
00485   name = prefix + "_Chi2Prob";
00486   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, 0., 1.); 
00487  
00488   name = prefix + "_nRecTrks";
00489   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 501, -0.5, 500.5);
00490 
00491   name = prefix + "_AverageTrackWeight";
00492   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 100, -0.1, 1.1);
00493 
00494   name = prefix + "_Mass";
00495   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 65, 0., 6.5);
00496 
00497   name = prefix + "_RecPt";
00498   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 2000, 0., 1000.);
00499 
00500   name = prefix + "_RecEta";
00501   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 200, -3., 3.);
00502 
00503   name = prefix + "_RecCharge";
00504   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 21, -0.5, 20.5);
00505 
00506   name = prefix + "_RecTrackPt";
00507   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 2000, 0., 1000.);
00508 
00509   name = prefix + "_RecTrackEta";
00510   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 200, -3., 3.);
00511 
00512   name = prefix + "_nSimTrks";
00513   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 501, -0.5, 500.5); 
00514 
00515   name = prefix + "_SimPt";
00516   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 2000, 0., 1000.);
00517 
00518   name = prefix + "_SimEta";
00519   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 200, -3., 3.);
00520 
00521   name = prefix + "_SimCharge";
00522   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 21, -0.5, 20.5);
00523   
00524   name = prefix + "_SimTrackPt";
00525   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 500, 0., 500.);
00526 
00527   name = prefix + "_SimTrackEta";
00528   TH1Index_[name] = fs_->make<TH1D>(name.c_str(), name.c_str(), 200, -3., 3.);
00529 
00530 
00531 }
00532 
00533 
00534 void SVTagInfoValidationAnalyzer::fillRecoToSim(std::string const & prefix, reco::Vertex const & vertex, TrackingVertexRef const & simVertex)
00535 {
00536  
00537   double pullx = (vertex.x() - simVertex->position().x())/vertex.xError();
00538   double pully = (vertex.y() - simVertex->position().y())/vertex.yError();
00539   double pullz = (vertex.z() - simVertex->position().z())/vertex.zError();
00540 
00541   double resx = vertex.x() - simVertex->position().x();
00542   double resy = vertex.y() - simVertex->position().y();
00543   double resz = vertex.z() - simVertex->position().z();
00544 
00545   double chi2norm = vertex.normalizedChi2();
00546   double chi2prob = ChiSquaredProbability( vertex.chi2(), vertex.ndof() ); 
00547 
00548  double sum_weight = 0.;
00549  double weight = 0.;
00550  double tracksize =  vertex.tracksSize();
00551  math::XYZVector momentum;
00552  math::XYZTLorentzVector sum;
00553  int charge = 0;
00554   double thePiMass = 0.13957;
00555   for ( reco::Vertex::trackRef_iterator recDaughter = vertex.tracks_begin() ; recDaughter != vertex.tracks_end(); ++recDaughter) {
00556 
00557     ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec;
00558 
00559     vec.SetPx( (**recDaughter).px() );
00560     vec.SetPy( (**recDaughter).py() );
00561     vec.SetPz( (**recDaughter).pz() );
00562     vec.SetM(thePiMass); 
00563 
00564     sum += vec;
00565   
00566     weight = vertex.trackWeight(*recDaughter); 
00567     sum_weight += weight;
00568 
00569     math::XYZVector p;
00570     p = (**recDaughter).momentum();
00571     momentum += p; 
00572 
00573     charge += (*recDaughter)->charge();
00574 
00575     TH1Index_[prefix + "_RecTrackPt"]->Fill( (*recDaughter)->pt() );
00576     TH1Index_[prefix + "_RecTrackEta"]->Fill( (*recDaughter)->eta() );      
00577   }//end loop to recDaughters
00578   //cout << "                   average sum of weights = " << sum_weight/tracksize << endl; 
00579   
00580   double mass = sum.M();
00581   double pt = sqrt( momentum.Perp2() );
00582   double eta = momentum.Eta();
00583 
00584   math::XYZVector simmomentum;
00585   int simcharge = 0;
00586   for ( TrackingVertex::tp_iterator simDaughter = simVertex->daughterTracks_begin(); simDaughter != simVertex->daughterTracks_end(); ++simDaughter ){
00587 
00588     math::XYZVector p;
00589     p = (**simDaughter).momentum();
00590     simmomentum += p; 
00591 
00592     simcharge += (*simDaughter)->charge();
00593 
00594     TH1Index_[prefix + "_SimTrackPt"]->Fill( (*simDaughter)->pt() );
00595     TH1Index_[prefix + "_SimTrackEta"]->Fill( (*simDaughter)->eta() );   
00596   }
00597 
00598   double simpt = sqrt( simmomentum.Perp2() );
00599   double simeta = simmomentum.Eta();
00600 
00601   //cout << "[fillRecoToSim]  vertex.tracksSize() = " << vertex.tracksSize() << " ; simVertex->nDaughterTracks() = " << simVertex->nDaughterTracks() << endl;  
00602 
00603   TH1Index_[prefix + "_nRecTrks"]->Fill( vertex.tracksSize() );
00604   TH1Index_[prefix + "_nSimTrks"]->Fill( simVertex->nDaughterTracks() );
00605   TH1Index_[prefix + "_Pullx"]->Fill(pullx);
00606   TH1Index_[prefix + "_Pully"]->Fill(pully);
00607   TH1Index_[prefix + "_Pullz"]->Fill(pullz);
00608   TH1Index_[prefix + "_Resx"]->Fill(resx);
00609   TH1Index_[prefix + "_Resy"]->Fill(resy);
00610   TH1Index_[prefix + "_Resz"]->Fill(resz);
00611   TH1Index_[prefix + "_AverageTrackWeight"]->Fill( sum_weight/tracksize );   
00612   TH1Index_[prefix + "_Chi2Norm"]->Fill(chi2norm);
00613   TH1Index_[prefix + "_Chi2Prob"]->Fill(chi2prob);
00614   TH1Index_[prefix + "_RecPt"]->Fill(pt);
00615   TH1Index_[prefix + "_RecEta"]->Fill(eta);
00616   TH1Index_[prefix + "_RecCharge"]->Fill(charge);
00617   TH1Index_[prefix + "_Mass"]->Fill(mass);
00618   TH1Index_[prefix + "_SimPt"]->Fill(simpt);
00619   TH1Index_[prefix + "_SimEta"]->Fill(simeta);
00620   TH1Index_[prefix + "_SimCharge"]->Fill(simcharge);  
00621 }
00622 
00623 
00624 void SVTagInfoValidationAnalyzer::fillSimToReco(std::string const & prefix, reco::VertexBaseRef const & vertex, TrackingVertexRef const & simVertex)
00625 {
00626 
00627   double pullx = (vertex->x() - simVertex->position().x())/vertex->xError();
00628   double pully = (vertex->y() - simVertex->position().y())/vertex->yError();
00629   double pullz = (vertex->z() - simVertex->position().z())/vertex->zError();
00630 
00631   double resx = vertex->x() - simVertex->position().x();
00632   double resy = vertex->y() - simVertex->position().y();
00633   double resz = vertex->z() - simVertex->position().z();
00634 
00635   double chi2norm = vertex->normalizedChi2();
00636   double chi2prob = ChiSquaredProbability(vertex->chi2(), vertex->ndof()); 
00637 
00638   double sum_weight = 0.;
00639   double weight = 0.;
00640   double tracksize =  vertex->tracksSize();
00641   math::XYZVector momentum;
00642   math::XYZTLorentzVector sum;
00643   int charge = 0;
00644   double thePiMass = 0.13957;
00645   for ( reco::Vertex::trackRef_iterator recDaughter = vertex->tracks_begin() ; recDaughter != vertex->tracks_end(); ++recDaughter ) {
00646 
00647     ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec;
00648 
00649     vec.SetPx( (**recDaughter).px() );
00650     vec.SetPy( (**recDaughter).py() );
00651     vec.SetPz( (**recDaughter).pz() );
00652     vec.SetM(thePiMass); 
00653 
00654     sum += vec;
00655   
00656     weight = vertex->trackWeight(*recDaughter); 
00657     sum_weight += weight;
00658 
00659     math::XYZVector p;
00660     p = (**recDaughter).momentum();
00661     momentum += p; 
00662 
00663     charge += (*recDaughter)->charge();
00664 
00665     TH1Index_[prefix + "_RecTrackPt"]->Fill( (*recDaughter)->pt() );
00666     TH1Index_[prefix + "_RecTrackEta"]->Fill( (*recDaughter)->eta() );  
00667   
00668   }
00669   //cout << "                   average sum of weights = " << sum_weight/tracksize << endl; 
00670 
00671   double mass = sum.M();
00672   double pt = sqrt( momentum.Perp2() );
00673   double eta = momentum.Eta();
00674 
00675    math::XYZVector simmomentum;
00676   int simcharge = 0;
00677   for ( TrackingVertex::tp_iterator simDaughter = simVertex->daughterTracks_begin(); simDaughter != simVertex->daughterTracks_end(); ++simDaughter ){
00678 
00679     math::XYZVector p;
00680     p = (**simDaughter).momentum();
00681     simmomentum += p; 
00682 
00683     simcharge += (*simDaughter)->charge();
00684 
00685     TH1Index_[prefix + "_SimTrackPt"]->Fill( (*simDaughter)->pt() );
00686     TH1Index_[prefix + "_SimTrackEta"]->Fill( (*simDaughter)->eta() );   
00687   }
00688 
00689   double simpt = sqrt( simmomentum.Perp2() );
00690   double simeta = simmomentum.Eta();
00691 
00692   //cout << "[fillSimToReco]  vertex->tracksSize() = " << vertex->tracksSize() << " ; simVertex->nDaughterTracks() = " << simVertex->nDaughterTracks() << endl;  
00693 
00694   TH1Index_[prefix + "_nRecTrks"]->Fill( vertex->tracksSize() );
00695   TH1Index_[prefix + "_nSimTrks"]->Fill( simVertex->nDaughterTracks() );
00696   TH1Index_[prefix + "_Pullx"]->Fill(pullx);
00697   TH1Index_[prefix + "_Pully"]->Fill(pully);
00698   TH1Index_[prefix + "_Pullz"]->Fill(pullz);
00699   TH1Index_[prefix + "_Resx"]->Fill(resx);
00700   TH1Index_[prefix + "_Resy"]->Fill(resy);
00701   TH1Index_[prefix + "_Resz"]->Fill(resz);
00702   TH1Index_[prefix + "_AverageTrackWeight"]->Fill( sum_weight/tracksize );   
00703   TH1Index_[prefix + "_Chi2Norm"]->Fill(chi2norm);
00704   TH1Index_[prefix + "_Chi2Prob"]->Fill(chi2prob);
00705   TH1Index_[prefix + "_RecPt"]->Fill(pt);
00706   TH1Index_[prefix + "_RecEta"]->Fill(eta);
00707   TH1Index_[prefix + "_RecCharge"]->Fill(charge);
00708   TH1Index_[prefix + "_Mass"]->Fill(mass);
00709   TH1Index_[prefix + "_SimPt"]->Fill(simpt);
00710   TH1Index_[prefix + "_SimEta"]->Fill(simeta);
00711   TH1Index_[prefix + "_SimCharge"]->Fill(simcharge);  
00712 
00713 }
00714 
00715 
00716 void 
00717 SVTagInfoValidationAnalyzer::endJob() {
00718 
00719   std::cout << std::endl;
00720   std::cout << " ====== Total Number of analyzed events: " << n_event << " ====== " << std::endl;
00721   std::cout << " ====== Total Number of R2S All:                         " << rs_total_nall << " ====== " << std::endl;
00722   std::cout << " ====== Total Number of R2S SecondaryVertex:             " << rs_total_nsv  << " ====== " << std::endl;
00723   std::cout << " ====== Total Number of R2S BWeakDecay:                  " << rs_total_nbv  << " ====== " << std::endl;
00724   std::cout << " ====== Total Number of R2S BWeakDecay::SecondaryVertex: " << rs_total_nbsv << " ====== " << std::endl;
00725   std::cout << " ====== Total Number of R2S CWeakDecay:                  " << rs_total_ncv  << " ====== " << std::endl;
00726   std::cout << " ====== Total Number of R2S Light:                       " << rs_total_nlv  << " ====== " << std::endl;
00727   std::cout << std::endl; 
00728   std::cout << " ====== Total Number of S2R All:                         " << sr_total_nall << " ====== " << std::endl;
00729   std::cout << " ====== Total Number of S2R SecondaryVertex:             " << sr_total_nsv  << " ====== " << std::endl;
00730   std::cout << " ====== Total Number of S2R BWeakDecay:                  " << sr_total_nbv  << " ====== " << std::endl;
00731   std::cout << " ====== Total Number of S2R BWeakDecay::SecondaryVertex: " << sr_total_nbsv << " ====== " << std::endl;
00732   std::cout << " ====== Total Number of S2R CWeakDecay:                  " << sr_total_ncv  << " ====== " << std::endl;
00733   std::cout << " ====== Total Number of S2R Light:                       " << sr_total_nlv  << " ====== " << std::endl;
00734   std::cout << std::endl; 
00735   std::cout << " ====== Total Number of Fake Vertices:              " << total_nfake   << " ====== " << std::endl;
00736  
00737 } 
00738 
00739 
00740 DEFINE_FWK_MODULE(SVTagInfoValidationAnalyzer);