]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ACORDE/AliACORDEQADataMakerRec.cxx
added verbosity to QA histograms (Yves)
[u/mrichter/AliRoot.git] / ACORDE / AliACORDEQADataMakerRec.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 //  Produces the data needed to calculate the quality assurance. 
17 //  All data must be mergeable objects.
18
19
20 //  Authors:
21 //
22 //  Luciano Diaz Gonzalez <luciano.diaz@nucleares.unam.mx> (ICN-UNAM)
23 //  Mario Rodriguez Cahuantzi <mrodrigu@mail.cern.ch> (FCFM-BUAP)
24 //  Arturo Fernandez Tellez <afernan@mail.cern.ch (FCFM-BUAP)
25 //
26 //  Created: June 13th 2008
27 //---
28 // Last Update: Aug. 27th 2008 --> Implementation to declare QA expert histogram
29
30
31 // --- ROOT system ---
32 #include <TClonesArray.h>
33 #include <TFile.h> 
34 #include <TH1F.h> 
35 #include <TDirectory.h>
36 // --- Standard library ---
37
38 // --- AliRoot header files ---
39 #include "AliESDEvent.h"
40 #include "AliLog.h"
41 #include "AliACORDEdigit.h" 
42 #include "AliACORDEhit.h"
43 #include "AliACORDEQADataMakerRec.h"
44 #include "AliQAChecker.h"
45 #include "AliACORDERawReader.h"
46 #include "AliACORDERawStream.h"
47 ClassImp(AliACORDEQADataMakerRec)
48            
49 //____________________________________________________________________________ 
50 AliACORDEQADataMakerRec::AliACORDEQADataMakerRec():AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kACORDE), "ACORDE Quality Assurance Data Maker")
51 {
52
53 }
54 //____________________________________________________________________________ 
55 AliACORDEQADataMakerRec::AliACORDEQADataMakerRec(const AliACORDEQADataMakerRec& qadm):AliQADataMakerRec() 
56 {
57   SetName((const char*)qadm.GetName()) ; 
58   SetTitle((const char*)qadm.GetTitle()); 
59 }
60 //__________________________________________________________________
61 AliACORDEQADataMakerRec& AliACORDEQADataMakerRec::operator = (const AliACORDEQADataMakerRec& qadm )
62 {
63   // Equal operator.
64   this->~AliACORDEQADataMakerRec();
65   new(this) AliACORDEQADataMakerRec(qadm);
66   return *this;
67 }
68 //____________________________________________________________________________
69 void AliACORDEQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
70 {
71   //Detector specific actions at end of cycle
72   // do the QA checking
73   AliQAChecker::Instance()->Run(AliQAv1::kACORDE, task, list) ;
74 }
75
76 //____________________________________________________________________________
77 void AliACORDEQADataMakerRec::StartOfDetectorCycle()
78 {
79   //Detector specific actions at start of cycle
80
81 }
82  
83 //____________________________________________________________________________ 
84 void AliACORDEQADataMakerRec::InitRaws()
85 {
86   // create Raw histograms in Raw subdir
87
88   const Bool_t expert   = kTRUE ; 
89   const Bool_t saveCorr = kTRUE ; 
90   const Bool_t image    = kTRUE ; 
91   
92   TH1D *fhACORDEBitPattern[4];
93   fhACORDEBitPattern[0] = new TH1D("ACORDERawDataSM","ACORDE-SingleMuon;Bit Pattern;Counts",60,1,60);//AcordeSingleMuon BitPattern
94   fhACORDEBitPattern[1] = new TH1D("ACORDERawDataMM","ACORDE-MultiMuon;Bit Pattern;Counts",60,1,60);//AcordeMultiMuon BitPattern
95   fhACORDEBitPattern[2] = new TH1D("ACORDERawDataSMM","ACORDE-SingleMuonMultiplicity;Multiplicity;Counts",60,1,60);//AcordeSingleMuon Multiplicity
96   fhACORDEBitPattern[3] = new TH1D("ACORDERawDataMMM","ACORDE-MultiMuonMultiplicity;Multiplicity;Counts",60,1,60);//AcordeMultiMuon Multiplicity
97   for(Int_t i=0;i<4;i++) 
98     Add2RawsList(fhACORDEBitPattern[i],i,!expert, image, !saveCorr);
99 }
100 //____________________________________________________________________________ 
101 void AliACORDEQADataMakerRec::InitDigits()
102 {
103   // create Digits histograms in Digits subdir
104   
105   const Bool_t expert   = kTRUE ; 
106   const Bool_t image    = kTRUE ; 
107   TH1F *    fhDigitsModule;
108   TString   modulename;
109   modulename = "hDigitsModule";
110   fhDigitsModule = new TH1F(modulename.Data(),"hDigitsModuleSingle;# of modules;Counts",60,0,60);
111   Add2DigitsList(fhDigitsModule,0,!expert,image);
112     
113 }
114
115 //____________________________________________________________________________ 
116
117 void AliACORDEQADataMakerRec::InitRecPoints()
118 {
119   // create cluster histograms in RecPoint subdir
120   // Not needed for ACORDE by now !!!
121 }
122
123 //____________________________________________________________________________
124 void AliACORDEQADataMakerRec::InitESDs()
125 {
126   //create ESDs histograms in ESDs subdir
127
128   const Bool_t expert   = kTRUE ; 
129   const Bool_t image    = kTRUE ; 
130   
131   TH1F *    fhESDsSingle;
132   TH1F *    fhESDsMulti;
133
134    TString   name;
135
136    name = "hESDsSingle";
137    fhESDsSingle = new TH1F(name.Data(),"hESDsSingle;??;??",60,0,60);
138    Add2ESDsList(fhESDsSingle,0,!expert,image);
139
140    name = "hESDsMulti";
141    fhESDsMulti = new TH1F(name.Data(),"hESDsMulti;??;??",60,0,60);
142    Add2ESDsList(fhESDsMulti,1,!expert,image);
143 }
144 //____________________________________________________________________________
145 void AliACORDEQADataMakerRec::MakeRaws(AliRawReader* rawReader)
146 {
147   //fills QA histos for RAW
148   rawReader->Reset();
149   AliACORDERawStream rawStream(rawReader);
150   size_t contSingle=0;
151   size_t contMulti=0;
152   UInt_t dy[4];
153
154   bool kroSingle[60],kroMulti[60];
155   UInt_t tmpDy;
156
157   for(Int_t m=0;m<60;m++) {kroSingle[m]=0;kroMulti[m]=0;}
158
159 if(rawStream.Next())
160 {
161         dy[0]=rawStream.GetWord(0);
162         dy[1]=rawStream.GetWord(1);
163         dy[2]=rawStream.GetWord(2);
164         dy[3]=rawStream.GetWord(3);
165         tmpDy=dy[0];
166         for(Int_t r=0;r<30;++r)
167         {
168                 kroSingle[r] = tmpDy & 1;
169                 tmpDy>>=1;
170         }
171         tmpDy=dy[1];
172         for(Int_t r=30;r<60;++r)
173         {
174                 kroSingle[r] = tmpDy & 1;
175                 tmpDy>>=1;
176         }
177         tmpDy=dy[2];
178         for(Int_t r=0;r<30;++r)
179         {
180                 kroMulti[r] = tmpDy & 1;
181                 tmpDy>>=1;
182         }
183         tmpDy=dy[3];
184         for(Int_t r=30;r<60;++r)
185         {
186                 kroMulti[r] = tmpDy & 1;
187                 tmpDy>>=1;
188         }
189         contSingle=0;
190         contMulti=0;
191         for(Int_t r=0;r<60;++r)
192         {
193                 if(kroSingle[r]==1)
194                 {
195                         GetRawsData(0)->Fill(r+1);
196                         contSingle=contSingle+1;
197                 }
198                 if(kroMulti[r]==1)
199                 {
200                         GetRawsData(1)->Fill(r+1);
201                         contMulti++;
202                 }
203
204         }GetRawsData(2)->Fill(contSingle);GetRawsData(3)->Fill(contMulti);
205 }
206 }
207 //____________________________________________________________________________
208 void AliACORDEQADataMakerRec::MakeDigits( TTree *digitsTree)
209 {
210   //fills QA histos for Digits
211   TClonesArray * digits = new TClonesArray("AliACORDEdigit",1000);
212   TBranch * branch = digitsTree->GetBranch("ACORDEdigit");
213   if (!branch) {
214     AliWarning("ACORDE branch in Digits Tree not found");
215   } else {
216     branch->SetAddress(&digits);
217     for(Int_t track = 0 ; track < branch->GetEntries() ; track++) {
218       branch->GetEntry(track);
219       for(Int_t idigit = 0 ; idigit < digits->GetEntriesFast() ; idigit++) {
220         AliACORDEdigit *AcoDigit = (AliACORDEdigit*) digits->UncheckedAt(idigit);
221         if (!AcoDigit) {
222           AliError("The unchecked digit doesn't exist");
223           continue ;
224         }
225         GetDigitsData(0)->Fill(AcoDigit->GetModule()-1);
226       }
227     }
228   }
229 }
230
231 //____________________________________________________________________________
232 void AliACORDEQADataMakerRec::MakeESDs(AliESDEvent * esd)
233 {
234   //fills QA histos for ESD
235
236         AliESDACORDE * fESDACORDE= esd->GetACORDEData();
237         Int_t *fACORDEMultiMuon =fESDACORDE->GetACORDEMultiMuon();
238         Int_t *fACORDESingleMuon=fESDACORDE->GetACORDESingleMuon();
239
240        for(int i=0;i<60;i++){
241             if(fACORDESingleMuon[i]==1)
242                GetESDsData(0) -> Fill(i);
243             if(fACORDEMultiMuon[i]==1)
244                GetESDsData(1) -> Fill(i);
245         }
246
247
248
249
250 }
251