ATLAS Offline Software
Loading...
Searching...
No Matches
DerivationFramework::ZeeVertexRefittingTool Class Reference

#include <ZeeVertexRefittingTool.h>

Inheritance diagram for DerivationFramework::ZeeVertexRefittingTool:
Collaboration diagram for DerivationFramework::ZeeVertexRefittingTool:

Public Member Functions

virtual StatusCode initialize () override final
virtual StatusCode finalize () override final
virtual StatusCode addBranches (const EventContext &ctx) const override final

Private Member Functions

StatusCode makeZeePairs (const xAOD::ElectronContainer *particles, std::vector< std::vector< unsigned int > > &ZeePairs) const

Private Attributes

Gaudi::Property< std::string > m_expression {this, "ObjectRequirements", "true"}
Gaudi::Property< float > m_massCut {this, "LowMassCut", 0.f}
SG::ReadHandleKey< xAOD::VertexContainerm_primaryVertexKey {this, "PVContainerName", "PrimaryVertices", "" }
SG::ReadHandleKey< xAOD::ElectronContainerm_electronKey { this, "ElectronContainerName", "Electrons", "" }
SG::WriteHandleKey< xAOD::VertexContainerm_refitpvKey {this, "RefittedPVContainerName", "HggPrimaryVertices", "" }
SG::ReadHandleKey< xAOD::EventInfom_eventInfoKey { this, "EventInfoKey", "EventInfo", "" }
Gaudi::Property< std::vector< unsigned int > > m_MCSamples {this, "MCSamples", {} }
ToolHandle< Analysis::PrimaryVertexRefitterm_pvrefitter {this, "PrimaryVertexRefitterTool", "Analysis::PrimaryVertexRefitter"}

Detailed Description

Definition at line 36 of file ZeeVertexRefittingTool.h.

Member Function Documentation

◆ addBranches()

StatusCode DerivationFramework::ZeeVertexRefittingTool::addBranches ( const EventContext & ctx) const
finaloverridevirtual

Definition at line 59 of file ZeeVertexRefittingTool.cxx.

