CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/RecoEgamma/EgammaTools/src/ConversionLikelihoodCalculator.cc

Go to the documentation of this file.
00001 #include "RecoEgamma/EgammaTools/interface/ConversionLikelihoodCalculator.h"
00002 
00003 
00004 ConversionLikelihoodCalculator::ConversionLikelihoodCalculator()
00005 {
00006    reader_ = new TMVA::Reader("!Color:Silent");
00007    
00008 //   std::cout << "Init Reader()" << std::endl;
00009 
00010    reader_->AddVariable("log(e_over_p)",        &log_e_over_p_);
00011    reader_->AddVariable("log(abs(cot_theta))",  &log_abs_cot_theta_);
00012    reader_->AddVariable("log(abs(delta_phi))",  &log_abs_delta_phi_);
00013    reader_->AddVariable("log(chi2_max_pt)", &log_chi2_max_pt_);
00014    reader_->AddVariable("log(chi2_min_pt)", &log_chi2_min_pt_);
00015 
00016 }
00017 
00018 void ConversionLikelihoodCalculator::setWeightsFile(const char * weightsFile)
00019 {
00020 //   std::cout << "Before BookMVA " << weightsFile << std::endl;
00021    reader_->BookMVA("Likelihood", weightsFile);
00022 //   std::cout << "After  BookMVA" << std::endl;
00023 }
00024 
00025 double ConversionLikelihoodCalculator::calculateLikelihood(reco::ConversionRef conversion)
00026 {
00027    if (conversion->nTracks() != 2) return -1.;
00028    
00029    log_e_over_p_ = log(conversion->EoverP());
00030 
00031    log_abs_cot_theta_ = log(fabs(conversion->pairCotThetaSeparation()));
00032 
00033    double delta_phi = conversion->tracks()[0]->innerMomentum().phi()-conversion->tracks()[1]->innerMomentum().phi();
00034    double pi = 3.14159265;
00035    // phi normalization
00036    while (delta_phi > pi) delta_phi -= 2*pi;
00037    while (delta_phi < -pi) delta_phi += 2*pi;
00038    log_abs_delta_phi_ = log(fabs(delta_phi));
00039 
00040    double chi2_1 = conversion->tracks()[0]->normalizedChi2();
00041    double pt_1 = conversion->tracks()[0]->pt();
00042 
00043    double chi2_2 = conversion->tracks()[1]->normalizedChi2();
00044    double pt_2 = conversion->tracks()[1]->pt();
00045 
00046    double chi2_max_pt=chi2_1;
00047    double chi2_min_pt=chi2_2;
00048 
00049    if (pt_2 > pt_1) {
00050       chi2_max_pt=chi2_2;
00051       chi2_min_pt=chi2_1;
00052    }
00053 
00054    log_chi2_max_pt_ = log(chi2_max_pt);
00055    log_chi2_min_pt_ = log(chi2_min_pt);
00056 
00057 //   std::cout << "log_e_over_p_ " << log_e_over_p_ << std::endl;
00058 //   std::cout << "log_abs_cot_theta_ " << log_abs_cot_theta_ << std::endl;
00059 //   std::cout << "log_abs_delta_phi_ " << log_abs_delta_phi_ << std::endl;
00060 //   std::cout << "log_chi2_max_pt_ " << log_chi2_max_pt_ << std::endl;
00061 //   std::cout << "log_chi2_min_pt_ " << log_chi2_min_pt_ << std::endl;
00062    std::vector<Float_t> inputVec;
00063    inputVec.push_back(log_e_over_p_ );
00064    inputVec.push_back(log_abs_cot_theta_ );
00065    inputVec.push_back(log_abs_delta_phi_ );
00066    inputVec.push_back(log_chi2_max_pt_ );
00067    inputVec.push_back(log_chi2_min_pt_ );
00068    float like = reader_->EvaluateMVA(inputVec,"Likelihood");
00069 //   std::cout << "reader_->EvaluateMVA(\"Likelihood\") " << reader_->EvaluateMVA(inputVec,"Likelihood") << std::endl;
00070 
00071    return like;
00072 }
00073 
00074 double ConversionLikelihoodCalculator::calculateLikelihood(reco::Conversion& conversion)
00075 {
00076 
00077    if (conversion.nTracks() != 2) return -1.;
00078    
00079    log_e_over_p_ = log(conversion.EoverP());
00080 
00081    log_abs_cot_theta_ = log(fabs(conversion.pairCotThetaSeparation()));
00082 
00083 
00084    double delta_phi = conversion.tracksPin()[0].phi()-conversion.tracksPin()[1].phi();
00085    double pi = 3.14159265;
00086    // phi normalization
00087    while (delta_phi > pi) delta_phi -= 2*pi;
00088    while (delta_phi < -pi) delta_phi += 2*pi;
00089    log_abs_delta_phi_ = log(fabs(delta_phi));
00090 
00091    double chi2_1 = conversion.tracks()[0]->normalizedChi2();
00092    double pt_1 = conversion.tracks()[0]->pt();
00093 
00094    double chi2_2 = conversion.tracks()[1]->normalizedChi2();
00095    double pt_2 = conversion.tracks()[1]->pt();
00096 
00097    double chi2_max_pt=chi2_1;
00098    double chi2_min_pt=chi2_2;
00099 
00100    if (pt_2 > pt_1) {
00101       chi2_max_pt=chi2_2;
00102       chi2_min_pt=chi2_1;
00103    }
00104 
00105    log_chi2_max_pt_ = log(chi2_max_pt);
00106    log_chi2_min_pt_ = log(chi2_min_pt);
00107 
00108 //   std::cout << "log_e_over_p_ " << log_e_over_p_ << std::endl;
00109 //   std::cout << "log_abs_cot_theta_ " << log_abs_cot_theta_ << std::endl;
00110 //   std::cout << "log_abs_delta_phi_ " << log_abs_delta_phi_ << std::endl;
00111 //   std::cout << "log_chi2_max_pt_ " << log_chi2_max_pt_ << std::endl;
00112 //   std::cout << "log_chi2_min_pt_ " << log_chi2_min_pt_ << std::endl;
00113    std::vector<Float_t> inputVec;
00114    inputVec.push_back(log_e_over_p_ );
00115    inputVec.push_back(log_abs_cot_theta_ );
00116    inputVec.push_back(log_abs_delta_phi_ );
00117    inputVec.push_back(log_chi2_max_pt_ );
00118    inputVec.push_back(log_chi2_min_pt_ );
00119    float like = reader_->EvaluateMVA(inputVec,"Likelihood");
00120 //   std::cout << "reader_->EvaluateMVA(\"Likelihood\") " << reader_->EvaluateMVA(inputVec,"Likelihood") << std::endl;
00121 
00122    return like;
00123 }
00124 
00125