ATLAS_2014_I1327229.cc
Go to the documentation of this file.
00001 // -*- C++ -*- 00002 #include "Rivet/Analysis.hh" 00003 #include "Rivet/Projections/FinalState.hh" 00004 #include "Rivet/Projections/ChargedFinalState.hh" 00005 #include "Rivet/Projections/VisibleFinalState.hh" 00006 #include "Rivet/Projections/VetoedFinalState.hh" 00007 #include "Rivet/Projections/IdentifiedFinalState.hh" 00008 #include "Rivet/Projections/UnstableFinalState.hh" 00009 #include "Rivet/Projections/FastJets.hh" 00010 00011 00012 namespace Rivet { 00013 00014 00015 class ATLAS_2014_I1327229 : public Analysis { 00016 public: 00017 00018 /// Constructor 00019 ATLAS_2014_I1327229() 00020 : Analysis("ATLAS_2014_I1327229") { } 00021 00022 00023 /// Book histograms and initialise projections before the run 00024 void init() { 00025 00026 // To calculate the acceptance without having the fiducial lepton efficiencies included, this part can be turned off 00027 _use_fiducial_lepton_efficiency = true; 00028 00029 // Random numbers for simulation of ATLAS detector reconstruction efficiency 00030 srand (160385); 00031 00032 // Read in all signal regions 00033 _signal_regions = getSignalRegions(); 00034 00035 // Set number of events per signal region to 0 00036 for (size_t i = 0; i < _signal_regions.size(); i++) 00037 _eventCountsPerSR[_signal_regions[i]] = 0.0; 00038 00039 // Final state including all charged and neutral particles 00040 const FinalState fs(-5.0, 5.0, 1*GeV); 00041 addProjection(fs, "FS"); 00042 00043 // Final state including all charged particles 00044 addProjection(ChargedFinalState(-2.5, 2.5, 1*GeV), "CFS"); 00045 00046 // Final state including all visible particles (to calculate MET, Jets etc.) 00047 addProjection(VisibleFinalState(-5.0,5.0),"VFS"); 00048 00049 // Final state including all AntiKt 04 Jets 00050 VetoedFinalState vfs; 00051 vfs.addVetoPairId(PID::MUON); 00052 addProjection(FastJets(vfs, FastJets::ANTIKT, 0.4), "AntiKtJets04"); 00053 00054 // Final state including all unstable particles (including taus) 00055 addProjection(UnstableFinalState(Cuts::abseta < 5.0 && Cuts::pT > 5*GeV),"UFS"); 00056 00057 // Final state including all electrons 00058 IdentifiedFinalState elecs(Cuts::abseta < 2.47 && Cuts::pT > 10*GeV); 00059 elecs.acceptIdPair(PID::ELECTRON); 00060 addProjection(elecs, "elecs"); 00061 00062 // Final state including all muons 00063 IdentifiedFinalState muons(Cuts::abseta < 2.5 && Cuts::pT > 10*GeV); 00064 muons.acceptIdPair(PID::MUON); 00065 addProjection(muons, "muons"); 00066 00067 00068 00069 /// Book histograms: 00070 _h_HTlep_all = bookHisto1D("HTlep_all", 30,0,3000); 00071 _h_HTjets_all = bookHisto1D("HTjets_all", 30,0,3000); 00072 _h_MET_all = bookHisto1D("MET_all", 30,0,1500); 00073 _h_Meff_all = bookHisto1D("Meff_all", 50,0,5000); 00074 _h_min_pT_all = bookHisto1D("min_pT_all", 50, 0, 2000); 00075 _h_mT_all = bookHisto1D("mT_all", 50, 0, 2000); 00076 00077 _h_e_n = bookHisto1D("e_n", 10, -0.5, 9.5); 00078 _h_mu_n = bookHisto1D("mu_n", 10, -0.5, 9.5); 00079 _h_tau_n = bookHisto1D("tau_n", 10, -0.5, 9.5); 00080 00081 _h_pt_1_3l = bookHisto1D("pt_1_3l", 100, 0, 2000); 00082 _h_pt_2_3l = bookHisto1D("pt_2_3l", 100, 0, 2000); 00083 _h_pt_3_3l = bookHisto1D("pt_3_3l", 100, 0, 2000); 00084 _h_pt_1_2ltau = bookHisto1D("pt_1_2ltau", 100, 0, 2000); 00085 _h_pt_2_2ltau = bookHisto1D("pt_2_2ltau", 100, 0, 2000); 00086 _h_pt_3_2ltau = bookHisto1D("pt_3_2ltau", 100, 0, 2000); 00087 00088 _h_excluded = bookHisto1D("excluded", 2, -0.5, 1.5); 00089 } 00090 00091 00092 /// Perform the per-event analysis 00093 void analyze(const Event& event) { 00094 00095 // Muons 00096 Particles muon_candidates; 00097 const Particles charged_tracks = applyProjection<ChargedFinalState>(event, "CFS").particles(); 00098 const Particles visible_particles = applyProjection<VisibleFinalState>(event, "VFS").particles(); 00099 foreach (const Particle& mu, applyProjection<IdentifiedFinalState>(event, "muons").particlesByPt() ) { 00100 00101 // Calculate pTCone30 variable (pT of all tracks within dR<0.3 - pT of muon itself) 00102 double pTinCone = -mu.pT(); 00103 foreach (const Particle& track, charged_tracks ) { 00104 if (deltaR(mu.momentum(),track.momentum()) < 0.3 ) 00105 pTinCone += track.pT(); 00106 } 00107 00108 // Calculate eTCone30 variable (pT of all visible particles within dR<0.3) 00109 double eTinCone = 0.; 00110 foreach (const Particle& visible_particle, visible_particles) { 00111 if (visible_particle.abspid() != PID::MUON && inRange(deltaR(mu.momentum(),visible_particle.momentum()), 0.1, 0.3)) 00112 eTinCone += visible_particle.pT(); 00113 } 00114 00115 // Apply reconstruction efficiency and simulate reconstruction 00116 int muon_id = 13; 00117 if (mu.hasAncestor(15) || mu.hasAncestor(-15)) muon_id = 14; 00118 const double eff = (_use_fiducial_lepton_efficiency) ? apply_reco_eff(muon_id,mu) : 1.0; 00119 const bool keep_muon = rand()/static_cast<double>(RAND_MAX)<=eff; 00120 00121 // Keep muon if pTCone30/pT < 0.15 and eTCone30/pT < 0.2 and reconstructed 00122 if (keep_muon && pTinCone/mu.pT() <= 0.1 && eTinCone/mu.pT() < 0.1) 00123 muon_candidates.push_back(mu); 00124 } 00125 00126 // Electrons 00127 Particles electron_candidates; 00128 foreach (const Particle& e, applyProjection<IdentifiedFinalState>(event, "elecs").particlesByPt() ) { 00129 // Neglect electrons in crack regions 00130 if (inRange(e.abseta(), 1.37, 1.52)) continue; 00131 00132 // Calculate pTCone30 variable (pT of all tracks within dR<0.3 - pT of electron itself) 00133 double pTinCone = -e.pT(); 00134 foreach (const Particle& track, charged_tracks) { 00135 if (deltaR(e.momentum(), track.momentum()) < 0.3 ) pTinCone += track.pT(); 00136 } 00137 00138 // Calculate eTCone30 variable (pT of all visible particles (except muons) within dR<0.3) 00139 double eTinCone = 0.; 00140 foreach (const Particle& visible_particle, visible_particles) { 00141 if (visible_particle.abspid() != PID::MUON && inRange(deltaR(e.momentum(),visible_particle.momentum()), 0.1, 0.3)) 00142 eTinCone += visible_particle.pT(); 00143 } 00144 00145 // Apply reconstruction efficiency and simulate reconstruction 00146 int elec_id = 11; 00147 if (e.hasAncestor(15) || e.hasAncestor(-15)) elec_id = 12; 00148 const double eff = (_use_fiducial_lepton_efficiency) ? apply_reco_eff(elec_id,e) : 1.0; 00149 const bool keep_elec = rand()/static_cast<double>(RAND_MAX)<=eff; 00150 00151 // Keep electron if pTCone30/pT < 0.13 and eTCone30/pT < 0.2 and reconstructed 00152 if (keep_elec && pTinCone/e.pT() <= 0.1 && eTinCone/e.pT() < 0.1) 00153 electron_candidates.push_back(e); 00154 } 00155 00156 00157 // Taus 00158 Particles tau_candidates; 00159 foreach (const Particle& tau, applyProjection<UnstableFinalState>(event, "UFS").particles() ) { 00160 // Only pick taus out of all unstable particles 00161 if ( tau.abspid() != PID::TAU) continue; 00162 // Check that tau has decayed into daughter particles 00163 if (tau.genParticle()->end_vertex() == 0) continue; 00164 // Calculate visible tau momentum using the tau neutrino momentum in the tau decay 00165 FourMomentum daughter_tau_neutrino_momentum = get_tau_neutrino_momentum(tau); 00166 Particle tau_vis = tau; 00167 tau_vis.setMomentum(tau.momentum()-daughter_tau_neutrino_momentum); 00168 // keep only taus in certain eta region and above 15 GeV of visible tau pT 00169 if ( tau_vis.pT()/GeV <= 15.0 || tau_vis.abseta() > 2.5) continue; 00170 00171 // Get prong number (number of tracks) in tau decay and check if tau decays leptonically 00172 unsigned int nprong = 0; 00173 bool lep_decaying_tau = false; 00174 get_prong_number(tau.genParticle(),nprong,lep_decaying_tau); 00175 00176 // Apply reconstruction efficiency and simulate reconstruction 00177 int tau_id = 15; 00178 if (nprong == 1) tau_id = 15; 00179 else if (nprong == 3) tau_id = 16; 00180 00181 00182 const double eff = (_use_fiducial_lepton_efficiency) ? apply_reco_eff(tau_id,tau_vis) : 1.0; 00183 const bool keep_tau = rand()/static_cast<double>(RAND_MAX)<=eff; 00184 00185 // Keep tau if nprong = 1, it decays hadronically and it is reconstructed 00186 if ( !lep_decaying_tau && nprong == 1 && keep_tau) tau_candidates.push_back(tau_vis); 00187 } 00188 00189 // Jets (all anti-kt R=0.4 jets with pT > 30 GeV and eta < 4.9 00190 Jets jet_candidates; 00191 foreach (const Jet& jet, applyProjection<FastJets>(event, "AntiKtJets04").jetsByPt(30.0*GeV) ) { 00192 if (jet.abseta() < 4.9 ) jet_candidates.push_back(jet); 00193 } 00194 00195 // ETmiss 00196 Particles vfs_particles = applyProjection<VisibleFinalState>(event, "VFS").particles(); 00197 FourMomentum pTmiss; 00198 foreach (const Particle& p, vfs_particles ) pTmiss -= p.momentum(); 00199 double eTmiss = pTmiss.pT()/GeV; 00200 00201 00202 // ------------------------- 00203 // Overlap removal 00204 00205 // electron - electron 00206 Particles electron_candidates_2; 00207 for(size_t ie = 0; ie < electron_candidates.size(); ++ie) { 00208 const Particle& e = electron_candidates[ie]; 00209 bool away = true; 00210 // If electron pair within dR < 0.1: remove electron with lower pT 00211 for(size_t ie2 = 0; ie2 < electron_candidates_2.size(); ++ie2) { 00212 if (deltaR(e.momentum(),electron_candidates_2[ie2].momentum()) < 0.1 ) { 00213 away = false; 00214 break; 00215 } 00216 } 00217 // If isolated keep it 00218 if ( away ) 00219 electron_candidates_2.push_back( e ); 00220 } 00221 00222 // jet - electron 00223 Jets recon_jets; 00224 foreach (const Jet& jet, jet_candidates) { 00225 bool away = true; 00226 // If jet within dR < 0.2 of electron: remove jet 00227 foreach (const Particle& e, electron_candidates_2) { 00228 if (deltaR(e.momentum(), jet.momentum()) < 0.2 ) { 00229 away = false; 00230 break; 00231 } 00232 } 00233 00234 // jet - tau 00235 if ( away ) { 00236 // If jet within dR < 0.2 of tau: remove jet 00237 foreach (const Particle& tau, tau_candidates) { 00238 if (deltaR(tau.momentum(), jet.momentum()) < 0.2 ) { 00239 away = false; 00240 break; 00241 } 00242 } 00243 } 00244 // If isolated keep it 00245 if ( away ) 00246 recon_jets.push_back( jet ); 00247 } 00248 00249 // electron - jet 00250 Particles recon_leptons, recon_e; 00251 for (size_t ie = 0; ie < electron_candidates_2.size(); ++ie) { 00252 const Particle& e = electron_candidates_2[ie]; 00253 // If electron within 0.2 < dR < 0.4 from any jets: remove electron 00254 bool away = true; 00255 foreach (const Jet& jet, recon_jets) { 00256 if (deltaR(e.momentum(), jet.momentum()) < 0.4 ) { 00257 away = false; 00258 break; 00259 } 00260 } 00261 // electron - muon 00262 // If electron within dR < 0.1 of a muon: remove electron 00263 if (away) { 00264 foreach (const Particle& mu, muon_candidates) { 00265 if (deltaR(mu.momentum(),e.momentum()) < 0.1) { 00266 away = false; 00267 break; 00268 } 00269 } 00270 } 00271 // If isolated keep it 00272 if ( away ) { 00273 recon_e.push_back( e ); 00274 recon_leptons.push_back( e ); 00275 } 00276 } 00277 00278 // tau - electron 00279 Particles recon_tau; 00280 foreach (const Particle& tau, tau_candidates) { 00281 bool away = true; 00282 // If tau within dR < 0.2 of an electron: remove tau 00283 foreach (const Particle & e, recon_e) { 00284 if (deltaR(tau.momentum(),e.momentum()) < 0.2 ) { 00285 away = false; 00286 break; 00287 } 00288 } 00289 // tau - muon 00290 // If tau within dR < 0.2 of a muon: remove tau 00291 if (away) { 00292 foreach (const Particle& mu, muon_candidates) { 00293 if (deltaR(tau.momentum(), mu.momentum()) < 0.2 ) { 00294 away = false; 00295 break; 00296 } 00297 } 00298 } 00299 // If isolated keep it 00300 if (away) recon_tau.push_back( tau ); 00301 } 00302 00303 // muon - jet 00304 Particles recon_mu, trigger_mu; 00305 // If muon within dR < 0.4 of a jet: remove muon 00306 foreach (const Particle& mu, muon_candidates ) { 00307 bool away = true; 00308 foreach (const Jet& jet, recon_jets) { 00309 if (deltaR(mu.momentum(), jet.momentum()) < 0.4 ) { 00310 away = false; 00311 break; 00312 } 00313 } 00314 if (away) { 00315 recon_mu.push_back( mu ); 00316 recon_leptons.push_back( mu ); 00317 if (mu.abseta() < 2.4) trigger_mu.push_back( mu ); 00318 } 00319 } 00320 00321 // End overlap removal 00322 // --------------------- 00323 00324 // Jet cleaning 00325 if (rand()/static_cast<double>(RAND_MAX) <= 0.42) { 00326 foreach (const Jet& jet, recon_jets ) { 00327 const double eta = jet.rapidity(); 00328 const double phi = jet.azimuthalAngle(MINUSPI_PLUSPI); 00329 if(jet.pT() > 25*GeV && inRange(eta,-0.1,1.5) && inRange(phi,-0.9,-0.5)) vetoEvent; 00330 } 00331 } 00332 00333 // Event selection 00334 // Require at least 3 charged tracks in event 00335 if (charged_tracks.size() < 3) vetoEvent; 00336 00337 // And at least one e/mu passing trigger 00338 if( !( !recon_e.empty() && recon_e[0].pT()>26.*GeV) && 00339 !( !trigger_mu.empty() && trigger_mu[0].pT()>26.*GeV) ) { 00340 MSG_DEBUG("Hardest lepton fails trigger"); 00341 vetoEvent; 00342 } 00343 00344 // And only accept events with at least 2 electrons and muons and at least 3 leptons in total 00345 if (recon_mu.size() + recon_e.size() + recon_tau.size() < 3 || recon_leptons.size() < 2) vetoEvent; 00346 00347 // Getting the event weight 00348 const double weight = event.weight(); 00349 00350 // Sort leptons by decreasing pT 00351 sortByPt(recon_leptons); 00352 sortByPt(recon_tau); 00353 00354 // Calculate HTlep, fill lepton pT histograms & store chosen combination of 3 leptons 00355 double HTlep = 0.; 00356 Particles chosen_leptons; 00357 if (recon_leptons.size() > 2) { 00358 _h_pt_1_3l->fill(recon_leptons[0].pT()/GeV, weight); 00359 _h_pt_2_3l->fill(recon_leptons[1].pT()/GeV, weight); 00360 _h_pt_3_3l->fill(recon_leptons[2].pT()/GeV, weight); 00361 HTlep = (recon_leptons[0].pT() + recon_leptons[1].pT() + recon_leptons[2].pT())/GeV; 00362 chosen_leptons.push_back( recon_leptons[0] ); 00363 chosen_leptons.push_back( recon_leptons[1] ); 00364 chosen_leptons.push_back( recon_leptons[2] ); 00365 } 00366 else { 00367 _h_pt_1_2ltau->fill(recon_leptons[0].pT()/GeV, weight); 00368 _h_pt_2_2ltau->fill(recon_leptons[1].pT()/GeV, weight); 00369 _h_pt_3_2ltau->fill(recon_tau[0].pT()/GeV, weight); 00370 HTlep = recon_leptons[0].pT()/GeV + recon_leptons[1].pT()/GeV + recon_tau[0].pT()/GeV; 00371 chosen_leptons.push_back( recon_leptons[0] ); 00372 chosen_leptons.push_back( recon_leptons[1] ); 00373 chosen_leptons.push_back( recon_tau[0] ); 00374 } 00375 00376 // Calculate mT and mTW variable 00377 Particles mT_leptons; 00378 Particles mTW_leptons; 00379 for (size_t i1 = 0; i1 < 3; i1 ++) { 00380 for (size_t i2 = i1+1; i2 < 3; i2 ++) { 00381 double OSSF_inv_mass = isOSSF_mass(chosen_leptons[i1],chosen_leptons[i2]); 00382 if (OSSF_inv_mass != 0.) { 00383 for (size_t i3 = 0; i3 < 3 ; i3 ++) { 00384 if (i3 != i2 && i3 != i1) { 00385 mT_leptons.push_back(chosen_leptons[i3]); 00386 if ( fabs(91.0 - OSSF_inv_mass) < 20. ) 00387 mTW_leptons.push_back(chosen_leptons[i3]); 00388 } 00389 } 00390 } 00391 else { 00392 mT_leptons.push_back(chosen_leptons[0]); 00393 mTW_leptons.push_back(chosen_leptons[0]); 00394 } 00395 } 00396 } 00397 00398 sortByPt(mT_leptons); 00399 sortByPt(mTW_leptons); 00400 00401 double mT = sqrt(2*pTmiss.pT()/GeV*mT_leptons[0].pT()/GeV*(1-cos(pTmiss.phi()-mT_leptons[0].phi()))); 00402 double mTW = sqrt(2*pTmiss.pT()/GeV*mTW_leptons[0].pT()/GeV*(1-cos(pTmiss.phi()-mTW_leptons[0].phi()))); 00403 00404 // Calculate Min pT variable 00405 double min_pT = chosen_leptons[2].pT()/GeV; 00406 00407 // Number of prompt e/mu and had taus 00408 _h_e_n->fill(recon_e.size(),weight); 00409 _h_mu_n->fill(recon_mu.size(),weight); 00410 _h_tau_n->fill(recon_tau.size(),weight); 00411 00412 // Calculate HTjets variable 00413 double HTjets = 0.; 00414 foreach (const Jet& jet, recon_jets) 00415 HTjets += jet.pT()/GeV; 00416 00417 // Calculate meff variable 00418 double meff = eTmiss + HTjets; 00419 Particles all_leptons; 00420 foreach (const Particle& e, recon_e ) { 00421 meff += e.pT()/GeV; 00422 all_leptons.push_back( e ); 00423 } 00424 foreach (const Particle& mu, recon_mu) { 00425 meff += mu.pT()/GeV; 00426 all_leptons.push_back( mu ); 00427 } 00428 foreach (const Particle& tau, recon_tau) { 00429 meff += tau.pT()/GeV; 00430 all_leptons.push_back( tau ); 00431 } 00432 00433 // Fill histograms of kinematic variables 00434 _h_HTlep_all->fill(HTlep,weight); 00435 _h_HTjets_all->fill(HTjets,weight); 00436 _h_MET_all->fill(eTmiss,weight); 00437 _h_Meff_all->fill(meff,weight); 00438 _h_min_pT_all->fill(min_pT,weight); 00439 _h_mT_all->fill(mT,weight); 00440 00441 // Determine signal region (3l / 2ltau , onZ / offZ OSSF / offZ no-OSSF) 00442 // 3l vs. 2ltau 00443 string basic_signal_region; 00444 if (recon_mu.size() + recon_e.size() > 2) 00445 basic_signal_region += "3l_"; 00446 else if ( (recon_mu.size() + recon_e.size() == 2) && (recon_tau.size() > 0)) 00447 basic_signal_region += "2ltau_"; 00448 00449 // Is there an OSSF pair or a three lepton combination with an invariant mass close to the Z mass 00450 int onZ = isonZ(chosen_leptons); 00451 if (onZ == 1) basic_signal_region += "onZ"; 00452 else if (onZ == 0) { 00453 bool OSSF = isOSSF(chosen_leptons); 00454 if (OSSF) basic_signal_region += "offZ_OSSF"; 00455 else basic_signal_region += "offZ_noOSSF"; 00456 } 00457 00458 // Check in which signal regions this event falls and adjust event counters 00459 // INFO: The b-jet signal regions of the paper are not included in this Rivet implementation 00460 fillEventCountsPerSR(basic_signal_region,onZ,HTlep,eTmiss,HTjets,meff,min_pT,mTW,weight); 00461 } 00462 00463 00464 /// Normalise histograms etc., after the run 00465 void finalize() { 00466 00467 // Normalize to an integrated luminosity of 1 fb-1 00468 double norm = crossSection()/femtobarn/sumOfWeights(); 00469 00470 string best_signal_region = ""; 00471 double ratio_best_SR = 0.; 00472 00473 // Loop over all signal regions and find signal region with best sensitivity (ratio signal events/visible cross-section) 00474 for (size_t i = 0; i < _signal_regions.size(); i++) { 00475 double signal_events = _eventCountsPerSR[_signal_regions[i]] * norm; 00476 // Use expected upper limits to find best signal region: 00477 double UL95 = getUpperLimit(_signal_regions[i],false); 00478 double ratio = signal_events / UL95; 00479 if (ratio > ratio_best_SR) { 00480 best_signal_region = _signal_regions.at(i); 00481 ratio_best_SR = ratio; 00482 } 00483 } 00484 00485 double signal_events_best_SR = _eventCountsPerSR[best_signal_region] * norm; 00486 double exp_UL_best_SR = getUpperLimit(best_signal_region, false); 00487 double obs_UL_best_SR = getUpperLimit(best_signal_region, true); 00488 00489 00490 // Print out result 00491 cout << "----------------------------------------------------------------------------------------" << endl; 00492 cout << "Number of total events: " << sumOfWeights() << endl; 00493 cout << "Best signal region: " << best_signal_region << endl; 00494 cout << "Normalized number of signal events in this best signal region (per fb-1): " << signal_events_best_SR << endl; 00495 cout << "Efficiency*Acceptance: " << _eventCountsPerSR[best_signal_region]/sumOfWeights() << endl; 00496 cout << "Cross-section [fb]: " << crossSection()/femtobarn << endl; 00497 cout << "Expected visible cross-section (per fb-1): " << exp_UL_best_SR << endl; 00498 cout << "Ratio (signal events / expected visible cross-section): " << ratio_best_SR << endl; 00499 cout << "Observed visible cross-section (per fb-1): " << obs_UL_best_SR << endl; 00500 cout << "Ratio (signal events / observed visible cross-section): " << signal_events_best_SR/obs_UL_best_SR << endl; 00501 cout << "----------------------------------------------------------------------------------------" << endl; 00502 00503 cout << "Using the EXPECTED limits (visible cross-section) of the analysis: " << endl; 00504 if (signal_events_best_SR > exp_UL_best_SR) { 00505 cout << "Since the number of signal events > the visible cross-section, this model/grid point is EXCLUDED with 95% C.L." << endl; 00506 _h_excluded->fill(1); 00507 } 00508 else { 00509 cout << "Since the number of signal events < the visible cross-section, this model/grid point is NOT EXCLUDED." << endl; 00510 _h_excluded->fill(0); 00511 } 00512 cout << "----------------------------------------------------------------------------------------" << endl; 00513 00514 cout << "Using the OBSERVED limits (visible cross-section) of the analysis: " << endl; 00515 if (signal_events_best_SR > obs_UL_best_SR) { 00516 cout << "Since the number of signal events > the visible cross-section, this model/grid point is EXCLUDED with 95% C.L." << endl; 00517 _h_excluded->fill(1); 00518 } 00519 else { 00520 cout << "Since the number of signal events < the visible cross-section, this model/grid point is NOT EXCLUDED." << endl; 00521 _h_excluded->fill(0); 00522 } 00523 cout << "----------------------------------------------------------------------------------------" << endl; 00524 cout << "INFO: The b-jet signal regions of the paper are not included in this Rivet implementation." << endl; 00525 cout << "----------------------------------------------------------------------------------------" << endl; 00526 00527 00528 /// Normalize to cross section 00529 00530 if (norm != 0) { 00531 scale(_h_HTlep_all, norm); 00532 scale(_h_HTjets_all, norm); 00533 scale(_h_MET_all, norm); 00534 scale(_h_Meff_all, norm); 00535 scale(_h_min_pT_all, norm); 00536 scale(_h_mT_all, norm); 00537 00538 scale(_h_pt_1_3l, norm); 00539 scale(_h_pt_2_3l, norm); 00540 scale(_h_pt_3_3l, norm); 00541 scale(_h_pt_1_2ltau, norm); 00542 scale(_h_pt_2_2ltau, norm); 00543 scale(_h_pt_3_2ltau, norm); 00544 00545 scale(_h_e_n, norm); 00546 scale(_h_mu_n, norm); 00547 scale(_h_tau_n, norm); 00548 00549 scale(_h_excluded, norm); 00550 } 00551 00552 } 00553 00554 00555 /// Helper functions 00556 //@{ 00557 /// Function giving a list of all signal regions 00558 vector<string> getSignalRegions() { 00559 00560 // List of basic signal regions 00561 vector<string> basic_signal_regions; 00562 basic_signal_regions.push_back("3l_offZ_OSSF"); 00563 basic_signal_regions.push_back("3l_offZ_noOSSF"); 00564 basic_signal_regions.push_back("3l_onZ"); 00565 basic_signal_regions.push_back("2ltau_offZ_OSSF"); 00566 basic_signal_regions.push_back("2ltau_offZ_noOSSF"); 00567 basic_signal_regions.push_back("2ltau_onZ"); 00568 00569 // List of kinematic variables 00570 vector<string> kinematic_variables; 00571 kinematic_variables.push_back("HTlep"); 00572 kinematic_variables.push_back("METStrong"); 00573 kinematic_variables.push_back("METWeak"); 00574 kinematic_variables.push_back("Meff"); 00575 kinematic_variables.push_back("MeffStrong"); 00576 kinematic_variables.push_back("MeffMt"); 00577 kinematic_variables.push_back("MinPt"); 00578 00579 vector<string> signal_regions; 00580 // Loop over all kinematic variables and basic signal regions 00581 for (size_t i0 = 0; i0 < kinematic_variables.size(); i0++) { 00582 for (size_t i1 = 0; i1 < basic_signal_regions.size(); i1++) { 00583 // Is signal region onZ? 00584 int onZ = (basic_signal_regions[i1].find("onZ") != string::npos) ? 1 : 0; 00585 // Get cut values for this kinematic variable 00586 vector<int> cut_values = getCutsPerSignalRegion(kinematic_variables[i0], onZ); 00587 // Loop over all cut values 00588 for (size_t i2 = 0; i2 < cut_values.size(); i2++) { 00589 // Push signal region into vector 00590 signal_regions.push_back( kinematic_variables[i0] + "_" + basic_signal_regions[i1] + "_cut_" + toString(cut_values[i2]) ); 00591 } 00592 } 00593 } 00594 return signal_regions; 00595 } 00596 00597 00598 00599 /// Function giving all cut values per kinematic variable 00600 vector<int> getCutsPerSignalRegion(const string& signal_region, int onZ = 0) { 00601 vector<int> cutValues; 00602 00603 // Cut values for HTlep 00604 if (signal_region.compare("HTlep") == 0) { 00605 cutValues.push_back(0); 00606 cutValues.push_back(200); 00607 cutValues.push_back(500); 00608 cutValues.push_back(800); 00609 } 00610 // Cut values for MinPt 00611 else if (signal_region.compare("MinPt") == 0) { 00612 cutValues.push_back(0); 00613 cutValues.push_back(50); 00614 cutValues.push_back(100); 00615 cutValues.push_back(150); 00616 } 00617 // Cut values for METStrong (HTjets > 150 GeV) and METWeak (HTjets < 150 GeV) 00618 else if (signal_region.compare("METStrong") == 0 || signal_region.compare("METWeak") == 0) { 00619 cutValues.push_back(0); 00620 cutValues.push_back(100); 00621 cutValues.push_back(200); 00622 cutValues.push_back(300); 00623 } 00624 // Cut values for Meff 00625 if (signal_region.compare("Meff") == 0) { 00626 cutValues.push_back(0); 00627 cutValues.push_back(600); 00628 cutValues.push_back(1000); 00629 cutValues.push_back(1500); 00630 } 00631 // Cut values for MeffStrong (MET > 100 GeV) 00632 if ((signal_region.compare("MeffStrong") == 0 || signal_region.compare("MeffMt") == 0) && onZ ==1) { 00633 cutValues.push_back(0); 00634 cutValues.push_back(600); 00635 cutValues.push_back(1200); 00636 } 00637 00638 return cutValues; 00639 } 00640 00641 /// function fills map _eventCountsPerSR by looping over all signal regions 00642 /// and looking if the event falls into this signal region 00643 void fillEventCountsPerSR(const string& basic_signal_region, int onZ, 00644 double HTlep, double eTmiss, double HTjets, 00645 double meff, double min_pT, double mTW, 00646 double weight) { 00647 00648 // Get cut values for HTlep, loop over them and add event if cut is passed 00649 vector<int> cut_values = getCutsPerSignalRegion("HTlep", onZ); 00650 for (size_t i = 0; i < cut_values.size(); i++) { 00651 if (HTlep > cut_values[i]) 00652 _eventCountsPerSR[("HTlep_" + basic_signal_region + "_cut_" + toString(cut_values[i]))] += weight; 00653 } 00654 00655 // Get cut values for MinPt, loop over them and add event if cut is passed 00656 cut_values = getCutsPerSignalRegion("MinPt", onZ); 00657 for (size_t i = 0; i < cut_values.size(); i++) { 00658 if (min_pT > cut_values[i]) 00659 _eventCountsPerSR[("MinPt_" + basic_signal_region + "_cut_" + toString(cut_values[i]))] += weight; 00660 } 00661 00662 // Get cut values for METStrong, loop over them and add event if cut is passed 00663 cut_values = getCutsPerSignalRegion("METStrong", onZ); 00664 for (size_t i = 0; i < cut_values.size(); i++) { 00665 if (eTmiss > cut_values[i] && HTjets > 150.) 00666 _eventCountsPerSR[("METStrong_" + basic_signal_region + "_cut_" + toString(cut_values[i]))] += weight; 00667 } 00668 00669 // Get cut values for METWeak, loop over them and add event if cut is passed 00670 cut_values = getCutsPerSignalRegion("METWeak", onZ); 00671 for (size_t i = 0; i < cut_values.size(); i++) { 00672 if (eTmiss > cut_values[i] && HTjets <= 150.) 00673 _eventCountsPerSR[("METWeak_" + basic_signal_region + "_cut_" + toString(cut_values[i]))] += weight; 00674 } 00675 00676 // Get cut values for Meff, loop over them and add event if cut is passed 00677 cut_values = getCutsPerSignalRegion("Meff", onZ); 00678 for (size_t i = 0; i < cut_values.size(); i++) { 00679 if (meff > cut_values[i]) 00680 _eventCountsPerSR[("Meff_" + basic_signal_region + "_cut_" + toString(cut_values[i]))] += weight; 00681 } 00682 00683 // Get cut values for MeffStrong, loop over them and add event if cut is passed 00684 cut_values = getCutsPerSignalRegion("MeffStrong", onZ); 00685 for (size_t i = 0; i < cut_values.size(); i++) { 00686 if (meff > cut_values[i] && eTmiss > 100.) 00687 _eventCountsPerSR[("MeffStrong_" + basic_signal_region + "_cut_" + toString(cut_values[i]))] += weight; 00688 } 00689 00690 // Get cut values for MeffMt, loop over them and add event if cut is passed 00691 cut_values = getCutsPerSignalRegion("MeffMt", onZ); 00692 for (size_t i = 0; i < cut_values.size(); i++) { 00693 if (meff > cut_values[i] && mTW > 100. && onZ == 1) 00694 _eventCountsPerSR[("MeffMt_" + basic_signal_region + "_cut_" + toString(cut_values[i]))] += weight; 00695 } 00696 00697 } 00698 00699 /// Function returning 4-momentum of daughter-particle if it is a tau neutrino 00700 FourMomentum get_tau_neutrino_momentum(const Particle& p) { 00701 assert(p.abspid() == PID::TAU); 00702 const GenVertex* dv = p.genParticle()->end_vertex(); 00703 assert(dv != NULL); 00704 // Loop over all daughter particles 00705 for (GenVertex::particles_out_const_iterator pp = dv->particles_out_const_begin(); pp != dv->particles_out_const_end(); ++pp) { 00706 if (abs((*pp)->pdg_id()) == PID::NU_TAU) return FourMomentum((*pp)->momentum()); 00707 } 00708 return FourMomentum(); 00709 } 00710 00711 /// Function calculating the prong number of taus 00712 void get_prong_number(const GenParticle* p, unsigned int& nprong, bool& lep_decaying_tau) { 00713 assert(p != NULL); 00714 const GenVertex* dv = p->end_vertex(); 00715 assert(dv != NULL); 00716 for (GenVertex::particles_out_const_iterator pp = dv->particles_out_const_begin(); pp != dv->particles_out_const_end(); ++pp) { 00717 // If they have status 1 and are charged they will produce a track and the prong number is +1 00718 if ((*pp)->status() == 1 ) { 00719 const int id = (*pp)->pdg_id(); 00720 if (Rivet::PID::charge(id) != 0 ) ++nprong; 00721 // Check if tau decays leptonically 00722 if (( abs(id) == PID::ELECTRON || abs(id) == PID::MUON || abs(id) == PID::TAU ) && abs(p->pdg_id()) == PID::TAU) lep_decaying_tau = true; 00723 } 00724 // If the status of the daughter particle is 2 it is unstable and the further decays are checked 00725 else if ((*pp)->status() == 2 ) { 00726 get_prong_number((*pp),nprong,lep_decaying_tau); 00727 } 00728 } 00729 } 00730 00731 /// Function giving fiducial lepton efficiency 00732 double apply_reco_eff(int flavor, const Particle& p) { 00733 float pt = p.pT()/GeV; 00734 float eta = p.eta(); 00735 00736 double eff = 0.; 00737 00738 if (flavor == 11) { // weight prompt electron -- now including data/MC ID SF in eff. 00739 double avgrate = 0.685; 00740 float wz_ele[] = {0.0256,0.522,0.607,0.654,0.708,0.737,0.761,0.784,0.815,0.835,0.851,0.841,0.898}; 00741 // float ewz_ele[] = {0.000257,0.00492,0.00524,0.00519,0.00396,0.00449,0.00538,0.00513,0.00773,0.00753,0.0209,0.0964,0.259}; 00742 int ibin = 0; 00743 if(pt > 10 && pt < 15) ibin = 0; 00744 if(pt > 15 && pt < 20) ibin = 1; 00745 if(pt > 20 && pt < 25) ibin = 2; 00746 if(pt > 25 && pt < 30) ibin = 3; 00747 if(pt > 30 && pt < 40) ibin = 4; 00748 if(pt > 40 && pt < 50) ibin = 5; 00749 if(pt > 50 && pt < 60) ibin = 6; 00750 if(pt > 60 && pt < 80) ibin = 7; 00751 if(pt > 80 && pt < 100) ibin = 8; 00752 if(pt > 100 && pt < 200) ibin = 9; 00753 if(pt > 200 && pt < 400) ibin = 10; 00754 if(pt > 400 && pt < 600) ibin = 11; 00755 if(pt > 600) ibin = 12; 00756 double eff_pt = 0.; 00757 eff_pt = wz_ele[ibin]; 00758 00759 eta = fabs(eta); 00760 00761 float wz_ele_eta[] = {0.65,0.714,0.722,0.689,0.635,0.615}; 00762 // float ewz_ele_eta[] = {0.00642,0.00355,0.00335,0.004,0.00368,0.00422}; 00763 ibin = 0; 00764 if(eta > 0 && eta < 0.1) ibin = 0; 00765 if(eta > 0.1 && eta < 0.5) ibin = 1; 00766 if(eta > 0.5 && eta < 1.0) ibin = 2; 00767 if(eta > 1.0 && eta < 1.5) ibin = 3; 00768 if(eta > 1.5 && eta < 2.0) ibin = 4; 00769 if(eta > 2.0 && eta < 2.5) ibin = 5; 00770 double eff_eta = 0.; 00771 eff_eta = wz_ele_eta[ibin]; 00772 00773 eff = (eff_pt * eff_eta) / avgrate; 00774 } 00775 00776 if (flavor == 12) { // weight electron from tau 00777 double avgrate = 0.476; 00778 float wz_ele[] = {0.00855,0.409,0.442,0.55,0.632,0.616,0.615,0.642,0.72,0.617}; 00779 // float ewz_ele[] = {0.000573,0.0291,0.0366,0.0352,0.0363,0.0474,0.0628,0.0709,0.125,0.109}; 00780 int ibin = 0; 00781 if(pt > 10 && pt < 15) ibin = 0; 00782 if(pt > 15 && pt < 20) ibin = 1; 00783 if(pt > 20 && pt < 25) ibin = 2; 00784 if(pt > 25 && pt < 30) ibin = 3; 00785 if(pt > 30 && pt < 40) ibin = 4; 00786 if(pt > 40 && pt < 50) ibin = 5; 00787 if(pt > 50 && pt < 60) ibin = 6; 00788 if(pt > 60 && pt < 80) ibin = 7; 00789 if(pt > 80 && pt < 100) ibin = 8; 00790 if(pt > 100) ibin = 9; 00791 double eff_pt = 0.; 00792 eff_pt = wz_ele[ibin]; 00793 00794 eta = fabs(eta); 00795 00796 float wz_ele_eta[] = {0.546,0.5,0.513,0.421,0.47,0.433}; 00797 //float ewz_ele_eta[] = {0.0566,0.0257,0.0263,0.0263,0.0303,0.0321}; 00798 ibin = 0; 00799 if(eta > 0 && eta < 0.1) ibin = 0; 00800 if(eta > 0.1 && eta < 0.5) ibin = 1; 00801 if(eta > 0.5 && eta < 1.0) ibin = 2; 00802 if(eta > 1.0 && eta < 1.5) ibin = 3; 00803 if(eta > 1.5 && eta < 2.0) ibin = 4; 00804 if(eta > 2.0 && eta < 2.5) ibin = 5; 00805 double eff_eta = 0.; 00806 eff_eta = wz_ele_eta[ibin]; 00807 00808 eff = (eff_pt * eff_eta) / avgrate; 00809 } 00810 00811 if (flavor == 13) { // weight prompt muon 00812 int ibin = 0; 00813 if(pt > 10 && pt < 15) ibin = 0; 00814 if(pt > 15 && pt < 20) ibin = 1; 00815 if(pt > 20 && pt < 25) ibin = 2; 00816 if(pt > 25 && pt < 30) ibin = 3; 00817 if(pt > 30 && pt < 40) ibin = 4; 00818 if(pt > 40 && pt < 50) ibin = 5; 00819 if(pt > 50 && pt < 60) ibin = 6; 00820 if(pt > 60 && pt < 80) ibin = 7; 00821 if(pt > 80 && pt < 100) ibin = 8; 00822 if(pt > 100 && pt < 200) ibin = 9; 00823 if(pt > 200 && pt < 400) ibin = 10; 00824 if(pt > 400) ibin = 11; 00825 if(fabs(eta) < 0.1) { 00826 float wz_mu[] = {0.00705,0.402,0.478,0.49,0.492,0.499,0.527,0.512,0.53,0.528,0.465,0.465}; 00827 //float ewz_mu[] = {0.000298,0.0154,0.017,0.0158,0.0114,0.0123,0.0155,0.0133,0.0196,0.0182,0.0414,0.0414}; 00828 double eff_pt = 0.; 00829 eff_pt = wz_mu[ibin]; 00830 eff = eff_pt; 00831 } 00832 if(fabs(eta) > 0.1) { 00833 float wz_mu[] = {0.0224,0.839,0.887,0.91,0.919,0.923,0.925,0.925,0.922,0.918,0.884,0.834}; 00834 //float ewz_mu[] = {0.000213,0.00753,0.0074,0.007,0.00496,0.00534,0.00632,0.00583,0.00849,0.00804,0.0224,0.0963}; 00835 double eff_pt = 0.; 00836 eff_pt = wz_mu[ibin]; 00837 eff = eff_pt; 00838 } 00839 } 00840 00841 if (flavor == 14) { // weight muon from tau 00842 int ibin = 0; 00843 if(pt > 10 && pt < 15) ibin = 0; 00844 if(pt > 15 && pt < 20) ibin = 1; 00845 if(pt > 20 && pt < 25) ibin = 2; 00846 if(pt > 25 && pt < 30) ibin = 3; 00847 if(pt > 30 && pt < 40) ibin = 4; 00848 if(pt > 40 && pt < 50) ibin = 5; 00849 if(pt > 50 && pt < 60) ibin = 6; 00850 if(pt > 60 && pt < 80) ibin = 7; 00851 if(pt > 80 && pt < 100) ibin = 8; 00852 if(pt > 100) ibin = 9; 00853 00854 if(fabs(eta) < 0.1) { 00855 float wz_mu[] = {0.0,0.664,0.124,0.133,0.527,0.283,0.495,0.25,0.5,0.331}; 00856 //float ewz_mu[] = {0.0,0.192,0.0437,0.0343,0.128,0.107,0.202,0.125,0.25,0.191}; 00857 double eff_pt = 0.; 00858 eff_pt = wz_mu[ibin]; 00859 eff = eff_pt; 00860 } 00861 if(fabs(eta) > 0.1) { 00862 float wz_mu[] = {0.0,0.617,0.655,0.676,0.705,0.738,0.712,0.783,0.646,0.745}; 00863 //float ewz_mu[] = {0.0,0.043,0.0564,0.0448,0.0405,0.0576,0.065,0.0825,0.102,0.132}; 00864 double eff_pt = 0.; 00865 eff_pt = wz_mu[ibin]; 00866 eff = eff_pt; 00867 } 00868 } 00869 00870 if (flavor == 15) { // weight hadronic tau 1p 00871 double avgrate = 0.16; 00872 float wz_tau1p[] = {0.0,0.0311,0.148,0.229,0.217,0.292,0.245,0.307,0.227,0.277}; 00873 //float ewz_tau1p[] = {0.0,0.00211,0.0117,0.0179,0.0134,0.0248,0.0264,0.0322,0.0331,0.0427}; 00874 int ibin = 0; 00875 if(pt > 10 && pt < 15) ibin = 0; 00876 if(pt > 15 && pt < 20) ibin = 1; 00877 if(pt > 20 && pt < 25) ibin = 2; 00878 if(pt > 25 && pt < 30) ibin = 3; 00879 if(pt > 30 && pt < 40) ibin = 4; 00880 if(pt > 40 && pt < 50) ibin = 5; 00881 if(pt > 50 && pt < 60) ibin = 6; 00882 if(pt > 60 && pt < 80) ibin = 7; 00883 if(pt > 80 && pt < 100) ibin = 8; 00884 if(pt > 100) ibin = 9; 00885 double eff_pt = 0.; 00886 eff_pt = wz_tau1p[ibin]; 00887 00888 float wz_tau1p_eta[] = {0.166,0.15,0.188,0.175,0.142,0.109}; 00889 //float ewz_tau1p_eta[] ={0.0166,0.00853,0.0097,0.00985,0.00949,0.00842}; 00890 ibin = 0; 00891 if(eta > 0.0 && eta < 0.1) ibin = 0; 00892 if(eta > 0.1 && eta < 0.5) ibin = 1; 00893 if(eta > 0.5 && eta < 1.0) ibin = 2; 00894 if(eta > 1.0 && eta < 1.5) ibin = 3; 00895 if(eta > 1.5 && eta < 2.0) ibin = 4; 00896 if(eta > 2.0 && eta < 2.5) ibin = 5; 00897 double eff_eta = 0.; 00898 eff_eta = wz_tau1p_eta[ibin]; 00899 00900 eff = (eff_pt * eff_eta) / avgrate; 00901 } 00902 00903 return eff; 00904 } 00905 00906 /// Function giving observed and expected upper limits (on the visible cross-section) 00907 double getUpperLimit(const string& signal_region, bool observed) { 00908 00909 map<string,double> upperLimitsObserved; 00910 map<string,double> upperLimitsExpected; 00911 00912 upperLimitsObserved["HTlep_3l_offZ_OSSF_cut_0"] = 2.435; 00913 upperLimitsObserved["HTlep_3l_offZ_OSSF_cut_200"] = 0.704; 00914 upperLimitsObserved["HTlep_3l_offZ_OSSF_cut_500"] = 0.182; 00915 upperLimitsObserved["HTlep_3l_offZ_OSSF_cut_800"] = 0.147; 00916 upperLimitsObserved["HTlep_2ltau_offZ_OSSF_cut_0"] = 13.901; 00917 upperLimitsObserved["HTlep_2ltau_offZ_OSSF_cut_200"] = 1.677; 00918 upperLimitsObserved["HTlep_2ltau_offZ_OSSF_cut_500"] = 0.141; 00919 upperLimitsObserved["HTlep_2ltau_offZ_OSSF_cut_800"] = 0.155; 00920 upperLimitsObserved["HTlep_3l_offZ_noOSSF_cut_0"] = 1.054; 00921 upperLimitsObserved["HTlep_3l_offZ_noOSSF_cut_200"] = 0.341; 00922 upperLimitsObserved["HTlep_3l_offZ_noOSSF_cut_500"] = 0.221; 00923 upperLimitsObserved["HTlep_3l_offZ_noOSSF_cut_800"] = 0.140; 00924 upperLimitsObserved["HTlep_2ltau_offZ_noOSSF_cut_0"] = 4.276; 00925 upperLimitsObserved["HTlep_2ltau_offZ_noOSSF_cut_200"] = 0.413; 00926 upperLimitsObserved["HTlep_2ltau_offZ_noOSSF_cut_500"] = 0.138; 00927 upperLimitsObserved["HTlep_2ltau_offZ_noOSSF_cut_800"] = 0.150; 00928 upperLimitsObserved["HTlep_3l_onZ_cut_0"] = 29.804; 00929 upperLimitsObserved["HTlep_3l_onZ_cut_200"] = 3.579; 00930 upperLimitsObserved["HTlep_3l_onZ_cut_500"] = 0.466; 00931 upperLimitsObserved["HTlep_3l_onZ_cut_800"] = 0.298; 00932 upperLimitsObserved["HTlep_2ltau_onZ_cut_0"] = 205.091; 00933 upperLimitsObserved["HTlep_2ltau_onZ_cut_200"] = 3.141; 00934 upperLimitsObserved["HTlep_2ltau_onZ_cut_500"] = 0.290; 00935 upperLimitsObserved["HTlep_2ltau_onZ_cut_800"] = 0.157; 00936 upperLimitsObserved["METStrong_3l_offZ_OSSF_cut_0"] = 1.111; 00937 upperLimitsObserved["METStrong_3l_offZ_OSSF_cut_100"] = 0.354; 00938 upperLimitsObserved["METStrong_3l_offZ_OSSF_cut_200"] = 0.236; 00939 upperLimitsObserved["METStrong_3l_offZ_OSSF_cut_300"] = 0.150; 00940 upperLimitsObserved["METStrong_2ltau_offZ_OSSF_cut_0"] = 1.881; 00941 upperLimitsObserved["METStrong_2ltau_offZ_OSSF_cut_100"] = 0.406; 00942 upperLimitsObserved["METStrong_2ltau_offZ_OSSF_cut_200"] = 0.194; 00943 upperLimitsObserved["METStrong_2ltau_offZ_OSSF_cut_300"] = 0.134; 00944 upperLimitsObserved["METStrong_3l_offZ_noOSSF_cut_0"] = 0.770; 00945 upperLimitsObserved["METStrong_3l_offZ_noOSSF_cut_100"] = 0.295; 00946 upperLimitsObserved["METStrong_3l_offZ_noOSSF_cut_200"] = 0.149; 00947 upperLimitsObserved["METStrong_3l_offZ_noOSSF_cut_300"] = 0.140; 00948 upperLimitsObserved["METStrong_2ltau_offZ_noOSSF_cut_0"] = 2.003; 00949 upperLimitsObserved["METStrong_2ltau_offZ_noOSSF_cut_100"] = 0.806; 00950 upperLimitsObserved["METStrong_2ltau_offZ_noOSSF_cut_200"] = 0.227; 00951 upperLimitsObserved["METStrong_2ltau_offZ_noOSSF_cut_300"] = 0.138; 00952 upperLimitsObserved["METStrong_3l_onZ_cut_0"] = 6.383; 00953 upperLimitsObserved["METStrong_3l_onZ_cut_100"] = 0.959; 00954 upperLimitsObserved["METStrong_3l_onZ_cut_200"] = 0.549; 00955 upperLimitsObserved["METStrong_3l_onZ_cut_300"] = 0.182; 00956 upperLimitsObserved["METStrong_2ltau_onZ_cut_0"] = 10.658; 00957 upperLimitsObserved["METStrong_2ltau_onZ_cut_100"] = 0.637; 00958 upperLimitsObserved["METStrong_2ltau_onZ_cut_200"] = 0.291; 00959 upperLimitsObserved["METStrong_2ltau_onZ_cut_300"] = 0.227; 00960 upperLimitsObserved["METWeak_3l_offZ_OSSF_cut_0"] = 1.802; 00961 upperLimitsObserved["METWeak_3l_offZ_OSSF_cut_100"] = 0.344; 00962 upperLimitsObserved["METWeak_3l_offZ_OSSF_cut_200"] = 0.189; 00963 upperLimitsObserved["METWeak_3l_offZ_OSSF_cut_300"] = 0.148; 00964 upperLimitsObserved["METWeak_2ltau_offZ_OSSF_cut_0"] = 12.321; 00965 upperLimitsObserved["METWeak_2ltau_offZ_OSSF_cut_100"] = 0.430; 00966 upperLimitsObserved["METWeak_2ltau_offZ_OSSF_cut_200"] = 0.137; 00967 upperLimitsObserved["METWeak_2ltau_offZ_OSSF_cut_300"] = 0.134; 00968 upperLimitsObserved["METWeak_3l_offZ_noOSSF_cut_0"] = 0.562; 00969 upperLimitsObserved["METWeak_3l_offZ_noOSSF_cut_100"] = 0.153; 00970 upperLimitsObserved["METWeak_3l_offZ_noOSSF_cut_200"] = 0.154; 00971 upperLimitsObserved["METWeak_3l_offZ_noOSSF_cut_300"] = 0.141; 00972 upperLimitsObserved["METWeak_2ltau_offZ_noOSSF_cut_0"] = 2.475; 00973 upperLimitsObserved["METWeak_2ltau_offZ_noOSSF_cut_100"] = 0.244; 00974 upperLimitsObserved["METWeak_2ltau_offZ_noOSSF_cut_200"] = 0.141; 00975 upperLimitsObserved["METWeak_2ltau_offZ_noOSSF_cut_300"] = 0.142; 00976 upperLimitsObserved["METWeak_3l_onZ_cut_0"] = 24.769; 00977 upperLimitsObserved["METWeak_3l_onZ_cut_100"] = 0.690; 00978 upperLimitsObserved["METWeak_3l_onZ_cut_200"] = 0.198; 00979 upperLimitsObserved["METWeak_3l_onZ_cut_300"] = 0.138; 00980 upperLimitsObserved["METWeak_2ltau_onZ_cut_0"] = 194.360; 00981 upperLimitsObserved["METWeak_2ltau_onZ_cut_100"] = 0.287; 00982 upperLimitsObserved["METWeak_2ltau_onZ_cut_200"] = 0.144; 00983 upperLimitsObserved["METWeak_2ltau_onZ_cut_300"] = 0.130; 00984 upperLimitsObserved["Meff_3l_offZ_OSSF_cut_0"] = 2.435; 00985 upperLimitsObserved["Meff_3l_offZ_OSSF_cut_600"] = 0.487; 00986 upperLimitsObserved["Meff_3l_offZ_OSSF_cut_1000"] = 0.156; 00987 upperLimitsObserved["Meff_3l_offZ_OSSF_cut_1500"] = 0.140; 00988 upperLimitsObserved["Meff_2ltau_offZ_OSSF_cut_0"] = 13.901; 00989 upperLimitsObserved["Meff_2ltau_offZ_OSSF_cut_600"] = 0.687; 00990 upperLimitsObserved["Meff_2ltau_offZ_OSSF_cut_1000"] = 0.224; 00991 upperLimitsObserved["Meff_2ltau_offZ_OSSF_cut_1500"] = 0.155; 00992 upperLimitsObserved["Meff_3l_offZ_noOSSF_cut_0"] = 1.054; 00993 upperLimitsObserved["Meff_3l_offZ_noOSSF_cut_600"] = 0.249; 00994 upperLimitsObserved["Meff_3l_offZ_noOSSF_cut_1000"] = 0.194; 00995 upperLimitsObserved["Meff_3l_offZ_noOSSF_cut_1500"] = 0.145; 00996 upperLimitsObserved["Meff_2ltau_offZ_noOSSF_cut_0"] = 4.276; 00997 upperLimitsObserved["Meff_2ltau_offZ_noOSSF_cut_600"] = 0.772; 00998 upperLimitsObserved["Meff_2ltau_offZ_noOSSF_cut_1000"] = 0.218; 00999 upperLimitsObserved["Meff_2ltau_offZ_noOSSF_cut_1500"] = 0.204; 01000 upperLimitsObserved["Meff_3l_onZ_cut_0"] = 29.804; 01001 upperLimitsObserved["Meff_3l_onZ_cut_600"] = 2.933; 01002 upperLimitsObserved["Meff_3l_onZ_cut_1000"] = 0.912; 01003 upperLimitsObserved["Meff_3l_onZ_cut_1500"] = 0.225; 01004 upperLimitsObserved["Meff_2ltau_onZ_cut_0"] = 205.091; 01005 upperLimitsObserved["Meff_2ltau_onZ_cut_600"] = 1.486; 01006 upperLimitsObserved["Meff_2ltau_onZ_cut_1000"] = 0.641; 01007 upperLimitsObserved["Meff_2ltau_onZ_cut_1500"] = 0.204; 01008 upperLimitsObserved["MeffStrong_3l_offZ_OSSF_cut_0"] = 0.479; 01009 upperLimitsObserved["MeffStrong_3l_offZ_OSSF_cut_600"] = 0.353; 01010 upperLimitsObserved["MeffStrong_3l_offZ_OSSF_cut_1200"] = 0.187; 01011 upperLimitsObserved["MeffStrong_2ltau_offZ_OSSF_cut_0"] = 0.617; 01012 upperLimitsObserved["MeffStrong_2ltau_offZ_OSSF_cut_600"] = 0.320; 01013 upperLimitsObserved["MeffStrong_2ltau_offZ_OSSF_cut_1200"] = 0.281; 01014 upperLimitsObserved["MeffStrong_3l_offZ_noOSSF_cut_0"] = 0.408; 01015 upperLimitsObserved["MeffStrong_3l_offZ_noOSSF_cut_600"] = 0.240; 01016 upperLimitsObserved["MeffStrong_3l_offZ_noOSSF_cut_1200"] = 0.150; 01017 upperLimitsObserved["MeffStrong_2ltau_offZ_noOSSF_cut_0"] = 0.774; 01018 upperLimitsObserved["MeffStrong_2ltau_offZ_noOSSF_cut_600"] = 0.417; 01019 upperLimitsObserved["MeffStrong_2ltau_offZ_noOSSF_cut_1200"] = 0.266; 01020 upperLimitsObserved["MeffStrong_3l_onZ_cut_0"] = 1.208; 01021 upperLimitsObserved["MeffStrong_3l_onZ_cut_600"] = 0.837; 01022 upperLimitsObserved["MeffStrong_3l_onZ_cut_1200"] = 0.269; 01023 upperLimitsObserved["MeffStrong_2ltau_onZ_cut_0"] = 0.605; 01024 upperLimitsObserved["MeffStrong_2ltau_onZ_cut_600"] = 0.420; 01025 upperLimitsObserved["MeffStrong_2ltau_onZ_cut_1200"] = 0.141; 01026 upperLimitsObserved["MeffMt_3l_onZ_cut_0"] = 1.832; 01027 upperLimitsObserved["MeffMt_3l_onZ_cut_600"] = 0.862; 01028 upperLimitsObserved["MeffMt_3l_onZ_cut_1200"] = 0.222; 01029 upperLimitsObserved["MeffMt_2ltau_onZ_cut_0"] = 1.309; 01030 upperLimitsObserved["MeffMt_2ltau_onZ_cut_600"] = 0.481; 01031 upperLimitsObserved["MeffMt_2ltau_onZ_cut_1200"] = 0.146; 01032 upperLimitsObserved["MinPt_3l_offZ_OSSF_cut_0"] = 2.435; 01033 upperLimitsObserved["MinPt_3l_offZ_OSSF_cut_50"] = 0.500; 01034 upperLimitsObserved["MinPt_3l_offZ_OSSF_cut_100"] = 0.203; 01035 upperLimitsObserved["MinPt_3l_offZ_OSSF_cut_150"] = 0.128; 01036 upperLimitsObserved["MinPt_2ltau_offZ_OSSF_cut_0"] = 13.901; 01037 upperLimitsObserved["MinPt_2ltau_offZ_OSSF_cut_50"] = 0.859; 01038 upperLimitsObserved["MinPt_2ltau_offZ_OSSF_cut_100"] = 0.158; 01039 upperLimitsObserved["MinPt_2ltau_offZ_OSSF_cut_150"] = 0.155; 01040 upperLimitsObserved["MinPt_3l_offZ_noOSSF_cut_0"] = 1.054; 01041 upperLimitsObserved["MinPt_3l_offZ_noOSSF_cut_50"] = 0.295; 01042 upperLimitsObserved["MinPt_3l_offZ_noOSSF_cut_100"] = 0.148; 01043 upperLimitsObserved["MinPt_3l_offZ_noOSSF_cut_150"] = 0.137; 01044 upperLimitsObserved["MinPt_2ltau_offZ_noOSSF_cut_0"] = 4.276; 01045 upperLimitsObserved["MinPt_2ltau_offZ_noOSSF_cut_50"] = 0.314; 01046 upperLimitsObserved["MinPt_2ltau_offZ_noOSSF_cut_100"] = 0.134; 01047 upperLimitsObserved["MinPt_2ltau_offZ_noOSSF_cut_150"] = 0.140; 01048 upperLimitsObserved["MinPt_3l_onZ_cut_0"] = 29.804; 01049 upperLimitsObserved["MinPt_3l_onZ_cut_50"] = 1.767; 01050 upperLimitsObserved["MinPt_3l_onZ_cut_100"] = 0.690; 01051 upperLimitsObserved["MinPt_3l_onZ_cut_150"] = 0.301; 01052 upperLimitsObserved["MinPt_2ltau_onZ_cut_0"] = 205.091; 01053 upperLimitsObserved["MinPt_2ltau_onZ_cut_50"] = 1.050; 01054 upperLimitsObserved["MinPt_2ltau_onZ_cut_100"] = 0.155; 01055 upperLimitsObserved["MinPt_2ltau_onZ_cut_150"] = 0.146; 01056 upperLimitsObserved["nbtag_3l_offZ_OSSF_cut_0"] = 2.435; 01057 upperLimitsObserved["nbtag_3l_offZ_OSSF_cut_1"] = 0.865; 01058 upperLimitsObserved["nbtag_3l_offZ_OSSF_cut_2"] = 0.474; 01059 upperLimitsObserved["nbtag_2ltau_offZ_OSSF_cut_0"] = 13.901; 01060 upperLimitsObserved["nbtag_2ltau_offZ_OSSF_cut_1"] = 1.566; 01061 upperLimitsObserved["nbtag_2ltau_offZ_OSSF_cut_2"] = 0.426; 01062 upperLimitsObserved["nbtag_3l_offZ_noOSSF_cut_0"] = 1.054; 01063 upperLimitsObserved["nbtag_3l_offZ_noOSSF_cut_1"] = 0.643; 01064 upperLimitsObserved["nbtag_3l_offZ_noOSSF_cut_2"] = 0.321; 01065 upperLimitsObserved["nbtag_2ltau_offZ_noOSSF_cut_0"] = 4.276; 01066 upperLimitsObserved["nbtag_2ltau_offZ_noOSSF_cut_1"] = 2.435; 01067 upperLimitsObserved["nbtag_2ltau_offZ_noOSSF_cut_2"] = 1.073; 01068 upperLimitsObserved["nbtag_3l_onZ_cut_0"] = 29.804; 01069 upperLimitsObserved["nbtag_3l_onZ_cut_1"] = 3.908; 01070 upperLimitsObserved["nbtag_3l_onZ_cut_2"] = 0.704; 01071 upperLimitsObserved["nbtag_2ltau_onZ_cut_0"] = 205.091; 01072 upperLimitsObserved["nbtag_2ltau_onZ_cut_1"] = 9.377; 01073 upperLimitsObserved["nbtag_2ltau_onZ_cut_2"] = 0.657; 01074 upperLimitsExpected["HTlep_3l_offZ_OSSF_cut_0"] = 2.893; 01075 upperLimitsExpected["HTlep_3l_offZ_OSSF_cut_200"] = 1.175; 01076 upperLimitsExpected["HTlep_3l_offZ_OSSF_cut_500"] = 0.265; 01077 upperLimitsExpected["HTlep_3l_offZ_OSSF_cut_800"] = 0.155; 01078 upperLimitsExpected["HTlep_2ltau_offZ_OSSF_cut_0"] = 14.293; 01079 upperLimitsExpected["HTlep_2ltau_offZ_OSSF_cut_200"] = 1.803; 01080 upperLimitsExpected["HTlep_2ltau_offZ_OSSF_cut_500"] = 0.159; 01081 upperLimitsExpected["HTlep_2ltau_offZ_OSSF_cut_800"] = 0.155; 01082 upperLimitsExpected["HTlep_3l_offZ_noOSSF_cut_0"] = 0.836; 01083 upperLimitsExpected["HTlep_3l_offZ_noOSSF_cut_200"] = 0.340; 01084 upperLimitsExpected["HTlep_3l_offZ_noOSSF_cut_500"] = 0.218; 01085 upperLimitsExpected["HTlep_3l_offZ_noOSSF_cut_800"] = 0.140; 01086 upperLimitsExpected["HTlep_2ltau_offZ_noOSSF_cut_0"] = 4.132; 01087 upperLimitsExpected["HTlep_2ltau_offZ_noOSSF_cut_200"] = 0.599; 01088 upperLimitsExpected["HTlep_2ltau_offZ_noOSSF_cut_500"] = 0.146; 01089 upperLimitsExpected["HTlep_2ltau_offZ_noOSSF_cut_800"] = 0.148; 01090 upperLimitsExpected["HTlep_3l_onZ_cut_0"] = 32.181; 01091 upperLimitsExpected["HTlep_3l_onZ_cut_200"] = 4.879; 01092 upperLimitsExpected["HTlep_3l_onZ_cut_500"] = 0.473; 01093 upperLimitsExpected["HTlep_3l_onZ_cut_800"] = 0.266; 01094 upperLimitsExpected["HTlep_2ltau_onZ_cut_0"] = 217.801; 01095 upperLimitsExpected["HTlep_2ltau_onZ_cut_200"] = 3.676; 01096 upperLimitsExpected["HTlep_2ltau_onZ_cut_500"] = 0.235; 01097 upperLimitsExpected["HTlep_2ltau_onZ_cut_800"] = 0.150; 01098 upperLimitsExpected["METStrong_3l_offZ_OSSF_cut_0"] = 1.196; 01099 upperLimitsExpected["METStrong_3l_offZ_OSSF_cut_100"] = 0.423; 01100 upperLimitsExpected["METStrong_3l_offZ_OSSF_cut_200"] = 0.208; 01101 upperLimitsExpected["METStrong_3l_offZ_OSSF_cut_300"] = 0.158; 01102 upperLimitsExpected["METStrong_2ltau_offZ_OSSF_cut_0"] = 2.158; 01103 upperLimitsExpected["METStrong_2ltau_offZ_OSSF_cut_100"] = 0.461; 01104 upperLimitsExpected["METStrong_2ltau_offZ_OSSF_cut_200"] = 0.186; 01105 upperLimitsExpected["METStrong_2ltau_offZ_OSSF_cut_300"] = 0.138; 01106 upperLimitsExpected["METStrong_3l_offZ_noOSSF_cut_0"] = 0.495; 01107 upperLimitsExpected["METStrong_3l_offZ_noOSSF_cut_100"] = 0.284; 01108 upperLimitsExpected["METStrong_3l_offZ_noOSSF_cut_200"] = 0.150; 01109 upperLimitsExpected["METStrong_3l_offZ_noOSSF_cut_300"] = 0.146; 01110 upperLimitsExpected["METStrong_2ltau_offZ_noOSSF_cut_0"] = 1.967; 01111 upperLimitsExpected["METStrong_2ltau_offZ_noOSSF_cut_100"] = 0.732; 01112 upperLimitsExpected["METStrong_2ltau_offZ_noOSSF_cut_200"] = 0.225; 01113 upperLimitsExpected["METStrong_2ltau_offZ_noOSSF_cut_300"] = 0.147; 01114 upperLimitsExpected["METStrong_3l_onZ_cut_0"] = 7.157; 01115 upperLimitsExpected["METStrong_3l_onZ_cut_100"] = 1.342; 01116 upperLimitsExpected["METStrong_3l_onZ_cut_200"] = 0.508; 01117 upperLimitsExpected["METStrong_3l_onZ_cut_300"] = 0.228; 01118 upperLimitsExpected["METStrong_2ltau_onZ_cut_0"] = 12.441; 01119 upperLimitsExpected["METStrong_2ltau_onZ_cut_100"] = 0.534; 01120 upperLimitsExpected["METStrong_2ltau_onZ_cut_200"] = 0.243; 01121 upperLimitsExpected["METStrong_2ltau_onZ_cut_300"] = 0.218; 01122 upperLimitsExpected["METWeak_3l_offZ_OSSF_cut_0"] = 2.199; 01123 upperLimitsExpected["METWeak_3l_offZ_OSSF_cut_100"] = 0.391; 01124 upperLimitsExpected["METWeak_3l_offZ_OSSF_cut_200"] = 0.177; 01125 upperLimitsExpected["METWeak_3l_offZ_OSSF_cut_300"] = 0.144; 01126 upperLimitsExpected["METWeak_2ltau_offZ_OSSF_cut_0"] = 12.431; 01127 upperLimitsExpected["METWeak_2ltau_offZ_OSSF_cut_100"] = 0.358; 01128 upperLimitsExpected["METWeak_2ltau_offZ_OSSF_cut_200"] = 0.150; 01129 upperLimitsExpected["METWeak_2ltau_offZ_OSSF_cut_300"] = 0.135; 01130 upperLimitsExpected["METWeak_3l_offZ_noOSSF_cut_0"] = 0.577; 01131 upperLimitsExpected["METWeak_3l_offZ_noOSSF_cut_100"] = 0.214; 01132 upperLimitsExpected["METWeak_3l_offZ_noOSSF_cut_200"] = 0.155; 01133 upperLimitsExpected["METWeak_3l_offZ_noOSSF_cut_300"] = 0.140; 01134 upperLimitsExpected["METWeak_2ltau_offZ_noOSSF_cut_0"] = 2.474; 01135 upperLimitsExpected["METWeak_2ltau_offZ_noOSSF_cut_100"] = 0.382; 01136 upperLimitsExpected["METWeak_2ltau_offZ_noOSSF_cut_200"] = 0.144; 01137 upperLimitsExpected["METWeak_2ltau_offZ_noOSSF_cut_300"] = 0.146; 01138 upperLimitsExpected["METWeak_3l_onZ_cut_0"] = 26.305; 01139 upperLimitsExpected["METWeak_3l_onZ_cut_100"] = 1.227; 01140 upperLimitsExpected["METWeak_3l_onZ_cut_200"] = 0.311; 01141 upperLimitsExpected["METWeak_3l_onZ_cut_300"] = 0.188; 01142 upperLimitsExpected["METWeak_2ltau_onZ_cut_0"] = 205.198; 01143 upperLimitsExpected["METWeak_2ltau_onZ_cut_100"] = 0.399; 01144 upperLimitsExpected["METWeak_2ltau_onZ_cut_200"] = 0.166; 01145 upperLimitsExpected["METWeak_2ltau_onZ_cut_300"] = 0.140; 01146 upperLimitsExpected["Meff_3l_offZ_OSSF_cut_0"] = 2.893; 01147 upperLimitsExpected["Meff_3l_offZ_OSSF_cut_600"] = 0.649; 01148 upperLimitsExpected["Meff_3l_offZ_OSSF_cut_1000"] = 0.252; 01149 upperLimitsExpected["Meff_3l_offZ_OSSF_cut_1500"] = 0.150; 01150 upperLimitsExpected["Meff_2ltau_offZ_OSSF_cut_0"] = 14.293; 01151 upperLimitsExpected["Meff_2ltau_offZ_OSSF_cut_600"] = 0.657; 01152 upperLimitsExpected["Meff_2ltau_offZ_OSSF_cut_1000"] = 0.226; 01153 upperLimitsExpected["Meff_2ltau_offZ_OSSF_cut_1500"] = 0.154; 01154 upperLimitsExpected["Meff_3l_offZ_noOSSF_cut_0"] = 0.836; 01155 upperLimitsExpected["Meff_3l_offZ_noOSSF_cut_600"] = 0.265; 01156 upperLimitsExpected["Meff_3l_offZ_noOSSF_cut_1000"] = 0.176; 01157 upperLimitsExpected["Meff_3l_offZ_noOSSF_cut_1500"] = 0.146; 01158 upperLimitsExpected["Meff_2ltau_offZ_noOSSF_cut_0"] = 4.132; 01159 upperLimitsExpected["Meff_2ltau_offZ_noOSSF_cut_600"] = 0.678; 01160 upperLimitsExpected["Meff_2ltau_offZ_noOSSF_cut_1000"] = 0.243; 01161 upperLimitsExpected["Meff_2ltau_offZ_noOSSF_cut_1500"] = 0.184; 01162 upperLimitsExpected["Meff_3l_onZ_cut_0"] = 32.181; 01163 upperLimitsExpected["Meff_3l_onZ_cut_600"] = 3.219; 01164 upperLimitsExpected["Meff_3l_onZ_cut_1000"] = 0.905; 01165 upperLimitsExpected["Meff_3l_onZ_cut_1500"] = 0.261; 01166 upperLimitsExpected["Meff_2ltau_onZ_cut_0"] = 217.801; 01167 upperLimitsExpected["Meff_2ltau_onZ_cut_600"] = 1.680; 01168 upperLimitsExpected["Meff_2ltau_onZ_cut_1000"] = 0.375; 01169 upperLimitsExpected["Meff_2ltau_onZ_cut_1500"] = 0.178; 01170 upperLimitsExpected["MeffStrong_3l_offZ_OSSF_cut_0"] = 0.571; 01171 upperLimitsExpected["MeffStrong_3l_offZ_OSSF_cut_600"] = 0.386; 01172 upperLimitsExpected["MeffStrong_3l_offZ_OSSF_cut_1200"] = 0.177; 01173 upperLimitsExpected["MeffStrong_2ltau_offZ_OSSF_cut_0"] = 0.605; 01174 upperLimitsExpected["MeffStrong_2ltau_offZ_OSSF_cut_600"] = 0.335; 01175 upperLimitsExpected["MeffStrong_2ltau_offZ_OSSF_cut_1200"] = 0.249; 01176 upperLimitsExpected["MeffStrong_3l_offZ_noOSSF_cut_0"] = 0.373; 01177 upperLimitsExpected["MeffStrong_3l_offZ_noOSSF_cut_600"] = 0.223; 01178 upperLimitsExpected["MeffStrong_3l_offZ_noOSSF_cut_1200"] = 0.150; 01179 upperLimitsExpected["MeffStrong_2ltau_offZ_noOSSF_cut_0"] = 0.873; 01180 upperLimitsExpected["MeffStrong_2ltau_offZ_noOSSF_cut_600"] = 0.428; 01181 upperLimitsExpected["MeffStrong_2ltau_offZ_noOSSF_cut_1200"] = 0.210; 01182 upperLimitsExpected["MeffStrong_3l_onZ_cut_0"] = 2.034; 01183 upperLimitsExpected["MeffStrong_3l_onZ_cut_600"] = 1.093; 01184 upperLimitsExpected["MeffStrong_3l_onZ_cut_1200"] = 0.293; 01185 upperLimitsExpected["MeffStrong_2ltau_onZ_cut_0"] = 0.690; 01186 upperLimitsExpected["MeffStrong_2ltau_onZ_cut_600"] = 0.392; 01187 upperLimitsExpected["MeffStrong_2ltau_onZ_cut_1200"] = 0.156; 01188 upperLimitsExpected["MeffMt_3l_onZ_cut_0"] = 2.483; 01189 upperLimitsExpected["MeffMt_3l_onZ_cut_600"] = 0.845; 01190 upperLimitsExpected["MeffMt_3l_onZ_cut_1200"] = 0.255; 01191 upperLimitsExpected["MeffMt_2ltau_onZ_cut_0"] = 1.448; 01192 upperLimitsExpected["MeffMt_2ltau_onZ_cut_600"] = 0.391; 01193 upperLimitsExpected["MeffMt_2ltau_onZ_cut_1200"] = 0.146; 01194 upperLimitsExpected["MinPt_3l_offZ_OSSF_cut_0"] = 2.893; 01195 upperLimitsExpected["MinPt_3l_offZ_OSSF_cut_50"] = 0.703; 01196 upperLimitsExpected["MinPt_3l_offZ_OSSF_cut_100"] = 0.207; 01197 upperLimitsExpected["MinPt_3l_offZ_OSSF_cut_150"] = 0.143; 01198 upperLimitsExpected["MinPt_2ltau_offZ_OSSF_cut_0"] = 14.293; 01199 upperLimitsExpected["MinPt_2ltau_offZ_OSSF_cut_50"] = 0.705; 01200 upperLimitsExpected["MinPt_2ltau_offZ_OSSF_cut_100"] = 0.149; 01201 upperLimitsExpected["MinPt_2ltau_offZ_OSSF_cut_150"] = 0.155; 01202 upperLimitsExpected["MinPt_3l_offZ_noOSSF_cut_0"] = 0.836; 01203 upperLimitsExpected["MinPt_3l_offZ_noOSSF_cut_50"] = 0.249; 01204 upperLimitsExpected["MinPt_3l_offZ_noOSSF_cut_100"] = 0.135; 01205 upperLimitsExpected["MinPt_3l_offZ_noOSSF_cut_150"] = 0.136; 01206 upperLimitsExpected["MinPt_2ltau_offZ_noOSSF_cut_0"] = 4.132; 01207 upperLimitsExpected["MinPt_2ltau_offZ_noOSSF_cut_50"] = 0.339; 01208 upperLimitsExpected["MinPt_2ltau_offZ_noOSSF_cut_100"] = 0.149; 01209 upperLimitsExpected["MinPt_2ltau_offZ_noOSSF_cut_150"] = 0.145; 01210 upperLimitsExpected["MinPt_3l_onZ_cut_0"] = 32.181; 01211 upperLimitsExpected["MinPt_3l_onZ_cut_50"] = 2.260; 01212 upperLimitsExpected["MinPt_3l_onZ_cut_100"] = 0.438; 01213 upperLimitsExpected["MinPt_3l_onZ_cut_150"] = 0.305; 01214 upperLimitsExpected["MinPt_2ltau_onZ_cut_0"] = 217.801; 01215 upperLimitsExpected["MinPt_2ltau_onZ_cut_50"] = 1.335; 01216 upperLimitsExpected["MinPt_2ltau_onZ_cut_100"] = 0.162; 01217 upperLimitsExpected["MinPt_2ltau_onZ_cut_150"] = 0.149; 01218 upperLimitsExpected["nbtag_3l_offZ_OSSF_cut_0"] = 2.893; 01219 upperLimitsExpected["nbtag_3l_offZ_OSSF_cut_1"] = 0.923; 01220 upperLimitsExpected["nbtag_3l_offZ_OSSF_cut_2"] = 0.452; 01221 upperLimitsExpected["nbtag_2ltau_offZ_OSSF_cut_0"] = 14.293; 01222 upperLimitsExpected["nbtag_2ltau_offZ_OSSF_cut_1"] = 1.774; 01223 upperLimitsExpected["nbtag_2ltau_offZ_OSSF_cut_2"] = 0.549; 01224 upperLimitsExpected["nbtag_3l_offZ_noOSSF_cut_0"] = 0.836; 01225 upperLimitsExpected["nbtag_3l_offZ_noOSSF_cut_1"] = 0.594; 01226 upperLimitsExpected["nbtag_3l_offZ_noOSSF_cut_2"] = 0.298; 01227 upperLimitsExpected["nbtag_2ltau_offZ_noOSSF_cut_0"] = 4.132; 01228 upperLimitsExpected["nbtag_2ltau_offZ_noOSSF_cut_1"] = 2.358; 01229 upperLimitsExpected["nbtag_2ltau_offZ_noOSSF_cut_2"] = 0.958; 01230 upperLimitsExpected["nbtag_3l_onZ_cut_0"] = 32.181; 01231 upperLimitsExpected["nbtag_3l_onZ_cut_1"] = 3.868; 01232 upperLimitsExpected["nbtag_3l_onZ_cut_2"] = 0.887; 01233 upperLimitsExpected["nbtag_2ltau_onZ_cut_0"] = 217.801; 01234 upperLimitsExpected["nbtag_2ltau_onZ_cut_1"] = 9.397; 01235 upperLimitsExpected["nbtag_2ltau_onZ_cut_2"] = 0.787; 01236 01237 01238 01239 if (observed) return upperLimitsObserved[signal_region]; 01240 else return upperLimitsExpected[signal_region]; 01241 } 01242 01243 01244 /// Function checking if there is an OSSF lepton pair or a combination of 3 leptons with an invariant mass close to the Z mass 01245 int isonZ (const Particles& particles) { 01246 int onZ = 0; 01247 double best_mass_2 = 999.; 01248 double best_mass_3 = 999.; 01249 01250 // Loop over all 2 particle combinations to find invariant mass of OSSF pair closest to Z mass 01251 foreach (const Particle& p1, particles) { 01252 foreach (const Particle& p2, particles) { 01253 double mass_difference_2_old = fabs(91.0 - best_mass_2); 01254 double mass_difference_2_new = fabs(91.0 - (p1.momentum() + p2.momentum()).mass()/GeV); 01255 01256 // If particle combination is OSSF pair calculate mass difference to Z mass 01257 if ((p1.pid()*p2.pid() == -121 || p1.pid()*p2.pid() == -169)) { 01258 01259 // Get invariant mass closest to Z mass 01260 if (mass_difference_2_new < mass_difference_2_old) 01261 best_mass_2 = (p1.momentum() + p2.momentum()).mass()/GeV; 01262 // In case there is an OSSF pair take also 3rd lepton into account (e.g. from FSR and photon to electron conversion) 01263 foreach (const Particle& p3 , particles ) { 01264 double mass_difference_3_old = fabs(91.0 - best_mass_3); 01265 double mass_difference_3_new = fabs(91.0 - (p1.momentum() + p2.momentum() + p3.momentum()).mass()/GeV); 01266 if (mass_difference_3_new < mass_difference_3_old) 01267 best_mass_3 = (p1.momentum() + p2.momentum() + p3.momentum()).mass()/GeV; 01268 } 01269 } 01270 } 01271 } 01272 01273 // Pick the minimum invariant mass of the best OSSF pair combination and the best 3 lepton combination 01274 double best_mass = min(best_mass_2,best_mass_3); 01275 // if this mass is in a 20 GeV window around the Z mass, the event is classified as onZ 01276 if ( fabs(91.0 - best_mass) < 20. ) onZ = 1; 01277 return onZ; 01278 } 01279 01280 /// function checking if two leptons are an OSSF lepton pair and giving out the invariant mass (0 if no OSSF pair) 01281 double isOSSF_mass (const Particle& p1, const Particle& p2) { 01282 double inv_mass = 0.; 01283 // Is particle combination OSSF pair? 01284 if ((p1.pid()*p2.pid() == -121 || p1.pid()*p2.pid() == -169)) { 01285 // Get invariant mass 01286 inv_mass = (p1.momentum() + p2.momentum()).mass()/GeV; 01287 } 01288 return inv_mass; 01289 } 01290 01291 /// Function checking if there is an OSSF lepton pair 01292 bool isOSSF (const Particles& particles) { 01293 for (size_t i1=0 ; i1 < 3 ; i1 ++) { 01294 for (size_t i2 = i1+1 ; i2 < 3 ; i2 ++) { 01295 if ((particles[i1].pid()*particles[i2].pid() == -121 || particles[i1].pid()*particles[i2].pid() == -169)) { 01296 return true; 01297 } 01298 } 01299 } 01300 return false; 01301 } 01302 01303 //@} 01304 01305 private: 01306 01307 /// Histograms 01308 //@{ 01309 Histo1DPtr _h_HTlep_all, _h_HTjets_all, _h_MET_all, _h_Meff_all, _h_min_pT_all, _h_mT_all; 01310 Histo1DPtr _h_pt_1_3l, _h_pt_2_3l, _h_pt_3_3l, _h_pt_1_2ltau, _h_pt_2_2ltau, _h_pt_3_2ltau; 01311 Histo1DPtr _h_e_n, _h_mu_n, _h_tau_n; 01312 Histo1DPtr _h_excluded; 01313 //@} 01314 01315 /// Fiducial efficiencies to model the effects of the ATLAS detector 01316 bool _use_fiducial_lepton_efficiency; 01317 01318 /// List of signal regions and event counts per signal region 01319 vector<string> _signal_regions; 01320 map<string, double> _eventCountsPerSR; 01321 01322 }; 01323 01324 01325 DECLARE_RIVET_PLUGIN(ATLAS_2014_I1327229); 01326 01327 } 01328 Generated on Wed Oct 7 2015 12:09:11 for The Rivet MC analysis system by ![]() |