]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - ESDCheck/AliFMDQATask.cxx
No effective C++ option for compilation of C files
[u/mrichter/AliRoot.git] / ESDCheck / AliFMDQATask.cxx
... / ...
CommitLineData
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 TFile * currentFile = (dynamic_cast<TChain *>(fChain))->GetFile() ;
128 if ( !((entry-1)%100) )
129 AliInfo(Form("%s ----> Processing event # %lld", currentFile->GetName(), entry)) ;
130
131 // ************************ FMD *************************************
132
133 AliESDFMD * fmd = fESD->GetFMDData() ;
134
135 fmd->CheckNeedUShort(currentFile);
136
137 Int_t nFMD1 = 0, nFMD2i = 0, nFMD2o = 0, nFMD3i = 0, nFMD3o = 0 ;
138
139 UShort_t detector = 1 ;
140 for(detector = 1 ; detector <= fmd->MaxDetectors() ; detector++) {
141 Char_t ring = 'O' ;
142 UShort_t sector ;
143 for(sector = 0 ;sector < fmd->MaxSectors() ; sector++) {
144 UShort_t strip ;
145 for(strip = 0 ; strip < fmd->MaxStrips(); strip++) {
146 if(fmd->Multiplicity(detector, ring, sector, strip) != AliESDFMD::kInvalidMult)
147 RingSelector(detector, ring, fmd->Multiplicity(detector, ring, sector, strip)) ;
148 if( (fmd->Multiplicity(detector, ring, sector, strip) == AliESDFMD::kInvalidMult) && detector == 2 )
149 nFMD2o++ ;
150 if( (fmd->Multiplicity(detector, ring, sector, strip)==AliESDFMD::kInvalidMult) && detector == 3 )
151 nFMD3o++ ;
152 }
153 }
154 ring='I';
155 for(sector = 0; sector < fmd->MaxSectors(); sector++) {
156 UShort_t strip ;
157 for(strip = 0 ; strip < fmd->MaxStrips() ; strip++) {
158 if(fmd->Multiplicity(detector, ring, sector, strip) != AliESDFMD::kInvalidMult)
159 RingSelector(detector, ring, fmd->Multiplicity(detector, ring, sector, strip));
160 if( (fmd->Multiplicity(detector, ring, sector, strip) == AliESDFMD::kInvalidMult) && detector == 1 )
161 nFMD1++;
162 if( (fmd->Multiplicity(detector, ring, sector, strip) == AliESDFMD::kInvalidMult) && detector == 2 )
163 nFMD2i++;
164 if( (fmd->Multiplicity(detector, ring, sector, strip) == AliESDFMD::kInvalidMult) && detector == 3 )
165 nFMD3i++;
166 }
167 }
168 }
169
170 if(nFMD1>100+10240)
171 AliWarning(Form("number of missing strips in FMD1i too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
172 if(nFMD2i>100+10240)
173 AliWarning(Form("number of missing strips in FMD2i too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
174 if(nFMD2o>100+10240)
175 AliWarning(Form("number of missing strips in FMD2o too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
176 if(nFMD3i>100+10240)
177 AliWarning(Form("number of missing strips in FMD3i too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
178 if(nFMD3o>100+10240)
179 AliWarning(Form("number of missing strips in FMD3o too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
180
181 PostData(0, fOutputContainer);
182}
183
184//______________________________________________________________________________
185void AliFMDQATask::Terminate(Option_t *)
186{
187 // Processing when the event loop is ended
188
189 TCanvas * cFMD1 = new TCanvas("cFMD1", "FMD ESD Test", 400, 10, 600, 700);
190 cFMD1->Divide(3, 2) ;
191
192 cFMD1->cd(1) ;;
193 fhFMD1i->Draw() ;
194
195 cFMD1->cd(2) ;;
196 fhFMD2i->Draw() ;
197
198 cFMD1->cd(3) ;;
199 fhFMD2o->Draw() ;
200
201 cFMD1->cd(4) ;;
202 fhFMD3i->Draw() ;
203
204 cFMD1->cd(5) ;;
205 fhFMD3o->Draw() ;
206
207 cFMD1->Print("FMD.eps") ;
208
209 Bool_t rv1i = TestHisto(fhFMD1i) ;
210 Bool_t rv2i = TestHisto(fhFMD2i) ;
211 Bool_t rv2o = TestHisto(fhFMD2o) ;
212 Bool_t rv3i = TestHisto(fhFMD2i) ;
213 Bool_t rv3o = TestHisto(fhFMD2o) ;
214
215 if ( !(rv1i * rv2i * rv2o * rv3i * rv3o) )
216 AliWarning("Possible problem in file !!! Check output!") ;
217
218
219
220 char line[1024] ;
221 sprintf(line, ".!tar -zcvf %s.tar.gz *.eps", GetName()) ;
222 gROOT->ProcessLine(line);
223 sprintf(line, ".!rm -fR *.eps");
224 gROOT->ProcessLine(line);
225
226 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!! \n", GetName())) ;
227}
228
229//______________________________________________________________________________
230void AliFMDQATask::RingSelector(const UShort_t detector, const Char_t ring, const Float_t mult) const
231{
232 // fill the histograms for each ring in each detector layer
233
234 if(ring == 'I' && detector == 1)
235 fhFMD1i->Fill(mult) ;
236 if(ring == 'I' && detector == 2)
237 fhFMD2i->Fill(mult) ;
238 if(ring == 'O' && detector == 2)
239 fhFMD2o ->Fill(mult) ;
240 if(ring == 'I' && detector == 3)
241 fhFMD3i ->Fill(mult) ;
242 if(ring == 'O' && detector == 3)
243 fhFMD3o ->Fill(mult) ;
244}
245//______________________________________________________________________________
246Bool_t AliFMDQATask::TestHisto(TH1D * hTest) const
247{
248 // analyses the histogram with a Landau function
249
250 Float_t chiMax = 3, chiLow=0.5 ;
251 Float_t chiSq ;
252 Float_t mpv ;
253 Int_t ndf ;
254
255 FitAll(hTest, chiSq, ndf, mpv, chiMax, chiLow);
256
257 if( (chiSq > chiMax || chiSq < chiLow) || mpv < 0.6 || mpv > 1 ) {
258 hTest->Rebin(2) ;
259 FitAll(hTest, chiSq, ndf, mpv, chiMax, chiLow) ;
260 }
261
262 Bool_t ret = kFALSE ;
263 Char_t * test = "not OK";
264 Char_t * test2 = "not OK";
265
266 if(chiSq < chiMax && chiSq > chiLow)
267 test = "OK" ;
268 if(mpv > 0.6 && mpv < 1)
269 test2 = "OK" ;
270
271 if(test == "OK" && test2 == "OK")
272 ret = kTRUE;
273
274 if(test == "not OK" || test2 == "not OK") {
275 AliWarning("Bad fit results") ;
276 printf("Detector : %s\n", hTest->GetName()) ;
277 printf("Landau fit Chi Square / NDF = %f / %d which is %s\n", chiSq*ndf, ndf, test) ;
278 printf("Landau fit MPV is: %f which is %s\n", mpv, test2) ;
279 }
280 return ret;
281}
282//______________________________________________________________________________
283void AliFMDQATask::FitAll(TH1D* hTest, Float_t &chiSq, Int_t &ndf, Float_t &mpv, Float_t chiMax, Float_t chiLow ) const
284{
285 // fit a histogram with a Landau distribution and returns chi square
286
287 Float_t fitMax = hTest->GetXaxis()->GetXmax() ;
288
289 TH1D hTmp = *hTest ;
290 hTmp.SetAxisRange(0.4,fitMax) ;
291 Int_t maxBin = hTmp.GetMaximumBin();
292 Float_t max = hTmp.GetBinCenter(maxBin);
293
294 hTest->Fit("landau", "QOI", "", max-0.3, fitMax) ;
295 TF1* fitfunc = hTest->GetFunction("landau") ;
296 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
297 mpv = fitfunc->GetParameter(1) ;
298 ndf = fitfunc->GetNDF() ;
299
300 if( ( chiSq > chiMax || chiSq < chiLow ) || ( mpv < 0.6 || mpv > 1 ) ) {
301 hTest->Fit("landau", "QOI", "", max-0.2, fitMax) ;
302 fitfunc = hTest->GetFunction("landau") ;
303 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
304 mpv = fitfunc->GetParameter(1) ;
305 ndf = fitfunc->GetNDF() ;
306 }
307 if( ( chiSq >chiMax || chiSq < chiLow ) || ( mpv < 0.6 || mpv > 1 ) ) {
308 hTest->Fit("landau", "QOI", "", max-0.1, fitMax) ;
309 fitfunc = hTest->GetFunction("landau") ;
310 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
311 mpv = fitfunc->GetParameter(1) ;
312 ndf = fitfunc->GetNDF() ;
313 }
314 if( ( chiSq > chiMax || chiSq <chiLow ) || ( mpv < 0.6 || mpv > 1 ) ) {
315 hTest->Fit("landau", "QOI", "", max, fitMax) ;
316 fitfunc = hTest->GetFunction("landau") ;
317 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
318 mpv = fitfunc->GetParameter(1) ;
319 ndf = fitfunc->GetNDF();
320 }
321 if( ( chiSq > chiMax || chiSq < chiLow ) || ( mpv < 0.6 ||mpv > 1 ) ) {
322 hTest->Fit("landau", "QOI", "", max+0.1, fitMax) ;
323 fitfunc = hTest->GetFunction("landau") ;
324 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
325 mpv = fitfunc->GetParameter(1) ;
326 ndf = fitfunc->GetNDF() ;
327 }
328}