CMS 3D CMS Logo

ConversionLikelihoodCalculator.cc
Go to the documentation of this file.
2 
4 
6  reader_ = std::make_unique<TMVA::Reader>("!Color:Silent");
7 
8  // std::cout << "Init Reader()" << std::endl;
9 
10  reader_->AddVariable("log(e_over_p)", &log_e_over_p_);
11  reader_->AddVariable("log(abs(cot_theta))", &log_abs_cot_theta_);
12  reader_->AddVariable("log(abs(delta_phi))", &log_abs_delta_phi_);
13  reader_->AddVariable("log(chi2_max_pt)", &log_chi2_max_pt_);
14  reader_->AddVariable("log(chi2_min_pt)", &log_chi2_min_pt_);
15 }
16 
18  // std::cout << "Before BookMVA " << weightsFile << std::endl;
19  reader_->BookMVA("Likelihood", weightsFile);
20  // std::cout << "After BookMVA" << std::endl;
21 }
22 
24  if (conversion->nTracks() != 2)
25  return -1.;
26 
27  log_e_over_p_ = log(conversion->EoverP());
28 
29  log_abs_cot_theta_ = log(fabs(conversion->pairCotThetaSeparation()));
30 
31  double delta_phi = conversion->tracks()[0]->innerMomentum().phi() - conversion->tracks()[1]->innerMomentum().phi();
32  double pi = 3.14159265;
33  // phi normalization
34  while (delta_phi > pi)
35  delta_phi -= 2 * pi;
36  while (delta_phi < -pi)
37  delta_phi += 2 * pi;
39 
40  double chi2_1 = conversion->tracks()[0]->normalizedChi2();
41  double pt_1 = conversion->tracks()[0]->pt();
42 
43  double chi2_2 = conversion->tracks()[1]->normalizedChi2();
44  double pt_2 = conversion->tracks()[1]->pt();
45 
46  double chi2_max_pt = chi2_1;
47  double chi2_min_pt = chi2_2;
48 
49  if (pt_2 > pt_1) {
50  chi2_max_pt = chi2_2;
51  chi2_min_pt = chi2_1;
52  }
53 
54  log_chi2_max_pt_ = log(chi2_max_pt);
55  log_chi2_min_pt_ = log(chi2_min_pt);
56 
57  // std::cout << "log_e_over_p_ " << log_e_over_p_ << std::endl;
58  // std::cout << "log_abs_cot_theta_ " << log_abs_cot_theta_ << std::endl;
59  // std::cout << "log_abs_delta_phi_ " << log_abs_delta_phi_ << std::endl;
60  // std::cout << "log_chi2_max_pt_ " << log_chi2_max_pt_ << std::endl;
61  // std::cout << "log_chi2_min_pt_ " << log_chi2_min_pt_ << std::endl;
62  std::vector<Float_t> inputVec;
63  inputVec.push_back(log_e_over_p_);
64  inputVec.push_back(log_abs_cot_theta_);
65  inputVec.push_back(log_abs_delta_phi_);
66  inputVec.push_back(log_chi2_max_pt_);
67  inputVec.push_back(log_chi2_min_pt_);
68  float like = reader_->EvaluateMVA(inputVec, "Likelihood");
69  // std::cout << "reader_->EvaluateMVA(\"Likelihood\") " << reader_->EvaluateMVA(inputVec,"Likelihood") << std::endl;
70 
71  return like;
72 }
73 
75  if (conversion.nTracks() != 2)
76  return -1.;
77 
78  log_e_over_p_ = log(conversion.EoverP());
79 
80  log_abs_cot_theta_ = log(fabs(conversion.pairCotThetaSeparation()));
81 
82  double delta_phi = conversion.tracksPin()[0].phi() - conversion.tracksPin()[1].phi();
83  double pi = 3.14159265;
84  // phi normalization
85  while (delta_phi > pi)
86  delta_phi -= 2 * pi;
87  while (delta_phi < -pi)
88  delta_phi += 2 * pi;
90 
91  double chi2_1 = conversion.tracks()[0]->normalizedChi2();
92  double pt_1 = conversion.tracks()[0]->pt();
93 
94  double chi2_2 = conversion.tracks()[1]->normalizedChi2();
95  double pt_2 = conversion.tracks()[1]->pt();
96 
97  double chi2_max_pt = chi2_1;
98  double chi2_min_pt = chi2_2;
99 
100  if (pt_2 > pt_1) {
101  chi2_max_pt = chi2_2;
102  chi2_min_pt = chi2_1;
103  }
104 
105  log_chi2_max_pt_ = log(chi2_max_pt);
106  log_chi2_min_pt_ = log(chi2_min_pt);
107 
108  // std::cout << "log_e_over_p_ " << log_e_over_p_ << std::endl;
109  // std::cout << "log_abs_cot_theta_ " << log_abs_cot_theta_ << std::endl;
110  // std::cout << "log_abs_delta_phi_ " << log_abs_delta_phi_ << std::endl;
111  // std::cout << "log_chi2_max_pt_ " << log_chi2_max_pt_ << std::endl;
112  // std::cout << "log_chi2_min_pt_ " << log_chi2_min_pt_ << std::endl;
113  std::vector<Float_t> inputVec;
114  inputVec.push_back(log_e_over_p_);
115  inputVec.push_back(log_abs_cot_theta_);
116  inputVec.push_back(log_abs_delta_phi_);
117  inputVec.push_back(log_chi2_max_pt_);
118  inputVec.push_back(log_chi2_min_pt_);
119  float like = reader_->EvaluateMVA(inputVec, "Likelihood");
120  // std::cout << "reader_->EvaluateMVA(\"Likelihood\") " << reader_->EvaluateMVA(inputVec,"Likelihood") << std::endl;
121 
122  return like;
123 }
ConversionLikelihoodCalculator::reader_
std::unique_ptr< TMVA::Reader > reader_
Definition: ConversionLikelihoodCalculator.h:18
reco::Conversion
Definition: Conversion.h:23
edm::conversion
void conversion(EventAux const &from, EventAuxiliary &to)
Definition: EventAux.cc:9
ConversionLikelihoodCalculator::ConversionLikelihoodCalculator
ConversionLikelihoodCalculator()
Definition: ConversionLikelihoodCalculator.cc:5
ConversionLikelihoodCalculator::log_chi2_max_pt_
float log_chi2_max_pt_
Definition: ConversionLikelihoodCalculator.h:22
ConversionLikelihoodCalculator::calculateLikelihood
double calculateLikelihood(reco::ConversionRef conversion)
Definition: ConversionLikelihoodCalculator.cc:23
edm::Ref< ConversionCollection >
pfClustersFromHGC3DClusters_cfi.weightsFile
weightsFile
Definition: pfClustersFromHGC3DClusters_cfi.py:19
ConversionLikelihoodCalculator::log_chi2_min_pt_
float log_chi2_min_pt_
Definition: ConversionLikelihoodCalculator.h:23
ConversionLikelihoodCalculator::log_e_over_p_
float log_e_over_p_
Definition: ConversionLikelihoodCalculator.h:19
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
HLT_FULL_cff.delta_phi
delta_phi
Definition: HLT_FULL_cff.py:11672
ConversionLikelihoodCalculator::log_abs_delta_phi_
float log_abs_delta_phi_
Definition: ConversionLikelihoodCalculator.h:21
pi
const Double_t pi
Definition: trackSplitPlot.h:36
ConversionLikelihoodCalculator::log_abs_cot_theta_
float log_abs_cot_theta_
Definition: ConversionLikelihoodCalculator.h:20
ConversionLikelihoodCalculator.h
ConversionLikelihoodCalculator::setWeightsFile
void setWeightsFile(const char *weightsFile)
Definition: ConversionLikelihoodCalculator.cc:17
Conversion.h