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 
176  std::string variables_name_[Variable::VarNumber] = {
177  "CosThetaCS", "DeltaEta", "EtaMinus", "EtaPlus", "PhiCS", "PhiMinus", "PhiPlus", "Pt"};
178 
179  int variables_bins_number_[Variable::VarNumber]; // = {20, 20, 12, 12, 20, 16, 16, 100};
180  double variables_min_[Variable::VarNumber]; // = {-1, -4.8, -2.4, -2.4, -M_PI / 2, -M_PI, -M_PI, 0};
181  double variables_max_[Variable::VarNumber]; // = {+1, +4.8, +2.4, +2.4, +M_PI / 2, +M_PI, +M_PI, 100};
182 
187 };
188 
189 //
190 // member functions
191 //
192 
193 // ------------ method called for each event ------------
195  using namespace edm;
196 
198 
199  DiMuonValid::LV LV_mother(0., 0., 0., 0.);
200  for (reco::TrackCollection::const_iterator track1 = tC.begin(); track1 != tC.end(); track1++) {
201  DiMuonValid::LV LV_track1(track1->px(),
202  track1->py(),
203  track1->pz(),
204  sqrt((track1->p() * track1->p()) + mu_mass2_)); //old 106
205 
206  for (reco::TrackCollection::const_iterator track2 = track1 + 1; track2 != tC.end(); track2++) {
207  if (track1->charge() == track2->charge()) {
208  continue;
209  } // only reconstruct opposite charge pair
210 
211  DiMuonValid::LV LV_track2(
212  track2->px(), track2->py(), track2->pz(), sqrt((track2->p() * track2->p()) + mu_mass2_));
213 
214  LV_mother = LV_track1 + LV_track2;
215  double mother_mass = LV_mother.M();
216  th1f_mass->Fill(mother_mass);
217  double mother_pt = LV_mother.Pt();
218 
219  int charge1 = track1->charge();
220  double etaMu1 = track1->eta();
221  double phiMu1 = track1->phi();
222  double ptMu1 = track1->pt();
223 
224  int charge2 = track2->charge();
225  double etaMu2 = track2->eta();
226  double phiMu2 = track2->phi();
227  double ptMu2 = track2->pt();
228 
229  if (charge1 < 0) { // use Mu+ for charge1, Mu- for charge2
230  std::swap(charge1, charge2);
231  std::swap(etaMu1, etaMu2);
232  std::swap(phiMu1, phiMu2);
233  std::swap(ptMu1, ptMu2);
234  }
235 
236  const auto& tplus = track1->charge() > 0 ? track1 : track2;
237  const auto& tminus = track1->charge() < 0 ? track1 : track2;
238  TLorentzVector p4_tplus(tplus->px(), tplus->py(), tplus->pz(), sqrt((tplus->p() * tplus->p()) + mu_mass2_));
239  TLorentzVector p4_tminus(tminus->px(), tminus->py(), tminus->pz(), sqrt((tminus->p() * tminus->p()) + mu_mass2_));
240  std::pair<TLorentzVector, TLorentzVector> tktk_p4 = std::make_pair(p4_tplus, p4_tminus);
241  InvMassInEtaBins.fillTH1Plots(mother_mass, tktk_p4);
242  InvMassVsPhiPlusInEtaBins.fillTH2Plots(mother_mass, phiMu1, tktk_p4);
243  InvMassVsPhiMinusInEtaBins.fillTH2Plots(mother_mass, phiMu2, tktk_p4);
244 
245  //eta cut
246  if (etaMu1 < pair_etaminpos_ || etaMu1 > pair_etamaxpos_ || etaMu2 < pair_etaminneg_ ||
247  etaMu2 > pair_etamaxneg_) {
248  continue;
249  }
250 
251  double delta_eta = etaMu1 - etaMu2;
252 
253  double muplus = 1.0 / sqrt(2.0) * (LV_track1.E() + LV_track1.Z());
254  double muminus = 1.0 / sqrt(2.0) * (LV_track1.E() - LV_track1.Z());
255  double mubarplus = 1.0 / sqrt(2.0) * (LV_track2.E() + LV_track2.Z());
256  double mubarminus = 1.0 / sqrt(2.0) * (LV_track2.E() - LV_track2.Z());
257  //double costheta = 2.0 / Q.mag() / sqrt(pow(Q.mag(), 2) + pow(Q.Pt(), 2)) * (muplus * mubarminus - muminus * mubarplus);
258  double costhetaCS = 2.0 / LV_mother.mag() / sqrt(pow(LV_mother.mag(), 2) + pow(LV_mother.Pt(), 2)) *
259  (muplus * mubarminus - muminus * mubarplus);
260 
261  InvMassVsCosThetaCSInEtaBins.fillTH2Plots(mother_mass, costhetaCS, tktk_p4);
262 
263  DiMuonValid::LV Pbeam(0., 0., eBeam_, eBeam_);
264  auto R = Pbeam.Vect().Cross(LV_mother.Vect());
265  auto Runit = R.Unit();
266  auto Qt = LV_mother.Vect();
267  Qt.SetZ(0);
268  auto Qtunit = Qt.Unit();
269  DiMuonValid::LV D(LV_track1 - LV_track2);
270  auto Dt = D.Vect();
271  Dt.SetZ(0);
272  double tanphi =
273  sqrt(pow(LV_mother.mag(), 2) + pow(LV_mother.Pt(), 2)) / LV_mother.mag() * Dt.Dot(Runit) / Dt.Dot(Qtunit);
274  double phiCS = atan(tanphi);
275 
276  if (mother_mass > pair_mass_min_ && mother_mass < pair_mass_max_) {
277  th2d_mass_variables_[Variable::CosThetaCS]->Fill(mother_mass, costhetaCS, 1);
278  th2d_mass_variables_[Variable::DeltaEta]->Fill(mother_mass, delta_eta, 1);
279  th2d_mass_variables_[Variable::EtaMinus]->Fill(mother_mass, etaMu2, 1);
280  th2d_mass_variables_[Variable::EtaPlus]->Fill(mother_mass, etaMu1, 1);
281  th2d_mass_variables_[Variable::PhiCS]->Fill(mother_mass, phiCS, 1);
282  th2d_mass_variables_[Variable::PhiMinus]->Fill(mother_mass, phiMu2, 1);
283  th2d_mass_variables_[Variable::PhiPlus]->Fill(mother_mass, phiMu1, 1);
284  th2d_mass_variables_[Variable::Pt]->Fill(mother_mass, mother_pt, 1);
285  }
286  }
287  }
288 }
289 
290 // ------------ method called once each job just before starting event loop ------------
293  if (compressionSettings_ > 0) {
294  fs->file().SetCompressionSettings(compressionSettings_);
295  }
296 
297  th1f_mass = fs->make<TH1F>("hMass", "mass;m_{#mu#mu} [GeV];events", 200, 0., 200.);
298 
299  for (int i = 0; i < Variable::VarNumber; i++) {
300  std::string th2d_name = fmt::sprintf("th2d_mass_%s", variables_name_[i].c_str());
301  th2d_mass_variables_[i] = fs->make<TH2D>(th2d_name.c_str(),
302  th2d_name.c_str(),
307  variables_min_[i],
308  variables_max_[i]);
309  }
310 
311  // Z-> mm mass in eta bins
312  TFileDirectory dirResMassEta = fs->mkdir("TkTkMassInEtaBins");
313  InvMassInEtaBins.bookSet(dirResMassEta, th1f_mass);
314 
315  // Z->mm mass vs cos phi Collins-Soper in eta bins
316  TFileDirectory dirResMassVsCosThetaCSInEta = fs->mkdir("TkTkMassVsCosThetaCSInEtaBins");
317  InvMassVsCosThetaCSInEtaBins.bookSet(dirResMassVsCosThetaCSInEta, th2d_mass_variables_[Variable::CosThetaCS]);
318 
319  // Z->mm mass vs phi minus in eta bins
320  TFileDirectory dirResMassVsPhiMinusInEta = fs->mkdir("TkTkMassVsPhiMinusInEtaBins");
321  InvMassVsPhiMinusInEtaBins.bookSet(dirResMassVsPhiMinusInEta, th2d_mass_variables_[Variable::PhiMinus]);
322 
323  // Z->mm mass vs phi plus in eta bins
324  TFileDirectory dirResMassVsPhiPlusInEta = fs->mkdir("TkTkMassVsPhiPlusInEtaBins");
325  InvMassVsPhiPlusInEtaBins.bookSet(dirResMassVsPhiPlusInEta, th2d_mass_variables_[Variable::PhiPlus]);
326 }
327 
328 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
331  desc.setComment("Validates alignment payloads by evaluating bias in Z->mm mass distributions");
332  desc.addUntracked<int>("compressionSettings", -1);
333 
334  desc.add<double>("eBeam", 3500.)->setComment("beam energy in GeV");
335  desc.add<std::string>("TkTag", "ALCARECOTkAlZMuMu");
336 
337  desc.add<double>("Pair_mass_min", 60);
338  desc.add<double>("Pair_mass_max", 120);
339  desc.add<int>("Pair_mass_nbins", 120);
340  desc.add<double>("Pair_etaminpos", -2.4);
341  desc.add<double>("Pair_etamaxpos", 2.4);
342  desc.add<double>("Pair_etaminneg", -2.4);
343  desc.add<double>("Pair_etamaxneg", 2.4);
344 
345  desc.add<double>("Variable_CosThetaCS_xmin", -1.);
346  desc.add<double>("Variable_CosThetaCS_xmax", 1.);
347  desc.add<int>("Variable_CosThetaCS_nbins", 20);
348 
349  desc.add<double>("Variable_DeltaEta_xmin", -4.8);
350  desc.add<double>("Variable_DeltaEta_xmax", 4.8);
351  desc.add<int>("Variable_DeltaEta_nbins", 20);
352 
353  desc.add<double>("Variable_EtaMinus_xmin", -2.4);
354  desc.add<double>("Variable_EtaMinus_xmax", 2.4);
355  desc.add<int>("Variable_EtaMinus_nbins", 12);
356 
357  desc.add<double>("Variable_EtaPlus_xmin", -2.4);
358  desc.add<double>("Variable_EtaPlus_xmax", 2.4);
359  desc.add<int>("Variable_EtaPlus_nbins", 12);
360 
361  desc.add<double>("Variable_PhiCS_xmin", -M_PI / 2);
362  desc.add<double>("Variable_PhiCS_xmax", M_PI / 2);
363  desc.add<int>("Variable_PhiCS_nbins", 20);
364 
365  desc.add<double>("Variable_PhiMinus_xmin", -M_PI);
366  desc.add<double>("Variable_PhiMinus_xmax", M_PI);
367  desc.add<int>("Variable_PhiMinus_nbins", 16);
368 
369  desc.add<double>("Variable_PhiPlus_xmin", -M_PI);
370  desc.add<double>("Variable_PhiPlus_xmax", M_PI);
371  desc.add<int>("Variable_PhiPlus_nbins", 16);
372 
373  desc.add<double>("Variable_PairPt_xmin", 0.);
374  desc.add<double>("Variable_PairPt_xmax", 100.);
375  desc.add<int>("Variable_PairPt_nbins", 100);
376 
377  descriptions.addWithDefaultLabel(desc);
378 }
379 
380 //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
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
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29