Belle II Software development
HadronPrep.h
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#pragma once
10
11#include <map>
12#include <vector>
13#include <iostream>
14#include <fstream>
15
16#include <TH1D.h>
17#include <TH2F.h>
18#include <TF1.h>
19#include <TCanvas.h>
20#include <TLine.h>
21#include <TLegend.h>
22#include<TGraphErrors.h>
23
24#include <calibration/CalibrationAlgorithm.h>
25#include <framework/database/DBObjPtr.h>
26#include <cdc/utilities/CDCDedxHadSat.h>
27#include <framework/gearbox/Const.h>
28
29#include <cdc/calibration/CDCdEdx/HadronBgPrep.h>
30
31namespace Belle2 {
40 class HadronPrep {
41
42 public:
43
47 HadronPrep();
48
52 virtual ~HadronPrep() {};
53
57 HadronPrep(int bgbins, double upperbg, double lowerbg, int cosbins, double uppercos, double lowercos, double cut);
58
62 void prepareSample(std::shared_ptr<TTree> hadron, TFile*& outfile, const std::string& suffix,
63 const std::string& pdg, bool ismakePlots, bool correct);
64
68 void defineHisto(std::vector<TH1F*>& htemp, const std::string& title, const std::string& pdg);
69
73 void plotDist(std::map<int, std::vector<TH1F*>>& hist, const std::string& sname, const std::string& pdg);
74
78 void plotDist(std::vector<TH2F*>& hist, const std::string& sname, const std::string& pdg);
79
83 void setPars(TFile*& outfile, std::map<int, std::vector<TH1F*>>& hdedx_bgcosth, const std::string& pdg);
84
88 void plotGraph(const std::string& sname, const std::string& pdg);
89
93 void clear();
94
95 private:
96
97 std::map<int, std::vector<double>> m_sumcos;
98 std::map<int, std::vector<double>> m_sumbg;
99 std::map<int, std::vector<double>> m_means;
100 std::map<int, std::vector<double>> m_errors;
102 std::map<int, std::vector<int>> m_sumsize;
104 double m_dedxmax = 0.0;
105 double m_dedxmin = 99999999.0;
108 double m_bgMin;
109 double m_bgMax;
112 double m_cosMin;
113 double m_cosMax;
115 double m_cut;
117 };
119} // namespace Belle2
Class to prepare sample for hadron saturation calibration.
Definition: HadronPrep.h:40
void prepareSample(std::shared_ptr< TTree > hadron, TFile *&outfile, const std::string &suffix, const std::string &pdg, bool ismakePlots, bool correct)
function to prepare sample for monitoring plots, bg curve fitting and sigma vs ionz fitting
Definition: HadronPrep.cc:34
double m_dedxmax
variables to set maximum dedx mean
Definition: HadronPrep.h:104
int m_cosBins
bins for dedx histogram
Definition: HadronPrep.h:111
std::map< int, std::vector< double > > m_errors
error variable
Definition: HadronPrep.h:100
double m_cosMax
max range of dedx
Definition: HadronPrep.h:113
std::map< int, std::vector< int > > m_sumsize
variables for size
Definition: HadronPrep.h:102
HadronPrep()
Constructor: Sets the description, the properties and the parameters of the algorithm.
Definition: HadronPrep.cc:12
double m_dedxmin
variables to set minimum dedx mean
Definition: HadronPrep.h:105
void plotDist(std::map< int, std::vector< TH1F * > > &hist, const std::string &sname, const std::string &pdg)
function to plot the map of histograms
Definition: HadronPrep.cc:177
int m_bgBins
bins for dedx histogram
Definition: HadronPrep.h:107
virtual ~HadronPrep()
Destructor.
Definition: HadronPrep.h:52
void setPars(TFile *&outfile, std::map< int, std::vector< TH1F * > > &hdedx_bgcosth, const std::string &pdg)
function to fill the parameters like mean and reso in the tree
Definition: HadronPrep.cc:247
double m_bgMax
max range of dedx
Definition: HadronPrep.h:109
void plotGraph(const std::string &sname, const std::string &pdg)
function to make graph dedx vs cos in different bg bins
Definition: HadronPrep.cc:304
void defineHisto(std::vector< TH1F * > &htemp, const std::string &title, const std::string &pdg)
function to define histograms
Definition: HadronPrep.cc:157
double m_bgMin
min range of dedx
Definition: HadronPrep.h:108
double m_cosMin
min range of dedx
Definition: HadronPrep.h:112
double m_cut
cut to clean protons
Definition: HadronPrep.h:115
void clear()
function to clear the variables
Definition: HadronPrep.cc:367
std::map< int, std::vector< double > > m_means
mean variable
Definition: HadronPrep.h:99
std::map< int, std::vector< double > > m_sumbg
variables to add bg values
Definition: HadronPrep.h:98
std::map< int, std::vector< double > > m_sumcos
variables to add cos values
Definition: HadronPrep.h:97
Abstract base class for different kinds of events.