]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MFT/AliMFTPlane.cxx
end-of-line normalization
[u/mrichter/AliRoot.git] / MFT / AliMFTPlane.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 //====================================================================================================================================================
17 //
18 //      Class for the description of the structure for the planes of the ALICE Muon Forward Tracker
19 //
20 //      Contact author: antonio.uras@cern.ch
21 //
22 //====================================================================================================================================================
23
24 #include "TNamed.h"
25 #include "THnSparse.h"
26 #include "TClonesArray.h"
27 #include "TAxis.h"
28 #include "TPave.h"
29 #include "TCanvas.h"
30 #include "TH2D.h"
31 #include "TEllipse.h"
32 #include "TMath.h"
33 #include "AliLog.h"
34 #include "AliMFTConstants.h"
35 #include "AliMFTPlane.h"
36
37 const Double_t AliMFTPlane::fRadiusMin           = AliMFTConstants::fRadiusMin;
38 const Double_t AliMFTPlane::fActiveSuperposition = AliMFTConstants::fActiveSuperposition;
39 const Double_t AliMFTPlane::fHeightActive        = AliMFTConstants::fHeightActive;
40 const Double_t AliMFTPlane::fHeightReadout       = AliMFTConstants::fHeightReadout;
41 const Double_t AliMFTPlane::fSupportExtMargin    = AliMFTConstants::fSupportExtMargin;
42
43 ClassImp(AliMFTPlane)
44
45 //====================================================================================================================================================
46
47 AliMFTPlane::AliMFTPlane():
48   TNamed(),
49   fPlaneNumber(-1),
50   fZCenter(0), 
51   fRMinSupport(0), 
52   fRMax(0),
53   fRMaxSupport(0),
54   fPixelSizeX(0), 
55   fPixelSizeY(0), 
56   fThicknessActive(0), 
57   fThicknessSupport(0), 
58   fThicknessReadout(0),
59   fZCenterActiveFront(0),
60   fZCenterActiveBack(0),
61   fEquivalentSilicon(0),
62   fEquivalentSiliconBeforeFront(0),
63   fEquivalentSiliconBeforeBack(0),
64   fActiveElements(0),
65   fReadoutElements(0),
66   fSupportElements(0),
67   fHasPixelRectangularPatternAlongY(kFALSE),
68   fPlaneIsOdd(kFALSE)
69 {
70
71   // default constructor
72
73 }
74
75 //====================================================================================================================================================
76
77 AliMFTPlane::AliMFTPlane(const Char_t *name, const Char_t *title):
78   TNamed(name, title),
79   fPlaneNumber(-1),
80   fZCenter(0), 
81   fRMinSupport(0), 
82   fRMax(0),
83   fRMaxSupport(0),
84   fPixelSizeX(0), 
85   fPixelSizeY(0), 
86   fThicknessActive(0), 
87   fThicknessSupport(0), 
88   fThicknessReadout(0),
89   fZCenterActiveFront(0),
90   fZCenterActiveBack(0),
91   fEquivalentSilicon(0),
92   fEquivalentSiliconBeforeFront(0),
93   fEquivalentSiliconBeforeBack(0),
94   fActiveElements(0),
95   fReadoutElements(0),
96   fSupportElements(0),
97   fHasPixelRectangularPatternAlongY(kFALSE),
98   fPlaneIsOdd(kFALSE)
99 {
100
101   // constructor
102   fActiveElements  = new TClonesArray("THnSparseC");
103   fReadoutElements = new TClonesArray("THnSparseC");
104   fSupportElements = new TClonesArray("THnSparseC");
105   fActiveElements->SetOwner(kTRUE);
106   fReadoutElements->SetOwner(kTRUE);
107   fSupportElements->SetOwner(kTRUE);
108   
109 }
110
111 //====================================================================================================================================================
112
113 AliMFTPlane::AliMFTPlane(const AliMFTPlane& plane):
114   TNamed(plane),
115   fPlaneNumber(plane.fPlaneNumber),
116   fZCenter(plane.fZCenter), 
117   fRMinSupport(plane.fRMinSupport), 
118   fRMax(plane.fRMax),
119   fRMaxSupport(plane.fRMaxSupport),
120   fPixelSizeX(plane.fPixelSizeX), 
121   fPixelSizeY(plane.fPixelSizeY), 
122   fThicknessActive(plane.fThicknessActive), 
123   fThicknessSupport(plane.fThicknessSupport), 
124   fThicknessReadout(plane.fThicknessReadout),
125   fZCenterActiveFront(plane.fZCenterActiveFront),
126   fZCenterActiveBack(plane.fZCenterActiveBack),
127   fEquivalentSilicon(plane.fEquivalentSilicon),
128   fEquivalentSiliconBeforeFront(plane.fEquivalentSiliconBeforeFront),
129   fEquivalentSiliconBeforeBack(plane.fEquivalentSiliconBeforeBack),
130   fActiveElements(0),
131   fReadoutElements(0),
132   fSupportElements(0),
133   fHasPixelRectangularPatternAlongY(plane.fHasPixelRectangularPatternAlongY),
134   fPlaneIsOdd(plane.fPlaneIsOdd)
135 {
136
137   // copy constructor
138   fActiveElements  = new TClonesArray(*(plane.fActiveElements));
139   fActiveElements  -> SetOwner(kTRUE);
140   fReadoutElements = new TClonesArray(*(plane.fReadoutElements));
141   fReadoutElements -> SetOwner(kTRUE);
142   fSupportElements = new TClonesArray(*(plane.fSupportElements));
143   fSupportElements -> SetOwner(kTRUE);
144
145         
146 }
147
148 //====================================================================================================================================================
149
150 AliMFTPlane::~AliMFTPlane() {
151
152   AliInfo("Delete AliMFTPlane");
153   if(fActiveElements) fActiveElements->Delete();
154   delete fActiveElements; 
155   if(fReadoutElements) fReadoutElements->Delete();
156   delete fReadoutElements; 
157   if(fSupportElements) fSupportElements->Delete();
158   delete fSupportElements; 
159
160 }
161
162 //====================================================================================================================================================
163
164 void AliMFTPlane::Clear(const Option_t* /*opt*/) {
165
166   AliInfo("Clear AliMFTPlane");
167   if(fActiveElements) fActiveElements->Delete();
168   delete fActiveElements; fActiveElements=NULL;
169   if(fReadoutElements) fReadoutElements->Delete();
170   delete fReadoutElements;  fReadoutElements=NULL; 
171   if(fSupportElements) fSupportElements->Delete();
172   delete fSupportElements;   fSupportElements=NULL;
173
174 }
175
176 //====================================================================================================================================================
177
178 AliMFTPlane& AliMFTPlane::operator=(const AliMFTPlane& plane) {
179
180   // Assignment operator
181   
182   // check assignement to self
183   if (this != &plane) {
184     
185     // base class assignement
186     TNamed::operator=(plane);
187     
188     // clear memory
189     Clear("");
190     
191     fPlaneNumber                      = plane.fPlaneNumber;
192     fZCenter                          = plane.fZCenter; 
193     fRMinSupport                      = plane.fRMinSupport; 
194     fRMax                             = plane.fRMax;
195     fRMaxSupport                      = plane.fRMaxSupport;
196     fPixelSizeX                       = plane.fPixelSizeX;
197     fPixelSizeY                       = plane.fPixelSizeY; 
198     fThicknessActive                  = plane.fThicknessActive; 
199     fThicknessSupport                 = plane.fThicknessSupport; 
200     fThicknessReadout                 = plane.fThicknessReadout;
201     fZCenterActiveFront               = plane.fZCenterActiveFront;
202     fZCenterActiveBack                = plane.fZCenterActiveBack;
203     fEquivalentSilicon                = plane.fEquivalentSilicon;
204     fEquivalentSiliconBeforeFront     = plane.fEquivalentSiliconBeforeFront;
205     fEquivalentSiliconBeforeBack      = plane.fEquivalentSiliconBeforeBack;
206     fActiveElements = new TClonesArray(*(plane.fActiveElements));
207     fActiveElements -> SetOwner(kTRUE);
208     fReadoutElements = new TClonesArray(*(plane.fReadoutElements));
209     fReadoutElements -> SetOwner(kTRUE);
210     fSupportElements = new TClonesArray(*(plane.fSupportElements));
211     fSupportElements -> SetOwner(kTRUE);
212     fHasPixelRectangularPatternAlongY = plane.fHasPixelRectangularPatternAlongY;
213     fPlaneIsOdd                       = plane.fPlaneIsOdd;
214
215   }
216   
217   return *this;
218   
219 }
220
221 //====================================================================================================================================================
222
223 Bool_t AliMFTPlane::Init(Int_t    planeNumber,
224                          Double_t zCenter, 
225                          Double_t rMin, 
226                          Double_t rMax, 
227                          Double_t pixelSizeX, 
228                          Double_t pixelSizeY, 
229                          Double_t thicknessActive, 
230                          Double_t thicknessSupport, 
231                          Double_t thicknessReadout,
232                          Bool_t   hasPixelRectangularPatternAlongY) {
233
234   AliDebug(1, Form("Initializing Plane Structure for Plane %s", GetName()));
235
236   fPlaneNumber      = planeNumber;
237   fZCenter          = zCenter;
238   fRMinSupport      = rMin;
239   fRMax             = rMax;
240   fPixelSizeX       = pixelSizeX;
241   fPixelSizeY       = pixelSizeY;
242   fThicknessActive  = thicknessActive;
243   fThicknessSupport = thicknessSupport;
244   fThicknessReadout = thicknessReadout;
245
246   fHasPixelRectangularPatternAlongY = hasPixelRectangularPatternAlongY;
247
248   fZCenterActiveFront = fZCenter - 0.5*fThicknessSupport - 0.5*fThicknessActive;
249   fZCenterActiveBack  = fZCenter + 0.5*fThicknessSupport + 0.5*fThicknessActive;
250
251 //   if (fRMinSupport <= fRadiusMin) fRMinSupport = fRadiusMin;
252 //   else {
253 //     fRMinSupport = fRadiusMin + (fHeightActive-fActiveSuperposition) * Int_t((fRMinSupport-fRadiusMin)/(fHeightActive-fActiveSuperposition));
254 //   }
255   
256   if (fRMax < fRMinSupport+fHeightActive) fRMax = fRMinSupport + fHeightActive;
257
258   Int_t nLaddersWithinPipe = Int_t(fRMinSupport/(fHeightActive-fActiveSuperposition));
259   if (fRMinSupport-nLaddersWithinPipe*(fHeightActive-fActiveSuperposition) > 0.5*(fHeightActive-2*fActiveSuperposition)) fPlaneIsOdd = kTRUE;
260   else fPlaneIsOdd = kFALSE;
261
262   fRMax = fRMinSupport + (fHeightActive-fActiveSuperposition) * 
263     (Int_t((fRMax-fRMinSupport-fHeightActive)/(fHeightActive-fActiveSuperposition))+1) + fHeightActive;
264
265   fRMaxSupport = TMath::Sqrt(fHeightActive*(2.*rMax-fHeightActive) + fRMax*fRMax) + fSupportExtMargin;
266    
267   return kTRUE;
268  
269 }
270
271 //====================================================================================================================================================
272
273 Bool_t AliMFTPlane::CreateStructure() {
274
275   Int_t nBins[3]={0};
276   Double_t minPosition[3]={0}, maxPosition[3]={0};
277   
278   // ------------------- det elements: active + readout ----------------------------------
279
280   // 1st Section : below and above the beam pipe
281
282   Double_t lowEdgeActive = -1.*fRMax;
283   Double_t supEdgeActive = lowEdgeActive + fHeightActive;
284   Double_t zMinFront = fZCenter - 0.5*fThicknessSupport - fThicknessActive;
285   Double_t zMinBack  = fZCenter + 0.5*fThicknessSupport;
286   Double_t zMin = 0.;
287   Bool_t isFront = kTRUE;
288   
289   while (lowEdgeActive < 0) {
290     
291     Double_t extLimitAtLowEdgeActive = TMath::Sqrt((fRMax-TMath::Abs(lowEdgeActive)) * TMath::Abs(2*fRMax - (fRMax-TMath::Abs(lowEdgeActive))));
292     Double_t extLimitAtSupEdgeActive = TMath::Sqrt((fRMax-TMath::Abs(supEdgeActive)) * TMath::Abs(2*fRMax - (fRMax-TMath::Abs(supEdgeActive))));
293
294     // creating new det element: active + readout
295     
296     Double_t extLimitDetElem = TMath::Max(extLimitAtLowEdgeActive, extLimitAtSupEdgeActive);
297     
298     if (supEdgeActive<-1.*fRMinSupport+0.01 || lowEdgeActive>1.*fRMinSupport-0.01) {     // single element covering the row
299       
300       nBins[0] = TMath::Nint(2.*extLimitDetElem/fPixelSizeX);
301       nBins[1] = TMath::Nint(fHeightActive/fPixelSizeY);
302       nBins[2] = 1;
303
304       // element below the pipe
305       
306       if (isFront) zMin = zMinFront;
307       else         zMin = zMinBack;
308
309       minPosition[0] = -1.*extLimitDetElem;
310       minPosition[1] = lowEdgeActive;
311       minPosition[2] = zMin;
312       
313       maxPosition[0] = +1.*extLimitDetElem;
314       maxPosition[1] = supEdgeActive;
315       maxPosition[2] = zMin+fThicknessActive; 
316       
317       new ((*fActiveElements)[fActiveElements->GetEntries()]) THnSparseC(Form("MFTActiveElemHist_%02d%03d", fPlaneNumber, fActiveElements->GetEntries()), 
318                                                                          Form("MFTActiveElemHist_%02d%03d", fPlaneNumber, fActiveElements->GetEntries()), 
319                                                                          3, nBins, minPosition, maxPosition);
320
321       minPosition[1] = lowEdgeActive-fHeightReadout;
322       maxPosition[1] = lowEdgeActive;
323       
324       new ((*fReadoutElements)[fReadoutElements->GetEntries()]) THnSparseC(Form("MFTReadoutElemHist_%02d%03d", fPlaneNumber, fReadoutElements->GetEntries()), 
325                                                                            Form("MFTReadoutElemHist_%02d%03d", fPlaneNumber, fReadoutElements->GetEntries()), 
326                                                                            3, nBins, minPosition, maxPosition);
327
328       // specular element above the pipe
329
330       if (fPlaneIsOdd) {
331         if (isFront) zMin = zMinBack;
332         else         zMin = zMinFront;
333       }
334
335       minPosition[0] = -1.*extLimitDetElem;
336       minPosition[1] = -1.*supEdgeActive;
337       minPosition[2] = zMin;
338       
339       maxPosition[0] = +1.*extLimitDetElem;
340       maxPosition[1] = -1.*lowEdgeActive;
341       maxPosition[2] = zMin+fThicknessActive; 
342       
343       new ((*fActiveElements)[fActiveElements->GetEntries()]) THnSparseC(Form("MFTActiveElemHist_%02d%03d", fPlaneNumber, fActiveElements->GetEntries()), 
344                                                                          Form("MFTActiveElemHist_%02d%03d", fPlaneNumber, fActiveElements->GetEntries()), 
345                                                                          3, nBins, minPosition, maxPosition);
346
347       minPosition[1] = -1.*lowEdgeActive;
348       maxPosition[1] = -1.*(lowEdgeActive-fHeightReadout);
349
350       new ((*fReadoutElements)[fReadoutElements->GetEntries()]) THnSparseC(Form("MFTReadoutElemHist_%02d%03d", fPlaneNumber, fReadoutElements->GetEntries()), 
351                                                                            Form("MFTReadoutElemHist_%02d%03d", fPlaneNumber, fReadoutElements->GetEntries()), 
352                                                                            3, nBins, minPosition, maxPosition);
353
354     }
355     
356     else {     // two elements covering the row
357       
358       Double_t intLimitAtLowEdge = 0., intLimitAtSupEdge = 0.;
359       if (fRMinSupport-TMath::Abs(lowEdgeActive)>0.) intLimitAtLowEdge = TMath::Sqrt((fRMinSupport-TMath::Abs(lowEdgeActive)) * TMath::Abs(2*fRMinSupport - (fRMinSupport-TMath::Abs(lowEdgeActive))));
360       if (fRMinSupport-TMath::Abs(supEdgeActive)>0.) intLimitAtSupEdge = TMath::Sqrt((fRMinSupport-TMath::Abs(supEdgeActive)) * TMath::Abs(2*fRMinSupport - (fRMinSupport-TMath::Abs(supEdgeActive))));
361       Double_t intLimitDetElem = TMath::Max(intLimitAtLowEdge, intLimitAtSupEdge);
362       
363       nBins[0] = TMath::Nint((extLimitDetElem-intLimitDetElem)/fPixelSizeX);
364       nBins[1] = TMath::Nint(fHeightActive/fPixelSizeY);
365       nBins[2] = 1;
366       
367       // left element: below the beam pipe
368       
369       if (isFront) zMin = zMinFront;
370       else         zMin = zMinBack;
371
372       minPosition[0] = -1.*extLimitDetElem;
373       minPosition[1] = lowEdgeActive;
374       minPosition[2] = zMin;
375       
376       maxPosition[0] = -1.*intLimitDetElem;
377       maxPosition[1] = supEdgeActive;
378       maxPosition[2] = zMin+fThicknessActive; 
379       
380       new ((*fActiveElements)[fActiveElements->GetEntries()]) THnSparseC(Form("MFTActiveElemHist_%02d%03d", fPlaneNumber, fActiveElements->GetEntries()), 
381                                                                          Form("MFTActiveElemHist_%02d%03d", fPlaneNumber, fActiveElements->GetEntries()), 
382                                                                          3, nBins, minPosition, maxPosition);   
383       
384       minPosition[1] = lowEdgeActive-fHeightReadout;
385       maxPosition[1] = lowEdgeActive;
386       
387       new ((*fReadoutElements)[fReadoutElements->GetEntries()]) THnSparseC(Form("MFTReadoutElemHist_%02d%03d", fPlaneNumber, fReadoutElements->GetEntries()), 
388                                                                            Form("MFTReadoutElemHist_%02d%03d", fPlaneNumber, fReadoutElements->GetEntries()), 
389                                                                            3, nBins, minPosition, maxPosition);
390
391       // left element: above the beam pipe
392       
393       if (supEdgeActive < 0.5*fHeightActive) {
394         
395         if (fPlaneIsOdd) {
396           if (isFront) zMin = zMinBack;
397           else         zMin = zMinFront;
398         }
399         
400         minPosition[0] = -1.*extLimitDetElem;
401         minPosition[1] = -1.*supEdgeActive;
402         minPosition[2] = zMin;
403         
404         maxPosition[0] = -1.*intLimitDetElem;
405         maxPosition[1] = -1.*lowEdgeActive;
406         maxPosition[2] = zMin+fThicknessActive; 
407         
408         new ((*fActiveElements)[fActiveElements->GetEntries()]) THnSparseC(Form("MFTActiveElemHist_%02d%03d", fPlaneNumber, fActiveElements->GetEntries()), 
409                                                                            Form("MFTActiveElemHist_%02d%03d", fPlaneNumber, fActiveElements->GetEntries()), 
410                                                                            3, nBins, minPosition, maxPosition); 
411         
412         minPosition[1] = -1.*lowEdgeActive;
413         maxPosition[1] = -1.*(lowEdgeActive-fHeightReadout);
414         
415         new ((*fReadoutElements)[fReadoutElements->GetEntries()]) THnSparseC(Form("MFTReadoutElemHist_%02d%03d", fPlaneNumber, fReadoutElements->GetEntries()), 
416                                                                              Form("MFTReadoutElemHist_%02d%03d", fPlaneNumber, fReadoutElements->GetEntries()), 
417                                                                              3, nBins, minPosition, maxPosition);
418       
419       }
420
421       // right element: below the beam pipe
422       
423       if (isFront) zMin = zMinFront;
424       else         zMin = zMinBack;
425
426       minPosition[0] = +1.*intLimitDetElem;
427       minPosition[1] = lowEdgeActive;
428       minPosition[2] = zMin;
429       
430       maxPosition[0] = +1.*extLimitDetElem;
431       maxPosition[1] = supEdgeActive;
432       maxPosition[2] = zMin+fThicknessActive; 
433       
434       new ((*fActiveElements)[fActiveElements->GetEntries()]) THnSparseC(Form("MFTActiveElemHist_%02d%03d", fPlaneNumber, fActiveElements->GetEntries()), 
435                                                                          Form("MFTActiveElemHist_%02d%03d", fPlaneNumber, fActiveElements->GetEntries()), 
436                                                                          3, nBins, minPosition, maxPosition);   
437       
438       minPosition[1] = lowEdgeActive-fHeightReadout;
439       maxPosition[1] = lowEdgeActive;
440
441       new ((*fReadoutElements)[fReadoutElements->GetEntries()]) THnSparseC(Form("MFTReadoutElemHist_%02d%03d", fPlaneNumber, fReadoutElements->GetEntries()), 
442                                                                            Form("MFTReadoutElemHist_%02d%03d", fPlaneNumber, fReadoutElements->GetEntries()), 
443                                                                            3, nBins, minPosition, maxPosition);
444
445       // right element: above the beam pipe
446       
447       if (supEdgeActive < 0.5*fHeightActive) {
448
449         if (fPlaneIsOdd) {
450           if (isFront) zMin = zMinBack;
451           else         zMin = zMinFront;
452         }
453         
454         minPosition[0] = +1.*intLimitDetElem;
455         minPosition[1] = -1.*supEdgeActive;
456         minPosition[2] = zMin;
457         
458         maxPosition[0] = +1.*extLimitDetElem;
459         maxPosition[1] = -1.*lowEdgeActive;
460         maxPosition[2] = zMin+fThicknessActive; 
461         
462         new ((*fActiveElements)[fActiveElements->GetEntries()]) THnSparseC(Form("MFTActiveElemHist_%02d%03d", fPlaneNumber, fActiveElements->GetEntries()), 
463                                                                            Form("MFTActiveElemHist_%02d%03d", fPlaneNumber, fActiveElements->GetEntries()), 
464                                                                            3, nBins, minPosition, maxPosition); 
465         
466         minPosition[1] = -1.*lowEdgeActive;
467         maxPosition[1] = -1.*(lowEdgeActive-fHeightReadout);
468         
469         new ((*fReadoutElements)[fReadoutElements->GetEntries()]) THnSparseC(Form("MFTReadoutElemHist_%02d%03d", fPlaneNumber, fReadoutElements->GetEntries()), 
470                                                                              Form("MFTReadoutElemHist_%02d%03d", fPlaneNumber, fReadoutElements->GetEntries()), 
471                                                                              3, nBins, minPosition, maxPosition);
472
473       }
474       
475     }
476     
477     lowEdgeActive += fHeightActive - fActiveSuperposition;
478     supEdgeActive = lowEdgeActive + fHeightActive;
479     isFront = !isFront;
480     
481   }
482   
483   // ------------------- support element -------------------------------------------------
484   
485   nBins[0] = 1;
486   nBins[1] = 1;
487   nBins[2] = 1;
488   
489   minPosition[0] = -1.*fRMaxSupport;
490   minPosition[1] = -1.*fRMaxSupport;
491   minPosition[2] = fZCenter - 0.5*fThicknessSupport;
492   
493   maxPosition[0] = +1.*fRMaxSupport;
494   maxPosition[1] = +1.*fRMaxSupport;
495   maxPosition[2] = fZCenter + 0.5*fThicknessSupport;
496   
497   new ((*fSupportElements)[fSupportElements->GetEntries()]) THnSparseC(Form("MFTSupportElemHist_%02d%03d", fPlaneNumber, fSupportElements->GetEntries()), 
498                                                                        Form("MFTSupportElemHist_%02d%03d", fPlaneNumber, fSupportElements->GetEntries()), 
499                                                                        3, nBins, minPosition, maxPosition);
500
501   // --------------------------------------------------------------------------------------
502
503   AliDebug(1, Form("Structure completed for MFT plane %s", GetName()));
504
505   return kTRUE;
506   
507 }
508
509 //====================================================================================================================================================
510
511 THnSparseC* AliMFTPlane::GetActiveElement(Int_t id) {
512
513   if (id<0 || id>=GetNActiveElements()) return NULL;
514   else return (THnSparseC*) fActiveElements->At(id);
515
516 }
517
518 //====================================================================================================================================================
519
520 THnSparseC* AliMFTPlane::GetReadoutElement(Int_t id) {
521
522   if (id<0 || id>=GetNReadoutElements()) return NULL;
523   else return (THnSparseC*) fReadoutElements->At(id);
524
525 }
526
527 //====================================================================================================================================================
528
529 THnSparseC* AliMFTPlane::GetSupportElement(Int_t id) {
530
531   if (id<0 || id>=GetNSupportElements()) return NULL;
532   else return (THnSparseC*) fSupportElements->At(id);
533
534 }
535
536 //====================================================================================================================================================
537
538 void AliMFTPlane::DrawPlane(Option_t *opt) {
539
540   // ------------------- "FRONT" option ------------------
541
542   if (!strcmp(opt, "front")) {
543
544     TCanvas *cnv = new TCanvas("cnv", GetName(), 900, 900);
545     cnv->Draw();
546
547     TH2D *h = new TH2D("tmp", GetName(), 
548                        1, 1.1*GetSupportElement(0)->GetAxis(0)->GetXmin(), 1.1*GetSupportElement(0)->GetAxis(0)->GetXmax(), 
549                        1, 1.1*GetSupportElement(0)->GetAxis(1)->GetXmin(), 1.1*GetSupportElement(0)->GetAxis(1)->GetXmax());
550     h->SetXTitle("x [cm]");
551     h->SetYTitle("y [cm]");
552     h->Draw();
553
554     AliInfo("Created hist");
555
556     TEllipse *supportExt = new TEllipse(0.0, 0.0, fRMaxSupport, fRMaxSupport);
557     TEllipse *supportInt = new TEllipse(0.0, 0.0, fRMinSupport, fRMinSupport);
558     supportExt->SetFillColor(kCyan-10);
559     supportExt -> Draw("same");
560     supportInt -> Draw("same");
561
562     for (Int_t iEl=0; iEl<GetNActiveElements(); iEl++) {
563       if (!IsFront(GetActiveElement(iEl))) continue;
564       TPave *pave = new TPave(GetActiveElement(iEl)->GetAxis(0)->GetXmin(), 
565                               GetActiveElement(iEl)->GetAxis(1)->GetXmin(), 
566                               GetActiveElement(iEl)->GetAxis(0)->GetXmax(), 
567                               GetActiveElement(iEl)->GetAxis(1)->GetXmax(), 1);
568       pave -> SetFillColor(kGreen);
569       pave -> Draw("same");
570     }
571
572     for (Int_t iEl=0; iEl<GetNReadoutElements(); iEl++) {
573       if (!IsFront(GetReadoutElement(iEl))) continue;
574       TPave *pave = new TPave(GetReadoutElement(iEl)->GetAxis(0)->GetXmin(), 
575                               GetReadoutElement(iEl)->GetAxis(1)->GetXmin(), 
576                               GetReadoutElement(iEl)->GetAxis(0)->GetXmax(), 
577                               GetReadoutElement(iEl)->GetAxis(1)->GetXmax(), 1);
578       pave -> SetFillColor(kRed);
579       pave -> Draw("same");
580     }
581
582   }
583     
584   // ------------------- "BACK" option ------------------
585
586   else if (!strcmp(opt, "back")) {
587
588     TCanvas *cnv = new TCanvas("cnv", GetName(), 900, 900);
589     cnv->Draw();
590     
591     TH2D *h = new TH2D("tmp", GetName(), 
592                        1, 1.1*GetSupportElement(0)->GetAxis(0)->GetXmin(), 1.1*GetSupportElement(0)->GetAxis(0)->GetXmax(), 
593                        1, 1.1*GetSupportElement(0)->GetAxis(1)->GetXmin(), 1.1*GetSupportElement(0)->GetAxis(1)->GetXmax());
594     h->SetXTitle("x [cm]");
595     h->SetYTitle("y [cm]");
596     h->Draw();
597
598     TEllipse *supportExt = new TEllipse(0.0, 0.0, fRMaxSupport, fRMaxSupport);
599     TEllipse *supportInt = new TEllipse(0.0, 0.0, fRMinSupport, fRMinSupport);
600     supportExt -> SetFillColor(kCyan-10);
601     supportExt -> Draw("same");
602     supportInt -> Draw("same");
603
604     for (Int_t iEl=0; iEl<GetNActiveElements(); iEl++) {
605       if (IsFront(GetActiveElement(iEl))) continue;
606       TPave *pave = new TPave(GetActiveElement(iEl)->GetAxis(0)->GetXmin(), 
607                               GetActiveElement(iEl)->GetAxis(1)->GetXmin(), 
608                               GetActiveElement(iEl)->GetAxis(0)->GetXmax(), 
609                               GetActiveElement(iEl)->GetAxis(1)->GetXmax(), 1);
610       pave -> SetFillColor(kGreen);
611       pave -> Draw("same");
612     }
613
614     for (Int_t iEl=0; iEl<GetNReadoutElements(); iEl++) {
615       if (IsFront(GetReadoutElement(iEl))) continue;
616       TPave *pave = new TPave(GetReadoutElement(iEl)->GetAxis(0)->GetXmin(), 
617                               GetReadoutElement(iEl)->GetAxis(1)->GetXmin(), 
618                               GetReadoutElement(iEl)->GetAxis(0)->GetXmax(), 
619                               GetReadoutElement(iEl)->GetAxis(1)->GetXmax(), 1);
620       pave -> SetFillColor(kRed);
621       pave -> Draw("same");
622     }
623
624   }
625
626   // ------------------- "BOTH" option ------------------
627
628   else if (!strcmp(opt, "both")) {
629
630     TCanvas *cnv = new TCanvas("cnv", GetName(), 900, 900);
631     cnv->Draw();
632
633     TH2D *h = new TH2D("tmp", GetName(), 
634                        1, 1.1*GetSupportElement(0)->GetAxis(0)->GetXmin(), 1.1*GetSupportElement(0)->GetAxis(0)->GetXmax(), 
635                        1, 1.1*GetSupportElement(0)->GetAxis(1)->GetXmin(), 1.1*GetSupportElement(0)->GetAxis(1)->GetXmax());
636     h->SetXTitle("x [cm]");
637     h->SetYTitle("y [cm]");
638     h->Draw();
639
640     TEllipse *supportExt = new TEllipse(0.0, 0.0, fRMaxSupport, fRMaxSupport);
641     TEllipse *supportInt = new TEllipse(0.0, 0.0, fRMinSupport, fRMinSupport);
642     supportExt -> SetFillColor(kCyan-10);
643     supportExt -> Draw("same");
644     supportInt -> Draw("same");
645
646     for (Int_t iEl=0; iEl<GetNActiveElements(); iEl++) {
647       if (IsFront(GetActiveElement(iEl)) && GetActiveElement(iEl)->GetAxis(0)->GetXmin()<0.) {
648         TPave *pave = new TPave(GetActiveElement(iEl)->GetAxis(0)->GetXmin(), 
649                                 GetActiveElement(iEl)->GetAxis(1)->GetXmin(), 
650                                 TMath::Min(GetActiveElement(iEl)->GetAxis(0)->GetXmax(), 0.),
651                                 GetActiveElement(iEl)->GetAxis(1)->GetXmax(), 1);
652         pave -> SetFillColor(kGreen);
653         pave -> Draw("same");
654       }
655       else if (!IsFront(GetActiveElement(iEl)) && GetActiveElement(iEl)->GetAxis(0)->GetXmax()>0.) {
656         TPave *pave = new TPave(TMath::Max(GetActiveElement(iEl)->GetAxis(0)->GetXmin(), 0.), 
657                                 GetActiveElement(iEl)->GetAxis(1)->GetXmin(), 
658                                 GetActiveElement(iEl)->GetAxis(0)->GetXmax(), 
659                                 GetActiveElement(iEl)->GetAxis(1)->GetXmax(), 1);
660         pave -> SetFillColor(kGreen);
661         pave -> Draw("same");
662       }
663     }
664     
665     for (Int_t iEl=0; iEl<GetNReadoutElements(); iEl++) {
666       if (IsFront(GetReadoutElement(iEl)) && GetReadoutElement(iEl)->GetAxis(0)->GetXmin()<0.) {
667         TPave *pave = new TPave(GetReadoutElement(iEl)->GetAxis(0)->GetXmin(), 
668                                 GetReadoutElement(iEl)->GetAxis(1)->GetXmin(), 
669                                 TMath::Min(GetReadoutElement(iEl)->GetAxis(0)->GetXmax(), 0.), 
670                                 GetReadoutElement(iEl)->GetAxis(1)->GetXmax(), 1);
671         pave -> SetFillColor(kRed);
672         pave -> Draw("same");
673       }
674       else if (!IsFront(GetReadoutElement(iEl)) && GetReadoutElement(iEl)->GetAxis(0)->GetXmax()>0.) {
675         TPave *pave = new TPave(TMath::Max(GetReadoutElement(iEl)->GetAxis(0)->GetXmin(), 0.),  
676                                 GetReadoutElement(iEl)->GetAxis(1)->GetXmin(), 
677                                 GetReadoutElement(iEl)->GetAxis(0)->GetXmax(), 
678                                 GetReadoutElement(iEl)->GetAxis(1)->GetXmax(), 1);
679         pave -> SetFillColor(kRed);
680         pave -> Draw("same");
681       }
682     }
683     
684   }
685
686   // ------------------- "PROFILE" option ------------------
687
688   else if (!strcmp(opt, "profile")) {
689
690     TCanvas *cnv = new TCanvas("cnv", GetName(), 300, 900);
691     cnv->Draw();
692
693     TH2D *h = new TH2D("tmp", GetName(), 
694                        1, fZCenter-0.5, fZCenter+0.5, 
695                        1, 1.1*GetSupportElement(0)->GetAxis(1)->GetXmin(), 1.1*GetSupportElement(0)->GetAxis(1)->GetXmax());
696     h->SetXTitle("z [cm]");
697     h->SetYTitle("y [cm]");
698     h->Draw();
699
700     TPave *supportExt = new TPave(GetSupportElement(0)->GetAxis(2)->GetXmin(), -fRMaxSupport, 
701                                   GetSupportElement(0)->GetAxis(2)->GetXmax(),  fRMaxSupport);
702     TPave *supportInt = new TPave(GetSupportElement(0)->GetAxis(2)->GetXmin(), -fRMinSupport, 
703                                   GetSupportElement(0)->GetAxis(2)->GetXmax(),  fRMinSupport);
704     supportExt -> SetFillColor(kCyan-10);
705     supportInt -> SetFillColor(kCyan-10);
706     supportExt -> SetBorderSize(1);
707     supportInt -> SetBorderSize(1);
708     supportExt -> Draw("same");
709     supportInt -> Draw("same");
710
711     for (Int_t iEl=0; iEl<GetNActiveElements(); iEl++) {
712       TPave * pave = 0;
713       if (IsFront(GetActiveElement(iEl))) {
714         pave = new TPave(GetActiveElement(iEl)->GetAxis(2)->GetXmax() - 
715                          5*(GetActiveElement(iEl)->GetAxis(2)->GetXmax()-GetActiveElement(iEl)->GetAxis(2)->GetXmin()), 
716                          GetActiveElement(iEl)->GetAxis(1)->GetXmin(), 
717                          GetActiveElement(iEl)->GetAxis(2)->GetXmax(), 
718                          GetActiveElement(iEl)->GetAxis(1)->GetXmax(), 1);
719       }
720       else {
721         pave = new TPave(GetActiveElement(iEl)->GetAxis(2)->GetXmin(), 
722                          GetActiveElement(iEl)->GetAxis(1)->GetXmin(), 
723                          GetActiveElement(iEl)->GetAxis(2)->GetXmin() + 
724                          5*(GetActiveElement(iEl)->GetAxis(2)->GetXmax()-GetActiveElement(iEl)->GetAxis(2)->GetXmin()), 
725                          GetActiveElement(iEl)->GetAxis(1)->GetXmax(), 1);
726       } 
727       pave -> SetFillColor(kGreen);
728       pave -> Draw("same");
729     }
730     
731     for (Int_t iEl=0; iEl<GetNReadoutElements(); iEl++) {
732       TPave *pave = 0;
733       if (IsFront(GetReadoutElement(iEl))) {
734         pave = new TPave(GetReadoutElement(iEl)->GetAxis(2)->GetXmax() - 
735                          5*(GetReadoutElement(iEl)->GetAxis(2)->GetXmax()-GetReadoutElement(iEl)->GetAxis(2)->GetXmin()), 
736                          GetReadoutElement(iEl)->GetAxis(1)->GetXmin(), 
737                          GetReadoutElement(iEl)->GetAxis(2)->GetXmax(), 
738                          GetReadoutElement(iEl)->GetAxis(1)->GetXmax(), 1);
739       }
740       else {
741         pave = new TPave(GetReadoutElement(iEl)->GetAxis(2)->GetXmin(), 
742                          GetReadoutElement(iEl)->GetAxis(1)->GetXmin(), 
743                          GetReadoutElement(iEl)->GetAxis(2)->GetXmin() + 
744                          5*(GetReadoutElement(iEl)->GetAxis(2)->GetXmax()-GetReadoutElement(iEl)->GetAxis(2)->GetXmin()), 
745                          GetReadoutElement(iEl)->GetAxis(1)->GetXmax(), 1);
746       } 
747       pave -> SetFillColor(kRed);
748       pave -> Draw("same");
749     }
750     
751   }
752
753 }
754
755 //====================================================================================================================================================
756
757 Int_t AliMFTPlane::GetNumberOfChips(Option_t *opt) {
758
759   Int_t nChips = 0;
760
761   if (!strcmp(opt, "front")) {
762     for (Int_t iEl=0; iEl<GetNActiveElements(); iEl++) {
763       if (!IsFront(GetActiveElement(iEl))) continue;
764       Double_t length = GetActiveElement(iEl)->GetAxis(0)->GetXmax() - GetActiveElement(iEl)->GetAxis(0)->GetXmin();
765       nChips += Int_t (length/AliMFTConstants::fWidthChip) + 1;
766     }
767   }
768
769   else if (!strcmp(opt, "back")) {
770     for (Int_t iEl=0; iEl<GetNActiveElements(); iEl++) {
771       if (IsFront(GetActiveElement(iEl))) continue;
772       Double_t length = GetActiveElement(iEl)->GetAxis(0)->GetXmax() - GetActiveElement(iEl)->GetAxis(0)->GetXmin();
773       nChips += Int_t (length/AliMFTConstants::fWidthChip) + 1;
774     }
775   }
776
777   return nChips;
778
779 }
780
781 //====================================================================================================================================================