Skip to content
Snippets Groups Projects
Commit 80efe6ba authored by Adam Edward Barton's avatar Adam Edward Barton
Browse files

Merge branch 'master-JetEtmiss-NTUP_PHYSVAL' into 'master'

Improve protection against missing four-vectors

See merge request atlas/athena!38175
parents e8d34e0e 46b4c455
No related branches found
No related tags found
No related merge requests found
// -*- c++ -*- // -*- c++ -*-
/* /*
Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
*/ */
#ifndef JETMONITORING_JETKINEMATICHISTOS_H #ifndef JETMONITORING_JETKINEMATICHISTOS_H
...@@ -52,7 +52,7 @@ protected: ...@@ -52,7 +52,7 @@ protected:
TProfile2D* m_averagePtEtaPhi; TProfile2D* m_averagePtEtaPhi;
TProfile2D* m_averageE_EtaPhi; TProfile2D* m_averageE_EtaPhi;
int m_jetScale; std::string m_jetScale;
bool m_doN; bool m_doN;
bool m_doM; bool m_doM;
......
...@@ -15,7 +15,7 @@ JetKinematicHistos::JetKinematicHistos(const std::string &t) : JetHistoBase(t) ...@@ -15,7 +15,7 @@ JetKinematicHistos::JetKinematicHistos(const std::string &t) : JetHistoBase(t)
,m_phi(0) ,m_phi(0)
,m_m(0) ,m_m(0)
,m_e(0) ,m_e(0)
, m_jetScale(xAOD::JetAssignedScaleMomentum) ,m_jetScale("JetAssignedScaleMomentum")
{ {
declareProperty("JetScale", m_jetScale); declareProperty("JetScale", m_jetScale);
declareProperty("PlotNJet", m_doN = false); declareProperty("PlotNJet", m_doN = false);
...@@ -34,12 +34,11 @@ int JetKinematicHistos::buildHistos(){ ...@@ -34,12 +34,11 @@ int JetKinematicHistos::buildHistos(){
// -------------- Modify histo names/titles -------- // -------------- Modify histo names/titles --------
// if we're plotting a defined scale, modify the Histo titles and names // if we're plotting a defined scale, modify the Histo titles and names
std::map<xAOD::JetScale, std::string> scale2str( { std::map<std::string, std::string> scale2str( {
{ xAOD::JetEMScaleMomentum , "EMScale" } , { "JetEMScaleMomentum" , "EMScale" } ,
{ xAOD::JetConstitScaleMomentum , "ConstitScale" } } ); { "JetConstitScaleMomentum" , "ConstitScale" } } );
TString scaleTag = scale2str[ (xAOD::JetScale) m_jetScale ] ; // defaults to "" TString scaleTag = scale2str[ m_jetScale ] ; // defaults to ""
TString prefixn = scaleTag; TString prefixn = scaleTag;
if(prefixn != "") prefixn +="_"; if(prefixn != "") prefixn +="_";
// -------------- // --------------
...@@ -107,28 +106,31 @@ int JetKinematicHistos::fillHistosFromContainer(const xAOD::JetContainer & cont) ...@@ -107,28 +106,31 @@ int JetKinematicHistos::fillHistosFromContainer(const xAOD::JetContainer & cont)
int JetKinematicHistos::fillHistosFromJet(const xAOD::Jet &j){ int JetKinematicHistos::fillHistosFromJet(const xAOD::Jet &j){
{ if(m_jetScale != "JetAssignedScaleMomentum" && !j.isAvailable<float>(m_jetScale+"_pt")){
// m_jetScale is a property of the base tool if(m_doNConstit) m_nConstit->Fill( j.numConstituents() );
const xAOD::JetFourMom_t p4 = j.jetP4( (xAOD::JetScale) m_jetScale); return 0;
m_pt->Fill( p4.Pt()*toGeV ); }
m_eta->Fill( p4.Eta() );
m_phi->Fill( p4.Phi() );
if (p4.Pt()*toGeV > 200.0){ // high eta
m_pt_high->Fill( p4.Pt()*toGeV );
m_eta_high->Fill( p4.Eta() );
if(m_doE) m_e_high->Fill( p4.E()*toGeV );
if(m_doM) m_m_high->Fill( p4.M()*toGeV );
if(m_doNConstit) m_nConstit_high->Fill( j.numConstituents() );
}
if(m_doE) m_e->Fill( p4.E()*toGeV ); // m_jetScale is a property of the base tool
if(m_doM) m_m->Fill( p4.M()*toGeV ); const xAOD::JetFourMom_t p4 = j.jetP4(m_jetScale);
m_pt->Fill( p4.Pt()*toGeV );
if(m_doOccupancy) m_occupancyEtaPhi->Fill( p4.Eta(), p4.Phi() ); m_eta->Fill( p4.Eta() );
if(m_doAveragePt) m_averagePtEtaPhi->Fill( p4.Eta(), p4.Phi() , p4.Pt()*toGeV); m_phi->Fill( p4.Phi() );
if(m_doAverageE) m_averageE_EtaPhi->Fill( p4.Eta(), p4.Phi() , p4.E()*toGeV); if (p4.Pt()*toGeV > 200.0){ // high eta
m_pt_high->Fill( p4.Pt()*toGeV );
m_eta_high->Fill( p4.Eta() );
if(m_doE) m_e_high->Fill( p4.E()*toGeV );
if(m_doM) m_m_high->Fill( p4.M()*toGeV );
if(m_doNConstit) m_nConstit_high->Fill( j.numConstituents() );
} }
if(m_doE) m_e->Fill( p4.E()*toGeV );
if(m_doM) m_m->Fill( p4.M()*toGeV );
if(m_doOccupancy) m_occupancyEtaPhi->Fill( p4.Eta(), p4.Phi() );
if(m_doAveragePt) m_averagePtEtaPhi->Fill( p4.Eta(), p4.Phi() , p4.Pt()*toGeV);
if(m_doAverageE) m_averageE_EtaPhi->Fill( p4.Eta(), p4.Phi() , p4.E()*toGeV);
if(m_doNConstit) m_nConstit->Fill( j.numConstituents() ); if(m_doNConstit) m_nConstit->Fill( j.numConstituents() );
return 0; return 0;
} }
......
...@@ -78,8 +78,8 @@ attributeHistoManager.buildKnownTools(compactSpecification) ...@@ -78,8 +78,8 @@ attributeHistoManager.buildKnownTools(compactSpecification)
# Jet histogramming tools # Jet histogramming tools
jhm.addTool( JetKinematicHistos("allkinematics",PlotOccupancy=True, PlotAveragePt=True, PlotNJet=True , PlotNConstit = True) ) jhm.addTool( JetKinematicHistos("allkinematics",PlotOccupancy=True, PlotAveragePt=True, PlotNJet=True , PlotNConstit = True) )
jhm.addTool( JetKinematicHistos("basickinematics") ) jhm.addTool( JetKinematicHistos("basickinematics") )
#jhm.addTool( JetKinematicHistos("basickinematics_emscale", JetScale=0) ) #jhm.addTool( JetKinematicHistos("basickinematics_emscale", JetScale="JetEMScaleMomentum") )
#jhm.addTool( JetKinematicHistos("basickinematics_constscale", JetScale=1) ) #jhm.addTool( JetKinematicHistos("basickinematics_constscale", JetScale="JetConstitScaleMomentum") )
jhm.addTool( LeadingJetsRelations("leadingjetrel", jhm.addTool( LeadingJetsRelations("leadingjetrel",
HistoDef = [ HistoDef = [
......
...@@ -182,8 +182,8 @@ attributeHistoManager.buildKnownTools(compactSpecification) ...@@ -182,8 +182,8 @@ attributeHistoManager.buildKnownTools(compactSpecification)
# Jet histogramming tools # Jet histogramming tools
jhm.addTool( JetKinematicHistos("allkinematics",PlotOccupancy=True, PlotAveragePt=True, PlotNJet=True , PlotNConstit = True) ) jhm.addTool( JetKinematicHistos("allkinematics",PlotOccupancy=True, PlotAveragePt=True, PlotNJet=True , PlotNConstit = True) )
jhm.addTool( JetKinematicHistos("basickinematics") ) jhm.addTool( JetKinematicHistos("basickinematics") )
jhm.addTool( JetKinematicHistos("basickinematics_emscale", JetScale=0) ) jhm.addTool( JetKinematicHistos("basickinematics_emscale", JetScale="JetEMScaleMomentum") )
jhm.addTool( JetKinematicHistos("basickinematics_constscale", JetScale=1) ) jhm.addTool( JetKinematicHistos("basickinematics_constscale", JetScale="JetConstitScaleMomentum") )
jhm.addTool( LeadingJetsRelations("leadingjetrel", jhm.addTool( LeadingJetsRelations("leadingjetrel",
HistoDef = [ HistoDef = [
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment