Skip to content
Snippets Groups Projects
Commit 3a5d7e19 authored by Elin Bergeaas Kuutmann's avatar Elin Bergeaas Kuutmann Committed by Johannes Elmsheuser
Browse files

HLT RoI and chains per signature in DQ monitoring for Run 3 [ATR-19976] [ATR-21172]

parent 9a8109fd
No related branches found
No related tags found
No related merge requests found
#
# Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
# Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
#
'''@file TrigHLTMonitorAlgorithm.py
......@@ -17,6 +17,7 @@ def TrigGeneralMonConfig(inputFlags):
log_trighlt = logging.getLogger( 'TrigGeneralMonitorAlgorithm' )
################################
# HLT general monitoring
......@@ -38,15 +39,62 @@ def TrigGeneralMonConfig(inputFlags):
# Edit properties of a algorithm
# Placeholder
# Add a generic monitoring tool (a "group" in old language).
#
####### Signature names
signature_names = []
signature_names.append("AllChains") #0
signature_names.append("Electrons") #1
signature_names.append("Gamma") #2
signature_names.append("Muons") #3
signature_names.append("MissingET") #4
signature_names.append("Taus") #5
signature_names.append("Jets") #6
signature_names.append("MinBias") #7
# Add a monitoring tool per signature (a "group" in old language).
hltGroup = helper.addGroup(
trigHLTMonAlg,
'TrigHLTMonitor',
'HLT/ResultMon/'
)
hltAllGroup = helper.addGroup(
trigHLTMonAlg,
'TrigHLTAllMonitor',
'HLT/ResultMon/AllChains'
)
hltEleGroup = helper.addGroup(
trigHLTMonAlg,
'TrigHLTEleMonitor',
'HLT/ResultMon/Electrons'
)
hltGamGroup = helper.addGroup(
trigHLTMonAlg,
'TrigHLTGamMonitor',
'HLT/ResultMon/Gamma'
)
hltMuoGroup = helper.addGroup(
trigHLTMonAlg,
'TrigHLTMuoMonitor',
'HLT/ResultMon/Muons'
)
hltMETGroup = helper.addGroup(
trigHLTMonAlg,
'TrigHLTMETMonitor',
'HLT/ResultMon/MissingET'
)
hltTauGroup = helper.addGroup(
trigHLTMonAlg,
'TrigHLTTauMonitor',
'HLT/ResultMon/Taus'
)
hltJetGroup = helper.addGroup(
trigHLTMonAlg,
'TrigHLTJetMonitor',
'HLT/ResultMon/Jets'
)
#No RoI map for minbias
# Configure histograms
......@@ -68,42 +116,90 @@ def TrigGeneralMonConfig(inputFlags):
from AthenaConfiguration.AutoConfigFlags import GetFileMD
#The HLT chains
### Set up counters ########################################
## The HLT chains
counter_i = 1
for chain_name in GetFileMD(inputFlags.Input.Files)['TriggerMenu']['HLTChains']:
counter_i = counter_i+1
log_trighlt.debug('HLT chain_name = %s',chain_name)
max_hlt_chains = counter_i-1
log_trighlt.info('max_hlt_chains = %i', max_hlt_chains)
log_trighlt.debug('max_hlt_chains = %i', max_hlt_chains)
#The L1 items
## The L1 items
counter_L1 = 1
for item_name in GetFileMD(inputFlags.Input.Files)['TriggerMenu']['L1Items']:
counter_L1 = counter_L1+1
log_trighlt.debug('L1 item_name = %s',item_name)
max_L1_items = counter_L1-1
log_trighlt.info('max_L1_items = %i', max_L1_items)
log_trighlt.debug('max_L1_items = %i', max_L1_items)
##### L1 summary histogram ################################
hltGroup.defineHistogram('L1Events',title='Events per Item at L1;;Events',
path='',xbins=max_L1_items,xmin=1,xmax=max_L1_items+1)
#### HLT summary histograms #########################
hltGroup.defineHistogram('HLT_AllChainsRAW',title='AllChains;;Events',
path='AllChains',xbins=max_hlt_chains,
xmin=1,xmax=max_hlt_chains+1)
#### HLT summary histograms for the signatures #############
hltGroup.defineHistogram('eta,phi;HLT_AllChainsRoIs',type='TH2F',title='AllChains;#eta;#phi',
path='AllChains',xbins=64,xmin=-3.2,xmax=3.2,
ybins=64,ymin=-3.2,ymax=3.2)
triggerstatus = ['RAW','PS'] #all chains or prescaled chains
for sig in signature_names: #loop over signatures
log_trighlt.debug("Signature %s",sig)
### Events
titlename = sig+";;"+"Events"
for trigstatus in triggerstatus: #loop over RAW/PS
histname = "HLT_"+sig+trigstatus
log_trighlt.debug('Histname = %s', histname)
hltGroup.defineHistogram(histname,title=titlename,
path=sig,xbins=max_hlt_chains+1,
xmin=1,xmax=max_hlt_chains+1)
#### Summary histograms for the signatures #################
hltGroup.defineHistogram('HLT_ElectronsRAW',title='Placeholder;X;Y',
path='Electrons',xbins=1000000,xmin=-0.5,xmax=999999.5)
### RoIs, one per signature
histname = "eta,phi;HLT_"+signature_names[0]+"RoIs" #AllChains
titlename = signature_names[0]+";#eta;#phi"
hltAllGroup.defineHistogram(histname,type='TH2F',title=titlename,
path='',xbins=64,xmin=-3.2,xmax=3.2,
ybins=64,ymin=-3.2,ymax=3.2)
histname = "eta,phi;HLT_"+signature_names[1]+"RoIs" #Electrons
titlename = signature_names[1]+";#eta;#phi"
hltEleGroup.defineHistogram(histname,type='TH2F',title=titlename,
path='',xbins=64,xmin=-3.2,xmax=3.2,
ybins=64,ymin=-3.2,ymax=3.2)
histname = "eta,phi;HLT_"+signature_names[2]+"RoIs" #Gamma
titlename = signature_names[2]+";#eta;#phi"
hltGamGroup.defineHistogram(histname,type='TH2F',title=titlename,
path='',xbins=64,xmin=-3.2,xmax=3.2,
ybins=64,ymin=-3.2,ymax=3.2)
histname = "eta,phi;HLT_"+signature_names[3]+"RoIs" #Muons
titlename = signature_names[3]+";#eta;#phi"
hltMuoGroup.defineHistogram(histname,type='TH2F',title=titlename,
path='',xbins=64,xmin=-3.2,xmax=3.2,
ybins=64,ymin=-3.2,ymax=3.2)
histname = "eta,phi;HLT_"+signature_names[4]+"RoIs" #MissingET
titlename = signature_names[4]+";#eta;#phi"
hltMETGroup.defineHistogram(histname,type='TH2F',title=titlename,
path='',xbins=64,xmin=-3.2,xmax=3.2,
ybins=64,ymin=-3.2,ymax=3.2)
histname = "eta,phi;HLT_"+signature_names[5]+"RoIs" #Taus
titlename = signature_names[5]+";#eta;#phi"
hltTauGroup.defineHistogram(histname,type='TH2F',title=titlename,
path='',xbins=64,xmin=-3.2,xmax=3.2,
ybins=64,ymin=-3.2,ymax=3.2)
histname = "eta,phi;HLT_"+signature_names[6]+"RoIs" #Jets
titlename = signature_names[6]+";#eta;#phi"
hltJetGroup.defineHistogram(histname,type='TH2F',title=titlename,
path='',xbins=64,xmin=-3.2,xmax=3.2,
ybins=64,ymin=-3.2,ymax=3.2)
#####################################
......
#
# Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
# Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
#
'''@file TrigHLTMonitorAlgorithm.py
......@@ -85,3 +85,4 @@ def TrigHLTMonTopConfig(inputFlags):
return result
/*
Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
*/
#include "TrigHLTMonitorAlgorithm.h"
......@@ -31,14 +31,21 @@ StatusCode TrigHLTMonitorAlgorithm::fillHistograms( const EventContext& ctx ) co
using namespace Monitored;
StatusCode sc = StatusCode::FAILURE;
//Fetch the tools
auto tool = getGroup("TrigHLTMonitor");
auto toolAll = getGroup("TrigHLTAllMonitor");
auto toolEle = getGroup("TrigHLTEleMonitor");
auto toolGam = getGroup("TrigHLTGamMonitor");
auto toolMuo = getGroup("TrigHLTMuoMonitor");
auto toolMET = getGroup("TrigHLTMETMonitor");
auto toolTau = getGroup("TrigHLTTauMonitor");
auto toolJet = getGroup("TrigHLTJetMonitor");
////////////////////////////////////////
// Declare the quantities which should be monitored
// NB! The variables and histograms defined here must match the ones in the py file exactly!
auto HLT_ElectronsRAW = Monitored::Scalar<int>("HLT_ElectronsRAW",0);
auto HLT_AllChainsRAW = Monitored::Scalar<std::string>("HLT_AllChainsRAW");
auto L1Events = Monitored::Scalar<std::string>("L1Events");
......@@ -52,11 +59,10 @@ StatusCode TrigHLTMonitorAlgorithm::fillHistograms( const EventContext& ctx ) co
// all others are the variables to be saved.
// Alternative fill method. Get the group yourself, and pass it to the fill function.
ATH_MSG_INFO("Filling L1Events histogram");
ATH_MSG_DEBUG("Filling L1Events histogram");
for(unsigned int it=0; it<nL1Items; ++it) {
if( L1items[it] != "" ) {
ATH_MSG_DEBUG("L1Item " << it << " " << L1items << m_trigDecTool->isPassed(L1items[it]));
ATH_MSG_DEBUG("L1Item " << it << " " << L1items[it] );
if(m_trigDecTool->isPassed(L1items[it])) {
/// Fill L1 histogram
......@@ -70,57 +76,140 @@ StatusCode TrigHLTMonitorAlgorithm::fillHistograms( const EventContext& ctx ) co
////////////////////////////////////
// HLT chain monitoring
ATH_MSG_DEBUG( "HLT chain monitoring" );
ATH_MSG_DEBUG( "Setting up the regex map..." );
std::map<std::string,std::string> streams;
streams.insert(std::make_pair("AllChains", "HLT_.*"));
streams.insert(std::make_pair("Electrons", "HLT_e[0-9]+.*"));
streams.insert(std::make_pair("Gamma", "HLT_g[0-9]+.*"));
streams.insert(std::make_pair("Muons", "HLT_[0-9]*mu[0-9]+.*"));
streams.insert(std::make_pair("MissingET", "HLT_(t|x)e[0-9]+.*"));
streams.insert(std::make_pair("Taus", "HLT_(tau[0-9]*|trk.*Tau).*"));
streams.insert(std::make_pair("Jets", "HLT_(FJ|j)[0-9]+.*"));
streams.insert(std::make_pair("MinBias", "HLT_mb.*"));
//// Set the values of the monitored variables for the event
HLT_ElectronsRAW = GetEventInfo(ctx)->runNumber(); //DUMMY; to be updated.
ATH_MSG_INFO( "Filling HLT_AllChains and RoI information" );
std::vector< std::string > chainNames = m_trigDecTool->getChainGroup("HLT_.*")->getListOfTriggers();
unsigned int nHLTChains = chainNames.size();
for(unsigned int ith=0; ith<nHLTChains; ++ith) {
if( chainNames[ith] != "" ) {
ATH_MSG_DEBUG("HLTChain " << ith << " " << chainNames << m_trigDecTool->isPassed(chainNames[ith]));
if(m_trigDecTool->isPassed(chainNames[ith])) {
ATH_MSG_DEBUG( " Chain " << chainNames[ith] << " IS passed");
/// Fill plain chains histogram
HLT_AllChainsRAW = chainNames[ith];
ATH_MSG_DEBUG( "Fill HLT_AllChainsRAW" ); // with " << HLT_AllChainsRAW );
fill(tool,HLT_AllChainsRAW);
ATH_MSG_DEBUG( "Iterating over the regex map...");
std::map<std::string,std::string>::const_iterator strItr;
for (strItr=streams.begin();strItr!=streams.end(); ++strItr){
std::string signaturename = strItr->first;
std::string thisregex = strItr->second;
std::string histname_raw = "HLT_"+signaturename+"RAW";
std::string histname_ps = "HLT_"+signaturename+"PS";
auto HLT_RAW = Monitored::Scalar<std::string>(histname_raw);
auto HLT_PS = Monitored::Scalar<std::string>(histname_ps);
ATH_MSG_DEBUG( "Filling HLT" << signaturename << " and RoI information for " << thisregex );
std::vector< std::string > chainNames = m_trigDecTool->getChainGroup(thisregex)->getListOfTriggers();
unsigned int nHLTChains = chainNames.size();
for(unsigned int ith=0; ith<nHLTChains; ++ith) {
if( chainNames[ith] != "" ) {
ATH_MSG_DEBUG("HLTChain " << ith << " " << chainNames[ith] );
if(m_trigDecTool->isPassed(chainNames[ith])) {
ATH_MSG_DEBUG( " Chain " << chainNames[ith] << " IS passed");
/// Fill plain chains histogram
HLT_RAW = chainNames[ith];
ATH_MSG_DEBUG( "Fill HLT_RAW for " << signaturename << " and " << chainNames[ith]);
fill(tool,HLT_RAW);
/// Fill RoIs histogram
ATH_MSG_DEBUG("Fill RoI histograms for chain " << chainNames[ith] );
std::vector<LinkInfo<TrigRoiDescriptorCollection>> fvec = m_trigDecTool->features<TrigRoiDescriptorCollection>(chainNames[ith], TrigDefs::Physics, "", TrigDefs::lastFeatureOfType, "initialRoI"); //initialRoiString()
//Loop over RoIs
for (const LinkInfo<TrigRoiDescriptorCollection>& li : fvec) {
if( li.isValid() ) {
auto phi = Monitored::Scalar("phi", 0.0);
auto eta = Monitored::Scalar("eta", 0.0);
auto HLT_AllChainsRoIs = Monitored::Group(tool, eta, phi); //filled when going out of scope
const TrigRoiDescriptor* roi = *(li.link).cptr();
eta = roi->eta();
phi = roi->phi();
ATH_MSG_DEBUG( "RoI: eta = " << eta << ", phi = " << phi );
}
else {
ATH_MSG_WARNING( "TrigRoiDescriptorCollection for chain " << chainNames[ith] << " is not valid");
//If the chain is prescaled
ATH_MSG_DEBUG( "Prescale: " << m_trigDecTool->getPrescale(chainNames[ith]) );
if(m_trigDecTool->getPrescale(chainNames[ith])!=1) {
//NB! Right now very few chains are prescaled, so this histogram is seldom filled
HLT_PS = chainNames[ith];
ATH_MSG_DEBUG( "Fill HLT_PS for " << signaturename << " and " << chainNames[ith]);
fill(tool,HLT_PS);
}
}//end for (const LinkInfo<TrigRoiDescriptorCollection>& li : fvec)
}// end if(m_trigDecTool->isPassed(chainNames[ith]))
}// end if( chainNames[ith] != "" )
}//end for(unsigned int ith=0; ith<nHLTChains; ++ith)
/// Fill RoIs histogram
ATH_MSG_DEBUG("Fill RoI histograms for chain " << chainNames[ith] );
std::vector<LinkInfo<TrigRoiDescriptorCollection>> fvec = m_trigDecTool->features<TrigRoiDescriptorCollection>(chainNames[ith], TrigDefs::Physics, "", TrigDefs::lastFeatureOfType, initialRoIString());
//Loop over RoIs
for (const LinkInfo<TrigRoiDescriptorCollection>& li : fvec) {
if( li.isValid() ) {
auto phi = Monitored::Scalar("phi",0.0);
auto eta = Monitored::Scalar("eta",0.0);
if(signaturename=="AllChains") {
ATH_MSG_DEBUG( "RoI: filling for " << signaturename );
auto HLT_RoIs = Monitored::Group(toolAll, eta, phi);
const TrigRoiDescriptor* roi = *(li.link).cptr();
eta = roi->eta();
phi = roi->phi();
ATH_MSG_DEBUG( "RoI: eta = " << eta << ", phi = " << phi );
}
//Check signatures
if(signaturename=="Electrons") {
ATH_MSG_DEBUG( "RoI: filling for " << signaturename );
auto HLT_RoIs = Monitored::Group(toolEle, eta, phi);
const TrigRoiDescriptor* roi = *(li.link).cptr();
eta = roi->eta();
phi = roi->phi();
ATH_MSG_DEBUG( "RoI: eta = " << eta << ", phi = " << phi );
}
else if(signaturename=="Gamma") {
ATH_MSG_DEBUG( "RoI: filling for " << signaturename );
auto HLT_RoIs = Monitored::Group(toolGam, eta, phi);
const TrigRoiDescriptor* roi = *(li.link).cptr();
eta = roi->eta();
phi = roi->phi();
ATH_MSG_DEBUG( "RoI: eta = " << eta << ", phi = " << phi );
}
else if(signaturename=="Muons") {
ATH_MSG_DEBUG( "RoI: filling for " << signaturename );
auto HLT_RoIs = Monitored::Group(toolMuo, eta, phi);
const TrigRoiDescriptor* roi = *(li.link).cptr();
eta = roi->eta();
phi = roi->phi();
ATH_MSG_DEBUG( "RoI: eta = " << eta << ", phi = " << phi );
}
else if(signaturename=="MissingET") {
ATH_MSG_DEBUG( "RoI: filling for " << signaturename );
auto HLT_RoIs = Monitored::Group(toolMET, eta, phi);
const TrigRoiDescriptor* roi = *(li.link).cptr();
eta = roi->eta();
phi = roi->phi();
ATH_MSG_DEBUG( "RoI: eta = " << eta << ", phi = " << phi );
}
else if(signaturename=="Taus") {
ATH_MSG_DEBUG( "RoI: filling for " << signaturename );
auto HLT_RoIs = Monitored::Group(toolTau, eta, phi);
const TrigRoiDescriptor* roi = *(li.link).cptr();
eta = roi->eta();
phi = roi->phi();
ATH_MSG_DEBUG( "RoI: eta = " << eta << ", phi = " << phi );
}
else if(signaturename=="Jets") {
ATH_MSG_DEBUG( "RoI: filling for " << signaturename );
auto HLT_RoIs = Monitored::Group(toolJet, eta, phi);
const TrigRoiDescriptor* roi = *(li.link).cptr();
eta = roi->eta();
phi = roi->phi();
ATH_MSG_DEBUG( "RoI: eta = " << eta << ", phi = " << phi );
}
}//end if(li.isValid())
else {
ATH_MSG_WARNING( "TrigRoiDescriptorCollection for chain " << chainNames[ith] << " is not valid");
}
}//end for (const LinkInfo<TrigRoiDescriptorCollection>& li : fvec)
}// end if(m_trigDecTool->isPassed(chainNames[ith]))
}// end if( chainNames[ith] != "" )
}//end for(unsigned int ith=0; ith<nHLTChains; ++ith)
}//end loop over streams
//Placeholder
ATH_MSG_DEBUG( "Fill HLT_ElectronsRAW with " << HLT_ElectronsRAW );
fill(tool,HLT_ElectronsRAW);
//////////////////////////////////////
......@@ -131,7 +220,7 @@ StatusCode TrigHLTMonitorAlgorithm::fillHistograms( const EventContext& ctx ) co
//////////////////////////////////////
// End
ATH_MSG_INFO( "Finalizing the TrigHLTMonitorAlgorithm..." );
ATH_MSG_DEBUG( "Finalizing the TrigHLTMonitorAlgorithm..." );
return StatusCode::SUCCESS;
}
......@@ -168,11 +257,11 @@ StatusCode TrigHLTMonitorAlgorithm::fillResultAndConsistencyHistograms( const Ev
//THIS IS NOT IMPLEMENTED YET IN Run3
//bskeys_1 = HLTResult->getConfigSuperMasterKey();
//bskeys_2 = HLTResult->getConfigPrescalesKey();
ATH_MSG_INFO("sc_hltResult should not be SUCCESS");
ATH_MSG_DEBUG("sc_hltResult should not be SUCCESS");
}
else {
//For now, this will always happen
ATH_MSG_INFO("No HLTResult monitored (because it has not been implemented yet)");
ATH_MSG_DEBUG("====> No HLTResult monitored (because it has not been implemented yet)");
//FIXME: When the HLTResult IS implemented this should be a WARNING // ebergeas Sept 2020
//ATH_MSG_WARNING("No HLTResult found"); //Prints a warning message in the log file
ConfigConsistency_HLT=7;
......
/*
Copyright (C) 2002-2020 CERN for the benefit of the ATLAS collaboration
Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
*/
#ifndef TRIGHLTMONITORING_TRIGHLTMONITORALGORITHM_H
......@@ -29,5 +29,6 @@ class TrigHLTMonitorAlgorithm : public AthMonitorAlgorithm {
ToolHandle<Trig::TrigDecisionTool> m_trigDecTool;
ServiceHandle< TrigConf::ITrigConfigSvc > m_trigConfigSvc{ this, "TrigConfigSvc", "" };
StatusCode fillResultAndConsistencyHistograms( const EventContext& ctx ) const;
};
#endif
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