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