]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONQADataMakerRec.cxx
Setting the reco-params for the QADataMakerRec objects. At least we managed to remove...
[u/mrichter/AliRoot.git] / MUON / AliMUONQADataMakerRec.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 // --- MUON header files ---
19 #include "AliMUONQADataMakerRec.h"
20
21 #include "AliMUON2DMap.h"
22 #include "AliMUONCluster.h"  
23 #include "AliMUONConstants.h"  
24 #include "AliMUONDDLTrigger.h"
25 #include "AliMUONDarcHeader.h"
26 #include "AliMUONDigitMaker.h"
27 #include "AliMUONLocalStruct.h"
28 #include "AliMUONLocalTrigger.h"
29 #include "AliMUONRawStreamTracker.h"
30 #include "AliMUONRawStreamTrigger.h"
31 #include "AliMUONRegHeader.h"
32 #include "AliMUONTrackerCalibratedDataMaker.h"
33 #include "AliMUONTriggerDisplay.h"
34 #include "AliMUONVCluster.h"
35 #include "AliMUONVClusterStore.h"
36 #include "AliMUONVDigit.h"
37 #include "AliMUONVDigitStore.h"
38 #include "AliMUONVTrackerData.h"
39 #include "AliMUONVTriggerStore.h"
40 #include "AliMpCDB.h"
41 #include "AliMpConstants.h"
42 #include "AliMpDDLStore.h"
43 #include "AliMpDEIterator.h"
44 #include "AliMpDEManager.h"
45 #include "AliMpLocalBoard.h"
46 #include "AliMpStationType.h"
47 #include "AliMpTriggerCrate.h"
48 #include "AliRawEventHeaderBase.h"
49
50 // --- AliRoot header files ---
51 #include "AliCDBManager.h"
52 #include "AliCDBStorage.h"
53 #include "AliESDEvent.h"
54 #include "AliESDMuonTrack.h"
55 #include "AliESDMuonCluster.h"
56 #include "AliLog.h"
57 #include "AliRawReader.h"
58 #include "AliQAChecker.h"
59 #include "AliCodeTimer.h"
60
61 // --- ROOT system ---
62 #include <TClonesArray.h>
63 #include <TFile.h> 
64 #include <TH1F.h> 
65 #include <TH1I.h> 
66 #include <TH2F.h>
67 #include <TH3F.h> 
68 #include <TLorentzVector.h>
69 #include <Riostream.h>
70
71 //-----------------------------------------------------------------------------
72 /// \class AliMUONQADataMakerRec
73 ///
74 /// MUON base class for quality assurance data (histo) maker
75 ///
76 /// \author C. Finck, D. Stocco, L. Aphecetche
77
78 /// \cond CLASSIMP
79 ClassImp(AliMUONQADataMakerRec)
80 /// \endcond
81            
82 //____________________________________________________________________________ 
83 AliMUONQADataMakerRec::AliMUONQADataMakerRec() : 
84 AliQADataMakerRec(AliQA::GetDetName(AliQA::kMUON), "MUON Quality Assurance Data Maker"),
85 fIsInitRaws(kFALSE),
86 fIsInitRecPointsTracker(kFALSE),
87 fIsInitRecPointsTrigger(kFALSE),
88 fIsInitESDs(kFALSE),
89 fDigitStore(0x0),
90 fTriggerStore(0x0),
91 fDigitMaker(0x0),
92 fClusterStore(0x0),
93 fTrackerDataMaker(0x0)
94 {
95     /// ctor
96         
97         Ctor();
98 }
99
100 //____________________________________________________________________________ 
101 void
102 AliMUONQADataMakerRec::Ctor()
103 {
104         /// Init some members
105         fDigitStore = AliMUONVDigitStore::Create("AliMUONDigitStoreV1");
106         fDigitMaker = new AliMUONDigitMaker(kTRUE);
107 }
108
109 //____________________________________________________________________________ 
110 AliMUONQADataMakerRec::AliMUONQADataMakerRec(const AliMUONQADataMakerRec& qadm) :
111 AliQADataMakerRec(qadm),
112 fIsInitRaws(kFALSE),
113 fIsInitRecPointsTracker(kFALSE),
114 fIsInitRecPointsTrigger(kFALSE),
115 fIsInitESDs(kFALSE),
116 fDigitStore(0x0),
117 fTriggerStore(0x0),
118 fDigitMaker(0x0),
119 fClusterStore(0x0),
120 fTrackerDataMaker(0x0)
121 {
122     ///copy ctor 
123     SetName((const char*)qadm.GetName()) ; 
124     SetTitle((const char*)qadm.GetTitle()); 
125
126         // Do not copy the digit store and digit maker, but create its own ones
127         
128         Ctor();
129         
130 }
131
132 //__________________________________________________________________
133 AliMUONQADataMakerRec& AliMUONQADataMakerRec::operator = (const AliMUONQADataMakerRec& qadm )
134 {
135   /// Assignment operator
136
137   // check assignment to self
138   if (this == &qadm) return *this;
139
140   this->~AliMUONQADataMakerRec();
141   new(this) AliMUONQADataMakerRec(qadm);
142   return *this;
143 }
144
145 //__________________________________________________________________
146 AliMUONQADataMakerRec::~AliMUONQADataMakerRec()
147 {
148     /// dtor
149   
150   AliCodeTimerAuto("");
151   
152   delete fDigitStore;
153   delete fTriggerStore;
154   delete fDigitMaker;
155         delete fClusterStore;
156         delete fTrackerDataMaker;
157 }
158
159 //____________________________________________________________________________ 
160 void AliMUONQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray* list)
161 {
162         ///Detector specific actions at end of cycle
163         
164         // Display trigger histos in a more user friendly way
165         DisplayTriggerInfo(task);
166         
167         if ( task == AliQA::kRAWS ) 
168         {
169                 TIter next(list);
170                 TObject* o;
171                 Bool_t alreadyThere(kFALSE);
172                 while ( ( o = next() ) && !alreadyThere )
173                 {
174                         TString classname(o->ClassName());
175                         if ( classname.Contains("TrackerData") ) alreadyThere = kTRUE;
176                 }
177                 if (!alreadyThere) list->AddAt(fTrackerDataMaker->Data(),(Int_t)kTrackerData);
178         }
179         
180         // do the QA checking
181         AliQAChecker::Instance()->Run(AliQA::kMUON, task, list) ;
182 }
183
184 //____________________________________________________________________________ 
185 void AliMUONQADataMakerRec::InitRaws()
186 {
187     /// create Raws histograms in Raws subdir
188         
189         if ( ! AliCDBManager::Instance()->GetDefaultStorage() )
190         {
191                 AliError("CDB default storage not set. Cannot work.");
192                 fIsInitRaws=kFALSE;
193         }
194         
195         TH3F* h3 = new TH3F("hTriggerScalersBendPlane", "Trigger scalers in bending plane",
196                                                                                         4, 10.5, 14.5,
197                                                                                         234, 0.5, 234.5,
198                                                                                         16, -0.5, 15.5);
199         h3->GetXaxis()->SetTitle("Chamber");
200         h3->GetYaxis()->SetTitle("Board");
201         h3->GetZaxis()->SetTitle("Strip");
202         Add2RawsList(h3, kTriggerScalersBP);
203         
204         TH3F* h4 = new TH3F("hTriggerScalersNonBendPlane", "Trigger scalers in non-bending plane",
205                                                                                         4, 10.5, 14.5,
206                                                                                         234, 0.5, 234.5,
207                                                                                         16, -0.5, 15.5);
208         h4->GetXaxis()->SetTitle("Chamber");
209         h4->GetYaxis()->SetTitle("Board");
210         h4->GetZaxis()->SetTitle("Strip");
211         Add2RawsList(h4, kTriggerScalersNBP);
212         
213         AliMUONTriggerDisplay triggerDisplay;
214         TString histoName, histoTitle;
215         for(Int_t iCath=0; iCath<AliMpConstants::NofCathodes(); iCath++){
216                 TString cathName = ( iCath==0 ) ? "BendPlane" : "NonBendPlane";
217                 for(Int_t iChamber=0; iChamber<AliMpConstants::NofTriggerChambers(); iChamber++){
218                         histoName = Form("hScalers%sChamber%i", cathName.Data(), 11+iChamber);
219                         histoTitle = Form("Chamber %i: Scalers %s", 11+iChamber, cathName.Data());
220                         TH2F* h5 = (TH2F*)triggerDisplay.GetEmptyDisplayHisto(histoName, AliMUONTriggerDisplay::kDisplayStrips, 
221                                                                                                                                                                                                                                                 iCath, iChamber, histoTitle);
222                         Add2RawsList(h5, kTriggerScalersDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
223                 }
224         }
225         
226         fIsInitRaws = kTRUE;
227 }
228
229 //____________________________________________________________________________ 
230 void AliMUONQADataMakerRec::InitRecPoints()
231 {
232         /// create Reconstructed Points histograms in RecPoints subdir
233         InitRecPointsTrigger();
234         InitRecPointsTracker();
235 }
236
237 //____________________________________________________________________________ 
238 void AliMUONQADataMakerRec::InitRecPointsTracker()
239 {
240         /// create Reconstructed Points histograms in RecPoints subdir for the
241         /// MUON tracker subsystem.
242
243         AliMpDEIterator it;
244         
245         it.First();
246         
247         Int_t ndes(0);
248         
249         while ( !it.IsDone())
250         {
251                 Int_t detElemId = it.CurrentDEId();
252                 
253                 it.Next();
254
255                 if ( AliMpDEManager::GetStationType(detElemId) != AliMp::kStationTrigger )
256                 {
257                         ndes = TMath::Max(ndes,detElemId);
258
259                         TH1* h = new TH1I(Form("hTrackerClusterMultiplicityForDE%04d",detElemId),
260                                                                                                 Form("Multiplicity of the clusters in detection element %d",detElemId),
261                                                                                                 100,0,100);
262                         
263                         h->GetXaxis()->SetTitle("Detection Element Id");
264                         
265                         Add2RecPointsList(h,kTrackerClusterMultiplicityPerDE+detElemId);
266                         
267                         h =  new TH1I(Form("hTrackerClusterChargeForDE%04d",detElemId),
268                                                                                 Form("Charge of the clusters in detection element %d",detElemId),
269                                                                                 100,0,1000);
270
271                         h->GetXaxis()->SetTitle("Detection Element Id");
272
273                         Add2RecPointsList(h,kTrackerClusterChargePerDE+detElemId);
274
275                 }
276
277         }
278
279         TH1* h = new TH1I("hTrackerNumberOfClustersPerDE","Number of clusters per detection element",
280                                                                                 ndes, -0.5, ndes - 0.5);
281
282         h->GetXaxis()->SetTitle("Detection Element Id");
283
284         Add2RecPointsList(h, kTrackerNumberOfClustersPerDE);
285
286         for ( Int_t i = 0; i < AliMpConstants::NofTrackingChambers(); ++i ) 
287         {
288                 TH1* h1 = new TH1I("hTrackerNumberOfClustersPerChamber","Number of clusters per chamber",AliMpConstants::NofTrackingChambers(),-0.5,AliMpConstants::NofTrackingChambers()-0.5);
289                 Add2RecPointsList(h1,kTrackerNumberOfClustersPerChamber);
290                 h1 = new TH1I(Form("hTrackerClusterMultiplicityForChamber%d",i),
291                                                                  Form("Cluster multiplicity for chamber %d",i),
292                                                                  100,0,100);
293                 Add2RecPointsList(h1,kTrackerClusterMultiplicityPerChamber+i);
294                 h1 = new TH1I(Form("hTrackerClusterChargeForChamber%d",i),
295                                                                  Form("Cluster charge for chamber %d",i),
296                                                                  100,0,1000);
297                 Add2RecPointsList(h1,kTrackerClusterChargePerChamber+i);
298         }
299         
300         fIsInitRecPointsTracker=kTRUE;
301 }
302
303 //____________________________________________________________________________ 
304 void AliMUONQADataMakerRec::InitRecPointsTrigger()
305 {
306         /// create Reconstructed Points histograms in RecPoints subdir for the
307         /// MUON Trigger subsystem.
308         
309     TH3F* h0 = new TH3F("hTriggerDigitsBendPlane", "Trigger digits in bending plane",
310                         4, 10.5, 14.5,
311                         234, 0.5, 234.5,
312                         16, -0.5, 15.5);
313     h0->GetXaxis()->SetTitle("Chamber");
314     h0->GetYaxis()->SetTitle("Board");
315     h0->GetZaxis()->SetTitle("Strip");
316     Add2RecPointsList(h0, kTriggerDigitsBendPlane);
317
318     TH3F* h1 = new TH3F("hTriggerDigitsNonBendPlane", "Trigger digits in non-bending plane",
319                         4, 10.5, 14.5,
320                         234, 0.5, 234.5,
321                         16, -0.5, 15.5);
322     h1->GetXaxis()->SetTitle("Chamber");
323     h1->GetYaxis()->SetTitle("Board");
324     h1->GetZaxis()->SetTitle("Strip");
325     Add2RecPointsList(h1, kTriggerDigitsNonBendPlane);
326
327     TH1F* h2 = new TH1F("hTriggeredBoards", "Triggered boards", 234, 0.5, 234.5);
328     Add2RecPointsList(h2, kTriggeredBoards);
329
330     AliMUONTriggerDisplay triggerDisplay;
331     TString histoName, histoTitle;
332     for(Int_t iCath=0; iCath<AliMpConstants::NofCathodes(); iCath++){
333       TString cathName = ( iCath==0 ) ? "BendPlane" : "NonBendPlane";
334       for(Int_t iChamber=0; iChamber<AliMpConstants::NofTriggerChambers(); iChamber++){
335         histoName = Form("hTriggerDigits%sChamber%i", cathName.Data(), 11+iChamber);
336         histoTitle = Form("Chamber %i: Fired pads %s", 11+iChamber, cathName.Data());
337         TH2F* h3 = (TH2F*)triggerDisplay.GetEmptyDisplayHisto(histoName, AliMUONTriggerDisplay::kDisplayStrips, 
338                                                               iCath, iChamber, histoTitle);
339         Add2RecPointsList(h3, kTriggerDigitsDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
340       }
341     }
342
343     TH2F* h4 = (TH2F*)triggerDisplay.GetEmptyDisplayHisto("hFiredBoardsDisplay", AliMUONTriggerDisplay::kDisplayBoards,
344                                                           0, 0, "Fired boards");
345     Add2RecPointsList(h4, kTriggerBoardsDisplay);
346         
347         fIsInitRecPointsTrigger = kTRUE;
348 }
349
350
351 //____________________________________________________________________________ 
352 void AliMUONQADataMakerRec::InitESDs()
353 {
354     ///create ESDs histograms in ESDs subdir
355   TH1F* h0 = new TH1F("hESDnTracks", "ESDs track number distribution", 30, 0., 30.);  
356   Add2ESDsList(h0, kESDnTracks);
357
358   TH1F* h1 = new TH1F("hESDMomentum", "ESDs P distribution", 300, 0., 300) ; 
359   Add2ESDsList(h1, kESDMomentum);
360
361   TH1F* h2 = new TH1F("hESDPt", "ESDs Pt distribution", 200, 0., 50) ; 
362   Add2ESDsList(h2, kESDPt);
363
364   TH1F* h3 = new TH1F("hESDRapidity", "ESDs rapidity distribution", 200, -4.5,-2.) ; 
365   Add2ESDsList(h3, kESDRapidity);
366
367   for (Int_t i = 0; i < AliMUONConstants::NTrackingCh(); ++i) 
368   {
369     TH2F* h4 = new TH2F(Form("%s%d", "hESDClusterHitMap", i), 
370                      Form("%s %d", "ESD Clusters hit distribution for chamber", i),
371                      100, -1*AliMUONConstants::Rmax(i/2), AliMUONConstants::Rmax(i/2),
372                      100, -1*AliMUONConstants::Rmax(i/2), AliMUONConstants::Rmax(i/2)); 
373     Add2ESDsList(h4, kESDClusterHitMap+i);
374   }
375   
376   fIsInitESDs =  kTRUE;
377 }
378
379 //____________________________________________________________________________
380 void AliMUONQADataMakerRec::MakeRaws(AliRawReader* rawReader)
381 {
382     /// make QA for rawdata
383
384     if ( ! fIsInitRaws ) {
385       AliWarningStream() 
386         << "Skipping function due to a failure in Init" << endl;
387       return;
388     }    
389
390   if ( rawReader->GetType() == AliRawEventHeaderBase::kPhysicsEvent ) 
391   {
392     rawReader->Reset();
393     MakeRawsTracker(rawReader);
394   }
395   
396   rawReader->Reset();    
397   MakeRawsTrigger(rawReader);
398 }
399
400 //____________________________________________________________________________
401 void AliMUONQADataMakerRec::MakeRawsTracker(AliRawReader* rawReader)
402 {
403         /// make QA for rawdata tracker
404   
405         if (!fTrackerDataMaker) 
406         {
407                 const Bool_t histogram(kFALSE);
408                 const Bool_t fastDecoder(kTRUE);
409     
410 //    fTrackerDataMaker = new AliMUONTrackerRawDataMaker(rawReader,histogram,fastDecoder,takeRawReaderOwnership);
411
412                 fTrackerDataMaker = new AliMUONTrackerCalibratedDataMaker(GetRecoParam(),
413                                                               AliCDBManager::Instance()->GetRun(),
414                                                               rawReader,
415                                                                                                                                                                                                                                                         AliCDBManager::Instance()->GetDefaultStorage()->GetURI(),
416                                                                                                                                                                                                                                                         "NOGAIN",
417                                                                                                                                                                                                                                                         histogram,
418                                                                                                                                                                                                                                                         0.0,0.0,
419                                                               fastDecoder);
420                 
421                 fTrackerDataMaker->Data()->DisableChannelLevel(); // to save up disk space, we only store starting at the manu level
422                 
423                 fTrackerDataMaker->SetRunning(kTRUE);
424         }
425         
426         ((AliMUONTrackerCalibratedDataMaker*)fTrackerDataMaker)->SetRawReader(rawReader);
427         
428         fTrackerDataMaker->ProcessEvent();
429 }
430
431 //____________________________________________________________________________
432 void AliMUONQADataMakerRec::MakeRawsTrigger(AliRawReader* rawReader)
433 {
434         /// make QA for rawdata trigger
435         
436     // Get trigger scalers
437
438     Int_t loCircuit=0;
439     AliMpCDB::LoadDDLStore();
440
441     AliMUONRawStreamTrigger rawStreamTrig(rawReader);
442     while (rawStreamTrig.NextDDL()) 
443     {
444       // If not a scaler event, do nothing
445       Bool_t scalerEvent =  rawReader->GetDataHeader()->GetL1TriggerMessage() & 0x1;
446       if(!scalerEvent) break;
447
448       AliMUONDDLTrigger* ddlTrigger = rawStreamTrig.GetDDLTrigger();
449       AliMUONDarcHeader* darcHeader = ddlTrigger->GetDarcHeader();
450
451       Int_t nReg = darcHeader->GetRegHeaderEntries();
452     
453       for(Int_t iReg = 0; iReg < nReg ;iReg++)
454       {   //reg loop
455
456         // crate info  
457         AliMpTriggerCrate* crate = AliMpDDLStore::Instance()->
458           GetTriggerCrate(rawStreamTrig.GetDDL(), iReg);
459
460         AliMUONRegHeader* regHeader =  darcHeader->GetRegHeaderEntry(iReg);
461
462         // loop over local structures
463         Int_t nLocal = regHeader->GetLocalEntries();
464         for(Int_t iLocal = 0; iLocal < nLocal; iLocal++) 
465         {
466           AliMUONLocalStruct* localStruct = regHeader->GetLocalEntry(iLocal);
467         
468           // if card exist
469           if (!localStruct) continue;
470           
471           loCircuit = crate->GetLocalBoardId(localStruct->GetId());
472
473           if ( !loCircuit ) continue; // empty slot
474
475           AliMpLocalBoard* localBoard = AliMpDDLStore::Instance()->GetLocalBoard(loCircuit, false);
476
477           // skip copy cards
478           if( !localBoard->IsNotified()) 
479             continue;
480
481           Int_t cathode = localStruct->GetComptXY()%2;
482
483           ERaw hindex = (cathode==0) ? kTriggerScalersBP : kTriggerScalersNBP;
484
485           // loop over strips
486           for (Int_t ibitxy = 0; ibitxy < 16; ++ibitxy) {
487             if(localStruct->GetXY1(ibitxy) > 0)
488               ((TH3F*)GetRawsData(hindex))->Fill(11+0, loCircuit, ibitxy, 2*localStruct->GetXY1(ibitxy));
489             if(localStruct->GetXY2(ibitxy) > 0)
490               ((TH3F*)GetRawsData(hindex))->Fill(11+1, loCircuit, ibitxy, 2*localStruct->GetXY2(ibitxy));
491             if(localStruct->GetXY3(ibitxy) > 0)
492               ((TH3F*)GetRawsData(hindex))->Fill(11+2, loCircuit, ibitxy, 2*localStruct->GetXY3(ibitxy));
493             if(localStruct->GetXY4(ibitxy) > 0)
494               ((TH3F*)GetRawsData(hindex))->Fill(11+3, loCircuit, ibitxy, 2*localStruct->GetXY4(ibitxy));
495           } // loop on strips
496         } // iLocal
497       } // iReg
498     } // NextDDL
499 }
500
501 //____________________________________________________________________________
502 void AliMUONQADataMakerRec::MakeRecPoints(TTree* clustersTree)
503 {
504         /// Fill histograms from treeR
505         
506         if (fIsInitRecPointsTracker) MakeRecPointsTracker(clustersTree);
507         if (fIsInitRecPointsTrigger) MakeRecPointsTrigger(clustersTree);
508 }
509
510 //____________________________________________________________________________
511 void AliMUONQADataMakerRec::MakeRecPointsTracker(TTree* clustersTree)
512 {
513         /// Fill histograms related to tracker clusters 
514         
515         // First things first : do we have clusters in the TreeR ?
516         // In "normal" production mode, it should be perfectly normal
517         // *not* to have them.
518         // But if for some reason we de-activated the combined tracking,
519         // then we have clusters in TreeR, so let's take that opportunity
520         // to QA them...
521         
522         if (!fClusterStore)
523         {
524                 AliCodeTimerAuto("ClusterStore creation");
525                 fClusterStore = AliMUONVClusterStore::Create(*clustersTree);
526                 if (!fClusterStore) 
527                 {
528                         fIsInitRecPointsTracker = kFALSE;
529                         return;
530                 }
531         }
532         
533         AliCodeTimerAuto("");
534         
535         fClusterStore->Connect(*clustersTree,kFALSE);
536         clustersTree->GetEvent(0);
537
538         TIter next(fClusterStore->CreateIterator());
539         AliMUONVCluster* cluster;
540         
541         while ( ( cluster = static_cast<AliMUONVCluster*>(next()) ) )
542         {
543                 Int_t detElemId = cluster->GetDetElemId();
544                 Int_t chamberId = AliMpDEManager::GetChamberId(detElemId);
545                 
546                 GetRecPointsData(kTrackerNumberOfClustersPerDE)->Fill(detElemId);
547                 GetRecPointsData(kTrackerClusterChargePerDE+detElemId)->Fill(cluster->GetCharge());
548                 GetRecPointsData(kTrackerClusterMultiplicityPerDE+detElemId)->Fill(cluster->GetNDigits());
549
550                 GetRecPointsData(kTrackerNumberOfClustersPerChamber)->Fill(chamberId);
551                 GetRecPointsData(kTrackerClusterChargePerChamber+chamberId)->Fill(cluster->GetCharge());
552                 GetRecPointsData(kTrackerClusterMultiplicityPerChamber+chamberId)->Fill(cluster->GetNDigits());
553                 
554         }
555         
556         fClusterStore->Clear();
557 }
558
559 //____________________________________________________________________________
560 void AliMUONQADataMakerRec::MakeRecPointsTrigger(TTree* clustersTree)
561 {
562         /// makes data from trigger response
563       
564     // Fired pads info
565     fDigitStore->Clear();
566
567     if (!fTriggerStore) fTriggerStore = AliMUONVTriggerStore::Create(*clustersTree);
568     fTriggerStore->Clear();
569     fTriggerStore->Connect(*clustersTree, false);
570     clustersTree->GetEvent(0);
571
572     AliMUONLocalTrigger* locTrg;
573     TIter nextLoc(fTriggerStore->CreateLocalIterator());
574
575     while ( ( locTrg = static_cast<AliMUONLocalTrigger*>(nextLoc()) ) ) 
576     {
577       if (locTrg->IsNull()) continue;
578    
579       TArrayS xyPattern[2];
580       locTrg->GetXPattern(xyPattern[0]);
581       locTrg->GetYPattern(xyPattern[1]);
582
583       Int_t nBoard = locTrg->LoCircuit();
584
585       Bool_t xTrig=locTrg->IsTrigX();
586       Bool_t yTrig=locTrg->IsTrigY();
587     
588       if (xTrig && yTrig)
589         ((TH1F*)GetRecPointsData(kTriggeredBoards))->Fill(nBoard);
590
591       fDigitMaker->TriggerDigits(nBoard, xyPattern, *fDigitStore);
592     }
593
594     TIter nextDigit(fDigitStore->CreateIterator());
595     AliMUONVDigit* mDigit;
596     while ( ( mDigit = static_cast<AliMUONVDigit*>(nextDigit()) ) )
597     {
598       Int_t detElemId = mDigit->DetElemId();
599       Int_t ch = detElemId/100;
600       Int_t localBoard = mDigit->ManuId();
601       Int_t channel = mDigit->ManuChannel();
602       Int_t cathode = mDigit->Cathode();
603       ERecPoints hindex 
604         = ( cathode == 0 ) ? kTriggerDigitsBendPlane : kTriggerDigitsNonBendPlane;
605       
606       ((TH3F*)GetRecPointsData(hindex))->Fill(ch, localBoard, channel);
607     }
608 }
609
610 //____________________________________________________________________________
611 void AliMUONQADataMakerRec::MakeESDs(AliESDEvent* esd)
612 {
613     /// make QA data from ESDs
614
615     if ( ! fIsInitESDs ) {
616       AliWarningStream() 
617         << "Skipping function due to a failure in Init" << endl;
618       return;
619     }    
620
621     TLorentzVector v1;
622
623     Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ; 
624     GetESDsData(0)->Fill(nTracks);
625
626     for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) {
627
628       AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
629       
630       // skip "ghosts"
631       if (!muonTrack->ContainTrackerData()) continue;
632       
633       muonTrack->LorentzP(v1);
634
635       GetESDsData(1)->Fill(v1.P());
636       GetESDsData(2)->Fill(v1.Pt());
637       GetESDsData(3)->Fill(v1.Rapidity());
638
639       TClonesArray clusters =  muonTrack->GetClusters();
640
641       for (Int_t iCluster = 0; iCluster <clusters.GetEntriesFast(); ++iCluster) {
642         AliESDMuonCluster* cluster = (AliESDMuonCluster*)clusters.At(iCluster);
643         GetESDsData(kESDClusterHitMap+cluster->GetChamberId())
644           ->Fill(cluster->GetX(), cluster->GetY());
645       }
646     }
647 }
648
649 //____________________________________________________________________________ 
650 void AliMUONQADataMakerRec::StartOfDetectorCycle()
651 {
652     /// Detector specific actions at start of cycle
653   
654 }
655
656 //____________________________________________________________________________ 
657 void AliMUONQADataMakerRec::DisplayTriggerInfo(AliQA::TASKINDEX_t task)
658 {
659   //
660   /// Display trigger information in a user-friendly way:
661   /// from local board and strip numbers to their position on chambers
662   //
663   if(task!=AliQA::kRECPOINTS && task!=AliQA::kRAWS) return;
664
665   AliMUONTriggerDisplay triggerDisplay;
666   
667   TH3F* histoStrips=0x0;
668   TH2F* histoDisplayStrips=0x0;
669
670   for (Int_t iCath = 0; iCath < AliMpConstants::NofCathodes(); iCath++)
671   {
672     if(task==AliQA::kRECPOINTS){
673       ERecPoints hindex 
674         = ( iCath == 0 ) ? kTriggerDigitsBendPlane : kTriggerDigitsNonBendPlane;
675       histoStrips = (TH3F*)GetRecPointsData(hindex);
676     }
677     else if(task==AliQA::kRAWS){
678       ERaw hindex 
679         = ( iCath == 0 ) ? kTriggerScalersBP : kTriggerScalersNBP;
680       histoStrips = (TH3F*)GetRawsData(hindex);
681       if(histoStrips->GetEntries()==0) return; // No scalers found
682     }
683     
684     for (Int_t iChamber = 0; iChamber < AliMpConstants::NofTriggerChambers(); iChamber++)
685     {
686       if(task==AliQA::kRECPOINTS){
687         histoDisplayStrips = (TH2F*)GetRecPointsData(kTriggerDigitsDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
688       }
689       else if(task==AliQA::kRAWS){
690         histoDisplayStrips = (TH2F*)GetRawsData(kTriggerScalersDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
691       }
692       Int_t bin = histoStrips->GetXaxis()->FindBin(11+iChamber);
693       histoStrips->GetXaxis()->SetRange(bin,bin);
694       TH2F* inputHisto = (TH2F*)histoStrips->Project3D("zy");
695       triggerDisplay.FillDisplayHistogram(inputHisto, histoDisplayStrips, AliMUONTriggerDisplay::kDisplayStrips, iCath, iChamber);
696     } // iChamber
697   } // iCath
698
699   if(task!=AliQA::kRECPOINTS) return;
700   TH1F* histoBoards = (TH1F*)GetRecPointsData(kTriggeredBoards);
701   TH2F* histoDisplayBoards = (TH2F*)GetRecPointsData(kTriggerBoardsDisplay);
702   triggerDisplay.FillDisplayHistogram(histoBoards, histoDisplayBoards, AliMUONTriggerDisplay::kDisplayBoards, 0, 0);
703 }