]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONQADataMakerRec.cxx
AliMUONBlockHeader, AliMUONRawWriter:
[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
17 // --- ROOT system ---
18 #include <TClonesArray.h>
19 #include <TFile.h> 
20 #include <TH1F.h> 
21 #include <TH1I.h> 
22 #include <TH2F.h>
23 #include <TH3F.h> 
24 #include <TLorentzVector.h>
25
26 // --- AliRoot header files ---
27 #include "AliESDEvent.h"
28 #include "AliMUONConstants.h"  
29 #include "AliLog.h"
30 #include "AliRawReader.h"
31 #include "AliQAChecker.h"
32 #include "AliMpBusPatch.h"
33 #include "AliMUONCluster.h"  
34 #include "AliMUONRawStreamTracker.h"
35 #include "AliMUONRawStreamTrigger.h"
36
37 #include "AliMUONVClusterStore.h"
38 #include "AliMUONVCluster.h"
39 #include "AliESDMuonTrack.h"
40 #include "AliESDMuonCluster.h"
41
42 #include "AliMUONDigitMaker.h"
43 #include "AliMUONVDigitStore.h"
44 #include "AliMUONVTriggerStore.h"
45 #include "AliMUONVDigit.h"
46 #include "AliMUONLocalTrigger.h"
47
48 #include "AliMUONQADataMakerRec.h"
49
50 //-----------------------------------------------------------------------------
51 /// \class AliMUONQADataMakerRec
52 ///
53 /// MUON base class for quality assurance data (histo) maker
54 ///
55 /// \author C. Finck
56
57 /// \cond CLASSIMP
58 ClassImp(AliMUONQADataMakerRec)
59 /// \endcond
60            
61 //____________________________________________________________________________ 
62 AliMUONQADataMakerRec::AliMUONQADataMakerRec() : 
63     AliQADataMakerRec(AliQA::GetDetName(AliQA::kMUON), "MUON Quality Assurance Data Maker"),
64     fDigitStore(0x0),
65     fTriggerStore(0x0),
66     fDigitMaker(0x0)
67 {
68     /// ctor
69   fDigitStore = AliMUONVDigitStore::Create("AliMUONDigitStoreV1");
70   fDigitMaker = new AliMUONDigitMaker(kTRUE,kFALSE);
71
72 }
73
74 //____________________________________________________________________________ 
75 AliMUONQADataMakerRec::AliMUONQADataMakerRec(const AliMUONQADataMakerRec& qadm) :
76     AliQADataMakerRec(qadm),
77     fDigitStore(0x0),
78     fTriggerStore(0x0),
79     fDigitMaker(0x0)
80 {
81     ///copy ctor 
82     SetName((const char*)qadm.GetName()) ; 
83     SetTitle((const char*)qadm.GetTitle()); 
84
85     // Do not copy the digit store and digit maker, but create its own ones
86     fDigitStore = AliMUONVDigitStore::Create("AliMUONDigitStoreV1");
87     fDigitMaker = new AliMUONDigitMaker(kTRUE,kFALSE);
88 }
89
90 //__________________________________________________________________
91 AliMUONQADataMakerRec& AliMUONQADataMakerRec::operator = (const AliMUONQADataMakerRec& qadm )
92 {
93   /// Assignment operator
94
95   // check assignment to self
96   if (this == &qadm) return *this;
97
98   this->~AliMUONQADataMakerRec();
99   new(this) AliMUONQADataMakerRec(qadm);
100   return *this;
101 }
102
103 //__________________________________________________________________
104 AliMUONQADataMakerRec::~AliMUONQADataMakerRec()
105 {
106     /// dtor
107     delete fDigitStore;
108     delete fTriggerStore;
109     delete fDigitMaker;
110 }
111
112 //____________________________________________________________________________ 
113 void AliMUONQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray* list)
114 {
115     ///Detector specific actions at end of cycle
116     // do the QA checking
117     AliQAChecker::Instance()->Run(AliQA::kMUON, task, list) ;  
118 }
119
120 //____________________________________________________________________________ 
121 void AliMUONQADataMakerRec::InitRaws()
122 {
123     /// create Raws histograms in Raws subdir
124     TH1I* h0 = new TH1I("hRawBusPatch", "buspatch distribution",  1932, 1, 1932); 
125     Add2RawsList(h0, kRawBusPatch);
126
127     TH1I* h1 = new TH1I("hRawCharge", "Charge distribution in rawdata", 4096, 0, 4095); 
128     Add2RawsList(h1, kRawCharge);
129                 
130                 for (Int_t iDDL = 0; iDDL < 20; ++iDDL) 
131                 {
132                         TH1F* h2 = new TH1F(Form("%s%d", "hRawBusPatchDDL", iDDL), Form("%s %d","RAW Buspatch distribution for DDL", iDDL), 50, 0, 49); 
133                         Add2RawsList(h2, kRawBuspatchDDL+iDDL);
134                 }
135
136 }
137
138 //____________________________________________________________________________ 
139 void AliMUONQADataMakerRec::InitRecPoints()
140 {
141     /// create Reconstructed Points histograms in RecPoints subdir
142     TH3F *h2 = new TH3F("hTriggerDigitsBendPlane", "Trigger digits in bending plane",
143                         4, -0.5, 4. - 0.5,
144                         18, -0.5, 18. - 0.5,
145                         7*64, -0.5, 7.*64. - 0.5);
146     Add2RecPointsList(h2, kTriggerDigitsBendPlane);
147
148     TH3F *h3 = new TH3F("hTriggerDigitsNonBendPlane", "Trigger digits in non-bending plane",
149                         4, -0.5, 4. - 0.5,
150                         18, -0.5, 18. - 0.5,
151                         112, -0.5, 112. - 0.5);
152     Add2RecPointsList(h3, kTriggerDigitsNonBendPlane);
153 }
154
155
156 //____________________________________________________________________________ 
157 void AliMUONQADataMakerRec::InitESDs()
158 {
159     ///create ESDs histograms in ESDs subdir
160   TH1F* h0 = new TH1F("hESDnTracks", "ESDs track number distribution", 30, 0., 30.);  
161   Add2ESDsList(h0, kESDnTracks);
162
163   TH1F* h1 = new TH1F("hESDMomentum", "ESDs P distribution", 300, 0., 300) ; 
164   Add2ESDsList(h1, kESDMomentum);
165
166   TH1F* h2 = new TH1F("hESDPt", "ESDs Pt distribution", 200, 0., 50) ; 
167   Add2ESDsList(h2, kESDPt);
168
169   TH1F* h3 = new TH1F("hESDRapidity", "ESDs rapidity distribution", 200, -4.5,-2.) ; 
170   Add2ESDsList(h3, kESDRapidity);
171
172   for (Int_t i = 0; i < AliMUONConstants::NTrackingCh(); ++i) 
173   {
174     TH2F* h4 = new TH2F(Form("%s%d", "hESDClusterHitMap", i), 
175                      Form("%s %d", "ESD Clusters hit distribution for chamber", i),
176                      100, -1*AliMUONConstants::Rmax(i/2), AliMUONConstants::Rmax(i/2),
177                      100, -1*AliMUONConstants::Rmax(i/2), AliMUONConstants::Rmax(i/2)); 
178     Add2ESDsList(h4, kESDClusterHitMap+i);
179   }
180 }
181
182 //____________________________________________________________________________
183 void AliMUONQADataMakerRec::MakeRaws(AliRawReader* rawReader)
184 {
185     /// make QA for rawdata
186     Int_t busPatchId;
187     UShort_t manuId;  
188     UChar_t channelId;
189     UShort_t charge;
190
191     rawReader->Reset();
192     AliMUONRawStreamTracker rawStream(rawReader);
193     rawStream.First();
194     while( rawStream.Next(busPatchId, manuId, channelId, charge) ) {
195   
196       GetRawsData(kRawBusPatch)->Fill(busPatchId);
197       GetRawsData(kRawCharge)->Fill(charge);
198                         Int_t iDDL = rawStream.GetCurentDDL();
199                         GetRawsData(kRawBuspatchDDL + iDDL)->Fill(AliMpBusPatch::GetLocalBusID(busPatchId, iDDL));
200                 
201                   
202     } // Next digit
203 }
204
205 //____________________________________________________________________________
206 void AliMUONQADataMakerRec::MakeRecPoints(TTree* clustersTree)
207 {
208   
209     /// makes data from trigger response
210       
211     // Fired pads info
212     fDigitStore->Clear();
213
214     if (!fTriggerStore) fTriggerStore = AliMUONVTriggerStore::Create(*clustersTree);
215     fTriggerStore->Clear();
216     fTriggerStore->Connect(*clustersTree, false);
217     clustersTree->GetEvent(0);
218
219     AliMUONLocalTrigger* locTrg;
220     TIter nextLoc(fTriggerStore->CreateLocalIterator());
221
222     while ( ( locTrg = static_cast<AliMUONLocalTrigger*>(nextLoc()) ) ) 
223     {
224       if (locTrg->IsNull()) continue;
225    
226       TArrayS xyPattern[2];
227       locTrg->GetXPattern(xyPattern[0]);
228       locTrg->GetYPattern(xyPattern[1]);
229
230       Int_t nBoard = locTrg->LoCircuit();
231        fDigitMaker->TriggerDigits(nBoard, xyPattern, *fDigitStore);
232     }
233
234     TIter nextDigit(fDigitStore->CreateIterator());
235     AliMUONVDigit* mDigit;
236     while ( ( mDigit = static_cast<AliMUONVDigit*>(nextDigit()) ) )
237     {
238       Int_t detElemId = mDigit->DetElemId();
239       Int_t ch = detElemId/100 - 11;
240       Int_t slat = detElemId%100;
241       Int_t cathode = mDigit->Cathode();
242       Int_t ix = mDigit->PadX();
243       Int_t iy = mDigit->PadY();
244       Int_t maxY = (cathode==0) ? 64 : 1;
245       Int_t currPair = ix*maxY + iy;
246       ERecPoints hindex 
247         = ( cathode == 0 ) ? kTriggerDigitsBendPlane : kTriggerDigitsNonBendPlane;
248       
249       ((TH3F*)GetRecPointsData(hindex))->Fill(ch, slat, currPair);
250     }
251 }
252
253 //____________________________________________________________________________
254 void AliMUONQADataMakerRec::MakeESDs(AliESDEvent* esd)
255 {
256     /// make QA data from ESDs
257
258     TLorentzVector v1;
259
260     Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ; 
261     GetESDsData(0)->Fill(nTracks);
262
263     for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) {
264
265       AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
266       muonTrack->LorentzP(v1);
267
268       GetESDsData(1)->Fill(v1.P());
269       GetESDsData(2)->Fill(v1.Pt());
270       GetESDsData(3)->Fill(v1.Rapidity());
271
272       TClonesArray clusters =  muonTrack->GetClusters();
273
274       for (Int_t iCluster = 0; iCluster <clusters.GetEntriesFast(); ++iCluster) {
275         AliESDMuonCluster* cluster = (AliESDMuonCluster*)clusters.At(iCluster);
276         GetESDsData(kESDClusterHitMap+cluster->GetChamberId())
277           ->Fill(cluster->GetX(), cluster->GetY());
278       }
279     }
280 }
281
282 //____________________________________________________________________________ 
283 void AliMUONQADataMakerRec::StartOfDetectorCycle()
284 {
285     /// Detector specific actions at start of cycle
286   
287 }