The Gaudi Framework  master (181af51f)
Loading...
Searching...
No Matches
THistRead.cpp
Go to the documentation of this file.
1/***********************************************************************************\
2* (c) Copyright 1998-2024 CERN for the benefit of the LHCb and ATLAS collaborations *
3* *
4* This software is distributed under the terms of the Apache version 2 licence, *
5* copied verbatim in the file "LICENSE". *
6* *
7* In applying this licence, CERN does not waive the privileges and immunities *
8* granted to it by virtue of its status as an Intergovernmental Organization *
9* or submit itself to any jurisdiction. *
10\***********************************************************************************/
11// Include files
12#include "THistRead.h"
13
17#include <math.h>
18
19#include <TDirectory.h>
20#include <TError.h>
21#include <TFile.h>
22#include <TH1F.h>
23#include <TH2F.h>
24#include <TH3F.h>
25#include <TKey.h>
26#include <TTree.h>
27
29
30//------------------------------------------------------------------------------
31THistRead::THistRead( const std::string& name, ISvcLocator* pSvcLocator )
32 : Algorithm( name, pSvcLocator )
33//------------------------------------------------------------------------------
34{}
35
36//------------------------------------------------------------------------------
38//------------------------------------------------------------------------------
39{
40 if ( m_ths.retrieve().isFailure() ) {
41 error() << "Couldn't get THistSvc" << endmsg;
43 }
44
45 // stream read1, 1D in "/xxx"
46 StatusCode sc1a = m_ths->regHist( "/read1/xxx/1Dgauss" );
47 TH1* h1( nullptr );
48 StatusCode sc1b = m_ths->getHist( "/read1/xxx/1Dgauss", h1 );
49 if ( sc1a.isFailure() || sc1b.isFailure() || h1 == nullptr ) {
50 error() << "Couldn't read gauss1d" << endmsg;
51 } else {
52 info() << h1->GetName() << ": " << h1->GetEntries() << endmsg;
53 }
54
55 // stream read2, 2D tree in "/"
56 StatusCode sc2a = m_ths->regHist( "/read2/2Dgauss" );
57 TH2* h2( nullptr );
58 StatusCode sc2b = m_ths->getHist( "/read2/2Dgauss", h2 );
59 if ( sc2a.isFailure() || sc2b.isFailure() || h2 == nullptr ) {
60 error() << "Couldn't read 2Dgauss" << endmsg;
61 } else {
62 info() << h2->GetName() << ": " << h2->GetEntries() << endmsg;
63 }
64
65 // 3D tree in "/"
66 StatusCode sc3a = m_ths->regHist( "/read2/3Dgauss" );
67 TH3* h3( nullptr );
68 StatusCode sc3b = m_ths->getHist( "/read2/3Dgauss", h3 );
69 if ( sc3a.isFailure() || sc3b.isFailure() || h3 == nullptr ) {
70 error() << "Couldn't read 3Dgauss" << endmsg;
71 } else {
72 info() << h3->GetName() << ": " << h3->GetEntries() << endmsg;
73 }
74
75 // Profile in "/"
76 StatusCode sc4a = m_ths->regHist( "/read2/profile" );
77 TH1* tp( nullptr );
78 StatusCode sc4b = m_ths->getHist( "/read2/profile", tp );
79 if ( sc4a.isFailure() || sc4b.isFailure() || tp == nullptr ) {
80 error() << "Couldn't read profile" << endmsg;
81 } else {
82 info() << tp->GetName() << ": " << tp->GetEntries() << endmsg;
83 }
84
85 // Tree with branches in "/trees/stuff"
86 StatusCode sc5a = m_ths->regTree( "/read2/trees/stuff/treename" );
87 TTree* tr( nullptr );
88 StatusCode sc5b = m_ths->getTree( "/read2/trees/stuff/treename", tr );
89 if ( sc5a.isFailure() || sc5b.isFailure() || tr == nullptr ) {
90 error() << "Couldn't read tree" << endmsg;
91 } else {
92 info() << tr->GetName() << ": " << tr->GetEntries() << endmsg;
93 }
95}
96
97//------------------------------------------------------------------------------
99//------------------------------------------------------------------------------
100{
101 return StatusCode::SUCCESS;
102}
103
104//------------------------------------------------------------------------------
106//------------------------------------------------------------------------------
107{
108 debug() << "Finalizing..." << endmsg;
109 return StatusCode::SUCCESS;
110}
MsgStream & endmsg(MsgStream &s)
MsgStream Modifier: endmsg. Calls the output method of the MsgStream.
Definition MsgStream.h:198
#define DECLARE_COMPONENT(type)
MsgStream & error() const
shortcut for the method msgStream(MSG::ERROR)
MsgStream & debug() const
shortcut for the method msgStream(MSG::DEBUG)
MsgStream & info() const
shortcut for the method msgStream(MSG::INFO)
Algorithm(std::string name, ISvcLocator *svcloc, std::string version=PACKAGE_VERSION)
Constructor.
Definition Algorithm.h:98
const std::string & name() const override
The identifying name of the algorithm object.
The ISvcLocator is the interface implemented by the Service Factory in the Application Manager to loc...
Definition ISvcLocator.h:42
This class is used for returning status codes from appropriate routines.
Definition StatusCode.h:64
bool isFailure() const
Definition StatusCode.h:129
constexpr static const auto SUCCESS
Definition StatusCode.h:99
constexpr static const auto FAILURE
Definition StatusCode.h:100
StatusCode initialize() override
Definition THistRead.cpp:37
StatusCode finalize() override
ServiceHandle< ITHistSvc > m_ths
Definition THistRead.h:33
THistRead(const std::string &name, ISvcLocator *pSvcLocator)
Definition THistRead.cpp:31
StatusCode execute() override
Definition THistRead.cpp:98
STL namespace.