]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONQADataMakerRec.cxx
- Put the code to fill ESD with trigger data back in the correct framework:
[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 "AliMUONCluster.h"  
22 #include "AliMUONRawStreamTracker.h"
23 #include "AliMUONRawStreamTrigger.h"
24 #include "AliMUONConstants.h"  
25 #include "AliMUONVClusterStore.h"
26 #include "AliMUONVCluster.h"
27
28 #include "AliMUONDigitMaker.h"
29 #include "AliMUONVDigitStore.h"
30 #include "AliMUONVTriggerStore.h"
31 #include "AliMUONVDigit.h"
32 #include "AliMUONLocalTrigger.h"
33
34 #include "AliMUONDDLTrigger.h"
35 #include "AliMUONRegHeader.h"
36 #include "AliMUONDarcHeader.h"
37 #include "AliMUONLocalStruct.h"
38
39 #include "AliMUONTriggerDisplay.h"
40
41 #include "AliMpDDLStore.h"
42 #include "AliMpConstants.h"
43 #include "AliMpBusPatch.h"
44 #include "AliMpTriggerCrate.h"
45 #include "AliMpLocalBoard.h"
46 #include "AliMpCDB.h"
47
48 // --- AliRoot header files ---
49 #include "AliESDEvent.h"
50 #include "AliESDMuonTrack.h"
51 #include "AliESDMuonCluster.h"
52 #include "AliLog.h"
53 #include "AliRawReader.h"
54 #include "AliQAChecker.h"
55
56 // --- ROOT system ---
57 #include <TClonesArray.h>
58 #include <TFile.h> 
59 #include <TH1F.h> 
60 #include <TH1I.h> 
61 #include <TH2F.h>
62 #include <TH3F.h> 
63 #include <TLorentzVector.h>
64 #include <Riostream.h>
65
66 //-----------------------------------------------------------------------------
67 /// \class AliMUONQADataMakerRec
68 ///
69 /// MUON base class for quality assurance data (histo) maker
70 ///
71 /// \author C. Finck
72
73 /// \cond CLASSIMP
74 ClassImp(AliMUONQADataMakerRec)
75 /// \endcond
76            
77 //____________________________________________________________________________ 
78 AliMUONQADataMakerRec::AliMUONQADataMakerRec() : 
79     AliQADataMakerRec(AliQA::GetDetName(AliQA::kMUON), "MUON Quality Assurance Data Maker"),
80     fIsInitRaws(kFALSE),
81     fIsInitRecPoints(kFALSE),
82     fIsInitESDs(kFALSE),
83     fDigitStore(0x0),
84     fTriggerStore(0x0),
85     fDigitMaker(0x0)
86 {
87     /// ctor
88   fDigitStore = AliMUONVDigitStore::Create("AliMUONDigitStoreV1");
89   fDigitMaker = new AliMUONDigitMaker(kTRUE);
90 }
91
92 //____________________________________________________________________________ 
93 AliMUONQADataMakerRec::AliMUONQADataMakerRec(const AliMUONQADataMakerRec& qadm) :
94     AliQADataMakerRec(qadm),
95     fIsInitRaws(kFALSE),
96     fIsInitRecPoints(kFALSE),
97     fIsInitESDs(kFALSE),
98     fDigitStore(0x0),
99     fTriggerStore(0x0),
100     fDigitMaker(0x0)
101 {
102     ///copy ctor 
103     SetName((const char*)qadm.GetName()) ; 
104     SetTitle((const char*)qadm.GetTitle()); 
105
106     // Do not copy the digit store and digit maker, but create its own ones
107     fDigitStore = AliMUONVDigitStore::Create("AliMUONDigitStoreV1");
108     fDigitMaker = new AliMUONDigitMaker(kTRUE);
109 }
110
111 //__________________________________________________________________
112 AliMUONQADataMakerRec& AliMUONQADataMakerRec::operator = (const AliMUONQADataMakerRec& qadm )
113 {
114   /// Assignment operator
115
116   // check assignment to self
117   if (this == &qadm) return *this;
118
119   this->~AliMUONQADataMakerRec();
120   new(this) AliMUONQADataMakerRec(qadm);
121   return *this;
122 }
123
124 //__________________________________________________________________
125 AliMUONQADataMakerRec::~AliMUONQADataMakerRec()
126 {
127     /// dtor
128     delete fDigitStore;
129     delete fTriggerStore;
130     delete fDigitMaker;
131 }
132
133 //____________________________________________________________________________ 
134 void AliMUONQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray* list)
135 {
136     ///Detector specific actions at end of cycle
137
138     // Display trigger histos in a more user friendly way
139     DisplayTriggerInfo(task);
140
141     // do the QA checking
142     AliQAChecker::Instance()->Run(AliQA::kMUON, task, list) ;
143 }
144
145 //____________________________________________________________________________ 
146 void AliMUONQADataMakerRec::InitRaws()
147 {
148     /// create Raws histograms in Raws subdir
149     TH1I* h0 = new TH1I("hRawBusPatch", "buspatch distribution",  1932, 1, 1932); 
150     Add2RawsList(h0, kRawBusPatch);
151
152     TH1I* h1 = new TH1I("hRawCharge", "Charge distribution in rawdata", 4096, 0, 4095); 
153     Add2RawsList(h1, kRawCharge);
154                 
155     for (Int_t iDDL = 0; iDDL < 20; ++iDDL) 
156     {
157       TH1F* h2 = new TH1F(Form("%s%d", "hRawBusPatchDDL", iDDL), Form("%s %d","RAW Buspatch distribution for DDL", iDDL), 50, 0, 49); 
158       Add2RawsList(h2, kRawBuspatchDDL+iDDL);
159     }
160
161     TH3F* h3 = new TH3F("hTriggerScalersBendPlane", "Trigger scalers in bending plane",
162                         4, 10.5, 14.5,
163                         234, 0.5, 234.5,
164                         16, -0.5, 15.5);
165     h3->GetXaxis()->SetTitle("Chamber");
166     h3->GetYaxis()->SetTitle("Board");
167     h3->GetZaxis()->SetTitle("Strip");
168     Add2RawsList(h3, kTriggerScalersBP);
169
170     TH3F* h4 = new TH3F("hTriggerScalersNonBendPlane", "Trigger scalers in non-bending plane",
171                         4, 10.5, 14.5,
172                         234, 0.5, 234.5,
173                         16, -0.5, 15.5);
174     h4->GetXaxis()->SetTitle("Chamber");
175     h4->GetYaxis()->SetTitle("Board");
176     h4->GetZaxis()->SetTitle("Strip");
177     Add2RawsList(h4, kTriggerScalersNBP);
178
179     AliMUONTriggerDisplay triggerDisplay;
180     TString histoName, histoTitle;
181     for(Int_t iCath=0; iCath<AliMpConstants::NofCathodes(); iCath++){
182       TString cathName = ( iCath==0 ) ? "BendPlane" : "NonBendPlane";
183       for(Int_t iChamber=0; iChamber<AliMpConstants::NofTriggerChambers(); iChamber++){
184         histoName = Form("hScalers%sChamber%i", cathName.Data(), 11+iChamber);
185         histoTitle = Form("Chamber %i: Scalers %s", 11+iChamber, cathName.Data());
186         TH2F* h5 = (TH2F*)triggerDisplay.GetEmptyDisplayHisto(histoName, AliMUONTriggerDisplay::kDisplayStrips, 
187                                                               iCath, iChamber, histoTitle);
188         Add2RawsList(h5, kTriggerScalersDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
189       }
190     }
191     
192     fIsInitRaws = kTRUE;
193 }
194
195 //____________________________________________________________________________ 
196 void AliMUONQADataMakerRec::InitRecPoints()
197 {
198     /// create Reconstructed Points histograms in RecPoints subdir
199     TH3F* h0 = new TH3F("hTriggerDigitsBendPlane", "Trigger digits in bending plane",
200                         4, 10.5, 14.5,
201                         234, 0.5, 234.5,
202                         16, -0.5, 15.5);
203     h0->GetXaxis()->SetTitle("Chamber");
204     h0->GetYaxis()->SetTitle("Board");
205     h0->GetZaxis()->SetTitle("Strip");
206     Add2RecPointsList(h0, kTriggerDigitsBendPlane);
207
208     TH3F* h1 = new TH3F("hTriggerDigitsNonBendPlane", "Trigger digits in non-bending plane",
209                         4, 10.5, 14.5,
210                         234, 0.5, 234.5,
211                         16, -0.5, 15.5);
212     h1->GetXaxis()->SetTitle("Chamber");
213     h1->GetYaxis()->SetTitle("Board");
214     h1->GetZaxis()->SetTitle("Strip");
215     Add2RecPointsList(h1, kTriggerDigitsNonBendPlane);
216
217     TH1F* h2 = new TH1F("hTriggeredBoards", "Triggered boards", 234, 0.5, 234.5);
218     Add2RecPointsList(h2, kTriggeredBoards);
219
220     AliMUONTriggerDisplay triggerDisplay;
221     TString histoName, histoTitle;
222     for(Int_t iCath=0; iCath<AliMpConstants::NofCathodes(); iCath++){
223       TString cathName = ( iCath==0 ) ? "BendPlane" : "NonBendPlane";
224       for(Int_t iChamber=0; iChamber<AliMpConstants::NofTriggerChambers(); iChamber++){
225         histoName = Form("hTriggerDigits%sChamber%i", cathName.Data(), 11+iChamber);
226         histoTitle = Form("Chamber %i: Fired pads %s", 11+iChamber, cathName.Data());
227         TH2F* h3 = (TH2F*)triggerDisplay.GetEmptyDisplayHisto(histoName, AliMUONTriggerDisplay::kDisplayStrips, 
228                                                               iCath, iChamber, histoTitle);
229         Add2RecPointsList(h3, kTriggerDigitsDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
230       }
231     }
232
233     TH2F* h4 = (TH2F*)triggerDisplay.GetEmptyDisplayHisto("hFiredBoardsDisplay", AliMUONTriggerDisplay::kDisplayBoards,
234                                                           0, 0, "Fired boards");
235     Add2RecPointsList(h4, kTriggerBoardsDisplay);
236
237     fIsInitRecPoints = kTRUE;
238 }
239
240
241 //____________________________________________________________________________ 
242 void AliMUONQADataMakerRec::InitESDs()
243 {
244     ///create ESDs histograms in ESDs subdir
245   TH1F* h0 = new TH1F("hESDnTracks", "ESDs track number distribution", 30, 0., 30.);  
246   Add2ESDsList(h0, kESDnTracks);
247
248   TH1F* h1 = new TH1F("hESDMomentum", "ESDs P distribution", 300, 0., 300) ; 
249   Add2ESDsList(h1, kESDMomentum);
250
251   TH1F* h2 = new TH1F("hESDPt", "ESDs Pt distribution", 200, 0., 50) ; 
252   Add2ESDsList(h2, kESDPt);
253
254   TH1F* h3 = new TH1F("hESDRapidity", "ESDs rapidity distribution", 200, -4.5,-2.) ; 
255   Add2ESDsList(h3, kESDRapidity);
256
257   for (Int_t i = 0; i < AliMUONConstants::NTrackingCh(); ++i) 
258   {
259     TH2F* h4 = new TH2F(Form("%s%d", "hESDClusterHitMap", i), 
260                      Form("%s %d", "ESD Clusters hit distribution for chamber", i),
261                      100, -1*AliMUONConstants::Rmax(i/2), AliMUONConstants::Rmax(i/2),
262                      100, -1*AliMUONConstants::Rmax(i/2), AliMUONConstants::Rmax(i/2)); 
263     Add2ESDsList(h4, kESDClusterHitMap+i);
264   }
265   
266   fIsInitESDs =  kTRUE;
267 }
268
269 //____________________________________________________________________________
270 void AliMUONQADataMakerRec::MakeRaws(AliRawReader* rawReader)
271 {
272     /// make QA for rawdata
273
274     if ( ! fIsInitRaws ) {
275       AliWarningStream() 
276         << "Skipping function due to a failure in Init" << endl;
277       return;
278     }    
279
280     Int_t busPatchId;
281     UShort_t manuId;  
282     UChar_t channelId;
283     UShort_t charge;
284
285     rawReader->Reset();
286     AliMUONRawStreamTracker rawStreamTrack(rawReader);
287     rawStreamTrack.First();
288     while( rawStreamTrack.Next(busPatchId, manuId, channelId, charge) ) {
289   
290       GetRawsData(kRawBusPatch)->Fill(busPatchId);
291       GetRawsData(kRawCharge)->Fill(charge);
292       Int_t iDDL = rawStreamTrack.GetCurentDDL();
293       GetRawsData(kRawBuspatchDDL + iDDL)->Fill(AliMpBusPatch::GetLocalBusID(busPatchId, iDDL));
294                 
295                   
296     } // Next digit
297
298
299     // Get trigger scalers
300
301     Int_t loCircuit=0;
302     AliMpCDB::LoadDDLStore();
303
304     AliMUONRawStreamTrigger rawStreamTrig(rawReader);
305     while (rawStreamTrig.NextDDL()) 
306     {
307       // If not a scaler event, do nothing
308       Bool_t scalerEvent =  rawReader->GetDataHeader()->GetL1TriggerMessage() & 0x1;
309       if(!scalerEvent) break;
310
311       AliMUONDDLTrigger* ddlTrigger = rawStreamTrig.GetDDLTrigger();
312       AliMUONDarcHeader* darcHeader = ddlTrigger->GetDarcHeader();
313
314       Int_t nReg = darcHeader->GetRegHeaderEntries();
315     
316       for(Int_t iReg = 0; iReg < nReg ;iReg++)
317       {   //reg loop
318
319         // crate info  
320         AliMpTriggerCrate* crate = AliMpDDLStore::Instance()->
321           GetTriggerCrate(rawStreamTrig.GetDDL(), iReg);
322
323         AliMUONRegHeader* regHeader =  darcHeader->GetRegHeaderEntry(iReg);
324
325         // loop over local structures
326         Int_t nLocal = regHeader->GetLocalEntries();
327         for(Int_t iLocal = 0; iLocal < nLocal; iLocal++) 
328         {
329           AliMUONLocalStruct* localStruct = regHeader->GetLocalEntry(iLocal);
330         
331           // if card exist
332           if (!localStruct) continue;
333           
334           loCircuit = crate->GetLocalBoardId(localStruct->GetId());
335
336           if ( !loCircuit ) continue; // empty slot
337
338           AliMpLocalBoard* localBoard = AliMpDDLStore::Instance()->GetLocalBoard(loCircuit, false);
339
340           // skip copy cards
341           if( !localBoard->IsNotified()) 
342             continue;
343
344           Int_t cathode = localStruct->GetComptXY()%2;
345
346           ERaw hindex = (cathode==0) ? kTriggerScalersBP : kTriggerScalersNBP;
347
348           // loop over strips
349           for (Int_t ibitxy = 0; ibitxy < 16; ++ibitxy) {
350             if(localStruct->GetXY1(ibitxy) > 0)
351               ((TH3F*)GetRawsData(hindex))->Fill(11+0, loCircuit, ibitxy, 2*localStruct->GetXY1(ibitxy));
352             if(localStruct->GetXY2(ibitxy) > 0)
353               ((TH3F*)GetRawsData(hindex))->Fill(11+1, loCircuit, ibitxy, 2*localStruct->GetXY2(ibitxy));
354             if(localStruct->GetXY3(ibitxy) > 0)
355               ((TH3F*)GetRawsData(hindex))->Fill(11+2, loCircuit, ibitxy, 2*localStruct->GetXY3(ibitxy));
356             if(localStruct->GetXY4(ibitxy) > 0)
357               ((TH3F*)GetRawsData(hindex))->Fill(11+3, loCircuit, ibitxy, 2*localStruct->GetXY4(ibitxy));
358           } // loop on strips
359         } // iLocal
360       } // iReg
361     } // NextDDL
362 }
363
364 //____________________________________________________________________________
365 void AliMUONQADataMakerRec::MakeRecPoints(TTree* clustersTree)
366 {
367   
368     /// makes data from trigger response
369       
370     if ( ! fIsInitRecPoints ) {
371       AliWarningStream() 
372         << "Skipping function due to a failure in Init" << endl;
373       return;
374     }    
375
376     // Fired pads info
377     fDigitStore->Clear();
378
379     if (!fTriggerStore) fTriggerStore = AliMUONVTriggerStore::Create(*clustersTree);
380     fTriggerStore->Clear();
381     fTriggerStore->Connect(*clustersTree, false);
382     clustersTree->GetEvent(0);
383
384     AliMUONLocalTrigger* locTrg;
385     TIter nextLoc(fTriggerStore->CreateLocalIterator());
386
387     while ( ( locTrg = static_cast<AliMUONLocalTrigger*>(nextLoc()) ) ) 
388     {
389       if (locTrg->IsNull()) continue;
390    
391       TArrayS xyPattern[2];
392       locTrg->GetXPattern(xyPattern[0]);
393       locTrg->GetYPattern(xyPattern[1]);
394
395       Int_t nBoard = locTrg->LoCircuit();
396
397       Bool_t xTrig=locTrg->IsTrigX();
398       Bool_t yTrig=locTrg->IsTrigY();
399     
400       if (xTrig && yTrig)
401         ((TH1F*)GetRecPointsData(kTriggeredBoards))->Fill(nBoard);
402
403       fDigitMaker->TriggerDigits(nBoard, xyPattern, *fDigitStore);
404     }
405
406     TIter nextDigit(fDigitStore->CreateIterator());
407     AliMUONVDigit* mDigit;
408     while ( ( mDigit = static_cast<AliMUONVDigit*>(nextDigit()) ) )
409     {
410       Int_t detElemId = mDigit->DetElemId();
411       Int_t ch = detElemId/100;
412       Int_t localBoard = mDigit->ManuId();
413       Int_t channel = mDigit->ManuChannel();
414       Int_t cathode = mDigit->Cathode();
415       ERecPoints hindex 
416         = ( cathode == 0 ) ? kTriggerDigitsBendPlane : kTriggerDigitsNonBendPlane;
417       
418       ((TH3F*)GetRecPointsData(hindex))->Fill(ch, localBoard, channel);
419     }
420 }
421
422 //____________________________________________________________________________
423 void AliMUONQADataMakerRec::MakeESDs(AliESDEvent* esd)
424 {
425     /// make QA data from ESDs
426
427     if ( ! fIsInitESDs ) {
428       AliWarningStream() 
429         << "Skipping function due to a failure in Init" << endl;
430       return;
431     }    
432
433     TLorentzVector v1;
434
435     Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ; 
436     GetESDsData(0)->Fill(nTracks);
437
438     for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) {
439
440       AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
441       
442       // skip "ghosts"
443       if (!muonTrack->ContainTrackerData()) continue;
444       
445       muonTrack->LorentzP(v1);
446
447       GetESDsData(1)->Fill(v1.P());
448       GetESDsData(2)->Fill(v1.Pt());
449       GetESDsData(3)->Fill(v1.Rapidity());
450
451       TClonesArray clusters =  muonTrack->GetClusters();
452
453       for (Int_t iCluster = 0; iCluster <clusters.GetEntriesFast(); ++iCluster) {
454         AliESDMuonCluster* cluster = (AliESDMuonCluster*)clusters.At(iCluster);
455         GetESDsData(kESDClusterHitMap+cluster->GetChamberId())
456           ->Fill(cluster->GetX(), cluster->GetY());
457       }
458     }
459 }
460
461 //____________________________________________________________________________ 
462 void AliMUONQADataMakerRec::StartOfDetectorCycle()
463 {
464     /// Detector specific actions at start of cycle
465   
466 }
467
468 //____________________________________________________________________________ 
469 void AliMUONQADataMakerRec::DisplayTriggerInfo(AliQA::TASKINDEX_t task)
470 {
471   //
472   /// Display trigger information in a user-friendly way:
473   /// from local board and strip numbers to their position on chambers
474   //
475   if(task!=AliQA::kRECPOINTS && task!=AliQA::kRAWS) return;
476
477   AliMUONTriggerDisplay triggerDisplay;
478   
479   TH3F* histoStrips=0x0;
480   TH2F* histoDisplayStrips=0x0;
481
482   for (Int_t iCath = 0; iCath < AliMpConstants::NofCathodes(); iCath++)
483   {
484     if(task==AliQA::kRECPOINTS){
485       ERecPoints hindex 
486         = ( iCath == 0 ) ? kTriggerDigitsBendPlane : kTriggerDigitsNonBendPlane;
487       histoStrips = (TH3F*)GetRecPointsData(hindex);
488     }
489     else if(task==AliQA::kRAWS){
490       ERaw hindex 
491         = ( iCath == 0 ) ? kTriggerScalersBP : kTriggerScalersNBP;
492       histoStrips = (TH3F*)GetRawsData(hindex);
493       if(histoStrips->GetEntries()==0) return; // No scalers found
494     }
495     
496     for (Int_t iChamber = 0; iChamber < AliMpConstants::NofTriggerChambers(); iChamber++)
497     {
498       if(task==AliQA::kRECPOINTS){
499         histoDisplayStrips = (TH2F*)GetRecPointsData(kTriggerDigitsDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
500       }
501       else if(task==AliQA::kRAWS){
502         histoDisplayStrips = (TH2F*)GetRawsData(kTriggerScalersDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
503       }
504       Int_t bin = histoStrips->GetXaxis()->FindBin(11+iChamber);
505       histoStrips->GetXaxis()->SetRange(bin,bin);
506       TH2F* inputHisto = (TH2F*)histoStrips->Project3D("zy");
507       triggerDisplay.FillDisplayHistogram(inputHisto, histoDisplayStrips, AliMUONTriggerDisplay::kDisplayStrips, iCath, iChamber);
508     } // iChamber
509   } // iCath
510
511   if(task!=AliQA::kRECPOINTS) return;
512   TH1F* histoBoards = (TH1F*)GetRecPointsData(kTriggeredBoards);
513   TH2F* histoDisplayBoards = (TH2F*)GetRecPointsData(kTriggerBoardsDisplay);
514   triggerDisplay.FillDisplayHistogram(histoBoards, histoDisplayBoards, AliMUONTriggerDisplay::kDisplayBoards, 0, 0);
515 }