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