Rivet Analyses Reference

H1_1999_I504022

Forward pi0 meson production at HERA
Experiment: H1 (HERA)
Inspire ID: 504022
Status: VALIDATED
Authors:
  • Keila Moral Figueroa
  • Hannes Jung
References:
  • Phys.Lett.B462:440-452,1999
  • DOI: 10.1016/S0370-2693(99)00906-5
  • arXiv: hep-ex/9907030
Beams: p+ e-, e- p+
Beam energies: (820.0, 27.5); (27.5, 820.0) GeV
Run details:
  • $\pi^0$ production in the forward region detected by H1 in DIS ep scattering events at HERA. Cuts:0.1 $< y <$ 0.6, 2 $< Q^2 <$ 70 $GeV^2$, $\pi^0$ transverse momentum in the hadronic CMS. $p^*_T$ $ > $ 2.5 GeV, 5$^o< \theta < 25^o$, where $\theta$ is the polar angle of the produced $\pi^0$,$x_{\pi^0} > 0.001$ where $ x_{\pi^0} = E_{\pi^0}/ E_{proton}$ and $E_{proton} = 820GeV$ is the proton beam energy. RAPGAP 26.7 configuration parameters: 10^6 events.

High transverse momentum $\pi^0$ mesons have been measured with the H1 detector at HERA in deep-inelastic ep scattering events at low Bjorken-x, down to $x \approx 4 \cdot 10^{-5}$ .The measurement is performed in a region of small angles with respect to the proton remnant in the laboratory frame of reference, namely the forward region, and corresponds to central rapidity in the centre of mass system of the virtual photon and proton. This region is expected to be particularly sensitive to QCD effects in hadronic final states. Differential cross-sections for inclusive $\pi^0$ meson production are presented as a function of Bjorken-x and the four-momentum transfer $Q^2$, and as a function of transverse momentum and pseudorapidity.

Source code: H1_1999_I504022.cc
  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/UnstableParticles.hh"

#include "Rivet/Projections/DISKinematics.hh"

namespace Rivet {


  /// @brief Forward pi0 meson production at HERA (H1)
  class H1_1999_I504022 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(H1_1999_I504022);


    /// @name Analysis methods
    ///@{

    /// Book histograms and initialise projections before the run
    void init() {

      // Initialise and register projections

      // The basic final-state projection:
 
      declare(FinalState(Cuts::abseta < 7 ), "FS");
      declare(DISKinematics(), "Kinematics");
      declare(UnstableParticles(), "UFS");


      // Book histograms
      // take binning from reference data using HEPData ID (digits in "d01-x01-y01" etc.)
      
      book(_h["p_T>2.5&x"], 1, 1, 1);
      book(_h["Q:2.0-4.5-eta"], 2, 1, 1);
      book(_h["Q:2.0-4.5-p_T"], 3, 1, 1);
      book(_h["Q:4.5-15.0-x"], 4, 1, 1);
      book(_h["Q:4.5-15.0-eta"], 5, 1, 1);
      book(_h["Q:4.5-15.0-p_T"], 6, 1, 1);
      book(_h["Q:15.0-70.0-x"], 7, 1, 1);
      book(_h["Q:15.0-70.0-eta"], 8, 1, 1);
      book(_h["Q:15.0-70.0-p_T"], 9, 1, 1);
      book(_h["p_T>2.5&Q"], 10, 1, 1);
      book(_h["p_T>3.5&x"], 11, 1, 1);
      book(_h["p_T>3.5&Q"], 12, 1, 1);
    }


