1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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
18 //*-- Hans Hjersing Dalsgaard
19 //////////////////////////////////////////////////////////////////////////////
27 #include "AliFMDQATask.h"
31 //______________________________________________________________________________
32 AliFMDQATask::AliFMDQATask(const char *name) :
33 AliAnalysisTask(name,""),
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()) ;
49 //______________________________________________________________________________
50 AliFMDQATask::~AliFMDQATask()
54 fOutputContainer->Clear() ;
55 delete fOutputContainer ;
65 //______________________________________________________________________________
66 void AliFMDQATask::Init(const Option_t*)
68 // Initialisation of branch container and histograms
70 AliInfo(Form("*** Initialization of %s", GetName())) ;
73 fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
75 AliError(Form("Input 0 for %s not found\n", GetName()));
80 // One should first check if the branch address was taken by some other task
81 char ** address = (char **)GetBranchAddress(0, "ESD") ;
83 fESD = (AliESD *)(*address) ;
85 fChain->SetBranchAddress("ESD", &fESD) ;
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") ;
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);
103 // create output container
105 fOutputContainer = new TObjArray(5) ;
106 fOutputContainer->SetName(GetName()) ;
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) ;
115 //______________________________________________________________________________
116 void AliFMDQATask::Exec(Option_t *)
118 // Processing of one event
120 Long64_t entry = fChain->GetReadEntry() ;
123 AliError("fESD is not connected to the input!") ;
127 if ( !((entry-1)%100) )
128 AliInfo(Form("%s ----> Processing event # %lld", (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ;
130 // ************************ FMD *************************************
132 AliESDFMD * fmd = fESD->GetFMDData() ;
133 fmd->CheckNeedUShort(fChain->GetFile());
135 Int_t nFMD1 = 0, nFMD2i = 0, nFMD2o = 0, nFMD3i = 0, nFMD3o = 0 ;
137 UShort_t detector = 1 ;
138 for(detector = 1 ; detector <= fmd->MaxDetectors() ; detector++) {
141 for(sector = 0 ;sector < fmd->MaxSectors() ; sector++) {
143 for(strip = 0 ; strip < fmd->MaxStrips(); strip++) {
144 if(fmd->Multiplicity(detector, ring, sector, strip) != AliESDFMD::kInvalidMult)
145 RingSelector(detector, ring, fmd->Multiplicity(detector, ring, sector, strip)) ;
146 if( (fmd->Multiplicity(detector, ring, sector, strip) == AliESDFMD::kInvalidMult) && detector == 2 )
148 if( (fmd->Multiplicity(detector, ring, sector, strip)==AliESDFMD::kInvalidMult) && detector == 3 )
153 for(sector = 0; sector < fmd->MaxSectors(); sector++) {
155 for(strip = 0 ; strip < fmd->MaxStrips() ; strip++) {
156 if(fmd->Multiplicity(detector, ring, sector, strip) != AliESDFMD::kInvalidMult)
157 RingSelector(detector, ring, fmd->Multiplicity(detector, ring, sector, strip));
158 if( (fmd->Multiplicity(detector, ring, sector, strip) == AliESDFMD::kInvalidMult) && detector == 1 )
160 if( (fmd->Multiplicity(detector, ring, sector, strip) == AliESDFMD::kInvalidMult) && detector == 2 )
162 if( (fmd->Multiplicity(detector, ring, sector, strip) == AliESDFMD::kInvalidMult) && detector == 3 )
169 AliWarning(Form("number of missing strips in FMD1i too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
171 AliWarning(Form("number of missing strips in FMD2i too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
173 AliWarning(Form("number of missing strips in FMD2o too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
175 AliWarning(Form("number of missing strips in FMD3i too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
177 AliWarning(Form("number of missing strips in FMD3o too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
179 PostData(0, fOutputContainer);
182 //______________________________________________________________________________
183 void AliFMDQATask::Terminate(Option_t *)
185 // Processing when the event loop is ended
187 TCanvas * cFMD1 = new TCanvas("cFMD1", "FMD ESD Test", 400, 10, 600, 700);
188 cFMD1->Divide(3, 2) ;
205 cFMD1->Print("FMD.eps") ;
207 Bool_t rv1i = TestHisto(fhFMD1i) ;
208 Bool_t rv2i = TestHisto(fhFMD2i) ;
209 Bool_t rv2o = TestHisto(fhFMD2o) ;
210 Bool_t rv3i = TestHisto(fhFMD2i) ;
211 Bool_t rv3o = TestHisto(fhFMD2o) ;
213 if ( !(rv1i * rv2i * rv2o * rv3i * rv3o) )
214 AliWarning("Possible problem in file !!! Check output!") ;
219 sprintf(line, ".!tar -zcvf %s.tar.gz *.eps", GetName()) ;
220 gROOT->ProcessLine(line);
221 sprintf(line, ".!rm -fR *.eps");
222 gROOT->ProcessLine(line);
224 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!! \n", GetName())) ;
227 //______________________________________________________________________________
228 void AliFMDQATask::RingSelector(const UShort_t detector, const Char_t ring, const Float_t mult) const
230 // fill the histograms for each ring in each detector layer
232 if(ring == 'I' && detector == 1)
233 fhFMD1i->Fill(mult) ;
234 if(ring == 'I' && detector == 2)
235 fhFMD2i->Fill(mult) ;
236 if(ring == 'O' && detector == 2)
237 fhFMD2o ->Fill(mult) ;
238 if(ring == 'I' && detector == 3)
239 fhFMD3i ->Fill(mult) ;
240 if(ring == 'O' && detector == 3)
241 fhFMD3o ->Fill(mult) ;
243 //______________________________________________________________________________
244 Bool_t AliFMDQATask::TestHisto(TH1D * hTest) const
246 // analyses the histogram with a Landau function
248 Float_t chiMax = 3, chiLow=0.5 ;
253 FitAll(hTest, chiSq, ndf, mpv, chiMax, chiLow);
255 if( (chiSq > chiMax || chiSq < chiLow) || mpv < 0.6 || mpv > 1 ) {
257 FitAll(hTest, chiSq, ndf, mpv, chiMax, chiLow) ;
260 Bool_t ret = kFALSE ;
261 Char_t * test = "not OK";
262 Char_t * test2 = "not OK";
264 if(chiSq < chiMax && chiSq > chiLow)
266 if(mpv > 0.6 && mpv < 1)
269 if(test == "OK" && test2 == "OK")
272 if(test == "not OK" || test2 == "not OK") {
273 AliWarning("Bad fit results") ;
274 printf("Detector : %s\n", hTest->GetName()) ;
275 printf("Landau fit Chi Square / NDF = %f / %d which is %s\n", chiSq*ndf, ndf, test) ;
276 printf("Landau fit MPV is: %f which is %s\n", mpv, test2) ;
280 //______________________________________________________________________________
281 void AliFMDQATask::FitAll(TH1D* hTest, Float_t &chiSq, Int_t &ndf, Float_t &mpv, Float_t chiMax, Float_t chiLow ) const
283 // fit a histogram with a Landau distribution and returns chi square
285 Float_t fitMax = hTest->GetXaxis()->GetXmax() ;
288 hTmp.SetAxisRange(0.4,fitMax) ;
289 Int_t maxBin = hTmp.GetMaximumBin();
290 Float_t max = hTmp.GetBinCenter(maxBin);
292 hTest->Fit("landau", "QOI", "", max-0.3, fitMax) ;
293 TF1* fitfunc = hTest->GetFunction("landau") ;
294 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
295 mpv = fitfunc->GetParameter(1) ;
296 ndf = fitfunc->GetNDF() ;
298 if( ( chiSq > chiMax || chiSq < chiLow ) || ( mpv < 0.6 || mpv > 1 ) ) {
299 hTest->Fit("landau", "QOI", "", max-0.2, fitMax) ;
300 fitfunc = hTest->GetFunction("landau") ;
301 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
302 mpv = fitfunc->GetParameter(1) ;
303 ndf = fitfunc->GetNDF() ;
305 if( ( chiSq >chiMax || chiSq < chiLow ) || ( mpv < 0.6 || mpv > 1 ) ) {
306 hTest->Fit("landau", "QOI", "", max-0.1, fitMax) ;
307 fitfunc = hTest->GetFunction("landau") ;
308 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
309 mpv = fitfunc->GetParameter(1) ;
310 ndf = fitfunc->GetNDF() ;
312 if( ( chiSq > chiMax || chiSq <chiLow ) || ( mpv < 0.6 || mpv > 1 ) ) {
313 hTest->Fit("landau", "QOI", "", max, fitMax) ;
314 fitfunc = hTest->GetFunction("landau") ;
315 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
316 mpv = fitfunc->GetParameter(1) ;
317 ndf = fitfunc->GetNDF();
319 if( ( chiSq > chiMax || chiSq < chiLow ) || ( mpv < 0.6 ||mpv > 1 ) ) {
320 hTest->Fit("landau", "QOI", "", max+0.1, fitMax) ;
321 fitfunc = hTest->GetFunction("landau") ;
322 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
323 mpv = fitfunc->GetParameter(1) ;
324 ndf = fitfunc->GetNDF() ;