60 {
61 // skip mc samples not included in the MCSamples list
62 SG::ReadHandle<xAOD::EventInfo> eventInfo (m_eventInfoKey, ctx);
63
64 if(eventInfo->eventType(xAOD::EventInfo::IS_SIMULATION)) {
65 bool skipSample = true;
66 for (auto mcid : m_MCSamples) {
67 if (mcid==eventInfo->mcChannelNumber()) {
68 skipSample = false;
69 break;
70 }
71 }
72 if (skipSample) return StatusCode::SUCCESS;
73 }
74
75 // check that container we want to write in SG does not yet exist
76 if (evtStore()->contains<xAOD::VertexContainer>(m_refitpvKey.key())) {
77 ATH_MSG_ERROR("Tool is attempting to write a StoreGate key " << m_refitpvKey.key() << " which already exists. Please use a different key");
78 return StatusCode::FAILURE;
79 }
80
81
82 SG::ReadHandle<xAOD::VertexContainer> pv_cont (m_primaryVertexKey, ctx);
83
84 xAOD::VertexContainer* refittedPVContainer = new xAOD::VertexContainer;
85 xAOD::VertexAuxContainer* refittedPVAuxContainer = new xAOD::VertexAuxContainer;
86 refittedPVContainer->setStore( refittedPVAuxContainer );
87
88 SG::WriteHandle<xAOD::VertexContainer> vertexContainer(m_refitpvKey, ctx);
89 ATH_CHECK(vertexContainer.recordNonConst(std::unique_ptr< xAOD::VertexContainer >(refittedPVContainer),
90 std::unique_ptr< xAOD::VertexAuxContainer >(refittedPVAuxContainer)));
91
92 const xAOD::Vertex* pv = nullptr;
93 for ( const auto *v : *pv_cont ) {
94 if (v->vertexType()==xAOD::VxType::PriVtx) {
95 pv = v;
96 break;
97 }
98 }
99 if (!pv) {
100 return StatusCode::SUCCESS;
101 }
102 ATH_MSG_DEBUG("Found PV");
103
104 // retrieve particle collections
105 SG::ReadHandle<xAOD::ElectronContainer> electrons (m_electronKey, ctx);
106
107 // create the vector which will hold the Zee pairs
108 std::vector< std::vector<unsigned int> > eepairs;
109 CHECK( makeZeePairs( &*electrons, eepairs ) );
110
111 ATH_MSG_DEBUG("ee pairs found: " << eepairs.size());
112
113 for (auto pair : eepairs) {
114 std::vector<const xAOD::TrackParticle*> tps = {
117 };
118
119 TLorentzVector v0, v1, egamVec;
120 if(electrons->at(pair[0])->caloCluster())
121 {
122 v0.SetPtEtaPhiM(electrons->at(pair[0])->e()/cosh(electrons->at(pair[0])->caloCluster()->etaBE(2)),
123 electrons->at(pair[0])->caloCluster()->etaBE(2),
124 electrons->at(pair[0])->caloCluster()->phiBE(2),
125 0.0);
126 }
127
128 if(electrons->at(pair[1])->caloCluster())
129 {
130 v1.SetPtEtaPhiM(electrons->at(pair[1])->e()/cosh(electrons->at(pair[1])->caloCluster()->etaBE(2)),
131 electrons->at(pair[1])->caloCluster()->etaBE(2),
132 electrons->at(pair[1])->caloCluster()->phiBE(2),
133 0.0);
134 }
135
136 egamVec = v0+v1;
137
138 ATH_MSG_DEBUG("Refitting PV for e tracks: " << tps[0] << " " << tps[1]);
139 xAOD::Vertex* pv_ref = m_pvrefitter->refitVertex(pv,tps);
140 if (pv_ref) {
141 refittedPVContainer->push_back(pv_ref);
142
143 int ipv = 0;
144 for (const xAOD::Vertex* v : *pv_cont) {
145 xAOD::Vertex * nv = new xAOD::Vertex();
146 nv->makePrivateStore(v);
147 if(ipv !=0 ) refittedPVContainer->push_back(nv);
148 ipv++;
149 }
150
151 for ( xAOD::Vertex *v : *refittedPVContainer )
152 {
153 float vert_sumpt = (log10(xAOD::PVHelpers::getVertexSumPt(v)));
154 float vert_sumpt2 = (log10(xAOD::PVHelpers::getVertexSumPt(v,2, false)));
155
156 TLorentzVector vtxmom = xAOD::PVHelpers::getVertexMomentum(v, true, "");
157 float vert_dphi = (fabs(vtxmom.DeltaPhi(egamVec)));
158 //fill vertex variables
159 vertices_sumPt(*v) = vert_sumpt;
160 vertices_sumPt2(*v) = vert_sumpt2;
161 vertices_dPhi(*v) = vert_dphi;
162 }
163
164
165 ATH_MSG_DEBUG("refitted PV nTP: " << pv_ref->nTrackParticles() << " -- " << pv->nTrackParticles());
166 ATH_MSG_DEBUG("refitted PV z: " << pv_ref->z() << " -- " << pv->z());
167
168 if (pv_ref->nTrackParticles() < pv->nTrackParticles()) {
169 sumPt2(*pv_ref) = xAOD::PVHelpers::getVertexSumPt(pv_ref, 2, false);
170 //set links to electrons, used only for matching, not for vertexing
171 std::vector<ElementLink<xAOD::TrackParticleContainer> > electronTrackLinks = {
172 electrons->at(pair[0])->trackParticleLink(),
173 electrons->at(pair[1])->trackParticleLink()
174 };
175 electronTrackLinksDecor(*pv_ref) = electronTrackLinks;
176 } else {
177 ATH_MSG_DEBUG("Electrons from pair not used in the refitting ");
178 }
179 }
180 else {
181 ATH_MSG_DEBUG("refitting failed");
182 }
183 }
184
185 ATH_MSG_DEBUG("Vertex container size: " << refittedPVContainer->size());
186
187 return StatusCode::SUCCESS;
188 }
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
#define CHECK(...)
Evaluate an expression and check for errors.
value_type push_back(value_type pElem)
Add an element to the end of the collection.
StatusCode makeZeePairs(const xAOD::ElectronContainer *particles, std::vector< std::vector< unsigned int > > &ZeePairs) const
SG::ReadHandleKey< xAOD::ElectronContainer > m_electronKey
SG::ReadHandleKey< xAOD::VertexContainer > m_primaryVertexKey
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfoKey
ToolHandle< Analysis::PrimaryVertexRefitter > m_pvrefitter
Gaudi::Property< std::vector< unsigned int > > m_MCSamples
SG::WriteHandleKey< xAOD::VertexContainer > m_refitpvKey
void makePrivateStore()
Create a new (empty) private store for this object.
@ IS_SIMULATION
true: simulation, false: data
float z() const
Returns the z position.
size_t nTrackParticles() const
Get the number of tracks associated with this vertex.
bool contains(const std::string &s, const std::string &regx)
does a string contain the substring
Definition hcg.cxx:114
static const SG::AuxElement::Decorator< float > sumPt2("sumPt2")
static const SG::AuxElement::Decorator< float > vertices_sumPt("vertices_sumPt")
static const SG::AuxElement::Decorator< std::vector< ElementLink< xAOD::TrackParticleContainer > > > electronTrackLinksDecor("ElectronTrackLinks")
static const SG::AuxElement::Decorator< float > vertices_sumPt2("vertices_sumPt2")
static const SG::AuxElement::Decorator< float > vertices_dPhi("vertices_dPhi")
const xAOD::TrackParticle * getOriginalTrackParticle(const xAOD::Electron *el)
Helper function for getting the "Original" Track Particle (i.e before GSF) via the electron.
float getVertexSumPt(const xAOD::Vertex *vertex, int power=1, bool useAux=true)
Loop over track particles associated with vertex and return scalar sum of pT^power in GeV (from auxda...
TLorentzVector getVertexMomentum(const xAOD::Vertex *vertex, bool useAux=true, const std::string &derivationPrefix="")
Return vector sum of tracks associated with vertex (from auxdata if available and useAux = true)
@ PriVtx
Primary vertex.
VertexAuxContainer_v1 VertexAuxContainer
Definition of the current jet auxiliary container.
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.

◆ finalize()

StatusCode DerivationFramework::ZeeVertexRefittingTool::finalize ( )
finaloverridevirtual

Definition at line 53 of file ZeeVertexRefittingTool.cxx.

54 {
55 ATH_CHECK( finalizeParser() );
56 return StatusCode::SUCCESS;
57 }

◆ initialize()

StatusCode DerivationFramework::ZeeVertexRefittingTool::initialize ( )
finaloverridevirtual

Definition at line 31 of file ZeeVertexRefittingTool.cxx.

32 {
33
34 ATH_CHECK( m_pvrefitter.retrieve() );
35
36 if (!m_expression.empty()) {
37 ATH_CHECK(initializeParser(m_expression));
38 } else {
39 ATH_CHECK(initializeParser("true"));
40 }
41 ATH_CHECK( m_eventInfoKey.initialize() );
42 ATH_CHECK( m_primaryVertexKey.initialize() );
43 ATH_CHECK( m_electronKey.initialize() );
44 ATH_CHECK( m_refitpvKey.initialize() );
45 if (m_refitpvKey.key().empty()) {
46 ATH_MSG_ERROR("No SG name provided for the output of ZeeVertexRefittingTool!");
47 return StatusCode::FAILURE;
48 }
49
50 return StatusCode::SUCCESS;
51 }

◆ makeZeePairs()

StatusCode DerivationFramework::ZeeVertexRefittingTool::makeZeePairs ( const xAOD::ElectronContainer * particles,
std::vector< std::vector< unsigned int > > & ZeePairs ) const
private

Definition at line 190 of file ZeeVertexRefittingTool.cxx.

191 {
192 if (particles->size()<2) return StatusCode::SUCCESS;
193
194 // flags for the result of selection for each electron
195 std::vector<int> isSelected = m_parser->evaluateAsVector();
196 unsigned int nEntries = isSelected.size();
197
198 // if there are no particles in one of the two lists to combine, just leave function
199 if (nEntries==0) return StatusCode::SUCCESS;
200
201 // check the sizes are compatible
202 if (particles->size() != nEntries ) {
203 ATH_MSG_ERROR("Branch sizes incompatible - returning zero");
204 return StatusCode::FAILURE;
205 }
206
207 // Double loop to get the opposite-charge pairs with m>50 GeV
208 for (unsigned int i=0; i<nEntries-1; ++i) {
209 if (isSelected[i]!=1) continue;
210 float qi = particles->at(i)->charge();
211
212 for (unsigned int j=i+1; j<nEntries; ++j) {
213 if (isSelected[j]!=1) continue;
214 //std::vector<int> tmpPair; tmpPair.clear();
215 float qj = particles->at(j)->charge();
216 // opposite charge
217 if (qi*qj>=0) continue;
218
219 //invariant mass
220 float mass = (particles->at(i)->p4()+particles->at(j)->p4()).M();
221 if (mass<m_massCut) continue;
222
223 ZeePairs.push_back( {i, j} );
224 }
225 }
226 return StatusCode::SUCCESS;
227 }

Member Data Documentation

◆ m_electronKey

SG::ReadHandleKey<xAOD::ElectronContainer> DerivationFramework::ZeeVertexRefittingTool::m_electronKey { this, "ElectronContainerName", "Electrons", "" }
private

Definition at line 50 of file ZeeVertexRefittingTool.h.

50{ this, "ElectronContainerName", "Electrons", "" };

◆ m_eventInfoKey

SG::ReadHandleKey<xAOD::EventInfo> DerivationFramework::ZeeVertexRefittingTool::m_eventInfoKey { this, "EventInfoKey", "EventInfo", "" }
private

Definition at line 52 of file ZeeVertexRefittingTool.h.

52{ this, "EventInfoKey", "EventInfo", "" };

◆ m_expression

Gaudi::Property<std::string> DerivationFramework::ZeeVertexRefittingTool::m_expression {this, "ObjectRequirements", "true"}
private

Definition at line 46 of file ZeeVertexRefittingTool.h.

46{this, "ObjectRequirements", "true"};

◆ m_massCut

Gaudi::Property<float> DerivationFramework::ZeeVertexRefittingTool::m_massCut {this, "LowMassCut", 0.f}
private

Definition at line 47 of file ZeeVertexRefittingTool.h.

47{this, "LowMassCut", 0.f};

◆ m_MCSamples

Gaudi::Property<std::vector<unsigned int> > DerivationFramework::ZeeVertexRefittingTool::m_MCSamples {this, "MCSamples", {} }
private

Definition at line 54 of file ZeeVertexRefittingTool.h.

54{this, "MCSamples", {} };

◆ m_primaryVertexKey

SG::ReadHandleKey<xAOD::VertexContainer> DerivationFramework::ZeeVertexRefittingTool::m_primaryVertexKey {this, "PVContainerName", "PrimaryVertices", "" }
private

Definition at line 49 of file ZeeVertexRefittingTool.h.

49{this, "PVContainerName", "PrimaryVertices", "" };

◆ m_pvrefitter

ToolHandle< Analysis::PrimaryVertexRefitter > DerivationFramework::ZeeVertexRefittingTool::m_pvrefitter {this, "PrimaryVertexRefitterTool", "Analysis::PrimaryVertexRefitter"}
private

Definition at line 56 of file ZeeVertexRefittingTool.h.

56{this, "PrimaryVertexRefitterTool", "Analysis::PrimaryVertexRefitter"};

◆ m_refitpvKey

SG::WriteHandleKey<xAOD::VertexContainer> DerivationFramework::ZeeVertexRefittingTool::m_refitpvKey {this, "RefittedPVContainerName", "HggPrimaryVertices", "" }
private

Definition at line 51 of file ZeeVertexRefittingTool.h.

51{this, "RefittedPVContainerName", "HggPrimaryVertices", "" };

The documentation for this class was generated from the following files: