]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDDigit.cxx
Bugs in clusterizing and recon fixed
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDDigit.cxx
1 //  **************************************************************************
2 //  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 //  *                                                                        *
4 //  * Author: The ALICE Off-line Project.                                    *
5 //  * Contributors are mentioned in the code where appropriate.              *
6 //  *                                                                        *
7 //  * Permission to use, copy, modify and distribute this software and its   *
8 //  * documentation strictly for non-commercial purposes is hereby granted   *
9 //  * without fee, provided that the above copyright notice appears in all   *
10 //  * copies and that both the copyright notice and this permission notice   *
11 //  * appear in the supporting documentation. The authors make no claims     *
12 //  * about the suitability of this software for any purpose. It is          *
13 //  * provided "as is" without express or implied warranty.                  *
14 //  **************************************************************************
15
16 #include "AliHMPIDDigit.h" //class header
17 #include "AliHMPIDHit.h"   //Hit2Sdi() 
18 #include <TClonesArray.h> //Hit2Sdi() 
19 #include <TCanvas.h>      //TestSeg()
20 #include <TLine.h>        //TestSeg()
21 #include <TLatex.h>       //TestSeg()
22 #include <TPolyLine.h>    //DrawPc() 
23 #include <TGStatusBar.h>  //Zoom()
24 #include <TRootCanvas.h>  //Zoom() 
25 ClassImp(AliHMPIDDigit)
26
27 /*
28 Preface: all geometrical information (like left-right sides) is reported as seen from electronic side.
29
30
31 The DDL file starts with common header which size and structure is standartized and mandatory for all detectors. 
32 The header contains among other words, so called Equipment ID word. This unique value for each D-RORC is calculated as detector ID << 8 + DDL index. 
33 For HMPID the detector ID is 6 (reffered in the code as kRichRawId) while DDL indexes are from 0 to 13.
34
35 Common header might be followed by the private one although  HMPID has no any private header, just uses the common one.
36
37 Single HMPID D-RORC (with 2 channels) serves a single chamber so that channel 0 serves left half (PCs 0-2-4) 
38                                                                               1 serves right half(PCs 1-3-5) 
39
40 So the LDC -chamber-ddl map is:
41 DDL index  0 -> ch 0 left -> DDL ID 0x600          DDL index  1 -> ch 1 right -> DDL ID 0x601 
42 DDL index  2 -> ch 1 left -> DDL ID 0x602          DDL index  3 -> ch 2 right -> DDL ID 0x603 
43 DDL index  4 -> ch 2 left -> DDL ID 0x604          DDL index  5 -> ch 3 right -> DDL ID 0x605 
44 DDL index  6 -> ch 3 left -> DDL ID 0x606          DDL index  7 -> ch 4 right -> DDL ID 0x607 
45 DDL index  8 -> ch 4 left -> DDL ID 0x608          DDL index  9 -> ch 5 right -> DDL ID 0x609 
46 DDL index 10 -> ch 5 left -> DDL ID 0x60a          DDL index 11 -> ch 6 right -> DDL ID 0x60b 
47 DDL index 12 -> ch 6 left -> DDL ID 0x60c          DDL index 13 -> ch 7 right -> DDL ID 0x60d 
48
49 HMPID FEE as seen by single D-RORC is composed from a number of DILOGIC chips organized in vertical stack of rows. 
50 Each DILOGIC chip serves 48 channels for the 8x6 pads Channels counted from 0 to 47.
51
52 The mapping inside DILOGIC chip has the following structure (see from electronics side):
53 pady
54
55 5  04 10 16 22 28 34 40 46                   due to repetition in column structure we may introduce per column map:   
56 4  02 08 14 20 26 32 38 44                   pady= 0 1 2 3 4 5          
57 3  00 06 12 18 24 30 36 42                   addr= 5 3 1 0 2 4
58 2  01 07 13 19 25 31 37 43                   or vice versa 
59 1  03 09 15 21 27 33 39 45                   addr= 0 1 2 3 4 5
60 0  05 11 17 23 29 35 41 47                   pady= 3 2 4 1 5 0  
61     
62     0  1  2  3  4  5  6  7  padx
63
64 10 DILOGIC chips composes so called "row" in horizontal direction (reffered in the code as kNdil), so the row is 80x6 pads structure. 
65 DILOGIC chips in the row are counted from right to left as seen from electronics side, from 1 to 10.
66 24 rows are piled up forming the whole FEE served by single D-RORC, so one DDL sees 80x144 pads separated in 3 photocathodes.
67 Rows are counted from 1 to 24 from top    to bottom for right  half of the chamber (PCs 1-3-5) as seen from electronics side, meaning even LDC number
68                           and from bottom to top    for left   half of the chamber (PCs 0-2-4) as seen from electronics side, meaning odd LDC number.
69
70 HMPID raw word is 32 bits with the structure:   
71    00000             rrrrr                      dddd                               aaaaaa                          qqqqqqqqqqqq 
72  5 bits zero  5 bits row number (1..24)  4 bits DILOGIC chip number (1..10) 6 bits DILOGIC address (0..47)  12 bits QDC value (0..4095)
73 */
74 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
75 void AliHMPIDDigit::Hit2Sdi(AliHMPIDHit *pHit,TClonesArray *pSdiLst)
76 {
77 // Creates a list of sdigits out of provided hit
78 // Arguments: pHit- hit
79 //   Returns: none
80   
81   Int_t iSdiCnt=pSdiLst->GetEntries(); //list of sdigits contains sdigits from previous ivocations of Hit2Sdi, do not override them
82   AliHMPIDDigit dig;
83   for(Int_t i=0;i<9;i++){                                      //affected pads loop
84     dig.Set(pHit,i); //c,q,tid,x,y   create tmp sdigit for pad i around hit position
85     if(dig.PadPcX()==-1) continue;
86     new((*pSdiLst)[iSdiCnt++]) AliHMPIDDigit(dig);
87   }
88 }//Hit2Sdi()
89 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
90 void AliHMPIDDigit::Print(Option_t*)const
91 {
92 // Print current digit  
93 // Arguments: option string not used
94 //   Returns: none    
95   UInt_t w32; Raw(w32);
96   Printf("(ch=%1i,pc=%1i,x=%2i,y=%2i) (%7.3f,%7.3f) Q=%8.3f TID=(%5i,%5i,%5i) ddl=%i raw=0x%x (r=%2i,d=%2i,a=%2i) %s",
97                A2C(fPad),A2P(fPad),A2X(fPad),A2Y(fPad),LorsX(),LorsY(), Q(),  fTracks[0],fTracks[1],fTracks[2],DdlIdx(),w32,Row(),Dilogic(),Addr(), (IsOverTh(Q()))?"":"below thr");
98 }
99 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
100 void AliHMPIDDigit::PrintSize()
101 {
102 // Print all segmentaion related sizes
103 // Arguments: none
104 //   Returns: none    
105   Printf("-->pad    =(%6.2f,%6.2f) cm dead zone %.2f cm\n"
106          "-->PC     =(%6.2f,%6.2f) cm (%3i,%3i) pads\n"
107          "-->all PCs=(%6.2f,%6.2f) cm (%3i,%3i) pads",
108                SizePadX(),SizePadY(),SizeDead(),
109                SizePcX() ,SizePcY() ,kPadPcX ,kPadPcY,
110                SizeAllX(),SizeAllY(),kPadAllX,kPadAllY);
111 }
112 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
113 void AliHMPIDDigit::DrawSeg()
114 {
115 // Draws the picture of segmentation   
116 // Arguments: none
117 //   Returns: none     
118   TCanvas *pC=new TCanvas("pads","View from electronics side, IP is behind the picture.",1000,900);pC->ToggleEventStatus();
119   gPad->AddExec("test","AliHMPIDDigit::DrawZoom()");
120   DrawPc();  
121 }//TestSeg()
122 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
123 void AliHMPIDDigit::DrawZoom()
124 {
125 // Show info about current cursur position in status bar of the canvas
126 // Arguments: none
127 //   Returns: none     
128   TCanvas *pC=(TCanvas*)gPad; 
129   TRootCanvas *pRC= (TRootCanvas*)pC->GetCanvasImp();
130   TGStatusBar *pBar=pRC->GetStatusBar();
131   pBar->SetParts(5);
132   Float_t x=gPad->AbsPixeltoX(gPad->GetEventX());
133   Float_t y=gPad->AbsPixeltoY(gPad->GetEventY());
134   AliHMPIDDigit dig;dig.Manual1(1,x,y); UInt_t w32=0; 
135   if(IsInDead(x,y))
136     pBar->SetText("Out of sensitive area",4);    
137   else{
138     Int_t ddl=dig.Raw(w32);
139     pBar->SetText(Form("(p%i,x%i,y%i) ddl=%i 0x%x (r%i,d%i,a%i) (%.2f,%.2f)",
140         dig.Pc(),dig.PadPcX(),dig.PadPcY(),
141         ddl,w32,
142         dig.Row(),dig.Dilogic(),dig.Addr(),
143         dig.LorsX(),dig.LorsY()                            ),4);
144   }    
145   if(gPad->GetEvent()==1){
146     new TCanvas("zoom",Form("Row %i DILOGIC %i",dig.Row(),dig.Dilogic()));  
147   }
148 }//Zoom()
149 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
150 void AliHMPIDDigit::DrawPc(Bool_t isFill) 
151
152 // Utility methode draws HMPID chamber PCs on event display.
153 // Arguments: none
154 //   Returns: none      
155 //  y6  ----------  ----------
156 //      |        |  |        |
157 //      |    4   |  |    5   |
158 //  y5  ----------  ----------
159 //
160 //  y4  ----------  ----------
161 //      |        |  |        |
162 //      |    2   |  |    3   |   view from electronics side
163 //  y3  ----------  ----------
164 //          
165 //  y2  ----------  ----------
166 //      |        |  |        |
167 //      |    0   |  |    1   |
168 //  y1  ----------  ----------
169 //      x1      x2  x3       x4
170   gPad->Range(-5,-5,SizeAllX()+5,SizeAllY()+5); 
171   Float_t x1=0,x2=SizePcX(),x3=SizePcX()+SizeDead(),                                                  x4=SizeAllX();
172   Float_t y1=0,y2=SizePcY(),y3=SizePcY()+SizeDead(),y4=2*SizePcY()+SizeDead(),y5=SizeAllY()-SizePcY(),y6=SizeAllY();
173
174   Float_t xL[5]={x1,x1,x2,x2,x1}; //clockwise
175   Float_t xR[5]={x3,x3,x4,x4,x3};  
176   Float_t yD[5]={y1,y2,y2,y1,y1};
177   Float_t yC[5]={y3,y4,y4,y3,y3};  
178   Float_t yU[5]={y5,y6,y6,y5,y5};
179     
180   TLatex txt; txt.SetTextSize(0.01);
181   Int_t iColLeft=29,iColRight=41;
182   TPolyLine *pc=0;  TLine *pL;
183   AliHMPIDDigit dig;
184   for(Int_t iPc=0;iPc<kPcAll;iPc++){
185     if(iPc==4) pc=new TPolyLine(5,xL,yU); if(iPc==5) pc=new TPolyLine(5,xR,yU); //draw PCs
186     if(iPc==2) pc=new TPolyLine(5,xL,yC); if(iPc==3) pc=new TPolyLine(5,xR,yC);
187     if(iPc==0) pc=new TPolyLine(5,xL,yD); if(iPc==1) pc=new TPolyLine(5,xR,yD);
188     (iPc%2)? pc->SetFillColor(iColLeft): pc->SetFillColor(iColRight);
189     if(isFill) pc->Draw("f"); else pc->Draw();
190     if(iPc%2) {dig.Manual2(0,iPc,79,25); txt.DrawText(dig.LorsX()+2,dig.LorsY(),Form("PC%i",dig.Pc()));}//print PC#    
191     
192     txt.SetTextAlign(32);
193     for(Int_t iRow=0;iRow<8 ;iRow++){//draw row lines (horizontal)
194       dig.Manual2(0,iPc,0,iRow*6);   //set digit to the left-down pad of this row
195       if(iPc%2) txt.DrawText(dig.LorsX()-1           ,dig.LorsY(),Form("%i",dig.PadPcY())); //print PadY#    
196                 txt.DrawText(dig.LorsX()-1+(iPc%2)*67,dig.LorsY()+2,Form("r%i",dig.Row())); //print Row#    
197       pL=new TLine(dig.LorsX()-0.5*SizePadX(),dig.LorsY()-0.5*SizePadY(),dig.LorsX()+SizePcX()-0.5*SizePadX(),dig.LorsY()-0.5*SizePadY()); 
198       if(iRow!=0) pL->Draw(); 
199     }//row loop  
200     
201     txt.SetTextAlign(13);
202     for(Int_t iDil=0;iDil<10;iDil++){//draw dilogic lines (vertical)
203       dig.Manual2(0,iPc,iDil*8,0);       //set this digit to the left-down pad of this dilogic        
204                            txt.DrawText(dig.LorsX()  ,dig.LorsY()-1,Form("%i",dig.PadPcX()));   //print PadX# 
205       if(iPc==4 || iPc==5) txt.DrawText(dig.LorsX()+2,dig.LorsY()+42,Form("d%i",dig.Dilogic())); //print Dilogic#    
206       pL=new TLine(dig.LorsX()-0.5*SizePadX(),dig.LorsY()-0.5*SizePadY(),dig.LorsX()-0.5*SizePadX(),dig.LorsY()+SizePcY()-0.5*SizePadY()); 
207       if(iDil!=0)pL->Draw();
208     }//dilogic loop        
209   }//PC loop      
210 }//DrawPc()
211 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
212 void AliHMPIDDigit::Test()
213 {
214   AliHMPIDDigit d1,d2; Int_t ddl;UInt_t w32;
215   for(Int_t i=0;i<10000;i++){
216     Int_t ch=Int_t(gRandom->Rndm()*7);
217     Int_t pc=Int_t(gRandom->Rndm()*6);
218     Int_t px=Int_t(gRandom->Rndm()*80);
219     Int_t py=Int_t(gRandom->Rndm()*48);
220     d1.Manual2(ch,pc,px,py);                
221     ddl=d1.Raw(w32);    d2.Raw(ddl,w32);
222     if(d1.Compare(&d2)) Printf("Problem!!!");
223   }
224   Printf("OK");
225 }//Test()
226 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
227