Heavy ion Analysis Libriares
Loading...
Searching...
No Matches
QAPlot.cxx
1/*
2 * QAPlot.cxx
3 *
4 * Created on: 11 paź 2020
5 * Author: Daniel Wielanek
6 * E-mail: daniel.wielanek@gmail.com
7 * Warsaw University of Technology, Faculty of Physics
8 */
9
10#include "QAPlot.h"
11
12#include "ComplexEvent.h"
13#include "ComplexTrack.h"
14#include "Cout.h"
15#include "DataFormatManager.h"
16#include "Event.h"
17#include "HistogramManager.h"
18#include "QAPlotAxis.h"
19#include "StdHist.h"
20#include "Track.h"
21
22#include <TH1.h>
23#include <TH2.h>
24#include <TH3.h>
25#include <TVector3.h>
26
27namespace Hal {
28
29 QAPlot::QAPlot() : fUpdate(ECutUpdate::kNo), fReport(nullptr) {}
30
31 QAPlot::QAPlot(ECutUpdate upd) : QAPlot("", upd) {}
32
33 QAPlot::~QAPlot() {
34 if (fReport) delete fReport;
35 }
36
37 QAPlot::QAPlot(const QAPlot& other) : TNamed(other), fReport(nullptr) {
38 if (other.fReport) fReport = new QAPlotReport(*other.fReport);
39 fUpdate = other.fUpdate;
40 for (int i = 0; i < 3; i++) {
41 fSettings[i] = other.fSettings[i];
42 }
43 }
44
45 QAPlot::QAPlot(TString name, ECutUpdate upd) : fUpdate(upd), fReport(nullptr) { SetName(name); }
46
47 Bool_t QAPlot::Init(Int_t task_id) {
48 const Event* ev = DataFormatManager::Instance()->GetFormat(task_id, EFormatDepth::kNonBuffered);
49 if (ev == nullptr) return kFALSE;
50 switch (fUpdate) {
51 case ECutUpdate::kEvent: {
52 for (int j = 0; j < 3; j++) {
53 for (unsigned int i = 0; i < fSettings[j].size(); i++) {
54 if (fSettings[j][i].IsCustom()) continue;
55 Int_t fieldIDx = fSettings[j][i].GetFillFlagX();
56
57 TString axisNameX = ev->GetFieldName(fieldIDx);
58 fSettings[j][i].SetXaxisName(axisNameX);
59
60 if (j > 0) {
61 Int_t fieldIDy = fSettings[j][i].GetFillFlagY();
62 TString axisNameY = ev->GetFieldName(fieldIDy);
63 fSettings[j][i].SetYaxisName(axisNameY);
64 }
65 if (j > 1) {
66 Int_t fieldIDz = fSettings[j][i].GetFillFlagZ();
67 TString axisNameZ = ev->GetFieldName(fieldIDz);
68 fSettings[j][i].SetZaxisName(axisNameZ);
69 }
70 }
71 }
72 } break;
73 case ECutUpdate::kTrack: {
74 Track* track = ev->GetNewTrack();
75 ComplexEvent* zEvent = static_cast<ComplexEvent*>(ev->GetNewEvent());
76 track->SetEvent(zEvent);
77 Track* re_track = nullptr;
78 Track* im_track = nullptr;
79 const ComplexEvent* zev = dynamic_cast<const ComplexEvent*>(ev);
80 if (zev) {
81 ComplexTrack* z_track = (ComplexTrack*) track;
82 re_track = zev->GetNewRealTrack();
83 re_track->SetEvent(zEvent->GetRealEvent());
84 im_track = zev->GetNewImgTrack();
85 im_track->SetEvent(zEvent->GetImgEvent());
86 z_track->SetImgTrack(im_track);
87 z_track->SetRealTrack(re_track);
88 }
89
90 for (int j = 0; j < 3; j++) {
91 for (unsigned int i = 0; i < fSettings[j].size(); i++) {
92 if (fSettings[j][i].IsCustom()) continue;
93 Int_t fieldIDx = fSettings[j][i].GetFillFlagX();
94
95 TString axisNameX = track->GetFieldName(fieldIDx);
96 fSettings[j][i].SetXaxisName(axisNameX);
97
98 if (j > 0) {
99 Int_t fieldIDy = fSettings[j][i].GetFillFlagY();
100 TString axisNameY = track->GetFieldName(fieldIDy);
101 fSettings[j][i].SetYaxisName(axisNameY);
102 }
103 if (j > 1) {
104 Int_t fieldIDz = fSettings[j][i].GetFillFlagZ();
105 TString axisNameZ = track->GetFieldName(fieldIDz);
106 fSettings[j][i].SetZaxisName(axisNameZ);
107 }
108 }
109 }
110 delete track;
111 delete zEvent;
112 if (re_track) delete re_track;
113 if (im_track) delete im_track;
114 } break;
115 default: break;
116 }
117
118 fReport = new QAPlotReport(GetName(), GetSize1D(), GetSize2D(), GetSize3D());
119 for (int j = 0; j < 3; j++) {
120 for (unsigned int i = 0; i < fSettings[j].size(); i++) {
121 fReport->SetFlag(i, fSettings[j][i].GetFlag(), j + 1);
122 }
123 }
124 if (GetSize1D() > 0) {
125 fReport->f1dHistos->Init(GetSize1D(), 1, 0, 1);
126 for (int i = 0; i < GetSize1D(); i++) {
127 fReport->f1dHistos->OverwriteAt(((TH1D*) fSettings[0][i].MakeHisto()), i);
128 }
129 }
130 if (GetSize2D() > 0) {
131 fReport->f2dHistos->Init(GetSize2D(), 1, 0, 1, 1, 0, 1);
132 for (int i = 0; i < GetSize2D(); i++) {
133 fReport->f2dHistos->OverwriteAt(((TH2D*) fSettings[1][i].MakeHisto()), i);
134 }
135 }
136 if (GetSize3D() > 0) {
137 fReport->f3dHistos->Init(GetSize3D(), 1, 0, 1, 1, 0, 1, 1, 0, 1);
138 for (int i = 0; i < GetSize3D(); i++)
139 fReport->f3dHistos->OverwriteAt(((TH3D*) fSettings[2][i].MakeHisto()), i);
140 }
141 return kTRUE;
142 }
143
144 void QAPlot::Fill(TObject* obj) {
145 switch (fUpdate) {
146 case ECutUpdate::kEvent: {
147 FillEvent(static_cast<Event*>(obj));
148 } break;
149 case ECutUpdate::kTrack: {
150 FillTrack(static_cast<Track*>(obj));
151 } break;
152 default: break;
153 }
154 }
155
156 QAPlot& QAPlot::operator=(const QAPlot& other) {
157 if (&other == this) return *this;
158 if (fReport) delete fReport;
159 fReport = nullptr;
160 if (other.fReport) fReport = new QAPlotReport(*other.fReport);
161 TNamed::operator=(other);
162 fUpdate = other.fUpdate;
163 for (int i = 0; i < 3; i++) {
164 fSettings[i] = other.fSettings[i];
165 }
166 return *this;
167 }
168
169 QAPlotReport* QAPlot::GetReport() const {
170 QAPlotReport* rep = new QAPlotReport(*this->fReport);
171 rep->SetName(this->GetName());
172 rep->SetOriginClass(this->ClassName());
173 return rep;
174 }
175
176 Int_t QAPlot::AddTH1(TString name, const QAPlotAxis& x, TString flag) {
177 Int_t old_size = fSettings[0].size();
178 fSettings[0].push_back(QAHistoSettings(1));
179 fSettings[0][old_size].SetName(name);
180 fSettings[0][old_size].SetTitle(name);
181 fSettings[0][old_size].SetFlag(flag);
182 fSettings[0][old_size].SetFillFlagX(x.GetFlag());
183 fSettings[0][old_size].SetXaxis(x.GetNBins(), x.GetMin(), x.GetMax());
184 return old_size;
185 }
186
187 Int_t QAPlot::AddTH2(TString name, const QAPlotAxis& x, const QAPlotAxis& y, TString flag) {
188 Int_t old_size = fSettings[1].size();
189 fSettings[1].push_back(QAHistoSettings(2));
190 fSettings[1][old_size].SetName(name);
191 fSettings[1][old_size].SetTitle(name);
192 fSettings[1][old_size].SetFlag(flag);
193 fSettings[1][old_size].SetFillFlagX(x.GetFlag());
194 fSettings[1][old_size].SetFillFlagY(y.GetFlag());
195 fSettings[1][old_size].SetXaxis(x.GetNBins(), x.GetMin(), x.GetMax());
196 fSettings[1][old_size].SetYaxis(y.GetNBins(), y.GetMin(), y.GetMax());
197 return old_size;
198 }
199
200 Int_t QAPlot::AddTH3(TString name, const QAPlotAxis& x, const QAPlotAxis& y, const QAPlotAxis& z, TString flag) {
201 Int_t old_size = fSettings[2].size();
202 fSettings[2].push_back(QAHistoSettings(3));
203 fSettings[2][old_size].SetName(name);
204 fSettings[2][old_size].SetTitle(name);
205 fSettings[2][old_size].SetFlag(flag);
206 fSettings[2][old_size].SetFillFlagX(x.GetFlag());
207 fSettings[2][old_size].SetFillFlagY(y.GetFlag());
208 fSettings[2][old_size].SetFillFlagZ(z.GetFlag());
209 fSettings[2][old_size].SetXaxis(x.GetNBins(), x.GetMin(), x.GetMax());
210 fSettings[2][old_size].SetYaxis(y.GetNBins(), y.GetMin(), y.GetMax());
211 fSettings[2][old_size].SetZaxis(z.GetNBins(), z.GetMin(), z.GetMax());
212 return old_size;
213 }
214
215 void QAPlot::FillTrack(Track* track) {
216 for (unsigned int i = 0; i < fSettings[0].size(); i++) {
217 if (fSettings[0][i].IsCustom()) {
218 FillTrackCustom1D(track, (TH1D*) Get1D(i), i);
219 } else {
220 Get1D(i)->Fill(track->GetFieldVal(fSettings[0][i].GetFillFlagX()));
221 }
222 }
223 for (unsigned int i = 0; i < fSettings[1].size(); i++) {
224 if (fSettings[1][i].IsCustom()) {
225 FillTrackCustom2D(track, (TH2D*) Get2D(i), i);
226 } else {
227 Get2D(i)->Fill(track->GetFieldVal(fSettings[1][i].GetFillFlagX()), track->GetFieldVal(fSettings[1][i].GetFillFlagY()));
228 }
229 }
230 for (unsigned int i = 0; i < fSettings[2].size(); i++) {
231 if (fSettings[2][i].IsCustom()) {
232 FillTrackCustom3D(track, (TH3D*) Get3D(i), i);
233 } else {
234 Get3D(i)->Fill(track->GetFieldVal(fSettings[2][i].GetFillFlagX()),
235 track->GetFieldVal(fSettings[2][i].GetFillFlagY()),
236 track->GetFieldVal(fSettings[2][i].GetFillFlagZ()));
237 }
238 }
239 }
240
241 void QAPlot::FillEvent(Event* event) {
242 for (int i = 0; i < GetSize1D(); i++) {
243 if (fSettings[0][i].IsCustom()) {
244 FillEventCustom1D(event, (TH1D*) Get1D(i), i);
245 } else {
246 Get1D(i)->Fill(event->GetFieldVal(fSettings[0][i].GetFillFlagX()));
247 }
248 }
249 for (int i = 0; i < GetSize2D(); i++) {
250 if (fSettings[1][i].IsCustom()) {
251 FillEventCustom2D(event, (TH2D*) Get2D(i), i);
252 } else {
253 Get2D(i)->Fill(event->GetFieldVal(fSettings[1][i].GetFillFlagX()), event->GetFieldVal(fSettings[1][i].GetFillFlagY()));
254 }
255 }
256 for (int i = 0; i < GetSize3D(); i++) {
257 if (fSettings[2][i].IsCustom()) {
258 FillEventCustom3D(event, (TH3D*) Get3D(i), i);
259 } else {
260 Get3D(i)->Fill(event->GetFieldVal(fSettings[2][i].GetFillFlagX()),
261 event->GetFieldVal(fSettings[2][i].GetFillFlagY()),
262 event->GetFieldVal(fSettings[2][i].GetFillFlagZ()));
263 }
264 }
265 }
266
267 void QAPlot::Print(Option_t* /*opt*/) const {
268 for (int i = 0; i < 3; i++) {
269 for (unsigned int j = 0; j < fSettings[i].size(); j++) {
270 fSettings[i][j].Print();
271 }
272 }
273 }
274
275} // namespace Hal
const Event * GetFormat(Int_t task_id, EFormatDepth format_depth=EFormatDepth::kAll) const
static DataFormatManager * Instance()
Track * GetNewTrack() const
Definition Event.cxx:204
virtual Event * GetNewEvent() const
Definition Event.cxx:215
virtual TString GetFieldName(Int_t fieldID) const
Definition Event.cxx:272
void Init(Int_t size, HistogramAxisConf **axisconf, TString title, Bool_t Sumw=kFALSE)
void OverwriteAt(T *t, Int_t pos)
void OverwriteAt(T *t, Int_t pos)
void Init(Int_t size, HistogramAxisConf **axisconf, TString title, Bool_t Sumw=kFALSE)
void Init(Int_t size, HistogramAxisConf **axisconf, TString title, Bool_t Sumw=kFALSE)
void OverwriteAt(T *t, Int_t pos)
void SetFlag(Int_t no, TString flag, Int_t dim)
virtual void FillEventCustom3D(Event *, TH3D *, Int_t)
Definition QAPlot.h:105
virtual void Fill(TObject *obj)
Definition QAPlot.cxx:144
virtual Bool_t Init(Int_t id_task=-1)
Definition QAPlot.cxx:47
Int_t AddTH1(TString name, const QAPlotAxis &x, TString flag="")
Definition QAPlot.cxx:176
Int_t AddTH2(TString name, const QAPlotAxis &x, const QAPlotAxis &y, TString flag="")
Definition QAPlot.cxx:187
Int_t AddTH3(TString name, const QAPlotAxis &x, const QAPlotAxis &y, const QAPlotAxis &z, TString flag="")
Definition QAPlot.cxx:200
virtual void FillEventCustom2D(Event *, TH2D *, Int_t)
Definition QAPlot.h:93
TH1D * Get1D(Int_t no) const
Definition QAPlot.h:55
virtual void FillTrackCustom1D(Track *, TH1D *, Int_t)
Definition QAPlot.h:75
TH3D * Get3D(Int_t no) const
Definition QAPlot.h:67
TH2D * Get2D(Int_t no) const
Definition QAPlot.h:61
Int_t GetSize3D() const
Definition QAPlot.h:49
virtual void FillTrackCustom2D(Track *, TH2D *, Int_t)
Definition QAPlot.h:87
Int_t GetSize2D() const
Definition QAPlot.h:44
virtual void FillEventCustom1D(Event *, TH1D *, Int_t)
Definition QAPlot.h:81
virtual void FillTrackCustom3D(Track *, TH3D *, Int_t)
Definition QAPlot.h:99
void SetEvent(Event *event)
Definition Track.h:337
virtual Float_t GetFieldVal(Int_t fieldID) const
Definition Track.cxx:206
virtual TString GetFieldName(Int_t fieldID) const
Definition Track.cxx:244