ATLAS Offline Software
Loading...
Searching...
No Matches
Reco_V0Finder.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
5// Reco_V0Finder.cxx, (c) ATLAS Detector software
7// Author: Adam Barton
8#include "Reco_V0Finder.h"
11
12namespace DerivationFramework {
13
14 Reco_V0Finder::Reco_V0Finder(const std::string& t,
15 const std::string& n,
16 const IInterface* p) :
17 base_class(t,n,p),
18 m_v0FinderTool("InDet::V0FinderTool", this)
19 {
20
21 // Declare user-defined properties
22 declareProperty("CheckVertexContainers", m_CollectionsToCheck);
23 declareProperty("V0FinderTool", m_v0FinderTool);
24 }
25
26 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
27
29 {
30
31 ATH_MSG_DEBUG("in initialize()");
32 // get the V0Finder tool
33 ATH_CHECK( m_v0FinderTool.retrieve());
34 ATH_CHECK( m_v0DecoTool.retrieve());
35
36 ATH_CHECK(m_vertexKey.initialize());
37 ATH_CHECK(m_v0Key.initialize());
38 ATH_CHECK(m_ksKey.initialize());
39 ATH_CHECK(m_laKey.initialize());
40 ATH_CHECK(m_lbKey.initialize());
41
42 return StatusCode::SUCCESS;
43
44 }
45
46
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * *
48
49 StatusCode Reco_V0Finder::addBranches(const EventContext& ctx) const
50 {
51
52 bool callV0Finder = false;
53 // Jpsi container and its auxilliary store
54 for(const auto &str : m_CollectionsToCheck){
55 const xAOD::VertexContainer* vertContainer = nullptr;
56 ATH_CHECK( evtStore()->retrieve(vertContainer, str) );
57 if(vertContainer->size() == 0) {
58 ATH_MSG_DEBUG("Container VertexContainer (" << str << ") is empty");
59 }else{
60 callV0Finder = true;
61 ATH_MSG_DEBUG("Container VertexContainer (" << str << ") has events N= " << vertContainer->size());
62 break;//No point checking other containers
63 }
64 }
65
66 // InDetV0 container and its auxilliary store
67 //---- Recording section: write the results to StoreGate ---//
69 if ( h_V0.record(std::make_unique<xAOD::VertexContainer>() ,std::make_unique<xAOD::VertexAuxContainer>()).isFailure()){
70 ATH_MSG_ERROR("Storegate record of v0Container failed.");
71 return StatusCode::FAILURE;
72 }
73
75 if ( h_Ks.record(std::make_unique<xAOD::VertexContainer>() ,std::make_unique<xAOD::VertexAuxContainer>()).isFailure()){
76 ATH_MSG_ERROR("Storegate record of ksContainer failed.");
77 return StatusCode::FAILURE;
78 }
79
81 if( h_La.record(std::make_unique<xAOD::VertexContainer>() ,std::make_unique<xAOD::VertexAuxContainer>()).isFailure()){
82 ATH_MSG_ERROR("Storegate record of laContainer failed.");
83 return StatusCode::FAILURE;
84
85 }
87 if(h_Lb.record(std::make_unique<xAOD::VertexContainer>() ,std::make_unique<xAOD::VertexAuxContainer>()).isFailure()){
88 ATH_MSG_ERROR("Storegate record of lbContainer failed.");
89 return StatusCode::FAILURE;
90 }
91
92 xAOD::VertexContainer* v0Container(h_V0.ptr());
93 xAOD::VertexContainer* ksContainer(h_Ks.ptr());
94 xAOD::VertexContainer* laContainer(h_La.ptr());
95 xAOD::VertexContainer* lbContainer(h_Lb.ptr());
96 // Call V0Finder
97 if (callV0Finder) {
98 const xAOD::Vertex * primaryVertex(0);
99 SG::ReadHandle<xAOD::VertexContainer> importedVxContainer( m_vertexKey, ctx );
100 ATH_CHECK(importedVxContainer.isValid());
101
102 if (importedVxContainer->size()==0){
103 ATH_MSG_WARNING("You have no primary vertices: " << importedVxContainer->size());
104 } else {
105 primaryVertex = (*importedVxContainer)[0];
106 }
107 ATH_CHECK( m_v0FinderTool->performSearch(h_V0.ptr(),
108 h_Ks.ptr(),
109 h_La.ptr(),
110 h_Lb.ptr(),
111 primaryVertex, importedVxContainer.cptr(), ctx));
112
113 ATH_MSG_DEBUG("Reco_V0Finder v0Container->size() " << v0Container->size());
114 ATH_MSG_DEBUG("Reco_V0Finder ksContainer->size() " << ksContainer->size());
115 ATH_MSG_DEBUG("Reco_V0Finder laContainer->size() " << laContainer->size());
116 ATH_MSG_DEBUG("Reco_V0Finder lbContainer->size() " << lbContainer->size());
117
118
119 ATH_CHECK(m_v0DecoTool->decorateV0(h_V0.ptr(), ctx));
120 ATH_CHECK(m_v0DecoTool->decorateks(h_Ks.ptr() ,ctx));
121 ATH_CHECK(m_v0DecoTool->decoratela(h_La.ptr(), ctx));
122 ATH_CHECK(m_v0DecoTool->decoratelb(h_Lb.ptr(), ctx));
123 }
124
125 return StatusCode::SUCCESS;
126 }
127}
128
129
130
131
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
size_type size() const noexcept
Returns the number of elements in the collection.
SG::WriteHandleKey< xAOD::VertexContainer > m_v0Key
std::vector< std::string > m_CollectionsToCheck
ToolHandle< InDet::V0MainDecorator > m_v0DecoTool
ToolHandle< InDet::InDetV0FinderTool > m_v0FinderTool
SG::WriteHandleKey< xAOD::VertexContainer > m_ksKey
SG::WriteHandleKey< xAOD::VertexContainer > m_lbKey
virtual StatusCode addBranches(const EventContext &ctx) const override
SG::ReadHandleKey< xAOD::VertexContainer > m_vertexKey
Reco_V0Finder(const std::string &t, const std::string &n, const IInterface *p)
SG::WriteHandleKey< xAOD::VertexContainer > m_laKey
virtual bool isValid() override final
Can the handle be successfully dereferenced?
const_pointer_type cptr()
Dereference the pointer.
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
pointer_type ptr()
Dereference the pointer.
THE reconstruction tool.
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.