]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG0/TPC/AliROCRawAnalysisSelector.cxx
added more histograms. changed track cut method. introduced func to save hist for...
[u/mrichter/AliRoot.git] / PWG0 / TPC / AliROCRawAnalysisSelector.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 /* $Id$ */
17
18
19 #include "AliROCRawAnalysisSelector.h"
20
21 #include <AliLog.h>
22
23 #include <AliRawEvent.h>
24 #include <AliRawReaderRoot.h>
25 #include <AliRawEventHeaderBase.h>
26 #include <AliTPCRawStream.h>
27 #include <AliTPCParamSR.h>
28
29 #include <TFile.h>
30 #include <TTree.h>
31 #include <TCanvas.h>
32 #include <TTimeStamp.h>
33
34 #include <TPC/AliTPCRawHistograms.h>
35
36
37 ClassImp(AliROCRawAnalysisSelector)
38
39 AliROCRawAnalysisSelector::AliROCRawAnalysisSelector() :
40   TSelector(),
41   fRawEvent(0),
42   fTree(0),
43   fParam(0)  
44 {
45   //
46   // Constructor. Initialization of pointers
47   //
48
49   for (Int_t i=0; i<kTPCSectors; i++)
50     fHistograms[i] = 0;
51
52   fParam = new AliTPCParamSR;
53 }
54
55 AliROCRawAnalysisSelector::~AliROCRawAnalysisSelector()
56 {
57   //
58   // Destructor
59   //
60 }
61
62 void AliROCRawAnalysisSelector::SlaveBegin(TTree* tree)
63 {
64   //
65
66   if (tree != 0)
67     Init(tree);
68
69
70 void AliROCRawAnalysisSelector::Init(TTree* tree)
71 {
72   // The Init() function is called when the selector needs to initialize
73   // a new tree or chain. Typically here the branch addresses of the tree
74   // will be set. It is normaly not necessary to make changes to the
75   // generated code, but the routine can be extended by the user if needed.
76   // Init() will be called many times when running with PROOF.
77
78   fTree = tree;
79
80   // Set branch address
81   if (tree) 
82   {
83     AliDebug(AliLog::kInfo, "INFO: Tree found");
84
85     tree->SetBranchAddress("rawevent", &fRawEvent);
86   }
87 }
88
89 Bool_t AliROCRawAnalysisSelector::Process(Long64_t entry)
90 {
91   //
92   //
93   //
94
95   fTree->GetTree()->GetEntry(entry);
96   
97   AliDebug(AliLog::kInfo, Form("Processing event %lld", entry));
98   
99   if (!fRawEvent)
100   {
101     AliDebug(AliLog::kError, "fRawEvent empty");
102     return kFALSE;
103   }
104   
105   AliRawReaderRoot* rawReader = new AliRawReaderRoot(fRawEvent);
106
107   const AliRawEventHeaderBase* eventHeader = dynamic_cast<const AliRawEventHeaderBase*> (rawReader->GetEventHeader());
108   if (eventHeader) 
109   {
110     eventHeader->Print();
111     
112     UInt_t timeStamp = eventHeader->Get("Timestamp");
113     UInt_t eventType = eventHeader->Get("Type");
114     
115     AliDebug(AliLog::kInfo, Form("Time stamp: %s, event type %d", TTimeStamp(timeStamp).AsString(), eventType));
116   }           
117   
118   AliTPCRawStream* tpcRawStream = new AliTPCRawStream(rawReader);
119      
120   const Int_t kNIS = fParam->GetNInnerSector();
121   const Int_t kNOS = fParam->GetNOuterSector();
122   const Int_t kNS = kNIS + kNOS;
123   
124   for (Int_t sector = 0; sector < kNS; sector++) 
125   {
126     AliDebug(AliLog::kInfo, Form("*** Looking at sector %d ***", sector));
127             
128     Int_t nRows = 0;
129     Int_t nDDLs = 0, indexDDL = 0;
130     
131     if (sector < kNIS) 
132     {
133       nRows = fParam->GetNRowLow();
134       nDDLs = 2;
135       indexDDL = sector * 2;
136     }
137     else 
138     {
139       nRows = fParam->GetNRowUp();
140       nDDLs = 4;
141       indexDDL = (sector-kNIS) * 4 + kNIS * 2;
142     }
143     
144     // Loas the raw data for corresponding DDLs
145     rawReader->Reset();
146     tpcRawStream->SetOldRCUFormat(kTRUE);
147     rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
148     
149     AliDebug(AliLog::kDebug, Form("Selected DDLs %d ... %d", indexDDL, indexDDL+nDDLs-1));
150     
151     Int_t count = 0;
152     
153     while (tpcRawStream->Next())
154     {
155           if (tpcRawStream->GetSector() != sector)
156       {
157             AliDebug(AliLog::kError, Form("Sector index mismatch ! Expected (%d), but got (%d) !",sector,tpcRawStream->GetSector()));
158             return kFALSE;
159       }
160       
161       if ((count++ % 100000) == 0)
162             AliDebug(AliLog::kDebug, Form("Found %d. digit in sector %d: row %d, pad %d, time %d, signal %d", count, 
163             tpcRawStream->GetSector(), tpcRawStream->GetRow(), tpcRawStream->GetPad(), tpcRawStream->GetTime(), tpcRawStream->GetSignal()));
164
165       if (!fHistograms[sector])
166       {
167         // not sure if this is still needed, should prevent creation of the histogram in the opened root file
168         gROOT->cd();
169         fHistograms[sector] = new AliTPCRawHistograms(sector);
170       }
171       
172       fHistograms[sector]->FillDigit(tpcRawStream);
173     }
174   }
175    
176   return kTRUE;
177 }
178
179 void AliROCRawAnalysisSelector::SlaveTerminate()
180 {
181   //
182   for (Int_t i=0; i<kTPCSectors; i++)
183    if (fHistograms[i])
184      fOutput->Add(fHistograms[i]);
185
186
187 void AliROCRawAnalysisSelector::Terminate()
188 {
189   AliDebug(AliLog::kInfo, "Terminate....");
190
191   // TODO read from output list for PROOF
192     
193   TFile* file = TFile::Open("rocRaw.root", "RECREATE");
194   
195   for (Int_t i=0; i<kTPCSectors; i++)
196     if (fHistograms[i])
197       fHistograms[i]->SaveHistograms();
198
199   file->Close();
200
201   for (Int_t i=0; i<kTPCSectors; i++)
202     if (fHistograms[i])
203       fHistograms[i]->DrawHistograms();
204