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