]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveDet/AliEveTPCSector2D.cxx
Move contents of EVE/Alieve to EVE/EveDet as most code will remain there.
[u/mrichter/AliRoot.git] / EVE / EveDet / AliEveTPCSector2D.cxx
1 // $Id$
2 // Main authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007
3
4 /**************************************************************************
5  * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
6  * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for          *
7  * full copyright notice.                                                 *
8  **************************************************************************/
9
10 #include "AliEveTPCSector2D.h"
11 #include "AliEveTPCSector3D.h"
12
13 #include <EveDet/AliEveTPCData.h>
14 #include <EveDet/AliEveTPCSectorData.h>
15
16 #include <TEveManager.h>
17
18 #include <AliTPCParam.h>
19
20 #include <TBuffer3D.h>
21 #include <TBuffer3DTypes.h>
22 #include <TVirtualPad.h>
23 #include <TVirtualViewer3D.h>
24
25 #include <TH1S.h>
26 #include <TH2S.h>
27 #include <TVirtualPad.h>
28
29 using namespace std;
30
31 //______________________________________________________________________________
32 // AliEveTPCSector2D
33 //
34 // Displays TPC raw-data in 2D.
35 //
36 // fShowMax: true  - display maximum value for given time interval
37 //           false - display integral
38 // fAverage: only available when fShowMax = false; divide by time window width
39 //
40 // fUseTexture: use OpenGL textures to display data (fast rendering,
41 //   updates take the same time)
42 //
43
44 ClassImp(AliEveTPCSector2D)
45
46 /******************************************************************************/
47
48 AliEveTPCSector2D::AliEveTPCSector2D(const Text_t* n, const Text_t* t) :
49   AliEveTPCSectorViz(n,t),
50
51   fShowMax (kTRUE),
52   fAverage (kFALSE),
53
54   fUseTexture (kTRUE),
55   fPickEmpty  (kFALSE),
56   fPickMode   (1)
57 {}
58
59 AliEveTPCSector2D::~AliEveTPCSector2D()
60 {}
61
62 /******************************************************************************/
63
64 void AliEveTPCSector2D::MakeSector3D()
65 {
66   AliEveTPCSector3D* s = new AliEveTPCSector3D;
67   s->SetDataSource(fTPCData);
68   s->SetSectorID(fSectorID);
69   s->SetAutoTrans(fAutoTrans);
70   gEve->AddElement(s, this);
71   gEve->Redraw3D();
72 }
73
74 /******************************************************************************/
75
76 void AliEveTPCSector2D::ComputeBBox()
77 {
78   const AliEveTPCSectorData::SegmentInfo&  iSeg = AliEveTPCSectorData::GetInnSeg();
79   const AliEveTPCSectorData::SegmentInfo& o2Seg = AliEveTPCSectorData::GetOut2Seg();
80
81   BBoxInit();
82   Float_t w = o2Seg.GetNMaxPads()*o2Seg.GetPadWidth()/2;
83   fBBox[0] = -w;
84   fBBox[1] =  w;
85   fBBox[2] =  iSeg.GetRLow();
86   fBBox[3] =  o2Seg.GetRLow() + o2Seg.GetNRows()*o2Seg.GetPadHeight();
87   fBBox[4] = -0.5; // Fake z-width to 1 cm.
88   fBBox[5] =  0.5;
89 }
90
91 /******************************************************************************/
92
93 void AliEveTPCSector2D::PadSelected(Int_t row, Int_t pad)
94 {
95   // Called when ctrl-mouse-left-click registered over a pad.
96
97   // EVE -> Std convention
98   Int_t sseg = fSectorID, srow = row;
99   if (row >= AliEveTPCSectorData::GetInnSeg().GetNRows()) {
100     sseg += 36;
101     srow -= AliEveTPCSectorData::GetInnSeg().GetNRows();
102   }
103   switch (fPickMode)
104     {
105     case 0: {
106       printf("AliEveTPCSector2DGL::ProcessSelection segment=%d, row=%d, pad=%d\n",
107              sseg, srow, pad);
108       break;
109     }
110     case 1: {
111       AliEveTPCSectorData* sectorData = GetSectorData();
112       if (sectorData == 0) return;
113       Int_t mint = fMinTime;
114       Int_t maxt = fMaxTime;
115       TH1S* h = new TH1S(Form("Seg%d_Row%d_TEvePad%d", sseg, srow, pad),
116                          Form("Segment %d, Row %d, TEvePad %d", sseg, srow, pad),
117                          maxt - mint +1 , mint, maxt);
118       h->SetXTitle("Time");
119       h->SetYTitle("ADC");
120       AliEveTPCSectorData::PadIterator i = sectorData->MakePadIterator(row, pad);
121       while (i.Next())
122         h->Fill(i.Time(), i.Signal());
123       h->Draw();
124       gPad->Modified();
125       gPad->Update();
126       break;
127     }
128     case 2: {
129       AliEveTPCSectorData* sectorData = GetSectorData();
130       if (sectorData == 0) return;
131       Int_t mint = fMinTime;
132       Int_t maxt = fMaxTime;
133       Int_t npad = AliEveTPCSectorData::GetNPadsInRow(row);
134       TH2S* h = new TH2S(Form("Seg%d_Row%d", sseg, srow),
135                          Form("Segment %d, Row %d", sseg, srow),
136                          maxt - mint +1 , mint, maxt,
137                          npad, 0, npad - 1);
138       h->SetXTitle("Time");
139       h->SetYTitle("TEvePad");
140       h->SetZTitle("ADC");
141       AliEveTPCSectorData::RowIterator i = sectorData->MakeRowIterator(row);
142       while (i.NextPad())
143         while (i.Next())
144           h->Fill(i.Time(), i.TEvePad(), i.Signal());
145       h->Draw();
146       gPad->Modified();
147       gPad->Update();
148       break;
149     }
150     } // switch
151 }
152
153 /******************************************************************************/
154
155 void AliEveTPCSector2D::Paint(Option_t* )
156 {
157   if(fRnrSelf == kFALSE)
158     return;
159
160   TBuffer3D buffer(TBuffer3DTypes::kGeneric);
161
162   // Section kCore
163   buffer.fID           = this;
164   buffer.fColor        = 1;
165   buffer.fTransparency = 0;
166   fHMTrans.SetBuffer3D(buffer);
167   buffer.SetSectionsValid(TBuffer3D::kCore);
168
169   Int_t reqSections = gPad->GetViewer3D()->AddObject(buffer);
170   if (reqSections == TBuffer3D::kNone) {
171     return;
172   }
173
174   Error("AliEveTPCSector2D::Paint", "only direct OpenGL rendering supported.");
175 }