CMS 3D CMS Logo

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