]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ESDCheck/AliT0QATask.cxx
Updated version.
[u/mrichter/AliRoot.git] / ESDCheck / AliT0QATask.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 /* $Id$ */
17
18 //______________________________________________________________________________
19 // An analysis task to check the T0 data in simulated data
20 //
21 //*-- Alla Maevskaya
22 //////////////////////////////////////////////////////////////////////////////
23
24 #include <TCanvas.h>
25 #include <TChain.h>
26 #include <TFile.h> 
27 #include <TH1F.h>
28 #include <TLegend.h> 
29 #include <TROOT.h>
30 #include <TString.h> 
31
32 #include "AliT0QATask.h" 
33 #include "AliESD.h" 
34 #include "AliLog.h"
35 #include "AliESDVertex.h" 
36
37 //______________________________________________________________________________
38 AliT0QATask::AliT0QATask(const char *name) : 
39   AliAnalysisTask(name,""),  
40   fChain(0),
41   fESD(0), 
42   fOutputContainer(0),
43   fhT01(0),
44   fhT02(0),
45   fhT03(0)
46 {
47   // Constructor.
48   // Input slot #0 works with an Ntuple
49   DefineInput(0, TChain::Class());
50   // Output slot #0 writes into a TH1 container
51   DefineOutput(0,  TObjArray::Class()) ; 
52 }
53
54 //______________________________________________________________________________
55 AliT0QATask::~AliT0QATask()
56 {
57   // dtor
58   fOutputContainer->Clear() ; 
59   delete fOutputContainer ;
60  
61   delete fhT01 ; 
62   delete fhT02 ;
63   delete fhT03 ; 
64 }
65
66 //______________________________________________________________________________
67 void AliT0QATask::ConnectInputData(const Option_t*)
68 {
69   // Initialisation of branch container and histograms 
70     
71   AliInfo(Form("*** Initialization of %s", GetName())) ; 
72   
73   // Get input data
74   fChain = dynamic_cast<TChain *>(GetInputData(0)) ;
75   if (!fChain) {
76     AliError(Form("Input 0 for %s not found\n", GetName()));
77     return ;
78   }
79   
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   } else {
85     fESD = new AliESD();
86     SetBranchAddress(0, "ESD", &fESD);
87   }
88 }
89   
90 //________________________________________________________________________
91 void AliT0QATask::CreateOutputObjects()
92 {
93   // create histograms   
94
95   OpenFile(0) ; 
96
97   fhT01 = new TH1F("hRealVertex", "Primary vertex", 100,   -20,    20);
98   fhT02 = new TH1F("hT0start",    "T0 start time",  100, 12400, 12600);
99   fhT03 = new TH1F("hT0vertex",   "T0vertex",       100,   -20,    20);
100
101
102   // create output container
103   
104   fOutputContainer = new TObjArray(3) ; 
105   fOutputContainer->SetName(GetName()) ; 
106
107   fOutputContainer->AddAt(fhT01,             0) ; 
108   fOutputContainer->AddAt(fhT02,             1) ; 
109   fOutputContainer->AddAt(fhT03,             2) ; 
110 }
111
112 //______________________________________________________________________________
113 void AliT0QATask::Exec(Option_t *) 
114 {
115   // Processing of one event
116     
117   Long64_t entry = fChain->GetReadEntry() ;
118   
119   if (!fESD) {
120     AliError("fESD is not connected to the input!") ; 
121     return ; 
122   }
123   
124   if ( !((entry-1)%100) ) 
125     AliInfo(Form("%s ----> Processing event # %lld",  (dynamic_cast<TChain *>(fChain))->GetFile()->GetName(), entry)) ; 
126   
127   // ************************  T0 *************************************
128   
129   const AliESDVertex * vertex = fESD->GetPrimaryVertex() ; 
130   Double_t posZ = vertex->GetZv() ; 
131   fhT01->Fill( posZ ) ;
132
133   fhT02->Fill( fESD->GetT0() ) ;
134   
135   fhT03->Fill( fESD->GetT0zVertex() / 2. ) ;
136   
137   PostData(0, fOutputContainer);
138   
139 }
140
141 //______________________________________________________________________________
142 void AliT0QATask::Terminate(Option_t *)
143 {
144   // Processing when the event loop is ended
145   fOutputContainer = (TObjArray*)GetOutputData(0);
146   fhT01 = (TH1F*)fOutputContainer->At(0);
147   fhT02 = (TH1F*)fOutputContainer->At(1);
148   fhT03 = (TH1F*)fOutputContainer->At(2);
149
150   Bool_t problem = kFALSE ; 
151   AliInfo(Form(" *** %s Report:", GetName())) ; 
152
153   Float_t mean = fhT02->GetMean();
154
155   printf ("mean time T0 ps %f\n", mean) ;
156
157   if ( mean > 12600 || mean < 12400 ) {
158    AliWarning (" !!!!!!!!!!-----events sample is WRONG - T0 unreal -------");  
159    problem = kTRUE ;
160   }
161   TCanvas  * cTO1 = new TCanvas("cT01", "T0 ESD Test", 400, 10, 600, 700) ;
162   cTO1->Divide(2, 2) ;
163
164   cTO1->cd(1) ; 
165   fhT01->Draw() ; 
166     
167   cTO1->cd(2) ; 
168   fhT02->Draw() ;
169
170   cTO1->cd(3) ; 
171   fhT03->Draw() ; 
172
173
174   cTO1->Print("T0.eps");
175
176   char line[1024] ; 
177   sprintf(line, ".!tar -zcf %s.tar.gz *.eps", GetName()) ; 
178   gROOT->ProcessLine(line);
179   sprintf(line, ".!rm -fR *.eps"); 
180   gROOT->ProcessLine(line);
181  
182   AliInfo(Form("!!! All the eps files are in %s.tar.gz !!!", GetName())) ;
183
184   TString report ; 
185   if(problem)
186     report="Problems found, please check!!!";  
187   else 
188     report="OK";
189
190   AliInfo(Form("*** %s Summary Report: %s\n",GetName(), report.Data())) ; 
191 }