]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCAnalysisTaskcalib.cxx
Updated alignment:
[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   RegisterDebugOutput();
157   
158 }
159
160
161 void AliTPCAnalysisTaskcalib::Process(AliESDEvent *event) {
162   //
163   // Process ESD event
164   //
165   AliTPCcalibBase *job=0;
166   Int_t njobs = fCalibJobs->GetEntriesFast();
167   for (Int_t i=0;i<njobs;i++){
168     job = (AliTPCcalibBase*)fCalibJobs->UncheckedAt(i);
169     if (job) {
170       job->UpdateEventInfo(event);
171       if (job->AcceptTrigger())
172         job->Process(event);
173     }
174   }
175 }
176
177 void AliTPCAnalysisTaskcalib::Process(AliTPCseed *track) {
178   //
179   // Process TPC track
180   //
181   AliTPCcalibBase *job=0;
182   Int_t njobs = fCalibJobs->GetEntriesFast();
183   for (Int_t i=0;i<njobs;i++){
184     job = (AliTPCcalibBase*)fCalibJobs->UncheckedAt(i);
185     if (job)  
186       if (job->AcceptTrigger())
187         job->Process(track);
188   }
189 }
190
191 void AliTPCAnalysisTaskcalib::Process(AliESDtrack *track, Int_t run) {
192   //
193   // Process ESD track
194   //
195   AliTPCcalibBase *job=0;
196   Int_t njobs = fCalibJobs->GetEntriesFast();
197   for (Int_t i=0;i<njobs;i++){
198     job = (AliTPCcalibBase*)fCalibJobs->UncheckedAt(i);
199     if (job) 
200       if (job->AcceptTrigger())
201         job->Process(track,run);
202   }
203 }
204
205 Long64_t AliTPCAnalysisTaskcalib::Merge(TCollection *li) {
206   TIterator *i=fCalibJobs->MakeIterator();
207   AliTPCcalibBase *job;
208   Long64_t n=0;
209   while ((job=dynamic_cast<AliTPCcalibBase*>(i->Next())))
210     n+=job->Merge(li);
211   return n;
212 }
213
214 void AliTPCAnalysisTaskcalib::Analyze() {
215   //
216   // Analyze the content of the task
217   //
218   AliTPCcalibBase *job=0;
219   Int_t njobs = fCalibJobs->GetEntriesFast();
220   for (Int_t i=0;i<njobs;i++){
221     job = (AliTPCcalibBase*)fCalibJobs->UncheckedAt(i);
222     if (job) job->Analyze();
223   }
224 }
225
226
227 void AliTPCAnalysisTaskcalib::RegisterDebugOutput(){
228   //
229   //
230   //
231   AliTPCcalibBase *job=0;
232   Int_t njobs = fCalibJobs->GetEntriesFast();
233   for (Int_t i=0;i<njobs;i++){
234     job = (AliTPCcalibBase*)fCalibJobs->UncheckedAt(i);
235     if (job) job->RegisterDebugOutput(fDebugOutputPath.Data());
236   }
237   TString dsName=GetName();
238   dsName+=".root";
239   TFile fff(dsName.Data(),"recreate");
240   fCalibJobs->Write("TPCCalib",TObject::kSingleKey);
241   fff.Close();
242   //
243   // store  - copy debug output to the destination position
244   // currently ONLY for local copy
245   TString dsName2=fDebugOutputPath.Data();
246   gSystem->MakeDirectory(dsName2.Data());
247   dsName2+=gSystem->HostName();
248   gSystem->MakeDirectory(dsName2.Data());
249   dsName2+="/";
250   TTimeStamp s;
251   dsName2+=Int_t(s.GetNanoSec());
252   dsName2+="/";
253   gSystem->MakeDirectory(dsName2.Data());
254   dsName2+=dsName;
255   AliInfo(Form("copy %s\t%s\n",dsName.Data(),dsName2.Data()));
256   printf("copy %s\t%s\n",dsName.Data(),dsName2.Data());
257   TFile::Cp(dsName.Data(),dsName2.Data());
258
259 }