hpstr
The Heavy Photon Search Toolkit for Reconstruction (hpstr) provides an interface to physics data from the HPS experiment saved in the LCIO format and converts it into an ROOT based format.
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
VtxHistoProcessor.cxx
Go to the documentation of this file.
1
7#include "VtxHistoProcessor.h"
8
9VtxHistoProcessor::VtxHistoProcessor(const std::string& name, Process& process)
10 : Processor(name, process) {
11 }
12
15
17
18 std::cout << "Configuring VtxHistoProcessor" << std::endl;
19 try
20 {
21 debug_ = parameters.getInteger("debug");
22 rebin_ = parameters.getInteger("rebin");
23 selections_ = parameters.getVString("selections");
24 projections_ = parameters.getVString("projections");
25 }
26 catch (std::runtime_error& error)
27 {
28 std::cout << error.what() << std::endl;
29 }
30}
31
32void VtxHistoProcessor::initialize(std::string inFilename, std::string outFilename) {
33 // Init Files
34 inF_ = new TFile(inFilename.c_str());
35 outF_ = new TFile(outFilename.c_str(),"RECREATE");
36
37 //To save the projections
39
40 //Get the vtx vs InvM and p
41
42 for (auto selection_ : selections_) {
43 for (auto projection_ : projections_) {
44 if (debug_) {
45 std::cout<<"Getting... " <<(std::string("/")+selection_+"/"+selection_+"_"+projection_).c_str() << std::endl;
46 }
47
48 TH2F* proj2d = (TH2F*) inF_->Get((std::string("/")+selection_+"/"+selection_+"_"+projection_).c_str());
49
50 if (proj2d) {
51 _histos2d[selection_+"_"+projection_] = (TH2F*) inF_->Get((std::string("/")+selection_+"/"+selection_+"_"+projection_).c_str());
52 _histos2d[selection_+"_"+projection_]->SetDirectory(0);
53 _histos2d[selection_+"_"+projection_]->RebinX(rebin_);
54 }
55 else
56 std::cout<<"histo::"<<(std::string("/")+selection_+"/"+selection_+"_"+projection_).c_str()<<" does not exists."<<std::endl;
57 }//projection
58 }//selections
59}
60
62
63 for (it2d_ it = _histos2d.begin(); it!=_histos2d.end();++it) {
64 _histos1d[it->first+"_mu"] = new TH1F(
65 (it->first+"_mu").c_str(),
66 (it->first+"_mu").c_str(),
67 it->second->GetXaxis()->GetNbins(),
68 it->second->GetXaxis()->GetXmin(),
69 it->second->GetXaxis()->GetXmax());
70
71 _histos1d[it->first+"_sigma"] = new TH1F(
72 (it->first+"_sigma").c_str(),
73 (it->first+"_sigma").c_str(),
74 it->second->GetXaxis()->GetNbins(),
75 it->second->GetXaxis()->GetXmin(),
76 it->second->GetXaxis()->GetXmax());
77
78 _histos1d[it->first+"_mu"]->Sumw2();
79 _histos1d[it->first+"_sigma"]->Sumw2();
80
81 std::cout<<"Fitting::"<<it->first<<std::endl;
82 HistogramHelpers::profileYwithIterativeGaussFit(it->second,_histos1d[it->first+"_mu"],_histos1d[it->first+"_sigma"],1,0);
83 }
84 return true;
85}
86
88
89 outF_->cd();
90
91 for (it1d_ it = _histos1d.begin(); it!=_histos1d.end(); ++ it) {
92 it->second->Write();
93 }
94
95 for (it2d_ it = _histos2d.begin(); it!=_histos2d.end(); ++ it) {
96 it->second->Write();
97 }
98
99 outF_->Close();
100 inF_->Close();
101
102 //Close the projection file
104
105 delete inF_;
106}
107
#define DECLARE_PROCESSOR(CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
Definition Processor.h:139
description
Base class for all event processing components.
Definition Processor.h:34
TFile * outF_
Definition Processor.h:125
Insert description here. more details.
int debug_
Debug Level.
VtxHistoProcessor(const std::string &name, Process &process)
Constructor.
virtual void configure(const ParameterSet &parameters)
description
virtual void finalize()
description
std::map< std::string, TH2F * >::iterator it2d_
std::map< std::string, TH2F * > _histos2d
Map storing the 2D histograms.
std::map< std::string, TH1F * >::iterator it1d_
virtual void initialize(std::string inFilename, std::string outFilename)
description
std::map< std::string, TH1F * > _histos1d
Map storing the 1D biases and resolutions.
TFile * inF_
description
int rebin_
Rebin factor.
std::vector< std::string > projections_
2D histos to project
virtual bool process()
description
std::vector< std::string > selections_
Selection folder.
void profileYwithIterativeGaussFit(TH2 *hist, TH1 *mu_graph, TH1 *sigma_graph, int num_bins=1, int m_PrintLevel=0)
description
void OpenProjectionFile()
description
void CloseProjectionFile()
description