ATLAS Offline Software
Loading...
Searching...
No Matches
InDetTrackSmearingToolTester.cxx
Go to the documentation of this file.
1 /*
2 Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3*/
4
5// Framework include(s):
7
8// EDM include(s):
10
11// Local include(s):
14#include <TH1.h>
15
16
17namespace InDet {
18 InDetTrackSmearingToolTester::InDetTrackSmearingToolTester( const std::string& name, ISvcLocator* svcLoc )
19 : AthHistogramAlgorithm( name, svcLoc ),
20 m_smearTool( "InDet::InDetTrackSystematicsTools/InDetTrackSmearingTool", this ){
21 declareProperty( "TrackIP", m_Track_IP = "InDetTrackParticles" );
22 declareProperty( "SystematicEffects", m_systematicsNames );
23 declareProperty( "InDetTrackSmearingTool", m_smearTool );
24
25 }
26
28
29 // Greet the user:
30 ATH_MSG_INFO( "Initialising" );
31 ATH_MSG_DEBUG( "InDetTrackSmearingTool = " << m_smearTool );
32
33 // Retrieve the tools:
34 ATH_CHECK( m_smearTool.retrieve() );
35
36 for (const auto& name : m_systematicsNames) {
37 for (const auto& systpair : InDet::TrackSystematicMap) {
38 if (name == systpair.second.name()) {
39 m_systActive.insert(systpair.second);
40 }
41 }
42 }
43 auto systCode = m_smearTool->applySystematicVariation( m_systActive );
44 if (systCode != StatusCode::SUCCESS) {
45 ATH_MSG_ERROR( "Failure to apply systematic variation." );
46 return StatusCode::FAILURE;
47 }
48
49 ATH_CHECK( book( TH1F("d0_B", "original d0", 100, -5.0, 5.0) ) );
50 ATH_CHECK( book( TH1F("z0_B", "original z0", 100, -200.0, 200.0) ) );
51 ATH_CHECK( book( TH1F("d0sm", "d0 after smearing", 100, -5.0, 5.0) ) );
52 ATH_CHECK( book( TH1F("z0sm", "z0 after smearing", 100, -200.0, 200.0) ) );
53 ATH_CHECK( book( TH1F("subtraction_d0", "subtraction_d0", 100, -0.10, 0.10) ) );
54 ATH_CHECK( book( TH1F("subtraction_z0", "subtraction_z0", 100,-0.50, 0.50) ) );
55
56 // Return gracefully:
57 return StatusCode::SUCCESS;
58 }
59
61
62 // Create a shallow container copy and then apply the smearingtool to impact parameters:
63 const xAOD::TrackParticleContainer *IDParticles = nullptr;
64 ATH_CHECK( evtStore()->retrieve( IDParticles , m_Track_IP ) );
65 std::pair< xAOD::TrackParticleContainer*, xAOD::ShallowAuxContainer* > IDParticles_shallowCopy = xAOD::shallowCopyContainer( *IDParticles );
66 for( xAOD::TrackParticle* track : *IDParticles_shallowCopy.first ) {
67 double d0_1=0.,d0_2=0.,z0_1=0.,z0_2=0.;
68 d0_1=track->d0();
69 z0_1=track->z0();
70 hist("d0_B")->Fill( d0_1 );
71 hist("z0_B")->Fill( z0_1 );
72 if (m_smearTool->applyCorrection(*track) == CP::CorrectionCode::Error) {
73 ATH_MSG_ERROR( "Could not apply correction." );
74 }
75 d0_2=track->d0();
76 z0_2=track->z0();
77 hist("d0sm")->Fill( d0_2 );
78 hist("z0sm")->Fill( z0_2 );
79 hist("subtraction_d0")->Fill( d0_2 - d0_1 );
80 hist("subtraction_z0")->Fill( z0_2 - z0_1 );
81 }
82 delete IDParticles_shallowCopy.first;
83 delete IDParticles_shallowCopy.second;
84
85 // Return gracefully:
86 return StatusCode::SUCCESS;
87
88 } // End of the execute()
89
90} // namespace InDet
#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)
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
AthHistogramAlgorithm(const std::string &name, ISvcLocator *pSvcLocator)
Constructor with parameters:
StatusCode book(const TH1 &hist, const std::string &tDir="", const std::string &stream="")
Simplify the booking and registering (into THistSvc) of histograms.
TH1 * hist(const std::string &histName, const std::string &tDir="", const std::string &stream="")
Simplify the retrieval of registered histograms of any type.
@ Error
Some error happened during the object correction.
virtual StatusCode execute()
Function executing the algorithm.
std::string m_Track_IP
StoreGate key for the track container to investigate//--->delete in future.
virtual StatusCode initialize()
Function initialising the algorithm.
InDetTrackSmearingToolTester(const std::string &name, ISvcLocator *svcLoc)
Regular Algorithm constructor.
ToolHandle< IInDetTrackSmearingTool > m_smearTool
Connection to the smearing tool.
Primary Vertex Finder.
std::pair< std::unique_ptr< T >, std::unique_ptr< ShallowAuxContainer > > shallowCopyContainer(const T &cont, const EventContext &ctx)
Function making a shallow copy of a constant container.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".