BmnRoot
Loading...
Searching...
No Matches
BmnHistToF700.cxx
Go to the documentation of this file.
1#include "BmnHistToF700.h"
2
3BmnHistToF700::BmnHistToF700(TString title, TString path)
4 : BmnHist()
5{
6 fTitle = title;
7 fName = title + "_cl";
8 const Int_t NTimeBins = 5000;
9 const Int_t TimeUpper = 100;
10 const Int_t TimeLower = -100;
11 const Int_t NAmpBins = 500;
12 const Int_t AmpUpper = 200;
13 fSelectedPlane = -1;
14 fSelectedStrip = -1;
15 fSelectedSide = -1;
16 TString name;
17 name = fTitle + "_Leading_Time";
18 histLeadingTime = new TH1D(name, name, NTimeBins, TimeLower, TimeUpper);
19 histLeadingTime->GetXaxis()->SetTitle("Time, ns");
20 histLeadingTime->GetYaxis()->SetTitle("Activations count");
21 name = fTitle + "_Leading_Time_Specific";
22 histLeadingTimeSpecific = new TH1D(name, name, NTimeBins, TimeLower, TimeUpper);
23 histLeadingTimeSpecific->GetXaxis()->SetTitle("Time, ns");
24 histLeadingTimeSpecific->GetYaxis()->SetTitle("Activations count");
25 name = fTitle + "_Amplitude";
26 histAmp = new TH1D(name, name, NAmpBins, 0, AmpUpper);
27 histAmp->GetXaxis()->SetTitle("Amplitude, ns");
28 histAmp->GetYaxis()->SetTitle("Activations count");
29 name = fTitle + "_Amplitude_Specific";
30 histAmpSpecific = new TH1D(name, name, NAmpBins, 0, AmpUpper);
31 histAmpSpecific->GetXaxis()->SetTitle("Amplitude, ns");
32 histAmpSpecific->GetYaxis()->SetTitle("Activations count");
33 name = fTitle + "_Strip";
34 histStrip = new TH1I(name, name, TOF2_MAX_CHAMBERS * TOF2_MAX_STRIPS_IN_CHAMBER, 0,
36 histStrip->GetXaxis()->SetTitle("Strip #");
37 histStrip->GetYaxis()->SetTitle("Activations count");
38 // name = fTitle + "_StripSimult";
39 // histStripSimult = new TH1I(name, name, TOF2_MAX_STRIPS_IN_CHAMBER, 0, TOF2_MAX_STRIPS_IN_CHAMBER);
40 // name = fTitle + "_State";
41 // histState = new TH2F(name, name, TOF2_MAX_STRIPS_IN_CHAMBER, 0, TOF2_MAX_STRIPS_IN_CHAMBER, 2, 0, 2);
42
43 histSimultaneous.SetDirectory(0);
44 histL->SetDirectory(0);
45 histR->SetDirectory(0);
46 fServer = NULL;
47 frecoTree = NULL;
48 Events = NULL;
49 name = fTitle + "CanvasTimes";
50 canTimes = new TCanvas(name, name, PAD_WIDTH * TOF_ROWS, PAD_HEIGHT * TOF_COLS);
51 canTimes->Divide(TOF_ROWS, TOF_COLS);
52 canTimesPads.resize(TOF_ROWS * TOF_COLS);
53 for (Int_t iPad = 0; iPad < TOF_ROWS * TOF_COLS; iPad++) {
54 PadInfo* p = new PadInfo();
55 canTimesPads[iPad] = p;
56 canTimes->GetPad(iPad + 1)->SetGrid();
57 }
58 canTimesPads[0]->current = histLeadingTime;
59 canTimesPads[1]->current = histAmp;
60 canTimesPads[2]->current = histLeadingTimeSpecific;
61 canTimesPads[3]->current = histAmpSpecific;
62 canTimesPads[4]->current = histStrip;
63 // canTimesPads[5]->current = histStripSimult;
64 for (size_t iPad = 0; iPad < canTimesPads.size(); iPad++) {
65 if (canTimesPads[iPad]->current) {
66 Names.push_back(canTimesPads[iPad]->current->GetName());
67 BmnHist::SetHistStyleTH1(canTimesPads[iPad]->current);
68 }
69 }
70 canTimesPads[4]->logy = true;
71}
72
74{
75 delete histL;
76 delete histR;
77 delete Events;
78 delete canTimes;
79 if (fDir)
80 return;
81 for (auto pad : canTimesPads)
82 delete pad;
83}
84
86{
87 const Double_t ch2ns = 0.02344;
88 TClonesArray* digits = fDigiArrays->tof700;
89 if (!digits)
90 digits = fDigiArrays->tof701;
91 if (!digits)
92 return;
93 histL->Reset();
94 histR->Reset();
95 histSimultaneous.Reset();
96 // histState->Reset();
97 Events->Clear();
98 for (Int_t digIndex = 0; digIndex < digits->GetEntriesFast(); digIndex++) {
99 BmnTof1Digit* td = (BmnTof1Digit*)digits->At(digIndex);
100 Int_t strip = td->GetStrip();
101 histLeadingTime->Fill(td->GetTime());
102 histAmp->Fill(td->GetAmplitude() * ch2ns);
103 histStrip->Fill(strip + td->GetPlane() * TOF2_MAX_STRIPS_IN_CHAMBER);
104 if (((td->GetPlane() == fSelectedPlane) || (fSelectedPlane < 0))
105 && ((td->GetStrip() == fSelectedStrip) || (fSelectedStrip < 0)))
106 {
107 histAmpSpecific->Fill(td->GetAmplitude() * ch2ns);
108 histLeadingTimeSpecific->Fill(td->GetTime());
109 }
110
111 // new ((*Events)[Events->GetEntriesFast()])
112 // BmnTof2Digit(td->GetPlane(), td->GetStrip(), td->GetTime(), td->GetAmplitude(),
113 // td->GetDiff());
114 // frecoTree->Fill();
115 }
116}
117
118void BmnHistToF700::Register(THttpServer* serv)
119{
120 fServer = serv;
121 fServer->Register("/", this);
122 TString path = "/" + fTitle + "/";
123 fServer->Register(path, canTimes);
124
125 TString cmd = "/" + fName + "/->SetRefRun(%arg1%)";
126 TString cmdTitle = path + "SetRefRun";
127 fServer->RegisterCommand(cmdTitle.Data(), cmd.Data(), "button;");
128 fServer->Restrict(cmdTitle.Data(), "visible=shift");
129 fServer->Restrict(cmdTitle.Data(), "allow=shift");
130 fServer->Restrict(cmdTitle.Data(), "deny=guest");
131 fServer->SetItemField(path.Data(), "_monitoring", "2000");
132 fServer->SetItemField(path.Data(), "_layout", "grid3x3");
133 cmdTitle = path + "ChangeSlection";
134 fServer->RegisterCommand(cmdTitle, TString("/") + fName.Data() + "/->SetSelection(%arg1%,%arg2%)", "button;");
135 fServer->Restrict(cmdTitle.Data(), "visible=shift");
136 fServer->Restrict(cmdTitle.Data(), "allow=shift");
137 fServer->Restrict(cmdTitle.Data(), "deny=guest");
138 cmdTitle = path + TString("Reset");
139 fServer->RegisterCommand(cmdTitle, TString("/") + fName.Data() + "/->Reset()", "button;");
140 fServer->Restrict(cmdTitle.Data(), "visible=shift");
141 fServer->Restrict(cmdTitle.Data(), "allow=shift");
142 fServer->Restrict(cmdTitle.Data(), "deny=guest");
143}
144
145void BmnHistToF700::SetDir(TFile* outFile, TTree* recoTree)
146{
147 frecoTree = recoTree;
148 fDir = NULL;
149 if (outFile != NULL)
150 fDir = outFile->mkdir(fTitle + "_hists");
151 histLeadingTime->SetDirectory(fDir);
152 histLeadingTimeSpecific->SetDirectory(fDir);
153 histAmp->SetDirectory(fDir);
154 histAmpSpecific->SetDirectory(fDir);
155 histStrip->SetDirectory(fDir);
156 // histStripSimult->SetDirectory(fDir);
157 // histState->SetDirectory(fDir);
158 if (Events != NULL)
159 delete Events;
160 Events = new TClonesArray("BmnTof2Digit");
161 if (frecoTree != NULL)
162 frecoTree->Branch(fTitle.Data(), &Events);
163}
164
165void BmnHistToF700::SetSelection(Int_t Plane, Int_t Strip)
166{
167 SetPlane(Plane);
168 SetStrip(Strip);
169 TString command;
170 if (fSelectedPlane >= 0)
171 command = Form("fPlane == %d", fSelectedPlane);
172 if (fSelectedStrip >= 0) {
173 if (command.Length() > 0)
174 command = command + " && ";
175 command = command + Form("fStrip == %d", fSelectedStrip);
176 }
177 if (fSelectedSide >= 0) {
178 if (command.Length() > 0)
179 command = command + " && ";
180 command = command + Form("fSide == %d", fSelectedSide);
181 }
182 histAmpSpecific->Reset();
183 histAmpSpecific->SetTitle("Amplitude For: " + command);
184 histLeadingTimeSpecific->Reset();
185 histLeadingTimeSpecific->SetTitle("Leading Time For: " + command);
186 if (frecoTree != NULL) {
187 TString direction = "fAmplitude>>" + TString(histAmpSpecific->GetName());
188 frecoTree->Draw(direction, command, "");
189 direction = "fTime>>" + TString(histLeadingTimeSpecific->GetName());
190 frecoTree->Draw(direction, command, "");
191 }
192}
193
195{
196 histLeadingTime->Reset();
197 histLeadingTimeSpecific->Reset();
198 histAmp->Reset();
199 histAmpSpecific->Reset();
200 histStrip->Reset();
201 // histStripSimult->Reset();
202 // histState->Reset();
203}
204
206{
207 BmnHist::DrawRef(canTimes, &canTimesPads);
208}
209
211{
212 if (refID != id) {
213 TString FileName = Form("bmn_run%04d_hist.root", id);
214 printf("SetRefRun: %s\n", FileName.Data());
215 refRunName = FileName;
216 refID = id;
217 BmnHist::LoadRefRun(refID, refPath + FileName, fTitle, canTimesPads, Names);
218 DrawBoth();
219 }
220 return kBMNSUCCESS;
221}
222
224{
225 for (auto pad : canTimesPads) {
226 if (pad->ref)
227 delete pad->ref;
228 pad->ref = NULL;
229 }
230 refID = 0;
231}
#define PAD_WIDTH
Definition BmnAdcQA.cxx:3
#define PAD_HEIGHT
Definition BmnAdcQA.cxx:4
BmnStatus
Definition BmnEnums.h:24
@ kBMNSUCCESS
Definition BmnEnums.h:25
BmnStatus SetRefRun(Int_t id)
void FillFromDigi(DigiArrays *fDigiArrays)
void SetSelection(Int_t Plane, Int_t Strip)
void SetPlane(Int_t v)
void SetStrip(Int_t v)
virtual ~BmnHistToF700()
void Register(THttpServer *serv)
void SetDir(TFile *outFile=NULL, TTree *recoTree=NULL)
BmnHistToF700(TString title="ToF700", TString path="")
TTree * frecoTree
Definition BmnHist.h:88
Int_t refID
Definition BmnHist.h:92
TString refPath
Definition BmnHist.h:90
TDirectory * fDir
Definition BmnHist.h:89
static void DrawRef(unique_ptr< TCanvas > &canGemStrip, vector< PadInfo * > *canGemStripPads)
Definition BmnHist.cxx:15
static void SetHistStyleTH1(TH1 *h)
Definition BmnHist.cxx:195
TString refRunName
Definition BmnHist.h:91
static BmnStatus LoadRefRun(Int_t refID, TString FullName, TString fTitle, vector< PadInfo * > canPads, vector< TString > Names)
Definition BmnHist.cxx:100
THttpServer * fServer
Definition BmnHist.h:87
Short_t GetPlane() const
Short_t GetStrip() const
Float_t GetTime() const
Float_t GetAmplitude() const
TClonesArray * tof700
Definition DigiArrays.h:130
TClonesArray * tof701
Definition DigiArrays.h:131
Storage for pad content and it's options.
Definition PadInfo.h:20
#define TOF2_MAX_STRIPS_IN_CHAMBER
#define TOF2_MAX_CHAMBERS
#define TOF_ROWS
#define TOF_COLS