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 **************************************************************************/
18 //_________________________________________________________________________
19 // An analysis task to check the FMD data in simulated data
21 //*-- Hans Hjersing Dalsgaard
22 //////////////////////////////////////////////////////////////////////////////
31 #include "AliFMDQATask.h"
35 //______________________________________________________________________________
36 AliFMDQATask::AliFMDQATask(const char *name) :
37 AliAnalysisTask(name,""),
47 // Input slot #0 works with an Ntuple
48 DefineInput(0, TChain::Class());
49 // Output slot #0 writes into a TH1 container
50 DefineOutput(0, TObjArray::Class()) ;
53 //______________________________________________________________________________
54 AliFMDQATask::~AliFMDQATask()
58 fOutputContainer->Clear() ;
59 delete fOutputContainer ;
69 //______________________________________________________________________________
70 void AliFMDQATask::ConnectInputData(const Option_t*)
72 // Initialisation of branch container and histograms
74 AliInfo(Form("*** Initialization of %s", GetName())) ;
77 fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
79 AliError(Form("Input 0 for %s not found\n", GetName()));
83 // One should first check if the branch address was taken by some other task
84 char ** address = (char **)GetBranchAddress(0, "ESD");
86 fESD = (AliESD*)(*address);
89 SetBranchAddress(0, "ESD", &fESD);
93 //________________________________________________________________________
94 void AliFMDQATask::CreateOutputObjects()
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 TFile * currentFile = (dynamic_cast<TChain *>(fChain))->GetFile() ;
128 if ( !((entry-1)%100) )
129 AliInfo(Form("%s ----> Processing event # %lld", currentFile->GetName(), entry)) ;
131 // ************************ FMD *************************************
133 AliESDFMD * fmd = fESD->GetFMDData() ;
135 fmd->CheckNeedUShort(currentFile);
137 Int_t nFMD1 = 0, nFMD2i = 0, nFMD2o = 0, nFMD3i = 0, nFMD3o = 0 ;
139 UShort_t detector = 1 ;
140 for(detector = 1 ; detector <= fmd->MaxDetectors() ; detector++) {
143 for(sector = 0 ;sector < fmd->MaxSectors() ; sector++) {
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 )
150 if( (fmd->Multiplicity(detector, ring, sector, strip)==AliESDFMD::kInvalidMult) && detector == 3 )
155 for(sector = 0; sector < fmd->MaxSectors(); sector++) {
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 )
162 if( (fmd->Multiplicity(detector, ring, sector, strip) == AliESDFMD::kInvalidMult) && detector == 2 )
164 if( (fmd->Multiplicity(detector, ring, sector, strip) == AliESDFMD::kInvalidMult) && detector == 3 )
171 AliWarning(Form("number of missing strips in FMD1i too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
173 AliWarning(Form("number of missing strips in FMD2i too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
175 AliWarning(Form("number of missing strips in FMD2o too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
177 AliWarning(Form("number of missing strips in FMD3i too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
179 AliWarning(Form("number of missing strips in FMD3o too high in event number %lld in file", entry, fChain->GetCurrentFile()->GetName())) ;
181 PostData(0, fOutputContainer);
184 //______________________________________________________________________________
185 void AliFMDQATask::Terminate(Option_t *)
187 // Processing when the event loop is ended
189 fOutputContainer = (TObjArray*)GetOutputData(0);
190 fhFMD1i = (TH1D*)fOutputContainer->At(0);
191 fhFMD2i = (TH1D*)fOutputContainer->At(1);
192 fhFMD2o = (TH1D*)fOutputContainer->At(2);
193 fhFMD3i = (TH1D*)fOutputContainer->At(3);
194 fhFMD3o = (TH1D*)fOutputContainer->At(4);
196 TCanvas * cFMD1 = new TCanvas("cFMD1", "FMD ESD Test", 400, 10, 600, 700);
197 cFMD1->Divide(3, 2) ;
214 cFMD1->Print("FMD.eps") ;
216 Bool_t rv1i = TestHisto(fhFMD1i) ;
217 Bool_t rv2i = TestHisto(fhFMD2i) ;
218 Bool_t rv2o = TestHisto(fhFMD2o) ;
219 Bool_t rv3i = TestHisto(fhFMD2i) ;
220 Bool_t rv3o = TestHisto(fhFMD2o) ;
222 if ( !(rv1i * rv2i * rv2o * rv3i * rv3o) )
223 AliWarning("Possible problem in file !!! Check output!") ;
228 sprintf(line, ".!tar -zcvf %s.tar.gz *.eps", GetName()) ;
229 gROOT->ProcessLine(line);
230 sprintf(line, ".!rm -fR *.eps");
231 gROOT->ProcessLine(line);
233 AliInfo(Form("!!! All the eps files are in %s.tar.gz !!! \n", GetName())) ;
236 //______________________________________________________________________________
237 void AliFMDQATask::RingSelector(const UShort_t detector, const Char_t ring, const Float_t mult) const
239 // fill the histograms for each ring in each detector layer
241 if(ring == 'I' && detector == 1)
242 fhFMD1i->Fill(mult) ;
243 if(ring == 'I' && detector == 2)
244 fhFMD2i->Fill(mult) ;
245 if(ring == 'O' && detector == 2)
246 fhFMD2o ->Fill(mult) ;
247 if(ring == 'I' && detector == 3)
248 fhFMD3i ->Fill(mult) ;
249 if(ring == 'O' && detector == 3)
250 fhFMD3o ->Fill(mult) ;
252 //______________________________________________________________________________
253 Bool_t AliFMDQATask::TestHisto(TH1D * hTest) const
255 // analyses the histogram with a Landau function
257 Float_t chiMax = 3, chiLow=0.5 ;
262 FitAll(hTest, chiSq, ndf, mpv, chiMax, chiLow);
264 if( (chiSq > chiMax || chiSq < chiLow) || mpv < 0.6 || mpv > 1 ) {
266 FitAll(hTest, chiSq, ndf, mpv, chiMax, chiLow) ;
269 Bool_t ret = kFALSE ;
270 Char_t * test = "not OK";
271 Char_t * test2 = "not OK";
273 if(chiSq < chiMax && chiSq > chiLow)
275 if(mpv > 0.6 && mpv < 1)
278 if(test == "OK" && test2 == "OK")
281 if(test == "not OK" || test2 == "not OK") {
282 AliWarning("Bad fit results") ;
283 printf("Detector : %s\n", hTest->GetName()) ;
284 printf("Landau fit Chi Square / NDF = %f / %d which is %s\n", chiSq*ndf, ndf, test) ;
285 printf("Landau fit MPV is: %f which is %s\n", mpv, test2) ;
289 //______________________________________________________________________________
290 void AliFMDQATask::FitAll(TH1D* hTest, Float_t &chiSq, Int_t &ndf, Float_t &mpv, Float_t chiMax, Float_t chiLow ) const
292 // fit a histogram with a Landau distribution and returns chi square
294 Float_t fitMax = hTest->GetXaxis()->GetXmax() ;
297 hTmp.SetAxisRange(0.4,fitMax) ;
298 Int_t maxBin = hTmp.GetMaximumBin();
299 Float_t max = hTmp.GetBinCenter(maxBin);
301 hTest->Fit("landau", "QOI", "", max-0.3, fitMax) ;
302 TF1* fitfunc = hTest->GetFunction("landau") ;
303 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
304 mpv = fitfunc->GetParameter(1) ;
305 ndf = fitfunc->GetNDF() ;
307 if( ( chiSq > chiMax || chiSq < chiLow ) || ( mpv < 0.6 || mpv > 1 ) ) {
308 hTest->Fit("landau", "QOI", "", max-0.2, fitMax) ;
309 fitfunc = hTest->GetFunction("landau") ;
310 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
311 mpv = fitfunc->GetParameter(1) ;
312 ndf = fitfunc->GetNDF() ;
314 if( ( chiSq >chiMax || chiSq < chiLow ) || ( mpv < 0.6 || mpv > 1 ) ) {
315 hTest->Fit("landau", "QOI", "", max-0.1, fitMax) ;
316 fitfunc = hTest->GetFunction("landau") ;
317 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
318 mpv = fitfunc->GetParameter(1) ;
319 ndf = fitfunc->GetNDF() ;
321 if( ( chiSq > chiMax || chiSq <chiLow ) || ( mpv < 0.6 || mpv > 1 ) ) {
322 hTest->Fit("landau", "QOI", "", max, fitMax) ;
323 fitfunc = hTest->GetFunction("landau") ;
324 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
325 mpv = fitfunc->GetParameter(1) ;
326 ndf = fitfunc->GetNDF();
328 if( ( chiSq > chiMax || chiSq < chiLow ) || ( mpv < 0.6 ||mpv > 1 ) ) {
329 hTest->Fit("landau", "QOI", "", max+0.1, fitMax) ;
330 fitfunc = hTest->GetFunction("landau") ;
331 chiSq = fitfunc->GetChisquare() / fitfunc->GetNDF() ;
332 mpv = fitfunc->GetParameter(1) ;
333 ndf = fitfunc->GetNDF() ;