]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTriggerChamberEff.cxx
Generation of generic AOD by the test script of MUON
[u/mrichter/AliRoot.git] / MUON / AliMUONTriggerChamberEff.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 /* $Id$ */
17
18 //-----------------------------------------------------------------------------
19 /// \class AliMUONTriggerChamberEff
20 /// Implementation of the trigger chamber efficiency determination from
21 /// data, and returns the
22 /// efficiencyCells.dat with the calculated efficiencies
23 ///
24 /// \author Diego Stocco (Torino)
25 //-----------------------------------------------------------------------------
26
27 #include "AliMUONTriggerChamberEff.h"
28 #include "AliMUONVDigit.h"
29 #include "AliMUONConstants.h"
30 #include "AliMUONTriggerTrack.h"
31 #include "AliMUONDigitMaker.h"
32 #include "AliMUONLocalTrigger.h"
33 #include "AliMUONGeometryTransformer.h"
34
35 #include "AliMUONTrack.h"
36 #include "AliMUONTrackParam.h"
37 #include "AliMUONTrackExtrap.h"
38
39 #include "AliMUONDigitStoreV1.h"
40 #include "AliMUONVDigitStore.h"
41 #include "AliMUONVTriggerStore.h"
42 #include "AliMUONVTriggerTrackStore.h"
43 #include "AliMUONVTrackStore.h"
44
45 #include "AliMpVSegmentation.h"
46 #include "AliMpSegmentation.h"
47 #include "AliMpPad.h"
48 #include "AliMpDEIterator.h"
49 #include "AliMpPlaneType.h"
50 #include "AliMpDEManager.h"
51 #include "AliMpConstants.h"
52
53 #include "AliMpDDLStore.h"
54 #include "AliMpDDL.h"
55 #include "AliMpTriggerCrate.h"
56 #include "AliMpLocalBoard.h"
57
58 #include "AliLog.h"
59
60 #include <Riostream.h>
61 #include <TFile.h>
62 #include <TH1F.h>
63 #include <TMath.h>
64
65 #include <TSeqCollection.h>
66 #include <TTree.h>
67 #include <TROOT.h>
68
69 #include <TStyle.h>
70 #include <TCanvas.h>
71 #include <TH2.h>
72 #include <TSystem.h>
73 #include <TPaveLabel.h>
74
75 #include <cassert>
76
77 /// \cond CLASSIMP
78 ClassImp(AliMUONTriggerChamberEff)
79 /// \endcond
80
81 //_____________________________________________________________________________
82 AliMUONTriggerChamberEff::AliMUONTriggerChamberEff()
83 : TObject(),
84   fTransformer(0x0),
85   fDigitMaker(0x0),
86   fReproduceTrigResponse(kFALSE),
87   fPrintInfo(kFALSE),
88   fWriteOnESD(kFALSE),
89   fDebugLevel(0),
90   fkMaxDistance(99999.)
91 {
92 /// Default constructor
93
94     CheckConstants();
95     ResetArrays();
96 }
97
98
99 //_____________________________________________________________________________
100 AliMUONTriggerChamberEff::AliMUONTriggerChamberEff(const AliMUONGeometryTransformer* transformer,
101                                                    const AliMUONDigitMaker* digitMaker,
102                                                    Bool_t writeOnESD)
103 : TObject(),
104   fTransformer(transformer),
105   fDigitMaker(digitMaker),
106   fReproduceTrigResponse(kFALSE),
107   fPrintInfo(kFALSE),
108   fWriteOnESD(writeOnESD),
109   fDebugLevel(0),
110   fkMaxDistance(99999.)
111 {
112 /// Standard constructor
113
114     CheckConstants();
115     ResetArrays();
116 }
117
118
119 //_____________________________________________________________________________
120 AliMUONTriggerChamberEff::~AliMUONTriggerChamberEff()
121 {
122 /// Destructor
123     Bool_t writeOnESD=fWriteOnESD;
124     fWriteOnESD=kFALSE;
125     if(writeOnESD) SaveInESDFile();
126 }
127
128
129 //_____________________________________________________________________________
130 AliMUONTriggerChamberEff::AliMUONTriggerChamberEff(const AliMUONTriggerChamberEff& other)
131     :TObject(other),
132      fTransformer(0x0),
133      fDigitMaker(0x0),
134      fReproduceTrigResponse(other.fReproduceTrigResponse),
135      fPrintInfo(other.fPrintInfo),
136      fWriteOnESD(other.fWriteOnESD),
137      fDebugLevel(other.fDebugLevel),
138      fkMaxDistance(other.fkMaxDistance),
139      fTrigger44(other.fTrigger44),
140      fTrigger34(other.fTrigger34)
141 {
142 /// Copy constructor
143
144     for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
145         fInefficientSlat[chCath] = other.fInefficientSlat[chCath];
146         fHitPerSlat[chCath] = other.fHitPerSlat[chCath];
147         fInefficientBoard[chCath] = other.fInefficientBoard[chCath];
148         fHitPerBoard[chCath] = other.fHitPerBoard[chCath];
149     }
150 }
151
152
153 //_____________________________________________________________________________
154 AliMUONTriggerChamberEff& AliMUONTriggerChamberEff::operator=(const AliMUONTriggerChamberEff& other)
155 {
156     /// Asignment operator
157     // check assignement to self
158     if (this == &other)
159         return *this;
160
161     // base class assignement
162     TObject::operator=(other);
163
164     fTransformer = 0x0;
165     fDigitMaker = 0x0;
166     fReproduceTrigResponse = other.fReproduceTrigResponse;
167     fPrintInfo = other.fPrintInfo;
168     fWriteOnESD = other.fWriteOnESD;
169     fDebugLevel = other.fDebugLevel;
170
171     fTrigger44 = other.fTrigger44;
172     fTrigger34 = other.fTrigger34;
173
174     for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
175         fInefficientSlat[chCath] = other.fInefficientSlat[chCath];
176         fHitPerSlat[chCath] = other.fHitPerSlat[chCath];
177         fInefficientBoard[chCath] = other.fInefficientBoard[chCath];
178         fHitPerBoard[chCath] = other.fHitPerBoard[chCath];
179     }
180     return *this;
181 }
182
183
184 //_____________________________________________________________________________
185 void AliMUONTriggerChamberEff::ResetArrays()
186 {
187     //
188     /// Sets the data member counters to 0.
189     //
190
191     fTrigger44.Set(fgkNcathodes);
192     fTrigger34.Set(fgkNplanes);
193
194     for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
195         fInefficientSlat[chCath].Set(fgkNslats);
196         fHitPerSlat[chCath].Set(fgkNslats);
197         fInefficientBoard[chCath].Set(AliMpConstants::NofLocalBoards());
198         fHitPerBoard[chCath].Set(AliMpConstants::NofLocalBoards());
199     }
200
201     fTrigger44.Reset();
202     fTrigger34.Reset();
203
204     for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
205         fInefficientSlat[chCath].Reset();
206         fHitPerSlat[chCath].Reset();
207         fInefficientBoard[chCath].Reset();
208         fHitPerBoard[chCath].Reset();
209     }
210 }
211
212
213 //______________________________________________________________________________
214 Bool_t 
215 AliMUONTriggerChamberEff::TriggerDigits(const AliMUONVTriggerStore& triggerStore,
216                                         AliMUONVDigitStore& digitStore) const
217 {
218     //
219     /// make (S)Digit for trigger
220     //
221     digitStore.Clear();
222
223     AliMUONLocalTrigger* locTrg;
224     TIter next(triggerStore.CreateLocalIterator());
225   
226     while ( ( locTrg = static_cast<AliMUONLocalTrigger*>(next()) ) ) 
227     {
228         if (locTrg->IsNull()) continue;
229    
230         TArrayS xyPattern[fgkNcathodes];
231         locTrg->GetXPattern(xyPattern[0]);
232         locTrg->GetYPattern(xyPattern[1]);
233     
234         Int_t nBoard = locTrg->LoCircuit();
235         fDigitMaker->TriggerDigits(nBoard, xyPattern, digitStore);
236     }
237     return kTRUE;
238 }
239
240
241 //_____________________________________________________________________________
242 void AliMUONTriggerChamberEff::InfoDigit(AliMUONVDigitStore& digitStore)
243 {
244     //
245     /// Prints information on digits (for debugging)
246     //
247     TIter next(digitStore.CreateIterator());
248     AliMUONVDigit* mDigit=0x0;
249     
250     while ( ( mDigit = static_cast<AliMUONVDigit*>(next()) ) )
251     {
252         mDigit->Print();
253     } // end digit loop
254     printf("\n");
255 }
256
257
258 //_____________________________________________________________________________
259 Int_t AliMUONTriggerChamberEff::MatchingPad(AliMUONVDigitStore& digitStore, Int_t &detElemId,
260                                             Float_t coor[2], Bool_t isMatch[2],
261                                             TArrayI nboard[2], TArrayF &zRealMatch, Float_t y11)
262 {
263     //
264     /// Check slat and board number of digit matching track
265     //
266
267     enum {kBending, kNonBending};
268
269     Float_t minMatchDist[fgkNcathodes];
270     Int_t padsInCheckArea[fgkNcathodes];
271
272     for(Int_t cath=0; cath<fgkNcathodes; cath++){
273         isMatch[cath] = kFALSE;
274         minMatchDist[cath] = fkMaxDistance/10.;
275         padsInCheckArea[cath] = 0;
276     }
277     Int_t iChamber = AliMpDEManager::GetChamberId(detElemId);
278     Int_t ch = iChamber-10;
279     Float_t oldDeltaZ = AliMUONConstants::DefaultChamberZ(iChamber) - AliMUONConstants::DefaultChamberZ(10);
280     Float_t y = coor[1];
281     Int_t iSlat = detElemId%100;
282     Int_t trigDigitBendPlane = -1;
283     Int_t foundDetElemId = detElemId;
284     Float_t foundZmatch=999.;
285     Float_t yCoorAtPadZ=999.;
286
287     TIter next(digitStore.CreateIterator());
288     AliMUONVDigit* mDigit;
289     Int_t idigit=0;
290     
291     while ( ( mDigit = static_cast<AliMUONVDigit*>(next()) ) )
292     {
293         idigit++;
294         Int_t currDetElemId = mDigit->DetElemId();
295         Int_t currCh = AliMpDEManager::GetChamberId(currDetElemId);
296         if(currCh!=iChamber) continue;
297         Int_t currSlat = currDetElemId%100;
298         Int_t slatDiff = TMath::Abs(currSlat-iSlat);
299         if(slatDiff>1 && slatDiff<17) continue; // Check neighbour slats
300         Int_t cathode = mDigit->Cathode();
301         Int_t ix = mDigit->PadX();
302         Int_t iy = mDigit->PadY();
303         Float_t xpad, ypad, zpad;
304         const AliMpVSegmentation* seg = AliMpSegmentation::Instance()
305             ->GetMpSegmentation(currDetElemId,AliMp::GetCathodType(cathode));
306
307         AliMpPad pad = seg->PadByIndices(AliMpIntPair(ix,iy),kTRUE);
308         Float_t xlocal1 = pad.Position().X();
309         Float_t ylocal1 = pad.Position().Y();
310         Float_t dpx = pad.Dimensions().X();
311         Float_t dpy = pad.Dimensions().Y();
312         fTransformer->Local2Global(currDetElemId, xlocal1, ylocal1, 0, xpad, ypad, zpad);
313         if(fDebugLevel>2)printf("DetElemId = %i\tCathode = %i\t(x,y) Pad = (%i,%i) = (%.2f,%.2f)\tDim = (%.2f,%.2f)\tTrack = (%.2f,%.2f)\n",currDetElemId,cathode,ix,iy,xpad,ypad,dpx,dpy,coor[0],coor[1]);
314         // searching track intersection with chambers (second approximation)
315         if(ch%2==1){
316             //if(iChamber%2==1){
317             Float_t deltaZ = zpad - zRealMatch[0];
318             y = (coor[1]-y11)*deltaZ/oldDeltaZ + y11;
319             if(fDebugLevel>=3 && TMath::Abs(y-coor[1])>0.1)printf("oldDeltaZ = %7.2f   newDeltaZ = %7.2f\toldY = %7.2f   new y = %7.2f\n",oldDeltaZ,deltaZ,coor[1],y);
320         }
321         Float_t matchDist = PadMatchTrack(xpad, ypad, dpx, dpy, coor[0], y, ch);
322         if(matchDist<fkMaxDistance/2.) padsInCheckArea[cathode]++;
323         if(matchDist>minMatchDist[cathode])continue;
324         isMatch[cathode] = kTRUE;
325         minMatchDist[cathode] = matchDist;
326         foundDetElemId = currDetElemId;
327         foundZmatch=zpad;
328         yCoorAtPadZ=y;
329         if(cathode==kBending) trigDigitBendPlane = idigit;
330         for (Int_t loc=0; loc<pad.GetNofLocations(); loc++){
331             AliMpIntPair location = pad.GetLocation(loc);
332             nboard[cathode][loc] = location.GetFirst();
333         }
334         for(Int_t loc=pad.GetNofLocations(); loc<fgkNlocations; loc++){
335             nboard[cathode][loc]=-1;
336         }
337     }
338
339     for(Int_t cath=0; cath<fgkNcathodes; cath++){
340         if(padsInCheckArea[cath]>2) {
341             if(fDebugLevel>=1) printf("padsInCheckArea[%i] = %i\n",cath,padsInCheckArea[cath]);
342             return -500;
343         }
344     }
345
346     if(isMatch[kBending] || isMatch[kNonBending]){
347         detElemId = foundDetElemId;
348         zRealMatch[ch] = foundZmatch;
349         coor[1] = yCoorAtPadZ;
350         if(fDebugLevel>2){
351             Int_t whichCathode=kBending;
352             if(!isMatch[kBending])whichCathode=kNonBending;
353         }
354     }
355     return trigDigitBendPlane;
356 }
357
358 //_____________________________________________________________________________
359 Float_t AliMUONTriggerChamberEff::PadMatchTrack(Float_t xPad, Float_t yPad,
360                                                 Float_t dpx, Float_t dpy, 
361                                                 Float_t xTrackAtPad, Float_t yTrackAtPad,
362                                                 Int_t chamber)
363 {
364     //
365     /// Decides if the digit belongs to the trigger track.
366     //
367
368     Float_t maxDist = 2.;//3. // cm
369     Float_t maxDistCheckArea = 6.; // cm
370
371     Float_t matchDist = fkMaxDistance;
372
373     Float_t deltaX = TMath::Abs(xPad-xTrackAtPad)-dpx;
374     Float_t deltaY = TMath::Abs(yPad-yTrackAtPad)-dpy;
375     Float_t maxDistX = maxDist;
376     Float_t maxDistY = maxDist;
377     
378     if(fReproduceTrigResponse){
379         maxDistX = dpx;
380         maxDistY = dpy;
381         deltaX = TMath::Abs(xPad-xTrackAtPad);
382         deltaY = TMath::Abs(yPad-yTrackAtPad);
383         if(dpx<dpy && chamber>=2) maxDistX = 3.*dpx;// Non-bending plane: check the +- 1 strip between stations
384         if(dpy<dpx && chamber%2) maxDistY = 3.*dpy;// bending plane: check the +- 1 strip between planes in the same station
385     }
386
387     if(deltaX<=maxDistX && deltaY<=maxDistY) matchDist = TMath::Max(deltaX, deltaY);
388     else if(deltaX<=maxDistCheckArea && deltaY<=maxDistCheckArea) matchDist = fkMaxDistance/5.;
389     return matchDist;
390 }
391
392
393 //_____________________________________________________________________________
394 void AliMUONTriggerChamberEff::CalculateEfficiency(Int_t trigger44, Int_t trigger34,
395                                                    Float_t &efficiency, Float_t &error,
396                                                    Bool_t failuresAsInput)
397 {
398     //
399     /// Returns the efficiency.
400     //
401
402     efficiency=-9.;
403     error=0.;
404     if(trigger34>0){
405         efficiency=(Double_t)trigger44/((Double_t)trigger34);
406         if(failuresAsInput)efficiency=1.-(Double_t)trigger44/((Double_t)trigger34);
407     }
408     Double_t q = TMath::Abs(1-efficiency);
409     if(efficiency<0)error=0.0;
410     else error = TMath::Sqrt(efficiency*q/((Double_t)trigger34));
411 }
412
413
414 //_____________________________________________________________________________
415 Int_t AliMUONTriggerChamberEff::DetElemIdFromPos(Float_t x, Float_t y, 
416                                                  Int_t chamber, Int_t cathode)
417 {
418     //
419     /// Given the (x,y) position in the chamber,
420     /// it returns the corresponding slat
421     //
422
423     Int_t resultingDetElemId = -1;
424     AliMpDEIterator it;
425     Float_t minDist = 999.;
426     for ( it.First(chamber-1); ! it.IsDone(); it.Next() ){
427         Int_t detElemId = it.CurrentDEId();
428         Int_t ich = detElemId/100-10;
429         Float_t tolerance=0.2*((Float_t)ich);
430         Float_t currDist=9999.;
431
432         const AliMpVSegmentation* seg = 
433             AliMpSegmentation::Instance()
434             ->GetMpSegmentation(detElemId,AliMp::GetCathodType(cathode));
435         if (!seg) continue;
436
437         Float_t deltax = seg->Dimensions().X();
438         Float_t deltay = seg->Dimensions().Y();
439         Float_t xlocal1 =  -deltax;
440         Float_t ylocal1 =  -deltay;
441         Float_t xlocal2 =  +deltax;
442         Float_t ylocal2 =  +deltay;
443         Float_t xg01, yg01, zg1, xg02, yg02, zg2;
444         fTransformer->Local2Global(detElemId, xlocal1, ylocal1, 0, xg01, yg01, zg1);
445         fTransformer->Local2Global(detElemId, xlocal2, ylocal2, 0, xg02, yg02, zg2);
446
447         Float_t xg1 = xg01, xg2 = xg02, yg1 = yg01, yg2 = yg02;
448
449         if(xg01>xg02){
450             xg1 = xg02;
451             xg2 = xg01;
452         }
453         if(yg01>yg02){
454             yg1 = yg02;
455             yg2 = yg01;
456         }
457
458         if(x>=xg1-tolerance && x<=xg2+tolerance && y>=yg1-tolerance && y<=yg2+tolerance){ // takes into account errors in extrapolation
459             if(y<yg1) currDist = yg1-y;
460             else if(y>yg2) currDist = y-yg2;
461             if(currDist<minDist) {
462                 resultingDetElemId = detElemId;
463                 minDist=currDist;
464                 continue;
465             }
466             resultingDetElemId = detElemId;
467             break;
468         }
469     } // loop on detElemId
470     return resultingDetElemId;
471 }
472
473
474 //_____________________________________________________________________________
475 void AliMUONTriggerChamberEff::LocalBoardFromPos(Float_t x, Float_t y,
476                                                  Int_t detElemId, Int_t cathode,
477                                                  Int_t localBoard[4])
478 {
479     //
480     /// Given the (x,y) position in the chamber,
481     /// it returns the corresponding local board
482     //
483
484     for(Int_t loc=0; loc<fgkNlocations; loc++){
485         localBoard[loc]=-1;
486     }
487     Float_t xl, yl, zl;
488     fTransformer->Global2Local(detElemId, x, y, 0, xl, yl, zl);
489     TVector2 pos(xl,yl);
490     const AliMpVSegmentation* seg = 
491         AliMpSegmentation::Instance()
492           ->GetMpSegmentation(detElemId,AliMp::GetCathodType(cathode));
493     if (seg){
494         AliMpPad pad = seg->PadByPosition(pos,kFALSE);
495         for (Int_t loc=0; loc<pad.GetNofLocations(); loc++){
496             AliMpIntPair location = pad.GetLocation(loc);
497             localBoard[loc] = location.GetFirst();
498         }
499     }
500 }
501
502
503 //_____________________________________________________________________________
504 void AliMUONTriggerChamberEff::EventChamberEff(const AliMUONVTriggerStore& triggerStore,
505                                                const AliMUONVTriggerTrackStore& trigTrackStore,
506                                                const AliMUONVTrackStore& trackStore)
507 {
508     //
509     /// Main method.
510     /// It loops over the the trigger rec. tracks in the event.
511     /// Then it search for matching digits around the track.
512     /// Finally it calculates the efficiency for each trigger board.
513     /// Files with calculated efficiency are placed in the user defined outputDir.
514     //
515
516     if(!fTransformer || ! fDigitMaker) {
517         AliError(Form("AliMUONGeometryTransformer or AliMUONDigitMaker not properly initialized!!"));
518         return;
519     }
520
521     enum {kBending, kNonBending};
522     Float_t rad2deg = 180./TMath::Pi();
523
524     Int_t chOrder[fgkNchambers] = {0,2,1,3};
525
526     TArrayF zRealMatch(fgkNchambers);
527     TArrayF correctFactor(fgkNcathodes);
528
529     Bool_t match[fgkNchambers][fgkNcathodes];
530     Bool_t matchPad[fgkNcathodes];
531     for(Int_t cath=0; cath<fgkNcathodes; cath++){
532         matchPad[cath] = kFALSE;
533     }
534
535     TArrayF zMeanChamber(fgkNchambers);
536     for(Int_t ch=0; ch<fgkNchambers; ch++){
537         zMeanChamber[ch] = AliMUONConstants::DefaultChamberZ(10+ch);
538     }
539
540     TArrayI digitPerTrack(fgkNcathodes);
541
542     Float_t trackIntersectCh[fgkNchambers][fgkNcathodes];
543
544     TArrayI triggeredDigits[2];
545     for(Int_t itrack=0; itrack<2; itrack++){
546         triggeredDigits[itrack].Set(fgkNchambers);
547         triggeredDigits[itrack].Reset(-1);
548     }
549
550     TArrayI trigScheme[fgkNcathodes];
551     TArrayI slatThatTriggered[fgkNcathodes];
552     for(Int_t cath=0; cath<fgkNcathodes; cath++){
553         trigScheme[cath].Set(fgkNchambers);
554         slatThatTriggered[cath].Set(fgkNchambers);
555     }
556
557     Int_t boardThatTriggered[fgkNchambers][fgkNcathodes][fgkNlocations];
558     TArrayI nboard[fgkNcathodes];
559     for(Int_t cath=0; cath<fgkNcathodes; cath++){
560         nboard[cath].Set(fgkNlocations);
561     }
562     Int_t ineffBoard[fgkNlocations];
563     for(Int_t loc=0; loc<fgkNlocations; loc++){
564         ineffBoard[loc] = -1;
565     }
566
567     AliMUONDigitStoreV1 digitStore;   
568     TriggerDigits(triggerStore,digitStore);
569
570     AliMUONTriggerTrack *recTrigTrack = 0x0;
571
572     TIter next(trigTrackStore.CreateIterator());
573     
574     while ( ( recTrigTrack = static_cast<AliMUONTriggerTrack*>(next()) ) )
575     {
576         if(!IsCleanTrack(recTrigTrack, trackStore)) {
577             if(fDebugLevel>=1) printf("\tTrack %p (%f, %f) don't match tracker track: rejected!\n",(void *)recTrigTrack,recTrigTrack->GetX11(),recTrigTrack->GetY11());
578             continue;
579         }
580
581         digitPerTrack.Reset();
582         for(Int_t ch=0; ch<fgkNchambers; ch++){
583             zRealMatch[ch] = zMeanChamber[ch];
584             for(Int_t cath=0; cath<fgkNcathodes; cath++){
585                 match[ch][cath]=kFALSE;
586                 //slatThatTriggered[ch][cath]=-1;
587                 //trigScheme[ch][cath] = 0;
588                 for(Int_t loc=0; loc<fgkNlocations; loc++){
589                     boardThatTriggered[ch][cath][loc]=-1;
590                 }
591             }
592         }
593
594         for(Int_t cath=0; cath<fgkNcathodes; cath++){
595             slatThatTriggered[cath].Reset(-1);
596             trigScheme[cath].Reset();
597         }
598
599         Bool_t isClearEvent = kTRUE;
600         Bool_t doubleCountTrack = kFALSE;
601
602         Float_t x11 = recTrigTrack->GetX11();// x position (info from non-bending plane)
603         Float_t y11 = recTrigTrack->GetY11();// y position (info from bending plane)
604         Float_t thetaX = recTrigTrack->GetThetax();
605         Float_t thetaY = recTrigTrack->GetThetay();
606
607         if(fDebugLevel>=3) printf("\tTrack = %p\npos from track: (x,y) = (%f, %f), (thetaX, thetaY) = (%f, %f)\n",(void *)recTrigTrack,x11,y11,thetaX*rad2deg,thetaY*rad2deg);
608
609         for(Int_t ch=0; ch<fgkNchambers; ch++) { // chamber loop
610             Int_t currCh = chOrder[ch];
611             if(fDebugLevel>=2)
612                 printf("zMeanChamber[%i] = %.2f\tzRealMatch[0] = %.2f\n",currCh,zMeanChamber[currCh],zRealMatch[0]);
613
614             for(Int_t cath=0; cath<fgkNcathodes; cath++){
615                 correctFactor[cath]=1.;
616             }
617             // calculate corrections to trigger track theta
618             if(ch>=1)correctFactor[kNonBending] = zMeanChamber[0]/zRealMatch[0];// corrects x position
619             if(ch>=2)correctFactor[kBending] = (zMeanChamber[2] - zMeanChamber[0]) / (zRealMatch[2] - zRealMatch[0]);// corrects y position
620
621             // searching track intersection with chambers (first approximation)
622             Float_t deltaZ = zMeanChamber[currCh] - zMeanChamber[0];
623             trackIntersectCh[currCh][0] = zMeanChamber[currCh] * TMath::Tan(thetaX) * correctFactor[kNonBending];// x position (info from non-bending plane) 
624             trackIntersectCh[currCh][1] = y11 + deltaZ * TMath::Tan(thetaY) * correctFactor[kBending];// y position (info from bending plane)
625             Int_t detElemIdFromTrack = DetElemIdFromPos(trackIntersectCh[currCh][0], trackIntersectCh[currCh][1], 11+currCh, 0);
626             if(detElemIdFromTrack<0) {
627                 if(fDebugLevel>1) printf("Warning: trigger track outside trigger chamber\n");
628                 continue;
629             }
630                 
631             triggeredDigits[1][currCh] = MatchingPad(digitStore, detElemIdFromTrack, trackIntersectCh[currCh], matchPad, nboard, zRealMatch, y11);
632
633             // if MatchingPad = -500 => too many digits matching pad =>
634             //                       => Event not clear => Reject track
635             if(triggeredDigits[1][currCh]<-100){
636                 isClearEvent = kFALSE;
637                 if(fDebugLevel>=1) printf("Warning: track = %p (%i) matches many pads. Rejected!\n",(void *)recTrigTrack, detElemIdFromTrack);
638                 break;
639             }
640
641             // deciding if digit matches track
642             Bool_t isDiffLocBoard = kFALSE;
643             if(fReproduceTrigResponse && ch>2){
644                 for(Int_t cath=0; cath<fgkNcathodes; cath++){
645                     if(boardThatTriggered[currCh][cath][0]>=0){
646                         if(boardThatTriggered[currCh][cath][0]!=boardThatTriggered[currCh-1][cath][0]) isDiffLocBoard = kTRUE;
647                     }
648                 }
649             }
650
651             if(isDiffLocBoard && fDebugLevel>=1)printf("\tDifferent local board\n");
652
653             for(Int_t cath=0; cath<fgkNcathodes; cath++){
654                 match[currCh][cath] = (matchPad[cath] && !isDiffLocBoard);
655                 if(!match[currCh][cath]) continue;
656                 digitPerTrack[cath]++;
657                 trigScheme[cath][currCh]++;
658                 slatThatTriggered[cath][currCh] = detElemIdFromTrack;
659                 for(Int_t loc=0; loc<fgkNlocations; loc++){
660                     boardThatTriggered[currCh][cath][loc] = nboard[cath][loc];
661                 }
662             }
663         } // end chamber loop
664
665         for(Int_t cath=0; cath<fgkNcathodes; cath++){
666             if(digitPerTrack[cath]<3)isClearEvent = kFALSE;
667             if(fDebugLevel>=1 && !isClearEvent)printf("Warning: found %i digits for trigger track cathode %i.\nRejecting event\n", digitPerTrack[cath],cath);
668         }
669
670         if(!isClearEvent && !fReproduceTrigResponse) continue;
671
672         Int_t commonDigits = 0;
673         for(Int_t ch=0; ch<fgkNchambers; ch++){
674             if(triggeredDigits[1][ch]==triggeredDigits[0][ch]) commonDigits++; // Compare with previous track
675             triggeredDigits[0][ch] = triggeredDigits[1][ch]; // Store this track parameters for comparison with next one
676         }
677         if(commonDigits>=2){
678             doubleCountTrack=kTRUE;
679         }
680
681         if(!doubleCountTrack || fReproduceTrigResponse){
682             for(Int_t cath=0; cath<fgkNcathodes; cath++){
683                 Int_t is44 = 1;
684                 Bool_t goodForSlatEff = kTRUE;
685                 Bool_t goodForBoardEff = kTRUE;
686                 Int_t ineffSlat = -1;
687                 Int_t ineffDetElId = -1;
688                 Int_t firstSlat = slatThatTriggered[cath][0]%100;
689                 if(firstSlat<0) firstSlat=slatThatTriggered[cath][1]%100;
690                 Int_t firstBoard = boardThatTriggered[0][kBending][0];
691                 if(firstBoard<0) firstBoard=boardThatTriggered[1][kBending][0];
692                 for(Int_t ch=0; ch<fgkNchambers; ch++){
693                     Bool_t isCurrChIneff = kFALSE;
694                     is44 *= trigScheme[cath][ch];
695                     Int_t currSlat = slatThatTriggered[cath][ch]%100;
696                     if(currSlat<0){
697                         ineffDetElId = DetElemIdFromPos(trackIntersectCh[ch][0], trackIntersectCh[ch][1], 11+ch, cath);
698                         currSlat = ineffDetElId%100;
699                         ineffSlat = currSlat;
700                         isCurrChIneff = kTRUE;
701                     }
702                     if(currSlat!=firstSlat)goodForSlatEff=kFALSE;
703                     Bool_t atLeastOneLoc=kFALSE;
704                     if(isCurrChIneff) LocalBoardFromPos(trackIntersectCh[ch][0], trackIntersectCh[ch][1], ineffDetElId, cath, ineffBoard);
705                     for(Int_t loc=0; loc<fgkNlocations; loc++){
706                         Int_t currBoard = boardThatTriggered[ch][cath][loc];
707                         if(isCurrChIneff) currBoard = ineffBoard[loc];
708                         if(currBoard==firstBoard){
709                             atLeastOneLoc=kTRUE;
710                             break;
711                         }
712                     }
713                     if(!atLeastOneLoc)goodForBoardEff=kFALSE;
714                 } // end chamber loop
715
716                 // Trigger 4/4
717                 if(is44==1){
718                     fTrigger44[cath]++;
719                     if(fDebugLevel>=1)printf("Trigger44[%i] = %i\n",cath,fTrigger44[cath]);
720                     if(goodForSlatEff){
721                         for(Int_t ch=0; ch<fgkNchambers; ch++){
722                             Int_t chCath = fgkNchambers*cath + ch;
723                             fHitPerSlat[chCath][firstSlat]++;
724                             if(fDebugLevel>=1)printf("Slat that triggered = %i\n",slatThatTriggered[cath][ch]);
725                             if(goodForBoardEff && firstBoard>0){
726                                 fHitPerBoard[chCath][firstBoard-1]++;
727                                 if(fDebugLevel>=1)printf("Board that triggered = %i\n",firstBoard);
728                             }
729                             else if(fDebugLevel>=1) printf("Track = %p: Particle crossed different boards: rejected!\n",(void *)recTrigTrack);
730                         }
731                     }
732                     else if(fDebugLevel>=1) printf("Track = %p: Particle crossed different slats: rejected!\n",(void *)recTrigTrack);
733                 }
734
735                 // Trigger 3/4
736                 if(ineffDetElId>0){
737                     Int_t ineffCh = ineffDetElId/100-11;
738                     Int_t chCath = fgkNchambers*cath + ineffCh;
739                     fTrigger34[chCath]++;
740                     if(fDebugLevel>=1) printf("Trigger34[%i] = %i\n",chCath,fTrigger34[chCath]);
741                     if(goodForSlatEff){
742                         if(fDebugLevel>=1) printf("Slat non efficient = %i\n",ineffDetElId);
743                         fInefficientSlat[chCath][ineffSlat]++;
744
745                         if(goodForBoardEff && firstBoard>0){
746                             if(fDebugLevel>=1) printf("Board non efficient = %i\n",firstBoard);
747                             fInefficientBoard[chCath][firstBoard-1]++;
748                         }
749                         else if(fDebugLevel>=1) printf("Track = %p: Particle crossed different boards: rejected!\n",(void *)recTrigTrack);
750                     }
751                     else if(fDebugLevel>=1) printf("Track = %p: Particle crossed different slats: rejected!\n",(void *)recTrigTrack);
752                 }
753             } // end loop on cathodes
754         }
755         else if(doubleCountTrack){
756             if(fDebugLevel>=1)
757                 printf("\n\tTrack = %p: \nDouble Count Track: Track rejected!\n",(void *)recTrigTrack);
758         }
759     } // end trigger tracks loop
760
761     if(fPrintInfo) InfoDigit(digitStore);
762 }
763
764 //_____________________________________________________________________________
765 void AliMUONTriggerChamberEff::WriteEfficiencyMap(const char* outputDir)
766 {
767     //
768     /// Writes information on calculated efficiency.
769     /// It writes: triggerChamberEff.root file containing efficiency histograms.
770     //
771
772     char *cathCode[fgkNcathodes] = {"bendPlane", "nonBendPlane"};
773
774     char outFileName[100];
775
776     sprintf(outFileName, "%s/MUON.TriggerEfficiencyMap.root",outputDir);
777     TFile *outputHistoFile = new TFile(outFileName,"RECREATE");
778     TDirectory *dir = gDirectory;
779
780     char *yAxisTitle = "trigger efficiency (a.u.)";
781     char *xAxisTitle = "chamber";
782
783     const Int_t kNumOfBoards = AliMpConstants::NofLocalBoards();
784
785     TH1F *histoChamber[fgkNcathodes];
786     TH1F *histoSlat[8];
787     TH1F *histoBoard[8];
788
789     // ADDED for check
790     enum {kAllChEff, kChNonEff, kNumOfHistoTypes};
791     char *histoTypeName[kNumOfHistoTypes] = {"CountInCh", "NonCountInCh"};
792     char *histoTypeTitle[kNumOfHistoTypes] = {"counted", "non counted"};
793     TH1F *histoCheckSlat[8][kNumOfHistoTypes];
794     TH1F *histoCheckBoard[8][kNumOfHistoTypes];
795     // end ADDED for check
796
797     char histoName[40];
798     char histoTitle[90];
799
800     for(Int_t cath=0; cath<fgkNcathodes; cath++){
801         sprintf(histoName, "%sChamberEff", cathCode[cath]);
802         sprintf(histoTitle, "Chamber efficiency %s", cathCode[cath]);
803         histoChamber[cath] = new TH1F(histoName, histoTitle, fgkNchambers, 11-0.5, 15-0.5);
804         histoChamber[cath]->SetXTitle(xAxisTitle);
805         histoChamber[cath]->SetYTitle(yAxisTitle);
806         histoChamber[cath]->GetXaxis()->SetNdivisions(fgkNchambers);
807         for(Int_t ch=0; ch<fgkNchambers; ch++){
808             Int_t chCath = fgkNchambers*cath + ch;
809             sprintf(histoName, "%sSlatEffChamber%i", cathCode[cath], 11+ch);
810             sprintf(histoTitle, "Chamber %i: slat efficiency %s", 11+ch, cathCode[cath]);
811             histoSlat[chCath] = new TH1F(histoName, histoTitle, fgkNslats, 0-0.5, fgkNslats-0.5);
812             histoSlat[chCath]->SetXTitle("slat");
813             histoSlat[chCath]->SetYTitle(yAxisTitle);
814             histoSlat[chCath]->GetXaxis()->SetNdivisions(fgkNslats);
815                 
816             sprintf(histoName, "%sBoardEffChamber%i", cathCode[cath], 11+ch);
817             sprintf(histoTitle, "Chamber %i: board efficiency %s", 11+ch, cathCode[cath]);
818             histoBoard[chCath] = new TH1F(histoName, histoTitle, kNumOfBoards, 1-0.5, kNumOfBoards+1.-0.5);
819             histoBoard[chCath]->SetXTitle("boards");
820             histoBoard[chCath]->SetYTitle(yAxisTitle);
821             histoBoard[chCath]->GetXaxis()->SetNdivisions(kNumOfBoards);
822
823             // ADDED for check
824             for(Int_t hType=0; hType<kNumOfHistoTypes; hType++){
825                 sprintf(histoName, "%sSlat%s%i", cathCode[cath], histoTypeName[hType], 11+ch);
826                 sprintf(histoTitle, "Chamber %i: slat %s %s", 11+ch, histoTypeTitle[hType], cathCode[cath]);
827                 histoCheckSlat[chCath][hType] = new TH1F(histoName, histoTitle, fgkNslats, 0-0.5, fgkNslats-0.5);
828                 histoCheckSlat[chCath][hType]->SetXTitle("slat");
829                 histoCheckSlat[chCath][hType]->SetYTitle(yAxisTitle);
830                 histoCheckSlat[chCath][hType]->GetXaxis()->SetNdivisions(fgkNslats);
831
832                 sprintf(histoName, "%sBoard%s%i", cathCode[cath], histoTypeName[hType], 11+ch);
833                 sprintf(histoTitle, "Chamber %i: board %s %s", 11+ch, histoTypeTitle[hType], cathCode[cath]);
834                 histoCheckBoard[chCath][hType] = new TH1F(histoName, histoTitle, kNumOfBoards, 1-0.5, kNumOfBoards+1.-0.5);
835                 histoCheckBoard[chCath][hType]->SetXTitle("boards");
836                 histoCheckBoard[chCath][hType]->SetYTitle(yAxisTitle);
837                 histoCheckBoard[chCath][hType]->GetXaxis()->SetNdivisions(kNumOfBoards);
838             }
839             // end ADDED for check
840         }
841     }
842
843     Float_t efficiency, efficiencyError;
844     Int_t bin;
845
846     for(Int_t cath=0; cath<fgkNcathodes; cath++){
847         for(Int_t ch=0; ch<fgkNchambers; ch++){
848             Int_t chCath = fgkNchambers*cath + ch;
849             for(Int_t slat=0; slat<fgkNslats; slat++){
850                 CalculateEfficiency(fHitPerSlat[chCath][slat], fHitPerSlat[chCath][slat]+fInefficientSlat[chCath][slat], efficiency, efficiencyError, kFALSE);
851                 bin = histoSlat[chCath]->FindBin(slat);
852                 histoSlat[chCath]->SetBinContent(bin, efficiency);
853                 histoSlat[chCath]->SetBinError(bin, efficiencyError);
854
855                 // ADDED for check
856                 histoCheckSlat[chCath][kAllChEff]->SetBinContent(bin, fHitPerSlat[chCath][slat]);
857                 histoCheckSlat[chCath][kChNonEff]->SetBinContent(bin, fInefficientSlat[chCath][slat]);
858             }
859             CalculateEfficiency(fTrigger44[cath], fTrigger34[chCath]+fTrigger44[cath], efficiency, efficiencyError, kFALSE);
860             bin = histoChamber[cath]->FindBin(11+ch);
861             histoChamber[cath]->SetBinContent(bin, efficiency);
862             histoChamber[cath]->SetBinError(bin, efficiencyError);
863
864             for(Int_t board=0; board<kNumOfBoards; board++){
865                 CalculateEfficiency(fHitPerBoard[chCath][board], fHitPerBoard[chCath][board]+fInefficientBoard[chCath][board], efficiency, efficiencyError, kFALSE);
866                 bin = histoBoard[chCath]->FindBin(board+1);
867                 histoBoard[chCath]->SetBinContent(bin, efficiency);
868                 histoBoard[chCath]->SetBinError(bin, efficiencyError);
869
870                 // ADDED for check
871                 histoCheckBoard[chCath][kAllChEff]->SetBinContent(bin, fHitPerBoard[chCath][board]);
872                 histoCheckBoard[chCath][kChNonEff]->SetBinContent(bin, fInefficientBoard[chCath][board]);
873             }
874         }
875     }
876
877 // write all histos
878     outputHistoFile->cd();
879     dir->GetList()->Write();
880     outputHistoFile->Close();
881 }
882
883
884 //_____________________________________________________________________________
885 void AliMUONTriggerChamberEff::WriteEfficiencyMapTxt(const char* outputDir)
886 {
887     //
888     /// Writes the calculated efficiency in the text file efficiencyCells.dat
889     ///
890     /// The file can be further put in $ALICE_ROOT/MUON/data
891     /// and used to run simulations with measured trigger chamber efficiencies.
892     //
893
894     Int_t effOutWidth=4;
895
896     Float_t efficiency, efficiencyError;
897
898     Int_t aCapo[] = {16, 38, 60, 76, 92, 108, 117, 133, 155, 177, 193, 209, 225, 234};
899
900     char filename[70];
901     sprintf(filename, "%s/efficiencyCells.dat", outputDir);
902     ofstream outFile(filename);
903     outFile << "localBoards" << endl;
904     for(Int_t ch=0; ch<fgkNchambers; ch++){
905         //Print information
906         outFile << "\n\ndetElemId:\t" << 11+ch;
907         outFile << "00";
908         for(Int_t cath=0; cath<fgkNcathodes; cath++){
909             outFile << "\n cathode:\t" << cath << endl;
910             Int_t chCath = fgkNchambers*cath + ch;
911             Int_t currLine=0;
912             for(Int_t board=0; board<AliMpConstants::NofLocalBoards(); board++){
913
914                 if(board==aCapo[currLine]){
915                     outFile << endl;
916                     currLine++;
917                 }
918                 CalculateEfficiency(fHitPerBoard[chCath][board], fHitPerBoard[chCath][board]+fInefficientBoard[chCath][board], efficiency, efficiencyError, kFALSE);
919                 outFile << " " << setw(effOutWidth) << efficiency;
920             }// loop on boards
921             outFile << endl;
922         }// loop on cathodes
923     }// loop on chambers    
924 }
925
926
927 //_____________________________________________________________________________
928 Bool_t AliMUONTriggerChamberEff::IsCleanTrack(AliMUONTriggerTrack *triggerTrack,
929                                               const AliMUONVTrackStore& trackStore)
930 {
931     //
932     /// Try to match track from tracking system with trigger track
933     //
934     const Double_t kDistSigma[3]={1,1,0.02}; // sigma of distributions (trigger-track) X,Y,slopeY
935     const Double_t kMaxChi2MatchTrigger = 16.0;
936   
937     AliMUONTrackParam trackParam; 
938
939     Double_t distTriggerTrack[3];
940     Double_t xTrack, yTrack, ySlopeTrack, chi2;
941   
942     AliMUONTrack* track;
943     TIter next(trackStore.CreateIterator());
944     
945     while ( ( track = static_cast<AliMUONTrack*>(next()) ) )
946     {
947         trackParam = *((AliMUONTrackParam*) (track->GetTrackParamAtHit()->Last()));
948         AliMUONTrackExtrap::ExtrapToZ(&trackParam, AliMUONConstants::DefaultChamberZ(10)); // extrap to 1st trigger chamber
949     
950         xTrack = trackParam.GetNonBendingCoor();
951         yTrack = trackParam.GetBendingCoor();
952         ySlopeTrack = trackParam.GetBendingSlope();
953   
954         distTriggerTrack[0] = (triggerTrack->GetX11()-xTrack)/kDistSigma[0];
955         distTriggerTrack[1] = (triggerTrack->GetY11()-yTrack)/kDistSigma[1];
956         distTriggerTrack[2] = (TMath::Tan(triggerTrack->GetThetay())-ySlopeTrack)/kDistSigma[2];
957         chi2 = 0.;
958         for (Int_t iVar = 0; iVar < 3; iVar++) chi2 += distTriggerTrack[iVar]*distTriggerTrack[iVar];
959         chi2 /= 3.; // Normalized Chi2: 3 degrees of freedom (X,Y,slopeY)
960         if (chi2 < kMaxChi2MatchTrigger) return kTRUE;
961     }
962
963     return kFALSE;
964 }
965
966
967 //_____________________________________________________________________________
968 void AliMUONTriggerChamberEff::SaveInESDFile()
969 {
970     //
971     /// Store AliMUONTriggerChamberEff in esd file
972     //
973     TDirectory *dir = gDirectory;
974     TFile *logFile = 0x0;
975     TSeqCollection *list = gROOT->GetListOfFiles();
976     Int_t n = list->GetEntries();
977     for(Int_t i=0; i<n; i++) {
978         logFile = (TFile*)list->At(i);
979         if (strstr(logFile->GetName(), "AliESDs.root")) break;
980     }
981     if(logFile){
982         TTree *esdTree = (TTree*)logFile->Get("esdTree");
983         if(esdTree){
984             if(!esdTree->GetUserInfo()->FindObject("AliMUONTriggerChamberEff")){
985                 AliInfo(Form("Adding AliMUONTrigChamberEff in %s",logFile->GetName()));
986                 esdTree->GetUserInfo()->Add(this->Clone());
987                 esdTree->Write("",TObject::kOverwrite);
988             }
989         }
990     }
991     dir->cd();
992 }
993
994
995 //_____________________________________________________________________________
996 void AliMUONTriggerChamberEff::DisplayEfficiency(Bool_t perSlat)
997 {
998     //
999     /// Display calculated efficiency.
1000     //
1001
1002     const Int_t kNumOfBoards = AliMpConstants::NofLocalBoards();
1003
1004     Int_t side, col, line, nbx, slat;
1005     Float_t xCenter, yCenter, zCenter, xWidth, yWidth;
1006     Int_t x1, y1, x2, y2, board=0;
1007     char name[8], text[200];
1008
1009     gStyle->SetPalette(1);
1010
1011     TString boardName[234];
1012
1013     // instanciate the elec. mapping class
1014     AliMpDDLStore* ddlStore = AliMpDDLStore::Instance();
1015
1016     // loop over the trigger DDL (Right: 20, Left: 21)
1017     for (Int_t iDDL = 20; iDDL <= 21; ++iDDL) {
1018
1019       // get ddl object
1020       AliMpDDL* ddl = ddlStore->GetDDL(iDDL);
1021
1022       Int_t nCrate = ddl->GetNofTriggerCrates();
1023     
1024       // loop over the number of crate in DDL
1025       for (Int_t index = 0; index < nCrate; ++index) {
1026
1027         // get crate object
1028         AliMpTriggerCrate* crate = ddlStore->GetTriggerCrate(iDDL, index);
1029
1030         Int_t nLocal = crate->GetNofLocalBoards();
1031
1032         for (Int_t iLocal = 0; iLocal  < nLocal; ++iLocal) {
1033
1034           // get local board Id from crate object
1035           board = crate->GetLocalBoardId(iLocal);
1036           if(board>kNumOfBoards || board<=0) continue;
1037
1038           AliMpLocalBoard* localBoard = ddlStore->GetLocalBoard(board);
1039           boardName[board-1] = localBoard->GetName();
1040         }
1041       }
1042     }
1043
1044     char *cathCode[fgkNcathodes] = {"bendPlane", "nonBendPlane"};
1045
1046     Float_t boardsX = 257.00;  // cm
1047     Float_t boardsY = 307.00;  // cm
1048
1049     TH2F *histo[8];
1050     TPaveLabel *boardLabel[8][234];
1051
1052     char histoName[40];
1053     char histoTitle[90];
1054
1055     Float_t efficiency, efficiencyError;
1056
1057     for(Int_t cath=0; cath<fgkNcathodes; cath++){
1058         for(Int_t ch=0; ch<fgkNchambers; ch++){
1059             Int_t chCath = fgkNchambers*cath + ch;
1060             sprintf(histoName, "%sChamber%i", cathCode[cath], 11+ch);
1061             sprintf(histoTitle, "Chamber %i: efficiency %s", 11+ch, cathCode[cath]);
1062             histo[chCath] = new TH2F(histoName, histoTitle, (Int_t)boardsX, -boardsX, boardsX, (Int_t)boardsY, -boardsY, boardsY);
1063             histo[chCath]->SetXTitle("X (cm)");
1064             histo[chCath]->SetYTitle("Y (cm)");
1065         }
1066     }
1067
1068     TString mapspath = gSystem->Getenv("ALICE_ROOT");
1069     mapspath.Append("/MUON/data");
1070
1071     sprintf(text,"%s/guimapp11.txt",mapspath.Data());
1072     FILE *fmap = fopen(text,"r");
1073
1074     for (Int_t ib = 0; ib < kNumOfBoards; ib++) {
1075         fscanf(fmap,"%d   %d   %d   %d   %f   %f   %f   %f   %f   %s   \n",&side,&col,&line,&nbx,&xCenter,&yCenter,&xWidth,&yWidth,&zCenter,&name[0]);
1076
1077         slat = (line-5)%fgkNslats;
1078         for(Int_t iBoard=0; iBoard<kNumOfBoards; iBoard++){
1079             if(!boardName[iBoard].Contains(name)) continue;
1080             board = iBoard;
1081             break;
1082         }
1083
1084         for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
1085             x1 = histo[chCath]->GetXaxis()->FindBin(xCenter-xWidth/2.)+1;
1086             y1 = histo[chCath]->GetYaxis()->FindBin(yCenter-yWidth/2.)+1;
1087             x2 = histo[chCath]->GetXaxis()->FindBin(xCenter+xWidth/2.)-1;
1088             y2 = histo[chCath]->GetYaxis()->FindBin(yCenter+yWidth/2.)-1;
1089             
1090             boardLabel[chCath][board] = new TPaveLabel(xCenter-xWidth/2., yCenter-yWidth/2., xCenter+xWidth/2., yCenter+yWidth/2., Form("%3d",board+1));
1091             boardLabel[chCath][board]->SetFillStyle(0);
1092             boardLabel[chCath][board]->SetBorderSize(0);
1093
1094             if(!perSlat){
1095                 CalculateEfficiency(fHitPerBoard[chCath][board], fHitPerBoard[chCath][board]+fInefficientBoard[chCath][board], efficiency, efficiencyError, kFALSE);
1096             }
1097             else{
1098                 CalculateEfficiency(fHitPerSlat[chCath][slat], fHitPerSlat[chCath][slat]+fInefficientSlat[chCath][slat], efficiency, efficiencyError, kFALSE);
1099             }
1100             
1101             for(Int_t binX=x1; binX<=x2; binX++){
1102                 for(Int_t binY=y1; binY<=y2; binY++){
1103                     histo[chCath]->SetBinContent(binX, binY, efficiency);
1104                     histo[chCath]->SetBinError(binX, binY, efficiencyError);
1105                 }
1106             }
1107         }
1108     }
1109
1110     TCanvas *can[8];
1111     for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
1112         sprintf(histoName, "%sCan", histo[chCath]->GetName());
1113         sprintf(histoTitle, "%s", histo[chCath]->GetTitle());
1114         can[chCath] = new TCanvas(histoName, histoTitle, 100+10*chCath, 10*chCath, 700, 700);
1115         can[chCath]->SetRightMargin(0.14);
1116         can[chCath]->SetLeftMargin(0.12);
1117         histo[chCath]->GetZaxis()->SetRangeUser(0.,1.);
1118         histo[chCath]->GetYaxis()->SetTitleOffset(1.4);
1119         histo[chCath]->SetStats(kFALSE);
1120         histo[chCath]->Draw("COLZ");
1121         for (Int_t board = 0; board < kNumOfBoards; board++) {
1122             boardLabel[chCath][board]->Draw("same");
1123         }
1124     }
1125 }
1126
1127 //_____________________________________________________________________________
1128 void AliMUONTriggerChamberEff::CheckConstants() const
1129 {
1130 /// Check consistence of redefined constants 
1131
1132   assert(fgkNcathodes == AliMpConstants::NofCathodes());    
1133   assert(fgkNchambers == AliMpConstants::NofTriggerChambers());    
1134   assert(fgkNplanes == AliMpConstants::NofTriggerChambers() * fgkNcathodes);    
1135 }