R3BROOT
R3B analysis software
Loading...
Searching...
No Matches
R3BNeulandHitHist.cxx
Go to the documentation of this file.
1/******************************************************************************
2 * Copyright (C) 2019 GSI Helmholtzzentrum für Schwerionenforschung GmbH *
3 * Copyright (C) 2019-2025 Members of R3B Collaboration *
4 * *
5 * This software is distributed under the terms of the *
6 * GNU General Public Licence (GPL) version 3, *
7 * copied verbatim in the file "LICENSE". *
8 * *
9 * In applying this license GSI does not waive the privileges and immunities *
10 * granted to it by virtue of its status as an Intergovernmental Organization *
11 * or submit itself to any jurisdiction. *
12 ******************************************************************************/
13
14// -----------------------------------------------------------------------------
15// ----- R3BNeulandHitHist -----
16// ----- Created 22-04-2014 by D.Kresan -----
17// -----------------------------------------------------------------------------
18
19#include "R3BNeulandHitHist.h"
20#include "FairLogger.h"
21#include "FairRootManager.h"
22#include "R3BLosHitData.h"
23#include "R3BNeulandHit.h"
24#include "TClonesArray.h"
25#include "TH1F.h"
26#include "TH2F.h"
27#include "TMath.h"
28
29const Double_t C_LIGHT = 29.9792458; // Speed of light [cm/ns]
30
32 : fNofBars(100)
33 , fFirstPlaneHorisontal(kFALSE)
34 , fnEvents(0)
35 , fLandDigi(NULL)
36 , fLosHit(NULL)
37{
38}
39
40R3BNeulandHitHist::R3BNeulandHitHist(const char* name, Int_t iVerbose)
41 : FairTask(name, iVerbose)
42 , fNofBars(100)
43 , fFirstPlaneHorisontal(kFALSE)
44 , fnEvents(0)
45 , fLandDigi(NULL)
46 , fLosHit(NULL)
47{
48}
49
51
53{
54 FairRootManager* fMan = FairRootManager::Instance();
55 if (!fMan)
56 {
57 LOG(fatal) << "FairRootManager not found";
58 }
59 fLandDigi = dynamic_cast<TClonesArray*>(fMan->GetObject("NeulandHits"));
60 fLosHit = dynamic_cast<TClonesArray*>(fMan->GetObject("LosHit"));
61 CreateHistos();
62
63 return kSUCCESS;
64}
65
66void R3BNeulandHitHist::Exec(Option_t* option)
67{
68 if (!fLosHit || fLosHit->GetEntriesFast() == 0)
69 return;
70
71 R3BLosHitData* losHit = dynamic_cast<R3BLosHitData*>(fLosHit->At(0));
72 // Double_t startTime = losHit->GetTime(); TODO
73 Double_t startTime = 0;
74 Int_t count = 0;
75 Double_t tmean = 0.;
76
77 if (fLandDigi && fLandDigi->GetEntriesFast() > 0)
78 {
79 Int_t nLandDigi = fLandDigi->GetEntriesFast();
80 R3BNeulandHit* digi;
81 Int_t barId;
82 Double_t time, tof, norm_tof;
83 Double_t qdc, qdcL, qdcR;
84 Double_t x, y, z;
85 Double_t l;
86 Double_t beta;
87
88 for (Int_t i = 0; i < nLandDigi; i++)
89 {
90 digi = dynamic_cast<R3BNeulandHit*>(fLandDigi->At(i));
91
92 barId = digi->GetPaddle();
93 time = digi->GetT();
94 qdc = digi->GetE();
95 qdcL = digi->GetQdcL();
96 qdcR = digi->GetQdcR();
97 x = digi->GetPosition().X();
98 y = digi->GetPosition().Y();
99 z = digi->GetPosition().Z();
100 l = TMath::Sqrt(x * x + y * y + z * z);
101
102 tmean += time;
103 count += 1;
104
105 fh_land_qdcbarid->Fill(barId, qdc);
106 fh_land_yx->Fill(x, y);
107 if (barId < 51)
108 {
109 fh_land_yx1->Fill(x, y);
110 }
111 else if (barId < 101)
112 {
113 fh_land_yx2->Fill(x, y);
114 }
115
116 fh_los_corr->Fill(time, startTime);
117
118 fh_land_barid->Fill(barId);
119 fh_land_timebarid->Fill(barId, time);
120 tof = time - startTime - fTimeOffset;
121 beta = l / tof;
122 norm_tof = tof - (l - 1300.) / beta;
123 fh_land_qdctof->Fill(norm_tof, qdc);
124 fh_land_tofbarid->Fill(barId, tof);
125 fh_land_lbarid->Fill(barId, l);
126 fh_land_ltime->Fill(tof, l);
127 fh_land_betabarid->Fill(barId, beta);
128 fh_land_tof->Fill(tof);
129 if (beta > 29.7 && beta < 30.2)
130 fh_land_norm_tof->Fill(norm_tof);
131
132 fh_land_beta->Fill(beta);
133 fh_land_qdc->Fill(qdc);
134 if (beta > 16 && beta < 26)
135 fh_land_qdc_cut->Fill(qdc);
136 }
137
138 fnEvents += 1;
139 if (0 == (fnEvents % 100000))
140 {
141 LOG(info) << "R3BNeulandHitHist : " << fnEvents << " events collected, start time=" << startTime;
142 if (count)
143 {
144 LOG(info) << ", mean time=" << (tmean / (Double_t)count);
145 }
146 LOG(info);
147 }
148 }
149}
150
151void R3BNeulandHitHist::FinishTask() { WriteHistos(); }
152
153void R3BNeulandHitHist::CreateHistos()
154{
155 fh_land_barid = new TH1F("h_land_barid", "Bar ID", fNofBars, (Double_t)fNofBars + 0.5, 100.5);
156 fh_land_qdcbarid =
157 new TH2F("h_land_qdcbarid", "QDC vs Bar ID", fNofBars, 0.5, (Double_t)fNofBars + 0.5, 2000, 0., 40.);
158 fh_land_timebarid =
159 new TH2F("h_land_timebarid", "Time vs Bar ID", fNofBars, 0.5, (Double_t)fNofBars + 0.5, 4000, 0., 2000.);
160 fh_land_tofbarid =
161 new TH2F("h_land_tofbarid", "ToF vs Bar ID", fNofBars, 0.5, (Double_t)fNofBars + 0.5, 1000, 0., 100.);
162 fh_land_tof = new TH1F("h_land_tof", "ToF", 1000, 0., 100.);
163 fh_land_qdctof = new TH2F("h_land_qdctof", "QDC vs ToF", 2000, -50., 150., 2000, 0., 40.);
164 fh_land_betabarid =
165 new TH2F("h_land_betabarid", "Velocity vs Bar ID", fNofBars, 0.5, (Double_t)fNofBars + 0.5, 500, 0., 50.);
166 fh_land_lbarid =
167 new TH2F("h_land_lbarid", "Length vs Bar ID", fNofBars, 0.5, (Double_t)fNofBars + 0.5, 1000, 700., 900.);
168 fh_land_ltime = new TH2F("h_land_ltime", "Length vs Bar ID", 4000, 0., 2000., 1000, 700., 900.);
169 fh_land_yx = new TH2F("h_land_yx", "Y vs X", 340, -170., 170., 340, -170., 170.);
170 if (fFirstPlaneHorisontal)
171 {
172 fh_land_yx1 = new TH2F("h_land_yx1", "Y vs X", 340, -170., 170., 68, -170., 170.);
173 fh_land_yx2 = new TH2F("h_land_yx2", "Y vs X", 68, -170., 170., 340, -170., 170.);
174 }
175 else
176 {
177 fh_land_yx1 = new TH2F("h_land_yx1", "Y vs X", 68, -170., 170., 340, -170., 170.);
178 fh_land_yx2 = new TH2F("h_land_yx2", "Y vs X", 340, -170., 170., 68, -170., 170.);
179 }
180 fh_land_beta = new TH1F("h_land_beta", "Velocity", 4500, -10., 35.);
181 fh_land_qdc = new TH1F("h_land_qdc", "QDC", 2000, 0., 40.);
182 fh_land_qdc_cut = new TH1F("h_land_qdc_cut", "QDC cut", 2000, 0., 40.);
183 fh_los_time = new TH1F("h_los_time", "LOS Time", 2000, -1000., -800.);
184 fh_los_corr = new TH2F("h_los_corr", "LOS vs NeuLAND", 2000, 700., 900., 2000, -1000., -800.);
185
186 fh_land_norm_tof = new TH1F("h_land_norm_tof", "norm. ToF", 1000, 40., 50.);
187}
188
189void R3BNeulandHitHist::WriteHistos()
190{
191 fh_land_barid->Write();
192 fh_land_qdcbarid->Write();
193 fh_land_tof->Write();
194 fh_land_qdctof->Write();
195 fh_land_timebarid->Write();
196 fh_land_tofbarid->Write();
197 fh_land_yx->Write();
198 fh_land_yx1->Write();
199 fh_land_yx2->Write();
200 fh_land_lbarid->Write();
201 fh_land_ltime->Write();
202 fh_land_beta->Write();
203 fh_land_qdc->Write();
204 fh_land_qdc_cut->Write();
205 fh_land_betabarid->Write();
206 fh_los_time->Write();
207 fh_los_corr->Write();
208 fh_land_norm_tof->Write();
209}
210
ClassImp(R3B::Neuland::Cal2HitPar)
const Double_t C_LIGHT
virtual void Exec(Option_t *option)
virtual void FinishTask()
virtual InitStatus Init()
auto GetQdcR() const -> double
auto GetPaddle() const -> int
auto GetT() const -> double
auto GetPosition() const -> TVector3
auto GetQdcL() const -> double
auto GetE() const -> double