ATLAS Offline Software
Loading...
Searching...
No Matches
HIEventSelectionTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3*/
4
9//#include "HIPileupTool/HIPileupTool.h"
11
12#include <TString.h>
13
14#include <iomanip>
15#include <iostream>
16
17namespace HI
18{
19
20 const std::string HIEventSelectionTool::m_VERTEX_CNAME="PrimaryVertices";
21 const std::string HIEventSelectionTool::m_MBTS_CNAME="MBTSContainer";
22 const std::string HIEventSelectionTool::m_FWD_CNAME="MBTSForwardEventInfo";
23 const std::string HIEventSelectionTool::m_ZDC_CNAME="ZdcModules";
24 const std::string HIEventSelectionTool::m_SHAPE_CNAME="CaloSums"; //HIEventShape
25
26 HIEventSelectionTool::HIEventSelectionTool(const std::string& name) : asg::AsgTool(name),
27 m_name(name),
28 m_accept("HIEventSelection"),
31 {
32 // MBTS timing cut
33 declareProperty("TimeCut",m_timeCut,"maximum DeltaT cut for MBTS timing cut");
34 // Vertex
35 declareProperty("VertexSumPtCut",m_vtx_min_sumPt,"Minimum sum pT from good tracks for a good vertex");
36 declareProperty("VertexNtrkCut",m_vtx_min_nTrk,"Minimum number good tracks for a good vertex");
37 declareProperty("RejectPileup",m_reject_pileup,"If true, cut on events w/ pileup applied");
38 declareProperty("ExclusiveMode",m_exclusive_mode,"returns at first failure (for speed), doesn't evaluate subsequent cuts");
39 declareProperty("VertexSelectionTool",m_vertexSelectionTool);
40 declareProperty("PileupTool",m_hiPileupTool);
41
42 ATH_MSG_DEBUG("Creating HIEventSelectionTool named" << m_name);
43 }
44
46 {
47 ATH_MSG_DEBUG(Form("Deleting HIEventSelectionTool named %s",m_name.c_str()));
48 }
49
51 {
52 if (m_timeCut > 0.) {
53 m_checkTime = true;
54 ATH_MSG_INFO( "Maximum delta t = " << m_timeCut );
55 m_accept.addCut("goodMBTSTime","good MBTS timing");
56 }
57
58 m_accept.addCut("hasPrimaryVertex","has at least one primary vertex");
59
60 m_accept.addCut("pileupVeto","reject events with topological features indicating pileup");
61
62 if (m_vertexSelectionTool.empty())
63 {
64#ifdef ROOTCORE
65 HI::HIVertexSelectionTool* vxSelTool = new HIVertexSelectionTool("HIVxSel");
66 ATH_CHECK( vxSelTool->initialize() );
67 m_vertexSelectionTool = ToolHandle<HI::IHIVertexSelectionTool>(vxSelTool);
68#endif
69 }
70 ATH_CHECK( m_vertexSelectionTool.retrieve() );
71
72 if (m_hiPileupTool.empty())
73 {
74#ifdef ROOTCORE
75 HI::HIPileupTool* m_hipileuptool = new HIPileupTool("HIpileup");
76 ATH_CHECK( m_hipileuptool->initialize() );
77 m_hiPileupTool = ToolHandle<HI::IHIPileupTool>(m_hipileuptool);
78#endif
79 }
80 ATH_CHECK( m_hiPileupTool.retrieve() );
81
82 return StatusCode::SUCCESS;
83 }
84
86{
87
88 ATH_MSG_INFO( "Finalizing event selection tool." );
89 return StatusCode::SUCCESS;
90
91}
92
94 {
95 asg::AcceptData acceptData ( &m_accept );
96 if(m_checkTime){
97 checkMBTSTime( acceptData );
98 }
99 checkVertex( acceptData );
100 return m_accept;
101 }
102
103 //R.Longo - 14/10/2019 - From checkMBTSTime to checkMBTSTime( asg::AcceptData& acceptData )
105 {
106 const xAOD::ForwardEventInfoContainer* forwardEventInfoContainer = 0;
107 if(evtStore()->retrieve(forwardEventInfoContainer,m_FWD_CNAME).isFailure())
108 {
109 ATH_MSG_ERROR("Could not retrieve ForwardEventInfo object w/ name " << m_FWD_CNAME);
110 return;
111 }
112
113 acceptData.setCutResult("goodMBTSTime",false);
114 const xAOD::ForwardEventInfo* forwardEventInfo = forwardEventInfoContainer->at(0);
115 int countA = forwardEventInfo->countA();
116 int countC = forwardEventInfo->countC();
117 float timeA = forwardEventInfo->timeA();
118 float timeC = forwardEventInfo->timeC();
119
120 if (countA==0||countC==0) return; // if either time undefined, return false;
121
122 float timeDiff = timeA - timeC;
123
124 if ( (std::abs(timeDiff)<m_timeCut) ) acceptData.setCutResult("goodMBTSTime",true);
125
126 return;
127 }
128
129
130 //R.Longo - 14/10/2019 - Introducting acceptData instead of m_accept
132 {
133 asg::AcceptData acceptData ( &m_accept );
134 std::stringstream ss;
135 for(unsigned int icut=0; icut < acceptData.getNCuts(); icut++)
136 {
137 std::string cname=acceptData.getCutName(icut);
138 ss << std::setw(15) << cname
139 << std::setw(10) << acceptData.getCutResult(cname)
140 << std::setw(40) << acceptData.getCutDescription(cname)
141 << std::endl;
142 }
143 return ss.str();
144 }
145
146//R.Longo - 14/10/2019 - From checkVertex to checkVertex( asg::AcceptData& acceptData )
148 {
149 unsigned int numPrimaryVertices=getNPrimaryVertices();
150 acceptData.setCutResult("hasPrimaryVertex", (numPrimaryVertices > 0) );
151
152 const xAOD::ZdcModuleContainer* zdcMod = 0;
153 if(evtStore()->retrieve(zdcMod,m_ZDC_CNAME).isFailure())
154 {
155 ATH_MSG_ERROR("Could not retrieve ZdcModule object w/ name " << m_ZDC_CNAME);
156 return;
157 }
158
159 const xAOD::HIEventShapeContainer* hiev = 0;
160 if(evtStore()->retrieve(hiev,m_SHAPE_CNAME).isFailure())
161 {
162 ATH_MSG_ERROR("Could not retrieve HIEventShape object w/ name " << m_SHAPE_CNAME);
163 return;
164 }
165
166 if (m_reject_pileup) acceptData.setCutResult("pileupVeto", m_hiPileupTool->is_pileup( *hiev, *zdcMod ));
167 else acceptData.setCutResult("pileupVeto", true );
168
169
170 return;
171 }
172
173
175 {
176 unsigned int nPriVx = 0;
177
178 const xAOD::VertexContainer* pv_container = 0;
179 if(evtStore()->retrieve(pv_container,m_VERTEX_CNAME).isFailure())
180 {
181 ATH_MSG_ERROR("Could not retrieve VertexContainer object w/ name " << m_VERTEX_CNAME);
182 return 0;
183 }
184 //if( pv_container->size() < 2) return 0;
185
186 for(const auto vert: *pv_container)
187 {
188 if (m_vertexSelectionTool->accept(*vert)) nPriVx++;
189 }
190
191 return nPriVx;
192
193 }
194
195} // HI
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
static Double_t ss
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
ServiceHandle< StoreGateSvc > & evtStore()
const T * at(size_type n) const
Access an element, as an rvalue.
static const std::string m_FWD_CNAME
virtual StatusCode initialize() override
static const std::string m_SHAPE_CNAME
ToolHandle< HI::IHIPileupTool > m_hiPileupTool
HIEventSelectionTool(const std::string &name="HIEventSelection")
virtual const asg::AcceptInfo & getAcceptInfo() const override
static const std::string m_MBTS_CNAME
void checkVertex(asg::AcceptData &acceptData) const
static const std::string m_ZDC_CNAME
ToolHandle< HI::IHIVertexSelectionTool > m_vertexSelectionTool
void checkMBTSTime(asg::AcceptData &acceptData) const
static const std::string m_VERTEX_CNAME
unsigned int getNPrimaryVertices() const
virtual StatusCode finalize() override
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
virtual StatusCode initialize() override
unsigned int getNCuts() const
Get the number of cuts defined.
Definition AcceptData.h:64
const std::string & getCutDescription(const std::string &cutName) const
Get the description of a cut, based on the cut name.
Definition AcceptData.h:85
void setCutResult(const std::string &cutName, bool cutResult)
Set the result of a cut, based on the cut name (safer)
Definition AcceptData.h:134
const std::string & getCutName(unsigned int cutPosition) const
Get the name of a cut, based on the cut position (slow, avoid usage)
Definition AcceptData.h:78
bool getCutResult(const std::string &cutName) const
Get the result of a cut, based on the cut name (safer)
Definition AcceptData.h:98
AsgTool(const std::string &name)
Constructor specifying the tool instance's name.
Definition AsgTool.cxx:58
unsigned short countC() const
unsigned short countA() const
ZdcModuleContainer_v1 ZdcModuleContainer
ForwardEventInfo_v1 ForwardEventInfo
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
ForwardEventInfoContainer_v1 ForwardEventInfoContainer
HIEventShapeContainer_v2 HIEventShapeContainer
Define the latest version of the container.