FairRoot/PandaRoot
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
PndMvdDigiPixelDraw Class Reference

#include <PndMvdDigiPixelDraw.h>

Inheritance diagram for PndMvdDigiPixelDraw:

Public Member Functions

 PndMvdDigiPixelDraw ()
 
 PndMvdDigiPixelDraw (const char *name, Int_t iVerbose=1)
 
virtual ~PndMvdDigiPixelDraw ()
 
void Exec (Option_t *option)
 
void Reset ()
 

Protected Member Functions

TVector3 GetVector (TObject *)
 
InitStatus Init ()
 
void SortDigis (TClonesArray *digis)
 
TEveBoxSet * CreateNewBoxSet (TString &name)
 
 ClassDef (PndMvdDigiPixelDraw, 1)
 

Protected Attributes

std::map< Int_t, TEveBoxSet * > fModules
 
PndGeoHandlingfGeoH
 
TClonesArray * fClusterCands
 
TClonesArray * fRecoHits
 
Double_t fPixelSize
 
Double_t fBoxSize
 
Bool_t fUseCluster
 

Detailed Description

Definition at line 23 of file PndMvdDigiPixelDraw.h.

Constructor & Destructor Documentation

PndMvdDigiPixelDraw::PndMvdDigiPixelDraw ( )

Definition at line 23 of file PndMvdDigiPixelDraw.cxx.

23  : fPixelSize(0.01), fBoxSize(1)
24 {}
PndMvdDigiPixelDraw::PndMvdDigiPixelDraw ( const char *  name,
Int_t  iVerbose = 1 
)

Definition at line 26 of file PndMvdDigiPixelDraw.cxx.

26  : FairBoxSetDraw(name, iVerbose), fPixelSize(0.01), fBoxSize(1)
27 {}
TString name
Int_t iVerbose
PndMvdDigiPixelDraw::~PndMvdDigiPixelDraw ( )
virtual

Definition at line 29 of file PndMvdDigiPixelDraw.cxx.

References fGeoH.

30 {
31  delete(fGeoH);
32 }
PndGeoHandling * fGeoH

Member Function Documentation

PndMvdDigiPixelDraw::ClassDef ( PndMvdDigiPixelDraw  ,
 
)
protected
TEveBoxSet * PndMvdDigiPixelDraw::CreateNewBoxSet ( TString name)
protected

Definition at line 128 of file PndMvdDigiPixelDraw.cxx.

References fBoxSize.

Referenced by Exec().

129 {
130  TEveBoxSet* bs = new TEveBoxSet(name);
131  bs->Reset(TEveBoxSet::kBT_AABoxFixedDim, kFALSE, 32);
132  bs->SetDefWidth(fBoxSize);
133  bs->SetDefHeight(fBoxSize);
134  bs->SetDefDepth(fBoxSize);
135 
136  return bs;
137 }
TString name
void PndMvdDigiPixelDraw::Exec ( Option_t *  option)

Definition at line 45 of file PndMvdDigiPixelDraw.cxx.

References PndSdsCalcFePixel::CalcSensorColRow(), col, CreateNewBoxSet(), fBoxSize, fClusterCands, fe, fGeoH, fModules, fPixelSize, fRecoHits, fUseCluster, PndSdsDigi::GetCharge(), PndSdsCluster::GetClusterSize(), PndSdsCluster::GetDigiIndex(), PndSdsDigi::GetFE(), PndGeoHandling::GetMatrixShortId(), PndGeoHandling::GetPath(), PndSdsDigiPixel::GetPixelColumn(), PndSdsDigiPixel::GetPixelRow(), PndGeoHandling::GetSensorDimensionsShortId(), PndSdsDigi::GetSensorID(), PndSdsHit::GetSensorID(), i, PndGeoHandling::MasterToLocalShortId(), p, Reset(), row, PndSdsPixelDigiPar::SetFECols(), PndSdsPixelDigiPar::SetFERows(), PndSdsPixelDigiPar::SetMaxFEperCol(), PndSdsPixelDigiPar::SetMaxFEperRow(), t, TString, X, Y, and Z.