    /// Perform the per-event analysis
    void analyze(const Event& event) {

      /// @todo Do the event by event analysis here

        const DISKinematics& dk = apply<DISKinematics>(event, "Kinematics");

        // Get the DIS kinematics
        double xbj  = dk.x();
        double ybj = dk.y();
        double Q2 = dk.Q2()/GeV2;
        
        // Q2 and inelasticity cuts
        //cout << " after xbj " << xbj << endl;
        if (!inRange(ybj, 0.1, 0.6)) vetoEvent;
        //cout << " after ybj " << ybj << endl;
        if (!inRange(Q2, 2.0*GeV2, 70.0*GeV2)) vetoEvent;
        //cout << " after Q2 " << Q2 << endl;

      const FinalState& fs = apply<FinalState>(event, "FS");
      const size_t numParticles = fs.particles().size();
      //cout << " Num all     final state particles " << numParticles << endl;
      //proton beam energy: 820 GeV
      double e_proton = dk.beamHadron().E()/GeV;
      // Even if we only generate hadronic events, we still need a cut on numCharged >= 2.
      if (numParticles < 2) {
        MSG_DEBUG("Failed leptonic event cut");
        vetoEvent;
      }

      // Extracting the pi0 
      const UnstableParticles& ufs = apply<UnstableFinalState>(event, "UFS");
      //Get the hadronic CMS kinematics
      const LorentzTransform hcmboost = dk.boostHCM();

      for (const Particle& p : ufs.particles(Cuts::pid==PID::PI0)) {
        //cout << " Pid = " << p.pid() << endl;
        //Get the LAB kinematics
        double theta = p.theta();
        double eta = p.momentum().pseudorapidity();
        // double pT = p.momentum().pT()/GeV; 
        //Boost hcm
        const FourMomentum hcmMom = hcmboost.transform(p.momentum());
        double pThcm =hcmMom.pT();

        double e_pi0 = p.E()/GeV;
        double x_pi0_proton = e_pi0/e_proton;

        // epi0/e_proton, theta and pThcm cuts
        if(x_pi0_proton < 0.01) continue;
        //cout << " theta " << theta << " in deg " << theta/degree << endl;
        if (!inRange(theta/degree, 5, 25)) continue;
        if(pThcm < 2.5*GeV) continue;
        
        //Three cuts for Q2:
        if (Q2 > 2.0*GeV2 && Q2 < 4.5*GeV2){
          _h["Q:2.0-4.5-eta"]->fill(eta);
          _h["Q:2.0-4.5-p_T"]->fill(pThcm);
        } 
        if (Q2 > 4.5*GeV2 && Q2 < 15.0*GeV2){
          _h["Q:4.5-15.0-x"]->fill(xbj);
          _h["Q:4.5-15.0-eta"]->fill(eta);
          _h["Q:4.5-15.0-p_T"]->fill(pThcm);
        }  
        if (Q2 > 15.0*GeV2 && Q2 < 70.0*GeV2){
          _h["Q:15.0-70.0-x"]->fill(xbj);
          _h["Q:15.0-70.0-eta"]->fill(eta); 
          _h["Q:15.0-70.0-p_T"]->fill(pThcm);
        }

        //Two cuts for p_T:
        if (pThcm > 2.5*GeV){
          _h["p_T>2.5&x"]->fill(xbj);
          _h["p_T>2.5&Q"]->fill(Q2);
        }
        if (pThcm > 3.5*GeV){
          _h["p_T>3.5&x"]->fill(xbj);
          _h["p_T>3.5&Q"]->fill(Q2);
        }
        
      }

    }


    /// Normalise histograms etc., after the run
    void finalize() {
      const double sn = crossSection()/nanobarn/sumW();
      const double sp = crossSection()/picobarn/sumW();
      scale(_h["p_T>2.5&x"], sn);
      scale(_h["Q:2.0-4.5-eta"], sp);
      scale(_h["Q:2.0-4.5-p_T"], sp);
      scale(_h["Q:4.5-15.0-x"], sn);
      scale(_h["Q:4.5-15.0-eta"], sp);
      scale(_h["Q:4.5-15.0-p_T"], sp);
      scale(_h["Q:15.0-70.0-x"], sn);
      scale(_h["Q:15.0-70.0-eta"], sp);
      scale(_h["Q:15.0-70.0-p_T"], sp);
      scale(_h["p_T>2.5&Q"], sp);
      scale(_h["p_T>3.5&x"], sn);
      scale(_h["p_T>3.5&Q"], sp);
    }

    ///@}


    /// @name Histograms
    ///@{
    map<string, Histo1DPtr> _h;
    map<string, Profile1DPtr> _p;
    map<string, CounterPtr> _c;
    ///@}


  };


  RIVET_DECLARE_PLUGIN(H1_1999_I504022);

}