]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCalibKrTask.cxx
- Small fix for CAF that messed-up user tasks after filters
[u/mrichter/AliRoot.git] / TPC / AliTPCCalibKrTask.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 // The AliTPCCalibKrTask class description (TPC Kr calibration).
18 // The AliTPCCalibKrTask loops over tree of TPC Kr clusters and fills AliTPCCalibKr 
19 // calibration component. 
20 // 
21 // As the input it requires the tree with reconstructed Kr clusters (AliTPCclusterKr objects). 
22 // The AliTPCCalibKr output calibration component contains an array of TH3F histograms which can be stored 
23 // in the ouptut file.
24 //
25
26 /*
27  
28 // Usage example:
29 //
30
31 // -- Load toolkit
32 gSystem->AddIncludePath("-I$ALICE_ROOT/TPC/macros");
33 gROOT->LoadMacro("$ALICE_ROOT/TPC/macros/AliXRDPROOFtoolkit.cxx+");
34 AliXRDPROOFtoolkit tool;
35
36 // -- Make chain of files
37 TChain * chain = tool.MakeChain("KrClusters.txt","Kr","",1000,0);
38
39 // -- Create the analysis manager
40 AliAnalysisManager *mgr = new AliAnalysisManager("testAnalysis");
41
42 // -- Calibration component 
43 AliTPCCalibKr *calibObj = new AliTPCCalibKr;
44
45 // -- Add task
46 AliTPCCalibKrTask *task = new AliTPCCalibKrTask;
47 task->SetInputChain(chain);
48 task->SetTPCCalibKr(calibObj);
49 mgr->AddTask(task);
50
51 // -- Attach input
52 cInput  = mgr->CreateContainer("cInput", TChain::Class(), AliAnalysisManager::kInputContainer);
53 mgr->ConnectInput(task, 0, cInput);
54
55 // -- Attach output
56 cOutput = mgr->CreateContainer("cOutput", TList::Class(), AliAnalysisManager::kOutputContainer,"outHistFile.root");
57 mgr->ConnectOutput(task, 0, cOutput);
58
59 // -- Run analysis
60 mgr->InitAnalysis();
61 mgr->PrintStatus();
62 mgr->StartAnalysis("local", chain);
63
64 */
65
66 //Root includes
67 #include <TH1F.h>
68 #include <TH2F.h>
69 #include <TH3F.h>
70 #include <TString.h>
71 #include <TMath.h>
72 #include <TF1.h>
73 #include <TRandom.h>
74 #include <TDirectory.h>
75 #include <TFile.h>
76 //AliRoot includes
77 #include "AliRawReader.h"
78 #include "AliRawReaderRoot.h"
79 #include "AliRawReaderDate.h"
80 #include "AliTPCRawStream.h"
81 #include "AliTPCCalROC.h"
82 #include "AliTPCCalPad.h"
83 #include "AliTPCROC.h"
84 #include "AliMathBase.h"
85 #include "TTreeStream.h"
86 #include "AliTPCRawStreamFast.h"
87
88 //date
89 #include "event.h"
90
91 //header file
92 #include "AliTPCCalibKr.h"
93 #include "AliTPCCalibKrTask.h"
94
95 Int_t AliTPCCalibKrTask::evtNumber = 0;
96
97 ClassImp(AliTPCCalibKrTask)
98
99 AliTPCCalibKrTask::AliTPCCalibKrTask(const char *name) : 
100   AliAnalysisTask(name,""),
101   fClustKr(0),
102   fTPCCalibKr(0),
103   fTree(0),
104   fOutput(0)
105 {
106   //
107   // default constructor
108   //
109
110   DefineInput(0, TChain::Class());
111   DefineOutput(0, TList::Class());
112 }
113
114 //_____________________________________________________________________
115 AliTPCCalibKrTask::~AliTPCCalibKrTask() 
116 {
117   //
118   // destructor
119   //
120   if(fClustKr) delete fClustKr; fClustKr = 0;
121 }
122
123 //_____________________________________________________________________
124 void AliTPCCalibKrTask::ConnectInputData(Option_t *)
125 {
126   // Connect input data
127   // Called once
128
129   fTree = dynamic_cast<TTree*> (GetInputData(0));
130
131   if(!fTree) { 
132    Printf("ERROR: Could not read chain from input");
133   }
134   else {
135    fTree->SetBranchStatus("*",1); 
136   }
137
138   // set branch address
139   if(!fTree->GetBranch("Cl.")) {
140     Printf("ERROR: Could not get Cl. branch from input");
141   } else {
142    fTree->GetBranch("Cl.")->SetAddress(&fClustKr);
143   }
144 }
145
146 //_____________________________________________________________________
147 void AliTPCCalibKrTask::CreateOutputObjects()
148 {
149   // create object to the output 
150   fOutput = new TList;
151   fOutput->SetOwner();
152
153   if(fTPCCalibKr) fOutput->Add(fTPCCalibKr);
154   //fTPCCalibKr = new AliTPCCalibKr;
155   //if(fTPCCalibKr) fOutput->Add(fTPCCalibKr);
156   else
157      Printf("WARNING: AliTPCCalibKr is not added to the output");
158 }
159
160 //_____________________________________________________________________
161 Bool_t AliTPCCalibKrTask::ReadEntry(Int_t evt)
162 {
163   // read entry 
164   Long64_t centry = fTree->LoadTree(evt);
165   if(centry < 0) return kFALSE;
166
167   if(!fTree->GetBranch("Cl.")) 
168   {
169     Printf("ERROR: Could not get Cl. branch from input");
170         return kFALSE;
171   } else {
172    fTree->GetBranch("Cl.")->SetAddress(&fClustKr);
173   }
174
175   fTree->GetEntry(evt);
176
177 return kTRUE;
178 }
179
180 //_____________________________________________________________________
181 void AliTPCCalibKrTask::Exec(Option_t *)
182 {
183   // Main loop
184   // Called for each event
185   
186   // read entry
187   if(fClustKr) delete fClustKr; fClustKr=0;
188   Bool_t status = ReadEntry(evtNumber);
189   if(status==kTRUE) 
190   {
191           // Process output objects
192       if(fClustKr) fTPCCalibKr->Process(fClustKr);
193   }
194  
195   if( !( evtNumber % 10000) ) {
196     cout << evtNumber << endl; }
197
198   evtNumber++;
199
200   // Post output data.
201   PostData(0, fOutput);
202 }
203
204 //_____________________________________________________________________
205 void AliTPCCalibKrTask::Terminate(Option_t *) 
206 {
207   // Called once at the end of the event loop
208   cout << "Terminate " << endl;
209
210   fOutput = dynamic_cast<TList*> (GetOutputData(0));
211   if (!fOutput) {
212     Printf("ERROR: fOutput not available");
213     return;
214   }
215   
216   fTPCCalibKr = dynamic_cast<AliTPCCalibKr*> (fOutput->FindObject("AliTPCCalibKr"));
217   if (!fTPCCalibKr) {
218     Printf("WARNING: AliTPCCalibKr not available");
219     return;
220   }
221
222   for(int i=0; i<72; ++i) {
223          if( fTPCCalibKr->IsCSide(i) == kTRUE )
224            printf("C side chamber: %d, 3D histo entries: %10.f \n",i,fTPCCalibKr->GetHistoKr(i)->GetEntries());
225
226          if( fTPCCalibKr->IsCSide(i) == kFALSE )
227            printf("A side chamber: %d, 3D histo entries: %10.f \n",i,fTPCCalibKr->GetHistoKr(i)->GetEntries());
228   }
229 }