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