]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCAnalysisTaskcalib.cxx
Adding the new calibration class for the multiplicity correction
[u/mrichter/AliRoot.git] / TPC / AliTPCAnalysisTaskcalib.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
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 // ANALYSIS task to perrorm TPC calibration                                  //
20
21 //                                                                           //
22 ///////////////////////////////////////////////////////////////////////////////
23 #include "AliTPCAnalysisTaskcalib.h"
24 #include "TChain.h"
25 #include "AliTPCcalibBase.h"
26 #include "AliESDEvent.h"
27 #include "AliESDfriend.h"
28 #include "AliESDtrack.h"
29 #include "AliESDfriendTrack.h"
30 #include "AliTPCseed.h"
31 #include "AliESDInputHandler.h"
32 #include "AliAnalysisManager.h"
33 #include "TFile.h"
34 #include "TSystem.h"
35 #include "TTimeStamp.h"
36
37 ClassImp(AliTPCAnalysisTaskcalib)
38
39
40 AliTPCAnalysisTaskcalib::AliTPCAnalysisTaskcalib()
41   :AliAnalysisTask(),
42    fCalibJobs(0),
43    fESD(0),
44    fESDfriend(0),
45    fDebugOutputPath("")
46 {
47   //
48   // default constructor
49   // 
50   
51 }
52
53
54 AliTPCAnalysisTaskcalib::AliTPCAnalysisTaskcalib(const char *name) 
55   :AliAnalysisTask(name,""),
56    fCalibJobs(0),
57    fESD(0),
58    fESDfriend(0),
59    fDebugOutputPath("")
60 {
61   //
62   // Constructor
63   //
64   DefineInput(0, TChain::Class());
65   DefineOutput(0, TObjArray::Class());
66   fCalibJobs = new TObjArray(0);
67   fCalibJobs->SetOwner(kTRUE);
68 }
69
70 AliTPCAnalysisTaskcalib::~AliTPCAnalysisTaskcalib() {
71   //
72   // destructor
73   //
74   printf("AliTPCAnalysisTaskcalib::~AliTPCAnalysisTaskcalib");
75   fCalibJobs->Delete();
76 }
77
78 void AliTPCAnalysisTaskcalib::Exec(Option_t *) {
79   //
80   // Exec function
81   // Loop over tracks and call  Process function
82   if (!fESD) {
83     //Printf("ERROR: fESD not available");
84     return;
85   }
86   fESDfriend=static_cast<AliESDfriend*>(fESD->FindListObject("AliESDfriend"));
87   if (!fESDfriend) {
88     //Printf("ERROR: fESDfriend not available");
89     return;
90   }
91   Int_t n=fESD->GetNumberOfTracks();
92   Process(fESD);
93   Int_t run = fESD->GetRunNumber();
94   for (Int_t i=0;i<n;++i) {
95     AliESDfriendTrack *friendTrack=fESDfriend->GetTrack(i);
96     AliESDtrack *track=fESD->GetTrack(i);
97     TObject *calibObject=0;
98     AliTPCseed *seed=0;
99     if (!friendTrack) continue;
100     for (Int_t j=0;(calibObject=friendTrack->GetCalibObject(j));++j)
101       if ((seed=dynamic_cast<AliTPCseed*>(calibObject)))
102         break;
103     if (track) Process(track, run);
104     if (seed)
105       Process(seed);
106   }
107   PostData(0,fCalibJobs);
108 }
109
110 void AliTPCAnalysisTaskcalib::ConnectInputData(Option_t *) {
111   //
112   //
113   //
114   TTree* tree=dynamic_cast<TTree*>(GetInputData(0));
115   if (!tree) {
116     //Printf("ERROR: Could not read chain from input slot 0");
117   } 
118   else {
119     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
120     if (!esdH) {
121       //Printf("ERROR: Could not get ESDInputHandler");
122     } 
123     else {
124       fESD = esdH->GetEvent();
125       //Printf("*** CONNECTED NEW EVENT ****");
126     }
127   }
128 }
129
130 void AliTPCAnalysisTaskcalib::CreateOutputObjects() {
131   //
132   //
133   //
134   //OpenFile(0, "RECREATE");
135 }
136 void AliTPCAnalysisTaskcalib::Terminate(Option_t */*option*/) {
137   //
138   // Terminate
139   //
140   AliTPCcalibBase *job=0;
141   Int_t njobs = fCalibJobs->GetEntriesFast();
142   for (Int_t i=0;i<njobs;i++){
143     job = (AliTPCcalibBase*)fCalibJobs->UncheckedAt(i);
144     if (job) job->Terminate();
145   }
146   
147 }
148
149 void AliTPCAnalysisTaskcalib::FinishTaskOutput()
150 {
151   //
152   // According description in AliAnalisysTask this method is call 
153   // on the slaves before sending data
154   //
155   Terminate("slave");
156   if(!fDebugOutputPath.IsNull()) { 
157     RegisterDebugOutput();
158   }
159   
160 }
161
162
163 void AliTPCAnalysisTaskcalib::Process(AliESDEvent *event) {
164   //
165   // Process ESD event
166   //
167   AliTPCcalibBase *job=0;
168   Int_t njobs = fCalibJobs->GetEntriesFast();
169   for (Int_t i=0;i<njobs;i++){
170     job = (AliTPCcalibBase*)fCalibJobs->UncheckedAt(i);
171     if (job) {
172       job->UpdateEventInfo(event);
173       if (job->AcceptTrigger())
174         job->Process(event);
175     }
176   }
177 }
178
179 void AliTPCAnalysisTaskcalib::Process(AliTPCseed *track) {
180   //
181   // Process TPC track
182   //
183   AliTPCcalibBase *job=0;
184   Int_t njobs = fCalibJobs->GetEntriesFast();
185   for (Int_t i=0;i<njobs;i++){
186     job = (AliTPCcalibBase*)fCalibJobs->UncheckedAt(i);
187     if (job)  
188       if (job->AcceptTrigger())
189         job->Process(track);
190   }
191 }
192
193 void AliTPCAnalysisTaskcalib::Process(AliESDtrack *track, Int_t run) {
194   //
195   // Process ESD track
196   //
197   AliTPCcalibBase *job=0;
198   Int_t njobs = fCalibJobs->GetEntriesFast();
199   for (Int_t i=0;i<njobs;i++){
200     job = (AliTPCcalibBase*)fCalibJobs->UncheckedAt(i);
201     if (job) 
202       if (job->AcceptTrigger())
203         job->Process(track,run);
204   }
205 }
206
207 Long64_t AliTPCAnalysisTaskcalib::Merge(TCollection *li) {
208   TIterator *i=fCalibJobs->MakeIterator();
209   AliTPCcalibBase *job;
210   Long64_t n=0;
211   while ((job=dynamic_cast<AliTPCcalibBase*>(i->Next())))
212     n+=job->Merge(li);
213   return n;
214 }
215
216 void AliTPCAnalysisTaskcalib::Analyze() {
217   //
218   // Analyze the content of the task
219   //
220   AliTPCcalibBase *job=0;
221   Int_t njobs = fCalibJobs->GetEntriesFast();
222   for (Int_t i=0;i<njobs;i++){
223     job = (AliTPCcalibBase*)fCalibJobs->UncheckedAt(i);
224     if (job) job->Analyze();
225   }
226 }
227
228
229 void AliTPCAnalysisTaskcalib::RegisterDebugOutput(){
230   //
231   //
232   //
233   AliTPCcalibBase *job=0;
234   Int_t njobs = fCalibJobs->GetEntriesFast();
235   for (Int_t i=0;i<njobs;i++){
236     job = (AliTPCcalibBase*)fCalibJobs->UncheckedAt(i);
237     if (job) job->RegisterDebugOutput(fDebugOutputPath.Data());
238   }
239   TString dsName=GetName();
240   dsName+=".root";
241   TFile fff(dsName.Data(),"recreate");
242   fCalibJobs->Write("TPCCalib",TObject::kSingleKey);
243   fff.Close();
244   //
245   // store  - copy debug output to the destination position
246   // currently ONLY for local copy
247   TString dsName2=fDebugOutputPath.Data();
248   gSystem->MakeDirectory(dsName2.Data());
249   dsName2+=gSystem->HostName();
250   gSystem->MakeDirectory(dsName2.Data());
251   dsName2+="/";
252   TTimeStamp s;
253   dsName2+=Int_t(s.GetNanoSec());
254   dsName2+="/";
255   gSystem->MakeDirectory(dsName2.Data());
256   dsName2+=dsName;
257   AliInfo(Form("copy %s\t%s\n",dsName.Data(),dsName2.Data()));
258   printf("copy %s\t%s\n",dsName.Data(),dsName2.Data());
259   TFile::Cp(dsName.Data(),dsName2.Data());
260
261 }