Belle II Software development
EclPainterPolar Class Reference

Painter for EclData, polar energy/event_count distribution. More...

#include <EclPainterPolar.h>

Inheritance diagram for EclPainterPolar:
EclPainter

Public Types

enum  Type {
  PHI ,
  THETA
}
 Type for polar histogram. More...
 

Public Member Functions

 EclPainterPolar (EclData *data, Type type)
 Constructor for EclPainter subclass.
 
 EclPainterPolar (const EclPainterPolar &other)
 Copy constructor.
 
 ~EclPainterPolar ()
 Destructor for EclPainter subclass.
 
EclPainterPolaroperator= (const EclPainterPolar &other)
 Assignment operator.
 
virtual void getInformation (int px, int py, MultilineWidget *panel) override
 Sets the information to be displayed in the provided MultilineWidget.
 
Type getType ()
 Return subtype of ECLPainterPolar.
 
virtual void Draw () override
 Redraw the canvas.
 
void setData (EclData *data)
 Set EclData to display in painter.
 
EclDatagetData ()
 Return currently displayed EclData.
 
const EclDatagetData () const
 Return currently displayed EclData.
 
void setMapper (ECL::ECLChannelMapper *mapper)
 Set ECLChannelMapper for CellID <-> (crate, shaper, chid) conversion.
 
ECL::ECLChannelMappergetMapper ()
 Return currently set ECLChannelMapper.
 
void setDisplayedSubsystem (EclData::EclSubsystem sys)
 Change between the displayed ECL subsystem (barrel, forward and backward endcaps).
 
EclData::EclSubsystem getDisplayedSubsystem ()
 Get currently displayed ECL subsystem.
 
TString getSubsystemTitle (EclData::EclSubsystem subsys)
 Return title of ECL subsystem to use in painter.
 
virtual EclPainterhandleClick (int px, int py)
 Some EclPainters can shift to another view upon click.
 
virtual void setXRange (int x1, int x2)
 Set XRange for histogram in EclPainter.
 

Protected Member Functions

void getNewRootObjectName (char *buf, int size)
 Make unique name for next root object.
 

Private Member Functions

void cloneFrom (const EclPainterPolar &other)
 Clone attributes from other EclPainterPolar.
 
int channelToSegId (int channel)
 Convert ECL channel id to id of the phi (theta) segment.
 
void setTitles ()
 Update titles of the histogram.
 

Private Attributes

Type m_type
 Type for polar histogram.
 
TH2F * m_hist
 Histogram that generates Z-axis.
 
TCrown ** m_segs
 Phi (or theta) segments of the ECL.
 
TText ** m_labels
 Labels for phi segments.
 
EclDatam_ecl_data
 Data to draw.
 
ECL::ECLChannelMapperm_mapper
 mapper for CellID <-> (crate, shaper, chid) conversion.
 
EclData::EclSubsystem displayed_subsys
 Identifier of displayed ECL subsystem.
 

Static Private Attributes

static int m_obj_counter = 0
 Counter to make unique names for new root objects.
 

Detailed Description

Painter for EclData, polar energy/event_count distribution.

Definition at line 26 of file EclPainterPolar.h.

Member Enumeration Documentation

◆ Type

enum Type

Type for polar histogram.

Definition at line 31 of file EclPainterPolar.h.

31{PHI, THETA};

Constructor & Destructor Documentation

◆ EclPainterPolar() [1/2]

Constructor for EclPainter subclass.

Definition at line 25 of file EclPainterPolar.cc.

