]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDDigit.cxx
CRT becomes ACORDE
[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 Any given LDC collects data from a number of D-RORC cards connected to this same LDC by separate DDLs. This data is stored in corresponding number
29 of DDL buffers.
30
31 Each DDL buffer corresponds to a single D-RORC. The data this buffer contains are hardware generated by D-RORC. The buffer starts with common header 
32 which size and structure is standartized and mandatory for all detectors.
33 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. 
34 For HMPID the detector ID is 6 (reffered in the code as kRichRawId) while DDL indexes are from 0 to 13.
35
36 Common header might be followed by the private one although  HMPID has no any private header, just uses the common one.
37
38 Single HMPID D-RORC serves one half of a chamber i.e. 3 photocathodes even LDC for left part( PCs 0-2-4) and odd LDC for right part(1-3-5) 
39 as it's seen from electronics side.
40
41 So the LDC -chamber-ddl map is:
42 DDL index  0 -> ch 1 left -> DDL ID 0x600          DDL index  1 -> ch 1 right -> DDL ID 0x601 
43 DDL index  2 -> ch 2 left -> DDL ID 0x602          DDL index  3 -> ch 2 right -> DDL ID 0x603 
44 DDL index  4 -> ch 3 left -> DDL ID 0x604          DDL index  5 -> ch 3 right -> DDL ID 0x605 
45 DDL index  6 -> ch 4 left -> DDL ID 0x606          DDL index  7 -> ch 4 right -> DDL ID 0x607 
46 DDL index  8 -> ch 5 left -> DDL ID 0x608          DDL index  9 -> ch 5 right -> DDL ID 0x609 
47 DDL index 10 -> ch 6 left -> DDL ID 0x60a          DDL index 11 -> ch 6 right -> DDL ID 0x60b 
48 DDL index 12 -> ch 7 left -> DDL ID 0x60c          DDL index 13 -> ch 7 right -> DDL ID 0x60d 
49
50 HMPID FEE as seen by single D-RORC is composed from a number of DILOGIC chips organized in vertical stack of rows. 
51 Each DILOGIC chip serves 48 channels for the 8x6 pads (reffered in the code as kDilX,kDilY). Channels counted from 0 to 47.
52
53 The mapping inside DILOGIC chip has the following structure (see from electronics side):
54
55 5  04 10 16 22 28 34 40 46
56 4  02 08 14 20 26 32 38 44
57 3  00 06 12 18 24 30 36 42
58 2  01 07 13 19 25 31 37 43
59 1  03 09 15 21 27 33 39 45
60 0  05 11 17 23 29 35 41 47
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 Float_t 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   Float_t  x=(pHit->LorsX() > SizePcX())? pHit->LorsX()-SizePcX()-SizeDead():pHit->LorsX();          //sagita is for PC (0-64) and not for chamber   
82   Float_t  qdcEle=34.06311+0.2337070*x+5.807476e-3*x*x-2.956471e-04*x*x*x+2.310001e-06*x*x*x*x;      //reparametrised from DiMauro
83   
84   Int_t iNele=Int_t(pHit->E()/26e-9);  if(iNele<1) iNele=1;                                          //number of electrons created by hit
85   Float_t qdcTot=0; for(Int_t i=1;i<=iNele;i++) qdcTot-=qdcEle*TMath::Log(gRandom->Rndm()+1e-6);     //total qdc fro hit, 1e-6 is a protection against 0 from rndm
86   
87   AliHMPIDDigit dd(1,1,1,pHit->LorsX(),pHit->LorsY());                                                //tmp digit to shift hit y to the nearest anod wire  
88   Float_t y=  (pHit->LorsY() > dd.LorsY()) ? dd.LorsY()+0.21 :  dd.LorsY()-0.21;
89   
90   Int_t iSdiCnt=pSdiLst->GetEntries();
91   for(Int_t i=0;i<9;i++){                                      //affected pads loop
92     AliHMPIDDigit dig(pHit->Ch(),qdcTot,pHit->Tid(),pHit->LorsX(),y,i); //c,q,tid,x,y   create tmp sdigit for pad i around hit position
93     if(dig.PadX()==-1) continue;
94     new((*pSdiLst)[iSdiCnt++]) AliHMPIDDigit(dig);
95   }
96   return qdcTot;
97 }//Hit2Sdi()
98 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
99 void AliHMPIDDigit::Print(Option_t*)const
100 {
101 // Print current digit  
102 // Arguments: option string not used
103 //   Returns: none    
104   Printf("pad=(%2i,%2i,%3i,%3i),pos=(%7.2f,%7.2f) QDC=%8.3f, TID=(%5i,%5i,%5i) raw r=%2i d=%2i a=%2i",
105                Ch(),Pc(),PadX(),PadY(),LorsX(),LorsY(), Q(),  fTracks[0],fTracks[1],fTracks[2], Row(), Dilogic(), Addr());
106 }
107 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
108 void AliHMPIDDigit::PrintSize()
109 {
110 // Print all segmentaion related sizes
111 // Arguments: none
112 //   Returns: none    
113   Printf("-->pad    =(%6.2f,%6.2f) cm dead zone %.2f cm\n"
114          "-->PC     =(%6.2f,%6.2f) cm\n"
115          "-->all PCs=(%6.2f,%6.2f) cm",
116                SizePadX(),SizePadY(),SizeDead(),SizePcX(),SizePcY(),SizeAllX(),SizeAllY());
117 }
118 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
119 void AliHMPIDDigit::TestSeg()
120 {
121 // Draws the picture of segmentation   
122 // Arguments: none
123 //   Returns: none     
124   TCanvas *pC=new TCanvas("pads","View from electronics side, IP is behind the picture.");pC->ToggleEventStatus();
125   gPad->AddExec("test","AliHMPIDDigit::Zoom()");
126   DrawPc();  
127 }//TestSeg()
128 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
129 void AliHMPIDDigit::Zoom()
130 {
131 // Show info about current cursur position in status bar of the canvas
132 // Arguments: none
133 //   Returns: none     
134   TCanvas *pC=(TCanvas*)gPad; 
135   TRootCanvas *pRC= (TRootCanvas*)pC->GetCanvasImp();
136   TGStatusBar *pBar=pRC->GetStatusBar();
137   pBar->SetParts(5);
138   Float_t x=gPad->AbsPixeltoX(gPad->GetEventX());
139   Float_t y=gPad->AbsPixeltoY(gPad->GetEventY());
140   AliHMPIDDigit dig(1,100,1,x,y);  
141   if(IsInDead(x,y))
142     pBar->SetText("Out of sensitive area",4);    
143   else
144     pBar->SetText(Form("p%i x%i y%i r%i d%i a%i (%.2f,%.2f)",dig.Pc(),dig.PadX(),dig.PadY(),dig.Row(),dig.Dilogic(),dig.Addr(),dig.LorsX(),dig.LorsY()),4);
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 t; t.SetTextSize(0.01); t.SetTextAlign(22);
181   Int_t iColLeft=29,iColRight=41;
182   TPolyLine *pc=0;  TLine *pL; Float_t x0=0,y0=0,wRow=kDilY*SizePadY(),wDil=kDilX*SizePadX();
183   for(Int_t iPc=0;iPc<kPcAll;iPc++){
184     if(iPc==4) {pc=new TPolyLine(5,xL,yU);t.DrawText(x1-3.5,y6-20,"PC4");x0=x1;y0=y5;} if(iPc==5) {pc=new TPolyLine(5,xR,yU);t.DrawText(x4+3.5,y6-20,"PC5");x0=x3;y0=y5;}
185     if(iPc==2) {pc=new TPolyLine(5,xL,yC);t.DrawText(x1-3.5,y4-20,"PC2");x0=x1;y0=y3;} if(iPc==3) {pc=new TPolyLine(5,xR,yC);t.DrawText(x4+3.5,y4-20,"PC3");x0=x3;y0=y3;}
186     if(iPc==0) {pc=new TPolyLine(5,xL,yD);t.DrawText(x1-3.5,y2-20,"PC0");x0=x1;y0=y1;} if(iPc==1) {pc=new TPolyLine(5,xR,yD);t.DrawText(x4+3.5,y2-20,"PC1");x0=x3;y0=y1;}    
187     (iPc%2)? pc->SetFillColor(iColLeft): pc->SetFillColor(iColRight);
188     if(isFill) pc->Draw("f"); else pc->Draw();
189     for(Int_t i=1;i<=8 ;i++){//draw row lines (horizontal)
190       Float_t y=y0+i*wRow;
191       Int_t row=i+iPc/2*8; if(iPc%2!=0) row=25-row; t.DrawText(x0-1,y -3,Form("r%i",row)); 
192       if(i==8) break;                                      //do not draw the last line of PC
193       pL=new TLine(x0,y,x0+SizePcX(),y);   pL->Draw(); 
194     }    
195     for(Int_t iDil=1;iDil<=10;iDil++){Float_t x=x0+iDil*wDil;t.DrawText(x -3,y0-1,Form("d%i",11-iDil));   if(iDil==10) break;    pL=new TLine(x,y0,x,y0+SizePcY());  pL->Draw();}  
196   }      
197 }//DrawPc()
198 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++