Rivet analyses referenceBABAR_2021_I1867843Dalitz plot analysis of $\eta_c\to K^+K^-\eta^\prime$, $\pi^+\pi^-\eta$ and $\pi^+\pi^-\eta^\prime$Experiment: BABAR (PEP-II) Inspire ID: 1867843 Status: VALIDATED NOHEPDATA Authors:
Beam energies: ANY Run details:
Measurement of the mass distributions in the decays $\eta_c\to K^+K^-\eta^\prime$, $\pi^+\pi^-\eta$ and $\pi^+\pi^-\eta^\prime$ by BaBar. The data were read from the plots in the paper and therefore for some points the error bars are the size of the point. Also the sideband background from the plots has been subtracted. It is also not clear that any resolution effects have been unfolded. Source code: BABAR_2021_I1867843.cc 1// -*- C++ -*-
2#include "Rivet/Analysis.hh"
3#include "Rivet/Projections/DecayedParticles.hh"
4#include "Rivet/Projections/UnstableParticles.hh"
5
6namespace Rivet {
7
8
9 /// @brief eta_c Dalitz decays
10 class BABAR_2021_I1867843 : public Analysis {
11 public:
12
13 /// Constructor
14 RIVET_DEFAULT_ANALYSIS_CTOR(BABAR_2021_I1867843);
15
16
17 /// @name Analysis methods
18 /// @{
19
20 /// Book histograms and initialise projections before the run
21 void init() {
22 // Initialise and register projections
23 UnstableParticles ufs = UnstableParticles(Cuts::pid== 441);
24 declare(ufs, "UFS");
25 DecayedParticles ETAC(ufs);
26 ETAC.addStable(PID::PI0);
27 ETAC.addStable(PID::K0S);
28 ETAC.addStable(PID::ETA);
29 ETAC.addStable(PID::ETAPRIME);
30 declare(ETAC,"ETAC");
31 // histograms
32 book(_h_KK ,1,1,1);
33 book(_h_etaPK ,1,1,2);
34 book(_h_pipi1 ,2,1,1);
35 book(_h_etaPpi,2,1,2);
36 book(_h_pipi2 ,3,1,1);
37 book(_h_etapi ,3,1,2);
38 book(_dalitz1, "dalitz1",50,2.,6.5,50,2.,6.5);
39 book(_dalitz2, "dalitz2",50,0.,9. ,50,0.,9. );
40 book(_dalitz3, "dalitz3",50,0.,8. ,50,0.,8. );
41 }
42
43 /// Perform the per-event analysis
44 void analyze(const Event& event) {
45 static const map<PdgId,unsigned int> & mode1 = { { 321,1}, {-321,1}, { 331,1}};
46 static const map<PdgId,unsigned int> & mode2 = { { 211,1}, {-211,1}, { 331,1}};
47 static const map<PdgId,unsigned int> & mode3 = { { 211,1}, {-211,1}, { 221,1}};
48 DecayedParticles ETAC = apply<DecayedParticles>(event, "ETAC");
49 // loop over particles
50 for(unsigned int ix=0;ix<ETAC.decaying().size();++ix) {
51 //K+ K- eta'
52 if (ETAC.modeMatches(ix,3,mode1)&&
53 ETAC.decaying()[ix].mass()>2.93 && ETAC.decaying()[ix].mass()<3.03) {
54 const Particle & Kp = ETAC.decayProducts()[ix].at( 321)[0];
55 const Particle & Km = ETAC.decayProducts()[ix].at(-321)[0];
56 const Particle & etaP = ETAC.decayProducts()[ix].at( 331)[0];
57 double mplus = (Kp.momentum()+etaP.momentum()).mass2();
58 double mminus = (Km.momentum()+etaP.momentum()).mass2();
59 double mKK = (Kp.momentum()+Km .momentum()).mass2();
60 _h_KK ->fill(sqrt(mKK));
61 _h_etaPK->fill(sqrt(mplus));
62 _h_etaPK->fill(sqrt(mminus));
63 _dalitz1->fill(mplus,mminus);
64 }
65 // pi+ pi- eta'
66 else if (ETAC.modeMatches(ix,3,mode2)&&
67 ETAC.decaying()[ix].mass()>2.93 && ETAC.decaying()[ix].mass()<3.03) {
68 const Particle & pip = ETAC.decayProducts()[ix].at( 211)[0];
69 const Particle & pim = ETAC.decayProducts()[ix].at(-211)[0];
70 const Particle & etaP = ETAC.decayProducts()[ix].at( 331)[0];
71 double mplus = (pip.momentum()+etaP.momentum()).mass2();
72 double mminus = (pim.momentum()+etaP.momentum()).mass2();
73 double mpipi = (pip.momentum()+pim .momentum()).mass2();
74 _h_pipi1 ->fill(sqrt(mpipi));
75 _h_etaPpi->fill(sqrt(mplus));
76 _h_etaPpi->fill(sqrt(mminus));
77 _dalitz2->fill(mplus,mminus);
78 }
79 // pi+ pi- eta
80 else if (ETAC.modeMatches(ix,3,mode3)&&
81 ETAC.decaying()[ix].mass()>2.92 && ETAC.decaying()[ix].mass()<3.02) {
82 const Particle & pip = ETAC.decayProducts()[ix].at( 211)[0];
83 const Particle & pim = ETAC.decayProducts()[ix].at(-211)[0];
84 const Particle & eta = ETAC.decayProducts()[ix].at( 221)[0];
85 double mplus = (pip.momentum()+eta.momentum()).mass2();
86 double mminus = (pim.momentum()+eta.momentum()).mass2();
87 double mpipi = (pip.momentum()+pim .momentum()).mass2();
88 _h_pipi2->fill(sqrt(mpipi));
89 _h_etapi->fill(sqrt(mplus));
90 _h_etapi->fill(sqrt(mminus));
91 _dalitz3->fill(mplus,mminus);
92 }
93 }
94 }
95
96
97 /// Normalise histograms etc., after the run
98 void finalize() {
99 normalize(_h_KK );
100 normalize(_h_etaPK );
101 normalize(_h_pipi1 );
102 normalize(_h_etaPpi);
103 normalize(_h_pipi2 );
104 normalize(_h_etapi );
105 normalize(_dalitz1 );
106 normalize(_dalitz2 );
107 normalize(_dalitz3 );
108 }
109
110 /// @}
111
112
113 /// @name Histograms
114 /// @{
115 Histo1DPtr _h_KK,_h_etaPK;
116 Histo1DPtr _h_pipi1,_h_etaPpi;
117 Histo1DPtr _h_pipi2,_h_etapi;
118 Histo2DPtr _dalitz1,_dalitz2,_dalitz3;
119 /// @}
120
121
122 };
123
124
125 RIVET_DECLARE_PLUGIN(BABAR_2021_I1867843);
126
127}
|