Rivet Analyses Reference

H1_1999_I481112

Measurement of D* meson cross-sections at HERA and determination of the gluon density in the proton using NLO QCD
Experiment: H1 (HERA)
Inspire ID: 481112
Status: VALIDATED
Authors:
  • Luca Marsili
  • Hannes Jung
References:Beams: e+ p+, p+ e+
Beam energies: (27.5, 820.0); (820.0, 27.5) GeV
    No run details listed

With the H1 detector at the $ep$ collider HERA, $D^*$ meson production cross sections have been measured in deep inelastic scattering with four-momentum transfers $Q^2>2$ GeV$^2$ and in photoproduction at energies around $W_{\gamma p} \sim 88$ GeV and 194 GeV. Next-to-Leading Order QCD calculations are found to describe the differential cross sections within theoretical and experimental uncertainties. Using these calculations, the NLO gluon momentum distribution in the proton, $x_g g(x_g)$, has been extracted in the momentum fraction range $7.5 \times 10^{-4}< x_g <4 \times 10^{-2}$ at average scales $mu^2 =25$ to 50 GeV$^2$. The gluon momentum fraction $x_g$ has been obtained from the measured kinematics of the scattered electron and the $D^*$ meson in the final state.

Source code: H1_1999_I481112.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
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
// -*- C++ -*-
#include "Rivet/Analysis.hh"
#include "Rivet/Projections/FinalState.hh"
#include "Rivet/Projections/FastJets.hh"
#include "Rivet/Projections/DISKinematics.hh"
#include "Rivet/Projections/UnstableParticles.hh"

namespace Rivet {


  /// @brief Measurement of D* meson cross-sections at HERA
  class H1_1999_I481112 : public Analysis {
  public:

    /// Constructor
    RIVET_DEFAULT_ANALYSIS_CTOR(H1_1999_I481112);


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

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

      // The basic final-state projection:
      // all final-state particles within
      // the given eta acceptance
      const FinalState fs(Cuts::abseta < 1.5);
      declare(fs, "fs");
      // The final-state particles declared above are clustered using FastJet with
    
      //Initialize quantities needed for cuts
      declare(DISKinematics(), "Kinematics");
      declare(UnstableParticles(), "DStars");
      

      Histo1DPtr dummy;
     
