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