wavecatcher-analysis
Loading...
Searching...
No Matches
CosmicsBox.cc
Go to the documentation of this file.
1#include "CosmicsBox.h"
2
25void CosmicsBox::Print_Phi_ew(vector<int> phi_chx, vector<float> ly_C0, vector<int> SiPMchannels, float windowmin, float windowmax, float maxfrom, float maxto, int nbins, bool corr, bool periodic) {
26 checkData();
27
28 // plotting uncorrected/corrected phi_ew - spectra ; first: map the phi_i to the channels, like Alex did, but with new channel positions
29 // match channel number to channel index (still very specific)
30 int sipmnum = SiPMchannels.size();
31 vector<int> ch_index; for (int i = 0; i < sipmnum; i++) ch_index.push_back(0); // initialize the ch_index vector
32 for (int i = 0; i < static_cast<int>(active_channels.size()); i++) {
33 for (int j = 0; j < sipmnum; j++) {
34 if (active_channels[i] == SiPMchannels[j]) ch_index[j] = i;
35 }
36 }
37
38 // compute correction factor
39 vector<float> ly_corr;
40 if (corr) {
41 float ly_av = 0;
42 // average lightyield of all channels * sipmnum
43 for (int i = 0; i < sipmnum; i++) ly_av += ly_C0[i];
44 // divide ly of one channel by (1/sipmnum)*sipmnum*average ly of all channels
45 for (int i = 0; i < sipmnum; i++) ly_corr.push_back(sipmnum * ly_C0[i] / ly_av);
46 }
47 else for (int i = 0; i < sipmnum; i++) ly_corr.push_back(1); // no correction
48
49 // print correction factors
50 for (int i = 0; i < sipmnum; i++) cout << "Correction factor for channel " << SiPMchannels[i] << ":" << ly_corr[i] << endl;
51
52 // initialize canvas + histograms
53 gStyle->SetOptStat("emr"); gStyle->SetOptFit(1111); //draws a box with some histogram parameters
54 TString his_name(Form("#phi_ew-spectrum_from_ch%d_to_ch%d", SiPMchannels.front(), SiPMchannels.back()));
55 TCanvas* hisc = new TCanvas(his_name, his_name, 600, 400);
56
57 double min_angle = -180, max_angle = 180; // 360 deg plot range
58 if (periodic) { min_angle = -540; max_angle = 540; } // 1080 deg plot range
59 TH1* his = new TH1F(his_name, his_name, nbins, min_angle, max_angle);
60
61 // loop through all events and compute phi_ew
62 float lightyield, anglevaluex, anglevaluey, phi_ew = 0;
63 for (int i = 0; i < nevents; i++) {
64 if (!skip_event[i]) {
65 for (int j = 0; j < sipmnum; j++) { //loop through all SiPM-channels
66 TH1F* hisly = Getwf(i * nchannels + ch_index[j]);
67 lightyield = GetPeakIntegral(hisly, windowmin, windowmax, maxfrom, maxto, 0); //lightyield as the integral around maximum
68 anglevaluex += cos(phi_chx[j] * TMath::Pi() / 180) * lightyield / ly_corr[j]; //x part of vectorial addition
69 anglevaluey += sin(phi_chx[j] * TMath::Pi() / 180) * lightyield / ly_corr[j]; //y part of vectorial addition
70 }
71
72 phi_ew = atan2(anglevaluey, anglevaluex) * 180 / TMath::Pi(); //vectorial addition for phi_ew --> fill in histo
73 his->Fill(phi_ew);
74 if (periodic) {
75 his->Fill(phi_ew + 360);
76 his->Fill(phi_ew - 360);
77 }
78 anglevaluex = 0, anglevaluey = 0; //reset the x and y parts
79 }
80 }
81
82 //make histogram fancy + printing
83 his->GetXaxis()->SetTitle("#phi_{ew} [#circ]");
84 his->GetYaxis()->SetTitle("#Entries");
85 his->Draw();
86
87 if (periodic) {
89 int n_par = 4;
90 TF1* f = new TF1("fitf", fitf, min_angle, max_angle, n_par);
91 f->SetLineColor(2);
92 f->SetNpx(1000);
93
94 double max = his->GetMaximum();
95 double min = his->GetMinimum();
96 his->GetXaxis()->SetRangeUser(-180, 180); //only look at one period for starting values
97 double phi_est = his->GetXaxis()->GetBinCenter(his->GetMaximumBin());
98 his->GetXaxis()->SetRangeUser(0, 0); //reset scale
99
100 f->SetParName(0, "A"); f->SetParameter(0, max - min); f->SetParLimits(0, 1, 1e9);
101 f->SetParName(1, "#bar{#phi_{ew}}"); f->SetParameter(1, phi_est); f->SetParLimits(1, -180, 180);
102 f->SetParName(2, "#sigma"); f->SetParameter(2, 40); f->SetParLimits(2, 5, 360);
103 f->SetParName(3, "offset"); f->SetParameter(3, min); f->SetParLimits(3, TMath::Min(min - 1, 0.1), 1e9);
104
105 TFitResultPtr fresults = his->Fit(f, "LRS");
106 fit_results.push_back(fresults);
107 }
108
109 hisc->Update();
110 root_out->WriteObject(hisc, "Phi_ew");
111 root_out->WriteObject(his, "Phi_ew_his");
112
113 string pdf_name = "phi_ew_spectrum";
114 pdf_name += corr ? "_corr.pdf" : "_uncorr.pdf";
115 hisc->SaveAs(pdf_name.c_str()); //write the histogram to a .pdf-file
116}
void Print_Phi_ew(vector< int >, vector< float >, vector< int >, float=1, float=1, float=100, float=140, int=400, bool=true, bool=false)
Angular distribution of passing particles in the cosmics setup.
Definition CosmicsBox.cc:25
Periodic gauss fit function for angular reconstruction of SiPM array signal (SHiP SBT WOM tube)
TH1F * Getwf(int)
Helper that returns the waveform histogram for a certain waveform number number.
Definition ReadRun.cc:2547
int nevents
Number of triggered events in data.
Definition ReadRun.h:263
vector< bool > skip_event
Stores the event numbers which should be skipped in the analysis.
Definition ReadRun.h:311
vector< int > active_channels
Stores the numbers of the active channels.
Definition ReadRun.h:290
void checkData() const
Primitive check to see if data has been loaded.
Definition ReadRun.h:113
float GetPeakIntegral(TH1F *, float, float, float, float, int=0)
Calculate the integral around a peak with several options explained in GetIntWindow().
Definition ReadRun.cc:1588
TFile * root_out
Stores results of analysis.
Definition ReadRun.h:338
int nchannels
Number of active channels in data.
Definition ReadRun.h:265
vector< TFitResultPtr > fit_results
Stores the fit results of PrintChargeSpectrum() for all channels and all function calls in ascending ...
Definition ReadRun.h:303