ATTPCROOT  0.3.0-alpha
A ROOT-based framework for analyzing data from active target detectors
AtTabPad.cxx
Go to the documentation of this file.
1 #include "AtTabPad.h"
2 
3 #include "AtAuxPad.h" // for AtAuxPad
4 #include "AtDataManip.h"
5 #include "AtEvent.h"
6 #include "AtHit.h" // for AtHit
7 #include "AtMap.h"
8 #include "AtPad.h"
9 #include "AtPadArray.h"
10 #include "AtPadBase.h" // for AtPadBase
11 #include "AtRawEvent.h"
12 #include "AtTabInfo.h"
13 #include "AtViewerManager.h"
14 
15 #include <FairLogger.h>
16 
17 #include <TCanvas.h>
18 #include <TF1.h>
19 #include <TH1.h> // for TH1D
20 
21 #include <iostream> // for operator<<, basic_ostream::operator<<
22 #include <memory> // for make_unique, allocator, unique_ptr
23 namespace DataHandling {
24 class AtSubject;
25 }
26 
27 constexpr auto cRED = "\033[1;31m";
28 constexpr auto cYELLOW = "\033[1;33m";
29 constexpr auto cNORMAL = "\033[0m";
30 constexpr auto cGREEN = "\033[1;32m";
31 constexpr auto cBLUE = "\033[1;34m";
32 
34 
35 AtTabPad::AtTabPad(int nRow, int nCol, TString name) : AtTabCanvas(name, nRow, nCol)
36 {
37  if (AtViewerManager::Instance() == nullptr)
38  throw "AtViewerManager must be initialized before creating tabs!";
39 
41  fPadNum->Attach(this);
42 }
44 {
45  fPadNum->Detach(this);
46 }
47 
49 {
50  std::cout << " ===== AtEventTabPad::Init =====" << std::endl;
51 
52  if (fTabName == "AtPad")
53  fTabName = TString::Format(fTabName + " %d", fTabId);
54 
55  auto man = AtViewerManager::Instance();
56  fTabInfo->AddAugment(std::make_unique<AtTabInfoFairRoot<AtEvent>>(man->GetEventBranch()));
57  fTabInfo->AddAugment(std::make_unique<AtTabInfoFairRoot<AtRawEvent>>(man->GetRawEventBranch()));
58 
59  std::cout << " AtEventTabPad::Init : Initialization complete! "
60  << "\n";
61 }
62 
63 void AtTabPad::MakeTab(TEveWindowSlot *slot)
64 {
66 }
67 
69 {
70  auto fRawEvent = GetFairRootInfo<AtRawEvent>();
71  if (fRawEvent == nullptr) {
72  LOG(debug) << "fRawEvent is nullptr for tab " << fTabId << "! Please set the raw event branch.";
73  return;
74  }
75 
76  // Redraw any Auxiliary channels
77  for (auto &[pos, toDraw] : fDrawMap) {
78  fCanvas->cd(pos + 1);
79  auto hist = toDraw.second;
80  if (toDraw.first == PadDrawType::kAuxPad) {
81  hist->Reset();
82 
83  auto auxPad = fRawEvent->GetAuxPad(fAugNames[pos]);
84  if (auxPad != nullptr)
85  DrawAdc(hist, *auxPad);
86  hist->Draw();
87  }
88  }
89  UpdateCvsPad();
90 }
91 
93 {
94  if (sub == fPadNum)
95  DrawPad();
96 }
97 
98 void AtTabPad::DrawPad()
99 {
100 
101  auto fRawEvent = GetFairRootInfo<AtRawEvent>();
102 
103  if (fRawEvent == nullptr) {
104  LOG(debug) << "fRawEvent is nullptr for tab " << fTabId << "! Please set the raw event branch.";
105  return;
106  }
107 
108  auto fPad = fRawEvent->GetPad(fPadNum->Get());
109 
110  for (auto &[pos, toDraw] : fDrawMap) {
111  if (toDraw.first == PadDrawType::kAuxPad)
112  continue;
113 
114  fCanvas->cd(pos + 1);
115  auto hist = toDraw.second;
116 
117  if (fPad == nullptr) {
118  hist->Reset();
119  hist->Draw();
120  continue;
121  }
122 
123  switch (toDraw.first) {
124  case PadDrawType::kADC: DrawAdc(hist, *fPad); break;
125  case PadDrawType::kRawADC: DrawRawAdc(hist, *fPad); break;
126  case PadDrawType::kArrAug: DrawArrayAug(hist, *fPad, fAugNames[pos]); break;
127  case PadDrawType::kFPN: DrawFPN(hist, *fPad); break;
129  auto auxPad = fRawEvent->GetAuxPad(fAugNames[pos]);
130  if (auxPad != nullptr)
131  DrawAdc(hist, *auxPad);
132  break;
133  }
134 
135  if (fDrawHits.find(pos) != fDrawHits.end())
136  DrawHit(*fPad, fDrawHits[pos]);
137  }
138 
139  UpdateCvsPad();
140 }
141 
142 void AtTabPad::DrawAdc(TH1D *hist, const AtPad &pad)
143 {
144  for (int i = 0; i < 512; i++) {
145  hist->SetBinContent(i + 1, pad.GetADC()[i]);
146  }
147  hist->Draw();
148 }
149 
150 void AtTabPad::DrawFPN(TH1D *hist, const AtPad &pad)
151 {
152 
153  auto fRawEvent = GetFairRootInfo<AtRawEvent>();
154  if (fRawEvent == nullptr) {
155  LOG(error) << "fRawEvent is nullptr for tab " << fTabId << "! Please set the raw event branch.";
156  return;
157  }
158  auto map = AtViewerManager::Instance()->GetMap();
159  auto fpnPad = fRawEvent->GetFpn(map->GetNearestFPN(pad.GetPadNum()));
160 
161  for (int i = 0; i < 512; i++) {
162  hist->SetBinContent(i + 1, fpnPad->GetRawADC()[i]);
163  }
164  hist->Draw();
165 }
166 
167 void AtTabPad::DrawRawAdc(TH1D *hist, const AtPad &pad)
168 {
169  for (int i = 0; i < 512; i++) {
170  hist->SetBinContent(i + 1, pad.GetRawADC()[i]);
171  }
172  hist->Draw();
173 }
174 
175 void AtTabPad::DrawArrayAug(TH1D *hist, const AtPad &pad, TString augName)
176 {
177  auto aug = dynamic_cast<const AtPadArray *>(pad.GetAugment(augName.Data()));
178  if (aug == nullptr)
179  return;
180 
181  for (int i = 0; i < 512; i++) {
182  hist->SetBinContent(i + 1, aug->GetArray()[i]);
183  }
184  hist->Draw();
185 }
186 
187 void AtTabPad::DrawHit(const AtPad &pad, TF1Vec &vec)
188 {
189  vec.clear();
190 
191  // Loop through all hits and
192  auto event = GetFairRootInfo<AtEvent>();
193  for (auto &hit : event->GetHits()) {
194  if (hit->GetPadNum() != pad.GetPadNum())
195  continue;
196 
197  LOG(debug) << "Drawing hit with charge " << hit->GetCharge();
198  LOG(debug) << hit->GetPosition().Z() << " " << hit->GetPositionSigma().Z();
199  auto func = AtTools::GetHitFunctionTB(*hit);
200  if (func != nullptr)
201  vec.push_back(std::move(func));
202  }
203 
204  for (auto &func : vec) {
205  LOG(debug) << "Drawing hit function";
206  func->Draw("SAME");
207  }
208 }
209 std::string AtTabPad::GetName(int pos, PadDrawType type)
210 {
211  switch (type) {
212  case PadDrawType::kADC: return "ADC";
213  case PadDrawType::kRawADC: return "Raw ADC";
214  case PadDrawType::kArrAug: return fAugNames[pos];
215  case PadDrawType::kFPN: return "Closest FPN";
216  case PadDrawType::kAuxPad: return fAugNames[pos];
217  }
218  return "";
219 }
220 void AtTabPad::SetDraw(Int_t pos, PadDrawType type)
221 {
222  auto name = TString::Format("padHist_Tab%i_%i", fTabId, pos);
223  TH1D *padHist = new TH1D(name, GetName(pos, type).c_str(), 512, 0, 512);
224  fDrawMap.emplace(pos, std::make_pair(type, padHist));
225 }
226 
227 void AtTabPad::DrawHits(int row, int col)
228 {
229  fDrawHits[row * fCols + col] = TF1Vec();
230 }
231 void AtTabPad::DrawADC(int row, int col)
232 {
233  SetDraw(row * fCols + col, PadDrawType::kADC);
234 }
235 
236 void AtTabPad::DrawFPN(int row, int col)
237 {
238  SetDraw(row * fCols + col, PadDrawType::kFPN);
239 }
240 
241 void AtTabPad::DrawRawADC(int row, int col)
242 {
243  SetDraw(row * fCols + col, PadDrawType::kRawADC);
244 }
245 
246 void AtTabPad::DrawArrayAug(TString augName, int row, int col)
247 {
248  fAugNames.emplace(row * fCols + col, augName);
249  SetDraw(row * fCols + col, PadDrawType::kArrAug);
250 }
251 
252 void AtTabPad::DrawAuxADC(TString auxName, int row, int col)
253 {
254  fAugNames.emplace(row * fCols + col, auxName);
255  SetDraw(row * fCols + col, PadDrawType::kAuxPad);
256 }
257 
258 void AtTabPad::UpdateCvsPad()
259 {
260  fCanvas->Modified();
261  fCanvas->Update();
262 }
AtTabPad::~AtTabPad
~AtTabPad()
Definition: AtTabPad.cxx:43
AtPad.h
AtTabCanvas::fCanvas
TCanvas * fCanvas
Definition: AtTabCanvas.h:22
cNORMAL
constexpr auto cNORMAL
Definition: AtTabPad.cxx:29
AtRawEvent.h
AtTabPad::PadDrawType::kArrAug
@ kArrAug
AtEvent.h
AtAuxPad.h
AtTabPad::PadDrawType::kADC
@ kADC
AtTabPad::PadDrawType::kFPN
@ kFPN
cRED
constexpr auto cRED
Definition: AtTabPad.cxx:27
AtTabPad::InitTab
void InitTab() override
Definition: AtTabPad.cxx:48
AtTabPad::fAugNames
std::unordered_map< Int_t, std::string > fAugNames
Let root handle hist memory.
Definition: AtTabPad.h:41
AtTabPad::fDrawHits
std::unordered_map< Int_t, TF1Vec > fDrawHits
Definition: AtTabPad.h:44
DataHandling::AtSubject::Detach
void Detach(AtObserver *observer)
Detach an observer to stop getting notified when this subject changes.
Definition: AtDataSubject.h:37
AtViewerManager::Instance
static AtViewerManager * Instance()
Definition: AtViewerManager.cxx:43
AtTabPad::DrawRawADC
void DrawRawADC(int row=0, int col=0)
Definition: AtTabPad.cxx:241
AtTabPad::AtTabPad
AtTabPad(int nRow=1, int nCol=1, TString name="AtPad")
Definition: AtTabPad.cxx:35
AtTabPad
Definition: AtTabPad.h:33
AtTabCanvas::MakeTab
void MakeTab(TEveWindowSlot *slot) override
Create the gui components of the tab in the passed window slot.
Definition: AtTabCanvas.cxx:8
AtTabBase::fTabId
Int_t fTabId
Definition: AtTabBase.h:30
AtDataManip.h
AtViewerManager::GetPadNum
DataHandling::AtPadNum & GetPadNum()
Definition: AtViewerManager.h:75
AtTabPad::Exec
void Exec() override
Called after the run's Exec() to update tab.
Definition: AtTabPad.cxx:68
AtTabPad::fPadNum
DataHandling::AtPadNum * fPadNum
Definition: AtTabPad.h:42
AtPadArray.h
cGREEN
constexpr auto cGREEN
Definition: AtTabPad.cxx:30
AtHit.h
AtViewerManager::GetMap
AtMap * GetMap()
Definition: AtViewerManager.h:66
AtTools::GetHitFunctionTB
std::unique_ptr< TF1 > GetHitFunctionTB(const AtHit &hit, const AtDigiPar *par=nullptr)
Get charge as a function of TB.
Definition: AtDataManip.cxx:74
AtTabPad::DrawFPN
void DrawFPN(int row=0, int col=0)
Definition: AtTabPad.cxx:236
AtTabPad::TF1Vec
std::vector< std::unique_ptr< TF1 > > TF1Vec
Definition: AtTabPad.h:36
AtTabPad::Update
void Update(DataHandling::AtSubject *sub) override
Definition: AtTabPad.cxx:92
DataHandling
Definition: AtDataObserver.h:4
AtTabPad::PadDrawType::kAuxPad
@ kAuxPad
AtPad::GetAugment
AtPadBase * GetAugment(std::string name)
Definition: AtPad.cxx:82
AtTabPad::DrawHits
void DrawHits(int row=0, int col=0)
Definition: AtTabPad.cxx:227
AtTabPad::DrawADC
void DrawADC(int row=0, int col=0)
Definition: AtTabPad.cxx:231
DataHandling::AtSimpleType::Get
T Get() const
Definition: AtDataSubject.h:57
AtPad::GetPadNum
Int_t GetPadNum() const
Definition: AtPad.h:96
AtPad::GetADC
const trace & GetADC() const
Definition: AtPad.cxx:97
AtTabCanvas
Abstract class for a tab composed of a single TCanvas.
Definition: AtTabCanvas.h:20
AtTabInfoFairRoot< AtEvent >
AtTabPad::MakeTab
void MakeTab(TEveWindowSlot *) override
Create the gui components of the tab in the passed window slot.
Definition: AtTabPad.cxx:63
AtTabPad::DrawArrayAug
void DrawArrayAug(TString augName, int row=0, int col=0)
Definition: AtTabPad.cxx:246
AtTabPad::PadDrawType::kRawADC
@ kRawADC
DataHandling::AtSubject::Attach
void Attach(AtObserver *observer)
Attach an observer to get notified when this subject changes.
Definition: AtDataSubject.cxx:12
AtTabPad.h
AtMap.h
DataHandling::AtSubject
Definition: AtDataSubject.h:24
AtPadArray
Holds an addition array of doubles for an AtPad.
Definition: AtPadArray.h:24
AtPad
Container class for AtPadBase objects.
Definition: AtPad.h:38
cYELLOW
constexpr auto cYELLOW
Definition: AtTabPad.cxx:28
AtViewerManager.h
AtTabPad::DrawAuxADC
void DrawAuxADC(TString auxName, int row=0, int col=0)
Definition: AtTabPad.cxx:252
ClassImp
ClassImp(AtTabPad)
AtPad::GetRawADC
const rawTrace & GetRawADC() const
Definition: AtPad.h:104
AtPadBase.h
AtTabBase::fTabName
TString fTabName
Definition: AtTabBase.h:31
cBLUE
constexpr auto cBLUE
Definition: AtTabPad.cxx:31
AtTabCanvas::fCols
Int_t fCols
Definition: AtTabCanvas.h:23
AtTabInfo.h
AtTabBase::fTabInfo
std::unique_ptr< AtTabInfo > fTabInfo
Definition: AtTabBase.h:33
AtTabPad::fDrawMap
std::unordered_map< Int_t, std::pair< PadDrawType, TH1D * > > fDrawMap
Definition: AtTabPad.h:40