rivet is hosted by Hepforge, IPPP Durham

Rivet analyses reference

BABAR_2010_I853279

Dalitz plot analysis of $D^0\to K^0_S\pi^+\pi^-$ and $D^0\to K^0_SK^+K^-$
Experiment: BABAR (PEP-II)
Inspire ID: 853279
Status: VALIDATED NOHEPDATA
Authors:
  • Peter Richardson
References:
  • Phys.Rev.Lett. 105 (2010) 081803
Beams: * *
Beam energies: ANY
Run details:
  • Any process producing D0

Measurement of Kinematic distributions in the decays $D^0\to K^0_S\pi^+\pi^-$ and $D^0\to K^0_SK^+K^-$ by BaBar. The data were read from the plots in the paper. Resolution/acceptance effects have been not unfolded. Given the agreement with the model in the paper this analysis should only be used for qualitative studies.

Source code: BABAR_2010_I853279.cc
  1// -*- C++ -*-
  2#include "Rivet/Analysis.hh"
  3#include "Rivet/Projections/UnstableParticles.hh"
  4#include "Rivet/Projections/DecayedParticles.hh"
  5
  6namespace Rivet {
  7
  8
  9  /// @brief D0 -> KS0 pi+pi- and K+K-
 10  class BABAR_2010_I853279 : public Analysis {
 11  public:
 12
 13    /// Constructor
 14    RIVET_DEFAULT_ANALYSIS_CTOR(BABAR_2010_I853279);
 15
 16
 17    /// @name Analysis methods
 18    /// @{
 19
 20    /// Book histograms and initialise projections before the run
 21    void init() {
 22
 23      // Initialise and register projections
 24      UnstableParticles ufs = UnstableParticles(Cuts::abspid==421);
 25      declare(ufs, "UFS");
 26      DecayedParticles D0(ufs);
 27      D0.addStable(PID::PI0);
 28      D0.addStable(PID::K0S);
 29      D0.addStable(PID::ETA);
 30      D0.addStable(PID::ETAPRIME);
 31      declare(D0, "D0");
 32      // Histograms
 33      book(_h_Kpim,1,1,1);
 34      book(_h_Kpip,1,1,2);
 35      book(_h_pipi,1,1,3);
 36      book(_dalitz[0], "dalitz1",50,0.3,3.2,50,0.3,3.2);
 37      book(_h_K0Km,1,1,4);
 38      book(_h_K0Kp,1,1,5);
 39      book(_h_KpKm,1,1,6);
 40      book(_dalitz[1], "dalitz2",50,0.9,1.9,50,0.9,1.9);
 41    }
 42
 43    /// Perform the per-event analysis
 44    void analyze(const Event& event) {
 45      // define the decay mode
 46      static const map<PdgId,unsigned int> & modePi  = { { 310,1}, { 211,1},{-211,1}};
 47      static const map<PdgId,unsigned int> & modeK   = { { 310,1}, { 321,1},{-321,1}};
 48      DecayedParticles D0 = apply<DecayedParticles>(event, "D0");
 49      // loop over particles
 50      for(unsigned int ix=0;ix<D0.decaying().size();++ix) {
 51	int sign = D0.decaying()[ix].pid()/421;
 52	// KS0 pi+pi-
 53	if (D0.modeMatches(ix,3,modePi) ) {
 54	  const Particles & pip= D0.decayProducts()[ix].at( sign*211);
 55	  const Particles & pim= D0.decayProducts()[ix].at(-sign*211);
 56	  const Particles & K0 = D0.decayProducts()[ix].at( 310);
 57	  double mminus = (pim[0].momentum()+K0[0].momentum() ).mass2();
 58	  double mplus  = (pip[0].momentum()+K0[0].momentum() ).mass2();
 59	  double mpipi  = (pip[0].momentum()+pim[0].momentum()).mass2();
 60	  _h_Kpip->fill(mplus);
 61	  _h_Kpim->fill(mminus);
 62	  _h_pipi->fill(mpipi);
 63	  _dalitz[0]->fill(mminus,mplus); 
 64	}
 65	else if (D0.modeMatches(ix,3,modeK) ) {
 66	  const Particles & Kp = D0.decayProducts()[ix].at( sign*321);
 67	  const Particles & Km = D0.decayProducts()[ix].at(-sign*321);
 68	  const Particles & K0 = D0.decayProducts()[ix].at( 310);
 69	  double mminus = (Km[0].momentum()+K0[0].momentum() ).mass2();
 70	  double mplus  = (Kp[0].momentum()+K0[0].momentum() ).mass2();
 71	  double mKK  = (Kp[0].momentum()+Km[0].momentum()).mass2();
 72	  _h_K0Kp->fill(mplus);
 73	  _h_K0Km->fill(mminus);
 74	  _h_KpKm->fill(mKK);
 75	}
 76      }
 77    }
 78
 79
 80    /// Normalise histograms etc., after the run
 81    void finalize() {
 82      normalize(_h_Kpim);
 83      normalize(_h_Kpip);
 84      normalize(_h_pipi);
 85      normalize(_dalitz[0]);
 86      normalize(_h_K0Km);
 87      normalize(_h_K0Kp);
 88      normalize(_h_KpKm);
 89      normalize(_dalitz[1]);
 90    }
 91
 92    /// @}
 93
 94
 95    /// @name Histograms
 96    /// @{
 97    Histo1DPtr _h_Kpim,_h_pipi,_h_Kpip;
 98    Histo1DPtr _h_K0Km,_h_KpKm,_h_K0Kp;
 99    Histo2DPtr _dalitz[2];
100    /// @}
101
102
103  };
104
105
106  RIVET_DECLARE_PLUGIN(BABAR_2010_I853279);
107
108}