]> git.uio.no Git - u/mrichter/AliRoot.git/blame - MUON/AliMUONQADataMakerRec.cxx
Splitting of the QA maker into simulation and reconstruction dependent parts (Yves)
[u/mrichter/AliRoot.git] / MUON / AliMUONQADataMakerRec.cxx
CommitLineData
04236e67 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 <TLorentzVector.h>
24
25// --- AliRoot header files ---
26#include "AliESDEvent.h"
27#include "AliLog.h"
28#include "AliRawReader.h"
29#include "AliQAChecker.h"
30#include "AliMUONCluster.h"
31#include "AliMUONRawStreamTracker.h"
32#include "AliMUONRawStreamTrigger.h"
33
34#include "AliMUONVClusterStore.h"
35#include "AliMUONVCluster.h"
36#include "AliESDMuonTrack.h"
37
38#include "AliMUONQADataMakerRec.h"
39
40//-----------------------------------------------------------------------------
41/// \class AliMUONQADataMakerRec
42///
43/// MUON base class for quality assurance data (histo) maker
44///
45/// \author C. Finck
46
47/// \cond CLASSIMP
48ClassImp(AliMUONQADataMakerRec)
49/// \endcond
50
51//____________________________________________________________________________
52AliMUONQADataMakerRec::AliMUONQADataMakerRec() :
53 AliQADataMakerRec(AliQA::GetDetName(AliQA::kMUON), "MUON Quality Assurance Data Maker"),
54 fClusterStore(0x0)
55{
56 /// ctor
57}
58
59//____________________________________________________________________________
60AliMUONQADataMakerRec::AliMUONQADataMakerRec(const AliMUONQADataMakerRec& qadm) :
61 AliQADataMakerRec()
62{
63 ///copy ctor
64 SetName((const char*)qadm.GetName()) ;
65 SetTitle((const char*)qadm.GetTitle());
66}
67
68//__________________________________________________________________
69AliMUONQADataMakerRec& AliMUONQADataMakerRec::operator = (const AliMUONQADataMakerRec& qadm )
70{
71 /// Equal operator.
72 this->~AliMUONQADataMakerRec();
73 new(this) AliMUONQADataMakerRec(qadm);
74 return *this;
75}
76
77//__________________________________________________________________
78AliMUONQADataMakerRec::~AliMUONQADataMakerRec()
79{
80 /// dtor
81 delete fClusterStore;
82}
83
84//____________________________________________________________________________
85void AliMUONQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX task, TObjArray* list)
86{
87 ///Detector specific actions at end of cycle
88 // do the QA checking
89 AliQAChecker::Instance()->Run(AliQA::kMUON, task, list) ;
90}
91
92//____________________________________________________________________________
93void AliMUONQADataMakerRec::InitRaws()
94{
95 /// create Raws histograms in Raws subdir
96 TH1I* h0 = new TH1I("hRawBusPatch", "buspatch distribution", 1932, 1, 1932);
97 Add2RawsList(h0, 0);
98
99 TH1I* h1 = new TH1I("hRawCharge", "Charge distribution in rawdata", 4096, 0, 4095);
100 Add2RawsList(h1, 1);
101
102}
103
104//____________________________________________________________________________
105void AliMUONQADataMakerRec::InitRecPoints()
106{
107 /// create Reconstructed Points histograms in RecPoints subdir
108 TH1F* h0 = new TH1F("hClusterCharge", "Clusters Charge distribution", 1000, 0., 4095.);
109 Add2RecPointsList(h0, 0);
110
111 TH1I* h1 = new TH1I("hClusterDetElem", "DetElemId distribution in Clusters ", 1000, 100., 1100.);
112 Add2RecPointsList(h1, 1);
113}
114
115
116//____________________________________________________________________________
117void AliMUONQADataMakerRec::InitESDs()
118{
119 ///create ESDs histograms in ESDs subdir
120 TH1F* h0 = new TH1F("hESDnTracks", "ESDs track number distribution", 30, 0., 30.);
121 Add2ESDsList(h0, 0);
122
123 TH1F* h1 = new TH1F("hESDMomentum", "ESDs P distribution", 300, 0., 300) ;
124 Add2ESDsList(h1, 1);
125
126 TH1F* h2 = new TH1F("hESDPt", "ESDs Pt distribution", 200, 0., 50) ;
127 Add2ESDsList(h2, 2) ;
128
129 TH1F* h3 = new TH1F("hESDRapidity", "ESDs rapidity distribution", 200, -4.5,-2.) ;
130 Add2ESDsList(h3, 3) ;
131}
132
133//____________________________________________________________________________
134void AliMUONQADataMakerRec::MakeRaws(AliRawReader* rawReader)
135{
136 /// make QA for rawdata
137 Int_t busPatchId;
138 UShort_t manuId;
139 UChar_t channelId;
140 UShort_t charge;
141
142 rawReader->Reset();
143 AliMUONRawStreamTracker rawStream(rawReader);
144 rawStream.First();
145 while( rawStream.Next(busPatchId, manuId, channelId, charge) ) {
146
147 GetRawsData(0)->Fill(busPatchId);
148 GetRawsData(1)->Fill(charge);
149
150 } // Next digit
151
152}
153
154//____________________________________________________________________________
155void AliMUONQADataMakerRec::MakeRecPoints(TTree* clustersTree)
156{
157
158 /// makes data from RecPoints
159 if (!fClusterStore)
160 fClusterStore = AliMUONVClusterStore::Create(*clustersTree);
161 fClusterStore->Connect(*clustersTree, false);
162 clustersTree->GetEvent(0);
163
164 TIter next(fClusterStore->CreateIterator());
165
166 AliMUONVCluster* clus = 0x0;
167
168 while ( ( clus = static_cast<AliMUONVCluster*>(next()) ) )
169 {
170 GetRecPointsData(0)->Fill(clus->GetCharge());
171 GetRecPointsData(1)->Fill(clus->GetDetElemId());
172 }
173}
174
175//____________________________________________________________________________
176void AliMUONQADataMakerRec::MakeESDs(AliESDEvent* esd)
177{
178 /// make QA data from ESDs
179
180 TLorentzVector v1;
181
182 Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ;
183 GetESDsData(0)->Fill(nTracks);
184
185 for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) {
186
187 AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
188 muonTrack->LorentzP(v1);
189
190 GetESDsData(1)->Fill(v1.P());
191 GetESDsData(2)->Fill(v1.Pt());
192 GetESDsData(3)->Fill(v1.Rapidity());
193 }
194}
195
196//____________________________________________________________________________
197void AliMUONQADataMakerRec::StartOfDetectorCycle()
198{
199 /// Detector specific actions at start of cycle
200
201}