rivet is hosted by Hepforge, IPPP Durham
ATLAS_2011_I919017.cc
Go to the documentation of this file.
00001 // -*- C++ -*-
00002 #include "Rivet/Analysis.hh"
00003 #include "Rivet/Projections/ChargedFinalState.hh"
00004 #include "Rivet/Projections/FastJets.hh"
00005 #include "Rivet/Math/MathUtils.hh"
00006 
00007 namespace Rivet {
00008 
00009 
00010   namespace {
00011 
00012     inline double calcz(const Jet& j, const Particle& p) {
00013       const double num = j.momentum().vector3().dot(p.momentum().vector3());
00014       const double den = j.momentum().vector3().mod2();
00015       return num/den;
00016     }
00017 
00018     inline double calcptrel(const Jet& j, const Particle& p) {
00019       const double num = j.momentum().vector3().cross(p.momentum().vector3()).mod();
00020       const double den = j.momentum().vector3().mod();
00021       return num/den;
00022     }
00023 
00024     inline double calcr(const Jet& j, const Particle& p) {
00025       return deltaR(j.rapidity(), j.momentum().phi(), p.rapidity(), p.momentum().phi());
00026     }
00027 
00028     // For annulus area kludge
00029     /// @todo Improve somehow... need normalisation *without* bin width factors!
00030     inline double calcrweight(const Jet& j, const Particle& p) {
00031       size_t nBins_r = 26;
00032       double bins_r[] = { 0.00, 0.01, 0.02, 0.03, 0.04, 0.05, 0.06, 0.07, 0.08, 0.09, 0.10,
00033                           0.12, 0.14, 0.16, 0.18, 0.20, 0.22, 0.24, 0.26, 0.28, 0.30,
00034                           0.35, 0.40, 0.45, 0.50, 0.55, 0.60 };
00035       double r = calcr(j,p);
00036       for (size_t bin = 0 ; bin < nBins_r ; bin++) {
00037         if (r < bins_r[bin+1]) {
00038           double up = bins_r[bin+1];
00039           double down = bins_r[bin];
00040           return ((up-down)/(M_PI*(up*up-down*down)));
00041         }
00042       }
00043       return 1.0;
00044     }
00045   }
00046 
00047 
00048   class ATLAS_2011_I919017 : public Analysis {
00049   public:
00050 
00051     /// @name Constructors etc.
00052     //@{
00053 
00054     /// Constructor
00055     ATLAS_2011_I919017()
00056       : Analysis("ATLAS_2011_I919017"),
00057         _sumofweights04(0), _sumofweights06(0),
00058         _numjets04_00_05_04_06(0), _numjets04_00_05_06_10(0), _numjets04_00_05_10_15(0), _numjets04_00_05_15_24(0), _numjets04_00_05_24_40(0),
00059         _numjets06_00_05_04_06(0), _numjets06_00_05_06_10(0), _numjets06_00_05_10_15(0), _numjets06_00_05_15_24(0), _numjets06_00_05_24_40(0),
00060         _numjets04_05_10_04_06(0), _numjets04_05_10_06_10(0), _numjets04_05_10_10_15(0), _numjets04_05_10_15_24(0), _numjets04_05_10_24_40(0),
00061         _numjets06_05_10_04_06(0), _numjets06_05_10_06_10(0), _numjets06_05_10_10_15(0), _numjets06_05_10_15_24(0), _numjets06_05_10_24_40(0),
00062         _numjets04_10_15_04_06(0), _numjets04_10_15_06_10(0), _numjets04_10_15_10_15(0), _numjets04_10_15_15_24(0), _numjets04_10_15_24_40(0),
00063         _numjets06_10_15_04_06(0), _numjets06_10_15_06_10(0), _numjets06_10_15_10_15(0), _numjets06_10_15_15_24(0), _numjets06_10_15_24_40(0),
00064         _numjets04_15_19_04_06(0), _numjets04_15_19_06_10(0), _numjets04_15_19_10_15(0), _numjets04_15_19_15_24(0), _numjets04_15_19_24_40(0),
00065         _numjets06_15_19_04_06(0), _numjets06_15_19_06_10(0), _numjets06_15_19_10_15(0), _numjets06_15_19_15_24(0), _numjets06_15_19_24_40(0),
00066         _numjets04_00_19_04_06(0), _numjets04_00_19_06_10(0), _numjets04_00_19_10_15(0), _numjets04_00_19_15_24(0), _numjets04_00_19_24_40(0),
00067         _numjets06_00_19_04_06(0), _numjets06_00_19_06_10(0), _numjets06_00_19_10_15(0), _numjets06_00_19_15_24(0), _numjets06_00_19_24_40(0)
00068     {    }
00069 
00070     //@}
00071 
00072 
00073   public:
00074 
00075     /// @name Analysis methods
00076     //@{
00077 
00078     /// Book histograms and initialise projections before the run
00079     void init() {
00080       ChargedFinalState cfs(-2.5, 2.5, 0.3*GeV);
00081       FastJets trkjets04(cfs, FastJets::ANTIKT, 0.4);
00082       FastJets trkjets06(cfs, FastJets::ANTIKT, 0.6);
00083       addProjection(trkjets04, "Jets04");
00084       addProjection(trkjets06, "Jets06");
00085 
00086       // Book histograms
00087       _h_pt04_00_05 = bookHisto1D(1, 1, 1);
00088       _h_pt06_00_05 = bookHisto1D(2, 1, 1);
00089       _h_N04_00_05_04_06 = bookHisto1D(1, 2, 1+5);
00090       _h_N06_00_05_04_06 = bookHisto1D(2, 2, 1+5);
00091       _h_N04_00_05_06_10 = bookHisto1D(1, 2, 2+5);
00092       _h_N06_00_05_06_10 = bookHisto1D(2, 2, 2+5);
00093       _h_N04_00_05_10_15 = bookHisto1D(1, 2, 3+5);
00094       _h_N06_00_05_10_15 = bookHisto1D(2, 2, 3+5);
00095       _h_N04_00_05_15_24 = bookHisto1D(1, 2, 4+5);
00096       _h_N06_00_05_15_24 = bookHisto1D(2, 2, 4+5);
00097       _h_N04_00_05_24_40 = bookHisto1D(1, 2, 5+5);
00098       _h_N06_00_05_24_40 = bookHisto1D(2, 2, 5+5);
00099       _h_z04_00_05_04_06 = bookHisto1D(1, 3, 1+5);
00100       _h_z06_00_05_04_06 = bookHisto1D(2, 3, 1+5);
00101       _h_z04_00_05_06_10 = bookHisto1D(1, 3, 2+5);
00102       _h_z06_00_05_06_10 = bookHisto1D(2, 3, 2+5);
00103       _h_z04_00_05_10_15 = bookHisto1D(1, 3, 3+5);
00104       _h_z06_00_05_10_15 = bookHisto1D(2, 3, 3+5);
00105       _h_z04_00_05_15_24 = bookHisto1D(1, 3, 4+5);
00106       _h_z06_00_05_15_24 = bookHisto1D(2, 3, 4+5);
00107       _h_z04_00_05_24_40 = bookHisto1D(1, 3, 5+5);
00108       _h_z06_00_05_24_40 = bookHisto1D(2, 3, 5+5);
00109       _h_ptrel04_00_05_04_06 = bookHisto1D(1, 4, 1+5);
00110       _h_ptrel06_00_05_04_06 = bookHisto1D(2, 4, 1+5);
00111       _h_ptrel04_00_05_06_10 = bookHisto1D(1, 4, 2+5);
00112       _h_ptrel06_00_05_06_10 = bookHisto1D(2, 4, 2+5);
00113       _h_ptrel04_00_05_10_15 = bookHisto1D(1, 4, 3+5);
00114       _h_ptrel06_00_05_10_15 = bookHisto1D(2, 4, 3+5);
00115       _h_ptrel04_00_05_15_24 = bookHisto1D(1, 4, 4+5);
00116       _h_ptrel06_00_05_15_24 = bookHisto1D(2, 4, 4+5);
00117       _h_ptrel04_00_05_24_40 = bookHisto1D(1, 4, 5+5);
00118       _h_ptrel06_00_05_24_40 = bookHisto1D(2, 4, 5+5);
00119       _h_rdA04_00_05_04_06 = bookHisto1D(1, 5, 1+5);
00120       _h_rdA06_00_05_04_06 = bookHisto1D(2, 5, 1+5);
00121       _h_rdA04_00_05_06_10 = bookHisto1D(1, 5, 2+5);
00122       _h_rdA06_00_05_06_10 = bookHisto1D(2, 5, 2+5);
00123       _h_rdA04_00_05_10_15 = bookHisto1D(1, 5, 3+5);
00124       _h_rdA06_00_05_10_15 = bookHisto1D(2, 5, 3+5);
00125       _h_rdA04_00_05_15_24 = bookHisto1D(1, 5, 4+5);
00126       _h_rdA06_00_05_15_24 = bookHisto1D(2, 5, 4+5);
00127       _h_rdA04_00_05_24_40 = bookHisto1D(1, 5, 5+5);
00128       _h_rdA06_00_05_24_40 = bookHisto1D(2, 5, 5+5);
00129 
00130       _h_pt04_05_10 = bookHisto1D(1, 1, 2);
00131       _h_pt06_05_10 = bookHisto1D(2, 1, 2);
00132       _h_N04_05_10_04_06 = bookHisto1D(1, 2, 1+10);
00133       _h_N06_05_10_04_06 = bookHisto1D(2, 2, 1+10);
00134       _h_N04_05_10_06_10 = bookHisto1D(1, 2, 2+10);
00135       _h_N06_05_10_06_10 = bookHisto1D(2, 2, 2+10);
00136       _h_N04_05_10_10_15 = bookHisto1D(1, 2, 3+10);
00137       _h_N06_05_10_10_15 = bookHisto1D(2, 2, 3+10);
00138       _h_N04_05_10_15_24 = bookHisto1D(1, 2, 4+10);
00139       _h_N06_05_10_15_24 = bookHisto1D(2, 2, 4+10);
00140       _h_N04_05_10_24_40 = bookHisto1D(1, 2, 5+10);
00141       _h_N06_05_10_24_40 = bookHisto1D(2, 2, 5+10);
00142       _h_z04_05_10_04_06 = bookHisto1D(1, 3, 1+10);
00143       _h_z06_05_10_04_06 = bookHisto1D(2, 3, 1+10);
00144       _h_z04_05_10_06_10 = bookHisto1D(1, 3, 2+10);
00145       _h_z06_05_10_06_10 = bookHisto1D(2, 3, 2+10);
00146       _h_z04_05_10_10_15 = bookHisto1D(1, 3, 3+10);
00147       _h_z06_05_10_10_15 = bookHisto1D(2, 3, 3+10);
00148       _h_z04_05_10_15_24 = bookHisto1D(1, 3, 4+10);
00149       _h_z06_05_10_15_24 = bookHisto1D(2, 3, 4+10);
00150       _h_z04_05_10_24_40 = bookHisto1D(1, 3, 5+10);
00151       _h_z06_05_10_24_40 = bookHisto1D(2, 3, 5+10);
00152       _h_ptrel04_05_10_04_06 = bookHisto1D(1, 4, 1+10);
00153       _h_ptrel06_05_10_04_06 = bookHisto1D(2, 4, 1+10);
00154       _h_ptrel04_05_10_06_10 = bookHisto1D(1, 4, 2+10);
00155       _h_ptrel06_05_10_06_10 = bookHisto1D(2, 4, 2+10);
00156       _h_ptrel04_05_10_10_15 = bookHisto1D(1, 4, 3+10);
00157       _h_ptrel06_05_10_10_15 = bookHisto1D(2, 4, 3+10);
00158       _h_ptrel04_05_10_15_24 = bookHisto1D(1, 4, 4+10);
00159       _h_ptrel06_05_10_15_24 = bookHisto1D(2, 4, 4+10);
00160       _h_ptrel04_05_10_24_40 = bookHisto1D(1, 4, 5+10);
00161       _h_ptrel06_05_10_24_40 = bookHisto1D(2, 4, 5+10);
00162       _h_rdA04_05_10_04_06 = bookHisto1D(1, 5, 1+10);
00163       _h_rdA06_05_10_04_06 = bookHisto1D(2, 5, 1+10);
00164       _h_rdA04_05_10_06_10 = bookHisto1D(1, 5, 2+10);
00165       _h_rdA06_05_10_06_10 = bookHisto1D(2, 5, 2+10);
00166       _h_rdA04_05_10_10_15 = bookHisto1D(1, 5, 3+10);
00167       _h_rdA06_05_10_10_15 = bookHisto1D(2, 5, 3+10);
00168       _h_rdA04_05_10_15_24 = bookHisto1D(1, 5, 4+10);
00169       _h_rdA06_05_10_15_24 = bookHisto1D(2, 5, 4+10);
00170       _h_rdA04_05_10_24_40 = bookHisto1D(1, 5, 5+10);
00171       _h_rdA06_05_10_24_40 = bookHisto1D(2, 5, 5+10);
00172 
00173       _h_pt04_10_15 = bookHisto1D(1, 1, 3);
00174       _h_pt06_10_15 = bookHisto1D(2, 1, 3);
00175       _h_N04_10_15_04_06 = bookHisto1D(1, 2, 1+15);
00176       _h_N06_10_15_04_06 = bookHisto1D(2, 2, 1+15);
00177       _h_N04_10_15_06_10 = bookHisto1D(1, 2, 2+15);
00178       _h_N06_10_15_06_10 = bookHisto1D(2, 2, 2+15);
00179       _h_N04_10_15_10_15 = bookHisto1D(1, 2, 3+15);
00180       _h_N06_10_15_10_15 = bookHisto1D(2, 2, 3+15);
00181       _h_N04_10_15_15_24 = bookHisto1D(1, 2, 4+15);
00182       _h_N06_10_15_15_24 = bookHisto1D(2, 2, 4+15);
00183       _h_N04_10_15_24_40 = bookHisto1D(1, 2, 5+15);
00184       _h_N06_10_15_24_40 = bookHisto1D(2, 2, 5+15);
00185       _h_z04_10_15_04_06 = bookHisto1D(1, 3, 1+15);
00186       _h_z06_10_15_04_06 = bookHisto1D(2, 3, 1+15);
00187       _h_z04_10_15_06_10 = bookHisto1D(1, 3, 2+15);
00188       _h_z06_10_15_06_10 = bookHisto1D(2, 3, 2+15);
00189       _h_z04_10_15_10_15 = bookHisto1D(1, 3, 3+15);
00190       _h_z06_10_15_10_15 = bookHisto1D(2, 3, 3+15);
00191       _h_z04_10_15_15_24 = bookHisto1D(1, 3, 4+15);
00192       _h_z06_10_15_15_24 = bookHisto1D(2, 3, 4+15);
00193       _h_z04_10_15_24_40 = bookHisto1D(1, 3, 5+15);
00194       _h_z06_10_15_24_40 = bookHisto1D(2, 3, 5+15);
00195       _h_ptrel04_10_15_04_06 = bookHisto1D(1, 4, 1+15);
00196       _h_ptrel06_10_15_04_06 = bookHisto1D(2, 4, 1+15);
00197       _h_ptrel04_10_15_06_10 = bookHisto1D(1, 4, 2+15);
00198       _h_ptrel06_10_15_06_10 = bookHisto1D(2, 4, 2+15);
00199       _h_ptrel04_10_15_10_15 = bookHisto1D(1, 4, 3+15);
00200       _h_ptrel06_10_15_10_15 = bookHisto1D(2, 4, 3+15);
00201       _h_ptrel04_10_15_15_24 = bookHisto1D(1, 4, 4+15);
00202       _h_ptrel06_10_15_15_24 = bookHisto1D(2, 4, 4+15);
00203       _h_ptrel04_10_15_24_40 = bookHisto1D(1, 4, 5+15);
00204       _h_ptrel06_10_15_24_40 = bookHisto1D(2, 4, 5+15);
00205       _h_rdA04_10_15_04_06 = bookHisto1D(1, 5, 1+15);
00206       _h_rdA06_10_15_04_06 = bookHisto1D(2, 5, 1+15);
00207       _h_rdA04_10_15_06_10 = bookHisto1D(1, 5, 2+15);
00208       _h_rdA06_10_15_06_10 = bookHisto1D(2, 5, 2+15);
00209       _h_rdA04_10_15_10_15 = bookHisto1D(1, 5, 3+15);
00210       _h_rdA06_10_15_10_15 = bookHisto1D(2, 5, 3+15);
00211       _h_rdA04_10_15_15_24 = bookHisto1D(1, 5, 4+15);
00212       _h_rdA06_10_15_15_24 = bookHisto1D(2, 5, 4+15);
00213       _h_rdA04_10_15_24_40 = bookHisto1D(1, 5, 5+15);
00214       _h_rdA06_10_15_24_40 = bookHisto1D(2, 5, 5+15);
00215 
00216       _h_pt04_15_19 = bookHisto1D(1, 1, 4);
00217       _h_pt06_15_19 = bookHisto1D(2, 1, 4);
00218       _h_N04_15_19_04_06 = bookHisto1D(1, 2, 1+20);
00219       _h_N06_15_19_04_06 = bookHisto1D(2, 2, 1+20);
00220       _h_N04_15_19_06_10 = bookHisto1D(1, 2, 2+20);
00221       _h_N06_15_19_06_10 = bookHisto1D(2, 2, 2+20);
00222       _h_N04_15_19_10_15 = bookHisto1D(1, 2, 3+20);
00223       _h_N06_15_19_10_15 = bookHisto1D(2, 2, 3+20);
00224       _h_N04_15_19_15_24 = bookHisto1D(1, 2, 4+20);
00225       _h_N06_15_19_15_24 = bookHisto1D(2, 2, 4+20);
00226       _h_N04_15_19_24_40 = bookHisto1D(1, 2, 5+20);
00227       _h_N06_15_19_24_40 = bookHisto1D(2, 2, 5+20);
00228       _h_z04_15_19_04_06 = bookHisto1D(1, 3, 1+20);
00229       _h_z06_15_19_04_06 = bookHisto1D(2, 3, 1+20);
00230       _h_z04_15_19_06_10 = bookHisto1D(1, 3, 2+20);
00231       _h_z06_15_19_06_10 = bookHisto1D(2, 3, 2+20);
00232       _h_z04_15_19_10_15 = bookHisto1D(1, 3, 3+20);
00233       _h_z06_15_19_10_15 = bookHisto1D(2, 3, 3+20);
00234       _h_z04_15_19_15_24 = bookHisto1D(1, 3, 4+20);
00235       _h_z06_15_19_15_24 = bookHisto1D(2, 3, 4+20);
00236       _h_z04_15_19_24_40 = bookHisto1D(1, 3, 5+20);
00237       _h_z06_15_19_24_40 = bookHisto1D(2, 3, 5+20);
00238       _h_ptrel04_15_19_04_06 = bookHisto1D(1, 4, 1+20);
00239       _h_ptrel06_15_19_04_06 = bookHisto1D(2, 4, 1+20);
00240       _h_ptrel04_15_19_06_10 = bookHisto1D(1, 4, 2+20);
00241       _h_ptrel06_15_19_06_10 = bookHisto1D(2, 4, 2+20);
00242       _h_ptrel04_15_19_10_15 = bookHisto1D(1, 4, 3+20);
00243       _h_ptrel06_15_19_10_15 = bookHisto1D(2, 4, 3+20);
00244       _h_ptrel04_15_19_15_24 = bookHisto1D(1, 4, 4+20);
00245       _h_ptrel06_15_19_15_24 = bookHisto1D(2, 4, 4+20);
00246       _h_ptrel04_15_19_24_40 = bookHisto1D(1, 4, 5+20);
00247       _h_ptrel06_15_19_24_40 = bookHisto1D(2, 4, 5+20);
00248       _h_rdA04_15_19_04_06 = bookHisto1D(1, 5, 1+20);
00249       _h_rdA06_15_19_04_06 = bookHisto1D(2, 5, 1+20);
00250       _h_rdA04_15_19_06_10 = bookHisto1D(1, 5, 2+20);
00251       _h_rdA06_15_19_06_10 = bookHisto1D(2, 5, 2+20);
00252       _h_rdA04_15_19_10_15 = bookHisto1D(1, 5, 3+20);
00253       _h_rdA06_15_19_10_15 = bookHisto1D(2, 5, 3+20);
00254       _h_rdA04_15_19_15_24 = bookHisto1D(1, 5, 4+20);
00255       _h_rdA06_15_19_15_24 = bookHisto1D(2, 5, 4+20);
00256       _h_rdA04_15_19_24_40 = bookHisto1D(1, 5, 5+20);
00257       _h_rdA06_15_19_24_40 = bookHisto1D(2, 5, 5+20);
00258 
00259       _h_N04_00_19_04_06 = bookHisto1D(1, 2, 1+0);
00260       _h_N06_00_19_04_06 = bookHisto1D(2, 2, 1+0);
00261       _h_N04_00_19_06_10 = bookHisto1D(1, 2, 2+0);
00262       _h_N06_00_19_06_10 = bookHisto1D(2, 2, 2+0);
00263       _h_N04_00_19_10_15 = bookHisto1D(1, 2, 3+0);
00264       _h_N06_00_19_10_15 = bookHisto1D(2, 2, 3+0);
00265       _h_N04_00_19_15_24 = bookHisto1D(1, 2, 4+0);
00266       _h_N06_00_19_15_24 = bookHisto1D(2, 2, 4+0);
00267       _h_N04_00_19_24_40 = bookHisto1D(1, 2, 5+0);
00268       _h_N06_00_19_24_40 = bookHisto1D(2, 2, 5+0);
00269       _h_z04_00_19_04_06 = bookHisto1D(1, 3, 1+0);
00270       _h_z06_00_19_04_06 = bookHisto1D(2, 3, 1+0);
00271       _h_z04_00_19_06_10 = bookHisto1D(1, 3, 2+0);
00272       _h_z06_00_19_06_10 = bookHisto1D(2, 3, 2+0);
00273       _h_z04_00_19_10_15 = bookHisto1D(1, 3, 3+0);
00274       _h_z06_00_19_10_15 = bookHisto1D(2, 3, 3+0);
00275       _h_z04_00_19_15_24 = bookHisto1D(1, 3, 4+0);
00276       _h_z06_00_19_15_24 = bookHisto1D(2, 3, 4+0);
00277       _h_z04_00_19_24_40 = bookHisto1D(1, 3, 5+0);
00278       _h_z06_00_19_24_40 = bookHisto1D(2, 3, 5+0);
00279       _h_ptrel04_00_19_04_06 = bookHisto1D(1, 4, 1+0);
00280       _h_ptrel06_00_19_04_06 = bookHisto1D(2, 4, 1+0);
00281       _h_ptrel04_00_19_06_10 = bookHisto1D(1, 4, 2+0);
00282       _h_ptrel06_00_19_06_10 = bookHisto1D(2, 4, 2+0);
00283       _h_ptrel04_00_19_10_15 = bookHisto1D(1, 4, 3+0);
00284       _h_ptrel06_00_19_10_15 = bookHisto1D(2, 4, 3+0);
00285       _h_ptrel04_00_19_15_24 = bookHisto1D(1, 4, 4+0);
00286       _h_ptrel06_00_19_15_24 = bookHisto1D(2, 4, 4+0);
00287       _h_ptrel04_00_19_24_40 = bookHisto1D(1, 4, 5+0);
00288       _h_ptrel06_00_19_24_40 = bookHisto1D(2, 4, 5+0);
00289       _h_rdA04_00_19_04_06 = bookHisto1D(1, 5, 1+0);
00290       _h_rdA06_00_19_04_06 = bookHisto1D(2, 5, 1+0);
00291       _h_rdA04_00_19_06_10 = bookHisto1D(1, 5, 2+0);
00292       _h_rdA06_00_19_06_10 = bookHisto1D(2, 5, 2+0);
00293       _h_rdA04_00_19_10_15 = bookHisto1D(1, 5, 3+0);
00294       _h_rdA06_00_19_10_15 = bookHisto1D(2, 5, 3+0);
00295       _h_rdA04_00_19_15_24 = bookHisto1D(1, 5, 4+0);
00296       _h_rdA06_00_19_15_24 = bookHisto1D(2, 5, 4+0);
00297       _h_rdA04_00_19_24_40 = bookHisto1D(1, 5, 5+0);
00298       _h_rdA06_00_19_24_40 = bookHisto1D(2, 5, 5+0);
00299     }
00300 
00301 
00302     /// Perform the per-event analysis
00303     void analyze(const Event& event) {
00304       const double weight = event.weight();
00305 
00306       const Jets& jets04 = applyProjection<JetAlg>(event, "Jets04").jets();
00307       if (!jets04.empty()) {
00308         _sumofweights04 += weight;
00309         foreach (const Jet& j, jets04) {
00310           const double jetpt = j.pT();
00311           if (fabs(j.rapidity()) < 0.5) {
00312             _h_pt04_00_05->fill(jetpt/GeV, weight);
00313             if (inRange(jetpt/GeV, 4., 6.)) {
00314               _numjets04_00_05_04_06 += weight;
00315               _h_N04_00_05_04_06->fill(j.particles().size(),weight);
00316               foreach (const Particle& p, j.particles()) {
00317                 _h_z04_00_05_04_06->fill(calcz(j,p),weight);
00318                 _h_ptrel04_00_05_04_06->fill(calcptrel(j,p),weight);
00319                 _h_rdA04_00_05_04_06->fill(calcr(j,p),weight*calcrweight(j,p));
00320               }
00321             }
00322             if (inRange(jetpt/GeV, 6., 10.)) {
00323               _numjets04_00_05_06_10 += weight;
00324               _h_N04_00_05_06_10->fill(j.particles().size(),weight);
00325               foreach (const Particle& p, j.particles()) {
00326                 _h_z04_00_05_06_10->fill(calcz(j,p),weight);
00327                 _h_ptrel04_00_05_06_10->fill(calcptrel(j,p),weight);
00328                 _h_rdA04_00_05_06_10->fill(calcr(j,p),weight*calcrweight(j,p));
00329               }
00330             }
00331             if (inRange(jetpt/GeV, 10., 15.)) {
00332               _numjets04_00_05_10_15 += weight;
00333               _h_N04_00_05_10_15->fill(j.particles().size(),weight);
00334               foreach (const Particle& p, j.particles()) {
00335                 _h_z04_00_05_10_15->fill(calcz(j,p),weight);
00336                 _h_ptrel04_00_05_10_15->fill(calcptrel(j,p),weight);
00337                 _h_rdA04_00_05_10_15->fill(calcr(j,p),weight*calcrweight(j,p));
00338               }
00339             }
00340             if (inRange(jetpt/GeV, 15., 24.)) {
00341               _numjets04_00_05_15_24 += weight;
00342               _h_N04_00_05_15_24->fill(j.particles().size(),weight);
00343               foreach (const Particle& p, j.particles()) {
00344                 _h_z04_00_05_15_24->fill(calcz(j,p),weight);
00345                 _h_ptrel04_00_05_15_24->fill(calcptrel(j,p),weight);
00346                 _h_rdA04_00_05_15_24->fill(calcr(j,p),weight*calcrweight(j,p));
00347               }
00348             }
00349             if (inRange(jetpt/GeV, 24., 40.)) {
00350               _numjets04_00_05_24_40 += weight;
00351               _h_N04_00_05_24_40->fill(j.particles().size(),weight);
00352               foreach (const Particle& p, j.particles()) {
00353                 _h_z04_00_05_24_40->fill(calcz(j,p),weight);
00354                 _h_ptrel04_00_05_24_40->fill(calcptrel(j,p),weight);
00355                 _h_rdA04_00_05_24_40->fill(calcr(j,p),weight*calcrweight(j,p));
00356               }
00357             }
00358           }
00359           if (fabs(j.rapidity()) > 0.5 && fabs(j.rapidity()) < 1.0) {
00360             _h_pt04_05_10->fill(jetpt/GeV, weight);
00361             if (inRange(jetpt/GeV, 4., 6.)) {
00362               _numjets04_05_10_04_06 += weight;
00363               _h_N04_05_10_04_06->fill(j.particles().size(),weight);
00364               foreach (const Particle& p, j.particles()) {
00365                 _h_z04_05_10_04_06->fill(calcz(j,p),weight);
00366                 _h_ptrel04_05_10_04_06->fill(calcptrel(j,p),weight);
00367                 _h_rdA04_05_10_04_06->fill(calcr(j,p),weight*calcrweight(j,p));
00368               }
00369             }
00370             if (inRange(jetpt/GeV, 6., 10.)) {
00371               _numjets04_05_10_06_10 += weight;
00372               _h_N04_05_10_06_10->fill(j.particles().size(),weight);
00373               foreach (const Particle& p, j.particles()) {
00374                 _h_z04_05_10_06_10->fill(calcz(j,p),weight);
00375                 _h_ptrel04_05_10_06_10->fill(calcptrel(j,p),weight);
00376                 _h_rdA04_05_10_06_10->fill(calcr(j,p),weight*calcrweight(j,p));
00377               }
00378             }
00379             if (inRange(jetpt/GeV, 10., 15.)) {
00380               _numjets04_05_10_10_15 += weight;
00381               _h_N04_05_10_10_15->fill(j.particles().size(),weight);
00382               foreach (const Particle& p, j.particles()) {
00383                 _h_z04_05_10_10_15->fill(calcz(j,p),weight);
00384                 _h_ptrel04_05_10_10_15->fill(calcptrel(j,p),weight);
00385                 _h_rdA04_05_10_10_15->fill(calcr(j,p),weight*calcrweight(j,p));
00386               }
00387             }
00388             if (inRange(jetpt/GeV, 15., 24.)) {
00389               _numjets04_05_10_15_24 += weight;
00390               _h_N04_05_10_15_24->fill(j.particles().size(),weight);
00391               foreach (const Particle& p, j.particles()) {
00392                 _h_z04_05_10_15_24->fill(calcz(j,p),weight);
00393                 _h_ptrel04_05_10_15_24->fill(calcptrel(j,p),weight);
00394                 _h_rdA04_05_10_15_24->fill(calcr(j,p),weight*calcrweight(j,p));
00395               }
00396             }
00397             if (inRange(jetpt/GeV, 24., 40.)) {
00398               _numjets04_05_10_24_40 += weight;
00399               _h_N04_05_10_24_40->fill(j.particles().size(),weight);
00400               foreach (const Particle& p, j.particles()) {
00401                 _h_z04_05_10_24_40->fill(calcz(j,p),weight);
00402                 _h_ptrel04_05_10_24_40->fill(calcptrel(j,p),weight);
00403                 _h_rdA04_05_10_24_40->fill(calcr(j,p),weight*calcrweight(j,p));
00404               }
00405             }
00406           }
00407           if (fabs(j.rapidity()) > 1.0 && fabs(j.rapidity()) < 1.5) {
00408             _h_pt04_10_15->fill(jetpt/GeV, weight);
00409             if (inRange(jetpt/GeV, 4., 6.)) {
00410               _numjets04_10_15_04_06 += weight;
00411               _h_N04_10_15_04_06->fill(j.particles().size(),weight);
00412               foreach (const Particle& p, j.particles()) {
00413                 _h_z04_10_15_04_06->fill(calcz(j,p),weight);
00414                 _h_ptrel04_10_15_04_06->fill(calcptrel(j,p),weight);
00415                 _h_rdA04_10_15_04_06->fill(calcr(j,p),weight*calcrweight(j,p));
00416               }
00417             }
00418             if (inRange(jetpt/GeV, 6., 10.)) {
00419               _numjets04_10_15_06_10 += weight;
00420               _h_N04_10_15_06_10->fill(j.particles().size(),weight);
00421               foreach (const Particle& p, j.particles()) {
00422                 _h_z04_10_15_06_10->fill(calcz(j,p),weight);
00423                 _h_ptrel04_10_15_06_10->fill(calcptrel(j,p),weight);
00424                 _h_rdA04_10_15_06_10->fill(calcr(j,p),weight*calcrweight(j,p));
00425               }
00426             }
00427             if (inRange(jetpt/GeV, 10., 15.)) {
00428               _numjets04_10_15_10_15 += weight;
00429               _h_N04_10_15_10_15->fill(j.particles().size(),weight);
00430               foreach (const Particle& p, j.particles()) {
00431                 _h_z04_10_15_10_15->fill(calcz(j,p),weight);
00432                 _h_ptrel04_10_15_10_15->fill(calcptrel(j,p),weight);
00433                 _h_rdA04_10_15_10_15->fill(calcr(j,p),weight*calcrweight(j,p));
00434               }
00435             }
00436             if (inRange(jetpt/GeV, 15., 24.)) {
00437               _numjets04_10_15_15_24 += weight;
00438               _h_N04_10_15_15_24->fill(j.particles().size(),weight);
00439               foreach (const Particle& p, j.particles()) {
00440                 _h_z04_10_15_15_24->fill(calcz(j,p),weight);
00441                 _h_ptrel04_10_15_15_24->fill(calcptrel(j,p),weight);
00442                 _h_rdA04_10_15_15_24->fill(calcr(j,p),weight*calcrweight(j,p));
00443               }
00444             }
00445             if (inRange(jetpt/GeV, 24., 40.)) {
00446               _numjets04_10_15_24_40 += weight;
00447               _h_N04_10_15_24_40->fill(j.particles().size(),weight);
00448               foreach (const Particle& p, j.particles()) {
00449                 _h_z04_10_15_24_40->fill(calcz(j,p),weight);
00450                 _h_ptrel04_10_15_24_40->fill(calcptrel(j,p),weight);
00451                 _h_rdA04_10_15_24_40->fill(calcr(j,p),weight*calcrweight(j,p));
00452               }
00453             }
00454           }
00455           if (fabs(j.rapidity()) > 1.5 && fabs(j.rapidity()) < 1.9) {
00456             _h_pt04_15_19->fill(jetpt/GeV, weight);
00457             if (inRange(jetpt/GeV, 4., 6.)) {
00458               _numjets04_15_19_04_06 += weight;
00459               _h_N04_15_19_04_06->fill(j.particles().size(),weight);
00460               foreach (const Particle& p, j.particles()) {
00461                 _h_z04_15_19_04_06->fill(calcz(j,p),weight);
00462                 _h_ptrel04_15_19_04_06->fill(calcptrel(j,p),weight);
00463                 _h_rdA04_15_19_04_06->fill(calcr(j,p),weight*calcrweight(j,p));
00464               }
00465             }
00466             if (inRange(jetpt/GeV, 6., 10.)) {
00467               _numjets04_15_19_06_10 += weight;
00468               _h_N04_15_19_06_10->fill(j.particles().size(),weight);
00469               foreach (const Particle& p, j.particles()) {
00470                 _h_z04_15_19_06_10->fill(calcz(j,p),weight);
00471                 _h_ptrel04_15_19_06_10->fill(calcptrel(j,p),weight);
00472                 _h_rdA04_15_19_06_10->fill(calcr(j,p),weight*calcrweight(j,p));
00473               }
00474             }
00475             if (inRange(jetpt/GeV, 10., 15.)) {
00476               _numjets04_15_19_10_15 += weight;
00477               _h_N04_15_19_10_15->fill(j.particles().size(),weight);
00478               foreach (const Particle& p, j.particles()) {
00479                 _h_z04_15_19_10_15->fill(calcz(j,p),weight);
00480                 _h_ptrel04_15_19_10_15->fill(calcptrel(j,p),weight);
00481                 _h_rdA04_15_19_10_15->fill(calcr(j,p),weight*calcrweight(j,p));
00482               }
00483             }
00484             if (inRange(jetpt/GeV, 15., 24.)) {
00485               _numjets04_15_19_15_24 += weight;
00486               _h_N04_15_19_15_24->fill(j.particles().size(),weight);
00487               foreach (const Particle& p, j.particles()) {
00488                 _h_z04_15_19_15_24->fill(calcz(j,p),weight);
00489                 _h_ptrel04_15_19_15_24->fill(calcptrel(j,p),weight);
00490                 _h_rdA04_15_19_15_24->fill(calcr(j,p),weight*calcrweight(j,p));
00491               }
00492             }
00493             if (inRange(jetpt/GeV, 24., 40.)) {
00494               _numjets04_15_19_24_40 += weight;
00495               _h_N04_15_19_24_40->fill(j.particles().size(),weight);
00496               foreach (const Particle& p, j.particles()) {
00497                 _h_z04_15_19_24_40->fill(calcz(j,p),weight);
00498                 _h_ptrel04_15_19_24_40->fill(calcptrel(j,p),weight);
00499                 _h_rdA04_15_19_24_40->fill(calcr(j,p),weight*calcrweight(j,p));
00500               }
00501             }
00502           } // 1.5 < rapidity < 1.9
00503           if (fabs(j.rapidity()) < 1.9) {
00504             if (inRange(jetpt/GeV, 4., 6.)) {
00505               _numjets04_00_19_04_06 += weight;
00506               _h_N04_00_19_04_06->fill(j.particles().size(),weight);
00507               foreach (const Particle& p, j.particles()) {
00508                 _h_z04_00_19_04_06->fill(calcz(j,p),weight);
00509                 _h_ptrel04_00_19_04_06->fill(calcptrel(j,p),weight);
00510                 _h_rdA04_00_19_04_06->fill(calcr(j,p),weight*calcrweight(j,p));
00511               }
00512             }
00513             if (inRange(jetpt/GeV, 6., 10.)) {
00514               _numjets04_00_19_06_10 += weight;
00515               _h_N04_00_19_06_10->fill(j.particles().size(),weight);
00516               foreach (const Particle& p, j.particles()) {
00517                 _h_z04_00_19_06_10->fill(calcz(j,p),weight);
00518                 _h_ptrel04_00_19_06_10->fill(calcptrel(j,p),weight);
00519                 _h_rdA04_00_19_06_10->fill(calcr(j,p),weight*calcrweight(j,p));
00520               }
00521             }
00522             if (inRange(jetpt/GeV, 10., 15.)) {
00523               _numjets04_00_19_10_15 += weight;
00524               _h_N04_00_19_10_15->fill(j.particles().size(),weight);
00525               foreach (const Particle& p, j.particles()) {
00526                 _h_z04_00_19_10_15->fill(calcz(j,p),weight);
00527                 _h_ptrel04_00_19_10_15->fill(calcptrel(j,p),weight);
00528                 _h_rdA04_00_19_10_15->fill(calcr(j,p),weight*calcrweight(j,p));
00529               }
00530             }
00531             if (inRange(jetpt/GeV, 15., 24.)) {
00532               _numjets04_00_19_15_24 += weight;
00533               _h_N04_00_19_15_24->fill(j.particles().size(),weight);
00534               foreach (const Particle& p, j.particles()) {
00535                 _h_z04_00_19_15_24->fill(calcz(j,p),weight);
00536                 _h_ptrel04_00_19_15_24->fill(calcptrel(j,p),weight);
00537                 _h_rdA04_00_19_15_24->fill(calcr(j,p),weight*calcrweight(j,p));
00538               }
00539             }
00540             if (inRange(jetpt/GeV, 24., 40.)) {
00541               _numjets04_00_19_24_40 += weight;
00542               _h_N04_00_19_24_40->fill(j.particles().size(),weight);
00543               foreach (const Particle& p, j.particles()) {
00544                 _h_z04_00_19_24_40->fill(calcz(j,p),weight);
00545                 _h_ptrel04_00_19_24_40->fill(calcptrel(j,p),weight);
00546                 _h_rdA04_00_19_24_40->fill(calcr(j,p),weight*calcrweight(j,p));
00547               }
00548             }
00549           } // 0.0 < rapidity < 1.9
00550         } // each jet
00551       } // jets04 not empty
00552 
00553       const Jets& jets06 = applyProjection<JetAlg>(event, "Jets06").jets();
00554       if (!jets06.empty()) {
00555         _sumofweights06 += weight;
00556         foreach (const Jet& j, jets06) {
00557           const double jetpt = j.pT();
00558           if (fabs(j.rapidity()) < 0.5) {
00559             _h_pt06_00_05->fill(jetpt/GeV, weight);
00560             if (inRange(jetpt/GeV, 4., 6.)) {
00561               _numjets06_00_05_04_06 += weight;
00562               _h_N06_00_05_04_06->fill(j.particles().size(),weight);
00563               foreach (const Particle& p, j.particles()) {
00564                 _h_z06_00_05_04_06->fill(calcz(j,p),weight);
00565                 _h_ptrel06_00_05_04_06->fill(calcptrel(j,p),weight);
00566                 _h_rdA06_00_05_04_06->fill(calcr(j,p),weight*calcrweight(j,p));
00567               }
00568             }
00569             if (inRange(jetpt/GeV, 6., 10.)) {
00570               _numjets06_00_05_06_10 += weight;
00571               _h_N06_00_05_06_10->fill(j.particles().size(),weight);
00572               foreach (const Particle& p, j.particles()) {
00573                 _h_z06_00_05_06_10->fill(calcz(j,p),weight);
00574                 _h_ptrel06_00_05_06_10->fill(calcptrel(j,p),weight);
00575                 _h_rdA06_00_05_06_10->fill(calcr(j,p),weight*calcrweight(j,p));
00576               }
00577             }
00578             if (inRange(jetpt/GeV, 10., 15.)) {
00579               _numjets06_00_05_10_15 += weight;
00580               _h_N06_00_05_10_15->fill(j.particles().size(),weight);
00581               foreach (const Particle& p, j.particles()) {
00582                 _h_z06_00_05_10_15->fill(calcz(j,p),weight);
00583                 _h_ptrel06_00_05_10_15->fill(calcptrel(j,p),weight);
00584                 _h_rdA06_00_05_10_15->fill(calcr(j,p),weight*calcrweight(j,p));
00585               }
00586             }
00587             if (inRange(jetpt/GeV, 15., 24.)) {
00588               _numjets06_00_05_15_24 += weight;
00589               _h_N06_00_05_15_24->fill(j.particles().size(),weight);
00590               foreach (const Particle& p, j.particles()) {
00591                 _h_z06_00_05_15_24->fill(calcz(j,p),weight);
00592                 _h_ptrel06_00_05_15_24->fill(calcptrel(j,p),weight);
00593                 _h_rdA06_00_05_15_24->fill(calcr(j,p),weight*calcrweight(j,p));
00594               }
00595             }
00596             if (inRange(jetpt/GeV, 24., 40.)) {
00597               _numjets06_00_05_24_40 += weight;
00598               _h_N06_00_05_24_40->fill(j.particles().size(),weight);
00599               foreach (const Particle& p, j.particles()) {
00600                 _h_z06_00_05_24_40->fill(calcz(j,p),weight);
00601                 _h_ptrel06_00_05_24_40->fill(calcptrel(j,p),weight);
00602                 _h_rdA06_00_05_24_40->fill(calcr(j,p),weight*calcrweight(j,p));
00603               }
00604             }
00605           }
00606           if (fabs(j.rapidity()) > 0.5 && fabs(j.rapidity()) < 1.0) {
00607             _h_pt06_05_10->fill(jetpt/GeV, weight);
00608             if (inRange(jetpt/GeV, 4., 6.)) {
00609               _numjets06_05_10_04_06 += weight;
00610               _h_N06_05_10_04_06->fill(j.particles().size(),weight);
00611               foreach (const Particle& p, j.particles()) {
00612                 _h_z06_05_10_04_06->fill(calcz(j,p),weight);
00613                 _h_ptrel06_05_10_04_06->fill(calcptrel(j,p),weight);
00614                 _h_rdA06_05_10_04_06->fill(calcr(j,p),weight*calcrweight(j,p));
00615               }
00616             }
00617             if (inRange(jetpt/GeV, 6., 10.)) {
00618               _numjets06_05_10_06_10 += weight;
00619               _h_N06_05_10_06_10->fill(j.particles().size(),weight);
00620               foreach (const Particle& p, j.particles()) {
00621                 _h_z06_05_10_06_10->fill(calcz(j,p),weight);
00622                 _h_ptrel06_05_10_06_10->fill(calcptrel(j,p),weight);
00623                 _h_rdA06_05_10_06_10->fill(calcr(j,p),weight*calcrweight(j,p));
00624               }
00625             }
00626             if (inRange(jetpt/GeV, 10., 15.)) {
00627               _numjets06_05_10_10_15 += weight;
00628               _h_N06_05_10_10_15->fill(j.particles().size(),weight);
00629               foreach (const Particle& p, j.particles()) {
00630                 _h_z06_05_10_10_15->fill(calcz(j,p),weight);
00631                 _h_ptrel06_05_10_10_15->fill(calcptrel(j,p),weight);
00632                 _h_rdA06_05_10_10_15->fill(calcr(j,p),weight*calcrweight(j,p));
00633               }
00634             }
00635             if (inRange(jetpt/GeV, 15., 24.)) {
00636               _numjets06_05_10_15_24 += weight;
00637               _h_N06_05_10_15_24->fill(j.particles().size(),weight);
00638               foreach (const Particle& p, j.particles()) {
00639                 _h_z06_05_10_15_24->fill(calcz(j,p),weight);
00640                 _h_ptrel06_05_10_15_24->fill(calcptrel(j,p),weight);
00641                 _h_rdA06_05_10_15_24->fill(calcr(j,p),weight*calcrweight(j,p));
00642               }
00643             }
00644             if (inRange(jetpt/GeV, 24., 40.)) {
00645               _numjets06_05_10_24_40 += weight;
00646               _h_N06_05_10_24_40->fill(j.particles().size(),weight);
00647               foreach (const Particle& p, j.particles()) {
00648                 _h_z06_05_10_24_40->fill(calcz(j,p),weight);
00649                 _h_ptrel06_05_10_24_40->fill(calcptrel(j,p),weight);
00650                 _h_rdA06_05_10_24_40->fill(calcr(j,p),weight*calcrweight(j,p));
00651               }
00652             }
00653           }
00654           if (fabs(j.rapidity()) > 1.0 && fabs(j.rapidity()) < 1.5) {
00655             _h_pt06_10_15->fill(jetpt/GeV, weight);
00656             if (inRange(jetpt/GeV, 4., 6.)) {
00657               _numjets06_10_15_04_06 += weight;
00658               _h_N06_10_15_04_06->fill(j.particles().size(),weight);
00659               foreach (const Particle& p, j.particles()) {
00660                 _h_z06_10_15_04_06->fill(calcz(j,p),weight);
00661                 _h_ptrel06_10_15_04_06->fill(calcptrel(j,p),weight);
00662                 _h_rdA06_10_15_04_06->fill(calcr(j,p),weight*calcrweight(j,p));
00663               }
00664             }
00665             if (inRange(jetpt/GeV, 6., 10.)) {
00666               _numjets06_10_15_06_10 += weight;
00667               _h_N06_10_15_06_10->fill(j.particles().size(),weight);
00668               foreach (const Particle& p, j.particles()) {
00669                 _h_z06_10_15_06_10->fill(calcz(j,p),weight);
00670                 _h_ptrel06_10_15_06_10->fill(calcptrel(j,p),weight);
00671                 _h_rdA06_10_15_06_10->fill(calcr(j,p),weight*calcrweight(j,p));
00672               }
00673             }
00674             if (inRange(jetpt/GeV, 10., 15.)) {
00675               _numjets06_10_15_10_15 += weight;
00676               _h_N06_10_15_10_15->fill(j.particles().size(),weight);
00677               foreach (const Particle& p, j.particles()) {
00678                 _h_z06_10_15_10_15->fill(calcz(j,p),weight);
00679                 _h_ptrel06_10_15_10_15->fill(calcptrel(j,p),weight);
00680                 _h_rdA06_10_15_10_15->fill(calcr(j,p),weight*calcrweight(j,p));
00681               }
00682             }
00683             if (inRange(jetpt/GeV, 15., 24.)) {
00684               _numjets06_10_15_15_24 += weight;
00685               _h_N06_10_15_15_24->fill(j.particles().size(),weight);
00686               foreach (const Particle& p, j.particles()) {
00687                 _h_z06_10_15_15_24->fill(calcz(j,p),weight);
00688                 _h_ptrel06_10_15_15_24->fill(calcptrel(j,p),weight);
00689                 _h_rdA06_10_15_15_24->fill(calcr(j,p),weight*calcrweight(j,p));
00690               }
00691             }
00692             if (inRange(jetpt/GeV, 24., 40.)) {
00693               _numjets06_10_15_24_40 += weight;
00694               _h_N06_10_15_24_40->fill(j.particles().size(),weight);
00695               foreach (const Particle& p, j.particles()) {
00696                 _h_z06_10_15_24_40->fill(calcz(j,p),weight);
00697                 _h_ptrel06_10_15_24_40->fill(calcptrel(j,p),weight);
00698                 _h_rdA06_10_15_24_40->fill(calcr(j,p),weight*calcrweight(j,p));
00699               }
00700             }
00701           }
00702           if (fabs(j.rapidity()) > 1.5 && fabs(j.rapidity()) < 1.9) {
00703             _h_pt06_15_19->fill(jetpt/GeV, weight);
00704             if (inRange(jetpt/GeV, 4., 6.)) {
00705               _numjets06_15_19_04_06 += weight;
00706               _h_N06_15_19_04_06->fill(j.particles().size(),weight);
00707               foreach (const Particle& p, j.particles()) {
00708                 _h_z06_15_19_04_06->fill(calcz(j,p),weight);
00709                 _h_ptrel06_15_19_04_06->fill(calcptrel(j,p),weight);
00710                 _h_rdA06_15_19_04_06->fill(calcr(j,p),weight*calcrweight(j,p));
00711               }
00712             }
00713             if (inRange(jetpt/GeV, 6., 10.)) {
00714               _numjets06_15_19_06_10 += weight;
00715               _h_N06_15_19_06_10->fill(j.particles().size(),weight);
00716               foreach (const Particle& p, j.particles()) {
00717                 _h_z06_15_19_06_10->fill(calcz(j,p),weight);
00718                 _h_ptrel06_15_19_06_10->fill(calcptrel(j,p),weight);
00719                 _h_rdA06_15_19_06_10->fill(calcr(j,p),weight*calcrweight(j,p));
00720               }
00721             }
00722             if (inRange(jetpt/GeV, 10., 15.)) {
00723               _numjets06_15_19_10_15 += weight;
00724               _h_N06_15_19_10_15->fill(j.particles().size(),weight);
00725               foreach (const Particle& p, j.particles()) {
00726                 _h_z06_15_19_10_15->fill(calcz(j,p),weight);
00727                 _h_ptrel06_15_19_10_15->fill(calcptrel(j,p),weight);
00728                 _h_rdA06_15_19_10_15->fill(calcr(j,p),weight*calcrweight(j,p));
00729               }
00730             }
00731             if (inRange(jetpt/GeV, 15., 24.)) {
00732               _numjets06_15_19_15_24 += weight;
00733               _h_N06_15_19_15_24->fill(j.particles().size(),weight);
00734               foreach (const Particle& p, j.particles()) {
00735                 _h_z06_15_19_15_24->fill(calcz(j,p),weight);
00736                 _h_ptrel06_15_19_15_24->fill(calcptrel(j,p),weight);
00737                 _h_rdA06_15_19_15_24->fill(calcr(j,p),weight*calcrweight(j,p));
00738               }
00739             }
00740             if (inRange(jetpt/GeV, 24., 40.)) {
00741               _numjets06_15_19_24_40 += weight;
00742               _h_N06_15_19_24_40->fill(j.particles().size(),weight);
00743               foreach (const Particle& p, j.particles()) {
00744                 _h_z06_15_19_24_40->fill(calcz(j,p),weight);
00745                 _h_ptrel06_15_19_24_40->fill(calcptrel(j,p),weight);
00746                 _h_rdA06_15_19_24_40->fill(calcr(j,p),weight*calcrweight(j,p));
00747               }
00748             }
00749           } // 1.5 < rapidity < 1.9
00750           if (fabs(j.rapidity()) < 1.9) {
00751             if (inRange(jetpt/GeV, 4., 6.)) {
00752               _numjets06_00_19_04_06 += weight;
00753               _h_N06_00_19_04_06->fill(j.particles().size(),weight);
00754               foreach (const Particle& p, j.particles()) {
00755                 _h_z06_00_19_04_06->fill(calcz(j,p),weight);
00756                 _h_ptrel06_00_19_04_06->fill(calcptrel(j,p),weight);
00757                 _h_rdA06_00_19_04_06->fill(calcr(j,p),weight*calcrweight(j,p));
00758               }
00759             }
00760             if (inRange(jetpt/GeV, 6., 10.)) {
00761               _numjets06_00_19_06_10 += weight;
00762               _h_N06_00_19_06_10->fill(j.particles().size(),weight);
00763               foreach (const Particle& p, j.particles()) {
00764                 _h_z06_00_19_06_10->fill(calcz(j,p),weight);
00765                 _h_ptrel06_00_19_06_10->fill(calcptrel(j,p),weight);
00766                 _h_rdA06_00_19_06_10->fill(calcr(j,p),weight*calcrweight(j,p));
00767               }
00768             }
00769             if (inRange(jetpt/GeV, 10., 15.)) {
00770               _numjets06_00_19_10_15 += weight;
00771               _h_N06_00_19_10_15->fill(j.particles().size(),weight);
00772               foreach (const Particle& p, j.particles()) {
00773                 _h_z06_00_19_10_15->fill(calcz(j,p),weight);
00774                 _h_ptrel06_00_19_10_15->fill(calcptrel(j,p),weight);
00775                 _h_rdA06_00_19_10_15->fill(calcr(j,p),weight*calcrweight(j,p));
00776               }
00777             }
00778             if (inRange(jetpt/GeV, 15., 24.)) {
00779               _numjets06_00_19_15_24 += weight;
00780               _h_N06_00_19_15_24->fill(j.particles().size(),weight);
00781               foreach (const Particle& p, j.particles()) {
00782                 _h_z06_00_19_15_24->fill(calcz(j,p),weight);
00783                 _h_ptrel06_00_19_15_24->fill(calcptrel(j,p),weight);
00784                 _h_rdA06_00_19_15_24->fill(calcr(j,p),weight*calcrweight(j,p));
00785               }
00786             }
00787             if (inRange(jetpt/GeV, 24., 40.)) {
00788               _numjets06_00_19_24_40 += weight;
00789               _h_N06_00_19_24_40->fill(j.particles().size(),weight);
00790               foreach (const Particle& p, j.particles()) {
00791                 _h_z06_00_19_24_40->fill(calcz(j,p),weight);
00792                 _h_ptrel06_00_19_24_40->fill(calcptrel(j,p),weight);
00793                 _h_rdA06_00_19_24_40->fill(calcr(j,p),weight*calcrweight(j,p));
00794               }
00795             }
00796           }
00797         } // each jet
00798       } // jets06 not empty
00799 
00800     } // end of event
00801 
00802 
00803     /// Normalise histograms etc., after the run
00804     void finalize() {
00805 
00806       // pT histos: d2sigma_jet / deta dpT
00807       const double xsec = crossSection()/microbarn;
00808       safeinvscale(_h_pt04_00_05, _sumofweights04*(2*0.5)/xsec);
00809       safeinvscale(_h_pt06_00_05, _sumofweights06*(2*0.5)/xsec);
00810       safeinvscale(_h_pt04_05_10, _sumofweights04*(2*0.5)/xsec);
00811       safeinvscale(_h_pt06_05_10, _sumofweights06*(2*0.5)/xsec);
00812       safeinvscale(_h_pt04_10_15, _sumofweights04*(2*0.5)/xsec);
00813       safeinvscale(_h_pt06_10_15, _sumofweights06*(2*0.5)/xsec);
00814       safeinvscale(_h_pt04_15_19, _sumofweights04*(2*0.4)/xsec);
00815       safeinvscale(_h_pt06_15_19, _sumofweights06*(2*0.4)/xsec);
00816 
00817       // N histos: 1/N_jet dN_jet / dN^{ch}_jet
00818       safeinvscale(_h_N04_00_05_04_06, _numjets04_00_05_04_06);
00819       safeinvscale(_h_N06_00_05_04_06, _numjets06_00_05_04_06);
00820       safeinvscale(_h_N04_00_05_06_10, _numjets04_00_05_06_10);
00821       safeinvscale(_h_N06_00_05_06_10, _numjets06_00_05_06_10);
00822       safeinvscale(_h_N04_00_05_10_15, _numjets04_00_05_10_15);
00823       safeinvscale(_h_N06_00_05_10_15, _numjets06_00_05_10_15);
00824       safeinvscale(_h_N04_00_05_15_24, _numjets04_00_05_15_24);
00825       safeinvscale(_h_N06_00_05_15_24, _numjets06_00_05_15_24);
00826       safeinvscale(_h_N04_00_05_24_40, _numjets04_00_05_24_40);
00827       safeinvscale(_h_N06_00_05_24_40, _numjets06_00_05_24_40);
00828       safeinvscale(_h_N04_05_10_04_06, _numjets04_05_10_04_06);
00829       safeinvscale(_h_N06_05_10_04_06, _numjets06_05_10_04_06);
00830       safeinvscale(_h_N04_05_10_06_10, _numjets04_05_10_06_10);
00831       safeinvscale(_h_N06_05_10_06_10, _numjets06_05_10_06_10);
00832       safeinvscale(_h_N04_05_10_10_15, _numjets04_05_10_10_15);
00833       safeinvscale(_h_N06_05_10_10_15, _numjets06_05_10_10_15);
00834       safeinvscale(_h_N04_05_10_15_24, _numjets04_05_10_15_24);
00835       safeinvscale(_h_N06_05_10_15_24, _numjets06_05_10_15_24);
00836       safeinvscale(_h_N04_05_10_24_40, _numjets04_05_10_24_40);
00837       safeinvscale(_h_N06_05_10_24_40, _numjets06_05_10_24_40);
00838       safeinvscale(_h_N04_10_15_04_06, _numjets04_10_15_04_06);
00839       safeinvscale(_h_N06_10_15_04_06, _numjets06_10_15_04_06);
00840       safeinvscale(_h_N04_10_15_06_10, _numjets04_10_15_06_10);
00841       safeinvscale(_h_N06_10_15_06_10, _numjets06_10_15_06_10);
00842       safeinvscale(_h_N04_10_15_10_15, _numjets04_10_15_10_15);
00843       safeinvscale(_h_N06_10_15_10_15, _numjets06_10_15_10_15);
00844       safeinvscale(_h_N04_10_15_15_24, _numjets04_10_15_15_24);
00845       safeinvscale(_h_N06_10_15_15_24, _numjets06_10_15_15_24);
00846       safeinvscale(_h_N04_10_15_24_40, _numjets04_10_15_24_40);
00847       safeinvscale(_h_N06_10_15_24_40, _numjets06_10_15_24_40);
00848       safeinvscale(_h_N04_15_19_04_06, _numjets04_15_19_04_06);
00849       safeinvscale(_h_N06_15_19_04_06, _numjets06_15_19_04_06);
00850       safeinvscale(_h_N04_15_19_06_10, _numjets04_15_19_06_10);
00851       safeinvscale(_h_N06_15_19_06_10, _numjets06_15_19_06_10);
00852       safeinvscale(_h_N04_15_19_10_15, _numjets04_15_19_10_15);
00853       safeinvscale(_h_N06_15_19_10_15, _numjets06_15_19_10_15);
00854       safeinvscale(_h_N04_15_19_15_24, _numjets04_15_19_15_24);
00855       safeinvscale(_h_N06_15_19_15_24, _numjets06_15_19_15_24);
00856       safeinvscale(_h_N04_15_19_24_40, _numjets04_15_19_24_40);
00857       safeinvscale(_h_N06_15_19_24_40, _numjets06_15_19_24_40);
00858       safeinvscale(_h_N04_00_19_04_06, _numjets04_00_19_04_06);
00859       safeinvscale(_h_N06_00_19_04_06, _numjets06_00_19_04_06);
00860       safeinvscale(_h_N04_00_19_06_10, _numjets04_00_19_06_10);
00861       safeinvscale(_h_N06_00_19_06_10, _numjets06_00_19_06_10);
00862       safeinvscale(_h_N04_00_19_10_15, _numjets04_00_19_10_15);
00863       safeinvscale(_h_N06_00_19_10_15, _numjets06_00_19_10_15);
00864       safeinvscale(_h_N04_00_19_15_24, _numjets04_00_19_15_24);
00865       safeinvscale(_h_N06_00_19_15_24, _numjets06_00_19_15_24);
00866       safeinvscale(_h_N04_00_19_24_40, _numjets04_00_19_24_40);
00867       safeinvscale(_h_N06_00_19_24_40, _numjets06_00_19_24_40);
00868 
00869       // z histos: 1/N_jet dN_track / dz_track
00870       safeinvscale(_h_z04_00_05_04_06, _numjets04_00_05_04_06);
00871       safeinvscale(_h_z06_00_05_04_06, _numjets06_00_05_04_06);
00872       safeinvscale(_h_z04_00_05_06_10, _numjets04_00_05_06_10);
00873       safeinvscale(_h_z06_00_05_06_10, _numjets06_00_05_06_10);
00874       safeinvscale(_h_z04_00_05_10_15, _numjets04_00_05_10_15);
00875       safeinvscale(_h_z06_00_05_10_15, _numjets06_00_05_10_15);
00876       safeinvscale(_h_z04_00_05_15_24, _numjets04_00_05_15_24);
00877       safeinvscale(_h_z06_00_05_15_24, _numjets06_00_05_15_24);
00878       safeinvscale(_h_z04_00_05_24_40, _numjets04_00_05_24_40);
00879       safeinvscale(_h_z06_00_05_24_40, _numjets06_00_05_24_40);
00880       safeinvscale(_h_z04_05_10_04_06, _numjets04_05_10_04_06);
00881       safeinvscale(_h_z06_05_10_04_06, _numjets06_05_10_04_06);
00882       safeinvscale(_h_z04_05_10_06_10, _numjets04_05_10_06_10);
00883       safeinvscale(_h_z06_05_10_06_10, _numjets06_05_10_06_10);
00884       safeinvscale(_h_z04_05_10_10_15, _numjets04_05_10_10_15);
00885       safeinvscale(_h_z06_05_10_10_15, _numjets06_05_10_10_15);
00886       safeinvscale(_h_z04_05_10_15_24, _numjets04_05_10_15_24);
00887       safeinvscale(_h_z06_05_10_15_24, _numjets06_05_10_15_24);
00888       safeinvscale(_h_z04_05_10_24_40, _numjets04_05_10_24_40);
00889       safeinvscale(_h_z06_05_10_24_40, _numjets06_05_10_24_40);
00890       safeinvscale(_h_z04_10_15_04_06, _numjets04_10_15_04_06);
00891       safeinvscale(_h_z06_10_15_04_06, _numjets06_10_15_04_06);
00892       safeinvscale(_h_z04_10_15_06_10, _numjets04_10_15_06_10);
00893       safeinvscale(_h_z06_10_15_06_10, _numjets06_10_15_06_10);
00894       safeinvscale(_h_z04_10_15_10_15, _numjets04_10_15_10_15);
00895       safeinvscale(_h_z06_10_15_10_15, _numjets06_10_15_10_15);
00896       safeinvscale(_h_z04_10_15_15_24, _numjets04_10_15_15_24);
00897       safeinvscale(_h_z06_10_15_15_24, _numjets06_10_15_15_24);
00898       safeinvscale(_h_z04_10_15_24_40, _numjets04_10_15_24_40);
00899       safeinvscale(_h_z06_10_15_24_40, _numjets06_10_15_24_40);
00900       safeinvscale(_h_z04_15_19_04_06, _numjets04_15_19_04_06);
00901       safeinvscale(_h_z06_15_19_04_06, _numjets06_15_19_04_06);
00902       safeinvscale(_h_z04_15_19_06_10, _numjets04_15_19_06_10);
00903       safeinvscale(_h_z06_15_19_06_10, _numjets06_15_19_06_10);
00904       safeinvscale(_h_z04_15_19_10_15, _numjets04_15_19_10_15);
00905       safeinvscale(_h_z06_15_19_10_15, _numjets06_15_19_10_15);
00906       safeinvscale(_h_z04_15_19_15_24, _numjets04_15_19_15_24);
00907       safeinvscale(_h_z06_15_19_15_24, _numjets06_15_19_15_24);
00908       safeinvscale(_h_z04_15_19_24_40, _numjets04_15_19_24_40);
00909       safeinvscale(_h_z06_15_19_24_40, _numjets06_15_19_24_40);
00910       safeinvscale(_h_z04_00_19_04_06, _numjets04_00_19_04_06);
00911       safeinvscale(_h_z06_00_19_04_06, _numjets06_00_19_04_06);
00912       safeinvscale(_h_z04_00_19_06_10, _numjets04_00_19_06_10);
00913       safeinvscale(_h_z06_00_19_06_10, _numjets06_00_19_06_10);
00914       safeinvscale(_h_z04_00_19_10_15, _numjets04_00_19_10_15);
00915       safeinvscale(_h_z06_00_19_10_15, _numjets06_00_19_10_15);
00916       safeinvscale(_h_z04_00_19_15_24, _numjets04_00_19_15_24);
00917       safeinvscale(_h_z06_00_19_15_24, _numjets06_00_19_15_24);
00918       safeinvscale(_h_z04_00_19_24_40, _numjets04_00_19_24_40);
00919       safeinvscale(_h_z06_00_19_24_40, _numjets06_00_19_24_40);
00920 
00921       // ptrel histos: 1/N_jet dN_track / dptrel_track
00922       safeinvscale(_h_ptrel04_00_05_04_06, _numjets04_00_05_04_06);
00923       safeinvscale(_h_ptrel06_00_05_04_06, _numjets06_00_05_04_06);
00924       safeinvscale(_h_ptrel04_00_05_06_10, _numjets04_00_05_06_10);
00925       safeinvscale(_h_ptrel06_00_05_06_10, _numjets06_00_05_06_10);
00926       safeinvscale(_h_ptrel04_00_05_10_15, _numjets04_00_05_10_15);
00927       safeinvscale(_h_ptrel06_00_05_10_15, _numjets06_00_05_10_15);
00928       safeinvscale(_h_ptrel04_00_05_15_24, _numjets04_00_05_15_24);
00929       safeinvscale(_h_ptrel06_00_05_15_24, _numjets06_00_05_15_24);
00930       safeinvscale(_h_ptrel04_00_05_24_40, _numjets04_00_05_24_40);
00931       safeinvscale(_h_ptrel06_00_05_24_40, _numjets06_00_05_24_40);
00932       safeinvscale(_h_ptrel04_05_10_04_06, _numjets04_05_10_04_06);
00933       safeinvscale(_h_ptrel06_05_10_04_06, _numjets06_05_10_04_06);
00934       safeinvscale(_h_ptrel04_05_10_06_10, _numjets04_05_10_06_10);
00935       safeinvscale(_h_ptrel06_05_10_06_10, _numjets06_05_10_06_10);
00936       safeinvscale(_h_ptrel04_05_10_10_15, _numjets04_05_10_10_15);
00937       safeinvscale(_h_ptrel06_05_10_10_15, _numjets06_05_10_10_15);
00938       safeinvscale(_h_ptrel04_05_10_15_24, _numjets04_05_10_15_24);
00939       safeinvscale(_h_ptrel06_05_10_15_24, _numjets06_05_10_15_24);
00940       safeinvscale(_h_ptrel04_05_10_24_40, _numjets04_05_10_24_40);
00941       safeinvscale(_h_ptrel06_05_10_24_40, _numjets06_05_10_24_40);
00942       safeinvscale(_h_ptrel04_10_15_04_06, _numjets04_10_15_04_06);
00943       safeinvscale(_h_ptrel06_10_15_04_06, _numjets06_10_15_04_06);
00944       safeinvscale(_h_ptrel04_10_15_06_10, _numjets04_10_15_06_10);
00945       safeinvscale(_h_ptrel06_10_15_06_10, _numjets06_10_15_06_10);
00946       safeinvscale(_h_ptrel04_10_15_10_15, _numjets04_10_15_10_15);
00947       safeinvscale(_h_ptrel06_10_15_10_15, _numjets06_10_15_10_15);
00948       safeinvscale(_h_ptrel04_10_15_15_24, _numjets04_10_15_15_24);
00949       safeinvscale(_h_ptrel06_10_15_15_24, _numjets06_10_15_15_24);
00950       safeinvscale(_h_ptrel04_10_15_24_40, _numjets04_10_15_24_40);
00951       safeinvscale(_h_ptrel06_10_15_24_40, _numjets06_10_15_24_40);
00952       safeinvscale(_h_ptrel04_15_19_04_06, _numjets04_15_19_04_06);
00953       safeinvscale(_h_ptrel06_15_19_04_06, _numjets06_15_19_04_06);
00954       safeinvscale(_h_ptrel04_15_19_06_10, _numjets04_15_19_06_10);
00955       safeinvscale(_h_ptrel06_15_19_06_10, _numjets06_15_19_06_10);
00956       safeinvscale(_h_ptrel04_15_19_10_15, _numjets04_15_19_10_15);
00957       safeinvscale(_h_ptrel06_15_19_10_15, _numjets06_15_19_10_15);
00958       safeinvscale(_h_ptrel04_15_19_15_24, _numjets04_15_19_15_24);
00959       safeinvscale(_h_ptrel06_15_19_15_24, _numjets06_15_19_15_24);
00960       safeinvscale(_h_ptrel04_15_19_24_40, _numjets04_15_19_24_40);
00961       safeinvscale(_h_ptrel06_15_19_24_40, _numjets06_15_19_24_40);
00962 
00963       safeinvscale(_h_ptrel04_00_19_04_06, _numjets04_00_19_04_06);
00964       safeinvscale(_h_ptrel06_00_19_04_06, _numjets06_00_19_04_06);
00965       safeinvscale(_h_ptrel04_00_19_06_10, _numjets04_00_19_06_10);
00966       safeinvscale(_h_ptrel06_00_19_06_10, _numjets06_00_19_06_10);
00967       safeinvscale(_h_ptrel04_00_19_10_15, _numjets04_00_19_10_15);
00968       safeinvscale(_h_ptrel06_00_19_10_15, _numjets06_00_19_10_15);
00969       safeinvscale(_h_ptrel04_00_19_15_24, _numjets04_00_19_15_24);
00970       safeinvscale(_h_ptrel06_00_19_15_24, _numjets06_00_19_15_24);
00971       safeinvscale(_h_ptrel04_00_19_24_40, _numjets04_00_19_24_40);
00972       safeinvscale(_h_ptrel06_00_19_24_40, _numjets06_00_19_24_40);
00973 
00974       // r histos: 1/N_jet dN_track / dA
00975       safeinvscale(_h_rdA04_00_05_04_06, _numjets04_00_05_04_06);
00976       safeinvscale(_h_rdA06_00_05_04_06, _numjets06_00_05_04_06);
00977       safeinvscale(_h_rdA04_00_05_06_10, _numjets04_00_05_06_10);
00978       safeinvscale(_h_rdA06_00_05_06_10, _numjets06_00_05_06_10);
00979       safeinvscale(_h_rdA04_00_05_10_15, _numjets04_00_05_10_15);
00980       safeinvscale(_h_rdA06_00_05_10_15, _numjets06_00_05_10_15);
00981       safeinvscale(_h_rdA04_00_05_15_24, _numjets04_00_05_15_24);
00982       safeinvscale(_h_rdA06_00_05_15_24, _numjets06_00_05_15_24);
00983       safeinvscale(_h_rdA04_00_05_24_40, _numjets04_00_05_24_40);
00984       safeinvscale(_h_rdA06_00_05_24_40, _numjets06_00_05_24_40);
00985       safeinvscale(_h_rdA04_05_10_04_06, _numjets04_05_10_04_06);
00986       safeinvscale(_h_rdA06_05_10_04_06, _numjets06_05_10_04_06);
00987       safeinvscale(_h_rdA04_05_10_06_10, _numjets04_05_10_06_10);
00988       safeinvscale(_h_rdA06_05_10_06_10, _numjets06_05_10_06_10);
00989       safeinvscale(_h_rdA04_05_10_10_15, _numjets04_05_10_10_15);
00990       safeinvscale(_h_rdA06_05_10_10_15, _numjets06_05_10_10_15);
00991       safeinvscale(_h_rdA04_05_10_15_24, _numjets04_05_10_15_24);
00992       safeinvscale(_h_rdA06_05_10_15_24, _numjets06_05_10_15_24);
00993       safeinvscale(_h_rdA04_05_10_24_40, _numjets04_05_10_24_40);
00994       safeinvscale(_h_rdA06_05_10_24_40, _numjets06_05_10_24_40);
00995       safeinvscale(_h_rdA04_10_15_04_06, _numjets04_10_15_04_06);
00996       safeinvscale(_h_rdA06_10_15_04_06, _numjets06_10_15_04_06);
00997       safeinvscale(_h_rdA04_10_15_06_10, _numjets04_10_15_06_10);
00998       safeinvscale(_h_rdA06_10_15_06_10, _numjets06_10_15_06_10);
00999       safeinvscale(_h_rdA04_10_15_10_15, _numjets04_10_15_10_15);
01000       safeinvscale(_h_rdA06_10_15_10_15, _numjets06_10_15_10_15);
01001       safeinvscale(_h_rdA04_10_15_15_24, _numjets04_10_15_15_24);
01002       safeinvscale(_h_rdA06_10_15_15_24, _numjets06_10_15_15_24);
01003       safeinvscale(_h_rdA04_10_15_24_40, _numjets04_10_15_24_40);
01004       safeinvscale(_h_rdA06_10_15_24_40, _numjets06_10_15_24_40);
01005       safeinvscale(_h_rdA04_15_19_04_06, _numjets04_15_19_04_06);
01006       safeinvscale(_h_rdA06_15_19_04_06, _numjets06_15_19_04_06);
01007       safeinvscale(_h_rdA04_15_19_06_10, _numjets04_15_19_06_10);
01008       safeinvscale(_h_rdA06_15_19_06_10, _numjets06_15_19_06_10);
01009       safeinvscale(_h_rdA04_15_19_10_15, _numjets04_15_19_10_15);
01010       safeinvscale(_h_rdA06_15_19_10_15, _numjets06_15_19_10_15);
01011       safeinvscale(_h_rdA04_15_19_15_24, _numjets04_15_19_15_24);
01012       safeinvscale(_h_rdA06_15_19_15_24, _numjets06_15_19_15_24);
01013       safeinvscale(_h_rdA04_15_19_24_40, _numjets04_15_19_24_40);
01014       safeinvscale(_h_rdA06_15_19_24_40, _numjets06_15_19_24_40);
01015 
01016       safeinvscale(_h_rdA04_00_19_04_06, _numjets04_00_19_04_06);
01017       safeinvscale(_h_rdA06_00_19_04_06, _numjets06_00_19_04_06);
01018       safeinvscale(_h_rdA04_00_19_06_10, _numjets04_00_19_06_10);
01019       safeinvscale(_h_rdA06_00_19_06_10, _numjets06_00_19_06_10);
01020       safeinvscale(_h_rdA04_00_19_10_15, _numjets04_00_19_10_15);
01021       safeinvscale(_h_rdA06_00_19_10_15, _numjets06_00_19_10_15);
01022       safeinvscale(_h_rdA04_00_19_15_24, _numjets04_00_19_15_24);
01023       safeinvscale(_h_rdA06_00_19_15_24, _numjets06_00_19_15_24);
01024       safeinvscale(_h_rdA04_00_19_24_40, _numjets04_00_19_24_40);
01025       safeinvscale(_h_rdA06_00_19_24_40, _numjets06_00_19_24_40);
01026     }
01027 
01028     //@}
01029 
01030 
01031   private:
01032 
01033     void safeinvscale(Histo1DPtr h, double denom) {
01034       if (denom != 0) {
01035         scale(h, 1.0/denom);
01036       } else {
01037         normalize(h, 0);
01038       }
01039     }
01040 
01041 
01042     /// Event weights
01043     double _sumofweights04, _sumofweights06;
01044 
01045 
01046     /// Jet counters
01047     double _numjets04_00_05_04_06, _numjets04_00_05_06_10, _numjets04_00_05_10_15, _numjets04_00_05_15_24, _numjets04_00_05_24_40;
01048     double _numjets06_00_05_04_06, _numjets06_00_05_06_10, _numjets06_00_05_10_15, _numjets06_00_05_15_24, _numjets06_00_05_24_40;
01049     double _numjets04_05_10_04_06, _numjets04_05_10_06_10, _numjets04_05_10_10_15, _numjets04_05_10_15_24, _numjets04_05_10_24_40;
01050     double _numjets06_05_10_04_06, _numjets06_05_10_06_10, _numjets06_05_10_10_15, _numjets06_05_10_15_24, _numjets06_05_10_24_40;
01051     double _numjets04_10_15_04_06, _numjets04_10_15_06_10, _numjets04_10_15_10_15, _numjets04_10_15_15_24, _numjets04_10_15_24_40;
01052     double _numjets06_10_15_04_06, _numjets06_10_15_06_10, _numjets06_10_15_10_15, _numjets06_10_15_15_24, _numjets06_10_15_24_40;
01053     double _numjets04_15_19_04_06, _numjets04_15_19_06_10, _numjets04_15_19_10_15, _numjets04_15_19_15_24, _numjets04_15_19_24_40;
01054     double _numjets06_15_19_04_06, _numjets06_15_19_06_10, _numjets06_15_19_10_15, _numjets06_15_19_15_24, _numjets06_15_19_24_40;
01055     double _numjets04_00_19_04_06, _numjets04_00_19_06_10, _numjets04_00_19_10_15, _numjets04_00_19_15_24, _numjets04_00_19_24_40;
01056     double _numjets06_00_19_04_06, _numjets06_00_19_06_10, _numjets06_00_19_10_15, _numjets06_00_19_15_24, _numjets06_00_19_24_40;
01057 
01058 
01059   private:
01060 
01061     /// @name Histograms
01062     //@{
01063 
01064     Histo1DPtr _h_pt04_00_05, _h_pt06_00_05;
01065     Histo1DPtr _h_N04_00_05_04_06, _h_N06_00_05_04_06;
01066     Histo1DPtr _h_N04_00_05_06_10, _h_N06_00_05_06_10;
01067     Histo1DPtr _h_N04_00_05_10_15, _h_N06_00_05_10_15;
01068     Histo1DPtr _h_N04_00_05_15_24, _h_N06_00_05_15_24;
01069     Histo1DPtr _h_N04_00_05_24_40, _h_N06_00_05_24_40;
01070     Histo1DPtr _h_z04_00_05_04_06, _h_z06_00_05_04_06;
01071     Histo1DPtr _h_z04_00_05_06_10, _h_z06_00_05_06_10;
01072     Histo1DPtr _h_z04_00_05_10_15, _h_z06_00_05_10_15;
01073     Histo1DPtr _h_z04_00_05_15_24, _h_z06_00_05_15_24;
01074     Histo1DPtr _h_z04_00_05_24_40, _h_z06_00_05_24_40;
01075     Histo1DPtr _h_ptrel04_00_05_04_06, _h_ptrel06_00_05_04_06;
01076     Histo1DPtr _h_ptrel04_00_05_06_10, _h_ptrel06_00_05_06_10;
01077     Histo1DPtr _h_ptrel04_00_05_10_15, _h_ptrel06_00_05_10_15;
01078     Histo1DPtr _h_ptrel04_00_05_15_24, _h_ptrel06_00_05_15_24;
01079     Histo1DPtr _h_ptrel04_00_05_24_40, _h_ptrel06_00_05_24_40;
01080     Histo1DPtr _h_rdA04_00_05_04_06, _h_rdA06_00_05_04_06;
01081     Histo1DPtr _h_rdA04_00_05_06_10, _h_rdA06_00_05_06_10;
01082     Histo1DPtr _h_rdA04_00_05_10_15, _h_rdA06_00_05_10_15;
01083     Histo1DPtr _h_rdA04_00_05_15_24, _h_rdA06_00_05_15_24;
01084     Histo1DPtr _h_rdA04_00_05_24_40, _h_rdA06_00_05_24_40;
01085 
01086     Histo1DPtr _h_pt04_05_10, _h_pt06_05_10;
01087     Histo1DPtr _h_N04_05_10_04_06, _h_N06_05_10_04_06;
01088     Histo1DPtr _h_N04_05_10_06_10, _h_N06_05_10_06_10;
01089     Histo1DPtr _h_N04_05_10_10_15, _h_N06_05_10_10_15;
01090     Histo1DPtr _h_N04_05_10_15_24, _h_N06_05_10_15_24;
01091     Histo1DPtr _h_N04_05_10_24_40, _h_N06_05_10_24_40;
01092     Histo1DPtr _h_z04_05_10_04_06, _h_z06_05_10_04_06;
01093     Histo1DPtr _h_z04_05_10_06_10, _h_z06_05_10_06_10;
01094     Histo1DPtr _h_z04_05_10_10_15, _h_z06_05_10_10_15;
01095     Histo1DPtr _h_z04_05_10_15_24, _h_z06_05_10_15_24;
01096     Histo1DPtr _h_z04_05_10_24_40, _h_z06_05_10_24_40;
01097     Histo1DPtr _h_ptrel04_05_10_04_06, _h_ptrel06_05_10_04_06;
01098     Histo1DPtr _h_ptrel04_05_10_06_10, _h_ptrel06_05_10_06_10;
01099     Histo1DPtr _h_ptrel04_05_10_10_15, _h_ptrel06_05_10_10_15;
01100     Histo1DPtr _h_ptrel04_05_10_15_24, _h_ptrel06_05_10_15_24;
01101     Histo1DPtr _h_ptrel04_05_10_24_40, _h_ptrel06_05_10_24_40;
01102     Histo1DPtr _h_rdA04_05_10_04_06, _h_rdA06_05_10_04_06;
01103     Histo1DPtr _h_rdA04_05_10_06_10, _h_rdA06_05_10_06_10;
01104     Histo1DPtr _h_rdA04_05_10_10_15, _h_rdA06_05_10_10_15;
01105     Histo1DPtr _h_rdA04_05_10_15_24, _h_rdA06_05_10_15_24;
01106     Histo1DPtr _h_rdA04_05_10_24_40, _h_rdA06_05_10_24_40;
01107 
01108     Histo1DPtr _h_pt04_10_15, _h_pt06_10_15;
01109     Histo1DPtr _h_N04_10_15_04_06, _h_N06_10_15_04_06;
01110     Histo1DPtr _h_N04_10_15_06_10, _h_N06_10_15_06_10;
01111     Histo1DPtr _h_N04_10_15_10_15, _h_N06_10_15_10_15;
01112     Histo1DPtr _h_N04_10_15_15_24, _h_N06_10_15_15_24;
01113     Histo1DPtr _h_N04_10_15_24_40, _h_N06_10_15_24_40;
01114     Histo1DPtr _h_z04_10_15_04_06, _h_z06_10_15_04_06;
01115     Histo1DPtr _h_z04_10_15_06_10, _h_z06_10_15_06_10;
01116     Histo1DPtr _h_z04_10_15_10_15, _h_z06_10_15_10_15;
01117     Histo1DPtr _h_z04_10_15_15_24, _h_z06_10_15_15_24;
01118     Histo1DPtr _h_z04_10_15_24_40, _h_z06_10_15_24_40;
01119     Histo1DPtr _h_ptrel04_10_15_04_06, _h_ptrel06_10_15_04_06;
01120     Histo1DPtr _h_ptrel04_10_15_06_10, _h_ptrel06_10_15_06_10;
01121     Histo1DPtr _h_ptrel04_10_15_10_15, _h_ptrel06_10_15_10_15;
01122     Histo1DPtr _h_ptrel04_10_15_15_24, _h_ptrel06_10_15_15_24;
01123     Histo1DPtr _h_ptrel04_10_15_24_40, _h_ptrel06_10_15_24_40;
01124     Histo1DPtr _h_rdA04_10_15_04_06, _h_rdA06_10_15_04_06;
01125     Histo1DPtr _h_rdA04_10_15_06_10, _h_rdA06_10_15_06_10;
01126     Histo1DPtr _h_rdA04_10_15_10_15, _h_rdA06_10_15_10_15;
01127     Histo1DPtr _h_rdA04_10_15_15_24, _h_rdA06_10_15_15_24;
01128     Histo1DPtr _h_rdA04_10_15_24_40, _h_rdA06_10_15_24_40;
01129 
01130     Histo1DPtr _h_pt04_15_19, _h_pt06_15_19;
01131     Histo1DPtr _h_N04_15_19_04_06, _h_N06_15_19_04_06;
01132     Histo1DPtr _h_N04_15_19_06_10, _h_N06_15_19_06_10;
01133     Histo1DPtr _h_N04_15_19_10_15, _h_N06_15_19_10_15;
01134     Histo1DPtr _h_N04_15_19_15_24, _h_N06_15_19_15_24;
01135     Histo1DPtr _h_N04_15_19_24_40, _h_N06_15_19_24_40;
01136     Histo1DPtr _h_z04_15_19_04_06, _h_z06_15_19_04_06;
01137     Histo1DPtr _h_z04_15_19_06_10, _h_z06_15_19_06_10;
01138     Histo1DPtr _h_z04_15_19_10_15, _h_z06_15_19_10_15;
01139     Histo1DPtr _h_z04_15_19_15_24, _h_z06_15_19_15_24;
01140     Histo1DPtr _h_z04_15_19_24_40, _h_z06_15_19_24_40;
01141     Histo1DPtr _h_ptrel04_15_19_04_06, _h_ptrel06_15_19_04_06;
01142     Histo1DPtr _h_ptrel04_15_19_06_10, _h_ptrel06_15_19_06_10;
01143     Histo1DPtr _h_ptrel04_15_19_10_15, _h_ptrel06_15_19_10_15;
01144     Histo1DPtr _h_ptrel04_15_19_15_24, _h_ptrel06_15_19_15_24;
01145     Histo1DPtr _h_ptrel04_15_19_24_40, _h_ptrel06_15_19_24_40;
01146     Histo1DPtr _h_rdA04_15_19_04_06, _h_rdA06_15_19_04_06;
01147     Histo1DPtr _h_rdA04_15_19_06_10, _h_rdA06_15_19_06_10;
01148     Histo1DPtr _h_rdA04_15_19_10_15, _h_rdA06_15_19_10_15;
01149     Histo1DPtr _h_rdA04_15_19_15_24, _h_rdA06_15_19_15_24;
01150     Histo1DPtr _h_rdA04_15_19_24_40, _h_rdA06_15_19_24_40;
01151 
01152     Histo1DPtr _h_N04_00_19_04_06, _h_N06_00_19_04_06;
01153     Histo1DPtr _h_N04_00_19_06_10, _h_N06_00_19_06_10;
01154     Histo1DPtr _h_N04_00_19_10_15, _h_N06_00_19_10_15;
01155     Histo1DPtr _h_N04_00_19_15_24, _h_N06_00_19_15_24;
01156     Histo1DPtr _h_N04_00_19_24_40, _h_N06_00_19_24_40;
01157     Histo1DPtr _h_z04_00_19_04_06, _h_z06_00_19_04_06;
01158     Histo1DPtr _h_z04_00_19_06_10, _h_z06_00_19_06_10;
01159     Histo1DPtr _h_z04_00_19_10_15, _h_z06_00_19_10_15;
01160     Histo1DPtr _h_z04_00_19_15_24, _h_z06_00_19_15_24;
01161     Histo1DPtr _h_z04_00_19_24_40, _h_z06_00_19_24_40;
01162     Histo1DPtr _h_ptrel04_00_19_04_06, _h_ptrel06_00_19_04_06;
01163     Histo1DPtr _h_ptrel04_00_19_06_10, _h_ptrel06_00_19_06_10;
01164     Histo1DPtr _h_ptrel04_00_19_10_15, _h_ptrel06_00_19_10_15;
01165     Histo1DPtr _h_ptrel04_00_19_15_24, _h_ptrel06_00_19_15_24;
01166     Histo1DPtr _h_ptrel04_00_19_24_40, _h_ptrel06_00_19_24_40;
01167     Histo1DPtr _h_rdA04_00_19_04_06, _h_rdA06_00_19_04_06;
01168     Histo1DPtr _h_rdA04_00_19_06_10, _h_rdA06_00_19_06_10;
01169     Histo1DPtr _h_rdA04_00_19_10_15, _h_rdA06_00_19_10_15;
01170     Histo1DPtr _h_rdA04_00_19_15_24, _h_rdA06_00_19_15_24;
01171     Histo1DPtr _h_rdA04_00_19_24_40, _h_rdA06_00_19_24_40;
01172 
01173     //@}
01174 
01175   };
01176 
01177 
01178   // This global object acts as a hook for the plugin system
01179   DECLARE_RIVET_PLUGIN(ATLAS_2011_I919017);
01180 
01181 }