wavecatcher-analysis
Loading...
Searching...
No Matches
PMT.cc
Go to the documentation of this file.
1#include "PMT.h"
2
9void PMT::PrintChargeSpectrumPMT(float windowlow, float windowhi, float start, float end, float rangestart, float rangeend, int nbins) {
10
11 PrintChargeSpectrumPMT_cnt++;
12
13 string ctitle("charge spectra PMT" + to_string(PrintChargeSpectrumPMT_cnt));
14 TCanvas* chargec = new TCanvas(ctitle.c_str(), ctitle.c_str(), 600, 400);
16
17 int current_canvas = 0;
18
19 for (int i = 0; i < nchannels; i++) {
20 if (PlotChannel(i)) {
21 current_canvas++;
22
23 TH1F* his;
24 his = ChargeSpectrum(i, windowlow, windowhi, start, end, rangestart, rangeend, nbins);
25 if (windowlow + windowhi > 0.) his->GetXaxis()->SetTitle("integral in mV#timesns");
26 else his->GetXaxis()->SetTitle("amplitude in mV");
27 chargec->cd(current_canvas);
28
29 his->GetYaxis()->SetTitle("#Entries");
30 his->GetXaxis()->SetTitle("integral in mV#timesns");
31 his->Draw();
32 stringstream allname; allname << his->GetEntries() << " entries";
33 his->SetTitle(allname.str().c_str());
34
35 TString name(Form("ChargeSpectrumPMT channel_%02d", active_channels[i]));
36 root_out->WriteObject(his, name.Data());
37
38 // fit sum of two gaussians
39 // estimate starting values of fit by fitting only one gauss
40 auto gauss = new TF1("gauss", "gaus", rangestart, rangeend);
41 TFitResultPtr fres_est = his->Fit(gauss, "QRSL");
42 // now use these fit results to fit the sum of two gauss
43 auto two_gauss = new TF1("two gaussians", "gaus(0)+gaus(3)", rangestart, rangeend); two_gauss->SetTitle("Sum of two gauss");
44 two_gauss->SetParameters(fres_est->Parameter(0) * .95, fres_est->Parameter(1) * .95, fres_est->Parameter(2) * .95, fres_est->Parameter(0) * .3, fres_est->Parameter(1) * 1.05, fres_est->Parameter(2) * .85); // factors are pretty much random
45 two_gauss->SetParName(0, "A_{pedestal}");
46 two_gauss->SetParName(1, "#mu_{pedestal}");
47 two_gauss->SetParName(2, "#sigma_{pedestal}"); two_gauss->SetParLimits(2, 1e-9, 1e3);
48 two_gauss->SetParName(3, "A_{SPE}");
49 two_gauss->SetParName(4, "#mu_{SPE}");
50 two_gauss->SetParName(5, "#sigma_{SPE}"); two_gauss->SetParLimits(5, 1e-9, 1e3);
51
52 if (!PrintChargeSpectrumPMT_pars.empty()) {
53 for (int j = 0; j < static_cast<int>(PrintChargeSpectrumPMT_pars.size()); j++) two_gauss->SetParameter(j, PrintChargeSpectrumPMT_pars[j]);
54 }
55
56 //two_gauss->SetLineColor(4);
57 TFitResultPtr fresults = his->Fit(two_gauss, "RSL");
58 fit_results.push_back(fresults);
59
60 two_gauss->Draw("same");
61
62 auto pedestal = new TF1("pedestal", "gaus", rangestart, rangeend); pedestal->SetTitle("pedestal");
63 pedestal->SetParameters(fresults->Parameter(0), fresults->Parameter(1), fresults->Parameter(2));
64 pedestal->SetLineColor(3);
65 pedestal->Draw("same");
66
67 auto pepeak = new TF1("pepeak", "gaus", rangestart, rangeend); pepeak->SetTitle("pepeak");
68 pepeak->SetParameters(fresults->Parameter(3), fresults->Parameter(4), fresults->Parameter(5));
69 pepeak->SetLineColor(4);
70 pepeak->Draw("same");
71
72 gPad->BuildLegend();
73 }
74 }
75
76 chargec->Update();
77 root_out->WriteObject(chargec, "ChargeSpectraPMT");
78}
79
87void PMT::PrintChargeSpectrumPMTthreshold(float windowlow, float windowhi, float rangestart, float rangeend, int nbins, double threshold, bool calculate_SiPM_DCR) {
88
89 PrintChargeSpectrumPMTthreshold_cnt++;
90
91 gStyle->SetOptStat(0); // 11 is title + entries
92
93 // show fraction of events above 0.5 pe charge = pedestal + gain/2
94 // dark count rate for SiPMs (currently only automated for fit function Fitf)
95 // need to call the SiPM fit function before this one for this functionality
96 bool use_fit_result_for_threshold = false;
97 if (threshold == 999) {
98 calculate_SiPM_DCR = true;
99 use_fit_result_for_threshold = true;
100 }
101
102 string unit(" mV");
103 string title("amplitude in mV"); // amplitude spectrum not good for fitting, will be biased
104 if (windowlow + windowhi > 0.) {
105 unit = " mV#timesns";
106 title = "integral in mV#timesns";
107 }
108 string ctitle("charge spectra PMT threshold" + to_string(PrintChargeSpectrumPMTthreshold_cnt));
109 TCanvas* chargec = new TCanvas(ctitle.c_str(), ctitle.c_str(), 600, 400);
111
112 int current_canvas = 0;
113 double threshold_bin_center = 0;
114
115 for (int i = 0; i < nchannels; i++) {
116 if (PlotChannel(i)) {
117
118 if (use_fit_result_for_threshold) threshold = fit_results[current_canvas]->Parameter(6) + fit_results[current_canvas]->Parameter(5) / 2.;
119
120 chargec->cd(++current_canvas);
121
122 TH1F* his;
123 his = ChargeSpectrum(i, windowlow, windowhi, rangestart, rangeend, rangestart, rangeend, nbins);
124 his->GetXaxis()->SetTitle(title.c_str());
125 his->GetYaxis()->SetTitle("#Entries");
126 his->Draw();
127 stringstream allname; allname << his->GetEntries() << " entries";
128 his->SetTitle(allname.str().c_str());
129
130 auto his_lo = (TH1F*)his->Clone();
131 his_lo->GetXaxis()->SetRange(his_lo->GetXaxis()->FindBin(rangestart), his_lo->GetXaxis()->FindBin(threshold));
132 his_lo->SetLineColor(2);
133 his_lo->SetFillColor(2);
134 his_lo->Draw("LF2 same");
135
136 stringstream lonamefrac;
137 stringstream lonamerate;
138 lonamefrac << 100. * his->Integral(his->GetXaxis()->FindBin(rangestart), his->GetXaxis()->FindBin(threshold)) / his->GetEntries() << "% <= " << threshold << unit;
139 lonamerate << "<0.5 pe=" << threshold << unit << " -> " << his->Integral(his->GetXaxis()->FindBin(rangestart), his->GetXaxis()->FindBin(threshold)) / his->GetEntries() / (1.e-3 * (rangeend - rangestart)) << " MHz";
140 cout << "\n" << lonamerate.str().c_str() << endl;
141 cout << "\n" << lonamefrac.str().c_str() << endl;
142 his_lo->SetTitle(lonamerate.str().c_str());
143 if (!calculate_SiPM_DCR) his_lo->SetTitle(lonamefrac.str().c_str());
144
145 auto his_hi = (TH1F*)his->Clone();
146 his_hi->GetXaxis()->SetRange(his_hi->GetXaxis()->FindBin(threshold), his_lo->GetXaxis()->FindBin(rangeend));
147 his_hi->SetLineColor(1);
148 his_hi->SetFillColor(1);
149 his_hi->Draw("LF2 same");
150
151 stringstream hinamefrac;
152 stringstream hinamerate;
153 hinamefrac << 100. * his->Integral(his->GetXaxis()->FindBin(threshold) + 1, his->GetXaxis()->FindBin(rangeend)) / his->GetEntries() << "% > " << threshold << unit;
154 hinamerate << ">0.5 pe=" << threshold << unit << " -> " << his->Integral(his->GetXaxis()->FindBin(threshold) + 1, his->GetXaxis()->FindBin(rangeend)) / his->GetEntries() / (1.e-3 * (rangeend - rangestart)) << " MHz";
155 cout << "\n" << hinamerate.str().c_str() << endl;
156 cout << "\n" << hinamefrac.str().c_str() << endl;
157 his_hi->SetTitle(hinamerate.str().c_str());
158 if (!calculate_SiPM_DCR) his_hi->SetTitle(hinamefrac.str().c_str());
159
160 gPad->BuildLegend();
161
162 threshold_bin_center = his->GetXaxis()->GetBinCenter(his->GetXaxis()->FindBin(threshold) + 1);
163 cout << "\n PMT charge spectrum is counting events above threshold from bin center >= " << threshold_bin_center << unit << " for a threshold setting of " << threshold << unit << "\n\n";
164 }
165 }
166
167 chargec->Update();
168 root_out->WriteObject(chargec, "ChargeSpectraPMTthreshold");
169}
static void SplitCanvas(TCanvas *&, vector< int >, vector< int >)
Helper to split canvas according to the number of channels to be plotted.
Definition Helpers.cc:97
vector< float > PrintChargeSpectrumPMT_pars
Starting values of the fit parameters for PrintChargeSpectrumPMT()
Definition PMT.h:24
void PrintChargeSpectrumPMTthreshold(float=0, float=0, float=0, float=300, int=750, double=4, bool=false)
Print "charge" spectrum with highlight a threshold.
Definition PMT.cc:87
void PrintChargeSpectrumPMT(float, float, float=0, float=300, float=-50, float=600, int=750)
"Charge" spectrum optimized for PMT signals
Definition PMT.cc:9
TH1F * ChargeSpectrum(int, float, float, float=0, float=300, float=-50, float=600, int=750)
Histogram of the "charge" spectrum for one channel.
Definition ReadRun.cc:1807
vector< int > active_channels
Stores the numbers of the active channels.
Definition ReadRun.h:290
TFile * root_out
Stores results of analysis.
Definition ReadRun.h:338
bool PlotChannel(int)
Check if a channel index should be plotted according to ReadRun::plot_active_channels.
Definition ReadRun.cc:2643
int nchannels
Number of active channels in data.
Definition ReadRun.h:265
vector< int > plot_active_channels
Stores the numbers of the active channels which should be plotted.
Definition ReadRun.h:295
vector< TFitResultPtr > fit_results
Stores the fit results of PrintChargeSpectrum() for all channels and all function calls in ascending ...
Definition ReadRun.h:303