ATTPCROOT  0.3.0-alpha
A ROOT-based framework for analyzing data from active target detectors
AtEventDrawTask.cxx
Go to the documentation of this file.
1 
7 #include "AtEventDrawTask.h"
8 // IWYU pragma: no_include <ext/alloc_traits.h>
9 #include "AtAuxPad.h" // for AtAuxPad
10 #include "AtBaseEvent.h" // for AtBaseEvent::AuxPadMap
11 #include "AtEvent.h" // for AtEvent, hitVector
12 #include "AtEventManager.h" // for AtEventManager
13 #include "AtFindVertex.h" //for vertex
14 #include "AtHit.h" // for AtHit
15 #include "AtHitCluster.h" // for AtHitCluster
16 #include "AtMap.h" // for AtMap
17 #include "AtPad.h" // for AtPad
18 #include "AtPadReference.h"
19 #include "AtPattern.h"
20 #include "AtPatternEvent.h" // for AtPatternEvent
21 #include "AtRawEvent.h" // for AtRawEvent, AuxPadMap
22 #include "AtTrack.h" // for AtTrack, operator<<
23 
24 #include <FairLogger.h> // for Logger, LOG
25 #include <FairRootManager.h> // for FairRootManager
26 
27 #include <Math/Point3D.h> // for PositionVector3D
28 #include <TAttMarker.h> // for kFullDotMedium
29 #include <TAxis.h> // for TAxis
30 #include <TCanvas.h> // for TCanvas
31 #include <TClonesArray.h> // for TClonesArray
32 #include <TColor.h> // for TColor
33 #include <TEveBoxSet.h> // for TEveBoxSet, TEveBoxSet::kBT_AABox
34 #include <TEveElement.h> // for TEveElement
35 #include <TEveLine.h> // for TEveLine
36 #include <TEveManager.h> // for TEveManager, gEve
37 #include <TEvePointSet.h> // for TEvePointSet
38 #include <TEveRGBAPalette.h> // for TEveRGBAPalette
39 #include <TEveTrans.h> // for TEveTrans
40 #include <TEveTreeTools.h> // for TEvePointSelectorConsumer, TEvePoint...
41 #include <TGraph.h> // for TGraph
42 #include <TH1.h> // for TH1D, TH1I, TH1F
43 #include <TH2.h> // for TH2F
44 #include <TH2Poly.h> // for TH2Poly
45 #include <TH3.h> // for TH3F
46 #include <TList.h> // for TList
47 #include <TNamed.h> // for TNamed
48 #include <TObject.h> // for TObject
49 #include <TPad.h> // for TPad
50 #include <TPaletteAxis.h> // for TPaletteAxis
51 #include <TROOT.h> // for TROOT, gROOT
52 #include <TRandom.h> // for TRandom
53 #include <TSeqCollection.h> // for TSeqCollection
54 #include <TStyle.h> // for TStyle, gStyle
55 #include <TVirtualPad.h> // for TVirtualPad, gPad
56 #include <TVirtualX.h> // for TVirtualX
57 
58 #include "S800Ana.h"
59 #include "S800Calc.h" // for S800Calc, CRDC, MultiHitTOF, IC
60 
61 #include <algorithm> // for max
62 #include <array> // for array
63 #include <cstdio> // for sprintf
64 #include <exception> // for exception
65 #include <iostream> // for cout
66 #include <map> // for operator!=, _Rb_tree_const_iterator
67 #include <string> // for allocator, operator<<
68 #include <utility> // for pair
69 
70 constexpr auto cRED = "\033[1;31m";
71 constexpr auto cYELLOW = "\033[1;33m";
72 constexpr auto cNORMAL = "\033[0m";
73 constexpr auto cGREEN = "\033[1;32m";
74 constexpr auto cBLUE = "\033[1;34m";
75 
76 using namespace std;
77 
79 
81  : fIs2DPlotRange(kFALSE), fUnpackHough(kFALSE), fEventArray(nullptr), fEventManager(nullptr), fRawevent(nullptr),
82  fDetmap(nullptr), fThreshold(0), fHitSet(nullptr), fCorrectedHitSet(nullptr), fhitBoxSet(nullptr),
83  fPadPlanePal(nullptr), fHitColor(kPink), fHitSize(1), fHitStyle(kFullDotMedium), fCvsPadPlane(nullptr),
84  fPadPlane(nullptr), fCvsPadWave(nullptr), fPadWave(nullptr), fCvsPadAll(nullptr), fCvsQEvent(nullptr),
85  fQEventHist(nullptr), fQEventHist_H(nullptr), fCvsHoughSpace(nullptr), fHoughSpace(nullptr),
86  fCvsRhoVariance(nullptr), fRhoVariance(nullptr), fCvsPhi(nullptr), fCvsMesh(nullptr), fMesh(nullptr),
87  fCvs3DHist(nullptr), f3DHist(nullptr), fCvsRad(nullptr), fRadVSTb(nullptr), fCvsTheta(nullptr), fTheta(nullptr),
88  fAtMapPtr(nullptr), fMinZ(0), fMaxZ(1344), fMinX(432), fMaxX(-432), f3DHitStyle(0), fMultiHit(0),
89  fSaveTextData(false), f3DThreshold(0), fRawEventBranchName("AtRawEvent"), fEventBranchName("AtEventH"),
90  fPatternEventBranchName("AtPatternEvent"), fVertexSize(1.5), fVertexStyle(kFullCircle), fCvsPID(nullptr),
91  fPID(nullptr), fCvsPID2(nullptr), fPID2(nullptr), fS800Calc(nullptr)
92 {
93 
94  Char_t padhistname[256];
95 
96  for (Int_t i = 0; i < fNumPads; i++) {
97  sprintf(padhistname, "pad_%d", i);
98  fPadAll[i] = new TH1I(padhistname, padhistname, 512, 0, 511);
99  }
100 
101  Char_t phihistname[256];
102 
103  for (Int_t i = 0; i < 5; i++) {
104  sprintf(phihistname, "PhiDistr_%d", i);
105  fPhiDistr[i] = new TH1D(phihistname, phihistname, 180.0, -180.0, 180.0);
106  if (i == 0)
107  fPhiDistr[i]->SetLineColor(kRed);
108  else if (i == 1)
109  fPhiDistr[i]->SetLineColor(kBlue);
110  else if (i == 2)
111  fPhiDistr[i]->SetLineColor(kGreen);
112  else if (i == 3)
113  fPhiDistr[i]->SetLineColor(kCyan);
114  else if (i == 4)
115  fPhiDistr[i]->SetLineColor(kMagenta);
116  fPhiDistr[i]->SetLineWidth(2);
117  fPhiDistr[i]->GetYaxis()->SetRangeUser(0., 20.);
118  }
119 
120  fIsRawData = kFALSE;
121 
122  fIniHit = new AtHit();
123  fIniHitRansac = new AtHit();
124  fTrackNum = 0;
125 
127 
128  fRGBAPalette = new TEveRGBAPalette(0, 4096);
129 }
130 
132 
134 {
135 
136  std::cout << " ===== AtEventDrawTask::Init =====" << std::endl;
137 
138  gROOT->Reset();
139  FairRootManager *ioMan = FairRootManager::Instance();
141 
142  if (fDetmap == nullptr) {
143  LOG(fatal) << "Map was never set using the function SetMap() in AtEventDrawTask!";
144  return kFATAL;
145  }
146 
147  fDetmap->SetName("fMap");
148  gROOT->GetListOfSpecials()->Add(fDetmap.get());
149 
150  fEventArray = dynamic_cast<TClonesArray *>(ioMan->GetObject(fEventBranchName));
151  if (fEventArray)
152  LOG(INFO) << cGREEN << "Event Array Found in branch " << fEventBranchName << "." << cNORMAL << std::endl;
153 
154  fCorrectedEventArray = dynamic_cast<TClonesArray *>(ioMan->GetObject(fCorrectedEventBranchName));
156  LOG(INFO) << cGREEN << "Corrected Event Array Found in branch " << fCorrectedEventBranchName << "." << cNORMAL
157  << std::endl;
158 
159  fRawEventArray = dynamic_cast<TClonesArray *>(ioMan->GetObject(fRawEventBranchName));
160  if (fRawEventArray) {
161  LOG(INFO) << cGREEN << "Raw Event Array Found in branch " << fRawEventBranchName << "." << cNORMAL << std::endl;
162  fIsRawData = kTRUE;
163  }
164 
165  fPatternEventArray = dynamic_cast<TClonesArray *>(ioMan->GetObject(fPatternEventBranchName));
166  if (fPatternEventArray)
167  LOG(INFO) << cGREEN << "Pattern Event Array Found in branch " << fPatternEventBranchName << "." << cNORMAL
168  << std::endl;
169 
170  fS800Calc = dynamic_cast<S800Calc *>(ioMan->GetObject("s800cal"));
171  if (fS800Calc)
172  LOG(INFO) << cGREEN << "S800Calc Found." << cNORMAL;
173 
174  gStyle->SetPalette(55);
175 
177  fCvsPadWave->SetName("fCvsPadWave");
178  gROOT->GetListOfSpecials()->Add(fCvsPadWave);
179  DrawPadWave();
180  fCvsPadPlane = fEventManager->GetCvsPadPlane(); // There is a problem if the pad plane is drawn first
181  fCvsPadPlane->ToggleEventStatus();
182  fCvsPadPlane->AddExec("ex", "AtEventDrawTask::SelectPad(\"fRawEvent\")");
183  DrawPadPlane();
185  DrawPadAll();
187  Draw3DHist();
188  fCvsQEvent = new TCanvas("fCvsQEvent", "fCvsQEvent");
189  DrawQEvent();
190  fCvsRhoVariance = new TCanvas("fCvsRhoVariance", "fCvsRhoVariance");
191  DrawRhoVariance();
193  DrawPhiReco();
195  DrawMesh();
197  DrawRad();
199  DrawTheta();
201  DrawThetaxPhi();
204  DrawMC();
205  /*fCvsQuadrant1 = fEventManager->GetCvsQuadrant1();
206  fCvsQuadrant2 = fEventManager->GetCvsQuadrant2();
207  fCvsQuadrant3 = fEventManager->GetCvsQuadrant3();
208  fCvsQuadrant4 = fEventManager->GetCvsQuadrant4();
209  DrawHoughSpaceProto();*/
211  DrawAux();
212 
214  DrawPID();
216  DrawPID2();
217 
218  S800Ana s800Ana;
219  s800Ana = fEventManager->GetS800Ana(); // MTDCRanges must be set before calling AtEventDrawTask::Intit()
220  fMTDCXfRange = s800Ana.GetMTDCXfRange();
221  fMTDCObjRange = s800Ana.GetMTDCObjRange();
222  fTofObjCorr = s800Ana.GetTofObjCorr();
223 
224  //******* NO CALLS TO TCANVAS BELOW THIS ONE
226  DrawHoughSpace();
227 
228  std::cout << " AtEventDrawTask::Init : Initialization complete! "
229  << "\n";
230  return kSUCCESS;
231 }
232 
233 void AtEventDrawTask::Exec(Option_t *option)
234 {
235  Reset();
236  ResetPadAll();
237  ResetPhiDistr();
238 
239  // if (fRawEventArray)
240  // DrawPadAll();
241 
242  if (fEventArray) {
243  DrawHitPoints();
244  }
245  if (fS800Calc) {
246  DrawS800();
247  }
248 
249  gEve->Redraw3D(kFALSE);
250 
251  UpdateCvsPadPlane();
252  UpdateCvsPadWave();
253  UpdateCvsPadAll();
254  UpdateCvsQEvent();
255  UpdateCvsRhoVariance();
256  UpdateCvsPhi();
257  UpdateCvsMesh();
258  UpdateCvs3DHist();
259  // UpdateCvsPID();
260  // UpdateCvsPID2();
262  UpdateCvsHoughSpace();
263  UpdateCvsRad();
264  UpdateCvsTheta();
265  UpdateCvsThetaxPhi();
266  UpdateCvsMC();
267  // UpdateCvsQuadrants();
268  }
269 }
270 void AtEventDrawTask::DrawAuxChannels()
271 {
272  for (auto &fAuxChannel : fAuxChannels)
273  fAuxChannel->Reset(nullptr);
274 
275  if (fIsRawData) {
276  fRawevent = dynamic_cast<AtRawEvent *>(fRawEventArray->At(0));
277  fRawevent->SetName("fRawEvent");
278  gROOT->GetListOfSpecials()->Add(fRawevent);
279 
280  // Draw aux channels
281  int numAux = 0;
282 
283  for (auto &padIt : fRawevent->GetAuxPads()) {
284  const AtAuxPad &pad = padIt.second;
285  if (numAux < 9) {
286  std::cout << cYELLOW << " Auxiliary Channel " << numAux << " - Name " << pad.GetAuxName() << cNORMAL
287  << std::endl;
288  auto adc = pad.GetADC();
289  for (int i = 0; i < 512; ++i)
290  fAuxChannels[numAux]->SetBinContent(i, adc[i]);
291  numAux++;
292  }
293  }
294  if (numAux + 1 == 9)
295  std::cout << cYELLOW << "Warning: More auxiliary channels than expected (max 9)" << cNORMAL << std::endl;
296  }
297 }
298 
299 void AtEventDrawTask::DrawRecoHits()
300 {
301  if (!fPatternEventArray)
302  return;
303 
304  auto *patternEvent = dynamic_cast<AtPatternEvent *>(fPatternEventArray->At(0));
305  if (!patternEvent) {
306  LOG(info) << "No pattern event";
307  return;
308  }
309 
310  auto TrackCand = patternEvent->GetTrackCand();
311  fHitLine.clear();
312  fVertex.clear();
313 
314  if (fDrawVertexFromLines) {
315  AtFindVertex findVtx(12);
316  // findVtx.FindVertexMultipleLines(TrackCand, 2);
317  findVtx.FindVertex(TrackCand, fMinTracksPerVertex);
318  std::vector<tracksFromVertex> tv;
319  tv = findVtx.GetTracksVertex();
320 
321  for (size_t i = 0; i < tv.size(); i++) {
322  fVertex.push_back(new TEvePointSet(Form("Vertex_%zu", i), 0, TEvePointSelectorConsumer::kTVT_XYZ));
323  fVertex[i]->SetMarkerColor(kViolet);
324  fVertex[i]->SetMarkerSize(fVertexSize);
325  fVertex[i]->SetMarkerStyle(fVertexStyle);
326  fVertex[i]->SetNextPoint(tv.at(i).vertex.X() / 10., tv.at(i).vertex.Y() / 10., tv.at(i).vertex.Z() / 10.);
327  }
328  }
329  for (Int_t i = 0; i < TrackCand.size(); i++) {
330 
331  Color_t trackColor = GetTrackColor(i) + 1;
332 
333  // Get the line to draw
334  AtTrack track = TrackCand.at(i);
335 
336  if (track.GetPattern() != nullptr && track.GetPattern()->GetEveElement() != nullptr) {
337  fHitLine.push_back(track.GetPattern()->GetEveElement());
338  fHitLine.back()->SetMainColor(trackColor);
339  fHitLine.back()->SetElementName(Form("line_%i", (int)fHitLine.size() - 1));
340  }
341 
342  std::vector<AtHit> trackHits = track.GetHitArrayObject();
343 
344  /* fHitSetTFHC.push_back(
345  std::make_unique<TEvePointSet>(Form("HitMC_%d", i), 0, TEvePointSelectorConsumer::kTVT_XYZ));
346  */
347  fHitSetTFHC.push_back(new TEvePointSet(Form("HitMC_%d", i), 0, TEvePointSelectorConsumer::kTVT_XYZ));
348  fHitSetTFHC[i]->SetMarkerColor(trackColor);
349  fHitSetTFHC[i]->SetMarkerSize(fHitSize);
350  fHitSetTFHC[i]->SetMarkerStyle(fHitStyle);
351 
352  for (auto &trackHit : trackHits) {
353  auto position = trackHit.GetPosition();
354  fHitSetTFHC[i]->SetNextPoint(position.X() / 10., position.Y() / 10., position.Z() / 10.);
355  }
356 
357  std::vector<AtHitCluster> *hitClusters = track.GetHitClusterArray();
358  // fHitClusterSet.push_back(std::make_unique<TEveBoxSet>(Form("HitCluster_%d", i)));
359  fHitClusterSet.push_back(new TEveBoxSet(Form("HitCluster_%d", i)));
360  fHitClusterSet[i]->Reset(TEveBoxSet::kBT_AABox, kFALSE, 64);
361  // fHitClusterSet[i]->SetPalette(fRGBAPalette);
362  // fHitClusterSet[i]->DigitValue(2000);
363 
364  if (hitClusters->size() > 0) {
365 
366  for (auto hitCluster : *hitClusters) {
367  auto clusPos = hitCluster.GetPosition();
368  double boxSize = 0.6;
369 
370  fHitClusterSet[i]->AddBox(clusPos.X() / 10. - boxSize / 2.0, clusPos.Y() / 10. - boxSize / 2.0,
371  clusPos.Z() / 10. - boxSize / 2.0, boxSize, boxSize, boxSize);
372  }
373  }
374 
375  fHitClusterSet[i]->UseSingleColor();
376  fHitClusterSet[i]->SetMainColor(trackColor);
377  fHitClusterSet[i]->SetMainTransparency(50);
378 
379  fHitClusterSet[i]->RefitPlex();
380  TEveTrans &tHitClusterBoxPos = fHitClusterSet[i]->RefMainTrans();
381  tHitClusterBoxPos.SetPos(0.0, 0.0, 0.0);
382  } // End loop over tracks
383 
384  // Add noise hits
385  // fHitSetTFHC.push_back(std::make_unique<TEvePointSet>("HitsNoise", 0, TEvePointSelectorConsumer::kTVT_XYZ));
386  fHitSetTFHC.push_back(new TEvePointSet("HitsNoise", 0, TEvePointSelectorConsumer::kTVT_XYZ));
387  fHitSetTFHC.back()->SetMainColor(kRed);
388  fHitSetTFHC.back()->SetMarkerSize(fHitSize);
389  fHitSetTFHC.back()->SetMarkerStyle(fHitStyle);
390 
391  for (auto &hit : patternEvent->GetNoiseHits()) {
392  auto position = hit->GetPosition();
393  fHitSetTFHC.back()->SetNextPoint(position.X() / 10., position.Y() / 10., position.Z() / 10.);
394  }
395 
396  std::cout << cRED << " Found " << TrackCand.size() << " track candidates " << cNORMAL << std::endl;
397 
398  for (auto &elem : fHitSetTFHC)
399  gEve->AddElement(elem);
400  for (auto &elem : fHitClusterSet)
401  gEve->AddElement(elem);
402  for (auto &elem : fHitLine)
403  gEve->AddElement(elem);
404  for (auto &elem : fVertex)
405  gEve->AddElement(elem);
406 }
407 
408 void AtEventDrawTask::DrawHitPoints()
409 {
410  DrawAuxChannels();
411  fMesh->Reset(nullptr);
412 
413  f3DHist->Reset(nullptr);
414  TRandom r(0);
415 
416  std::ofstream dumpEvent;
417  dumpEvent.open("event.dat");
418 
419  fQEventHist_H->Reset(nullptr);
420  auto *event = dynamic_cast<AtEvent *>(fEventArray->At(0));
421 
422  Double_t Qevent = event->GetEventCharge();
423  Double_t RhoVariance = event->GetRhoVariance();
424  auto MeshArray = event->GetMesh();
425  Int_t eventID = event->GetEventID();
426  TString TSevt = " Event ID : ";
427  TString TSpad = " Pad ID : ";
428  dumpEvent << TSevt << eventID << std::endl;
429 
430  if (fEventManager->GetEraseQEvent()) {
431  fQEventHist->Reset();
432  fRhoVariance->Reset();
433  }
434 
435  fQEventHist->Fill(Qevent);
436  fQEventHist_H->Fill(Qevent);
437  fRhoVariance->Fill(RhoVariance);
438 
439  for (Int_t i = 0; i < 512; i++) {
440 
441  fMesh->SetBinContent(i, MeshArray[i]);
442  }
443 
445  DrawRecoHits();
446 
448 
449  fhitBoxSet = new TEveBoxSet("hitBox");
450  fhitBoxSet->Reset(TEveBoxSet::kBT_AABox, kTRUE, 64);
451 
452  Int_t nHits = event->GetNumHits();
453  fHitSet = new TEvePointSet("Hit", nHits, TEvePointSelectorConsumer::kTVT_XYZ);
454  fHitSet->SetOwnIds(kTRUE);
455  fHitSet->SetMarkerColor(fHitColor);
456  fHitSet->SetMarkerSize(fHitSize);
457  fHitSet->SetMarkerStyle(fHitStyle);
458  std::cout << cBLUE << " Number of hits : " << nHits << cNORMAL << std::endl;
459 
460  for (Int_t iHit = 0; iHit < nHits; iHit++) {
461 
462  AtHit hit = *event->GetHits().at(iHit);
463  Int_t PadNumHit = hit.GetPadNum();
464  Int_t PadMultHit = event->GetHitPadMult(PadNumHit);
465 
466  if (hit.GetCharge() < fThreshold)
467  continue;
468  if (PadMultHit > fMultiHit)
469  continue;
470  auto position = hit.GetPosition();
471 
473  fHitSet->SetMarkerColor(fHitColor);
474  fHitSet->SetNextPoint(position.X() / 10., position.Y() / 10., position.Z() / 10.); // Convert into cm
475  fHitSet->SetPointId(new TNamed(Form("Hit %d", iHit), ""));
476  fPadPlane->Fill(position.X(), position.Y(), hit.GetCharge());
477  }
478 
479  if (fIsRawData) {
480  AtPad *RawPad = fRawevent->GetPad(PadNumHit);
481  if (RawPad != nullptr) {
482  auto adc = RawPad->GetADC();
483  for (Int_t i = 0; i < 512; i++) {
484 
486  if (adc[i] > f3DThreshold)
487  f3DHist->Fill(position.X() / 10., position.Y() / 10., i, adc[i]);
488  }
489  }
490  }
491 
492  if (fSaveTextData) {
494  dumpEvent << position.X() << " " << position.Y() << " " << position.Z() << " " << hit.GetTimeStamp() << " "
495  << hit.GetCharge() << std::endl;
496  }
497  }
498 
499  if (fCorrectedEventArray != nullptr) {
500  std::cout << "Adding corrected hits" << std::endl;
501  auto *eventCorr = dynamic_cast<AtEvent *>(fCorrectedEventArray->At(0));
502  fCorrectedHitSet = new TEvePointSet("Hit2", nHits, TEvePointSelectorConsumer::kTVT_XYZ);
503  fCorrectedHitSet->SetOwnIds(kTRUE);
504  fCorrectedHitSet->SetMarkerColor(kBlue);
505  fCorrectedHitSet->SetMarkerSize(fHitSize);
506  fCorrectedHitSet->SetMarkerStyle(fHitStyle);
507 
508  for (Int_t iHit = 0; iHit < nHits; iHit++) {
509  if (eventCorr == nullptr) {
510  std::cout << "Corrected event was empty!" << std::endl;
511  break;
512  }
513  AtHit hit = *eventCorr->GetHits().at(iHit);
514  auto position = hit.GetPosition();
515  if (hit.GetCharge() < fThreshold)
516  continue;
517  fCorrectedHitSet->SetNextPoint(position.X() / 10., position.Y() / 10., position.Z() / 10.); // Convert into cm
518  fCorrectedHitSet->SetPointId(new TNamed(Form("HitCorr %d", iHit), ""));
519  }
520  }
521 
523 
524  fPadPlane->Draw("zcol");
525  gPad->Update();
526  fPadPlanePal = dynamic_cast<TPaletteAxis *>(fPadPlane->GetListOfFunctions()->FindObject("palette"));
527 
528  for (Int_t iHit = 0; iHit < nHits; iHit++) {
529 
530  AtHit hit = *event->GetHits().at(iHit);
531  auto position = hit.GetPosition();
532 
533  if (f3DHitStyle == 0) {
534 
535  Float_t HitBoxYDim = hit.GetCharge() * 0.001;
536  Float_t HitBoxZDim = 0.05;
537  Float_t HitBoxXDim = 0.05;
538 
540  fhitBoxSet->AddBox(position.X() / 10. - HitBoxXDim / 2.0, position.Y() / 10.,
541  position.Z() / 10. - HitBoxZDim / 2.0, HitBoxXDim, HitBoxYDim,
542  HitBoxZDim); // This coordinates are x,y,z in our system
543  }
544 
545  } else if (f3DHitStyle == 1) {
546 
547  Float_t HitBoxYDim = hit.GetCharge() * 0.0002;
548  Float_t HitBoxZDim = hit.GetCharge() * 0.0002;
549  Float_t HitBoxXDim = hit.GetCharge() * 0.0002;
550 
552  fhitBoxSet->AddBox(position.X() / 10. - HitBoxXDim / 2.0, position.Y() / 10. - HitBoxYDim / 2.0,
553  position.Z() / 10. - HitBoxZDim / 2.0, HitBoxXDim, HitBoxYDim,
554  HitBoxZDim); // This coordinates are x,y,z in our system
555  }
556  }
557 
558  Float_t xrgb = 255, yrgb = 0, zrgb = 0;
559  if (fPadPlanePal) {
560 
561  Int_t cHit = fPadPlanePal->GetValueColor(hit.GetCharge());
562  TColor *hitBoxColor = gROOT->GetColor(cHit);
563  hitBoxColor->GetRGB(xrgb, yrgb, zrgb);
564 
565  // std::cout<<" xrgb : "<<xrgb<<std::endl;
566  // std::cout<<" yrgb : "<<yrgb<<std::endl;
567  // std::cout<<" zrgb : "<<zrgb<<std::endl;
568  // std::cout<<fPadPlanePal<<std::endl;
569  }
570 
571  fhitBoxSet->DigitColor(xrgb * 255, yrgb * 255, zrgb * 255, 0);
572  }
573 
575 
576  fhitBoxSet->RefitPlex();
577  TEveTrans &tHitBoxPos = fhitBoxSet->RefMainTrans();
578  tHitBoxPos.SetPos(0.0, 0.0, 0.0);
579 
580  // for(Int_t i=0;i<hitSphereArray.size();i++) gEve->AddElement(hitSphereArray[i]);
581 
582  if (fIsRawData) {
583  Int_t nPads = fRawevent->GetNumPads();
584  std::cout << "Num of pads : " << nPads << std::endl;
585 
586  if (fEventManager->GetDrawAllPad()) {
587 
588  std::cout << "Setting pads for Draw All Pads." << std::endl;
589 
590  for (Int_t iPad = 0; iPad < nPads; iPad++) {
591 
592  AtPad *fPad = fRawevent->GetPad(iPad);
593  // std::cout<<"Pad num : "<<iPad<<" Is Valid? : "<<fPad->GetValidPad()
594  // << std::endl;
595 
596  if (!fPad->GetValidPad())
597  continue;
598 
599  auto rawadc = fPad->GetRawADC();
600  auto adc = fPad->GetADC();
601  // dumpEvent<<TSpad<<fPad->GetPadNum()<<std::endl;
602 
603  for (Int_t j = 0; j < 512; j++)
604  fPadAll[iPad]->SetBinContent(j, adc[j]);
605  }
606  }
607  }
608 
609  // Adding raw data points
611 
612  if (fCorrectedHitSet)
613  gEve->AddElement(fCorrectedHitSet);
614  gEve->AddElement(fHitSet);
615  gEve->AddElement(fhitBoxSet);
616  }
617 
618  dumpEvent.close();
619 }
620 
621 void AtEventDrawTask::DrawS800()
622 {
623  // fS800Calc = dynamic_cast<S800Calc*> (fS800CalcArray->At(0));
624 
625  if (fS800Calc->GetIsInCut()) {
626  S800Ana s800Ana;
627  s800Ana.SetMTDCXfRange(fMTDCXfRange);
628  s800Ana.SetMTDCObjRange(fMTDCObjRange);
629  s800Ana.SetTofObjCorr(fTofObjCorr);
630 
631  s800Ana.Calc(fS800Calc);
632 
633  if (s800Ana.GetObjCorr_ToF() != -999)
634  fPID->Fill(s800Ana.GetObjCorr_ToF(), s800Ana.GetXfObj_ToF());
635  if (s800Ana.GetObjCorr_ToF() != -999)
636  fPID2->Fill(s800Ana.GetObjCorr_ToF(), s800Ana.GetICSum_E());
637  }
638 }
639 
641 {
642  if (fHitSet) {
643  fHitSet->Reset();
644  gEve->RemoveElement(fHitSet, fEventManager);
645  }
646 
647  if (fCorrectedHitSet) {
648  fCorrectedHitSet->Reset();
649  gEve->RemoveElement(fCorrectedHitSet, fEventManager);
650  }
651 
652  if (fhitBoxSet) {
653  fhitBoxSet->Reset();
654  gEve->RemoveElement(fhitBoxSet, fEventManager);
655  }
656 
658 
659  if (fPatternEventArray) {
660  for (auto elem : fHitClusterSet)
661  gEve->RemoveElement(elem, fEventManager);
662  for (auto elem : fHitSetTFHC)
663  gEve->RemoveElement(elem, fEventManager);
664  for (auto elem : fHitLine)
665  gEve->RemoveElement(elem, fEventManager);
666  for (auto elem : fVertex)
667  gEve->RemoveElement(elem, fEventManager);
668  fHitClusterSet.clear();
669  fHitSetTFHC.clear();
670  fHitLine.clear();
671  fVertex.clear();
672  }
673 
674  } // Draw Minimization
675 
676  if (fPadPlane != nullptr)
677  fPadPlane->Reset(nullptr);
678 }
679 
680 void AtEventDrawTask::DrawPadPlane()
681 {
682  if (fPadPlane) {
683  fPadPlane->Reset(nullptr);
684  return;
685  }
686 
687  fDetmap->GeneratePadPlane();
688  fPadPlane = fDetmap->GetPadPlane();
689  fCvsPadPlane->cd();
690  // fPadPlane -> Draw("COLZ L0"); //0 == bin lines adre not drawn
691  fPadPlane->Draw("COL L0");
692  fPadPlane->SetMinimum(1.0);
693  gStyle->SetOptStat(0);
694  gStyle->SetPalette(103);
695  gPad->Update();
696 }
697 
698 void AtEventDrawTask::DrawPadWave()
699 {
700  fPadWave = new TH1I("fPadWave", "fPadWave", 512, 0, 511);
701  gROOT->GetListOfSpecials()->Add(fPadWave);
702  fCvsPadWave->cd();
703  fPadWave->Draw();
704 }
705 
706 void AtEventDrawTask::DrawPadAll()
707 {
708 
709  fCvsPadAll->cd();
710 
711  std::cout << "Starting to draw pads" << std::endl;
712  for (auto &i : fPadAll) {
713  i->GetYaxis()->SetRangeUser(0, 2500);
714  // TODO: make it pad number independent / retrieve the quadrant info
715  i->Draw("SAME");
716  }
717  std::cout << "Finished drawing." << std::endl;
718 }
719 
720 void AtEventDrawTask::DrawQEvent()
721 {
722 
723  fQEventHist = new TH1D("fQEventHist", "fQEventHist", 300, 0., 2000000.);
724  fQEventHist_H = new TH1D("fQEventHist_H", "fQEventHist_H", 300, 0., 2000000.);
725  fQEventHist_H->SetLineColor(kRed);
726  fQEventHist_H->SetFillStyle(1);
727  fCvsQEvent->cd();
728  fQEventHist->Draw();
729  fQEventHist_H->Draw("SAMES");
730 }
731 
732 void AtEventDrawTask::DrawRhoVariance()
733 {
734 
735  fCvsRhoVariance->cd();
736  fRhoVariance = new TH1D("fRhoVariance", "fRhoVariance", 4000, 0., 1000000.);
737  fRhoVariance->Draw();
738  fRhoVariance->SetLineColor(kRed);
739 }
740 
741 void AtEventDrawTask::DrawHoughSpace()
742 {
743  fCvsHoughSpace->cd();
744  fHoughSpace = new TH2F("HistHoughXY", "HistHoughXY", 100, 0, 3.15, 500, -1000, 1000);
745  fHoughSpace->Draw("colz");
746 }
747 
748 void AtEventDrawTask::DrawPhiReco()
749 {
750  fCvsPhi->cd();
751  // fPhiDistr = new TH1D("PhiDist","PhiDist",90.0,0.0,90.0);
752  for (auto &i : fPhiDistr) {
753  i->Draw("SAME");
754  }
755 }
756 
757 void AtEventDrawTask::DrawMesh()
758 {
759 
760  fCvsMesh->cd();
761  fMesh = new TH1F("Mesh", "Mesh", 512, 0, 511);
762  fMesh->Draw();
763 }
764 
765 void AtEventDrawTask::Draw3DHist()
766 {
767 
768  fCvs3DHist->cd();
769  f3DHist = new TH3F("gl3DHist", "gl3DHist", 50, -25.0, 25.0, 50, -25.0, 25.0, 50, 0, 512);
770  gStyle->SetPalette(55);
771  // gStyle->SetCanvasPreferGL(kTRUE);
772 
773  f3DHist->SetFillColor(2);
774  f3DHist->Draw("box");
775  // f3DHist -> Draw("glbox3");
776  // f3DHist -> Draw("glcol"); //TODO: Not working, strange behavior
777 }
778 
779 void AtEventDrawTask::DrawRad()
780 {
781 
782  fCvsRad->cd();
783  fRadVSTb = new TH2F("RadVSTb", "RadVSTb", 100, 0, 512, 100, 0, 250);
784  fRadVSTb->SetMarkerStyle(22);
785  fRadVSTb->SetMarkerColor(kRed);
786  fRadVSTb->Draw();
787 }
788 
789 void AtEventDrawTask::DrawTheta()
790 {
791 
792  fCvsTheta->cd();
793  // fTheta = new TH1F("Theta","Theta",512,0,511);
794  fTheta = new TH2F("Theta", "Theta", 512, 0, 511, 500, 0, 2.0);
795  fTheta->SetMarkerStyle(22);
796  fTheta->SetMarkerColor(kRed);
797  // gStyle->SetErrorX(0);
798  fTheta->Draw("");
799 }
800 
801 void AtEventDrawTask::DrawThetaxPhi()
802 {
803 
804  fCvsThetaxPhi->cd();
805  fThetaxPhi = new TH2F("ThetaxPhi", "ThetaxPhi", 512, 0, 511, 100, -1000, 1000);
806  fThetaxPhi->SetMarkerStyle(22);
807  fThetaxPhi->SetMarkerColor(kRed);
808 
809  fThetaxPhi_Ini_RANSAC = new TH2F("ThetaxPhi_Ini_RANSAC", "ThetaxPhi_Ini_RANSAC", 512, 0, 511, 100, -1000, 1000);
810  fThetaxPhi_Ini_RANSAC->SetMarkerStyle(20);
811  fThetaxPhi_Ini_RANSAC->SetMarkerColor(kGreen);
812 
813  fThetaxPhi_Ini = new TH2F("ThetaxPhi_Ini", "ThetaxPhi_Ini", 512, 0, 511, 100, -1000, 1000);
814  fThetaxPhi_Ini->SetMarkerStyle(23);
815  fThetaxPhi_Ini->SetMarkerColor(kBlue);
816  fThetaxPhi_Ini->SetMarkerSize(1.0);
817 
818  fThetaxPhi->Draw("");
819  fThetaxPhi_Ini->Draw("SAMES");
820  fThetaxPhi_Ini_RANSAC->Draw("SAMES");
821 }
822 
823 void AtEventDrawTask::DrawMC()
824 {
825 
826  fCvsMC_XY->cd();
827 
828  fMC_XY_exp = new TGraph();
829  fMC_XY_exp->SetPoint(1, 0, 0);
830  fMC_XY_exp->SetMarkerStyle(20);
831  fMC_XY_exp->SetMarkerSize(1.0);
832  fMC_XY_exp->SetMarkerColor(kBlack);
833  fMC_XY_exp->Draw("AP");
834 
835  fMC_XY = new TGraph();
836  fMC_XY->SetPoint(1, 0, 0);
837  fMC_XY->SetMarkerStyle(20);
838  fMC_XY->SetMarkerSize(1.0);
839  fMC_XY->SetMarkerColor(kRed);
840  fMC_XY->Draw("P");
841 
842  fMC_XY_back = new TGraph();
843  fMC_XY_back->SetPoint(1, 0, 0);
844  fMC_XY_back->SetMarkerStyle(22);
845  fMC_XY_back->SetMarkerSize(1.0);
846  fMC_XY_back->SetMarkerColor(6);
847  fMC_XY_back->Draw("P");
848 
849  fMC_XY_int = new TGraph();
850  fMC_XY_int->SetPoint(1, 0, 0);
851  fMC_XY_int->SetMarkerStyle(22);
852  fMC_XY_int->SetMarkerSize(1.0);
853  fMC_XY_int->SetMarkerColor(8);
854  fMC_XY_int->Draw("P");
855 
856  fCvsMC_Z->cd();
857 
858  fMC_ZX = new TGraph();
859  fMC_ZX->SetPoint(1, 0, 0);
860  fMC_ZX->SetMarkerStyle(20);
861  fMC_ZX->SetMarkerSize(1.0);
862  fMC_ZX->SetMarkerColor(kRed);
863  fMC_ZX->Draw("AP");
864 
865  fMC_ZY = new TGraph();
866  fMC_ZY->SetPoint(1, 0, 0);
867  fMC_ZY->SetMarkerStyle(20);
868  fMC_ZY->SetMarkerSize(1.0);
869  fMC_ZY->SetMarkerColor(kBlack);
870  fMC_ZY->Draw("P");
871 
872  fMC_ZX_back = new TGraph();
873  fMC_ZX_back->SetPoint(1, 0, 0);
874  fMC_ZX_back->SetMarkerStyle(22);
875  fMC_ZX_back->SetMarkerSize(1.0);
876  fMC_ZX_back->SetMarkerColor(6);
877  fMC_ZX_back->Draw("P");
878 
879  fMC_ZY_back = new TGraph();
880  fMC_ZY_back->SetPoint(1, 0, 0);
881  fMC_ZY_back->SetMarkerStyle(22);
882  fMC_ZY_back->SetMarkerSize(1.0);
883  fMC_ZY_back->SetMarkerColor(6);
884  fMC_ZY_back->Draw("P");
885 
886  fMC_ZX_int = new TGraph();
887  fMC_ZX_int->SetPoint(1, 0, 0);
888  fMC_ZX_int->SetMarkerStyle(22);
889  fMC_ZX_int->SetMarkerSize(1.0);
890  fMC_ZX_int->SetMarkerColor(kRed);
891  fMC_ZX_int->Draw("P");
892 
893  fMC_ZY_int = new TGraph();
894  fMC_ZY_int->SetPoint(1, 0, 0);
895  fMC_ZY_int->SetMarkerStyle(22);
896  fMC_ZY_int->SetMarkerSize(1.0);
897  fMC_ZY_int->SetMarkerColor(kBlack);
898  fMC_ZY_int->Draw("P");
899 }
900 void AtEventDrawTask::DrawAux()
901 {
902  fCvsAux->Divide(3, 3);
903  for (int i = 0; i < 9; ++i) {
904  fAuxChannels[i] = new TH1F(Form("AuxCh%i", i), Form("AuxChannel%i", i), 512, 0, 511);
905  fCvsAux->cd(i + 1);
906  fAuxChannels[i]->Draw();
907  }
908 }
909 
910 void AtEventDrawTask::DrawPID()
911 {
912 
913  fCvsPID->cd();
914  // fPID = new TH2F("PID","PID",3000,-250,500,2000,0,500);
915  fPID = new TH2F("PID", "PID", 200, -150, 50, 300, 150, 450);
916  fPID->Draw("colz");
917 }
918 
919 void AtEventDrawTask::DrawPID2()
920 {
921 
922  fCvsPID2->cd();
923  // fPID2 = new TH2F("PID2","PID2",3000,-250,500,2000,0,500);
924  fPID2 = new TH2F("PID2", "PID2", 200, -150, 50, 300, 150, 450);
925  fPID2->Draw("colz");
926 }
927 
928 void AtEventDrawTask::UpdateCvsPID()
929 {
930  fCvsPID->Modified();
931  fCvsPID->Update();
932 }
933 
934 void AtEventDrawTask::UpdateCvsPID2()
935 {
936  fCvsPID2->Modified();
937  fCvsPID2->Update();
938 }
939 
940 void AtEventDrawTask::UpdateCvsAux()
941 {
942 
943  TPad *Pad_1 = dynamic_cast<TPad *>(fCvsAux->GetPad(1));
944  Pad_1->Modified();
945  Pad_1->Update();
946  TPad *Pad_2 = dynamic_cast<TPad *>(fCvsAux->GetPad(2));
947  Pad_2->Modified();
948  Pad_2->Update();
949  TPad *Pad_3 = dynamic_cast<TPad *>(fCvsAux->GetPad(3));
950  Pad_3->Modified();
951  Pad_3->Update();
952  TPad *Pad_4 = dynamic_cast<TPad *>(fCvsAux->GetPad(4));
953  Pad_4->Modified();
954  Pad_4->Update();
955  TPad *Pad_5 = dynamic_cast<TPad *>(fCvsAux->GetPad(5));
956  Pad_5->Modified();
957  Pad_5->Update();
958  TPad *Pad_6 = dynamic_cast<TPad *>(fCvsAux->GetPad(6));
959  Pad_6->Modified();
960  Pad_6->Update();
961  TPad *Pad_7 = dynamic_cast<TPad *>(fCvsAux->GetPad(7));
962  Pad_7->Modified();
963  Pad_7->Update();
964  TPad *Pad_8 = dynamic_cast<TPad *>(fCvsAux->GetPad(8));
965  Pad_8->Modified();
966  Pad_8->Update();
967  TPad *Pad_9 = dynamic_cast<TPad *>(fCvsAux->GetPad(9));
968  Pad_9->Modified();
969  Pad_9->Update();
970  fCvsAux->Modified();
971  fCvsAux->Update();
972 }
973 void AtEventDrawTask::UpdateCvsPadPlane()
974 {
975  fHoughSpace->Draw("colz");
976  fCvsPadPlane->Modified();
977  fCvsPadPlane->Update();
978 
979  /*if (paxis) {
980  if(fIs2DPlotRange) {
981  paxis -> SetX1NDC(0.940);
982  paxis -> SetX2NDC(0.955);
983  paxis -> SetLabelSize(0.08);
984  paxis -> GetAxis() -> SetTickSize(0.01);
985  } else {
986  paxis -> SetX1NDC(0.905);
987  paxis -> SetX2NDC(0.94);
988  }
989 
990  fCvsPadPlane -> Modified();
991  fCvsPadPlane -> Update();
992  }*/
993 }
994 
995 void AtEventDrawTask::UpdateCvsPadWave()
996 {
997  fCvsPadWave->Modified();
998  fCvsPadWave->Update();
999 
1000  // TPaletteAxis *paxis
1001  // = (TPaletteAxis *) fPadPlane->GetListOfFunctions()->FindObject("palette");
1002 }
1003 
1004 void AtEventDrawTask::UpdateCvsPadAll()
1005 {
1006  fCvsPadAll->Modified();
1007  fCvsPadAll->Update();
1008 
1009  // TPaletteAxis *paxis
1010  // = (TPaletteAxis *) fPadPlane->GetListOfFunctions()->FindObject("palette");
1011 }
1012 
1013 void AtEventDrawTask::UpdateCvsQEvent()
1014 {
1015  fCvsQEvent->Modified();
1016  fCvsQEvent->Update();
1017 }
1018 
1019 void AtEventDrawTask::UpdateCvsRhoVariance()
1020 {
1021  fCvsRhoVariance->Modified();
1022  fCvsRhoVariance->Update();
1023 }
1024 
1025 void AtEventDrawTask::UpdateCvsHoughSpace()
1026 {
1027  fCvsHoughSpace->Modified();
1028  fCvsHoughSpace->Update();
1029 }
1030 
1031 void AtEventDrawTask::UpdateCvsPhi()
1032 {
1033  // if(fPhiDistr!=NULL)fPhiDistr->Draw();
1034  fCvsPhi->Modified();
1035  fCvsPhi->Update();
1036 }
1037 
1038 void AtEventDrawTask::UpdateCvsMesh()
1039 {
1040 
1041  fCvsMesh->Modified();
1042  fCvsMesh->Update();
1043 }
1044 
1045 void AtEventDrawTask::UpdateCvs3DHist()
1046 {
1047 
1048  fCvs3DHist->Modified();
1049  fCvs3DHist->Update();
1050 }
1051 
1052 void AtEventDrawTask::UpdateCvsRad()
1053 {
1054 
1055  fCvsRad->Modified();
1056  fCvsRad->Update();
1057 }
1058 
1059 void AtEventDrawTask::UpdateCvsTheta()
1060 {
1061 
1062  fCvsTheta->Modified();
1063  fCvsTheta->Update();
1064 }
1065 
1066 void AtEventDrawTask::UpdateCvsThetaxPhi()
1067 {
1068 
1069  fCvsThetaxPhi->Modified();
1070  fCvsThetaxPhi->Update();
1071 }
1072 
1073 void AtEventDrawTask::UpdateCvsQuadrants()
1074 {
1075  fCvsQuadrant1->Modified();
1076  fCvsQuadrant1->Update();
1077  fCvsQuadrant2->Modified();
1078  fCvsQuadrant2->Update();
1079  fCvsQuadrant3->Modified();
1080  fCvsQuadrant3->Update();
1081  fCvsQuadrant4->Modified();
1082  fCvsQuadrant4->Update();
1083 }
1084 
1085 void AtEventDrawTask::UpdateCvsMC()
1086 {
1087  fMC_XY_exp->GetXaxis()->SetRangeUser(-300.0, 300);
1088  fMC_XY_exp->GetYaxis()->SetRangeUser(-300.0, 300);
1089  fMC_XY->GetXaxis()->SetRangeUser(-300.0, 300);
1090  fMC_XY->GetYaxis()->SetRangeUser(-300.0, 300);
1091  fMC_XY_int->GetXaxis()->SetRangeUser(-300.0, 300);
1092  fMC_XY_int->GetYaxis()->SetRangeUser(-300.0, 300);
1093  fMC_XY_back->GetXaxis()->SetRangeUser(-300.0, 300);
1094  fMC_XY_back->GetYaxis()->SetRangeUser(-300.0, 300);
1095  fCvsMC_XY->Modified();
1096  fCvsMC_XY->Update();
1097 
1098  fMC_ZX_int->GetXaxis()->SetRangeUser(0, 1000);
1099  fMC_ZX_int->GetYaxis()->SetRangeUser(-300.0, 300);
1100  fMC_ZY_int->GetXaxis()->SetRangeUser(0, 1000);
1101  fMC_ZY_int->GetYaxis()->SetRangeUser(-300.0, 300);
1102  fMC_ZX->GetXaxis()->SetRangeUser(0, 1000);
1103  fMC_ZX->GetYaxis()->SetRangeUser(-300.0, 300);
1104  fMC_ZY->GetXaxis()->SetRangeUser(0, 1000);
1105  fMC_ZY->GetYaxis()->SetRangeUser(-300.0, 300);
1106  fMC_ZX_back->GetXaxis()->SetRangeUser(0, 1000);
1107  fMC_ZX_back->GetYaxis()->SetRangeUser(-300.0, 300);
1108  fMC_ZY_back->GetXaxis()->SetRangeUser(0, 1000);
1109  fMC_ZY_back->GetYaxis()->SetRangeUser(-300.0, 300);
1110  fCvsMC_Z->Modified();
1111  fCvsMC_Z->Update();
1112 }
1113 
1114 void AtEventDrawTask::SetHitAttributes(Color_t color, Size_t size, Style_t style)
1115 {
1116  fHitColor = color;
1117  fHitSize = size;
1118  fHitStyle = style;
1119 }
1120 
1121 /*void
1122  AtEventDrawTask::SetHitClusterAttributes(Color_t color, Size_t size, Style_t style)
1123  {
1124  fHitClusterColor = color;
1125  fHitClusterSize = size;
1126  fHitClusterStyle = style;
1127  }*/
1128 
1129 /*void
1130  AtEventDrawTask::SetRiemannAttributes(Color_t color, Size_t size, Style_t style)
1131  {
1132  fRiemannColor = color;
1133  fRiemannSize = size;
1134  fRiemannStyle = style;
1135  }*/
1136 
1137 void AtEventDrawTask::SelectPad(const char *rawevt)
1138 {
1139 
1140  try {
1141  int event = gPad->GetEvent();
1142  if (event != 11)
1143  return; // may be comment this line
1144  TObject *select = gPad->GetSelected();
1145  if (!select)
1146  return;
1147  if (select->InheritsFrom(TH2Poly::Class())) {
1148  auto *h = dynamic_cast<TH2Poly *>(select);
1149  gPad->GetCanvas()->FeedbackMode(kTRUE);
1150  AtRawEvent *tRawEvent = nullptr;
1151  tRawEvent = dynamic_cast<AtRawEvent *>(gROOT->GetListOfSpecials()->FindObject(rawevt));
1152  if (tRawEvent == nullptr) {
1153  std::cout << " = AtEventDrawTask::SelectPad NULL pointer for the AtRawEvent! Please select an event first "
1154  << std::endl;
1155  return;
1156  }
1157 
1158  int pyold = gPad->GetUniqueID();
1159  int px = gPad->GetEventX();
1160  int py = gPad->GetEventY();
1161  float uxmin = gPad->GetUxmin();
1162  float uxmax = gPad->GetUxmax();
1163  int pxmin = gPad->XtoAbsPixel(uxmin);
1164  int pxmax = gPad->XtoAbsPixel(uxmax);
1165  if (pyold)
1166  TVirtualX::Instance()->DrawLine(pxmin, pyold, pxmax, pyold);
1167  TVirtualX::Instance()->DrawLine(pxmin, py, pxmax, py);
1168  gPad->SetUniqueID(py);
1169  Float_t upx = gPad->AbsPixeltoX(px);
1170  Float_t upy = gPad->AbsPixeltoY(py);
1171  Double_t x = gPad->PadtoX(upx);
1172  Double_t y = gPad->PadtoY(upy);
1173  Int_t bin = h->FindBin(x, y);
1174  const char *bin_name = h->GetBinName(bin);
1175  // std::cout<<" X : "<<x<<" Y: "<<y<<std::endl;
1176  // std::cout<<bin_name<<std::endl;
1177  std::cout << " ==========================" << std::endl;
1178  std::cout << " Bin number selected : " << bin << " Bin name :" << bin_name << std::endl;
1179 
1180  AtMap *tmap = nullptr;
1181  tmap = dynamic_cast<AtMap *>(gROOT->GetListOfSpecials()->FindObject("fMap"));
1182  Int_t tPadNum = tmap->BinToPad(bin);
1183 
1184  std::cout << " Bin : " << bin << " to Pad : " << tPadNum << std::endl;
1185  std::cout << " Electronic mapping: " << tmap->GetPadRef(tPadNum) << std::endl;
1186  AtPad *tPad = tRawEvent->GetPad(tPadNum);
1187 
1188  if (tPad == nullptr)
1189  return;
1190 
1191  std::cout << " Event ID (Select Pad) : " << tRawEvent->GetEventID() << std::endl;
1192  std::cout << " Raw Event Pad Num " << tPad->GetPadNum() << std::endl;
1193  std::cout << std::endl;
1194 
1195  TH1I *tPadWave = nullptr;
1196  tPadWave = dynamic_cast<TH1I *>(gROOT->GetListOfSpecials()->FindObject("fPadWave"));
1197  auto rawadc = tPad->GetRawADC();
1198  auto adc = tPad->GetADC();
1199  if (tPadWave == nullptr) {
1200  std::cout << " = AtEventDrawTask::SelectPad NULL pointer for the TH1I! Please enable SetPersistance for "
1201  "Unpacking task or select an event first "
1202  << std::endl;
1203  return;
1204  }
1205  tPadWave->Reset();
1206  // tPadWaveSub->Reset();
1207  for (Int_t i = 0; i < 512; i++) {
1208 
1209  // tPadWave->SetBinContent(i,rawadc[i]);
1210  tPadWave->SetBinContent(i, adc[i]);
1211  // tPadWaveSub->SetBinContent(i,adc[i]);
1212  }
1213 
1214  TCanvas *tCvsPadWave = nullptr;
1215  tCvsPadWave = dynamic_cast<TCanvas *>(gROOT->GetListOfSpecials()->FindObject("fCvsPadWave"));
1216  if (tCvsPadWave == nullptr) {
1217  std::cout << " = AtEventDrawTask::SelectPad NULL pointer for the TCanvas! Please select an event first "
1218  << std::endl;
1219  return;
1220  }
1221  tCvsPadWave->cd();
1222  tPadWave->Draw();
1223  // tPadWaveSub->Draw("SAME");
1224  tCvsPadWave->Update();
1225  }
1226 
1227  } catch (const std::exception &e) {
1228 
1229  std::cout << cRED << " Exception caught in Select Pad " << e.what() << cNORMAL << "\n";
1230  }
1231 }
1232 
1233 void AtEventDrawTask::DrawWave(Int_t PadNum)
1234 {
1235 
1236  // Bool_t IsValid=kFALSE;
1237  // AtPad *pad = fRawevent->GetPad(0);
1238  // AtPad *pad= fRawevent->GetPad(PadNum,IsValid);
1239  // std::cout<<" Raw Event Pad Num "<<pad->GetPadNum()<<" Is Valid? : "<<IsValidPad<<std::endl;
1240 }
1241 
1242 void AtEventDrawTask::ResetPadAll()
1243 {
1244  for (auto &i : fPadAll) {
1245  i->Reset(nullptr);
1246  }
1247 }
1248 
1249 void AtEventDrawTask::ResetPhiDistr()
1250 {
1251 
1252  for (auto &i : fPhiDistr) {
1253  i->Reset(nullptr);
1254  }
1255 }
1256 
1258 {
1259  f3DHitStyle = 0;
1260 }
1261 
1263 {
1264  f3DHitStyle = 1;
1265 }
1266 
1268 {
1269  fMultiHit = hitMax;
1270 }
1271 
1273 {
1274  fSaveTextData = kTRUE;
1275 }
1276 
1277 EColor AtEventDrawTask::GetTrackColor(int i)
1278 {
1279  std::vector<EColor> colors = {kAzure, kOrange, kViolet, kTeal, kMagenta, kBlue, kViolet, kYellow, kCyan, kAzure};
1280  if (i < 10) {
1281  return colors.at(i);
1282  } else
1283  return kAzure;
1284 }
1285 
1286 void AtEventDrawTask::SetRawEventBranch(TString branchName)
1287 {
1288  fRawEventBranchName = branchName;
1289 }
1290 
1291 void AtEventDrawTask::SetEventBranch(TString branchName)
1292 {
1293  fEventBranchName = branchName;
1294 }
AtEventDrawTask::fRawevent
AtRawEvent * fRawevent
Definition: AtEventDrawTask.h:58
AtEventDrawTask::fEventArray
TClonesArray * fEventArray
Definition: AtEventDrawTask.h:52
AtMap
Definition: AtMap.h:33
AtEvent::GetEventCharge
Double_t GetEventCharge() const
Definition: AtEvent.h:109
AtPad.h
AtAuxPad::GetAuxName
std::string GetAuxName() const
Definition: AtAuxPad.h:38
AtEventDrawTask::fMC_ZY_int
TGraph * fMC_ZY_int
Definition: AtEventDrawTask.h:127
AtEventDrawTask::fMC_ZX_int
TGraph * fMC_ZX_int
Definition: AtEventDrawTask.h:124
AtRawEvent.h
AtEventDrawTask.h
cYELLOW
constexpr auto cYELLOW
Definition: AtEventDrawTask.cxx:71
AtEventDrawTask::fQEventHist
TH1D * fQEventHist
Definition: AtEventDrawTask.h:85
AtEventDrawTask::fCvsHoughSpace
TCanvas * fCvsHoughSpace
Definition: AtEventDrawTask.h:87
AtBaseEvent.h
AtEventManager::GetEraseQEvent
Bool_t GetEraseQEvent()
Definition: AtEventManager.h:147
AtPatternEvent
Definition: AtPatternEvent.h:19
AtEventManager::GetCvsMC_Z
TCanvas * GetCvsMC_Z()
Definition: AtEventManager.h:129
AtEventDrawTask::fMC_ZY_back
TGraph * fMC_ZY_back
Definition: AtEventDrawTask.h:128
AtEventDrawTask::fTrackNum
Int_t fTrackNum
Definition: AtEventDrawTask.h:149
AtEvent.h
AtEventDrawTask::fPhiDistr
TH1D * fPhiDistr[5]
Definition: AtEventDrawTask.h:92
AtRawEvent::GetPad
AtPad * GetPad(Int_t padNum)
Definition: AtRawEvent.h:124
AtEventDrawTask::fHitSize
Size_t fHitSize
Definition: AtEventDrawTask.h:73
AtEventDrawTask::fCvsPadWave
TCanvas * fCvsPadWave
Definition: AtEventDrawTask.h:80
AtEventDrawTask::Exec
void Exec(Option_t *option)
Definition: AtEventDrawTask.cxx:233
cGREEN
constexpr auto cGREEN
Definition: AtEventDrawTask.cxx:73
AtEventDrawTask::fThreshold
Int_t fThreshold
Definition: AtEventDrawTask.h:62
AtEventDrawTask::fCorrectedEventArray
TClonesArray * fCorrectedEventArray
Definition: AtEventDrawTask.h:53
AtEventManager::GetDrawReconstruction
Bool_t GetDrawReconstruction()
Definition: AtEventManager.h:146
AtEventDrawTask::fCvsThetaxPhi
TCanvas * fCvsThetaxPhi
Definition: AtEventDrawTask.h:101
AtAuxPad.h
AtPadReference.h
S800Ana::SetTofObjCorr
void SetTofObjCorr(std::vector< Double_t > vec)
Definition: S800Ana.cxx:98
AtEventDrawTask::fRadVSTb
TH2F * fRadVSTb
Definition: AtEventDrawTask.h:98
ClassImp
ClassImp(AtEventDrawTask)
AtEventDrawTask::fCvsQuadrant4
TCanvas * fCvsQuadrant4
Definition: AtEventDrawTask.h:109
AtPad::GetValidPad
Bool_t GetValidPad() const
Definition: AtPad.h:97
AtBaseEvent::GetAuxPads
const AuxPadMap & GetAuxPads() const
Definition: AtBaseEvent.h:79
AtHitCluster.h
AtEventDrawTask::fPatternEventBranchName
TString fPatternEventBranchName
Definition: AtEventDrawTask.h:50
AtEventDrawTask::fAuxChannels
TH1F * fAuxChannels[9]
Definition: AtEventDrawTask.h:111
S800Ana::GetXfObj_ToF
Double_t GetXfObj_ToF()
Definition: S800Ana.cxx:127
AtEventDrawTask::fMesh
TH1F * fMesh
Definition: AtEventDrawTask.h:94
cNORMAL
constexpr auto cNORMAL
Definition: AtEventDrawTask.cxx:72
AtEventDrawTask::fRawEventArray
TClonesArray * fRawEventArray
Definition: AtEventDrawTask.h:54
S800Ana::SetMTDCXfRange
void SetMTDCXfRange(std::vector< Double_t > vec)
Definition: S800Ana.cxx:106
S800Calc
Definition: S800Calc.h:455
AtEventDrawTask::fIsRawData
Bool_t fIsRawData
Definition: AtEventDrawTask.h:140
AtEventDrawTask::fPID
TH2F * fPID
Definition: AtEventDrawTask.h:165
AtEventDrawTask::fCvsPID2
TCanvas * fCvsPID2
Definition: AtEventDrawTask.h:166
S800Calc::GetIsInCut
Bool_t GetIsInCut()
Definition: S800Calc.h:487
AtPattern.h
AtEventDrawTask::~AtEventDrawTask
~AtEventDrawTask()
AtEventDrawTask::fEventBranchName
TString fEventBranchName
Definition: AtEventDrawTask.h:48
AtEventDrawTask::fCvs3DHist
TCanvas * fCvs3DHist
Definition: AtEventDrawTask.h:95
AtEventManager::GetCvsRad
TCanvas * GetCvsRad()
Definition: AtEventManager.h:121
AtEventManager::GetCvsPhi
TCanvas * GetCvsPhi()
Definition: AtEventManager.h:118
AtEventDrawTask::fQEventHist_H
TH1D * fQEventHist_H
Definition: AtEventDrawTask.h:86
AtEventDrawTask::fCvsQuadrant3
TCanvas * fCvsQuadrant3
Definition: AtEventDrawTask.h:107
S800Ana::GetMTDCObjRange
std::vector< Double_t > GetMTDCObjRange()
Definition: S800Ana.cxx:115
AtEventDrawTask::fCvsMC_XY
TCanvas * fCvsMC_XY
Definition: AtEventDrawTask.h:117
AtEventDrawTask::fHitSetTFHC
std::vector< TEvePointSet * > fHitSetTFHC
Definition: AtEventDrawTask.h:157
S800Ana::GetObjCorr_ToF
Double_t GetObjCorr_ToF()
Definition: S800Ana.cxx:131
AtEventDrawTask::fCvsPhi
TCanvas * fCvsPhi
Definition: AtEventDrawTask.h:91
AtEventDrawTask::fPadWave
TH1I * fPadWave
Definition: AtEventDrawTask.h:81
AtTrack::GetHitArrayObject
std::vector< AtHit > GetHitArrayObject()
Definition: AtTrack.h:82
AtEventDrawTask::fCvsQuadrant1
TCanvas * fCvsQuadrant1
Definition: AtEventDrawTask.h:103
AtEventManager::GetCvsPadAll
TCanvas * GetCvsPadAll()
Definition: AtEventManager.h:115
AtEventDrawTask::DrawWave
void DrawWave(Int_t PadNum)
Definition: AtEventDrawTask.cxx:1233
AtEventDrawTask::fNumPads
static const Int_t fNumPads
Definition: AtEventDrawTask.h:45
AtEvent
Definition: AtEvent.h:22
AtEventDrawTask::fRawEventBranchName
TString fRawEventBranchName
Definition: AtEventDrawTask.h:47
AtEventDrawTask::fCvsRhoVariance
TCanvas * fCvsRhoVariance
Definition: AtEventDrawTask.h:89
AtEventDrawTask::SetHitAttributes
void SetHitAttributes(Color_t, Size_t, Style_t)
Definition: AtEventDrawTask.cxx:1114
AtEventDrawTask::f3DHitStyle
Int_t f3DHitStyle
Definition: AtEventDrawTask.h:135
AtEventDrawTask::SetSaveTextData
void SetSaveTextData()
Definition: AtEventDrawTask.cxx:1272
AtEventDrawTask::fMinTracksPerVertex
Int_t fMinTracksPerVertex
Definition: AtEventDrawTask.h:151
AtEventManager::GetCvsPadPlane
TCanvas * GetCvsPadPlane()
Definition: AtEventManager.h:113
AtRawEvent
Definition: AtRawEvent.h:34
AtEventDrawTask::fUnpackHough
Bool_t fUnpackHough
Definition: AtEventDrawTask.h:44
AtEventDrawTask::fMC_ZX_back
TGraph * fMC_ZX_back
Definition: AtEventDrawTask.h:125
AtEventDrawTask::fCvsTheta
TCanvas * fCvsTheta
Definition: AtEventDrawTask.h:99
AtEventManager::Get3DThreshold
Float_t Get3DThreshold()
Definition: AtEventManager.h:153
AtEventDrawTask::fhitBoxSet
TEveBoxSet * fhitBoxSet
Definition: AtEventDrawTask.h:68
AtEventDrawTask::fPadAll
TH1I * fPadAll[fNumPads]
Definition: AtEventDrawTask.h:83
AtEventDrawTask::fTheta
TH2F * fTheta
Definition: AtEventDrawTask.h:100
AtEventDrawTask::fPadPlanePal
TPaletteAxis * fPadPlanePal
Definition: AtEventDrawTask.h:70
AtEventDrawTask::fCorrectedHitSet
TEvePointSet * fCorrectedHitSet
Definition: AtEventDrawTask.h:66
S800Ana.h
AtEventDrawTask::fCvsPadPlane
TCanvas * fCvsPadPlane
Definition: AtEventDrawTask.h:78
AtHit::GetTimeStamp
Int_t GetTimeStamp() const
Definition: AtHit.h:86
AtEventDrawTask::fPID2
TH2F * fPID2
Definition: AtEventDrawTask.h:167
AtEventDrawTask::fRhoVariance
TH1D * fRhoVariance
Definition: AtEventDrawTask.h:90
AtHit::GetPosition
const XYZPoint & GetPosition() const
Definition: AtHit.h:79
AtEventDrawTask::Set3DHitStyleBar
void Set3DHitStyleBar()
Definition: AtEventDrawTask.cxx:1257
AtEventManager::GetCvsPadWave
TCanvas * GetCvsPadWave()
Definition: AtEventManager.h:114
AtHit.h
AtEventDrawTask::fCvsPID
TCanvas * fCvsPID
Definition: AtEventDrawTask.h:164
AtEventManager::Instance
static AtEventManager * Instance()
Definition: AtEventManager.cxx:63
AtEventDrawTask::fHitSet
TEvePointSet * fHitSet
Definition: AtEventDrawTask.h:65
AtEventDrawTask::fHitClusterSet
std::vector< TEveBoxSet * > fHitClusterSet
Definition: AtEventDrawTask.h:158
AtTrack::GetHitClusterArray
std::vector< AtHitCluster > * GetHitClusterArray()
Definition: AtTrack.h:90
AtPatternEvent::GetTrackCand
std::vector< AtTrack > & GetTrackCand()
Definition: AtPatternEvent.h:56
AtEventManager::GetCvsAux
TCanvas * GetCvsAux()
Definition: AtEventManager.h:130
AtEventDrawTask::fCvsQuadrant2
TCanvas * fCvsQuadrant2
Definition: AtEventDrawTask.h:105
AtEventDrawTask::Reset
void Reset()
Definition: AtEventDrawTask.cxx:640
y
const double * y
Definition: lmcurve.cxx:20
AtTrack
Definition: AtTrack.h:25
AtEventDrawTask::fRGBAPalette
TEveRGBAPalette * fRGBAPalette
Definition: AtEventDrawTask.h:162
AtEventDrawTask
Definition: AtEventDrawTask.h:41
AtEventDrawTask::fCorrectedEventBranchName
TString fCorrectedEventBranchName
Definition: AtEventDrawTask.h:49
AtEventDrawTask::fCvsPadAll
TCanvas * fCvsPadAll
Definition: AtEventDrawTask.h:82
S800Calc.h
S800Ana::GetTofObjCorr
std::vector< Double_t > GetTofObjCorr()
Definition: S800Ana.cxx:111
cBLUE
constexpr auto cBLUE
Definition: AtEventDrawTask.cxx:74
AtTrack.h
AtEventDrawTask::fMC_XY
TGraph * fMC_XY
Definition: AtEventDrawTask.h:118
AtPatternEvent.h
AtEventDrawTask::fS800Calc
S800Calc * fS800Calc
Definition: AtEventDrawTask.h:168
AtEventDrawTask::f3DHist
TH3F * f3DHist
Definition: AtEventDrawTask.h:96
S800Ana::GetMTDCXfRange
std::vector< Double_t > GetMTDCXfRange()
Definition: S800Ana.cxx:119
AtEventManager::GetCvsThetaxPhi
TCanvas * GetCvsThetaxPhi()
Definition: AtEventManager.h:123
AtEventDrawTask::fVertex
std::vector< TEvePointSet * > fVertex
Definition: AtEventDrawTask.h:160
AtEventDrawTask::fMC_XY_back
TGraph * fMC_XY_back
Definition: AtEventDrawTask.h:121
AtEventManager.h
AtEventDrawTask::SetEventBranch
void SetEventBranch(TString branchName)
Definition: AtEventDrawTask.cxx:1291
AtEventDrawTask::fSaveTextData
Bool_t fSaveTextData
Definition: AtEventDrawTask.h:137
AtEventDrawTask::fMC_ZY
TGraph * fMC_ZY
Definition: AtEventDrawTask.h:126
S800Ana::GetICSum_E
Double_t GetICSum_E()
Definition: S800Ana.cxx:135
AtEventDrawTask::SelectPad
static void SelectPad(const char *rawevt)
Definition: AtEventDrawTask.cxx:1137
AtPad::GetPadNum
Int_t GetPadNum() const
Definition: AtPad.h:96
AtEventDrawTask::fVertexStyle
Size_t fVertexStyle
Definition: AtEventDrawTask.h:76
AtPad::GetADC
const trace & GetADC() const
Definition: AtPad.cxx:97
AtEventManager::GetCvsMesh
TCanvas * GetCvsMesh()
Definition: AtEventManager.h:119
AtTrack::GetPattern
const AtPatterns::AtPattern * GetPattern() const
Definition: AtTrack.h:84
AtEventDrawTask::fMultiHit
Int_t fMultiHit
Definition: AtEventDrawTask.h:136
AtBaseEvent::GetEventID
ULong_t GetEventID() const
Definition: AtBaseEvent.h:67
AtEventDrawTask::fThetaxPhi_Ini
TH2F * fThetaxPhi_Ini
Definition: AtEventDrawTask.h:114
AtEventDrawTask::fVertexSize
Size_t fVertexSize
Definition: AtEventDrawTask.h:75
AtEventDrawTask::fDrawVertexFromLines
Bool_t fDrawVertexFromLines
Definition: AtEventDrawTask.h:141
AtEventManager::GetCvsTheta
TCanvas * GetCvsTheta()
Definition: AtEventManager.h:122
AtEventDrawTask::fCvsRad
TCanvas * fCvsRad
Definition: AtEventDrawTask.h:97
AtEventDrawTask::SetRawEventBranch
void SetRawEventBranch(TString branchName)
Definition: AtEventDrawTask.cxx:1286
AtEventDrawTask::fThetaxPhi_Ini_RANSAC
TH2F * fThetaxPhi_Ini_RANSAC
Definition: AtEventDrawTask.h:115
cRED
constexpr auto cRED
Event display task.
Definition: AtEventDrawTask.cxx:70
AtEventDrawTask::fCvsQEvent
TCanvas * fCvsQEvent
Definition: AtEventDrawTask.h:84
AtEventManager::GetDrawAllPad
Bool_t GetDrawAllPad()
Definition: AtEventManager.h:145
AtEventManager::GetCvsHoughSpace
TCanvas * GetCvsHoughSpace()
Definition: AtEventManager.h:117
AtEventDrawTask::AtEventDrawTask
AtEventDrawTask()
Definition: AtEventDrawTask.cxx:80
AtMap.h
AtRawEvent::GetNumPads
Int_t GetNumPads() const
Definition: AtRawEvent.h:122
AtEventManager::GetCvs3DHist
TCanvas * GetCvs3DHist()
Definition: AtEventManager.h:120
AtEventDrawTask::Set3DHitStyleBox
void Set3DHitStyleBox()
Definition: AtEventDrawTask.cxx:1262
AtEventDrawTask::fThetaxPhi
TH2F * fThetaxPhi
Definition: AtEventDrawTask.h:102
AtEventManager::GetCvsPID2
TCanvas * GetCvsPID2()
Definition: AtEventManager.h:132
AtMap::BinToPad
virtual Int_t BinToPad(Int_t binval)=0
AtEventDrawTask::fMC_XY_exp
TGraph * fMC_XY_exp
Definition: AtEventDrawTask.h:119
AtEventManager::GetToggleCorrData
Bool_t GetToggleCorrData()
Definition: AtEventManager.h:154
AtPad
Container class for AtPadBase objects.
Definition: AtPad.h:38
AtEventDrawTask::fIniHit
AtHit const * fIniHit
Definition: AtEventDrawTask.h:143
AtEventDrawTask::f3DThreshold
Float_t f3DThreshold
Definition: AtEventDrawTask.h:138
AtPatterns::AtPattern::GetEveElement
virtual TEveElement * GetEveElement() const =0
Get visual representation of pattern.
AtEventDrawTask::fCvsMC_Z
TCanvas * fCvsMC_Z
Definition: AtEventDrawTask.h:122
AtEventDrawTask::fHitLine
std::vector< TEveElement * > fHitLine
Definition: AtEventDrawTask.h:159
AtEventDrawTask::fPatternEventArray
TClonesArray * fPatternEventArray
Definition: AtEventDrawTask.h:55
AtHit::GetPadNum
Int_t GetPadNum() const
Definition: AtHit.h:83
AtEventDrawTask::fCvsAux
TCanvas * fCvsAux
Definition: AtEventDrawTask.h:112
AtFindVertex
Definition: AtFindVertex.h:27
AtPad::GetRawADC
const rawTrace & GetRawADC() const
Definition: AtPad.h:104
AtEventDrawTask::fIniHitRansac
AtHit const * fIniHitRansac
Definition: AtEventDrawTask.h:144
S800Ana::Calc
void Calc(S800Calc *s800calc)
Definition: S800Ana.cxx:151
AtEventDrawTask::Init
InitStatus Init()
Definition: AtEventDrawTask.cxx:133
S800Ana
Definition: S800Ana.h:19
S800Ana::SetMTDCObjRange
void SetMTDCObjRange(std::vector< Double_t > vec)
Definition: S800Ana.cxx:102
AtAuxPad
Definition: AtAuxPad.h:25
AtEventManager::GetCvsPID
TCanvas * GetCvsPID()
Definition: AtEventManager.h:131
AtEventDrawTask::fMC_XY_int
TGraph * fMC_XY_int
Definition: AtEventDrawTask.h:120
AtMap::GetPadRef
AtPadReference GetPadRef(int padNum) const
Definition: AtMap.cxx:281
AtEventDrawTask::fHoughSpace
TH2F * fHoughSpace
Definition: AtEventDrawTask.h:88
AtEventManager::GetS800Ana
S800Ana GetS800Ana()
Definition: AtEventManager.h:134
AtEventDrawTask::fHitStyle
Style_t fHitStyle
Definition: AtEventDrawTask.h:74
AtEventDrawTask::SetMultiHit
void SetMultiHit(Int_t hitMax)
Definition: AtEventDrawTask.cxx:1267
AtEventDrawTask::fHitColor
Color_t fHitColor
Definition: AtEventDrawTask.h:72
AtEventDrawTask::fEventManager
AtEventManager * fEventManager
Definition: AtEventDrawTask.h:57
AtEventDrawTask::fPadPlane
TH2Poly * fPadPlane
Definition: AtEventDrawTask.h:79
AtHit::GetCharge
Double_t GetCharge() const
Definition: AtHit.h:82
AtEventManager::GetCvsMC_XY
TCanvas * GetCvsMC_XY()
Definition: AtEventManager.h:128
AtHit
Point in space with charge.
Definition: AtHit.h:27
AtEventDrawTask::fMC_ZX
TGraph * fMC_ZX
Definition: AtEventDrawTask.h:123
AtEventDrawTask::fCvsMesh
TCanvas * fCvsMesh
Definition: AtEventDrawTask.h:93
AtEventDrawTask::fDetmap
std::shared_ptr< AtMap > fDetmap
Definition: AtEventDrawTask.h:60