25 :
26 EclPainter(data)
27{
28 m_type = type;
29
30 char obj_name[255];
31 getNewRootObjectName(obj_name, 255);
32 m_hist = new TH2F(obj_name, "title",
33 60, 0.0, 1.0, 60, 0.0, 1.0);
34
35 Double_t deg2rad = TMath::Pi() / 180;
36
37 m_segs = new TCrown*[36];
38 m_labels = new TText*[36];
39 char label_txt[32];
40 for (int i = 0; i < 36; i++) {
41 m_segs[i] = new TCrown(0.5, 0.5, 0.3, 0.4, (i - 9) * 10, (i - 8) * 10);
42
43 float x = 0.475 + 0.45 * TMath::Cos(deg2rad * (i - 9) * 10);
44 float y = 0.48 + 0.44 * TMath::Sin(deg2rad * (i - 9) * 10);
45 snprintf(label_txt, 32, "%d", i * 10);
46 m_labels[i] = new TText(x, y, label_txt);
47 m_labels[i]->SetTextSize(0.03);
48 }
49}
Type m_type
Type for polar histogram.
TText ** m_labels
Labels for phi segments.
TCrown ** m_segs
Phi (or theta) segments of the ECL.
TH2F * m_hist
Histogram that generates Z-axis.
Painter for EclData, parent class, created with EclPainterFactory.
Definition: EclPainter.h:29
void getNewRootObjectName(char *buf, int size)
Make unique name for next root object.
Definition: EclPainter.cc:88

◆ EclPainterPolar() [2/2]

EclPainterPolar ( const EclPainterPolar other)
inline

Copy constructor.

Definition at line 38 of file EclPainterPolar.h.

38: EclPainter(other) { cloneFrom(other); }
void cloneFrom(const EclPainterPolar &other)
Clone attributes from other EclPainterPolar.
EclPainter(EclData *data)
Default constructor.
Definition: EclPainter.cc:24

◆ ~EclPainterPolar()

Destructor for EclPainter subclass.

Definition at line 51 of file EclPainterPolar.cc.

52{
53 delete m_hist;
54}

Member Function Documentation

◆ channelToSegId()

int channelToSegId ( int  channel)
private

Convert ECL channel id to id of the phi (theta) segment.

Definition at line 66 of file EclPainterPolar.cc.

67{
68 if (m_type == PHI)
69 return getData()->getPhiId(ch) / 4;
70 else if (m_type == THETA) {
71 int theta_id = getData()->getThetaId(ch);
72
73 if (theta_id < 23)
74 return 3 + theta_id * 12 / 23;
75 else
76 return 21 + (theta_id - 23) * 12 / 23;
77 }
78
79 return 0;
80}
int getThetaId(int ch)
ECL CellId -> theta_id.
Definition: EclData.cc:288
int getPhiId(int ch)
ECL CellId -> phi_id.
Definition: EclData.cc:278
EclData * getData()
Return currently displayed EclData.
Definition: EclPainter.h:46

◆ cloneFrom()

void cloneFrom ( const EclPainterPolar other)
private

Clone attributes from other EclPainterPolar.

Definition at line 56 of file EclPainterPolar.cc.

57{
58 m_type = other.m_type;
59 m_hist = new TH2F(*other.m_hist);
60 m_segs = new TCrown*[36];
61 for (int i = 0; i < 36; i++) { m_segs[i] = other.m_segs[i]; }
62 m_labels = new TText*[36];
63 for (int i = 0; i < 36; i++) { m_labels[i] = other.m_labels[i]; }
64}

◆ Draw()

void Draw ( )
overridevirtual

Redraw the canvas.

Implements EclPainter.

Definition at line 108 of file EclPainterPolar.cc.

