ATTPCROOT  0.3.0-alpha
A ROOT-based framework for analyzing data from active target detectors
AtTabFF.cxx
Go to the documentation of this file.
1 #include "AtTabFF.h"
2 
3 #include "AtContainerManip.h" // for GetPointerVector
4 #include "AtE12014.h"
5 #include "AtMCResult.h"
6 #include "AtViewerManager.h"
7 #include "AtViewerManagerSubject.h" // for AtBranch, AtTreeEntry
8 
9 #include <FairLogger.h> // for Logger, LOG
10 #include <FairRootManager.h> // for FairRootManager
11 
12 #include <TCanvas.h>
13 #include <TClonesArray.h> // for TClonesArray
14 #include <TH1.h> // for TH1F
15 #include <THStack.h> // for THStack
16 #include <TObject.h> // for TObject
17 #include <TString.h>
18 
19 #include <ostream> // for endl
20 namespace DataHandling {
21 class AtSubject;
22 }
23 
24 AtTabFF::AtTabFF(DataHandling::AtBranch &fissionBranch, bool plotADC)
25  : AtTabCanvas("FF", 2, 2), fEvent(AtViewerManager::Instance()->GetEventBranch()),
26  fRawEvent(AtViewerManager::Instance()->GetRawEventBranch()), fFissionEvent(fissionBranch),
27  fEntry(AtViewerManager::Instance()->GetCurrentEntry())
28 {
29  fStacks[0] = std::make_unique<THStack>("hFF1", "FF 1");
30  fStacks[2] = std::make_unique<THStack>("hFF2", "FF 2");
31  fStacks[1] = std::make_unique<THStack>("hExp", "Experimental dE/dX curves");
32  if (plotADC)
33  fStacks[3] = std::make_unique<THStack>("hSim", "Experimental ADC Sum");
34  else
35  fStacks[3] = std::make_unique<THStack>("hSim", "Simulated dE/dX curves");
36 
37  std::array<Color_t, 2> colors = {9, 31};
38  for (int i = 0; i < 2; ++i) {
39  fSimdQdZ[i] = std::make_unique<TH1F>(TString::Format("sim%d", i), TString::Format("Sim FF %d", i), 512, 0, 512);
40  fExpdQdZ[i] = std::make_unique<TH1F>(TString::Format("exp%d", i), TString::Format("Exp FF %d", i), 512, 0, 512);
41  fExpADCSum[i] =
42  std::make_unique<TH1F>(TString::Format("expADC%d", i), TString::Format("Exp FF ADC Sum %d", i), 512, 0, 512);
43  ;
44  fSimdQdZ[i]->SetDirectory(nullptr);
45  fExpdQdZ[i]->SetDirectory(nullptr);
46  fExpADCSum[i]->SetDirectory(nullptr);
47 
48  fSimdQdZ[i]->SetLineColor(kRed);
49  fExpdQdZ[i]->SetLineColor(kBlue);
50  fExpADCSum[i]->SetLineColor(kBlue);
51 
52  // FF i should be added to i*2 (they're the first coloumn
53  fStacks[i * 2]->Add(fSimdQdZ[i].get());
54  fStacks[i * 2]->Add(fExpdQdZ[i].get());
55 
56  fStacks[1]->Add(fExpdQdZ[i].get());
57  if (plotADC)
58  fStacks[3]->Add(fExpADCSum[i].get());
59  else
60  fStacks[3]->Add(fSimdQdZ[i].get());
61  }
62 
63  fEntry.Attach(this);
65 }
66 
68 {
69  fEntry.Detach(this);
71 }
72 
74 {
75  UpdateEvent();
76  DrawCanvas();
77 }
78 
79 std::vector<AtHit *> AtTabFF::GetFragmentHits(AtEvent *event)
80 {
81 
82  return {};
83 }
84 
86 {
87  if (fEvent.Get() == nullptr)
88  return;
89 
90  for (int i = 0; i < 2; ++i) {
91  std::vector<double> exp;
92  std::vector<double> sim;
93  std::vector<double> adcSum;
94  auto tbMin = E12014::fTBMin;
95  E12014::fTBMin = 0;
96  int goodHits = E12014::FillHitSums(exp, sim, fFissionEvent->GetFragHits(i),
98  nullptr, &adcSum, fRawEvent.Get());
99  LOG(info) << "Good hits " << goodHits;
100  /*
101  using namespace SampleConsensus;
102  AtSampleConsensus sac(Estimators::kRANSAC, AtPatterns::PatternType::kLine, RandomSample::SampleMethod::kUniform);
103  auto patEvent = sac.Solve(fEvent.Get());
104 
105  E12014::FillSimSum(sim, ContainerManip::GetPointerVector(patEvent.GetTrackCand()[1].GetHitArray()));
106  */
107  E12014::fTBMin = tbMin;
111  auto fResultArray = dynamic_cast<TClonesArray *>(FairRootManager::Instance()->GetObject("AtMCResult"));
112  if (fResultArray) {
113  auto result = dynamic_cast<MCFitter::AtMCResult *>(fResultArray->At(0));
114  fSimdQdZ[i]->Scale(result->fParameters["Amp"]);
115  if (i == 0) {
116  LOG(info) << "Z1: " << result->fParameters["Z0"] << " A1: " << result->fParameters["A0"];
117  LOG(info) << "Z2: " << result->fParameters["Z1"] << " A2: " << result->fParameters["A1"];
118  LOG(info) << "ObjQ: " << result->fParameters["ObjQ"] << " ObjPos: " << result->fParameters["ObjPos"]
119  << " Amp: " << result->fParameters["Amp"];
120  bool inCut = result->fParameters["ObjQ"] < 15 && result->fParameters["Amp"] < 0.5;
121  LOG(info) << std::endl << (inCut ? "In Cut" : "Out Cut") << std::endl;
122  }
123  }
124  }
125 }
126 
128 {
129  // Loop through each stack. If its an even
130  for (int i = 0; i < 4; ++i) {
131  fCanvas->cd(i + 1);
132  // We are drawing a sim vs experiment
133  if (i % 2 == 0) {
134  fSimdQdZ[i / 2]->SetLineColor(kRed);
135  fExpdQdZ[i / 2]->SetLineColor(kBlue);
136  fStacks[i]->Draw("nostack;hist");
137  } else if (i == 1) { // Drawing experimental stack
138  // fExpdQdZ[0]->SetLineColor(9);
139  // fExpdQdZ[1]->SetLineColor(31);
140  fStacks[i]->Draw("nostack;hist");
141  } else if (i == 3) { // Drawing simulation stack
142  // fSimdQdZ[0]->SetLineColor(9);
143  // fSimdQdZ[1]->SetLineColor(31);
144  fStacks[i]->Draw("nostack;hist");
145  }
146  }
147  UpdateCanvas();
148 }
AtTabCanvas::fCanvas
TCanvas * fCanvas
Definition: AtTabCanvas.h:22
AtTabFF::fEvent
AtTabInfoFairRoot< AtEvent > fEvent
Definition: AtTabFF.h:45
ClassImp
ClassImp(AtTabFF)
AtTabFF::Update
void Update(DataHandling::AtSubject *sub) override
Definition: AtTabFF.cxx:73
DataHandling::AtBranch
Subject for the branch in the FairRoot tree.
Definition: AtViewerManagerSubject.h:38
AtTabFF::fFissionEvent
AtTabInfoFairRoot< AtFissionEvent > fFissionEvent
Definition: AtTabFF.h:47
AtTabFF::fExpdQdZ
std::array< TH1Ptr, 2 > fExpdQdZ
Definition: AtTabFF.h:51
AtTabFF::DrawCanvas
void DrawCanvas()
Definition: AtTabFF.cxx:127
AtMCResult.h
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
AtTabFF
Definition: AtTabFF.h:37
AtEvent
Definition: AtEvent.h:22
AtE12014.h
AtTabFF::fRawEvent
AtTabInfoFairRoot< AtRawEvent > fRawEvent
Definition: AtTabFF.h:46
AtTabInfoFairRoot::Get
T * Get()
Definition: AtTabInfo.h:102
AtViewerManagerSubject.h
E12014::fTBMin
static int fTBMin
Definition: AtE12014.h:25
AtTabFF::UpdateEvent
void UpdateEvent()
Definition: AtTabFF.cxx:85
AtTabFF::AtTabFF
AtTabFF(DataHandling::AtBranch &fissionBranch, bool plotADC=false)
Definition: AtTabFF.cxx:24
MCFitter::AtMCResult
Definition: AtMCResult.h:18
AtTabFF::fEntry
DataHandling::AtTreeEntry & fEntry
Definition: AtTabFF.h:48
AtEvent::GetHits
const HitVector & GetHits() const
Definition: AtEvent.h:116
AtViewerManager::GetEventBranch
DataHandling::AtBranch & GetEventBranch()
Definition: AtViewerManager.h:72
DataHandling
Definition: AtDataObserver.h:4
ContainerManip::GetPointerVector
std::vector< T * > GetPointerVector(const std::vector< T > &vec)
Definition: AtContainerManip.h:105
ContainerManip::SetHistFromData
void SetHistFromData(TH1 &hist, const T &data, double xMin=0, double xMax=0)
Use contents of data to set the bin contents of hist.
Definition: AtContainerManip.h:24
AtTabFF::GetFragmentHits
std::vector< AtHit * > GetFragmentHits(AtEvent *event)
Definition: AtTabFF.cxx:79
AtTabCanvas
Abstract class for a tab composed of a single TCanvas.
Definition: AtTabCanvas.h:20
AtViewerManager
Definition: AtViewerManager.h:33
AtContainerManip.h
AtTabFF::fExpADCSum
std::array< TH1Ptr, 2 > fExpADCSum
Definition: AtTabFF.h:52
DataHandling::AtSubject::Attach
void Attach(AtObserver *observer)
Attach an observer to get notified when this subject changes.
Definition: AtDataSubject.cxx:12
DataHandling::AtSubject
Definition: AtDataSubject.h:24
AtTabFF::~AtTabFF
~AtTabFF()
Definition: AtTabFF.cxx:67
AtTabCanvas::UpdateCanvas
void UpdateCanvas()
Definition: AtTabCanvas.cxx:24
E12014::FillHitSums
static int FillHitSums(std::vector< double > &exp, std::vector< double > &sim, const std::vector< AtHit * > &expHits, const std::vector< AtHit * > &simHits, int threshold=0, float saturationThreshold=std::numeric_limits< float >::max(), const AtDigiPar *par=nullptr, std::vector< double > *expADC=nullptr, AtRawEvent *expEvent=nullptr)
Definition: AtE12014.cxx:102
AtViewerManager.h
AtTabFF::fSimdQdZ
std::array< TH1Ptr, 2 > fSimdQdZ
Definition: AtTabFF.h:50
AtTabFF::fStacks
std::array< THStackPtr, 4 > fStacks
Definition: AtTabFF.h:54
E12014::fSatThreshold
static double fSatThreshold
Definition: AtE12014.h:27
AtFissionEvent::GetFragHits
std::vector< AtHit * > GetFragHits() const
Definition: AtFissionEvent.cxx:86
AtTabFF.h