]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ESDCheck/AliFMDQATask.cxx
Add AliAnalysisGoodies
[u/mrichter/AliRoot.git] / ESDCheck / AliFMDQATask.cxx
CommitLineData
1dfe075f 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// An analysis task to check the FMD data in simulated data
17//
18//*-- Hans Hjersing Dalsgaard
19//////////////////////////////////////////////////////////////////////////////
20
21#include <TChain.h>
22#include <TFile.h>
23#include <TCanvas.h>
24#include <TH1D.h>
25#include <TF1.h>
26
27#include "AliFMDQATask.h"
28#include "AliESD.h"
29#include "AliLog.h"
30
31//______________________________________________________________________________
32AliFMDQATask::AliFMDQATask(const char *name) :
33 AliAnalysisTask(name,""),
34 fChain(0),
35 fESD(0),
36 fhFMD1i(0),
37 fhFMD2i(0),
38 fhFMD2o(0),
39 fhFMD3i(0),
40 fhFMD3o(0)
41{
42 // Constructor.
43 // Input slot #0 works with an Ntuple
44 DefineInput(0, TChain::Class());
45 // Output slot #0 writes into a TH1 container
46 DefineOutput(0, TObjArray::Class()) ;
47}
48
49//______________________________________________________________________________
50AliFMDQATask::~AliFMDQATask()
51{
52 // dtor
53
54 fOutputContainer->Clear() ;
55 delete fOutputContainer ;
56
57 delete fhFMD1i ;
58 delete fhFMD2i ;
59 delete fhFMD2o ;
60 delete fhFMD3i ;
61 delete fhFMD3o ;
62
63}
64
65//______________________________________________________________________________
66void AliFMDQATask::Init(const Option_t*)
67{
68 // Initialisation of branch container and histograms
69
70 AliInfo(Form("*** Initialization of %s", GetName())) ;
71
72 // Get input data
73 fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
74 if (!fChain) {
75 AliError(Form("Input 0 for %s not found\n", GetName()));
76 return ;
77 }
78
79 if (!fESD) {
80 // One should first check if the branch address was taken by some other task
81 char ** address = (char **)GetBranchAddress(0, "ESD") ;
82 if (address)
83 fESD = (AliESD *)(*address) ;
84 if (!fESD)
85 fChain->SetBranchAddress("ESD", &fESD) ;
86 }
87 // The output objects will be written to
88 TDirectory * cdir = gDirectory ;
89 // Open a file for output #0
90 char outputName[1024] ;
91 sprintf(outputName, "%s.root", GetName() ) ;
92 OpenFile(0, outputName , "RECREATE") ;
93 if (cdir)
94 cdir->cd() ;
95
96 // create histograms
97 fhFMD1i = new TH1D("FMD1i", "FMD1i", 100, -0.5, 3);
98 fhFMD2i = new TH1D("FMD2i", "FMD2i", 100, -0.5, 3);
99 fhFMD2o = new TH1D("FMD2o", "FMD2o", 100, -0.5, 3);
100 fhFMD3i = new TH1D("FMD3i", "FMD3i", 100, -0.5, 3);
101 fhFMD3o = new TH1D("FMD3o", "FMD3o", 100, -0.5, 3);
102
103 // create output container
104
105 fOutputContainer = new TObjArray(5) ;
106 fOutputContainer->SetName(GetName()) ;
107
108 fOutputContainer->AddAt(fhFMD1i, 0) ;
109 fOutputContainer->AddAt(fhFMD2i, 1) ;
110 fOutputContainer->AddAt(fhFMD2o, 2) ;
111 fOutputContainer->AddAt(fhFMD3i, 3) ;
112 fOutputContainer->AddAt(fhFMD3o, 4) ;
113}
114
115//______________________________________________________________________________
116void AliFMDQATask::Exec(Option_t *)
117{
118 // Processing of one event
119
120 Long64_t entry = fChain->GetReadEntry() ;
121
122 if (!fESD) {
123 AliError("fESD is not connected to the input!") ;
124 return ;
125 }
126
127 if ( !((entry-1)%100) )
128 AliInfo(Form("%s ----> Processing event # %lld", (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ;
129
130 // ************************ FMD *************************************
131
132 AliESDFMD * fmd = fESD->GetFMDData() ;
133
134 Int_t nFMD1 = 0, nFMD2i = 0, nFMD2o = 0, nFMD3i = 0, nFMD3o = 0 ;
135
136 UShort_t detector = 1 ;
137 for(detector = 1 ; detector <= fmd->MaxDetectors() ; detector++) {
138 Char_t ring = 'O' ;
139 UShort_t sector ;
140 for(sector = 0 ;sector < fmd->MaxSectors() ; sector++) {
141 UShort_t strip ;
142 for(strip = 0 ; strip < fmd->MaxStrips(); strip++) {
143 if(fmd->Multiplicity(detector, ring, sector, strip) != AliESDFMD::kInvalidMult)
144 RingSelector(detector, ring, fmd->Multiplicity(detector, ring, sector, strip)) ;
145 if( (fmd->Multiplicity(detector, ring, sector, strip) == AliESDFMD::kInvalidMult) && detector == 2 )
146 nFMD2o++ ;
147 if( (fmd->Multiplicity(detector, ring, sector, strip)==AliESDFMD::kInvalidMult) && detector == 3 )
148 nFMD3o++ ;
149 }
150 }
151 ring='I';
152 for(sector = 0; sector < fmd->MaxSectors(); sector++) {
153 UShort_t strip ;
154 for(strip = 0 ; strip < fmd->MaxStrips() ; strip++) {
155 if(fmd->Multiplicity(detector, ring, sector, strip) != AliESDFMD::kInvalidMult)
156 RingSelector(detector, ring, fmd->Multiplicity(detector, ring, sector, strip));
157 if( (fmd->Multiplicity(detector, ring, sector, strip) == AliESDFMD::kInvalidMult) && detector == 1 )
158 nFMD1++;
159 if( (fmd->Multiplicity(detector, ring, sector, strip) == AliESDFMD::kInvalidMult) && detector == 2 )
160 nFMD2i++;
161 if( (fmd->Multiplicity(detector, ring, sector, strip) == AliESDFMD::kInvalidMult) && detector == 3 )
162 nFMD3i++;
163 }
164 }
165 }
166
167 if(nFMD1>100+10240)
168 AliWarning(Form("number of missing strips in FMD1i too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
169 if(nFMD2i>100+10240)
170 AliWarning(Form("number of missing strips in FMD2i too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
171 if(nFMD2o>100+10240)
172 AliWarning(Form("number of missing strips in FMD2o too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
173 if(nFMD3i>100+10240)
174 AliWarning(Form("number of missing strips in FMD3i too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
175 if(nFMD3o>100+10240)
176 AliWarning(Form("number of missing strips in FMD3o too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
177
178 PostData(0, fOutputContainer);
179}
180
181//______________________________________________________________________________
182void AliFMDQATask::Terminate(Option_t *)
183{
184 // Processing when the event loop is ended
185
186 TCanvas * cFMD1 = new TCanvas("cFMD1", "FMD ESD Test", 400, 10, 600, 700);
187 cFMD1->Divide(3, 2) ;
188
189 cFMD1->cd(1) ;;
190 fhFMD1i->Draw() ;
191
192 cFMD1->cd(2) ;;
193 fhFMD2i->Draw() ;
194
195 cFMD1->cd(3) ;;
196 fhFMD2o->Draw() ;
197
198 cFMD1->cd(4) ;;
199 fhFMD3i->Draw() ;
200
201 cFMD1->cd(5) ;;
202 fhFMD3o->Draw() ;
203
204 cFMD1->Print("FMD.eps") ;
205
206 Bool_t rv1i = TestHisto(fhFMD1i) ;
207 Bool_t rv2i = TestHisto(fhFMD2i) ;
208 Bool_t rv2o = TestHisto(fhFMD2o) ;
209 Bool_t rv3i = TestHisto(fhFMD2i) ;
210 Bool_t rv3o = TestHisto(fhFMD2o) ;
211
212 if ( !(rv1i * rv2i * rv2o * rv3i * rv3o) )
213 AliWarning("Possible problem in file !!! Check output!") ;
214
215
216
217 char line[1024] ;
218 sprintf(line, ".!tar -zcvf %s.tar.gz *.eps", GetName()) ;
219 gROOT->ProcessLine(line);
220 sprintf(line, ".!rm -fR *.eps");
221 gROOT->ProcessLine(line);
222
223 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!! \n", GetName())) ;
224}
225
226//______________________________________________________________________________
227void AliFMDQATask::RingSelector(const UShort_t detector, const Char_t ring, const Float_t mult) const
228{
229 // fill the histograms for each ring in each detector layer
230
231 if(ring == 'I' && detector == 1)
232 fhFMD1i->Fill(mult) ;
233 if(ring == 'I' && detector == 2)
234 fhFMD2i->Fill(mult) ;
235 if(ring == 'O' && detector == 2)
236 fhFMD2o ->Fill(mult) ;
237 if(ring == 'I' && detector == 3)
238 fhFMD3i ->Fill(mult) ;
239 if(ring == 'O' && detector == 3)
240 fhFMD3o ->Fill(mult) ;
241}
242//______________________________________________________________________________
243Bool_t AliFMDQATask::TestHisto(TH1D * hTest) const
244{
245 // analyses the histogram with a Landau function
246
247 Float_t chiMax = 3, chiLow=0.5 ;
248 Float_t chiSq ;
249 Float_t mpv ;
250 Int_t ndf ;
251
252 FitAll(hTest, chiSq, ndf, mpv, chiMax, chiLow);
253
254 if( (chiSq > chiMax || chiSq < chiLow) || mpv < 0.6 || mpv > 1 ) {
255 hTest->Rebin(2) ;
256 FitAll(hTest, chiSq, ndf, mpv, chiMax, chiLow) ;
257 }
258
259 Bool_t ret = kFALSE ;
260 Char_t * test = "not OK";
261 Char_t * test2 = "not OK";
262
263 if(chiSq < chiMax && chiSq > chiLow)
264 test = "OK" ;
265 if(mpv > 0.6 && mpv < 1)
266 test2 = "OK" ;
267
268 if(test == "OK" && test2 == "OK")
269 ret = kTRUE;
270
271 if(test == "not OK" || test2 == "not OK") {
272 AliWarning("Bad fit results") ;
273 printf("Detector : %s\n", hTest->GetName()) ;
274 printf("Landau fit Chi Square / NDF = %f / %d which is %s\n", chiSq*ndf, ndf, test) ;
275 printf("Landau fit MPV is: %f which is %s\n", mpv, test2) ;
276 }
277 return ret;
278}
279//______________________________________________________________________________
280void AliFMDQATask::FitAll(TH1D* hTest, Float_t &chiSq, Int_t &ndf, Float_t &mpv, Float_t chiMax, Float_t chiLow ) const
281{
282 // fit a histogram with a Landau distribution and returns chi square
283
284 Float_t fitMax = hTest->GetXaxis()->GetXmax() ;
285
286 TH1D hTmp = *hTest ;
287 hTmp.SetAxisRange(0.4,fitMax) ;
288 Int_t maxBin = hTmp.GetMaximumBin();
289 Float_t max = hTmp.GetBinCenter(maxBin);
290
291 hTest->Fit("landau", "QOI", "", max-0.3, fitMax) ;
292 TF1* fitfunc = hTest->GetFunction("landau") ;
293 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
294 mpv = fitfunc->GetParameter(1) ;
295 ndf = fitfunc->GetNDF() ;
296
297 if( ( chiSq > chiMax || chiSq < chiLow ) || ( mpv < 0.6 || mpv > 1 ) ) {
298 hTest->Fit("landau", "QOI", "", max-0.2, fitMax) ;
299 fitfunc = hTest->GetFunction("landau") ;
300 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
301 mpv = fitfunc->GetParameter(1) ;
302 ndf = fitfunc->GetNDF() ;
303 }
304 if( ( chiSq >chiMax || chiSq < chiLow ) || ( mpv < 0.6 || mpv > 1 ) ) {
305 hTest->Fit("landau", "QOI", "", max-0.1, fitMax) ;
306 fitfunc = hTest->GetFunction("landau") ;
307 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
308 mpv = fitfunc->GetParameter(1) ;
309 ndf = fitfunc->GetNDF() ;
310 }
311 if( ( chiSq > chiMax || chiSq <chiLow ) || ( mpv < 0.6 || mpv > 1 ) ) {
312 hTest->Fit("landau", "QOI", "", max, fitMax) ;
313 fitfunc = hTest->GetFunction("landau") ;
314 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
315 mpv = fitfunc->GetParameter(1) ;
316 ndf = fitfunc->GetNDF();
317 }
318 if( ( chiSq > chiMax || chiSq < chiLow ) || ( mpv < 0.6 ||mpv > 1 ) ) {
319 hTest->Fit("landau", "QOI", "", max+0.1, fitMax) ;
320 fitfunc = hTest->GetFunction("landau") ;
321 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
322 mpv = fitfunc->GetParameter(1) ;
323 ndf = fitfunc->GetNDF() ;
324 }
325}