109{
110 setTitles();
111
112 EclData* data = getData();
113 const int* ev_counts = data->getEventCounts();
114 const float* energy_sums = data->getEnergySums();
115
116 float seg_val[36];
117 TCrown** segs = m_segs;
118 for (int i = 0; i < 36; i++)
119 seg_val[i] = 0;
120
121 for (int i = 1; i <= getData()->getCrystalCount(); i++) {
122 if (!data->isCrystalInSubsystem(i, getDisplayedSubsystem())) continue;
123 int id = channelToSegId(i);
124 if (GetMode())
125 seg_val[id] += energy_sums[i];
126 else
127 seg_val[id] += ev_counts[i];
128 }
129
130 float max = 0;
131 float min = seg_val[0];
132 for (int i = 0; i < 36; i++) {
133 if (max < seg_val[i])
134 max = seg_val[i];
135 if (min > seg_val[i])
136 min = seg_val[i];
137 }
138
139 Double_t r[5] = { 0.00, 0.00, 0.87, 1.00, 0.51 };
140 Double_t g[5] = { 0.00, 0.81, 1.00, 0.20, 0.00 };
141 Double_t b[5] = { 0.51, 1.00, 0.12, 0.00, 0.00 };
142 Double_t stop[5] = { 0.00, 0.34, 0.61, 0.84, 1.00 };
143
144 int palette = TColor::CreateGradientColorTable(5, stop, r, g, b, 37);
145 m_hist->Reset();
146 m_hist->Fill(0.5, 0.05, 0.1);
147 m_hist->SetMaximum(max);
148 m_hist->SetMinimum(min);
149 m_hist->Draw("COLZ");
150
151 for (int i = 0; i < 36; i++) {
152 float val = 36.0 * TMath::Log(1 + seg_val[i]) / TMath::Log(1 + max);
153 segs[i]->SetFillColor(palette + val);
154 segs[i]->Draw(/*"SAME"*/);
155 m_labels[i]->Draw();
156 }
157}
This class contains data for ECLSimHit's and provides several relevant conversion functions for bette...
Definition: EclData.h:31
static int getCrystalCount()
Get number of crystals in ECL.
Definition: EclData.cc:148
void setTitles()
Update titles of the histogram.
int channelToSegId(int channel)
Convert ECL channel id to id of the phi (theta) segment.
EclData::EclSubsystem getDisplayedSubsystem()
Get currently displayed ECL subsystem.
Definition: EclPainter.cc:51
int GetMode()
Returns current displayed mode (0 shows event count, 1 shows total energy)
Definition: geometry.cc:19

◆ getData() [1/2]

EclData * getData ( )
inlineinherited

Return currently displayed EclData.

Definition at line 46 of file EclPainter.h.

46{ return m_ecl_data; }
EclData * m_ecl_data
Data to draw.
Definition: EclPainter.h:111

◆ getData() [2/2]

const EclData * getData ( ) const
inlineinherited

Return currently displayed EclData.

Definition at line 48 of file EclPainter.h.

48{ return m_ecl_data; }

◆ getDisplayedSubsystem()

EclData::EclSubsystem getDisplayedSubsystem ( )
inherited

Get currently displayed ECL subsystem.

Definition at line 51 of file EclPainter.cc.

52{
53 return displayed_subsys;
54}
EclData::EclSubsystem displayed_subsys
Identifier of displayed ECL subsystem.
Definition: EclPainter.h:116

◆ getInformation()

void getInformation ( int  px,
int  py,
MultilineWidget panel 
)
overridevirtual

Sets the information to be displayed in the provided MultilineWidget.

Parameters
pxX coordinate of mouse cursor.
pyY coordinate of mouse cursor.
panelMultilineWidget to display the information

Reimplemented from EclPainter.

Definition at line 98 of file EclPainterPolar.cc.

99{
100 EclPainter::getInformation(px, py, panel);
101}
virtual void getInformation(int px, int py, MultilineWidget *panel)
Sets the information to be displayed in the provided MultilineWidget.
Definition: EclPainter.cc:72

◆ getMapper()

ECL::ECLChannelMapper * getMapper ( )
inherited

Return currently set ECLChannelMapper.

Definition at line 41 of file EclPainter.cc.

42{
43 return m_mapper;
44}
ECL::ECLChannelMapper * m_mapper
mapper for CellID <-> (crate, shaper, chid) conversion.
Definition: EclPainter.h:113

◆ getNewRootObjectName()

void getNewRootObjectName ( char *  buf,
int  size 
)
protectedinherited

Make unique name for next root object.

Definition at line 88 of file EclPainter.cc.

89{
90 snprintf(buf, n, "ECL DATA_%d", m_obj_counter++);
91}
static int m_obj_counter
Counter to make unique names for new root objects.
Definition: EclPainter.h:109

◆ getSubsystemTitle()

TString getSubsystemTitle ( EclData::EclSubsystem  subsys)
inherited

Return title of ECL subsystem to use in painter.

Definition at line 56 of file EclPainter.cc.

