]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSQADataMakerSim.cxx
Introducing event specie in QA (Yves)
[u/mrichter/AliRoot.git] / PHOS / AliPHOSQADataMakerSim.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 /* $Id$ */
18
19 /*
20   Produces the data needed to calculate the quality assurance. 
21   All data must be mergeable objects.
22   Y. Schutz CERN July 2007
23 */
24
25 // --- ROOT system ---
26 #include <TClonesArray.h>
27 #include <TFile.h> 
28 #include <TH1F.h> 
29 #include <TH1I.h> 
30 #include <TH2F.h> 
31 #include <TTree.h>
32
33 // --- Standard library ---
34
35 // --- AliRoot header files ---
36 #include "AliESDCaloCluster.h"
37 #include "AliLog.h"
38 #include "AliPHOSDigit.h"
39 #include "AliPHOSHit.h"
40 #include "AliPHOSQADataMakerSim.h"
41 #include "AliQAChecker.h"
42
43 ClassImp(AliPHOSQADataMakerSim)
44            
45 //____________________________________________________________________________ 
46 AliPHOSQADataMakerSim::AliPHOSQADataMakerSim() : 
47   AliQADataMakerSim(AliQA::GetDetName(AliQA::kPHOS), "PHOS Quality Assurance Data Maker"),
48   fHits(0x0)
49 {
50   // ctor
51   fHits = new TClonesArray("AliPHOSHit", 1000);
52 }
53
54 //____________________________________________________________________________ 
55 AliPHOSQADataMakerSim::AliPHOSQADataMakerSim(const AliPHOSQADataMakerSim& qadm) :
56   AliQADataMakerSim(),
57   fHits(0x0)
58 {
59   //copy ctor 
60   SetName((const char*)qadm.GetName()) ; 
61   SetTitle((const char*)qadm.GetTitle()); 
62   fHits = new TClonesArray("AliPHOSHit", 1000);
63 }
64
65 //__________________________________________________________________
66 AliPHOSQADataMakerSim& AliPHOSQADataMakerSim::operator = (const AliPHOSQADataMakerSim& qadm )
67 {
68   // Assign operator.
69   this->~AliPHOSQADataMakerSim();
70   new(this) AliPHOSQADataMakerSim(qadm);
71   return *this;
72 }
73  
74 //____________________________________________________________________________ 
75 void AliPHOSQADataMakerSim::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray ** list)
76 {
77   //Detector specific actions at end of cycle
78   // do the QA checking
79   AliQAChecker::Instance()->Run(AliQA::kPHOS, task, list) ;  
80 }
81
82 //____________________________________________________________________________ 
83 void AliPHOSQADataMakerSim::InitHits()
84 {
85   // create Hits histograms in Hits subdir
86   Bool_t expert   = kTRUE ; 
87   TH1F * h0 = new TH1F("hPhosHits",    "Hits energy distribution in PHOS",       100, 0., 100.) ; 
88   h0->Sumw2() ;
89   Add2HitsList(h0, kHits, !expert) ;
90   TH1I * h1 = new TH1I("hPhosHitsMul", "Hits multiplicity distribution in PHOS", 500, 0., 10000) ; 
91   h1->Sumw2() ;
92   Add2HitsList(h1, kHitsMul, !expert) ;
93   
94 }
95
96 //____________________________________________________________________________ 
97 void AliPHOSQADataMakerSim::InitDigits()
98 {
99   // create Digits histograms in Digits subdir
100   Bool_t expert   = kTRUE ; 
101   TH1I * h0 = new TH1I("hPhosDigits",    "Digits amplitude distribution in PHOS",    500, 0, 1000) ; 
102   h0->Sumw2() ;
103   Add2DigitsList(h0, kDigits, !expert) ;
104   TH1I * h1 = new TH1I("hPhosDigitsMul", "Digits multiplicity distribution in PHOS", 2000, 0, 10000) ; 
105   h1->Sumw2() ;
106   Add2DigitsList(h1, kDigitsMul, !expert) ;
107 }
108
109 //____________________________________________________________________________ 
110 void AliPHOSQADataMakerSim::InitSDigits()
111 {
112   // create SDigits histograms in SDigits subdir
113   Bool_t expert   = kTRUE ; 
114   TH1F * h0 = new TH1F("hPhosSDigits",    "SDigits energy distribution in PHOS",       500, 0., 1000.) ; 
115   h0->Sumw2() ;
116   Add2SDigitsList(h0, kSDigits, !expert) ;
117   TH1I * h1 = new TH1I("hPhosSDigitsMul", "SDigits multiplicity distribution in PHOS", 500, 0,  1000) ; 
118   h1->Sumw2() ;
119   Add2SDigitsList(h1, kSDigitsMul, !expert) ;
120 }
121
122 //____________________________________________________________________________
123 void AliPHOSQADataMakerSim::MakeHits()
124 {
125   //make QA data from Hits
126   
127   TIter next(fHits) ; 
128   AliPHOSHit * hit ; 
129   while ( (hit = dynamic_cast<AliPHOSHit *>(next())) ) {
130     GetHitsData(kHits)->Fill( hit->GetEnergy()) ;
131   }
132 }
133
134 //____________________________________________________________________________
135 void AliPHOSQADataMakerSim::MakeHits(TTree * hitTree)
136 {
137   // make QA data from Hit Tree
138   
139   TBranch * branch = hitTree->GetBranch("PHOS") ;
140   if ( ! branch ) {
141     AliWarning("PHOS branch in Hit Tree not found") ; 
142   } else {
143     Int_t nHits = 0;
144     branch->SetAddress(&fHits) ;
145     for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
146       branch->GetEntry(ientry) ;
147       nHits += fHits->GetEntriesFast();
148       MakeHits() ; 
149       fHits->Clear();
150     }   
151     GetHitsData(1)->Fill(nHits) ;
152   }
153 }
154
155 //____________________________________________________________________________
156 void AliPHOSQADataMakerSim::MakeDigits(TClonesArray * digits)
157 {
158   // makes data from Digits
159
160     GetDigitsData(1)->Fill(digits->GetEntriesFast()) ; 
161     TIter next(digits) ; 
162     AliPHOSDigit * digit ; 
163     while ( (digit = dynamic_cast<AliPHOSDigit *>(next())) ) {
164       GetDigitsData(kDigits)->Fill( digit->GetEnergy()) ;
165     }  
166 }
167
168 //____________________________________________________________________________
169 void AliPHOSQADataMakerSim::MakeDigits(TTree * digitTree)
170 {
171         // makes data from Digit Tree
172         TClonesArray * digits = new TClonesArray("AliPHOSDigit", 1000) ; 
173
174         TBranch * branch = digitTree->GetBranch("PHOS") ;
175         if ( ! branch ) {
176                 AliWarning("PHOS branch in Digit Tree not found") ; 
177         } else {
178                 branch->SetAddress(&digits) ;
179                 branch->GetEntry(0) ; 
180                 MakeDigits(digits) ; 
181         }
182 }
183
184 //____________________________________________________________________________
185 void AliPHOSQADataMakerSim::MakeSDigits(TClonesArray * sdigits)
186 {
187   // makes data from SDigits
188   
189         GetSDigitsData(1)->Fill(sdigits->GetEntriesFast()) ; 
190     TIter next(sdigits) ; 
191     AliPHOSDigit * sdigit ; 
192     while ( (sdigit = dynamic_cast<AliPHOSDigit *>(next())) ) {
193       GetSDigitsData(kSDigits)->Fill( sdigit->GetEnergy()) ;
194     } 
195 }
196
197 //____________________________________________________________________________
198 void AliPHOSQADataMakerSim::MakeSDigits(TTree * sdigitTree)
199 {
200         // makes data from SDigit Tree
201         TClonesArray * sdigits = new TClonesArray("AliPHOSDigit", 1000) ; 
202
203         TBranch * branch = sdigitTree->GetBranch("PHOS") ;
204         if ( ! branch ) {
205                 AliWarning("PHOS branch in SDigit Tree not found") ; 
206         } else {
207                 branch->SetAddress(&sdigits) ;
208                 branch->GetEntry(0) ;
209                 MakeSDigits(sdigits) ; 
210         }
211 }
212
213 //____________________________________________________________________________ 
214 void AliPHOSQADataMakerSim::StartOfDetectorCycle()
215 {
216   //Detector specific actions at start of cycle
217   
218 }