46 {
47  PndSdsPixelDigiPar digipar;
48  digipar.SetFECols(104);
49  digipar.SetFERows(104);
50  digipar.SetMaxFEperCol(10);
51  digipar.SetMaxFEperRow(10);
52  PndSdsCalcFePixel calc(digipar);
53  Int_t col, row, fe;
54  Reset();
55  if (fUseCluster == kTRUE){
56  for (int hitNr = 0; hitNr < fRecoHits->GetEntriesFast(); hitNr++){
57  PndSdsHit* myRecoHit = (PndSdsHit*)fRecoHits->At(hitNr);
58  TVector3 recoVector;
59  myRecoHit->Position(recoVector);
60  TVector3 recoLocal = fGeoH->MasterToLocalShortId(recoVector, myRecoHit->GetSensorID());
61  TVector3 SensDim = fGeoH->GetSensorDimensionsShortId(myRecoHit->GetSensorID());
62  TString detName = fGeoH->GetPath(myRecoHit->GetSensorID());
63  Int_t sensorId = myRecoHit->GetSensorID();
64  if (fModules[sensorId] == 0)
65  fModules[sensorId] = CreateNewBoxSet(detName);
66 
67  PndSdsCluster* myCluster = (PndSdsCluster*)fClusterCands->At(myRecoHit->GetRefIndex());
68  for (int clusterNr = 0; clusterNr < myCluster->GetClusterSize(); clusterNr++){
69  PndSdsDigiPixel* p = (PndSdsDigiPixel*)fList->At(myCluster->GetDigiIndex(clusterNr));
70  col = p->GetPixelColumn();
71  row = p->GetPixelRow();
72  fe = p->GetFE();
73  std::cout << "Fe, Col, Row: " << fe << " " << col << "/" << row << std::endl;
74  calc.CalcSensorColRow(col, row, p->GetFE());
75  std::cout << "Col, Row " << p->GetSensorID() << " : " << col << "/" << row << std::endl;
76  TVector3 digiLocal(((col+0.5) - (SensDim.X()*100))*fPixelSize, ((row+0.5)-(SensDim.Y()*100))*fPixelSize, 0);
77  std::cout << "DigiLocal: " << digiLocal.X() << " " << digiLocal.Y() << " " << digiLocal.Z() << std::endl;
78  std::cout << "RecoLocal: " << recoLocal.X() << " " << recoLocal.Y() << " " << recoLocal.Z() << std::endl;
79  TVector3 diff = digiLocal - recoLocal;
80  std::cout << "Diff: " << diff.X() << " " << diff.Y() << " " << diff.Z() << std::endl;
81  TVector3 final = diff * (fBoxSize/fPixelSize);
82  final += recoLocal;
83  std::cout << "Final: " << final.X() << " " << final.Y() << " " << final.Z() << std::endl;
84 
85  fModules[p->GetSensorID()]->AddBox(final.X(), final.Y(), final.Z());
86  fModules[p->GetSensorID()]->DigitValue(p->GetCharge());
87  }
88  }
89  }
90  else {
91  for (Int_t i=0; i<fList->GetEntriesFast(); ++i) {
92  PndSdsDigiPixel* p = (PndSdsDigiPixel*)fList->At(i);
93  TString detName = fGeoH->GetPath(p->GetSensorID());
94  Int_t sensorId = p->GetSensorID();
95 
96  if (fModules[sensorId] == 0)
97  fModules[sensorId] = CreateNewBoxSet(detName);
98 
99  col = p->GetPixelColumn();
100  row = p->GetPixelRow();
101  fe = p->GetFE();
102  std::cout << "Fe, Col, Row: " << fe << " " << col << "/" << row << std::endl;
103  calc.CalcSensorColRow(col, row, p->GetFE());
104  std::cout << "Col, Row " << p->GetSensorID() << " : " << col << "/" << row << std::endl;
105  TVector3 SensDim = fGeoH->GetSensorDimensionsShortId(p->GetSensorID());
106  std::cout << "SensDim: " << SensDim.X() << " " << SensDim.Y() << " " << SensDim.Z() << std::endl;
107  fModules[p->GetSensorID()]->AddBox(((col+0.5) - (SensDim.X()*100))*0.01, ((row+0.5)-(SensDim.Y()*100))*0.01, 0.01);
108  fModules[p->GetSensorID()]->DigitValue(p->GetCharge());
109  }
110  }
111 
112  for (boxSetMapIter it = fModules.begin(); it != fModules.end(); it++){
113  TGeoHMatrix testMatrix = *(fGeoH->GetMatrixShortId(it->first));
114  //TGeoHMatrix invMatrix = testMatrix.Inverse();
115  //Double_t scale[] = {0.01,0.01,0.01}; //Dim is 0.1 mm //[R.K. 01/2017] unused variable
116  //testMatrix.SetScale(scale);
117  TEveTrans& t = it->second->RefMainTrans();
118  t.SetFrom(testMatrix);
119  //TEveElement* el = (TEveElement*)it->second; //[R.K. 01/2017] unused variable
120  TEveElement* man = (TEveElement*)fEventManager;
121 // gEve->AddElement(it->second, fEventManager);
122  gEve->AddElement(it->second, man);
123  }
124 
125  gEve->Redraw3D(kFALSE);
126 }
int row
Definition: anaLmdDigi.C:67
Int_t GetPixelRow() const
TClonesArray * fRecoHits
Int_t GetClusterSize() const
Definition: PndSdsCluster.h:39
Int_t GetSensorID() const
Definition: PndSdsDigi.h:59
std::map< Int_t, TEveBoxSet * >::iterator boxSetMapIter
Int_t i
Definition: run_full.C:25
Class to store the Digis which belong to one cluster This class holds the information which Digi belo...
Definition: PndSdsCluster.h:19
Int_t GetPixelColumn() const
TEveBoxSet * CreateNewBoxSet(TString &name)
int col
Definition: anaLmdDigi.C:67
Double_t GetCharge() const
Definition: PndSdsDigi.h:60
std::map< Int_t, TEveBoxSet * > fModules
Int_t GetFE() const
Definition: PndSdsDigi.h:57
double Y
Definition: anaLmdDigi.C:68
TString GetPath(Int_t shortID)
for a given shortID the path is returned
Double_t p
Definition: anasim.C:58
TVector3 GetSensorDimensionsShortId(Int_t shortId)
TGeoHMatrix * GetMatrixShortId(Int_t shortId)
Int_t GetDigiIndex(Int_t i) const
Definition: PndSdsCluster.h:40
TVector3 MasterToLocalShortId(const TVector3 &master, const Int_t &shortId)
TClonesArray * fClusterCands
PndGeoHandling * fGeoH
double X
Definition: anaLmdDigi.C:68
int fe
Definition: anaLmdDigi.C:67
double Z
Definition: anaLmdDigi.C:68
TTree * t
Definition: bump_analys.C:13
Data class to store the digi output of a pixel module.
void SetMaxFEperRow(Int_t x)
void SetMaxFEperCol(Int_t x)
Int_t GetSensorID() const
Definition: PndSdsHit.h:90
Class to calculate the position of digis on a front-end from the digis on a sensor.
void SetFECols(Int_t x)
void SetFERows(Int_t x)
Digitization Parameter Class for SDS-Pixel part.
TVector3 PndMvdDigiPixelDraw::GetVector ( TObject *  )
inlineprotected

Definition at line 35 of file PndMvdDigiPixelDraw.h.

35 {return TVector3();}; // obj //[R.K.03/2017] unused variable(s)
InitStatus PndMvdDigiPixelDraw::Init ( )
protected

Definition at line 34 of file PndMvdDigiPixelDraw.cxx.

References fClusterCands, fGeoH, fRecoHits, fUseCluster, Init(), and PndGeoHandling::Instance().

35 {
38  fClusterCands = (TClonesArray *)fManager->GetObject("MVDPixelClusterCand");
39  fRecoHits = (TClonesArray *)fManager->GetObject("MVDHitsPixel");
40 
41  if(fClusterCands != 0 && fRecoHits != 0) fUseCluster = kTRUE;
42  return kSUCCESS;
43 }
TClonesArray * fRecoHits
static PndGeoHandling * Instance()
TClonesArray * fClusterCands
PndGeoHandling * fGeoH
fRun Init()
Definition: NHitsPerEvent.C:20
void PndMvdDigiPixelDraw::Reset ( )

Definition at line 139 of file PndMvdDigiPixelDraw.cxx.

References fModules.

Referenced by Exec().

140 {
141  for (boxSetMapIter it = fModules.begin(); it != fModules.end(); it++){
142  //TEveElement* el = (TEveElement*)it->second; //[R.K. 01/2017] unused variable
143  TEveElement* man = (TEveElement*)fEventManager;
144 // gEve->AddElement(it->second, fEventManager);
145  gEve->AddElement(it->second, man);
146  }
147  fModules.clear();
148 }
std::map< Int_t, TEveBoxSet * >::iterator boxSetMapIter
std::map< Int_t, TEveBoxSet * > fModules
void PndMvdDigiPixelDraw::SortDigis ( TClonesArray *  digis)
protected

Member Data Documentation

Double_t PndMvdDigiPixelDraw::fBoxSize
protected

Definition at line 48 of file PndMvdDigiPixelDraw.h.

Referenced by CreateNewBoxSet(), and Exec().

TClonesArray* PndMvdDigiPixelDraw::fClusterCands
protected

Definition at line 44 of file PndMvdDigiPixelDraw.h.

Referenced by Exec(), and Init().

PndGeoHandling* PndMvdDigiPixelDraw::fGeoH
protected

Definition at line 43 of file PndMvdDigiPixelDraw.h.

Referenced by Exec(), Init(), and ~PndMvdDigiPixelDraw().

std::map<Int_t, TEveBoxSet* > PndMvdDigiPixelDraw::fModules
protected

Definition at line 42 of file PndMvdDigiPixelDraw.h.

Referenced by Exec(), and Reset().

Double_t PndMvdDigiPixelDraw::fPixelSize
protected

Definition at line 47 of file PndMvdDigiPixelDraw.h.

Referenced by Exec().

TClonesArray* PndMvdDigiPixelDraw::fRecoHits
protected

Definition at line 45 of file PndMvdDigiPixelDraw.h.

Referenced by Exec(), and Init().

Bool_t PndMvdDigiPixelDraw::fUseCluster
protected

Definition at line 50 of file PndMvdDigiPixelDraw.h.

Referenced by Exec(), and Init().


The documentation for this class was generated from the following files: