ATTPCROOT  0.3.0-alpha
A ROOT-based framework for analyzing data from active target detectors
AtTabFission.cxx
Go to the documentation of this file.
1 #include "AtTabFission.h"
2 
3 #include "AtFissionEvent.h" // for AtFissionEvent
4 #include "AtPattern.h" // for AtPattern
5 #include "AtTabInfo.h" // for AtTabInfoFairRoot, AtTabInfo
6 #include "AtTrack.h" // for AtTrack
7 
8 #include <FairLogger.h> // for LOG
9 
10 #include <TAttMarker.h> // for TAttMarker
11 #include <TEveElement.h> // for TEveElement
12 #include <TEveEventManager.h> // for TEveEventManager
13 #include <TEveManager.h> // for TEveManager, gEve
14 #include <TEvePointSet.h> // for TEvePointSet
15 
16 #include <array> // for array
17 #include <utility> // for move
18 namespace DataHandling {
19 class AtSubject;
20 }
21 
23 
25  : AtTabMain(),
26  fCorrHitSet({std::make_unique<TEvePointSet>("Corrected Frag 0"),
27  std::make_unique<TEvePointSet>("Corrected Frag 1"), std::make_unique<TEvePointSet>("Corrected Beam")})
28 {
29 }
31 {
33  gEve->AddEvent(fEveFissionEvent.get());
34 
35  fFissionEventBranch.SetBranchName("AtFissionEvent");
36  auto fissionInfo = std::make_unique<AtTabInfoFairRoot<AtFissionEvent>>(fFissionEventBranch);
37  fTabInfo->AddAugment(std::move(fissionInfo));
38 
39  fUncorrHitSet->SetDestroyOnZeroRefCnt(false);
40  for (auto &set : fCorrHitSet)
41  set->SetDestroyOnZeroRefCnt(false);
42 }
43 
45 {
46  if (sub == fEntry)
47  UpdateFissionElements();
48 
49  AtTabMain::Update(sub);
50 }
51 
52 void AtTabFission::UpdateFissionElements()
53 {
54  if (fEveFissionEvent == nullptr)
55  return;
56 
58  auto fissionEvent = GetFairRootInfo<AtFissionEvent>();
59 
60  if (fissionEvent == nullptr) {
61  LOG(debug) << "Cannot update fission event: no event availible.";
62  return;
63  }
64 
65  // Grab the pattern
66  try {
67  auto &track = fissionEvent->GetYTrack();
68 
69  fEveFissionEvent->RemoveElements();
70 
71  // Add the corrected hits to the viewer
72  for (int i = 0; i < 2; ++i) {
73  auto hitSet = fCorrHitSet[i].get();
74  SetPointsFromHits(*hitSet, fissionEvent->GetFragHitsCorr(i));
75  fHitAttr.Copy(*hitSet);
76  hitSet->SetMarkerColor(GetTrackColor(i));
77  fEveFissionEvent->AddElement(hitSet);
78  }
79  auto hitSet = fCorrHitSet[2].get();
80  fHitAttr.Copy(*hitSet);
81  SetPointsFromHits(*hitSet, fissionEvent->GetBeamHitsCorr());
82  hitSet->SetMarkerColor(GetTrackColor(2));
83  fEveFissionEvent->AddElement(hitSet);
84 
85  // Add the uncorrected hits to the viewer
86  hitSet = fUncorrHitSet.get();
87  fHitAttr.Copy(*hitSet);
88  SetPointsFromHits(*hitSet, fissionEvent->GetFragHits());
89  fEveFissionEvent->AddElement(fUncorrHitSet.get());
90 
91  // Add the pattern to the viewer
92  auto pattern = track.GetPattern()->GetEveElement();
93  pattern->SetDestroyOnZeroRefCnt(false);
94  pattern->SetMainColor(kRed);
95  pattern->FindChild("beam")->SetRnrState(false);
96  fEveFissionEvent->AddElement(pattern);
97  } catch (...) {
98  fEveFissionEvent->RemoveElements();
99  }
100 }
101 
103 {
104  fEveEvent->SetRnrState(false);
105  fEvePatternEvent->SetRnrState(false);
106  fEveFissionEvent->SetRnrState(true);
107  fUncorrHitSet->SetRnrState(false);
108 }
AtTabMain
Main tab in viewer for 3D and pad selection.
Definition: AtTabMain.h:34
AtTabMain::fEntry
DataHandling::AtTreeEntry * fEntry
Definition: AtTabMain.h:61
AtFissionEvent.h
AtTabFission.h
AtPattern.h
ClassImp
ClassImp(AtTabFission)
DataHandling::AtSubject::Notify
void Notify()
Notify all attached subjects that something changed.
Definition: AtDataSubject.cxx:7
AtTabMain::SetPointsFromHits
void SetPointsFromHits(TEvePointSet &hitSet, const std::vector< std::unique_ptr< AtHit >> &hits)
Definition: AtTabMain.cxx:290
AtTabFission::fUncorrHitSet
TEvePointSetPtr fUncorrHitSet
Definition: AtTabFission.h:26
AtTabMain::GetTrackColor
Color_t GetTrackColor(int i)
Definition: AtTabMain.cxx:129
AtTabMain::InitTab
void InitTab() override
Definition: AtTabMain.cxx:92
DataHandling
Definition: AtDataObserver.h:4
AtTrack.h
AtTabFission::fEveFissionEvent
TEveEventManagerPtr fEveFissionEvent
Definition: AtTabFission.h:24
AtTabFission::fCorrHitSet
std::array< TEvePointSetPtr, 3 > fCorrHitSet
Definition: AtTabFission.h:27
AtTabFission
Tab for hangling fission events in viewer.
Definition: AtTabFission.h:22
AtTabFission::UpdateRenderState
void UpdateRenderState() override
Definition: AtTabFission.cxx:102
DataHandling::AtSubject
Definition: AtDataSubject.h:24
AtTabMain::fEvePatternEvent
TEveEventManagerPtr fEvePatternEvent
Definition: AtTabMain.h:42
AtTabFission::fFissionEventBranch
DataHandling::AtBranch fFissionEventBranch
Definition: AtTabFission.h:29
AtTabFission::AtTabFission
AtTabFission()
Definition: AtTabFission.cxx:24
DataHandling::AtBranch::SetBranchName
void SetBranchName(const TString &name)
Will notify on change.
Definition: AtViewerManagerSubject.cxx:17
AtTabMain::fHitAttr
TAttMarker fHitAttr
Definition: AtTabMain.h:50
AtTabMain::Update
void Update(DataHandling::AtSubject *sub) override
Definition: AtTabMain.cxx:139
AtTabFission::InitTab
void InitTab() override
Definition: AtTabFission.cxx:30
AtTabInfo.h
AtTabFission::Update
void Update(DataHandling::AtSubject *sub) override
Definition: AtTabFission.cxx:44
AtTabBase::fTabInfo
std::unique_ptr< AtTabInfo > fTabInfo
Definition: AtTabBase.h:33
AtTabMain::fEveEvent
TEveEventManagerPtr fEveEvent
Definition: AtTabMain.h:39