CMS 3D CMS Logo

DiMuonValidation.cc
Go to the documentation of this file.
1 // system includes
2 #include <memory>
3 #include <cmath>
4 #include <fmt/printf.h>
5 
6 // user include files
23 
24 // ROOT includes
25 #include "TH1D.h"
26 #include "TH2D.h"
27 #include "TH3D.h"
28 
29 namespace DiMuonValid {
31 }
32 //
33 // class declaration
34 //
35 class DiMuonValidation : public edm::one::EDAnalyzer<edm::one::SharedResources> {
36 public:
38  : eBeam_(pset.getParameter<double>("eBeam")),
39  compressionSettings_(pset.getUntrackedParameter<int>("compressionSettings", -1)),
40  TkTag_(pset.getParameter<std::string>("TkTag")),
41  pair_mass_min_(pset.getParameter<double>("Pair_mass_min")),
42  pair_mass_max_(pset.getParameter<double>("Pair_mass_max")),
43  pair_mass_nbins_(pset.getParameter<int>("Pair_mass_nbins")),
44  pair_etaminpos_(pset.getParameter<double>("Pair_etaminpos")),
45  pair_etamaxpos_(pset.getParameter<double>("Pair_etamaxpos")),
46  pair_etaminneg_(pset.getParameter<double>("Pair_etaminneg")),
47  pair_etamaxneg_(pset.getParameter<double>("Pair_etamaxneg")),
48  variable_CosThetaCS_xmin_(pset.getParameter<double>("Variable_CosThetaCS_xmin")),
49  variable_CosThetaCS_xmax_(pset.getParameter<double>("Variable_CosThetaCS_xmax")),
50  variable_CosThetaCS_nbins_(pset.getParameter<int>("Variable_CosThetaCS_nbins")),
51  variable_DeltaEta_xmin_(pset.getParameter<double>("Variable_DeltaEta_xmin")),
52  variable_DeltaEta_xmax_(pset.getParameter<double>("Variable_DeltaEta_xmax")),
53  variable_DeltaEta_nbins_(pset.getParameter<int>("Variable_DeltaEta_nbins")),
54  variable_EtaMinus_xmin_(pset.getParameter<double>("Variable_EtaMinus_xmin")),
55  variable_EtaMinus_xmax_(pset.getParameter<double>("Variable_EtaMinus_xmax")),
56  variable_EtaMinus_nbins_(pset.getParameter<int>("Variable_EtaMinus_nbins")),
57  variable_EtaPlus_xmin_(pset.getParameter<double>("Variable_EtaPlus_xmin")),
58  variable_EtaPlus_xmax_(pset.getParameter<double>("Variable_EtaPlus_xmax")),
59  variable_EtaPlus_nbins_(pset.getParameter<int>("Variable_EtaPlus_nbins")),
60  variable_PhiCS_xmin_(pset.getParameter<double>("Variable_PhiCS_xmin")),
61  variable_PhiCS_xmax_(pset.getParameter<double>("Variable_PhiCS_xmax")),
62  variable_PhiCS_nbins_(pset.getParameter<int>("Variable_PhiCS_nbins")),
63  variable_PhiMinus_xmin_(pset.getParameter<double>("Variable_PhiMinus_xmin")),
64  variable_PhiMinus_xmax_(pset.getParameter<double>("Variable_PhiMinus_xmax")),
65  variable_PhiMinus_nbins_(pset.getParameter<int>("Variable_PhiMinus_nbins")),
66  variable_PhiPlus_xmin_(pset.getParameter<double>("Variable_PhiPlus_xmin")),
67  variable_PhiPlus_xmax_(pset.getParameter<double>("Variable_PhiPlus_xmax")),
68  variable_PhiPlus_nbins_(pset.getParameter<int>("Variable_PhiPlus_nbins")),
69  variable_PairPt_xmin_(pset.getParameter<double>("Variable_PairPt_xmin")),
70  variable_PairPt_xmax_(pset.getParameter<double>("Variable_PairPt_xmax")),
71  variable_PairPt_nbins_(pset.getParameter<int>("Variable_PairPt_nbins")) {
72  usesResource(TFileService::kSharedResource);
73  theTrackCollectionToken_ = consumes<reco::TrackCollection>(TkTag_);
74 
75  variables_min_[Variable::CosThetaCS] = variable_CosThetaCS_xmin_;
77  variables_min_[Variable::EtaMinus] = variable_EtaMinus_xmin_;
78  variables_min_[Variable::EtaPlus] = variable_EtaPlus_xmin_;
79  variables_min_[Variable::PhiCS] = variable_PhiCS_xmin_;
80  variables_min_[Variable::PhiMinus] = variable_PhiMinus_xmin_;
81  variables_min_[Variable::PhiPlus] = variable_PhiPlus_xmin_;
82  variables_min_[Variable::Pt] = variable_PairPt_xmin_;
83 
84  variables_max_[Variable::CosThetaCS] = variable_CosThetaCS_xmax_;
86  variables_max_[Variable::EtaMinus] = variable_EtaMinus_xmax_;
87  variables_max_[Variable::EtaPlus] = variable_EtaPlus_xmax_;
88  variables_max_[Variable::PhiCS] = variable_PhiCS_xmax_;
89  variables_max_[Variable::PhiMinus] = variable_PhiMinus_xmax_;
90  variables_max_[Variable::PhiPlus] = variable_PhiPlus_xmax_;
91  variables_max_[Variable::Pt] = variable_PairPt_xmax_;
92 
93  variables_bins_number_[Variable::CosThetaCS] = variable_CosThetaCS_nbins_;
101  }
102 
103  ~DiMuonValidation() override = default;
104 
105  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
106 
107  // enumeration to contain the variables to plot
108  enum Variable {
110  DeltaEta = 1,
111  EtaMinus = 2,
112  EtaPlus = 3,
113  PhiCS = 4,
114  PhiMinus = 5,
115  PhiPlus = 6,
116  Pt = 7,
118  };
119 
120  static constexpr double mu_mass2_ = 0.105658 * 0.105658; //The invariant mass of muon 105.658MeV
121 
122 private:
123  void beginJob() override;
124  void analyze(const edm::Event&, const edm::EventSetup&) override;
125 
126  // ----------member data ---------------------------
127  float eBeam_;
130 
138 
142 
146 
150 
154 
158 
162 
166 
170 
172 
173  TH1F* th1f_mass;
174  TH2D* th2d_mass_variables_[Variable::VarNumber]; // actual histograms
175  TH3D* th3d_mass_vs_eta_phi_plus_; // 3D histogram for scatter plot vs eta / phi (mu+)
176  TH3D* th3d_mass_vs_eta_phi_minus_; // 3D histogram for scatter plot vs eta / phi (mu-)
177 
178  std::string variables_name_[Variable::VarNumber] = {
179  "CosThetaCS", "DeltaEta", "EtaMinus", "EtaPlus", "PhiCS", "PhiMinus", "PhiPlus", "Pt"};
180 
181  std::string variables_title_[Variable::VarNumber] = {"Cos#theta_{CS}",
182  "#Delta#eta(#mu^{-},#mu^{+})",
183  "#mu^{-} #eta",
184  "#mu^{+} #eta",
185  "#phi_{CS} [rad]",
186  "#mu^{-} #phi [rad]",
187  "#mu^{+} #phi [rad]",
188  "p_{T} [GeV]"};
189 
190  int variables_bins_number_[Variable::VarNumber]; // = {20, 20, 12, 12, 20, 16, 16, 100};
191  double variables_min_[Variable::VarNumber]; // = {-1, -4.8, -2.4, -2.4, -M_PI / 2, -M_PI, -M_PI, 0};
192  double variables_max_[Variable::VarNumber]; // = {+1, +4.8, +2.4, +2.4, +M_PI / 2, +M_PI, +M_PI, 100};
193 
198 };
199 
200 //
201 // member functions
202 //
203 
204 // ------------ method called for each event ------------
206  using namespace edm;
207 
209 
210  DiMuonValid::LV LV_mother(0., 0., 0., 0.);
211  for (reco::TrackCollection::const_iterator track1 = tC.begin(); track1 != tC.end(); track1++) {
212  DiMuonValid::LV LV_track1(track1->px(),
213  track1->py(),
214  track1->pz(),
215  sqrt((track1->p() * track1->p()) + mu_mass2_)); //old 106
216 
217  for (reco::TrackCollection::const_iterator track2 = track1 + 1; track2 != tC.end(); track2++) {
218  if (track1->charge() == track2->charge()) {
219  continue;
220  } // only reconstruct opposite charge pair
221 
222  DiMuonValid::LV LV_track2(
223  track2->px(), track2->py(), track2->pz(), sqrt((track2->p() * track2->p()) + mu_mass2_));
224 
225  LV_mother = LV_track1 + LV_track2;
226  double mother_mass = LV_mother.M();
227  th1f_mass->Fill(mother_mass);
228  double mother_pt = LV_mother.Pt();
229 
230  int charge1 = track1->charge();
231  double etaMu1 = track1->eta();
232  double phiMu1 = track1->phi();
233  double ptMu1 = track1->pt();
234 
235  int charge2 = track2->charge();
236  double etaMu2 = track2->eta();
237  double phiMu2 = track2->phi();
238  double ptMu2 = track2->pt();
239 
240  if (charge1 < 0) { // use Mu+ for charge1, Mu- for charge2
241  std::swap(charge1, charge2);
242  std::swap(etaMu1, etaMu2);
243  std::swap(phiMu1, phiMu2);
244  std::swap(ptMu1, ptMu2);
245  }
246 
247  const auto& tplus = track1->charge() > 0 ? track1 : track2;
248  const auto& tminus = track1->charge() < 0 ? track1 : track2;
249  TLorentzVector p4_tplus(tplus->px(), tplus->py(), tplus->pz(), sqrt((tplus->p() * tplus->p()) + mu_mass2_));
250  TLorentzVector p4_tminus(tminus->px(), tminus->py(), tminus->pz(), sqrt((tminus->p() * tminus->p()) + mu_mass2_));
251  std::pair<TLorentzVector, TLorentzVector> tktk_p4 = std::make_pair(p4_tplus, p4_tminus);
252  InvMassInEtaBins.fillTH1Plots(mother_mass, tktk_p4);
253  InvMassVsPhiPlusInEtaBins.fillTH2Plots(mother_mass, phiMu1, tktk_p4);
254  InvMassVsPhiMinusInEtaBins.fillTH2Plots(mother_mass, phiMu2, tktk_p4);
255 
256  //eta cut
257  if (etaMu1 < pair_etaminpos_ || etaMu1 > pair_etamaxpos_ || etaMu2 < pair_etaminneg_ ||
258  etaMu2 > pair_etamaxneg_) {
259  continue;
260  }
261 
262  double delta_eta = etaMu1 - etaMu2;
263 
264  double muplus = 1.0 / sqrt(2.0) * (LV_track1.E() + LV_track1.Z());
265  double muminus = 1.0 / sqrt(2.0) * (LV_track1.E() - LV_track1.Z());
266  double mubarplus = 1.0 / sqrt(2.0) * (LV_track2.E() + LV_track2.Z());
267  double mubarminus = 1.0 / sqrt(2.0) * (LV_track2.E() - LV_track2.Z());
268  //double costheta = 2.0 / Q.mag() / sqrt(pow(Q.mag(), 2) + pow(Q.Pt(), 2)) * (muplus * mubarminus - muminus * mubarplus);
269  double costhetaCS = 2.0 / LV_mother.mag() / sqrt(pow(LV_mother.mag(), 2) + pow(LV_mother.Pt(), 2)) *
270  (muplus * mubarminus - muminus * mubarplus);
271 
272  InvMassVsCosThetaCSInEtaBins.fillTH2Plots(mother_mass, costhetaCS, tktk_p4);
273 
274  DiMuonValid::LV Pbeam(0., 0., eBeam_, eBeam_);
275  auto R = Pbeam.Vect().Cross(LV_mother.Vect());
276  auto Runit = R.Unit();
277  auto Qt = LV_mother.Vect();
278  Qt.SetZ(0);
279  auto Qtunit = Qt.Unit();
280  DiMuonValid::LV D(LV_track1 - LV_track2);
281  auto Dt = D.Vect();
282  Dt.SetZ(0);
283  double tanphi =
284  sqrt(pow(LV_mother.mag(), 2) + pow(LV_mother.Pt(), 2)) / LV_mother.mag() * Dt.Dot(Runit) / Dt.Dot(Qtunit);
285  double phiCS = atan(tanphi);
286 
287  if (mother_mass > pair_mass_min_ && mother_mass < pair_mass_max_) {
288  th2d_mass_variables_[Variable::CosThetaCS]->Fill(mother_mass, costhetaCS, 1);
289  th2d_mass_variables_[Variable::DeltaEta]->Fill(mother_mass, delta_eta, 1);
290  th2d_mass_variables_[Variable::EtaMinus]->Fill(mother_mass, etaMu2, 1);
291  th2d_mass_variables_[Variable::EtaPlus]->Fill(mother_mass, etaMu1, 1);
292  th2d_mass_variables_[Variable::PhiCS]->Fill(mother_mass, phiCS, 1);
293  th2d_mass_variables_[Variable::PhiMinus]->Fill(mother_mass, phiMu2, 1);
294  th2d_mass_variables_[Variable::PhiPlus]->Fill(mother_mass, phiMu1, 1);
295  th2d_mass_variables_[Variable::Pt]->Fill(mother_mass, mother_pt, 1);
296 
297  th3d_mass_vs_eta_phi_plus_->Fill(mother_mass, etaMu1, phiMu1);
298  th3d_mass_vs_eta_phi_minus_->Fill(mother_mass, etaMu2, phiMu2);
299  }
300  }
301  }
302 }
303 
304 // ------------ method called once each job just before starting event loop ------------
307  if (compressionSettings_ > 0) {
308  fs->file().SetCompressionSettings(compressionSettings_);
309  }
310 
311  th1f_mass = fs->make<TH1F>("hMass", "mass;m_{#mu#mu} [GeV];events", 200, 0., 200.);
312 
313  for (int i = 0; i < Variable::VarNumber; i++) {
314  std::string th2d_name = fmt::sprintf("th2d_mass_%s", variables_name_[i].c_str());
316  fs->make<TH2D>(th2d_name.c_str(),
317  fmt::format("{};M_{{#mu^{{-}}#mu^{{+}}}} [GeV];{}", th2d_name, variables_title_[i]).c_str(),
322  variables_min_[i],
323  variables_max_[i]);
324  }
325 
326  // 3D histogram for eta/phi map (mu+)
328  fs->make<TH3D>("th3d_mass_vs_eta_phi_plus",
329  "th3d_mass_vs_eta_phi_plus;M_{#mu^{-}#mu^{+}} [GeV];#mu^{+} #eta;#mu^{+} #phi [rad]",
333  variables_bins_number_[Variable::EtaPlus],
334  variables_min_[Variable::EtaPlus],
335  variables_max_[Variable::EtaPlus],
336  variables_bins_number_[Variable::PhiPlus],
337  variables_min_[Variable::PhiPlus],
338  variables_max_[Variable::PhiPlus]);
339 
340  // 3D histogram for eta/phi map (mu+)
342  fs->make<TH3D>("th3d_mass_vs_eta_phi_minus",
343  "th3d_mass_vs_eta_phi_minus;M_{#mu^{-}#mu^{+}} [GeV];#mu^{-} #eta;#mu^{-} #phi [rad]",
347  variables_bins_number_[Variable::EtaMinus],
348  variables_min_[Variable::EtaMinus],
349  variables_max_[Variable::EtaMinus],
350  variables_bins_number_[Variable::PhiMinus],
351  variables_min_[Variable::PhiMinus],
352  variables_max_[Variable::PhiMinus]);
353 
354  // Z-> mm mass in eta bins
355  TFileDirectory dirResMassEta = fs->mkdir("TkTkMassInEtaBins");
356  InvMassInEtaBins.bookSet(dirResMassEta, th1f_mass);
357 
358  // Z->mm mass vs cos phi Collins-Soper in eta bins
359  TFileDirectory dirResMassVsCosThetaCSInEta = fs->mkdir("TkTkMassVsCosThetaCSInEtaBins");
360  InvMassVsCosThetaCSInEtaBins.bookSet(dirResMassVsCosThetaCSInEta, th2d_mass_variables_[Variable::CosThetaCS]);
361 
362  // Z->mm mass vs phi minus in eta bins
363  TFileDirectory dirResMassVsPhiMinusInEta = fs->mkdir("TkTkMassVsPhiMinusInEtaBins");
364  InvMassVsPhiMinusInEtaBins.bookSet(dirResMassVsPhiMinusInEta, th2d_mass_variables_[Variable::PhiMinus]);
365 
366  // Z->mm mass vs phi plus in eta bins
367  TFileDirectory dirResMassVsPhiPlusInEta = fs->mkdir("TkTkMassVsPhiPlusInEtaBins");
368  InvMassVsPhiPlusInEtaBins.bookSet(dirResMassVsPhiPlusInEta, th2d_mass_variables_[Variable::PhiPlus]);
369 }
370 
371 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
374  desc.setComment("Validates alignment payloads by evaluating bias in Z->mm mass distributions");
375  desc.addUntracked<int>("compressionSettings", -1);
376 
377  desc.add<double>("eBeam", 3500.)->setComment("beam energy in GeV");
378  desc.add<std::string>("TkTag", "ALCARECOTkAlZMuMu");
379 
380  desc.add<double>("Pair_mass_min", 60.);
381  desc.add<double>("Pair_mass_max", 120.);
382  desc.add<int>("Pair_mass_nbins", 120);
383  desc.add<double>("Pair_etaminpos", -2.4);
384  desc.add<double>("Pair_etamaxpos", 2.4);
385  desc.add<double>("Pair_etaminneg", -2.4);
386  desc.add<double>("Pair_etamaxneg", 2.4);
387 
388  desc.add<double>("Variable_CosThetaCS_xmin", -1.);
389  desc.add<double>("Variable_CosThetaCS_xmax", 1.);
390  desc.add<int>("Variable_CosThetaCS_nbins", 20);
391 
392  desc.add<double>("Variable_DeltaEta_xmin", -4.8);
393  desc.add<double>("Variable_DeltaEta_xmax", 4.8);
394  desc.add<int>("Variable_DeltaEta_nbins", 20);
395 
396  desc.add<double>("Variable_EtaMinus_xmin", -2.4);
397  desc.add<double>("Variable_EtaMinus_xmax", 2.4);
398  desc.add<int>("Variable_EtaMinus_nbins", 12);
399 
400  desc.add<double>("Variable_EtaPlus_xmin", -2.4);
401  desc.add<double>("Variable_EtaPlus_xmax", 2.4);
402  desc.add<int>("Variable_EtaPlus_nbins", 12);
403 
404  desc.add<double>("Variable_PhiCS_xmin", -M_PI / 2);
405  desc.add<double>("Variable_PhiCS_xmax", M_PI / 2);
406  desc.add<int>("Variable_PhiCS_nbins", 20);
407 
408  desc.add<double>("Variable_PhiMinus_xmin", -M_PI);
409  desc.add<double>("Variable_PhiMinus_xmax", M_PI);
410  desc.add<int>("Variable_PhiMinus_nbins", 16);
411 
412  desc.add<double>("Variable_PhiPlus_xmin", -M_PI);
413  desc.add<double>("Variable_PhiPlus_xmax", M_PI);
414  desc.add<int>("Variable_PhiPlus_nbins", 16);
415 
416  desc.add<double>("Variable_PairPt_xmin", 0.);
417  desc.add<double>("Variable_PairPt_xmax", 100.);
418  desc.add<int>("Variable_PairPt_nbins", 100);
419 
420  descriptions.addWithDefaultLabel(desc);
421 }
422 
423 //define this as a plug-in
static const std::string kSharedResource
Definition: TFileService.h:76
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
DiLeptonHelp::PlotsVsDiLeptonRegion InvMassVsPhiMinusInEtaBins
DiMuonValidation(const edm::ParameterSet &pset)
static constexpr double mu_mass2_
void fillTH1Plots(const float val, const std::pair< TLorentzVector, TLorentzVector > &momenta)
double variable_CosThetaCS_xmax_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
DiLeptonHelp::PlotsVsDiLeptonRegion InvMassVsCosThetaCSInEtaBins
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:112
double variable_CosThetaCS_xmin_
DiLeptonHelp::PlotsVsDiLeptonRegion InvMassInEtaBins
int iEvent
Definition: GenABIO.cc:224
TH2D * th2d_mass_variables_[Variable::VarNumber]
void bookSet(const TFileDirectory &fs, const TH1 *histo)
T sqrt(T t)
Definition: SSEVec.h:19
std::string variables_name_[Variable::VarNumber]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
#define M_PI
TH3D * th3d_mass_vs_eta_phi_minus_
TH3D * th3d_mass_vs_eta_phi_plus_
void beginJob() override
double variables_min_[Variable::VarNumber]
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:141
DiLeptonHelp::PlotsVsDiLeptonRegion InvMassVsPhiPlusInEtaBins
double variables_max_[Variable::VarNumber]
HLT enums.
void analyze(const edm::Event &, const edm::EventSetup &) override
int variables_bins_number_[Variable::VarNumber]
reco::Particle::LorentzVector LV
edm::EDGetTokenT< reco::TrackCollection > theTrackCollectionToken_
~DiMuonValidation() override=default
void fillTH2Plots(const float valX, const float valY, const std::pair< TLorentzVector, TLorentzVector > &momenta)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
std::string variables_title_[Variable::VarNumber]
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29