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