]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveDet/AliEveTRDData.cxx
Move contents of EVE/Alieve to EVE/EveDet as most code will remain there.
[u/mrichter/AliRoot.git] / EVE / EveDet / AliEveTRDData.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 #include "AliEveTRDData.h"
10 #include "AliEveTRDModuleImp.h"
11
12 #include "AliLog.h"
13 #include "AliTRDhit.h"
14 #include "AliTRDcluster.h"
15 #include "AliTRDcalibDB.h"
16 #include "AliTRDpadPlane.h"
17 #include "AliTRDgeometry.h"
18 #include "AliTRDdigitsManager.h"
19
20 using namespace std;
21
22 ClassImp(AliEveTRDHits)
23 ClassImp(AliEveTRDDigits)
24 ClassImp(AliEveTRDClusters)
25
26 ///////////////////////////////////////////////////////////
27 /////////////   AliEveTRDDigits             /////////////////////
28 ///////////////////////////////////////////////////////////
29
30 //________________________________________________________
31 AliEveTRDDigits::AliEveTRDDigits(AliEveTRDChamber *p): TEveQuadSet("digits", ""), fParent(p)
32 {}
33
34 //________________________________________________________
35 void    AliEveTRDDigits::SetData(AliTRDdigitsManager *digits)
36 {
37
38         fData.Allocate(fParent->rowMax, fParent->colMax, fParent->timeMax);
39 //      digits->Expand();
40         for (Int_t  row = 0;  row <  fParent->rowMax;  row++)
41                 for (Int_t  col = 0;  col <  fParent->colMax;  col++)
42                         for (Int_t time = 0; time < fParent->timeMax; time++) {
43                                 if(digits->GetDigitAmp(row, col, time, fParent->GetID()) < 0) continue;
44                                 fData.SetDataUnchecked(row, col, time, digits->GetDigitAmp(row, col, time, fParent->GetID()));
45         }
46 }
47
48 //________________________________________________________
49 void AliEveTRDDigits::ComputeRepresentation()
50 {
51   // Calculate digits representation according to user settings. The
52   // user can set the following parameters:
53   // - digits scale (log/lin)
54   // - digits threshold
55   // - digits apparence (quads/boxes)
56
57   TEveQuadSet::Reset(TEveQuadSet::kQT_FreeQuad, kTRUE, 64);
58   // MT fBoxes.fBoxes.clear();
59
60   Double_t colSize, rowSize, scale;
61   Double_t x, y, z;
62
63   Int_t charge;
64   Float_t t0;
65   Float_t timeBinSize;
66
67   AliTRDcalibDB* calibration = AliTRDcalibDB::Instance();
68   Double_t cloc[4][3], cglo[3];
69   Int_t color, dimension;
70   fData.Expand();
71   for (Int_t  row = 0;  row <  fParent->rowMax;  row++) {
72     rowSize = .5 * fParent->fPadPlane->GetRowSize(row);
73     z = fParent->fPadPlane->GetRowPos(row) - rowSize;
74
75     for (Int_t  col = 0;  col <  fParent->colMax;  col++) {
76       colSize = .5 * fParent->fPadPlane->GetColSize(col);
77       y = fParent->fPadPlane->GetColPos(col) - colSize;
78       t0 = calibration->GetT0(fParent->fDet, col, row);
79       timeBinSize = calibration->GetVdrift(fParent->fDet, col, row)/fParent->samplingFrequency;
80
81       for (Int_t time = 0; time < fParent->timeMax; time++) {
82         charge = fData.GetDataUnchecked(row, col, time);
83         if (charge < fParent->GetDigitsThreshold()) continue;
84
85         x = fParent->fX0 - (time+0.5-t0)*timeBinSize;
86         scale = fParent->GetDigitsLog() ? TMath::Log(float(charge))/TMath::Log(1024.) : charge/1024.;
87         color  = 50+int(scale*50.);
88
89         cloc[0][2] = z - rowSize * scale;
90         cloc[0][1] = y - colSize * scale;
91         cloc[0][0] = x;
92
93         cloc[1][2] = z - rowSize * scale;
94         cloc[1][1] = y + colSize * scale;
95         cloc[1][0] = x;
96
97         cloc[2][2] = z + rowSize * scale;
98         cloc[2][1] = y + colSize * scale;
99         cloc[2][0] = x;
100
101         cloc[3][2] = z + rowSize * scale;
102         cloc[3][1] = y - colSize * scale;
103         cloc[3][0] = x;
104
105         Float_t* p = 0;
106         if( fParent->GetDigitsBox()){
107           // MT fBoxes.fBoxes.push_back(Box());
108           // MT fBoxes.fBoxes.back().color[0] = (UChar_t)color;
109           // MT fBoxes.fBoxes.back().color[1] = (UChar_t)color;
110           // MT fBoxes.fBoxes.back().color[2] = (UChar_t)color;
111           // MT fBoxes.fBoxes.back().color[3] = (UChar_t)color;
112           // MT p = fBoxes.fBoxes.back().vertices;
113           dimension = 2;
114         } else {
115           AddQuad((Float_t*)0);
116           QuadColor(color);
117           p = ((QFreeQuad_t*) fLastDigit)->fVertices;
118           dimension = 1;
119         }
120
121         for(int id=0; id<dimension; id++)
122           for (Int_t ic = 0; ic < 4; ic++) {
123             cloc[ic][0] -= .5 * id * timeBinSize;
124             fParent->fGeo->RotateBack(fParent->fDet,cloc[ic],cglo);
125             p[0] = cglo[0]; p[1] = cglo[1]; p[2] = cglo[2];
126             p+=3;
127           }
128       }  // end time loop
129     }  // end col loop
130   }  // end row loop
131   fData.Compress(1);
132 }
133
134 //________________________________________________________
135 void AliEveTRDDigits::Paint(Option_t *option)
136 {
137         if(fParent->GetDigitsBox()) fBoxes.Paint(option);
138         else TEveQuadSet::Paint(option);
139 }
140
141 //________________________________________________________
142 void AliEveTRDDigits::Reset()
143 {
144         TEveQuadSet::Reset(TEveQuadSet::kQT_FreeQuad, kTRUE, 64);
145         // MT fBoxes.fBoxes.clear();
146         fData.Reset();
147 }
148
149 ///////////////////////////////////////////////////////////
150 /////////////   AliEveTRDHits               /////////////////////
151 ///////////////////////////////////////////////////////////
152
153 //________________________________________________________
154 AliEveTRDHits::AliEveTRDHits(AliEveTRDChamber *p):TEvePointSet("hits", 20), fParent(p)
155 {}
156
157 //________________________________________________________
158 void AliEveTRDHits::PointSelected(Int_t n)
159 {
160         fParent->SpawnEditor();
161         AliTRDhit *h = dynamic_cast<AliTRDhit*>(GetPointId(n));
162         printf("\nDetector             : %d\n", h->GetDetector());
163         printf("Region of production : %c\n", h->FromAmplification() ? 'A' : 'D');
164         printf("TR photon            : %s\n", h->FromTRphoton() ? "Yes" : "No");
165         printf("Charge               : %d\n", h->GetCharge());
166         printf("MC track label       : %d\n", h->GetTrack());
167         printf("Time from collision  : %f\n", h->GetTime());
168 }
169
170
171 ///////////////////////////////////////////////////////////
172 /////////////   AliEveTRDHits               /////////////////////
173 ///////////////////////////////////////////////////////////
174
175 //________________________________________________________
176 AliEveTRDClusters::AliEveTRDClusters(AliEveTRDChamber *p):AliEveTRDHits(p)
177 {}
178
179 //________________________________________________________
180 void AliEveTRDClusters::PointSelected(Int_t n)
181 {
182         fParent->SpawnEditor();
183         AliTRDcluster *c = dynamic_cast<AliTRDcluster*>(GetPointId(n));
184         printf("\nDetector             : %d\n", c->GetDetector());
185         printf("Charge               : %f\n", c->GetQ());
186         printf("Sum S                : %4.0f\n", c->GetSumS());
187         printf("Time bin             : %d\n", c->GetLocalTimeBin());
188         printf("Signals              : ");
189         Short_t *cSignals = c->GetSignals();
190         for(Int_t ipad=0; ipad<7; ipad++) printf("%d ", cSignals[ipad]); printf("\n");
191         printf("Central pad          : %d\n", c->GetPadCol());
192         printf("MC track labels      : ");
193         for(Int_t itrk=0; itrk<3; itrk++) printf("%d ", c->GetLabel(itrk)); printf("\n");
194 // Bool_t       AliCluster::GetGlobalCov(Float_t* cov) const
195 // Bool_t       AliCluster::GetGlobalXYZ(Float_t* xyz) const
196 // Float_t      AliCluster::GetSigmaY2() const
197 // Float_t      AliCluster::GetSigmaYZ() const
198 // Float_t      AliCluster::GetSigmaZ2() const
199 }
200
201 ///////////////////////////////////////////////////////////
202 /////////////   AliEveTRDHitsEditor         /////////////////////
203 ///////////////////////////////////////////////////////////
204 AliEveTRDHitsEditor::AliEveTRDHitsEditor(const TGWindow* p, Int_t width, Int_t height, UInt_t options, Pixel_t back) : TGedFrame(p, width, height, options, back)
205 {
206         MakeTitle("TRD Hits");
207
208 }
209
210 AliEveTRDHitsEditor::~AliEveTRDHitsEditor()
211 {}
212
213 void AliEveTRDHitsEditor::SetModel(TObject* obj)
214 {
215         fM = dynamic_cast<AliEveTRDHits*>(obj);
216
217 //      Float_t x, y, z;
218 //      for(int ihit=0; ihit<fM->GetN(); ihit++){
219 //              fM->GetPoint(ihit, x, y, z);
220 //              printf("%3d : x=%6.3f y=%6.3f z=%6.3f\n", ihit, x, y, z);
221 //      }
222 }
223
224 ///////////////////////////////////////////////////////////
225 /////////////   AliEveTRDDigitsEditor       /////////////////////
226 ///////////////////////////////////////////////////////////
227 AliEveTRDDigitsEditor::AliEveTRDDigitsEditor(const TGWindow* p, Int_t width, Int_t height, UInt_t options, Pixel_t back) : TGedFrame(p, width, height, options, back)
228 {
229         MakeTitle("TRD Digits");
230
231 }
232
233 AliEveTRDDigitsEditor::~AliEveTRDDigitsEditor()
234 {}
235
236 void AliEveTRDDigitsEditor::SetModel(TObject* obj)
237 {
238         fM = dynamic_cast<AliEveTRDDigits*>(obj);
239         fM->fParent->SpawnEditor();
240
241 //      printf("Chamber %d", fM->fParent->GetID());
242 //      for (Int_t  row = 0;  row <  fM->fParent->GetRowMax();  row++)
243 //              for (Int_t  col = 0;  col <  fM->fParent->GetColMax();  col++)
244 //                      for (Int_t time = 0; time < fM->fParent->GetTimeMax(); time++) {
245 //                              printf("\tA(%d %d %d) = %d\n", row, col, time, fM->fData.GetDataUnchecked(row, col, time));
246 //                      }
247 }