57{
58 switch (subsys) {
59 case EclData::BARR:
60 return TString("Barrel");
61 case EclData::FORW:
62 return TString("Forward endcap");
63 case EclData::BACKW:
64 return TString("Backward endcap");
65 case EclData::ALL:
66 return TString("Full ECL");
67 default:
68 return TString();
69 }
70}

◆ getType()

EclPainterPolar::Type getType ( )

Return subtype of ECLPainterPolar.

Definition at line 103 of file EclPainterPolar.cc.

104{
105 return m_type;
106}

◆ handleClick()

EclPainter * handleClick ( int  px,
int  py 
)
virtualinherited

Some EclPainters can shift to another view upon click.

(For example, clicking on crate reveals histogram of shapers in that crate)

Returns
EclPainter with new perspective/range.

Reimplemented in EclPainter1D.

Definition at line 79 of file EclPainter.cc.

80{
81 return nullptr;
82}

◆ operator=()

EclPainterPolar & operator= ( const EclPainterPolar other)
inline

Assignment operator.

Definition at line 45 of file EclPainterPolar.h.

45{ cloneFrom(other); return *this; }

◆ setData()

void setData ( EclData data)
inlineinherited

Set EclData to display in painter.


Definition at line 44 of file EclPainter.h.

44{ m_ecl_data = data; }

◆ setDisplayedSubsystem()

void setDisplayedSubsystem ( EclData::EclSubsystem  sys)
inherited

Change between the displayed ECL subsystem (barrel, forward and backward endcaps).

Definition at line 46 of file EclPainter.cc.

47{
48 displayed_subsys = sys;
49}

◆ setMapper()

void setMapper ( ECL::ECLChannelMapper mapper)
inherited

Set ECLChannelMapper for CellID <-> (crate, shaper, chid) conversion.

Definition at line 36 of file EclPainter.cc.

37{
38 m_mapper = mapper;
39}

◆ setTitles()

void setTitles ( )
private

Update titles of the histogram.

Definition at line 82 of file EclPainterPolar.cc.

83{
84 const char* name[2][2] = {
85 {"Events per phi", "Events per theta"},
86 {"Energy per phi", "Energy per theta"}
87 };
88 const char* zname[2] = {"Event count", "Energy (MeV)"};
89
90 TString title = TString(name[GetMode()][m_type]) + " (" +
92
93 m_hist->SetTitle(title);
94 m_hist->SetZTitle(zname[GetMode()]);
95 m_hist->GetZaxis()->CenterTitle();
96}
TString getSubsystemTitle(EclData::EclSubsystem subsys)
Return title of ECL subsystem to use in painter.
Definition: EclPainter.cc:56

◆ setXRange()

void setXRange ( int  x1,
int  x2 
)
virtualinherited

Set XRange for histogram in EclPainter.

Reimplemented in EclPainter1D.

Definition at line 84 of file EclPainter.cc.

85{
86}

Member Data Documentation

◆ displayed_subsys

EclData::EclSubsystem displayed_subsys
privateinherited

Identifier of displayed ECL subsystem.

Definition at line 116 of file EclPainter.h.

◆ m_ecl_data

EclData* m_ecl_data
privateinherited

Data to draw.

Definition at line 111 of file EclPainter.h.

◆ m_hist

TH2F* m_hist
private

Histogram that generates Z-axis.

Definition at line 51 of file EclPainterPolar.h.

◆ m_labels

TText** m_labels
private

Labels for phi segments.

Definition at line 55 of file EclPainterPolar.h.

◆ m_mapper

ECL::ECLChannelMapper* m_mapper
privateinherited

mapper for CellID <-> (crate, shaper, chid) conversion.

Definition at line 113 of file EclPainter.h.

◆ m_obj_counter

int m_obj_counter = 0
staticprivateinherited

Counter to make unique names for new root objects.

Definition at line 109 of file EclPainter.h.

◆ m_segs

TCrown** m_segs
private

Phi (or theta) segments of the ECL.

Definition at line 53 of file EclPainterPolar.h.

◆ m_type

Type m_type
private

Type for polar histogram.

Definition at line 49 of file EclPainterPolar.h.


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