ATLAS_2012_I1204447.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 namespace Rivet { 00012 00013 00014 class ATLAS_2012_I1204447 : public Analysis { 00015 public: 00016 00017 /// Constructor 00018 ATLAS_2012_I1204447() 00019 : Analysis("ATLAS_2012_I1204447") { } 00020 00021 00022 /// Book histograms and initialise projections before the run 00023 void init() { 00024 00025 // To calculate the acceptance without having the fiducial lepton efficiencies included, this part can be turned off 00026 _use_fiducial_lepton_efficiency = true; 00027 00028 // Random numbers for simulation of ATLAS detector reconstruction efficiency 00029 srand(160385); 00030 00031 // Read in all signal regions 00032 _signal_regions = getSignalRegions(); 00033 00034 // Set number of events per signal region to 0 00035 for (size_t i = 0; i < _signal_regions.size(); i++) 00036 _eventCountsPerSR[_signal_regions[i]] = 0.0; 00037 00038 // Final state including all charged and neutral particles 00039 const FinalState fs(-5.0, 5.0, 1*GeV); 00040 addProjection(fs, "FS"); 00041 00042 // Final state including all charged particles 00043 addProjection(ChargedFinalState(Cuts::abseta < 2.5 && Cuts::pT > 1*GeV), "CFS"); 00044 00045 // Final state including all visible particles (to calculate MET, Jets etc.) 00046 addProjection(VisibleFinalState(Cuts::abseta < 5.0), "VFS"); 00047 00048 // Final state including all AntiKt 04 Jets 00049 VetoedFinalState vfs; 00050 vfs.addVetoPairId(PID::MUON); 00051 addProjection(FastJets(vfs, FastJets::ANTIKT, 0.4), "AntiKtJets04"); 00052 00053 // Final state including all unstable particles (including taus) 00054 addProjection(UnstableFinalState(Cuts::abseta < 5.0 && Cuts::pT > 5*GeV), "UFS"); 00055 00056 // Final state including all electrons 00057 IdentifiedFinalState elecs(Cuts::abseta < 2.47 && Cuts::pT > 10*GeV); 00058 elecs.acceptIdPair(PID::ELECTRON); 00059 addProjection(elecs, "elecs"); 00060 00061 // Final state including all muons 00062 IdentifiedFinalState muons(Cuts::abseta < 2.5 && Cuts::pT > 10*GeV); 00063 muons.acceptIdPair(PID::MUON); 00064 addProjection(muons, "muons"); 00065 00066 // Book histograms 00067 _h_HTlep_all = bookHisto1D("HTlep_all" , 30, 0, 1500); 00068 _h_HTjets_all = bookHisto1D("HTjets_all", 30, 0, 1500); 00069 _h_MET_all = bookHisto1D("MET_all" , 20, 0, 1000); 00070 _h_Meff_all = bookHisto1D("Meff_all" , 30, 0, 3000); 00071 00072 _h_e_n = bookHisto1D("e_n" , 10, -0.5, 9.5); 00073 _h_mu_n = bookHisto1D("mu_n" , 10, -0.5, 9.5); 00074 _h_tau_n = bookHisto1D("tau_n", 10, -0.5, 9.5); 00075 00076 _h_pt_1_3l = bookHisto1D("pt_1_3l", 100, 0, 2000); 00077 _h_pt_2_3l = bookHisto1D("pt_2_3l", 100, 0, 2000); 00078 _h_pt_3_3l = bookHisto1D("pt_3_3l", 100, 0, 2000); 00079 _h_pt_1_2ltau = bookHisto1D("pt_1_2ltau", 100, 0, 2000); 00080 _h_pt_2_2ltau = bookHisto1D("pt_2_2ltau", 100, 0, 2000); 00081 _h_pt_3_2ltau = bookHisto1D("pt_3_2ltau", 100, 0, 2000); 00082 00083 _h_excluded = bookHisto1D("excluded", 2, -0.5, 1.5); 00084 } 00085 00086 00087 /// Perform the per-event analysis 00088 void analyze(const Event& event) { 00089 00090 // Muons 00091 Particles muon_candidates; 00092 const Particles charged_tracks = applyProjection<ChargedFinalState>(event, "CFS").particles(); 00093 const Particles visible_particles = applyProjection<VisibleFinalState>(event, "VFS").particles(); 00094 foreach (const Particle& mu, applyProjection<IdentifiedFinalState>(event, "muons").particlesByPt()) { 00095 // Calculate pTCone30 variable (pT of all tracks within dR<0.3 - pT of muon itself) 00096 double pTinCone = -mu.pT(); 00097 foreach (const Particle& track, charged_tracks) { 00098 if (deltaR(mu.momentum(), track.momentum()) < 0.3) 00099 pTinCone += track.pT(); 00100 } 00101 00102 // Calculate eTCone30 variable (pT of all visible particles within dR<0.3) 00103 double eTinCone = 0.; 00104 foreach (const Particle& visible_particle, visible_particles) { 00105 if (visible_particle.abspid() != PID::MUON && inRange(deltaR(mu.momentum(), visible_particle.momentum()), 0.1, 0.3)) 00106 eTinCone += visible_particle.pT(); 00107 } 00108 00109 // Apply reconstruction efficiency and simulate reco 00110 int muon_id = 13; 00111 if ( mu.hasAncestor(15) || mu.hasAncestor(-15)) muon_id = 14; 00112 const double eff = (_use_fiducial_lepton_efficiency) ? apply_reco_eff(muon_id, mu) : 1.0; 00113 const bool keep_muon = rand()/static_cast<double>(RAND_MAX) <= eff; 00114 00115 // Keep muon if pTCone30/pT < 0.15 and eTCone30/pT < 0.2 and reconstructed 00116 if (keep_muon && pTinCone/mu.pT() <= 0.15 && eTinCone/mu.pT() < 0.2) 00117 muon_candidates.push_back(mu); 00118 } 00119 00120 00121 // Electrons 00122 Particles electron_candidates; 00123 foreach (const Particle& e, applyProjection<IdentifiedFinalState>(event, "elecs").particlesByPt()) { 00124 // Neglect electrons in crack regions 00125 if (inRange(e.abseta(), 1.37, 1.52)) continue; 00126 00127 // Calculate pTCone30 variable (pT of all tracks within dR<0.3 - pT of electron itself) 00128 double pTinCone = -e.pT(); 00129 foreach (const Particle& track, charged_tracks) { 00130 if (deltaR(e.momentum(), track.momentum()) < 0.3) pTinCone += track.pT(); 00131 } 00132 00133 // Calculate eTCone30 variable (pT of all visible particles (except muons) within dR<0.3) 00134 double eTinCone = 0.; 00135 foreach (const Particle& visible_particle, visible_particles) { 00136 if (visible_particle.abspid() != PID::MUON && inRange(deltaR(e.momentum(), visible_particle.momentum()), 0.1, 0.3)) 00137 eTinCone += visible_particle.pT(); 00138 } 00139 00140 // Apply reconstruction efficiency and simulate reco 00141 int elec_id = 11; 00142 if (e.hasAncestor(15) || e.hasAncestor(-15)) elec_id = 12; 00143 const double eff = (_use_fiducial_lepton_efficiency) ? apply_reco_eff(elec_id, e) : 1.0; 00144 const bool keep_elec = rand()/static_cast<double>(RAND_MAX) <= eff; 00145 00146 // Keep electron if pTCone30/pT < 0.13 and eTCone30/pT < 0.2 and reconstructed 00147 if (keep_elec && pTinCone/e.pT() <= 0.13 && eTinCone/e.pT() < 0.2) 00148 electron_candidates.push_back(e); 00149 } 00150 00151 00152 // Taus 00153 /// @todo This could benefit from a tau finder projection 00154 Particles tau_candidates; 00155 foreach (const Particle& tau, applyProjection<UnstableFinalState>(event, "UFS").particlesByPt()) { 00156 // Only pick taus out of all unstable particles 00157 if (tau.abspid() != PID::TAU) continue; 00158 00159 // Check that tau has decayed into daughter particles 00160 /// @todo Huh? Unstable taus with no decay vtx? Can use Particle.isStable()? But why in this situation? 00161 if (tau.genParticle()->end_vertex() == 0) continue; 00162 00163 // Calculate visible tau pT from pT of tau neutrino in tau decay for pT and |eta| cuts 00164 FourMomentum daughter_tau_neutrino_momentum = get_tau_neutrino_mom(tau); 00165 Particle tau_vis = tau; 00166 tau_vis.setMomentum(tau.momentum()-daughter_tau_neutrino_momentum); 00167 // keep only taus in certain eta region and above 15 GeV of visible tau pT 00168 if ( tau_vis.pT() <= 15.0*GeV || tau_vis.abseta() > 2.5) continue; 00169 00170 // Get prong number (number of tracks) in tau decay and check if tau decays leptonically 00171 unsigned int nprong = 0; 00172 bool lep_decaying_tau = false; 00173 get_prong_number(tau.genParticle(), nprong, lep_decaying_tau); 00174 00175 // Apply reconstruction efficiency 00176 int tau_id = 15; 00177 if (nprong == 1) tau_id = 15; 00178 else if (nprong == 3) tau_id = 16; 00179 00180 // Get fiducial lepton efficiency simulate reco efficiency 00181 const double eff = (_use_fiducial_lepton_efficiency) ? apply_reco_eff(tau_id, tau_vis) : 1.0; 00182 const bool keep_tau = rand()/static_cast<double>(RAND_MAX) <= eff; 00183 00184 // Keep tau if nprong = 1, it decays hadronically, and it's reconstructed by the detector 00185 if ( !lep_decaying_tau && nprong == 1 && keep_tau) tau_candidates.push_back(tau_vis); 00186 } 00187 00188 00189 // Jets (all anti-kt R=0.4 jets with pT > 25 GeV and eta < 4.9) 00190 Jets jet_candidates; 00191 foreach (const Jet& jet, applyProjection<FastJets>(event, "AntiKtJets04").jetsByPt(25*GeV)) { 00192 if (jet.abseta() < 4.9) jet_candidates.push_back(jet); 00193 } 00194 00195 00196 // ETmiss 00197 Particles vfs_particles = applyProjection<VisibleFinalState>(event, "VFS").particles(); 00198 FourMomentum pTmiss; 00199 foreach (const Particle& p, vfs_particles) pTmiss -= p.momentum(); 00200 double eTmiss = pTmiss.pT()/GeV; 00201 00202 00203 //------------------ 00204 // Overlap removal 00205 00206 // electron - electron 00207 Particles electron_candidates_2; 00208 for (size_t ie = 0; ie < electron_candidates.size(); ++ie) { 00209 const Particle & e = electron_candidates[ie]; 00210 bool away = true; 00211 // If electron pair within dR < 0.1: remove electron with lower pT 00212 for (size_t ie2=0; ie2 < electron_candidates_2.size(); ++ie2) { 00213 if ( deltaR( e.momentum(), electron_candidates_2[ie2].momentum()) < 0.1 ) { 00214 away = false; 00215 break; 00216 } 00217 } 00218 // If isolated keep it 00219 if ( away ) 00220 electron_candidates_2.push_back( e ); 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 // jet - tau 00234 if (away) { 00235 // If jet within dR < 0.2 of tau: remove jet 00236 foreach (const Particle& tau, tau_candidates) { 00237 if (deltaR(tau.momentum(), jet.momentum()) < 0.2) { 00238 away = false; 00239 break; 00240 } 00241 } 00242 } 00243 // If isolated keep it 00244 if ( away ) 00245 recon_jets.push_back( jet ); 00246 } 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 += e; 00274 recon_leptons += e; 00275 } 00276 } 00277 00278 00279 // tau - electron 00280 Particles recon_tau; 00281 foreach ( const Particle& tau, tau_candidates ) { 00282 bool away = true; 00283 // If tau within dR < 0.2 of an electron: remove tau 00284 foreach ( const Particle& e, recon_e ) { 00285 if (deltaR( tau.momentum(), e.momentum()) < 0.2) { 00286 away = false; 00287 break; 00288 } 00289 } 00290 // tau - muon 00291 // If tau within dR < 0.2 of a muon: remove tau 00292 if (away) { 00293 foreach (const Particle& mu, muon_candidates) { 00294 if (deltaR(tau.momentum(), mu.momentum()) < 0.2) { 00295 away = false; 00296 break; 00297 } 00298 } 00299 } 00300 // If isolated keep it 00301 if (away) recon_tau.push_back( tau ); 00302 } 00303 00304 // Muon - jet isolation 00305 Particles recon_mu, trigger_mu; 00306 // If muon within dR < 0.4 of a jet, remove muon 00307 foreach (const Particle& mu, muon_candidates) { 00308 bool away = true; 00309 foreach (const Jet& jet, recon_jets) { 00310 if ( deltaR( mu.momentum(), jet.momentum()) < 0.4 ) { 00311 away = false; 00312 break; 00313 } 00314 } 00315 if (away) { 00316 recon_mu.push_back( mu ); 00317 recon_leptons.push_back( mu ); 00318 if (mu.abseta() < 2.4) trigger_mu.push_back( mu ); 00319 } 00320 } 00321 00322 // End overlap removal 00323 //------------------ 00324 00325 00326 // Jet cleaning 00327 if (rand()/static_cast<double>(RAND_MAX) <= 0.42) { 00328 foreach (const Jet& jet, recon_jets) { 00329 const double eta = jet.rapidity(); 00330 const double phi = jet.azimuthalAngle(MINUSPI_PLUSPI); 00331 if (jet.pT() > 25*GeV && inRange(eta, -0.1, 1.5) && inRange(phi, -0.9, -0.5)) vetoEvent; 00332 } 00333 } 00334 00335 00336 // Post-isolation event cuts 00337 // Require at least 3 charged tracks in event 00338 if (charged_tracks.size() < 3) vetoEvent; 00339 00340 // And at least one e/mu passing trigger 00341 if (!( !recon_e .empty() && recon_e[0] .pT() > 25*GeV) && 00342 !( !trigger_mu.empty() && trigger_mu[0].pT() > 25*GeV) ) { 00343 MSG_DEBUG("Hardest lepton fails trigger"); 00344 vetoEvent; 00345 } 00346 00347 // And only accept events with at least 2 electrons and muons and at least 3 leptons in total 00348 if (recon_mu.size() + recon_e.size() + recon_tau.size() < 3 || recon_leptons.size() < 2) vetoEvent; 00349 00350 // Now it's worth getting the event weight 00351 const double weight = event.weight(); 00352 00353 // Sort leptons by decreasing pT 00354 sortByPt(recon_leptons); 00355 sortByPt(recon_tau); 00356 00357 // Calculate HTlep, fill lepton pT histograms & store chosen combination of 3 leptons 00358 double HTlep = 0.; 00359 Particles chosen_leptons; 00360 if ( recon_leptons.size() > 2 ) { 00361 _h_pt_1_3l->fill(recon_leptons[0].perp()/GeV, weight); 00362 _h_pt_2_3l->fill(recon_leptons[1].perp()/GeV, weight); 00363 _h_pt_3_3l->fill(recon_leptons[2].perp()/GeV, weight); 00364 HTlep = (recon_leptons[0].pT() + recon_leptons[1].pT() + recon_leptons[2].pT())/GeV; 00365 chosen_leptons.push_back( recon_leptons[0] ); 00366 chosen_leptons.push_back( recon_leptons[1] ); 00367 chosen_leptons.push_back( recon_leptons[2] ); 00368 } 00369 else { 00370 _h_pt_1_2ltau->fill(recon_leptons[0].perp()/GeV, weight); 00371 _h_pt_2_2ltau->fill(recon_leptons[1].perp()/GeV, weight); 00372 _h_pt_3_2ltau->fill(recon_tau[0].perp()/GeV, weight); 00373 HTlep = (recon_leptons[0].pT() + recon_leptons[1].pT() + recon_tau[0].pT())/GeV ; 00374 chosen_leptons.push_back( recon_leptons[0] ); 00375 chosen_leptons.push_back( recon_leptons[1] ); 00376 chosen_leptons.push_back( recon_tau[0] ); 00377 } 00378 00379 // Number of prompt e/mu and had taus 00380 _h_e_n ->fill(recon_e.size() , weight); 00381 _h_mu_n ->fill(recon_mu.size() , weight); 00382 _h_tau_n->fill(recon_tau.size(), weight); 00383 00384 // Calculate HTjets 00385 double HTjets = 0.; 00386 foreach ( const Jet & jet, recon_jets ) 00387 HTjets += jet.perp()/GeV; 00388 00389 // Calculate meff 00390 double meff = eTmiss + HTjets; 00391 Particles all_leptons; 00392 foreach ( const Particle & e , recon_e ) { 00393 meff += e.perp()/GeV; 00394 all_leptons.push_back( e ); 00395 } 00396 foreach ( const Particle & mu, recon_mu ) { 00397 meff += mu.perp()/GeV; 00398 all_leptons.push_back( mu ); 00399 } 00400 foreach ( const Particle & tau, recon_tau ) { 00401 meff += tau.perp()/GeV; 00402 all_leptons.push_back( tau ); 00403 } 00404 00405 // Fill histogram of kinematic variables 00406 _h_HTlep_all ->fill(HTlep , weight); 00407 _h_HTjets_all->fill(HTjets, weight); 00408 _h_MET_all ->fill(eTmiss, weight); 00409 _h_Meff_all ->fill(meff , weight); 00410 00411 // Determine signal region (3l/2ltau, onZ/offZ) 00412 string basic_signal_region; 00413 if ( recon_mu.size() + recon_e.size() > 2 ) 00414 basic_signal_region += "3l_"; 00415 else if ( (recon_mu.size() + recon_e.size() == 2) && (recon_tau.size() > 0)) 00416 basic_signal_region += "2ltau_"; 00417 // Is there an OSSF pair or a three lepton combination with an invariant mass close to the Z mass 00418 int onZ = isonZ(chosen_leptons); 00419 if (onZ == 1) basic_signal_region += "onZ"; 00420 else if (onZ == 0) basic_signal_region += "offZ"; 00421 // Check in which signal regions this event falls and adjust event counters 00422 fillEventCountsPerSR(basic_signal_region, onZ, HTlep, eTmiss, HTjets, meff, weight); 00423 } 00424 00425 00426 /// Normalise histograms etc., after the run 00427 void finalize() { 00428 00429 // Normalize to an integrated luminosity of 1 fb-1 00430 double norm = crossSection()/femtobarn/sumOfWeights(); 00431 string best_signal_region = ""; 00432 double ratio_best_SR = 0.; 00433 00434 // Loop over all signal regions and find signal region with best sensitivity (ratio signal events/visible cross-section) 00435 for (size_t i = 0; i < _signal_regions.size(); i++) { 00436 double signal_events = _eventCountsPerSR[_signal_regions[i]] * norm; 00437 // Use expected upper limits to find best signal region 00438 double UL95 = getUpperLimit(_signal_regions[i], false); 00439 double ratio = signal_events / UL95; 00440 if (ratio > ratio_best_SR) { 00441 best_signal_region = _signal_regions[i]; 00442 ratio_best_SR = ratio; 00443 } 00444 } 00445 00446 double signal_events_best_SR = _eventCountsPerSR[best_signal_region] * norm; 00447 double exp_UL_best_SR = getUpperLimit(best_signal_region, false); 00448 double obs_UL_best_SR = getUpperLimit(best_signal_region, true); 00449 00450 // Print out result 00451 cout << "----------------------------------------------------------------------------------------" << endl; 00452 cout << "Best signal region: " << best_signal_region << endl; 00453 cout << "Normalized number of signal events in this best signal region (per fb-1): " << signal_events_best_SR << endl; 00454 cout << "Efficiency*Acceptance: " << _eventCountsPerSR[best_signal_region]/sumOfWeights() << endl; 00455 cout << "Cross-section [fb]: " << crossSection()/femtobarn << endl; 00456 cout << "Expected visible cross-section (per fb-1): " << exp_UL_best_SR << endl; 00457 cout << "Ratio (signal events / expected visible cross-section): " << ratio_best_SR << endl; 00458 cout << "Observed visible cross-section (per fb-1): " << obs_UL_best_SR << endl; 00459 cout << "Ratio (signal events / observed visible cross-section): " << signal_events_best_SR/obs_UL_best_SR << endl; 00460 cout << "----------------------------------------------------------------------------------------" << endl; 00461 00462 cout << "Using the EXPECTED limits (visible cross-section) of the analysis: " << endl; 00463 if (signal_events_best_SR > exp_UL_best_SR) { 00464 cout << "Since the number of signal events > the visible cross-section, this model/grid point is EXCLUDED with 95% CL." << endl; 00465 _h_excluded->fill(1); 00466 } 00467 else { 00468 cout << "Since the number of signal events < the visible cross-section, this model/grid point is NOT EXCLUDED." << endl; 00469 _h_excluded->fill(0); 00470 } 00471 cout << "----------------------------------------------------------------------------------------" << endl; 00472 00473 cout << "Using the OBSERVED limits (visible cross-section) of the analysis: " << endl; 00474 if (signal_events_best_SR > obs_UL_best_SR) { 00475 cout << "Since the number of signal events > the visible cross-section, this model/grid point is EXCLUDED with 95% CL." << endl; 00476 _h_excluded->fill(1); 00477 } 00478 else { 00479 cout << "Since the number of signal events < the visible cross-section, this model/grid point is NOT EXCLUDED." << endl; 00480 _h_excluded->fill(0); 00481 } 00482 cout << "----------------------------------------------------------------------------------------" << endl; 00483 00484 00485 // Normalize to cross section 00486 if (norm != 0) { 00487 scale(_h_HTlep_all, norm); 00488 scale(_h_HTjets_all, norm); 00489 scale(_h_MET_all, norm); 00490 scale(_h_Meff_all, norm); 00491 00492 scale(_h_pt_1_3l, norm); 00493 scale(_h_pt_2_3l, norm); 00494 scale(_h_pt_3_3l, norm); 00495 scale(_h_pt_1_2ltau, norm); 00496 scale(_h_pt_2_2ltau, norm); 00497 scale(_h_pt_3_2ltau, norm); 00498 00499 scale(_h_e_n, norm); 00500 scale(_h_mu_n, norm); 00501 scale(_h_tau_n, norm); 00502 00503 scale(_h_excluded, signal_events_best_SR); 00504 } 00505 } 00506 00507 00508 /// Helper functions 00509 //@{ 00510 00511 /// Function giving a list of all signal regions 00512 vector<string> getSignalRegions() { 00513 00514 // List of basic signal regions 00515 vector<string> basic_signal_regions; 00516 basic_signal_regions.push_back("3l_offZ"); 00517 basic_signal_regions.push_back("3l_onZ"); 00518 basic_signal_regions.push_back("2ltau_offZ"); 00519 basic_signal_regions.push_back("2ltau_onZ"); 00520 00521 // List of kinematic variables 00522 vector<string> kinematic_variables; 00523 kinematic_variables.push_back("HTlep"); 00524 kinematic_variables.push_back("METStrong"); 00525 kinematic_variables.push_back("METWeak"); 00526 kinematic_variables.push_back("Meff"); 00527 kinematic_variables.push_back("MeffStrong"); 00528 00529 vector<string> signal_regions; 00530 // Loop over all kinematic variables and basic signal regions 00531 for (size_t i0 = 0; i0 < kinematic_variables.size(); i0++) { 00532 for (size_t i1 = 0; i1 < basic_signal_regions.size(); i1++) { 00533 // Is signal region onZ? 00534 int onZ = (basic_signal_regions[i1].find("onZ") != string::npos) ? 1 : 0; 00535 // Get cut values for this kinematic variable 00536 vector<int> cut_values = getCutsPerSignalRegion(kinematic_variables[i0], onZ); 00537 // Loop over all cut values 00538 for (size_t i2 = 0; i2 < cut_values.size(); i2++) { 00539 // push signal region into vector 00540 signal_regions.push_back( (kinematic_variables[i0] + "_" + basic_signal_regions[i1] + "_cut_" + toString(i2)) ); 00541 } 00542 } 00543 } 00544 return signal_regions; 00545 } 00546 00547 00548 /// Function giving all cut vales per kinematic variable (taking onZ for MET into account) 00549 vector<int> getCutsPerSignalRegion(const string& signal_region, int onZ=0) { 00550 vector<int> cutValues; 00551 00552 // Cut values for HTlep 00553 if (signal_region.compare("HTlep") == 0) { 00554 cutValues.push_back(0); 00555 cutValues.push_back(100); 00556 cutValues.push_back(150); 00557 cutValues.push_back(200); 00558 cutValues.push_back(300); 00559 } 00560 // Cut values for METStrong (HTjets > 100 GeV) and METWeak (HTjets < 100 GeV) 00561 else if (signal_region.compare("METStrong") == 0 || signal_region.compare("METWeak") == 0) { 00562 if (onZ == 0) cutValues.push_back(0); 00563 else if (onZ == 1) cutValues.push_back(20); 00564 cutValues.push_back(50); 00565 cutValues.push_back(75); 00566 } 00567 // Cut values for Meff and MeffStrong (MET > 75 GeV) 00568 if (signal_region.compare("Meff") == 0 || signal_region.compare("MeffStrong") == 0) { 00569 cutValues.push_back(0); 00570 cutValues.push_back(150); 00571 cutValues.push_back(300); 00572 cutValues.push_back(500); 00573 } 00574 00575 return cutValues; 00576 } 00577 00578 00579 /// function fills map EventCountsPerSR by looping over all signal regions 00580 /// and looking if the event falls into this signal region 00581 void fillEventCountsPerSR(const string& basic_signal_region, int onZ, 00582 double HTlep, double eTmiss, 00583 double HTjets, double meff, 00584 double weight) { 00585 00586 // Get cut values for HTlep, loop over them and add event if cut is passed 00587 vector<int> cut_values = getCutsPerSignalRegion("HTlep", onZ); 00588 for (size_t i = 0; i < cut_values.size(); i++) { 00589 if (HTlep > cut_values[i]) 00590 _eventCountsPerSR[("HTlep_" + basic_signal_region + "_cut_" + toString(cut_values[i]))] += weight; 00591 } 00592 00593 // Get cut values for METStrong, loop over them and add event if cut is passed 00594 cut_values = getCutsPerSignalRegion("METStrong", onZ); 00595 for (size_t i = 0; i < cut_values.size(); i++) { 00596 if (eTmiss > cut_values[i] && HTjets > 100.) 00597 _eventCountsPerSR[("METStrong_" + basic_signal_region + "_cut_" + toString(cut_values[i]))] += weight; 00598 } 00599 00600 // Get cut values for METWeak, loop over them and add event if cut is passed 00601 cut_values = getCutsPerSignalRegion("METWeak", onZ); 00602 for (size_t i = 0; i < cut_values.size(); i++) { 00603 if (eTmiss > cut_values[i] && HTjets <= 100.) 00604 _eventCountsPerSR[("METWeak_" + basic_signal_region + "_cut_" + toString(cut_values[i]))] += weight; 00605 } 00606 00607 // Get cut values for Meff, loop over them and add event if cut is passed 00608 cut_values = getCutsPerSignalRegion("Meff", onZ); 00609 for (size_t i = 0; i < cut_values.size(); i++) { 00610 if (meff > cut_values[i]) 00611 _eventCountsPerSR[("Meff_" + basic_signal_region + "_cut_" + toString(cut_values[i]))] += weight; 00612 } 00613 00614 // Get cut values for MeffStrong, loop over them and add event if cut is passed 00615 cut_values = getCutsPerSignalRegion("MeffStrong", onZ); 00616 for (size_t i = 0; i < cut_values.size(); i++) { 00617 if (meff > cut_values[i] && eTmiss > 75.) 00618 _eventCountsPerSR[("MeffStrong_" + basic_signal_region + "_cut_" + toString(cut_values[i]))] += weight; 00619 } 00620 } 00621 00622 00623 /// Function returning 4-vector of daughter-particle if it is a tau neutrino 00624 /// @todo Move to TauFinder and make less HepMC-ish 00625 FourMomentum get_tau_neutrino_mom(const Particle& p) { 00626 assert(p.abspid() == PID::TAU); 00627 const GenVertex* dv = p.genParticle()->end_vertex(); 00628 assert(dv != NULL); 00629 for (GenVertex::particles_out_const_iterator pp = dv->particles_out_const_begin(); pp != dv->particles_out_const_end(); ++pp) { 00630 if (abs((*pp)->pdg_id()) == PID::NU_TAU) return FourMomentum((*pp)->momentum()); 00631 } 00632 return FourMomentum(); 00633 } 00634 00635 00636 /// Function calculating the prong number of taus 00637 /// @todo Move to TauFinder and make less HepMC-ish 00638 void get_prong_number(const GenParticle* p, unsigned int& nprong, bool& lep_decaying_tau) { 00639 assert(p != NULL); 00640 //const int tau_barcode = p->barcode(); 00641 const GenVertex* dv = p->end_vertex(); 00642 assert(dv != NULL); 00643 for (GenVertex::particles_out_const_iterator pp = dv->particles_out_const_begin(); pp != dv->particles_out_const_end(); ++pp) { 00644 // If they have status 1 and are charged they will produce a track and the prong number is +1 00645 if ((*pp)->status() == 1 ) { 00646 const int id = (*pp)->pdg_id(); 00647 if (Rivet::PID::charge(id) != 0 ) ++nprong; 00648 // Check if tau decays leptonically 00649 // @todo Can a tau decay include a tau in its decay daughters?! 00650 if ((abs(id) == PID::ELECTRON || abs(id) == PID::MUON || abs(id) == PID::TAU) && abs(p->pdg_id()) == PID::TAU) lep_decaying_tau = true; 00651 } 00652 // If the status of the daughter particle is 2 it is unstable and the further decays are checked 00653 else if ((*pp)->status() == 2 ) { 00654 get_prong_number(*pp, nprong, lep_decaying_tau); 00655 } 00656 } 00657 } 00658 00659 00660 /// Function giving fiducial lepton efficiency 00661 double apply_reco_eff(int flavor, const Particle& p) { 00662 float pt = p.pT()/GeV; 00663 float eta = p.eta(); 00664 00665 double eff = 0.; 00666 //double err = 0.; 00667 00668 if (flavor == 11) { // weight prompt electron -- now including data/MC ID SF in eff. 00669 //float rho = 0.820; 00670 float p0 = 7.34; float p1 = 0.8977; 00671 //float ep0= 0.5 ; float ep1= 0.0087; 00672 eff = p1 - p0/pt; 00673 00674 //double err0 = ep0/pt; // d(eff)/dp0 00675 //double err1 = ep1; // d(eff)/dp1 00676 //err = sqrt(err0*err0 + err1*err1 - 2*rho*err0*err1); 00677 00678 double avgrate = 0.6867; 00679 float wz_ele_eta[] = {0.588717,0.603674,0.666135,0.747493,0.762202,0.675051,0.751606,0.745569,0.665333,0.610432,0.592693,}; 00680 //float ewz_ele_eta[] ={0.00292902,0.002476,0.00241209,0.00182319,0.00194339,0.00299785,0.00197339,0.00182004,0.00241793,0.00245997,0.00290394,}; 00681 int ibin = 3; 00682 00683 if (eta >= -2.5 && eta < -2.0) ibin = 0; 00684 if (eta >= -2.0 && eta < -1.5) ibin = 1; 00685 if (eta >= -1.5 && eta < -1.0) ibin = 2; 00686 if (eta >= -1.0 && eta < -0.5) ibin = 3; 00687 if (eta >= -0.5 && eta < -0.1) ibin = 4; 00688 if (eta >= -0.1 && eta < 0.1) ibin = 5; 00689 if (eta >= 0.1 && eta < 0.5) ibin = 6; 00690 if (eta >= 0.5 && eta < 1.0) ibin = 7; 00691 if (eta >= 1.0 && eta < 1.5) ibin = 8; 00692 if (eta >= 1.5 && eta < 2.0) ibin = 9; 00693 if (eta >= 2.0 && eta < 2.5) ibin = 10; 00694 00695 double eff_eta = wz_ele_eta[ibin]; 00696 //double err_eta = ewz_ele_eta[ibin]; 00697 00698 eff = (eff*eff_eta)/avgrate; 00699 } 00700 00701 if (flavor == 12) { // weight electron from tau 00702 //float rho = 0.884; 00703 float p0 = 6.799; float p1 = 0.842; 00704 //float ep0= 0.664; float ep1= 0.016; 00705 eff = p1 - p0/pt; 00706 00707 //double err0 = ep0/pt; // d(eff)/dp0 00708 //double err1 = ep1; // d(eff)/dp1 00709 //err = sqrt(err0*err0 + err1*err1 - 2*rho*err0*err1); 00710 00711 double avgrate = 0.5319; 00712 float wz_elet_eta[] = {0.468945,0.465953,0.489545,0.58709,0.59669,0.515829,0.59284,0.575828,0.498181,0.463536,0.481738,}; 00713 //float ewz_elet_eta[] ={0.00933795,0.00780868,0.00792679,0.00642083,0.00692652,0.0101568,0.00698452,0.00643524,0.0080002,0.00776238,0.0094699,}; 00714 int ibin = 3; 00715 00716 if (eta >= -2.5 && eta < -2.0) ibin = 0; 00717 if (eta >= -2.0 && eta < -1.5) ibin = 1; 00718 if (eta >= -1.5 && eta < -1.0) ibin = 2; 00719 if (eta >= -1.0 && eta < -0.5) ibin = 3; 00720 if (eta >= -0.5 && eta < -0.1) ibin = 4; 00721 if (eta >= -0.1 && eta < 0.1) ibin = 5; 00722 if (eta >= 0.1 && eta < 0.5) ibin = 6; 00723 if (eta >= 0.5 && eta < 1.0) ibin = 7; 00724 if (eta >= 1.0 && eta < 1.5) ibin = 8; 00725 if (eta >= 1.5 && eta < 2.0) ibin = 9; 00726 if (eta >= 2.0 && eta < 2.5) ibin = 10; 00727 00728 double eff_eta = wz_elet_eta[ibin]; 00729 //double err_eta = ewz_elet_eta[ibin]; 00730 00731 eff = (eff*eff_eta)/avgrate; 00732 00733 } 00734 00735 if (flavor == 13) {// weight prompt muon 00736 00737 //if eta>0.1 00738 float p0 = -18.21; float p1 = 14.83; float p2 = 0.9312; 00739 //float ep0= 5.06; float ep1= 1.9; float ep2=0.00069; 00740 00741 if ( fabs(eta) < 0.1) { 00742 p0 = 7.459; p1 = 2.615; p2 = 0.5138; 00743 //ep0 = 10.4; ep1 = 4.934; ep2 = 0.0034; 00744 } 00745 00746 double arg = ( pt-p0 )/( 2.*p1 ) ; 00747 eff = 0.5 * p2 * (1.+erf(arg)); 00748 //err = 0.1*eff; 00749 } 00750 00751 if (flavor == 14) {// weight muon from tau 00752 00753 if (fabs(eta) < 0.1) { 00754 float p0 = -1.756; float p1 = 12.38; float p2 = 0.4441; 00755 //float ep0= 10.39; float ep1= 7.9; float ep2=0.022; 00756 double arg = ( pt-p0 )/( 2.*p1 ) ; 00757 eff = 0.5 * p2 * (1.+erf(arg)); 00758 //err = 0.1*eff; 00759 } 00760 else { 00761 float p0 = 2.102; float p1 = 0.8293; 00762 //float ep0= 0.271; float ep1= 0.0083; 00763 eff = p1 - p0/pt; 00764 //double err0 = ep0/pt; // d(eff)/dp0 00765 //double err1 = ep1; // d(eff)/dp1 00766 //err = sqrt(err0*err0 + err1*err1 - 2*rho*err0*err1); 00767 } 00768 } 00769 00770 if (flavor == 15) {// weight hadronic tau 1p 00771 00772 float wz_tau1p[] = {0.0249278,0.146978,0.225049,0.229212,0.21519,0.206152,0.201559,0.197917,0.209249,0.228336,0.193548,}; 00773 //float ewz_tau1p[] ={0.00178577,0.00425252,0.00535052,0.00592126,0.00484684,0.00612941,0.00792099,0.0083006,0.0138307,0.015568,0.0501751,}; 00774 int ibin = 0; 00775 if (pt > 15) ibin = 1; 00776 if (pt > 20) ibin = 2; 00777 if (pt > 25) ibin = 3; 00778 if (pt > 30) ibin = 4; 00779 if (pt > 40) ibin = 5; 00780 if (pt > 50) ibin = 6; 00781 if (pt > 60) ibin = 7; 00782 if (pt > 80) ibin = 8; 00783 if (pt > 100) ibin = 9; 00784 if (pt > 200) ibin = 10; 00785 00786 eff = wz_tau1p[ibin]; 00787 //err = ewz_tau1p[ibin]; 00788 00789 00790 double avgrate = 0.1718; 00791 float wz_tau1p_eta[] = {0.162132,0.176393,0.139619,0.178813,0.185144,0.210027,0.203937,0.178688,0.137034,0.164216,0.163713,}; 00792 //float ewz_tau1p_eta[] ={0.00706705,0.00617989,0.00506798,0.00525172,0.00581865,0.00865675,0.00599245,0.00529877,0.00506368,0.00617025,0.00726219,}; 00793 00794 ibin = 3; 00795 if (eta >= -2.5 && eta < -2.0) ibin = 0; 00796 if (eta >= -2.0 && eta < -1.5) ibin = 1; 00797 if (eta >= -1.5 && eta < -1.0) ibin = 2; 00798 if (eta >= -1.0 && eta < -0.5) ibin = 3; 00799 if (eta >= -0.5 && eta < -0.1) ibin = 4; 00800 if (eta >= -0.1 && eta < 0.1) ibin = 5; 00801 if (eta >= 0.1 && eta < 0.5) ibin = 6; 00802 if (eta >= 0.5 && eta < 1.0) ibin = 7; 00803 if (eta >= 1.0 && eta < 1.5) ibin = 8; 00804 if (eta >= 1.5 && eta < 2.0) ibin = 9; 00805 if (eta >= 2.0 && eta < 2.5) ibin = 10; 00806 00807 double eff_eta = wz_tau1p_eta[ibin]; 00808 //double err_eta = ewz_tau1p_eta[ibin]; 00809 00810 eff = (eff*eff_eta)/avgrate; 00811 } 00812 00813 if (flavor == 16) { //weight hadronic tau 3p 00814 00815 float wz_tau3p[] = {0.000587199,0.00247181,0.0013031,0.00280112,}; 00816 //float ewz_tau3p[] ={0.000415091,0.000617187,0.000582385,0.00197792,}; 00817 00818 int ibin = 0; 00819 if (pt > 15) ibin = 1; 00820 if (pt > 20) ibin = 2; 00821 if (pt > 40) ibin = 3; 00822 if (pt > 80) ibin = 4; 00823 00824 eff = wz_tau3p[ibin]; 00825 //err = ewz_tau3p[ibin]; 00826 } 00827 00828 return eff; 00829 } 00830 00831 00832 /// Function giving observed upper limit (visible cross-section) 00833 double getUpperLimit(const string& signal_region, bool observed) { 00834 map<string,double> upperLimitsObserved; 00835 upperLimitsObserved["HTlep_3l_offZ_cut_0"] = 11.; 00836 upperLimitsObserved["HTlep_3l_offZ_cut_100"] = 8.7; 00837 upperLimitsObserved["HTlep_3l_offZ_cut_150"] = 4.0; 00838 upperLimitsObserved["HTlep_3l_offZ_cut_200"] = 4.4; 00839 upperLimitsObserved["HTlep_3l_offZ_cut_300"] = 1.6; 00840 upperLimitsObserved["HTlep_2ltau_offZ_cut_0"] = 25.; 00841 upperLimitsObserved["HTlep_2ltau_offZ_cut_100"] = 14.; 00842 upperLimitsObserved["HTlep_2ltau_offZ_cut_150"] = 6.1; 00843 upperLimitsObserved["HTlep_2ltau_offZ_cut_200"] = 3.3; 00844 upperLimitsObserved["HTlep_2ltau_offZ_cut_300"] = 1.2; 00845 upperLimitsObserved["HTlep_3l_onZ_cut_0"] = 48.; 00846 upperLimitsObserved["HTlep_3l_onZ_cut_100"] = 38.; 00847 upperLimitsObserved["HTlep_3l_onZ_cut_150"] = 14.; 00848 upperLimitsObserved["HTlep_3l_onZ_cut_200"] = 7.2; 00849 upperLimitsObserved["HTlep_3l_onZ_cut_300"] = 4.5; 00850 upperLimitsObserved["HTlep_2ltau_onZ_cut_0"] = 85.; 00851 upperLimitsObserved["HTlep_2ltau_onZ_cut_100"] = 53.; 00852 upperLimitsObserved["HTlep_2ltau_onZ_cut_150"] = 11.0; 00853 upperLimitsObserved["HTlep_2ltau_onZ_cut_200"] = 5.2; 00854 upperLimitsObserved["HTlep_2ltau_onZ_cut_300"] = 3.0; 00855 upperLimitsObserved["METStrong_3l_offZ_cut_0"] = 2.6; 00856 upperLimitsObserved["METStrong_3l_offZ_cut_50"] = 2.1; 00857 upperLimitsObserved["METStrong_3l_offZ_cut_75"] = 2.1; 00858 upperLimitsObserved["METStrong_2ltau_offZ_cut_0"] = 4.2; 00859 upperLimitsObserved["METStrong_2ltau_offZ_cut_50"] = 3.1; 00860 upperLimitsObserved["METStrong_2ltau_offZ_cut_75"] = 2.6; 00861 upperLimitsObserved["METStrong_3l_onZ_cut_20"] = 11.0; 00862 upperLimitsObserved["METStrong_3l_onZ_cut_50"] = 6.4; 00863 upperLimitsObserved["METStrong_3l_onZ_cut_75"] = 5.1; 00864 upperLimitsObserved["METStrong_2ltau_onZ_cut_20"] = 5.9; 00865 upperLimitsObserved["METStrong_2ltau_onZ_cut_50"] = 3.4; 00866 upperLimitsObserved["METStrong_2ltau_onZ_cut_75"] = 1.2; 00867 upperLimitsObserved["METWeak_3l_offZ_cut_0"] = 11.; 00868 upperLimitsObserved["METWeak_3l_offZ_cut_50"] = 5.3; 00869 upperLimitsObserved["METWeak_3l_offZ_cut_75"] = 3.1; 00870 upperLimitsObserved["METWeak_2ltau_offZ_cut_0"] = 23.; 00871 upperLimitsObserved["METWeak_2ltau_offZ_cut_50"] = 4.3; 00872 upperLimitsObserved["METWeak_2ltau_offZ_cut_75"] = 3.1; 00873 upperLimitsObserved["METWeak_3l_onZ_cut_20"] = 41.; 00874 upperLimitsObserved["METWeak_3l_onZ_cut_50"] = 16.; 00875 upperLimitsObserved["METWeak_3l_onZ_cut_75"] = 8.0; 00876 upperLimitsObserved["METWeak_2ltau_onZ_cut_20"] = 80.; 00877 upperLimitsObserved["METWeak_2ltau_onZ_cut_50"] = 4.4; 00878 upperLimitsObserved["METWeak_2ltau_onZ_cut_75"] = 1.8; 00879 upperLimitsObserved["Meff_3l_offZ_cut_0"] = 11.; 00880 upperLimitsObserved["Meff_3l_offZ_cut_150"] = 8.1; 00881 upperLimitsObserved["Meff_3l_offZ_cut_300"] = 3.1; 00882 upperLimitsObserved["Meff_3l_offZ_cut_500"] = 2.1; 00883 upperLimitsObserved["Meff_2ltau_offZ_cut_0"] = 25.; 00884 upperLimitsObserved["Meff_2ltau_offZ_cut_150"] = 12.; 00885 upperLimitsObserved["Meff_2ltau_offZ_cut_300"] = 3.9; 00886 upperLimitsObserved["Meff_2ltau_offZ_cut_500"] = 2.2; 00887 upperLimitsObserved["Meff_3l_onZ_cut_0"] = 48.; 00888 upperLimitsObserved["Meff_3l_onZ_cut_150"] = 37.; 00889 upperLimitsObserved["Meff_3l_onZ_cut_300"] = 11.; 00890 upperLimitsObserved["Meff_3l_onZ_cut_500"] = 4.8; 00891 upperLimitsObserved["Meff_2ltau_onZ_cut_0"] = 85.; 00892 upperLimitsObserved["Meff_2ltau_onZ_cut_150"] = 28.; 00893 upperLimitsObserved["Meff_2ltau_onZ_cut_300"] = 5.9; 00894 upperLimitsObserved["Meff_2ltau_onZ_cut_500"] = 1.9; 00895 upperLimitsObserved["MeffStrong_3l_offZ_cut_0"] = 3.8; 00896 upperLimitsObserved["MeffStrong_3l_offZ_cut_150"] = 3.8; 00897 upperLimitsObserved["MeffStrong_3l_offZ_cut_300"] = 2.8; 00898 upperLimitsObserved["MeffStrong_3l_offZ_cut_500"] = 2.1; 00899 upperLimitsObserved["MeffStrong_2ltau_offZ_cut_0"] = 3.9; 00900 upperLimitsObserved["MeffStrong_2ltau_offZ_cut_150"] = 4.0; 00901 upperLimitsObserved["MeffStrong_2ltau_offZ_cut_300"] = 2.9; 00902 upperLimitsObserved["MeffStrong_2ltau_offZ_cut_500"] = 1.5; 00903 upperLimitsObserved["MeffStrong_3l_onZ_cut_0"] = 10.0; 00904 upperLimitsObserved["MeffStrong_3l_onZ_cut_150"] = 10.0; 00905 upperLimitsObserved["MeffStrong_3l_onZ_cut_300"] = 6.8; 00906 upperLimitsObserved["MeffStrong_3l_onZ_cut_500"] = 3.9; 00907 upperLimitsObserved["MeffStrong_2ltau_onZ_cut_0"] = 1.6; 00908 upperLimitsObserved["MeffStrong_2ltau_onZ_cut_150"] = 1.4; 00909 upperLimitsObserved["MeffStrong_2ltau_onZ_cut_300"] = 1.5; 00910 upperLimitsObserved["MeffStrong_2ltau_onZ_cut_500"] = 0.9; 00911 00912 // Expected upper limits are also given but not used in this analysis 00913 map<string,double> upperLimitsExpected; 00914 upperLimitsExpected["HTlep_3l_offZ_cut_0"] = 11.; 00915 upperLimitsExpected["HTlep_3l_offZ_cut_100"] = 8.5; 00916 upperLimitsExpected["HTlep_3l_offZ_cut_150"] = 4.6; 00917 upperLimitsExpected["HTlep_3l_offZ_cut_200"] = 3.6; 00918 upperLimitsExpected["HTlep_3l_offZ_cut_300"] = 1.9; 00919 upperLimitsExpected["HTlep_2ltau_offZ_cut_0"] = 23.; 00920 upperLimitsExpected["HTlep_2ltau_offZ_cut_100"] = 14.; 00921 upperLimitsExpected["HTlep_2ltau_offZ_cut_150"] = 6.4; 00922 upperLimitsExpected["HTlep_2ltau_offZ_cut_200"] = 3.6; 00923 upperLimitsExpected["HTlep_2ltau_offZ_cut_300"] = 1.5; 00924 upperLimitsExpected["HTlep_3l_onZ_cut_0"] = 33.; 00925 upperLimitsExpected["HTlep_3l_onZ_cut_100"] = 25.; 00926 upperLimitsExpected["HTlep_3l_onZ_cut_150"] = 12.; 00927 upperLimitsExpected["HTlep_3l_onZ_cut_200"] = 6.5; 00928 upperLimitsExpected["HTlep_3l_onZ_cut_300"] = 3.1; 00929 upperLimitsExpected["HTlep_2ltau_onZ_cut_0"] = 94.; 00930 upperLimitsExpected["HTlep_2ltau_onZ_cut_100"] = 61.; 00931 upperLimitsExpected["HTlep_2ltau_onZ_cut_150"] = 9.9; 00932 upperLimitsExpected["HTlep_2ltau_onZ_cut_200"] = 4.5; 00933 upperLimitsExpected["HTlep_2ltau_onZ_cut_300"] = 1.9; 00934 upperLimitsExpected["METStrong_3l_offZ_cut_0"] = 3.1; 00935 upperLimitsExpected["METStrong_3l_offZ_cut_50"] = 2.4; 00936 upperLimitsExpected["METStrong_3l_offZ_cut_75"] = 2.3; 00937 upperLimitsExpected["METStrong_2ltau_offZ_cut_0"] = 4.8; 00938 upperLimitsExpected["METStrong_2ltau_offZ_cut_50"] = 3.3; 00939 upperLimitsExpected["METStrong_2ltau_offZ_cut_75"] = 2.1; 00940 upperLimitsExpected["METStrong_3l_onZ_cut_20"] = 8.7; 00941 upperLimitsExpected["METStrong_3l_onZ_cut_50"] = 4.9; 00942 upperLimitsExpected["METStrong_3l_onZ_cut_75"] = 3.8; 00943 upperLimitsExpected["METStrong_2ltau_onZ_cut_20"] = 7.3; 00944 upperLimitsExpected["METStrong_2ltau_onZ_cut_50"] = 2.8; 00945 upperLimitsExpected["METStrong_2ltau_onZ_cut_75"] = 1.5; 00946 upperLimitsExpected["METWeak_3l_offZ_cut_0"] = 10.; 00947 upperLimitsExpected["METWeak_3l_offZ_cut_50"] = 4.7; 00948 upperLimitsExpected["METWeak_3l_offZ_cut_75"] = 3.0; 00949 upperLimitsExpected["METWeak_2ltau_offZ_cut_0"] = 21.; 00950 upperLimitsExpected["METWeak_2ltau_offZ_cut_50"] = 4.0; 00951 upperLimitsExpected["METWeak_2ltau_offZ_cut_75"] = 2.6; 00952 upperLimitsExpected["METWeak_3l_onZ_cut_20"] = 30.; 00953 upperLimitsExpected["METWeak_3l_onZ_cut_50"] = 10.; 00954 upperLimitsExpected["METWeak_3l_onZ_cut_75"] = 5.4; 00955 upperLimitsExpected["METWeak_2ltau_onZ_cut_20"] = 88.; 00956 upperLimitsExpected["METWeak_2ltau_onZ_cut_50"] = 5.5; 00957 upperLimitsExpected["METWeak_2ltau_onZ_cut_75"] = 2.2; 00958 upperLimitsExpected["Meff_3l_offZ_cut_0"] = 11.; 00959 upperLimitsExpected["Meff_3l_offZ_cut_150"] = 8.8; 00960 upperLimitsExpected["Meff_3l_offZ_cut_300"] = 3.7; 00961 upperLimitsExpected["Meff_3l_offZ_cut_500"] = 2.1; 00962 upperLimitsExpected["Meff_2ltau_offZ_cut_0"] = 23.; 00963 upperLimitsExpected["Meff_2ltau_offZ_cut_150"] = 13.; 00964 upperLimitsExpected["Meff_2ltau_offZ_cut_300"] = 4.9; 00965 upperLimitsExpected["Meff_2ltau_offZ_cut_500"] = 2.4; 00966 upperLimitsExpected["Meff_3l_onZ_cut_0"] = 33.; 00967 upperLimitsExpected["Meff_3l_onZ_cut_150"] = 25.; 00968 upperLimitsExpected["Meff_3l_onZ_cut_300"] = 9.; 00969 upperLimitsExpected["Meff_3l_onZ_cut_500"] = 3.9; 00970 upperLimitsExpected["Meff_2ltau_onZ_cut_0"] = 94.; 00971 upperLimitsExpected["Meff_2ltau_onZ_cut_150"] = 35.; 00972 upperLimitsExpected["Meff_2ltau_onZ_cut_300"] = 6.8; 00973 upperLimitsExpected["Meff_2ltau_onZ_cut_500"] = 2.5; 00974 upperLimitsExpected["MeffStrong_3l_offZ_cut_0"] = 3.9; 00975 upperLimitsExpected["MeffStrong_3l_offZ_cut_150"] = 3.9; 00976 upperLimitsExpected["MeffStrong_3l_offZ_cut_300"] = 3.0; 00977 upperLimitsExpected["MeffStrong_3l_offZ_cut_500"] = 2.0; 00978 upperLimitsExpected["MeffStrong_2ltau_offZ_cut_0"] = 3.8; 00979 upperLimitsExpected["MeffStrong_2ltau_offZ_cut_150"] = 3.9; 00980 upperLimitsExpected["MeffStrong_2ltau_offZ_cut_300"] = 3.1; 00981 upperLimitsExpected["MeffStrong_2ltau_offZ_cut_500"] = 1.6; 00982 upperLimitsExpected["MeffStrong_3l_onZ_cut_0"] = 6.9; 00983 upperLimitsExpected["MeffStrong_3l_onZ_cut_150"] = 7.1; 00984 upperLimitsExpected["MeffStrong_3l_onZ_cut_300"] = 4.9; 00985 upperLimitsExpected["MeffStrong_3l_onZ_cut_500"] = 3.0; 00986 upperLimitsExpected["MeffStrong_2ltau_onZ_cut_0"] = 2.4; 00987 upperLimitsExpected["MeffStrong_2ltau_onZ_cut_150"] = 2.5; 00988 upperLimitsExpected["MeffStrong_2ltau_onZ_cut_300"] = 2.0; 00989 upperLimitsExpected["MeffStrong_2ltau_onZ_cut_500"] = 1.1; 00990 00991 if (observed) return upperLimitsObserved[signal_region]; 00992 else return upperLimitsExpected[signal_region]; 00993 } 00994 00995 00996 /// Function checking if there is an OSSF lepton pair or a combination of 3 leptons with an invariant mass close to the Z mass 00997 /// @todo Should the reference Z mass be 91.2? 00998 int isonZ (const Particles& particles) { 00999 int onZ = 0; 01000 double best_mass_2 = 999.; 01001 double best_mass_3 = 999.; 01002 01003 // Loop over all 2 particle combinations to find invariant mass of OSSF pair closest to Z mass 01004 foreach ( const Particle& p1, particles ) { 01005 foreach ( const Particle& p2, particles ) { 01006 double mass_difference_2_old = fabs(91.0 - best_mass_2); 01007 double mass_difference_2_new = fabs(91.0 - (p1.momentum() + p2.momentum()).mass()/GeV); 01008 01009 // If particle combination is OSSF pair calculate mass difference to Z mass 01010 if ( (p1.pid()*p2.pid() == -121 || p1.pid()*p2.pid() == -169) ) { 01011 01012 // Get invariant mass closest to Z mass 01013 if (mass_difference_2_new < mass_difference_2_old) 01014 best_mass_2 = (p1.momentum() + p2.momentum()).mass()/GeV; 01015 01016 // In case there is an OSSF pair take also 3rd lepton into account (e.g. from FSR and photon to electron conversion) 01017 foreach ( const Particle & p3 , particles ) { 01018 double mass_difference_3_old = fabs(91.0 - best_mass_3); 01019 double mass_difference_3_new = fabs(91.0 - (p1.momentum() + p2.momentum() + p3.momentum()).mass()/GeV); 01020 if (mass_difference_3_new < mass_difference_3_old) 01021 best_mass_3 = (p1.momentum() + p2.momentum() + p3.momentum()).mass()/GeV; 01022 } 01023 } 01024 } 01025 } 01026 01027 // Pick the minimum invariant mass of the best OSSF pair combination and the best 3 lepton combination 01028 // If this mass is in a 20 GeV window around the Z mass, the event is classified as onZ 01029 double best_mass = min(best_mass_2, best_mass_3); 01030 if (fabs(91.0 - best_mass) < 20) onZ = 1; 01031 return onZ; 01032 } 01033 01034 //@} 01035 01036 01037 private: 01038 01039 /// Histograms 01040 //@{ 01041 Histo1DPtr _h_HTlep_all, _h_HTjets_all, _h_MET_all, _h_Meff_all; 01042 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; 01043 Histo1DPtr _h_e_n, _h_mu_n, _h_tau_n; 01044 Histo1DPtr _h_excluded; 01045 //@} 01046 01047 /// Fiducial efficiencies to model the effects of the ATLAS detector 01048 bool _use_fiducial_lepton_efficiency; 01049 01050 /// List of signal regions and event counts per signal region 01051 vector<string> _signal_regions; 01052 map<string, double> _eventCountsPerSR; 01053 01054 }; 01055 01056 01057 01058 DECLARE_RIVET_PLUGIN(ATLAS_2012_I1204447); 01059 01060 } Generated on Wed Oct 7 2015 12:09:10 for The Rivet MC analysis system by ![]() |