]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/AliHMPIDDigit.cxx
Comics run analizer
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDDigit.cxx
CommitLineData
d3da6dc4 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()
25ClassImp(AliHMPIDDigit)
26
27/*
da08475b 28Preface: all geometrical information (like left-right sides) is reported as seen from electronic side.
d3da6dc4 29
da08475b 30
31The DDL file starts with common header which size and structure is standartized and mandatory for all detectors.
d3da6dc4 32The 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.
33For HMPID the detector ID is 6 (reffered in the code as kRichRawId) while DDL indexes are from 0 to 13.
34
35Common header might be followed by the private one although HMPID has no any private header, just uses the common one.
36
da08475b 37Single 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)
d3da6dc4 39
40So the LDC -chamber-ddl map is:
da08475b 41DDL index 0 -> ch 0 left -> DDL ID 0x600 DDL index 1 -> ch 1 right -> DDL ID 0x601
42DDL index 2 -> ch 1 left -> DDL ID 0x602 DDL index 3 -> ch 2 right -> DDL ID 0x603
43DDL index 4 -> ch 2 left -> DDL ID 0x604 DDL index 5 -> ch 3 right -> DDL ID 0x605
44DDL index 6 -> ch 3 left -> DDL ID 0x606 DDL index 7 -> ch 4 right -> DDL ID 0x607
45DDL index 8 -> ch 4 left -> DDL ID 0x608 DDL index 9 -> ch 5 right -> DDL ID 0x609
46DDL index 10 -> ch 5 left -> DDL ID 0x60a DDL index 11 -> ch 6 right -> DDL ID 0x60b
47DDL index 12 -> ch 6 left -> DDL ID 0x60c DDL index 13 -> ch 7 right -> DDL ID 0x60d
d3da6dc4 48
49HMPID FEE as seen by single D-RORC is composed from a number of DILOGIC chips organized in vertical stack of rows.
da08475b 50Each DILOGIC chip serves 48 channels for the 8x6 pads Channels counted from 0 to 47.
d3da6dc4 51
52The mapping inside DILOGIC chip has the following structure (see from electronics side):
da08475b 53pady
d3da6dc4 54
da08475b 555 04 10 16 22 28 34 40 46 due to repetition in column structure we may introduce per column map:
564 02 08 14 20 26 32 38 44 pady= 0 1 2 3 4 5
573 00 06 12 18 24 30 36 42 addr= 5 3 1 0 2 4
582 01 07 13 19 25 31 37 43 or vice versa
591 03 09 15 21 27 33 39 45 addr= 0 1 2 3 4 5
600 05 11 17 23 29 35 41 47 pady= 3 2 4 1 5 0
d3da6dc4 61
62 0 1 2 3 4 5 6 7 padx
63
6410 DILOGIC chips composes so called "row" in horizontal direction (reffered in the code as kNdil), so the row is 80x6 pads structure.
65DILOGIC chips in the row are counted from right to left as seen from electronics side, from 1 to 10.
6624 rows are piled up forming the whole FEE served by single D-RORC, so one DDL sees 80x144 pads separated in 3 photocathodes.
67Rows 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
70HMPID 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//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
75Float_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
da08475b 93 if(dig.PadPcX()==-1) continue;
d3da6dc4 94 new((*pSdiLst)[iSdiCnt++]) AliHMPIDDigit(dig);
95 }
96 return qdcTot;
97}//Hit2Sdi()
98//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
99void AliHMPIDDigit::Print(Option_t*)const
100{
101// Print current digit
102// Arguments: option string not used
103// Returns: none
da08475b 104 UInt_t w32; Raw(w32);
105 Printf("(ch=%1i,pc=%1i,x=%2i,y=%2i),pos=(%7.2f,%7.2f) Q=%8.3f TID=(%5i,%5i,%5i) ddl=%i raw=0x%x (r=%2i,d=%2i,a=%2i)",
106 A2C(fPad),A2P(fPad),A2X(fPad),A2Y(fPad),LorsX(),LorsY(), Q(), fTracks[0],fTracks[1],fTracks[2],DdlIdx(),w32,Row(),Dilogic(),Addr());
d3da6dc4 107}
108//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
109void AliHMPIDDigit::PrintSize()
110{
111// Print all segmentaion related sizes
112// Arguments: none
113// Returns: none
114 Printf("-->pad =(%6.2f,%6.2f) cm dead zone %.2f cm\n"
da08475b 115 "-->PC =(%6.2f,%6.2f) cm (%3i,%3i) pads\n"
116 "-->all PCs=(%6.2f,%6.2f) cm (%3i,%3i) pads",
117 SizePadX(),SizePadY(),SizeDead(),
118 SizePcX() ,SizePcY() ,kPadPcX ,kPadPcY,
119 SizeAllX(),SizeAllY(),kPadAllX,kPadAllY);
d3da6dc4 120}
121//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
da08475b 122void AliHMPIDDigit::DrawSeg()
d3da6dc4 123{
124// Draws the picture of segmentation
125// Arguments: none
126// Returns: none
da08475b 127 TCanvas *pC=new TCanvas("pads","View from electronics side, IP is behind the picture.",1000,900);pC->ToggleEventStatus();
128 gPad->AddExec("test","AliHMPIDDigit::DrawZoom()");
d3da6dc4 129 DrawPc();
130}//TestSeg()
131//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
da08475b 132void AliHMPIDDigit::DrawZoom()
d3da6dc4 133{
134// Show info about current cursur position in status bar of the canvas
135// Arguments: none
136// Returns: none
137 TCanvas *pC=(TCanvas*)gPad;
138 TRootCanvas *pRC= (TRootCanvas*)pC->GetCanvasImp();
139 TGStatusBar *pBar=pRC->GetStatusBar();
140 pBar->SetParts(5);
141 Float_t x=gPad->AbsPixeltoX(gPad->GetEventX());
142 Float_t y=gPad->AbsPixeltoY(gPad->GetEventY());
da08475b 143 AliHMPIDDigit dig(1,0,1,x,y); UInt_t w32=0;
d3da6dc4 144 if(IsInDead(x,y))
145 pBar->SetText("Out of sensitive area",4);
da08475b 146 else{
147 Int_t ddl=dig.Raw(w32);
148 pBar->SetText(Form("(p%i,x%i,y%i) ddl=%i 0x%x (r%i,d%i,a%i) (%.2f,%.2f)",
149 dig.Pc(),dig.PadPcX(),dig.PadPcY(),
150 ddl,w32,
151 dig.Row(),dig.Dilogic(),dig.Addr(),
152 dig.LorsX(),dig.LorsY() ),4);
153 }
d3da6dc4 154 if(gPad->GetEvent()==1){
155 new TCanvas("zoom",Form("Row %i DILOGIC %i",dig.Row(),dig.Dilogic()));
156 }
157}//Zoom()
158//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
159void AliHMPIDDigit::DrawPc(Bool_t isFill)
160{
161// Utility methode draws HMPID chamber PCs on event display.
162// Arguments: none
163// Returns: none
164// y6 ---------- ----------
165// | | | |
166// | 4 | | 5 |
167// y5 ---------- ----------
168//
169// y4 ---------- ----------
170// | | | |
171// | 2 | | 3 | view from electronics side
172// y3 ---------- ----------
173//
174// y2 ---------- ----------
175// | | | |
176// | 0 | | 1 |
177// y1 ---------- ----------
178// x1 x2 x3 x4
179 gPad->Range(-5,-5,SizeAllX()+5,SizeAllY()+5);
180 Float_t x1=0,x2=SizePcX(),x3=SizePcX()+SizeDead(), x4=SizeAllX();
181 Float_t y1=0,y2=SizePcY(),y3=SizePcY()+SizeDead(),y4=2*SizePcY()+SizeDead(),y5=SizeAllY()-SizePcY(),y6=SizeAllY();
182
183 Float_t xL[5]={x1,x1,x2,x2,x1}; //clockwise
184 Float_t xR[5]={x3,x3,x4,x4,x3};
185 Float_t yD[5]={y1,y2,y2,y1,y1};
186 Float_t yC[5]={y3,y4,y4,y3,y3};
187 Float_t yU[5]={y5,y6,y6,y5,y5};
188
da08475b 189 TLatex txt; txt.SetTextSize(0.01);
d3da6dc4 190 Int_t iColLeft=29,iColRight=41;
da08475b 191 TPolyLine *pc=0; TLine *pL;
192 AliHMPIDDigit dig;
d3da6dc4 193 for(Int_t iPc=0;iPc<kPcAll;iPc++){
da08475b 194 if(iPc==4) pc=new TPolyLine(5,xL,yU); if(iPc==5) pc=new TPolyLine(5,xR,yU); //draw PCs
195 if(iPc==2) pc=new TPolyLine(5,xL,yC); if(iPc==3) pc=new TPolyLine(5,xR,yC);
196 if(iPc==0) pc=new TPolyLine(5,xL,yD); if(iPc==1) pc=new TPolyLine(5,xR,yD);
d3da6dc4 197 (iPc%2)? pc->SetFillColor(iColLeft): pc->SetFillColor(iColRight);
198 if(isFill) pc->Draw("f"); else pc->Draw();
da08475b 199 if(iPc%2) {dig.Set(0,iPc,79,25); txt.DrawText(dig.LorsX()+2,dig.LorsY(),Form("PC%i",dig.Pc()));}//print PC#
200
201 txt.SetTextAlign(32);
202 for(Int_t iRow=0;iRow<8 ;iRow++){//draw row lines (horizontal)
203 dig.Set(0,iPc,0,iRow*6); //set digit to the left-down pad of this row
204 if(iPc%2) txt.DrawText(dig.LorsX()-1 ,dig.LorsY(),Form("%i",dig.PadPcY())); //print PadY#
205 txt.DrawText(dig.LorsX()-1+(iPc%2)*67,dig.LorsY()+2,Form("r%i",dig.Row())); //print Row#
206 pL=new TLine(dig.LorsX()-0.5*SizePadX(),dig.LorsY()-0.5*SizePadY(),dig.LorsX()+SizePcX()-0.5*SizePadX(),dig.LorsY()-0.5*SizePadY());
207 if(iRow!=0) pL->Draw();
208 }//row loop
209
210 txt.SetTextAlign(13);
211 for(Int_t iDil=0;iDil<10;iDil++){//draw dilogic lines (vertical)
212 dig.Set(0,iPc,iDil*8,0); //set this digit to the left-down pad of this dilogic
213 txt.DrawText(dig.LorsX() ,dig.LorsY()-1,Form("%i",dig.PadPcX())); //print PadX#
214 if(iPc==4 || iPc==5) txt.DrawText(dig.LorsX()+2,dig.LorsY()+42,Form("d%i",dig.Dilogic())); //print Dilogic#
215 pL=new TLine(dig.LorsX()-0.5*SizePadX(),dig.LorsY()-0.5*SizePadY(),dig.LorsX()-0.5*SizePadX(),dig.LorsY()+SizePcY()-0.5*SizePadY());
216 if(iDil!=0)pL->Draw();
217 }//dilogic loop
218 }//PC loop
d3da6dc4 219}//DrawPc()
220//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
da08475b 221void AliHMPIDDigit::Test()
222{
223 AliHMPIDDigit d1,d2; Int_t ddl;UInt_t w32;
224 for(Int_t i=0;i<10000;i++){
225 Int_t ch=Int_t(gRandom->Rndm()*7);
226 Int_t pc=Int_t(gRandom->Rndm()*6);
227 Int_t px=Int_t(gRandom->Rndm()*80);
228 Int_t py=Int_t(gRandom->Rndm()*48);
229 d1.Set(ch,pc,px,py);
230 ddl=d1.Raw(w32); d2.Raw(ddl,w32);
231 if(d1.Compare(&d2)) Printf("Problem!!!");
232 }
233 Printf("OK");
234}//Test()
235//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
236