]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALQADataMakerRec.cxx
correct coding violations
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALQADataMakerRec.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 Based on the QA code for PHOS written by Yves Schutz July 2007
17
18 Authors:  J.Klay (Cal Poly) May 2008
19           S. Salur LBL April 2008
20
21 */
22
23 // --- ROOT system ---
24 #include <TClonesArray.h>
25 #include <TFile.h> 
26 #include <TH1F.h> 
27 #include <TH1I.h> 
28 #include <TH2F.h> 
29
30 // --- Standard library ---
31
32 // --- AliRoot header files ---
33 #include "AliESDCaloCluster.h"
34 #include "AliESDCaloCells.h"
35 #include "AliESDEvent.h"
36 #include "AliLog.h"
37 #include "AliEMCALQADataMakerRec.h"
38 #include "AliQAChecker.h"
39 #include "AliEMCALRecPoint.h" 
40 #include "AliEMCALRawUtils.h"
41 #include "AliEMCALReconstructor.h"
42 #include "AliEMCALRecParam.h"
43 #include "AliRawReader.h"
44
45 ClassImp(AliEMCALQADataMakerRec)
46            
47 //____________________________________________________________________________ 
48   AliEMCALQADataMakerRec::AliEMCALQADataMakerRec() : 
49   AliQADataMakerRec(AliQA::GetDetName(AliQA::kEMCAL), "EMCAL Quality Assurance Data Maker")
50 {
51   // ctor
52 }
53
54 //____________________________________________________________________________ 
55 AliEMCALQADataMakerRec::AliEMCALQADataMakerRec(const AliEMCALQADataMakerRec& qadm) :
56   AliQADataMakerRec()
57 {
58   //copy ctor 
59   SetName((const char*)qadm.GetName()) ; 
60   SetTitle((const char*)qadm.GetTitle()); 
61 }
62
63 //__________________________________________________________________
64 AliEMCALQADataMakerRec& AliEMCALQADataMakerRec::operator = (const AliEMCALQADataMakerRec& qadm )
65 {
66   // Equal operator.
67   this->~AliEMCALQADataMakerRec();
68   new(this) AliEMCALQADataMakerRec(qadm);
69   return *this;
70 }
71  
72 //____________________________________________________________________________ 
73 void AliEMCALQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
74 {
75   //Detector specific actions at end of cycle
76   // do the QA checking
77   AliQAChecker::Instance()->Run(AliQA::kEMCAL, task, list) ;  
78 }
79
80 //____________________________________________________________________________ 
81 void AliEMCALQADataMakerRec::InitESDs()
82 {
83   //Create histograms to controll ESD
84   TH1F * h1 = new TH1F("hESDCaloClusterE",  "ESDs CaloCluster energy in EMCAL",    200, 0., 20.) ; 
85   h1->Sumw2() ;
86   Add2ESDsList(h1, kESDCaloClusE)  ;                                                     
87
88   TH1I * h2 = new TH1I("hESDCaloClusterM", "ESDs CaloCluster multiplicity in EMCAL", 100, 0,  100) ; 
89   h2->Sumw2() ;
90   Add2ESDsList(h2, kESDCaloClusM)  ;
91
92   TH1F * h3 = new TH1F("hESDCaloCellA",  "ESDs CaloCell amplitude in EMCAL",    500, 0., 250.) ; 
93   h3->Sumw2() ;
94   Add2ESDsList(h3, kESDCaloCellA)  ;  
95  
96   TH1I * h4 = new TH1I("hESDCaloCellM", "ESDs CaloCell multiplicity in EMCAL", 200, 0,  1000) ; 
97   h4->Sumw2() ;
98   Add2ESDsList(h4, kESDCaloCellM) ;
99         
100 }
101
102 //____________________________________________________________________________ 
103 void AliEMCALQADataMakerRec::InitRecPoints()
104 {
105   // create Reconstructed Points histograms in RecPoints subdir
106   TH1F* h0 = new TH1F("hEMCALRpE","EMCAL RecPoint energies",200, 0.,20.); //GeV
107   h0->Sumw2();
108   Add2RecPointsList(h0,kRecPE);
109
110   TH1I* h1 = new TH1I("hEMCALRpM","EMCAL RecPoint multiplicities",100,0,100);
111   h1->Sumw2();
112   Add2RecPointsList(h1,kRecPM);
113
114   TH1I* h2 = new TH1I("hEMCALRpDigM","EMCAL RecPoint Digit Multiplicities",20,0,20);
115   h2->Sumw2();
116   Add2RecPointsList(h2,kRecPDigM);
117
118 }
119
120 //____________________________________________________________________________ 
121 void AliEMCALQADataMakerRec::InitRaws()
122 {
123   // create Raws histograms in Raws subdir
124   //these need more thought
125   /*
126   TH1I * h0 = new TH1I("hLowEmcalSupermodules",    "Low Gain digits in EMCAL supermodules",       12, 0, 12) ;
127   h0->Sumw2() ;
128   Add2RawsList(h0, kNsmodLG) ;
129   TH1I * h1 = new TH1I("hHighEmcalSupermodules",   "High Gain Digits in EMCAL supermodules",       12, 0, 12) ;
130   h1->Sumw2() ;
131   Add2RawsList(h1, kNsmodHG) ;
132
133   TH1F * h2 = new TH1F("hLowEmcalRawtime", "Low Gain Time of raw digits in EMCAL", 500, -50., 200.) ;
134   h2->Sumw2() ;
135   Add2RawsList(h2, kLGtime) ;
136   TH1F * h3 = new TH1F("hHighEmcalRawtime", "High Gain Time of raw digits in EMCAL", 500, -50., 200.) ;
137   h3->Sumw2() ;
138   Add2RawsList(h3, kHGtime) ;
139
140   TH1F * h4 = new TH1F("hLowEmcalRawEnergy", "Low Gain Energy of raw digits in EMCAL", 500, 0., 1000.) ;
141   h4->Sumw2() ;
142   Add2RawsList(h4, kSpecLG) ;
143   TH1F * h5 = new TH1F("hHighEmcalRawEnergy", "High Gain Energy of raw digits in EMCAL",500,0., 1000.) ;
144   h5->Sumw2() ;
145   Add2RawsList(h5, kSpecHG) ;
146
147   TH1I * h6 = new TH1I("hLowNtot", "Low Gain Total Number of raw digits in EMCAL", 500, 0, 10000) ;
148   h6->Sumw2() ;
149   Add2RawsList(h6, kNtotLG) ;
150   TH1I * h7 = new TH1I("hHighNtot", "High Gain Total Number of raw digits in EMCAL",500,0, 10000) ;
151   h7->Sumw2() ;
152   Add2RawsList(h7, kNtotHG) ;
153
154   TH1F * h8 = new TH1F("hLowEtot", "Low Gain Total Energy of raw digits in EMCAL", 500, 0., 5000.) ;
155   h8->Sumw2() ;
156   Add2RawsList(h8, kEtotLG) ;
157   TH1F * h9 = new TH1F("hHighEtot", "High Gain Total Energy of raw digits in EMCAL",500,0., 100000.) ;
158   h9->Sumw2() ;
159   Add2RawsList(h9, kEtotHG) ;
160   */
161   
162 }
163
164 //____________________________________________________________________________
165 void AliEMCALQADataMakerRec::MakeESDs(AliESDEvent * esd)
166 {
167   // make QA data from ESDs
168
169   Int_t nTot = 0 ; 
170   for ( Int_t index = 0; index < esd->GetNumberOfCaloClusters() ; index++ ) {
171     AliESDCaloCluster * clu = esd->GetCaloCluster(index) ;
172     if( clu->IsEMCAL() ) {
173       GetESDsData(kESDCaloClusE)->Fill(clu->E()) ;
174       nTot++ ;
175     } 
176   }
177   GetESDsData(kESDCaloClusM)->Fill(nTot) ;
178
179   //fill calo cells
180   AliESDCaloCells* cells = esd->GetEMCALCells();
181   GetESDsData(kESDCaloCellM)->Fill(cells->GetNumberOfCells()) ;
182
183   for ( Int_t index = 0; index < cells->GetNumberOfCells() ; index++ ) {
184     GetESDsData(kESDCaloCellA)->Fill(cells->GetAmplitude(index)) ;
185   }
186
187 }
188
189 //____________________________________________________________________________
190 void AliEMCALQADataMakerRec::MakeRaws(AliRawReader* rawReader)
191 {
192   //Fill prepared histograms with Raw digit properties
193
194   //Raw histogram filling not yet implemented
195   //
196   //Need to figure out how to get the info we want without having to
197   //actually run Raw2Digits twice.
198   //I suspect what we actually want is a raw digits method, not a true
199   //emcal raw data method, but this doesn't seem to be allowed in
200   //AliQADataMakerRec.h
201
202         rawReader->Reset() ;
203 }
204
205 //____________________________________________________________________________
206 void AliEMCALQADataMakerRec::MakeRecPoints(TTree * clustersTree)
207 {
208   // makes data from RecPoints
209   TBranch *emcbranch = clustersTree->GetBranch("EMCALECARP");
210   if (!emcbranch) { 
211     AliError("can't get the branch with the EMCAL clusters !");
212     return;
213   }
214   TObjArray * emcrecpoints = new TObjArray(100) ;
215   emcbranch->SetAddress(&emcrecpoints);
216   emcbranch->GetEntry(0);
217   
218   GetRecPointsData(kRecPM)->Fill(emcrecpoints->GetEntriesFast()) ; 
219   TIter next(emcrecpoints) ; 
220   AliEMCALRecPoint * rp ; 
221   while ( (rp = dynamic_cast<AliEMCALRecPoint *>(next())) ) {
222     GetRecPointsData(kRecPE)->Fill( rp->GetEnergy()) ;
223     GetRecPointsData(kRecPDigM)->Fill(rp->GetMultiplicity());
224   }
225   emcrecpoints->Delete();
226   delete emcrecpoints;
227   
228 }
229
230 //____________________________________________________________________________ 
231 void AliEMCALQADataMakerRec::StartOfDetectorCycle()
232 {
233   //Detector specific actions at start of cycle
234   
235 }