ATTPCROOT  0.3.0-alpha
A ROOT-based framework for analyzing data from active target detectors
AtTabMain.cxx
Go to the documentation of this file.
1 #include "AtTabMain.h"
2 
3 #include "AtContainerManip.h" // for GetPointerVector
4 #include "AtEvent.h" // for AtEvent, AtEvent::HitVector
5 #include "AtHit.h" // for AtHit, AtHit::XYZPoint
6 #include "AtMap.h" // for AtMap
7 #include "AtPad.h" // for AtPad
8 #include "AtPadReference.h" // for operator<<
9 #include "AtPattern.h" // for AtPattern
10 #include "AtPatternEvent.h" // for AtPatternEvent
11 #include "AtRawEvent.h" // for AtRawEvent
12 #include "AtTabInfo.h" // for AtTabInfoFairRoot, AtTabInfoBase
13 #include "AtTrack.h" // for AtTrack
14 #include "AtViewerManager.h" // for AtViewerManager
15 
16 #include <FairLogger.h> // for LOG, Logger
17 
18 #include <Math/Point3D.h> // for PositionVector3D
19 #include <TAttMarker.h> // for TAttMarker
20 #include <TCanvas.h> // for TCanvas
21 #include <TEveBrowser.h> // for TEveBrowser
22 #include <TEveEventManager.h> // for TEveEventManager
23 #include <TEveGeoNode.h> // for TEveGeoTopNode
24 #include <TEveManager.h> // for TEveManager, gEve
25 #include <TEvePointSet.h> // for TEvePointSet
26 #include <TEveViewer.h> // for TEveViewer
27 #include <TEveWindow.h> // for TEveWindowPack, TEveWindowSlot, TEv...
28 #include <TGLCamera.h> // for TGLCamera
29 #include <TGLViewer.h> // for TGLViewer
30 #include <TGTab.h> // for TGTab
31 #include <TGeoManager.h> // for gGeoManager, TGeoManager
32 #include <TGeoVolume.h> // for TGeoVolume
33 #include <TH1.h> // for TH1I
34 #include <TH2Poly.h> // for TH2Poly
35 #include <TNamed.h> // for TNamed
36 #include <TObject.h> // for TObject
37 #include <TRootEmbeddedCanvas.h> // for TRootEmbeddedCanvas
38 #include <TString.h> // for Form, TString
39 #include <TStyle.h> // for TStyle, gStyle
40 #include <TVirtualPad.h> // for TVirtualPad, gPad
41 #include <TVirtualX.h> // for TVirtualX, gVirtualX
42 
43 #include <algorithm> // for max
44 #include <array> // for array
45 #include <cstdio> // for sprintf
46 #include <iostream> // for operator<<, endl, basic_ostream
47 #include <utility> // for move
48 // IWYU pragma: no_include <ext/alloc_traits.h>
49 namespace DataHandling {
50 class AtSubject;
51 }
52 class TGeoNode; // lines 45-45
53 
54 constexpr auto cRED = "\033[1;31m";
55 constexpr auto cYELLOW = "\033[1;33m";
56 constexpr auto cNORMAL = "\033[0m";
57 constexpr auto cGREEN = "\033[1;32m";
58 constexpr auto cBLUE = "\033[1;34m";
59 
61 
63 {
64  if (AtViewerManager::Instance() == nullptr)
65  throw "AtViewerManager must be initialized before creating tabs!";
66 
68  fPadNum->Attach(this);
69 
71  fEventBranch->Attach(this);
72 
74  fRawEventBranch->Attach(this);
75 
78 
80  fEntry->Attach(this);
81 };
82 
84 {
85  fPadNum->Detach(this);
86  fEventBranch->Detach(this);
87  fRawEventBranch->Detach(this);
89  fEntry->Detach(this);
90 }
91 
93 {
94  std::cout << " ===== AtTabMain::Init =====" << std::endl;
95 
96  gEve->AddEvent(fEveEvent.get());
97  fEveEvent->AddElement(fHitSet.get());
98 
99  gEve->AddEvent(fEvePatternEvent.get());
100 
101  fTabInfo->AddAugment(std::make_unique<AtTabInfoFairRoot<AtEvent>>(*fEventBranch));
102  fTabInfo->AddAugment(std::make_unique<AtTabInfoFairRoot<AtRawEvent>>(*fRawEventBranch));
103  fTabInfo->AddAugment(std::make_unique<AtTabInfoFairRoot<AtPatternEvent>>(*fPatternEventBranch));
104 
105  gStyle->SetPalette(55);
106 
107  std::cout << " AtTabMain::Init : Initialization complete! "
108  << "\n";
109 }
110 
111 void AtTabMain::ExpandNumPatterns(int num)
112 {
113 
114  // Expand vector so it's large enough for all of the patterns in the event
115  if (fPatternHitSets.size() < num)
116  LOG(info) << "Expanding number of patterns to " << num << " from " << fPatternHitSets.size() << std::endl;
117  while (fPatternHitSets.size() < num) {
118  int trackID = fPatternHitSets.size();
119 
120  auto trackSet = std::make_unique<TEvePointSet>(TString::Format("Track_%d", trackID));
121  trackSet->SetDestroyOnZeroRefCnt(false);
122  fHitAttr.Copy(*trackSet);
123  trackSet->SetMarkerColor(GetTrackColor(fPatternHitSets.size()));
124 
125  fPatternHitSets.push_back(std::move(trackSet));
126  }
127 }
128 
130 {
131  std::vector<Color_t> colors = {kBlue - 7, kGreen - 8, kOrange, kViolet, kYellow, kTeal - 6,
132  kMagenta + 1, kBlue, kViolet, kYellow, kCyan};
133  if (i < colors.size()) {
134  return colors.at(i);
135  } else
136  return kAzure;
137 }
138 
140 {
141  // If we should update the stuff that depends on the AtEvent
142  if (sub == fEventBranch || sub == fEntry) {
143  UpdatePadPlane();
144  UpdateEventElements();
145  }
146  if (sub == fPatternEventBranch || sub == fEntry) {
147  UpdatePatternEventElements();
148  }
149  if (sub == fRawEventBranch || sub == fEntry || sub == fPadNum) {
150  DrawWave(fPadNum->Get());
151  }
152 
153  // If we should update the 3D display
154  if (sub == fEventBranch || sub == fPatternEventBranch || sub == fEntry) {
155  gEve->Redraw3D(false); // false -> don't reset camera
156  }
157 }
158 
159 void AtTabMain::DrawPadPlane()
160 {
161  if (fPadPlane) {
162  fPadPlane->Reset(nullptr);
163  return;
164  }
165 
168  fPadPlane->SetBit(TH1::kNoTitle);
169  fCvsPadPlane->cd();
170  fPadPlane->Draw("COL L0");
171  fPadPlane->SetMinimum(1.0);
172  gStyle->SetOptStat(0);
173  gStyle->SetPalette(103);
174  gPad->Update();
175 }
176 
177 void AtTabMain::DrawPadWave()
178 {
179  char name[20];
180  sprintf(name, "fPadWave_DT%i", fTabId);
181  fPadWave = new TH1I(name, name, 512, 0, 511);
182  fPadWave->SetBit(TH1::kNoTitle);
183  fCvsPadWave->cd();
184  fPadWave->Draw();
185 }
186 
187 void AtTabMain::DumpEvent(std::string fileName)
188 {
189  auto fEvent = GetFairRootInfo<AtEvent>();
190  if (fEvent == nullptr) {
191  std::cout << "CRITICAL ERROR: Event missing for TabMain. Aborting draw." << std::endl;
192  return;
193  }
194 
195  std::ofstream dumpEvent;
196  dumpEvent.open(fileName);
197  dumpEvent << " Event ID : " << fEvent->GetEventID() << std::endl;
198 
199  for (auto &hit : fEvent->GetHits())
200  dumpEvent << hit->GetPosition().X() << "," << hit->GetPosition().Y() << "," << hit->GetPosition().Z() << ","
201  << hit->GetCharge() << std::endl;
202 
203  dumpEvent.close();
204 }
205 
207 {
208  fEveEvent->SetRnrState(true);
209  fEvePatternEvent->SetRnrState(false);
210 }
211 
212 void AtTabMain::UpdatePatternEventElements()
213 {
214  if (fEvePatternEvent == nullptr)
215  return;
216 
217  auto fPatternEvent = GetFairRootInfo<AtPatternEvent>();
218  if (fPatternEvent == nullptr) {
219  LOG(debug) << "Cannot update AtPatternEvent elements: no event availible";
220  return;
221  }
222 
223  // Make sure we have enough TEve elements to draw all the tracks
224  auto &tracks = fPatternEvent->GetTrackCand();
225  ExpandNumPatterns(tracks.size());
226 
227  // Remove all the elements, and re-add them
228  fEvePatternEvent->RemoveElements();
229  for (int i = 0; i < tracks.size(); ++i) {
230  if (tracks[i].GetPattern() == nullptr)
231  continue;
232 
233  // Update the hit points and re-add them to the event
234  auto hitSet = fPatternHitSets.at(i).get();
235  fHitAttr.Copy(*hitSet);
236  hitSet->SetMarkerColor(GetTrackColor(i));
237  SetPointsFromTrack(*hitSet, tracks[i]);
238  fEvePatternEvent->AddElement(hitSet);
239 
240  // Get the pattern and add it to the event
241  auto pattern = tracks[i].GetPattern()->GetEveElement();
242  pattern->SetDestroyOnZeroRefCnt(false);
243  pattern->SetMainColor(GetTrackColor(i));
244  fEvePatternEvent->AddElement(pattern);
245  }
246 }
247 
248 void AtTabMain::UpdateEventElements()
249 {
250  auto fEvent = GetFairRootInfo<AtEvent>();
251  if (fEvent == nullptr) {
252  LOG(debug) << "Cannot update AtEvent elements: no event availible";
253  return;
254  }
255 
256  auto &hits = fEvent->GetHits();
257  LOG(info) << cBLUE << " Number of hits : " << hits.size() << " in " << fEvent->GetEventID() << cNORMAL;
258 
259  SetPointsFromHits(*fHitSet, hits);
260 }
261 
262 void AtTabMain::UpdatePadPlane()
263 {
264 
265  if (fPadPlane)
266  fPadPlane->Reset(nullptr);
267  else
268  return;
269 
270  LOG(debug) << "Updating pad plane ";
271 
272  auto fEvent = GetFairRootInfo<AtEvent>();
273  if (fEvent == nullptr) {
274  LOG(debug) << "Cannot fill pad plane histogram: no event availible";
275  return;
276  }
277  auto &hits = fEvent->GetHits();
278 
279  for (auto &hit : hits) {
280  int padMultiHit = GetFairRootInfo<AtEvent>()->GetHitPadMult(hit->GetPadNum());
281  if (hit->GetCharge() < fThreshold || padMultiHit > fMaxHitMulti)
282  continue;
283  auto position = hit->GetPosition();
284  fPadPlane->Fill(position.X(), position.Y(), hit->GetCharge());
285  }
286 
287  fCvsPadPlane->Modified();
288  fCvsPadPlane->Update();
289 }
290 void AtTabMain::SetPointsFromHits(TEvePointSet &hitSet, const std::vector<std::unique_ptr<AtHit>> &hits)
291 {
293 }
294 
295 void AtTabMain::SetPointsFromHits(TEvePointSet &hitSet, const std::vector<AtHit *> &hits)
296 {
297  Int_t nHits = hits.size();
298 
299  hitSet.Reset(nHits);
300  hitSet.SetOwnIds(true);
301  fHitAttr.Copy(hitSet); // Copy attributes from fHitAttr into hitSet.
302 
303  for (Int_t iHit = 0; iHit < nHits; iHit++) {
304 
305  auto &hit = *hits.at(iHit);
306  Int_t PadMultHit = 0;
307  if (GetFairRootInfo<AtEvent>())
308  PadMultHit = GetFairRootInfo<AtEvent>()->GetHitPadMult(hit.GetPadNum());
309 
310  if (hit.GetCharge() < fThreshold || PadMultHit > fMaxHitMulti)
311  continue;
312 
313  auto position = hit.GetPosition();
314 
315  hitSet.SetNextPoint(position.X() / 10., position.Y() / 10., position.Z() / 10.); // Convert into cm
316  hitSet.SetPointId(new TNamed(Form("Hit %d", iHit), ""));
317  }
318 
319  gEve->ElementChanged(&hitSet);
320 }
321 
322 void AtTabMain::SetPointsFromTrack(TEvePointSet &hitSet, const AtTrack &track)
323 {
324  Int_t nHits = track.GetHitArray().size();
325 
326  hitSet.Reset(nHits);
327  hitSet.SetOwnIds(true);
328 
329  for (Int_t i = 0; i < nHits; i++) {
330 
331  auto &hit = *track.GetHitArray()[i];
332 
333  if (hit.GetCharge() < fThreshold)
334  continue;
335 
336  auto position = hit.GetPosition();
337  hitSet.SetNextPoint(position.X() / 10., position.Y() / 10., position.Z() / 10.); // Convert into cm
338  hitSet.SetPointId(new TNamed(Form("Hit %d", i), ""));
339  }
340 
341  gEve->ElementChanged(&hitSet);
342 }
343 
344 bool AtTabMain::DrawWave(Int_t PadNum)
345 {
346  fPadWave->Reset();
347  auto fRawEvent = GetFairRootInfo<AtRawEvent>();
348 
349  // std::cout << "checking fRawEvent" << std::endl;
350  if (fRawEvent == nullptr) {
351  return false;
352  }
353  // std::cout << "fRawEvent is not nullptr" << std::endl;
354  AtPad *fPad = fRawEvent->GetPad(PadNum);
355  if (fPad == nullptr) {
356  LOG(error) << "Pad in event is null!";
357  return false;
358  } else
359  LOG(debug) << "Drawing pad " << fPad->GetPadNum();
360 
361  auto adc = fPad->GetADC();
362 
363  for (Int_t i = 0; i < 512; i++) {
364  fPadWave->SetBinContent(i, adc[i]);
365  }
366 
367  fCvsPadWave->cd();
368  fPadWave->Draw();
369  fCvsPadWave->Modified();
370  fCvsPadWave->Update();
371 
372  return true;
373 }
374 
375 void AtTabMain::MakeTab(TEveWindowSlot *slot)
376 {
377  TEveWindowPack *pack = nullptr;
378 
379  // 3D
380  pack = slot->MakePack();
381  pack->SetElementName("Main");
382  pack->SetHorizontal();
383  // pack->SetVertical();
384  pack->SetShowTitleBar(kFALSE);
385 
386  pack->NewSlot()->MakeCurrent();
387  TEveViewer *view3D = gEve->SpawnNewViewer("3D View", "");
388  view3D->AddScene(gEve->GetGlobalScene());
389  view3D->AddScene(gEve->GetEventScene());
390  // }
391 
392  slot = pack->NewSlot();
393  TEveWindowPack *pack2 = slot->MakePack();
394  pack2->SetShowTitleBar(kFALSE);
395  pack2->SetVertical();
396  slot = pack2->NewSlot();
397  slot->StartEmbedding();
398  fCvsPadWave = new TCanvas("AtPad Canvas");
399  fCvsPadWave->ToggleEditor();
400  slot->StopEmbedding();
401 
402  // Pad Plane
403  slot = pack2->NewSlotWithWeight(1.5);
404  auto *ecvs = new TRootEmbeddedCanvas();
405  TEveWindowFrame *frame = slot->MakeFrame(ecvs);
406  frame->SetElementName("AtTPC Pad Plane");
407  pack->GetEveFrame()->SetShowTitleBar(kFALSE);
408  fCvsPadPlane = ecvs->GetCanvas();
409  fCvsPadPlane->AddExec("ex", "AtTabMain::SelectPad()");
410 
411  fCvsPadWave->SetName(TString::Format("fCvsPadWave_DT%i", fTabId));
412  DrawPadWave();
413 
414  fCvsPadPlane->ToggleEventStatus();
415  DrawPadPlane();
416 
417  if (gGeoManager) {
418  TGeoNode *geoNode = gGeoManager->GetTopNode();
419  Int_t option = 1;
420  Int_t level = 3;
421  Int_t nNodes = 10000;
422  auto *topNode = new TEveGeoTopNode(gGeoManager, geoNode, option, level, nNodes);
423  gEve->AddGlobalElement(topNode);
424 
425  Int_t transparency = 80;
426  gGeoManager->GetVolume("drift_volume")->SetTransparency(transparency);
427  gEve->FullRedraw3D(kTRUE);
428  }
429 
430  gEve->GetBrowser()->GetTabRight()->SetTab(1);
431 
432  gEve->Redraw3D(true, true);
433 
434  TGLViewer *dfViewer = gEve->GetDefaultGLViewer(); // Is this doing anything?
435  dfViewer->CurrentCamera().RotateRad(-.7, 0.5);
436  dfViewer->DoDraw();
438 }
439 
441 {
442 
443  // Return earlyt if this was not triggered by a left mouse click.
444  if (gPad->GetEvent() != 11)
445  return;
446  auto *h = dynamic_cast<TH2Poly *>(gPad->GetSelected());
447  if (!h)
448  return;
449 
450  gPad->GetCanvas()->FeedbackMode(true);
451 
452  int pyold = gPad->GetUniqueID();
453  int px = gPad->GetEventX();
454  int py = gPad->GetEventY();
455  float uxmin = gPad->GetUxmin();
456  float uxmax = gPad->GetUxmax();
457  int pxmin = gPad->XtoAbsPixel(uxmin);
458  int pxmax = gPad->XtoAbsPixel(uxmax);
459  if (pyold)
460  gVirtualX->DrawLine(pxmin, pyold, pxmax, pyold);
461  gVirtualX->DrawLine(pxmin, py, pxmax, py);
462  gPad->SetUniqueID(py);
463  Float_t upx = gPad->AbsPixeltoX(px);
464  Float_t upy = gPad->AbsPixeltoY(py);
465  Double_t x = gPad->PadtoX(upx);
466  Double_t y = gPad->PadtoY(upy);
467  Int_t bin = h->FindBin(x, y);
468  const char *bin_name = h->GetBinName(bin);
469  LOG(debug) << " ==========================";
470  LOG(debug) << " Bin number selected : " << bin << " Bin name :" << bin_name;
471 
473  if (tmap == nullptr) {
474  LOG(fatal) << "AtMap not set! Pass a valid map to the constructor of AtViewerManager!";
475  } else {
476  Int_t tPadNum = tmap->BinToPad(bin);
477  LOG(debug) << " Bin : " << bin << " to Pad : " << tPadNum;
478  LOG(debug) << " Electronic mapping: " << tmap->GetPadRef(tPadNum);
479  LOG(debug) << " Raw Event Pad Num " << tPadNum;
481  }
482 }
AtMap
Definition: AtMap.h:33
AtPad.h
AtTabMain::fRawEventBranch
DataHandling::AtBranch * fRawEventBranch
Definition: AtTabMain.h:59
AtRawEvent.h
cNORMAL
constexpr auto cNORMAL
Definition: AtTabMain.cxx:56
AtTabMain::UpdateRenderState
virtual void UpdateRenderState()
Definition: AtTabMain.cxx:206
AtTabMain::DumpEvent
void DumpEvent(std::string file)
Definition: AtTabMain.cxx:187
AtMap::GetPadPlane
TH2Poly * GetPadPlane()
Definition: AtMap.cxx:70
AtEvent.h
AtTabMain
Main tab in viewer for 3D and pad selection.
Definition: AtTabMain.h:34
AtTabMain::fEntry
DataHandling::AtTreeEntry * fEntry
Definition: AtTabMain.h:61
AtTabMain::AtTabMain
AtTabMain()
Definition: AtTabMain.cxx:62
AtPadReference.h
AtTabMain::fEventBranch
DataHandling::AtBranch * fEventBranch
Definition: AtTabMain.h:58
AtTabMain::fPatternEventBranch
DataHandling::AtBranch * fPatternEventBranch
Definition: AtTabMain.h:60
AtPattern.h
AtTabMain::fMaxHitMulti
Int_t fMaxHitMulti
Definition: AtTabMain.h:48
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
cRED
constexpr auto cRED
Definition: AtTabMain.cxx:54
AtViewerManager::GetRawEventBranch
DataHandling::AtBranch & GetRawEventBranch()
Definition: AtViewerManager.h:71
ClassImp
ClassImp(AtTabMain)
AtTabMain::fPatternHitSets
std::vector< TEvePointSetPtr > fPatternHitSets
Definition: AtTabMain.h:44
AtTabMain::fCvsPadPlane
TCanvas * fCvsPadPlane
Definition: AtTabMain.h:52
AtViewerManager::GetPatternEventBranch
DataHandling::AtBranch & GetPatternEventBranch()
Definition: AtViewerManager.h:73
AtTabMain::SetPointsFromTrack
void SetPointsFromTrack(TEvePointSet &hitSet, const AtTrack &track)
Definition: AtTabMain.cxx:322
AtViewerManager::GetCurrentEntry
DataHandling::AtTreeEntry & GetCurrentEntry()
Definition: AtViewerManager.h:74
AtTabMain::fPadPlane
TH2Poly * fPadPlane
Definition: AtTabMain.h:53
AtTabMain.h
AtTabMain::MakeTab
void MakeTab(TEveWindowSlot *slot) override
Create the gui components of the tab in the passed window slot.
Definition: AtTabMain.cxx:375
AtTabBase::fTabId
Int_t fTabId
Definition: AtTabBase.h:30
AtTabMain::SetPointsFromHits
void SetPointsFromHits(TEvePointSet &hitSet, const std::vector< std::unique_ptr< AtHit >> &hits)
Definition: AtTabMain.cxx:290
AtViewerManager::GetPadNum
DataHandling::AtPadNum & GetPadNum()
Definition: AtViewerManager.h:75
AtTabMain::fCvsPadWave
TCanvas * fCvsPadWave
Definition: AtTabMain.h:55
AtHit.h
AtViewerManager::GetMap
AtMap * GetMap()
Definition: AtViewerManager.h:66
AtTabMain::GetTrackColor
Color_t GetTrackColor(int i)
Definition: AtTabMain.cxx:129
AtTabMain::~AtTabMain
~AtTabMain()
Definition: AtTabMain.cxx:83
y
const double * y
Definition: lmcurve.cxx:20
AtTrack
Definition: AtTrack.h:25
AtTabMain::InitTab
void InitTab() override
Definition: AtTabMain.cxx:92
cBLUE
constexpr auto cBLUE
Definition: AtTabMain.cxx:58
AtViewerManager::GetEventBranch
DataHandling::AtBranch & GetEventBranch()
Definition: AtViewerManager.h:72
DataHandling
Definition: AtDataObserver.h:4
AtTrack.h
AtPatternEvent.h
ContainerManip::GetPointerVector
std::vector< T * > GetPointerVector(const std::vector< T > &vec)
Definition: AtContainerManip.h:105
AtTabBase
Definition: AtTabBase.h:27
cGREEN
constexpr auto cGREEN
Definition: AtTabMain.cxx:57
cYELLOW
constexpr auto cYELLOW
Definition: AtTabMain.cxx:55
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
AtTabInfoFairRoot< AtEvent >
AtContainerManip.h
AtTabMain::fHitSet
TEvePointSetPtr fHitSet
Definition: AtTabMain.h:40
DataHandling::AtSubject::Attach
void Attach(AtObserver *observer)
Attach an observer to get notified when this subject changes.
Definition: AtDataSubject.cxx:12
AtMap.h
DataHandling::AtSubject
Definition: AtDataSubject.h:24
AtTabMain::fEvePatternEvent
TEveEventManagerPtr fEvePatternEvent
Definition: AtTabMain.h:42
AtMap::BinToPad
virtual Int_t BinToPad(Int_t binval)=0
DataHandling::AtSimpleType::Set
void Set(T data, bool notify=true)
Definition: AtDataSubject.h:58
AtPad
Container class for AtPadBase objects.
Definition: AtPad.h:38
AtViewerManager.h
AtTabMain::fThreshold
Int_t fThreshold
Definition: AtTabMain.h:47
AtMap::GeneratePadPlane
virtual void GeneratePadPlane()=0
AtTabMain::SelectPad
static void SelectPad()
Definition: AtTabMain.cxx:440
AtTabMain::fPadWave
TH1I * fPadWave
Definition: AtTabMain.h:56
AtMap::GetPadRef
AtPadReference GetPadRef(int padNum) const
Definition: AtMap.cxx:281
AtTabMain::fHitAttr
TAttMarker fHitAttr
Definition: AtTabMain.h:50
AtTabMain::Update
void Update(DataHandling::AtSubject *sub) override
Definition: AtTabMain.cxx:139
AtTrack::GetHitArray
HitVector & GetHitArray()
Definition: AtTrack.h:79
AtTabInfo.h
AtTabBase::fTabInfo
std::unique_ptr< AtTabInfo > fTabInfo
Definition: AtTabBase.h:33
AtTabMain::fEveEvent
TEveEventManagerPtr fEveEvent
Definition: AtTabMain.h:39
AtTabMain::fPadNum
DataHandling::AtPadNum * fPadNum
Definition: AtTabMain.h:57