      book(_h["211"],2,1,1);
      book(_h["311"],3,1,1);
      book(_h["411"],4,1,1);
      book(_h["511"],5,1,1);
      book(_h["611"],6,1,1);
      book(_h["rap194"], 7,1,1);
      book(_h["pt194"], 8,1,1);
      book(_h["rap88"], 9,1,1);
      book(_h["pt88"], 10,1,1);
      _hpt.add(2.5, 3.5, book(dummy, 11,1,1));
      _hpt.add(3.5, 5.0, book(dummy, 11,1,2));
      _hpt.add(5.0, 10.0, book(dummy, 11,1,3));
      book(_h["1211"],12,1,1);
      book(_h["1212"],12,1,2);
      book(_h["1311"],13,1,1);
    }
    
    
    /// Perform the per-event analysis
    void analyze(const Event& event) {
     
      
      const DISKinematics& dk = applyProjection<DISKinematics>(event, "Kinematics");
      
      
      bool isDIS = false;
      bool ETAG44 = false;
      bool ETAG33 = false;

      const double y = dk.y();
      const double Q2 = dk.Q2();

      if(Q2 > 2 && Q2 <100 &&  y > 0.05 && y < 0.7) isDIS = true;
      if(Q2 < 0.009 && y > 0.02 && y <0.32  ) ETAG44 = true;
      if(Q2 < 0.01 &&  y > 0.29 && y <0.62 ) ETAG33 = true;
      if(isDIS == false && ETAG44 == false && ETAG33 == false) vetoEvent;
      
      
      //Creating array of D*
      Particles unstables ;
      if (isDIS ) { unstables = apply<ParticleFinder>(event, "DStars").particles(Cuts::pT > 1.5*GeV && Cuts::absetaIn(0,1.5));}
      else { unstables = apply<ParticleFinder>(event, "DStars").particles(Cuts::pT > 2.0*GeV && Cuts::absrapIn(0,1.5));}
      const Particles dstars = filter_select(unstables, [](const Particle& p){return p.abspid() == PID::DSPLUS;});
      
      if(dstars.empty() ) vetoEvent;
      MSG_DEBUG("D*" << dstars.size());
      
      const Particle& dstar = dstars.front();
      // boosting the system
      const LorentzTransform hcmboost = dk.boostHCM();
      const FourMomentum hcmMom = hcmboost.transform(dstar.momentum());

      //discriminate between dis and photoprod, and between ETA33 and ETA 44

      //kinematics quantities
      const double E = dstar.E();
      const double p_z = dstar.pz();
      //std::cout<<"y: "<<y<<endl;
      if(y<0.02) vetoEvent;
      const double m2 = 2.25; // charm mass^2
      const double E_e = dk.beamLepton().E();
      const double z = (E - p_z)/(2*y*E_e);
      
     
      const double M2 = (1.44*hcmMom.pT2() + m2)/(z*(1-z));
      const double x_g = (M2 + Q2)/(y*dk.s());
      // std::cout<<"s:  "<<dk.s()<<endl;
      // std::cout<<"M2: "<<M2<<endl;
      //  std::cout<<"z: "<<z<<endl;
      
      const double y_capp = dstar.rapidity();
      const double W = sqrt(dk.W2());
      //std::cout<<"x_g: "<<x_g<<endl;
      
      //perform the cuts
      if(isDIS == true){
         _h["211"]->fill(dstar.pT());
         _h["411"]->fill(dstar.eta());
         _h["511"]->fill(Q2);
         _h["611"]->fill(log10(x_g));
		//boosting to the hcm frame
         _h["311"]->fill(hcmMom.pT());
      }
      if(ETAG33 == true && abs(y_capp) < 1.5 && dstar.pT() > 2.5*GeV  ) {
		
         _h["rap194"]->fill(y_capp);
         _h["pt194"] ->fill(dstar.pT());
         _hpt.fill(dstar.pT(), y_capp);	
	
	if( W > 173 &&  W<273){
	   _h["1211"]->fill(log10(x_g));
	  
	}
	if( W > 130 &&  W<230){
	   _h["1212"]->fill(log10(x_g));
	  
	}
	
	
      }
      if(ETAG44 == true && abs(y_capp) < 1.5 && dstar.pT()> 2) { 

	   _h["pt88"]->fill(dstar.pT());
	   _h["rap88"]->fill(y_capp);
	   _h["1311"]->fill(log10(x_g));
	
      }
    }

    /// Normalise histograms etc., after the run
    void finalize() {
      // conversion factors from ep to gamma p xsections (as given in publication)
      const double F_etag33 = 0.0128;
      const double F_etag44 = 0.0838;

      double norm = crossSection()/nanobarn/sumW() ;
      scale( _h["211"], norm);
      
      scale(_h["311"], norm);  
      scale(_h["411"], norm);
      scale(_h["511"], norm);
      scale(_h["611"], norm);
      double norm_mub = crossSection()/microbarn/sumW();

      scale(_h["rap194"],norm_mub/F_etag33 );
      scale(_h["pt194"], norm_mub/F_etag33 );
      scale(_h["rap88"], norm_mub/F_etag44 );
      scale(_h["pt88"],  norm_mub/F_etag44 );
         
      _hpt.scale( norm/F_etag33, this);
      scale(_h["1211"], norm_mub/F_etag33);
      scale(_h["1212"], norm_mub/F_etag33);
      scale(_h["1311"], norm_mub/F_etag44);
    }

    ///@}


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


  };


  RIVET_DECLARE_PLUGIN(H1_1999_I481112);

}