]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCalibKrTask.cxx
TPC DQM update (Ruben, Peter C.)
[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 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 "AliTPCRawStream.h"
89 #include "AliTPCCalROC.h"
90 #include "AliTPCCalPad.h"
91 #include "AliTPCROC.h"
92 #include "AliMathBase.h"
93 #include "TTreeStream.h"
94
95 //date
96 #include "event.h"
97
98 //header file
99 #include "AliTPCCalibKr.h"
100 #include "AliTPCCalibKrTask.h"
101
102 Int_t AliTPCCalibKrTask::fEvtNumber = 0;
103
104 ClassImp(AliTPCCalibKrTask)
105
106 AliTPCCalibKrTask::AliTPCCalibKrTask(const char *name) : 
107   AliAnalysisTask(name,""),
108   fClustKr(0),
109   fTPCCalibKr(0),
110   fTree(0),
111   fOutput(0)
112 {
113   //
114   // default constructor
115   //
116
117   DefineInput(0, TChain::Class());
118   DefineOutput(0, TList::Class());
119 }
120
121 //_____________________________________________________________________
122 AliTPCCalibKrTask::~AliTPCCalibKrTask() 
123 {
124   //
125   // destructor
126   //
127   if(fOutput) fOutput->Delete();
128   delete fOutput; fOutput = 0;
129 }
130
131 //_____________________________________________________________________
132 void AliTPCCalibKrTask::ConnectInputData(Option_t *)
133 {
134   // Connect input data
135   // Called once
136
137   fTree = dynamic_cast<TTree*> (GetInputData(0));
138
139   if(!fTree) { 
140    Printf("ERROR: Could not read chain from input");
141    return;
142   }
143   else {
144    fTree->SetBranchStatus("*",1); 
145    fTree->SetBranchStatus("Cl.fCluster",0); 
146   }
147
148   // set branch address
149   if(!fTree->GetBranch("Cl.")) {
150     Printf("ERROR: Could not get Cl. branch from input");
151   } else {
152    fTree->GetBranch("Cl.")->SetAddress(&fClustKr);
153   }
154 }
155
156 //_____________________________________________________________________
157 void AliTPCCalibKrTask::CreateOutputObjects()
158 {
159   // create object to the output 
160   fOutput = new TList;
161   fOutput->SetOwner(); // is owner of the fTPCCalibKr objects
162
163   if(fTPCCalibKr) fOutput->Add(fTPCCalibKr);
164   //fTPCCalibKr = new AliTPCCalibKr;
165   //if(fTPCCalibKr) fOutput->Add(fTPCCalibKr);
166   else
167      Printf("WARNING: AliTPCCalibKr is not added to the output");
168 }
169
170 //_____________________________________________________________________
171 Bool_t AliTPCCalibKrTask::ReadEntry(Int_t evt)
172 {
173   // read entry 
174   Long64_t centry = fTree->LoadTree(evt);
175   if(centry < 0) return kFALSE;
176
177   if(!fTree->GetBranch("Cl.")) 
178   {
179     Printf("ERROR: Could not get Cl. branch from input");
180         return kFALSE;
181   } else {
182    fTree->GetBranch("Cl.")->SetAddress(&fClustKr);
183   }
184
185   fTree->GetEntry(evt);
186
187 return kTRUE;
188 }
189
190 //_____________________________________________________________________
191 void AliTPCCalibKrTask::Exec(Option_t *)
192 {
193   // Main loop
194   // Called for each event
195   
196   // read entry
197   if(fClustKr) delete fClustKr; fClustKr=0;
198   Bool_t status = ReadEntry(fEvtNumber);
199   if(status==kTRUE) 
200   {
201           // Process output objects
202       if(fClustKr) fTPCCalibKr->Process(fClustKr);
203   }
204  
205   if( !( fEvtNumber % 100000) ) {
206     cout << fEvtNumber << endl; }
207
208   fEvtNumber++;
209
210   // Post output data.
211   PostData(0, fOutput);
212 }
213
214 //_____________________________________________________________________
215 void AliTPCCalibKrTask::Terminate(Option_t *) 
216 {
217   // Called once at the end of the event loop
218   cout << "Terminate " << endl;
219
220   fOutput = dynamic_cast<TList*> (GetOutputData(0));
221   if (!fOutput) {
222     Printf("ERROR: fOutput not available");
223     return;
224   }
225   
226   fTPCCalibKr = dynamic_cast<AliTPCCalibKr*> (fOutput->FindObject("AliTPCCalibKr"));
227   if (!fTPCCalibKr) {
228     Printf("WARNING: AliTPCCalibKr not available");
229     return;
230   }
231
232   for(int i=0; i<72; ++i) {
233          if( fTPCCalibKr->IsCSide(i) == kTRUE )
234            printf("C side chamber: %d, 3D histo entries: %10.f \n",i,fTPCCalibKr->GetHistoKr(i)->GetEntries());
235
236          if( fTPCCalibKr->IsCSide(i) == kFALSE )
237            printf("A side chamber: %d, 3D histo entries: %10.f \n",i,fTPCCalibKr->GetHistoKr(i)->GetEntries());
238   }
239 }