rivet is hosted by Hepforge, IPPP Durham
ATLAS_2014_I1306615.cc
Go to the documentation of this file.
00001 /*
00002  * Rivet routine for H->yy differential cross-sections measurement
00003  * arXiv: http://arxiv.org/abs/ARXIV:1407.4222
00004  * HepData: http://hepdata.cedar.ac.uk/view/ins1306615
00005  * Author: Michaela Queitsch-Maitland <michaela.queitsch-maitland@cern.ch>
00006  */
00007 
00008 // -*- C++ -*-
00009 #include "Rivet/Analysis.hh"
00010 #include "Rivet/Projections/DressedLeptons.hh"
00011 #include "Rivet/Projections/VetoedFinalState.hh"
00012 #include "Rivet/Projections/FastJets.hh"
00013 
00014 namespace Rivet {
00015 
00016   class ATLAS_2014_I1306615 : public Analysis {
00017   public:
00018 
00019     // Constructor
00020     ATLAS_2014_I1306615()
00021       : Analysis("ATLAS_2014_I1306615")
00022     {    }
00023 
00024   public:
00025 
00026     // Book histograms and initialise projections before the run
00027     void init() {
00028 
00029       // Final state
00030       // All particles within |eta| < 5.0
00031       const FinalState FS(Cuts::abseta<5.0);
00032       addProjection(FS,"FS");
00033 
00034       // Project photons with pT > 25 GeV and |eta| < 2.37
00035       IdentifiedFinalState ph_FS(Cuts::abseta<2.37 && Cuts::pT>25.0*GeV);
00036       ph_FS.acceptIdPair(PID::PHOTON);
00037       addProjection(ph_FS, "PH_FS");
00038 
00039       // Project photons for dressing
00040       IdentifiedFinalState ph_dressing_FS(FS);
00041       ph_dressing_FS.acceptIdPair(PID::PHOTON);
00042 
00043       // Project bare electrons
00044       IdentifiedFinalState el_bare_FS(FS);
00045       el_bare_FS.acceptIdPair(PID::ELECTRON);
00046       addProjection(el_bare_FS,"el_bare_FS");
00047 
00048       // Project dressed electrons with pT > 15 GeV and |eta| < 2.47
00049       DressedLeptons el_dressed_FS(ph_dressing_FS, el_bare_FS, 0.1, Cuts::abseta < 2.47 && Cuts::pT > 15*GeV, true, false);
00050       addProjection(el_dressed_FS,"EL_DRESSED_FS");
00051 
00052       // Project bare muons
00053       IdentifiedFinalState mu_bare_FS(FS);
00054       mu_bare_FS.acceptIdPair(PID::MUON);
00055 
00056       // Project dressed muons with pT > 15 GeV and |eta| < 2.47
00057       //DressedLeptons mu_dressed_FS(ph_dressing_FS, mu_bare_FS, 0.1, true, -2.47, 2.47, 15.0*GeV, false);
00058       DressedLeptons mu_dressed_FS(ph_dressing_FS, mu_bare_FS, 0.1, Cuts::abseta < 2.47 && Cuts::pT > 15*GeV, true, false);
00059       addProjection(mu_dressed_FS,"MU_DRESSED_FS");
00060 
00061       // Final state excluding muons and neutrinos (for jet building and photon isolation)
00062       VetoedFinalState veto_mu_nu_FS(FS);
00063       veto_mu_nu_FS.vetoNeutrinos();
00064       veto_mu_nu_FS.addVetoPairId(PID::MUON);
00065       addProjection(veto_mu_nu_FS, "VETO_MU_NU_FS");
00066 
00067       // Build the anti-kT R=0.4 jets, using FinalState particles (vetoing muons and neutrinos)
00068       FastJets jets(veto_mu_nu_FS, FastJets::ANTIKT, 0.4);
00069       addProjection(jets, "JETS");
00070 
00071       // Book histograms
00072       // 1D distributions
00073       _h_pT_yy         = bookHisto1D(1,1,1);
00074       _h_y_yy          = bookHisto1D(2,1,1);
00075       _h_Njets30       = bookHisto1D(3,1,1);
00076       _h_Njets50       = bookHisto1D(4,1,1);
00077       _h_pT_j1         = bookHisto1D(5,1,1);
00078       _h_y_j1          = bookHisto1D(6,1,1);
00079       _h_HT            = bookHisto1D(7,1,1);
00080       _h_pT_j2         = bookHisto1D(8,1,1);
00081       _h_Dy_jj         = bookHisto1D(9,1,1);
00082       _h_Dphi_yy_jj    = bookHisto1D(10,1,1);
00083       _h_cosTS_CS      = bookHisto1D(11,1,1);
00084       _h_cosTS_CS_5bin = bookHisto1D(12,1,1);
00085       _h_Dphi_jj       = bookHisto1D(13,1,1);
00086       _h_pTt_yy        = bookHisto1D(14,1,1);
00087       _h_Dy_yy         = bookHisto1D(15,1,1);
00088       _h_tau_jet       = bookHisto1D(16,1,1);
00089       _h_sum_tau_jet   = bookHisto1D(17,1,1);
00090       _h_y_j2          = bookHisto1D(18,1,1);
00091       _h_pT_j3         = bookHisto1D(19,1,1);
00092       _h_m_jj          = bookHisto1D(20,1,1);
00093       _h_pT_yy_jj      = bookHisto1D(21,1,1);
00094 
00095       // 2D distributions of cosTS_CS x pT_yy
00096       _h_cosTS_pTyy_low  = bookHisto1D(22,1,1);
00097       _h_cosTS_pTyy_high = bookHisto1D(22,1,2);
00098       _h_cosTS_pTyy_rest = bookHisto1D(22,1,3);
00099 
00100       // 2D distributions of Njets x pT_yy
00101       _h_pTyy_Njets0 = bookHisto1D(23,1,1);
00102       _h_pTyy_Njets1 = bookHisto1D(23,1,2);
00103       _h_pTyy_Njets2 = bookHisto1D(23,1,3);
00104 
00105       _h_pTj1_excl = bookHisto1D(24,1,1);
00106 
00107       // Fiducial regions
00108       _h_fidXSecs = bookHisto1D(29,1,1);
00109     }
00110 
00111     // Perform the per-event analysis
00112     void analyze(const Event& event) {
00113 
00114       const double weight = event.weight();
00115       _weight = weight;
00116 
00117       // Get final state particles
00118       const ParticleVector& FS_ptcls          = applyProjection<FinalState>(event, "FS").particles();
00119       const ParticleVector& ptcls_veto_mu_nu  = applyProjection<VetoedFinalState>(event, "VETO_MU_NU_FS").particles();
00120       const ParticleVector& photons           = applyProjection<IdentifiedFinalState>(event, "PH_FS").particlesByPt();
00121       const vector<DressedLepton>& el_dressed = applyProjection<DressedLeptons>(event, "EL_DRESSED_FS").dressedLeptons();
00122       const vector<DressedLepton>& mu_dressed = applyProjection<DressedLeptons>(event, "MU_DRESSED_FS").dressedLeptons();
00123 
00124       // For isolation calculation
00125       float dR_iso    = 0.4;
00126       float ETcut_iso = 14.0;
00127       FourMomentum ET_iso;
00128 
00129       // Fiducial selection: pT > 25 GeV, |eta| < 2.37 and isolation (in cone deltaR = 0.4) is < 14 GeV
00130       vector<const Particle*> fid_photons;
00131       foreach (const Particle& ph, photons) {
00132 
00133     // Veto photons from hadron or tau decay
00134     if ( fromHadronDecay(ph) ) continue;
00135 
00136     // Calculate isolation
00137     ET_iso = - ph.momentum();
00138     // Loop over fs truth particles (excluding muons and neutrinos)
00139     foreach (const Particle& p, ptcls_veto_mu_nu) {
00140       // Check if the truth particle is in a cone of 0.4
00141       if ( deltaR(ph.momentum(), p.momentum()) < dR_iso )
00142         ET_iso += p.momentum();
00143     }
00144 
00145     // Check isolation
00146     if ( ET_iso.Et() > ETcut_iso ) continue;
00147 
00148     // Fill vector of photons passing fiducial selection
00149     fid_photons.push_back(&ph);
00150       }
00151 
00152       if(fid_photons.size() < 2) vetoEvent;
00153 
00154       const FourMomentum& y1 = fid_photons[0]->momentum();
00155       const FourMomentum& y2 = fid_photons[1]->momentum();
00156 
00157       double m_yy = (y1 + y2).mass();
00158 
00159       // Relative pT cuts
00160       if ( y1.pT() < 0.35 * m_yy || y2.pT() < 0.25 * m_yy ) vetoEvent;
00161 
00162       // Mass window cut
00163       if ( m_yy < 105 || m_yy > 160 ) vetoEvent;
00164 
00165       // -------------------------------------------- //
00166       // Passed diphoton baseline fiducial selection! //
00167       // -------------------------------------------- //
00168 
00169       // Electron selection
00170       vector<const Particle*> good_el;
00171       foreach(const DressedLepton& els, el_dressed)  {
00172 
00173     const Particle& el = els.constituentLepton();
00174     if ( el.momentum().pT()        < 15   ) continue;
00175     if ( fabs(el.momentum().eta()) > 2.47 ) continue;
00176     if ( deltaR(el.momentum(), y1) < 0.4  ) continue;
00177     if ( deltaR(el.momentum(), y2) < 0.4  ) continue;
00178     if ( fromHadronDecay(el)              ) continue; // Veto electrons from hadron or tau decay
00179     good_el.push_back(&el);
00180       }
00181 
00182       // Muon selection
00183       vector<const Particle*> good_mu;
00184       foreach(const DressedLepton& mus, mu_dressed)  {
00185 
00186     const Particle& mu = mus.constituentLepton();
00187     if ( mu.momentum().pT()        < 15   ) continue;
00188     if ( fabs(mu.momentum().eta()) > 2.47 ) continue;
00189     if ( deltaR(mu.momentum(), y1) < 0.4  ) continue;
00190     if ( deltaR(mu.momentum(), y2) < 0.4  ) continue;
00191     if ( fromHadronDecay(mu)              ) continue; // Veto muons from hadron or tau decay
00192     good_mu.push_back(&mu);
00193       }
00194 
00195       // Find prompt, invisible particles for missing ET calculation
00196       // Based on VisibleFinalState projection
00197       FourMomentum invisible(0,0,0,0);
00198       foreach (const Particle& p, FS_ptcls) {
00199 
00200     // Veto non-prompt particles (from hadron or tau decay)
00201         if ( fromHadronDecay(p) ) continue;
00202     // Charged particles are visible
00203     if ( PID::threeCharge( p.pid() ) != 0 ) continue;
00204     // Neutral hadrons are visible
00205     if ( PID::isHadron( p.pid() ) ) continue;
00206         // Photons are visible
00207     if ( p.pid() == PID::PHOTON ) continue;
00208         // Gluons are visible (for parton level analyses)
00209     if ( p.pid() == PID::GLUON ) continue;
00210     // Everything else is invisible
00211     invisible += p.momentum();
00212       }
00213       double MET = invisible.Et();
00214 
00215       // Jet selection
00216       // Get jets with pT > 25 GeV and |rapidity| < 4.4
00217       //const Jets& jets = applyProjection<FastJets>(event, "JETS").jetsByPt(25.0*GeV, MAXDOUBLE, -4.4, 4.4, RAPIDITY);
00218       const Jets& jets = applyProjection<FastJets>(event, "JETS").jetsByPt(Cuts::pT>25.0*GeV && Cuts::absrap <4.4);
00219 
00220       vector<const Jet*> jets_25;
00221       vector<const Jet*> jets_30;
00222       vector<const Jet*> jets_50;
00223 
00224       foreach (const Jet& jet, jets) {
00225 
00226     bool passOverlap = true;
00227     // Overlap with leading photons
00228     if ( deltaR(y1, jet.momentum()) < 0.4 ) passOverlap = false;
00229     if ( deltaR(y2, jet.momentum()) < 0.4 ) passOverlap = false;
00230 
00231     // Overlap with good electrons
00232     foreach (const Particle* el, good_el)
00233       if ( deltaR(el->momentum(), jet.momentum()) < 0.2 ) passOverlap = false;
00234 
00235     if ( ! passOverlap ) continue;
00236 
00237     if ( fabs(jet.momentum().eta()) < 2.4 || ( fabs(jet.momentum().eta()) > 2.4 && jet.momentum().pT() > 30 ) ) jets_25.push_back(&jet);
00238     if ( jet.momentum().pT() > 30 ) jets_30.push_back(&jet);
00239     if ( jet.momentum().pT() > 50 ) jets_50.push_back(&jet);
00240       }
00241 
00242       // Fiducial regions
00243       _h_fidXSecs->fill(1,_weight);
00244       if ( jets_30.size() >= 1 ) _h_fidXSecs->fill(2, _weight);
00245       if ( jets_30.size() >= 2 ) _h_fidXSecs->fill(3, _weight);
00246       if ( jets_30.size() >= 3 ) _h_fidXSecs->fill(4, _weight);
00247       if ( jets_30.size() >= 2 && passVBFCuts(y1 + y2, jets_30.at(0)->momentum(), jets_30.at(1)->momentum()) ) _h_fidXSecs->fill(5, _weight);
00248       if ( (good_el.size() + good_mu.size()) > 0 ) _h_fidXSecs->fill(6, _weight);
00249       if ( MET > 80 ) _h_fidXSecs->fill(7, _weight);
00250 
00251       // Fill histograms
00252       // Inclusive variables
00253       _pT_yy    = (y1 + y2).pT();
00254       _y_yy     = fabs( (y1 + y2).rapidity() );
00255       _cosTS_CS = cosTS_CS(y1, y2);
00256       _pTt_yy   = pTt(y1, y2);
00257       _Dy_yy    = fabs( deltaRap(y1, y2) );
00258 
00259       _Njets30 = jets_30.size() > 3 ? 3 : jets_30.size();
00260       _Njets50 = jets_50.size() > 3 ? 3 : jets_50.size();
00261       _h_Njets30->fill(_Njets30, _weight);
00262       _h_Njets50->fill(_Njets50, _weight);
00263 
00264       _pT_j1 = jets_30.size() > 0 ? jets_30.at(0)->momentum().pT() : 0.;
00265       _pT_j2 = jets_30.size() > 1 ? jets_30.at(1)->momentum().pT() : 0.;
00266       _pT_j3 = jets_30.size() > 2 ? jets_30.at(2)->momentum().pT() : 0.;
00267 
00268       _HT = 0.0;
00269       foreach (const Jet* jet, jets_30)
00270     _HT += jet->momentum().pT();
00271 
00272       _tau_jet     = tau_jet_max(y1 + y2, jets_25);
00273       _sum_tau_jet = sum_tau_jet(y1 + y2, jets_25);
00274 
00275       _h_pT_yy        ->fill(_pT_yy    ,_weight);
00276       _h_y_yy         ->fill(_y_yy     ,_weight);
00277       _h_pT_j1        ->fill(_pT_j1    ,_weight);
00278       _h_cosTS_CS     ->fill(_cosTS_CS ,_weight);
00279       _h_cosTS_CS_5bin->fill(_cosTS_CS ,_weight);
00280       _h_HT           ->fill(_HT       ,_weight);
00281       _h_pTt_yy       ->fill(_pTt_yy   ,_weight);
00282       _h_Dy_yy        ->fill(_Dy_yy    ,_weight);
00283       _h_tau_jet      ->fill(_tau_jet  ,_weight);
00284       _h_sum_tau_jet  ->fill(_sum_tau_jet,_weight);
00285 
00286       // >=1 jet variables
00287       if ( jets_30.size() >= 1 ) {
00288     FourMomentum j1 = jets_30[0]->momentum();
00289     _y_j1 = fabs( j1.rapidity() );
00290 
00291     _h_pT_j2->fill(_pT_j2 ,_weight);
00292     _h_y_j1 ->fill(_y_j1  ,_weight);
00293       }
00294 
00295       // >=2 jet variables
00296       if ( jets_30.size() >= 2 ) {
00297     FourMomentum j1 = jets_30[0]->momentum();
00298     FourMomentum j2 = jets_30[1]->momentum();
00299 
00300     _Dy_jj      = fabs( deltaRap(j1, j2) );
00301     _Dphi_jj    = fabs( deltaPhi(j1, j2) );
00302     _Dphi_yy_jj = fabs( deltaPhi(y1 + y2, j1 + j2) );
00303     _m_jj       = (j1 + j2).mass();
00304     _pT_yy_jj   = (y1 + y2 + j1 + j2).pT();
00305     _y_j2       = fabs( j2.rapidity() );
00306 
00307     _h_Dy_jj      ->fill(_Dy_jj     ,_weight);
00308     _h_Dphi_jj    ->fill(_Dphi_jj   ,_weight);
00309     _h_Dphi_yy_jj ->fill(_Dphi_yy_jj,_weight);
00310     _h_m_jj       ->fill(_m_jj      ,_weight);
00311     _h_pT_yy_jj   ->fill(_pT_yy_jj  ,_weight);
00312     _h_pT_j3      ->fill(_pT_j3     ,_weight);
00313     _h_y_j2       ->fill(_y_j2      ,_weight);
00314       }
00315 
00316       // 2D distributions of cosTS_CS x pT_yy
00317       if ( _pT_yy < 80 )
00318     _h_cosTS_pTyy_low->fill(_cosTS_CS, _weight);
00319       else if ( _pT_yy > 80 && _pT_yy < 200 )
00320     _h_cosTS_pTyy_high->fill(_cosTS_CS,_weight);
00321       else if ( _pT_yy > 200 )
00322     _h_cosTS_pTyy_rest->fill(_cosTS_CS,_weight);
00323 
00324       // 2D distributions of pT_yy x Njets
00325       if ( _Njets30 == 0 )
00326     _h_pTyy_Njets0->fill(_pT_yy, _weight);
00327       else if ( _Njets30 == 1 )
00328     _h_pTyy_Njets1->fill(_pT_yy, _weight);
00329       else if ( _Njets30 >= 2 )
00330     _h_pTyy_Njets2->fill(_pT_yy, _weight);
00331 
00332       if ( _Njets30 == 1 ) _h_pTj1_excl->fill(_pT_j1, _weight);
00333 
00334     }
00335 
00336     // Normalise histograms after the run
00337     void finalize() {
00338 
00339       const double xs = crossSectionPerEvent()/femtobarn;
00340 
00341       scale(_h_pT_yy, xs);
00342       scale(_h_y_yy, xs);
00343       scale(_h_pT_j1, xs);
00344       scale(_h_y_j1, xs);
00345       scale(_h_HT, xs);
00346       scale(_h_pT_j2, xs);
00347       scale(_h_Dy_jj, xs);
00348       scale(_h_Dphi_yy_jj, xs);
00349       scale(_h_cosTS_CS, xs);
00350       scale(_h_cosTS_CS_5bin, xs);
00351       scale(_h_Dphi_jj, xs);
00352       scale(_h_pTt_yy, xs);
00353       scale(_h_Dy_yy, xs);
00354       scale(_h_tau_jet, xs);
00355       scale(_h_sum_tau_jet, xs);
00356       scale(_h_y_j2, xs);
00357       scale(_h_pT_j3, xs);
00358       scale(_h_m_jj, xs);
00359       scale(_h_pT_yy_jj, xs);
00360       scale(_h_cosTS_pTyy_low, xs);
00361       scale(_h_cosTS_pTyy_high, xs);
00362       scale(_h_cosTS_pTyy_rest, xs);
00363       scale(_h_pTyy_Njets0, xs);
00364       scale(_h_pTyy_Njets1, xs);
00365       scale(_h_pTyy_Njets2, xs);
00366       scale(_h_pTj1_excl, xs);
00367       scale(_h_Njets30, xs);
00368       scale(_h_Njets50, xs);
00369       scale(_h_fidXSecs, xs);
00370     }
00371 
00372     // Trace event record to see if particle came from a hadron (or a tau from a hadron decay)
00373     // Based on fromDecay() function
00374     bool fromHadronDecay(const Particle& p ) {
00375       const GenVertex* prodVtx = p.genParticle()->production_vertex();
00376       if (prodVtx == NULL) return false;
00377       foreach (const GenParticle* ancestor, particles(prodVtx, HepMC::ancestors)) {
00378         const PdgId pid = ancestor->pdg_id();
00379         if (ancestor->status() == 2 && PID::isHadron(pid)) return true;
00380         if (ancestor->status() == 2 && (abs(pid) == PID::TAU && fromHadronDecay(ancestor))) return true;
00381       }
00382       return false;
00383     }
00384 
00385     // VBF-enhanced dijet topology selection cuts
00386     bool passVBFCuts(const FourMomentum &H, const FourMomentum &j1, const FourMomentum &j2) {
00387       return ( fabs(deltaRap(j1, j2)) > 2.8 && (j1 + j2).mass() > 400 && fabs(deltaPhi(H, j1 + j2)) > 2.6 );
00388     }
00389 
00390     // Cosine of the decay angle in the Collins-Soper frame
00391     double cosTS_CS(const FourMomentum &y1, const FourMomentum &y2) {
00392       return fabs( ( (y1.E() + y1.pz())* (y2.E() - y2.pz()) - (y1.E() - y1.pz()) * (y2.E() + y2.pz()) )
00393            / ((y1 + y2).mass() * sqrt(pow((y1 + y2).mass(), 2) + pow((y1 + y2).pt(), 2)) ) );
00394     }
00395 
00396     // Diphoton pT along thrust axis
00397     double pTt(const FourMomentum &y1, const FourMomentum &y2) {
00398       return fabs(y1.px() * y2.py() - y2.px() * y1.py()) / (y1 - y2).pT()*2;
00399     }
00400 
00401     // Tau of jet  (see paper for description)
00402     // tau_jet = mT/(2*cosh(y*)), where mT = pT (+) m, and y* = rapidty in Higgs rest frame
00403     double tau_jet( const FourMomentum &H, const FourMomentum &jet ) {
00404       return sqrt( pow(jet.pT(),2) + pow(jet.mass(),2) ) / (2.0 * cosh( jet.rapidity() - H.rapidity() ) );
00405     }
00406 
00407     // Maximal (leading) tau_jet (see paper for description)
00408     double tau_jet_max( const FourMomentum &H, const vector<const Jet*> jets, double tau_jet_cut = 8. ) {
00409       double max_tj = 0;
00410       for (size_t i=0; i < jets.size(); ++i) {
00411     FourMomentum jet = jets[i]->momentum();
00412     if ( tau_jet(H, jet) > tau_jet_cut )
00413       max_tj = max( tau_jet(H, jet), max_tj );
00414       }
00415       return max_tj;
00416     }
00417 
00418     // Scalar sum of tau for all jets (see paper for description)
00419     double sum_tau_jet( const FourMomentum &H, const vector<const Jet*> jets, double tau_jet_cut = 8. ) {
00420       double sum_tj = 0;
00421       for (size_t i=0; i < jets.size(); ++i) {
00422     FourMomentum jet = jets[i]->momentum();
00423     if ( tau_jet(H, jet) > tau_jet_cut )
00424       sum_tj += tau_jet(H, jet);
00425       }
00426       return sum_tj;
00427     }
00428 
00429   private:
00430 
00431     Histo1DPtr _h_pT_yy;
00432     Histo1DPtr _h_y_yy;
00433     Histo1DPtr _h_Njets30;
00434     Histo1DPtr _h_Njets50;
00435     Histo1DPtr _h_pT_j1;
00436     Histo1DPtr _h_y_j1;
00437     Histo1DPtr _h_HT;
00438     Histo1DPtr _h_pT_j2;
00439     Histo1DPtr _h_Dy_jj;
00440     Histo1DPtr _h_Dphi_yy_jj;
00441     Histo1DPtr _h_cosTS_CS;
00442     Histo1DPtr _h_cosTS_CS_5bin;
00443     Histo1DPtr _h_Dphi_jj;
00444     Histo1DPtr _h_pTt_yy;
00445     Histo1DPtr _h_Dy_yy;
00446     Histo1DPtr _h_tau_jet;
00447     Histo1DPtr _h_sum_tau_jet;
00448     Histo1DPtr _h_y_j2;
00449     Histo1DPtr _h_pT_j3;
00450     Histo1DPtr _h_m_jj;
00451     Histo1DPtr _h_pT_yy_jj;
00452     Histo1DPtr _h_cosTS_pTyy_low;
00453     Histo1DPtr _h_cosTS_pTyy_high;
00454     Histo1DPtr _h_cosTS_pTyy_rest;
00455     Histo1DPtr _h_pTyy_Njets0;
00456     Histo1DPtr _h_pTyy_Njets1;
00457     Histo1DPtr _h_pTyy_Njets2;
00458     Histo1DPtr _h_pTj1_excl;
00459     Histo1DPtr _h_fidXSecs;
00460 
00461     double _weight;
00462     int _Njets30;
00463     int _Njets50;
00464     double _pT_yy;
00465     double _y_yy;
00466     double _cosTS_CS;
00467     double _pT_j1;
00468     double _m_jj;
00469     double _y_j1;
00470     double _HT;
00471     double _pT_j2;
00472     double _y_j2;
00473     double _Dphi_yy_jj;
00474     double _pT_yy_jj;
00475     double _Dphi_jj;
00476     double _Dy_jj;
00477     double _pT_j3;
00478     double _pTt_yy;
00479     double _Dy_yy;
00480     double _tau_jet;
00481     double _sum_tau_jet;
00482   };
00483 
00484   // The hook for the plugin system
00485   DECLARE_RIVET_PLUGIN(ATLAS_2014_I1306615);
00486 
00487 }