]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSQADataMakerSim.cxx
Updated histogram limits (PHOS energy)
[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, 0, !expert) ;
90   TH1I * h1 = new TH1I("hPhosHitsMul", "Hits multiplicity distribution in PHOS", 500, 0., 10000) ; 
91   h1->Sumw2() ;
92   Add2HitsList(h1, 1, !expert) ;
93 }
94
95 //____________________________________________________________________________ 
96 void AliPHOSQADataMakerSim::InitDigits()
97 {
98   // create Digits histograms in Digits subdir
99   Bool_t expert   = kTRUE ; 
100   TH1I * h0 = new TH1I("hPhosDigits",    "Digits amplitude distribution in PHOS",    500, 0, 5000) ; 
101   h0->Sumw2() ;
102   Add2DigitsList(h0, 0, !expert) ;
103   TH1I * h1 = new TH1I("hPhosDigitsMul", "Digits multiplicity distribution in PHOS", 500, 0, 1000) ; 
104   h1->Sumw2() ;
105   Add2DigitsList(h1, 1, !expert) ;
106 }
107
108 //____________________________________________________________________________ 
109 void AliPHOSQADataMakerSim::InitSDigits()
110 {
111   // create SDigits histograms in SDigits subdir
112   Bool_t expert   = kTRUE ; 
113   TH1F * h0 = new TH1F("hPhosSDigits",    "SDigits energy distribution in PHOS",       100, 0., 100.) ; 
114   h0->Sumw2() ;
115   Add2SDigitsList(h0, 0, !expert) ;
116   TH1I * h1 = new TH1I("hPhosSDigitsMul", "SDigits multiplicity distribution in PHOS", 500, 0,  10000) ; 
117   h1->Sumw2() ;
118   Add2SDigitsList(h1, 1, !expert) ;
119 }
120
121 //____________________________________________________________________________
122 void AliPHOSQADataMakerSim::MakeHits()
123 {
124   //make QA data from Hits
125   
126   TIter next(fHits) ; 
127   AliPHOSHit * hit ; 
128   while ( (hit = dynamic_cast<AliPHOSHit *>(next())) ) {
129     GetHitsData(0)->Fill( hit->GetEnergy()) ;
130   }
131 }
132
133 //____________________________________________________________________________
134 void AliPHOSQADataMakerSim::MakeHits(TTree * hitTree)
135 {
136   // make QA data from Hit Tree
137   
138   TBranch * branch = hitTree->GetBranch("PHOS") ;
139   if ( ! branch ) {
140     AliWarning("PHOS branch in Hit Tree not found") ; 
141   } else {
142     Int_t nHits = 0;
143     branch->SetAddress(&fHits) ;
144     for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
145       branch->GetEntry(ientry) ;
146       nHits += fHits->GetEntriesFast();
147       MakeHits() ; 
148       fHits->Clear();
149     }   
150     GetHitsData(1)->Fill(nHits) ;
151   }
152 }
153
154 //____________________________________________________________________________
155 void AliPHOSQADataMakerSim::MakeDigits(TClonesArray * digits)
156 {
157   // makes data from Digits
158
159     GetDigitsData(1)->Fill(digits->GetEntriesFast()) ; 
160     TIter next(digits) ; 
161     AliPHOSDigit * digit ; 
162     while ( (digit = dynamic_cast<AliPHOSDigit *>(next())) ) {
163       GetDigitsData(0)->Fill( digit->GetEnergy()) ;
164     }  
165 }
166
167 //____________________________________________________________________________
168 void AliPHOSQADataMakerSim::MakeDigits(TTree * digitTree)
169 {
170         // makes data from Digit Tree
171         TClonesArray * digits = new TClonesArray("AliPHOSDigit", 1000) ; 
172
173         TBranch * branch = digitTree->GetBranch("PHOS") ;
174         if ( ! branch ) {
175                 AliWarning("PHOS branch in Digit Tree not found") ; 
176         } else {
177                 branch->SetAddress(&digits) ;
178                 branch->GetEntry(0) ; 
179                 MakeDigits(digits) ; 
180         }
181 }
182
183 //____________________________________________________________________________
184 void AliPHOSQADataMakerSim::MakeSDigits(TClonesArray * sdigits)
185 {
186   // makes data from SDigits
187   
188         GetSDigitsData(1)->Fill(sdigits->GetEntriesFast()) ; 
189     TIter next(sdigits) ; 
190     AliPHOSDigit * sdigit ; 
191     while ( (sdigit = dynamic_cast<AliPHOSDigit *>(next())) ) {
192       GetSDigitsData(0)->Fill( sdigit->GetEnergy()) ;
193     } 
194 }
195
196 //____________________________________________________________________________
197 void AliPHOSQADataMakerSim::MakeSDigits(TTree * sdigitTree)
198 {
199         // makes data from SDigit Tree
200         TClonesArray * sdigits = new TClonesArray("AliPHOSDigit", 1000) ; 
201
202         TBranch * branch = sdigitTree->GetBranch("PHOS") ;
203         if ( ! branch ) {
204                 AliWarning("PHOS branch in SDigit Tree not found") ; 
205         } else {
206                 branch->SetAddress(&sdigits) ;
207                 branch->GetEntry(0) ;
208                 MakeSDigits(sdigits) ; 
209         }
210 }
211
212 //____________________________________________________________________________ 
213 void AliPHOSQADataMakerSim::StartOfDetectorCycle()
214 {
215   //Detector specific actions at start of cycle
216   
217 }