1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
19 #include "AliROCRawAnalysisSelector.h"
23 #include <AliRawEvent.h>
24 #include <AliRawReaderRoot.h>
25 #include <AliRawEventHeaderBase.h>
26 #include <AliTPCRawStream.h>
27 #include <AliTPCParamSR.h>
32 #include <TTimeStamp.h>
34 #include <TPC/AliTPCRawHistograms.h>
37 ClassImp(AliROCRawAnalysisSelector)
39 AliROCRawAnalysisSelector::AliROCRawAnalysisSelector() :
46 // Constructor. Initialization of pointers
49 for (Int_t i=0; i<kTPCSectors; i++)
52 fParam = new AliTPCParamSR;
55 AliROCRawAnalysisSelector::~AliROCRawAnalysisSelector()
62 void AliROCRawAnalysisSelector::SlaveBegin(TTree* tree)
70 void AliROCRawAnalysisSelector::Init(TTree* tree)
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.
83 AliDebug(AliLog::kInfo, "INFO: Tree found");
85 tree->SetBranchAddress("rawevent", &fRawEvent);
89 Bool_t AliROCRawAnalysisSelector::Process(Long64_t entry)
95 fTree->GetTree()->GetEntry(entry);
97 AliDebug(AliLog::kInfo, Form("Processing event %lld", entry));
101 AliDebug(AliLog::kError, "fRawEvent empty");
105 AliRawReaderRoot* rawReader = new AliRawReaderRoot(fRawEvent);
107 const AliRawEventHeaderBase* eventHeader = dynamic_cast<const AliRawEventHeaderBase*> (rawReader->GetEventHeader());
110 eventHeader->Print();
112 UInt_t timeStamp = eventHeader->Get("Timestamp");
113 UInt_t eventType = eventHeader->Get("Type");
115 AliDebug(AliLog::kInfo, Form("Time stamp: %s, event type %d", TTimeStamp(timeStamp).AsString(), eventType));
118 AliTPCRawStream* tpcRawStream = new AliTPCRawStream(rawReader);
120 const Int_t kNIS = fParam->GetNInnerSector();
121 const Int_t kNOS = fParam->GetNOuterSector();
122 const Int_t kNS = kNIS + kNOS;
124 for (Int_t sector = 0; sector < kNS; sector++)
126 AliDebug(AliLog::kInfo, Form("*** Looking at sector %d ***", sector));
129 Int_t nDDLs = 0, indexDDL = 0;
133 nRows = fParam->GetNRowLow();
135 indexDDL = sector * 2;
139 nRows = fParam->GetNRowUp();
141 indexDDL = (sector-kNIS) * 4 + kNIS * 2;
144 // Loas the raw data for corresponding DDLs
146 tpcRawStream->SetOldRCUFormat(kTRUE);
147 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
149 AliDebug(AliLog::kDebug, Form("Selected DDLs %d ... %d", indexDDL, indexDDL+nDDLs-1));
153 while (tpcRawStream->Next())
155 if (tpcRawStream->GetSector() != sector)
157 AliDebug(AliLog::kError, Form("Sector index mismatch ! Expected (%d), but got (%d) !",sector,tpcRawStream->GetSector()));
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()));
165 if (!fHistograms[sector])
167 // not sure if this is still needed, should prevent creation of the histogram in the opened root file
169 fHistograms[sector] = new AliTPCRawHistograms(sector);
172 fHistograms[sector]->FillDigit(tpcRawStream);
179 void AliROCRawAnalysisSelector::SlaveTerminate()
182 for (Int_t i=0; i<kTPCSectors; i++)
184 fOutput->Add(fHistograms[i]);
187 void AliROCRawAnalysisSelector::Terminate()
189 AliDebug(AliLog::kInfo, "Terminate....");
191 // TODO read from output list for PROOF
193 TFile* file = TFile::Open("rocRaw.root", "RECREATE");
195 for (Int_t i=0; i<kTPCSectors; i++)
197 fHistograms[i]->SaveHistograms();
201 for (Int_t i=0; i<kTPCSectors; i++)
203 fHistograms[i]->DrawHistograms();