1 //**************************************************************************
2 //* This file is property of and copyright by the ALICE HLT Project *
3 //* ALICE Experiment at CERN, All rights reserved. *
5 //* Primary Authors: Sylwester Radomski radomski@physi.uni-heidelberg.de *
6 //* for The ALICE HLT Project. *
8 //* Permission to use, copy, modify and distribute this software and its *
9 //* documentation strictly for non-commercial purposes is hereby granted *
10 //* without fee, provided that the above copyright notice appears in all *
11 //* copies and that both the copyright notice and this permission notice *
12 //* appear in the supporting documentation. The authors make no claims *
13 //* about the suitability of this software for any purpose. It is *
14 //* provided "as is" without express or implied warranty. *
15 //**************************************************************************
17 /** @file AliHLTTRDClusterHistoComponent.cxx
18 @author Sylwester Radomski
19 @brief Component for ploting charge in clusters
28 #include "TObjString.h"
29 #include "TClonesArray.h"
32 #include "AliHLTTRDClusterHistoComponent.h"
33 #include "AliHLTTRDDefinitions.h"
34 #include "AliTRDcluster.h"
35 #include "AliCDBEntry.h"
36 #include "AliCDBManager.h"
37 #include "AliHLTTRDUtils.h"
39 //#include "AliHLTTRD.h"
43 /** ROOT macro for the implementation of ROOT specific class methods */
44 ClassImp(AliHLTTRDClusterHistoComponent)
46 AliHLTTRDClusterHistoComponent::AliHLTTRDClusterHistoComponent()
58 // see header file for class documentation
60 // refer to README to build package
62 // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
66 AliHLTTRDClusterHistoComponent::~AliHLTTRDClusterHistoComponent()
68 // see header file for class documentation
71 // Public functions to implement AliHLTComponent's interface.
72 // These functions are required for the registration process
74 const char* AliHLTTRDClusterHistoComponent::GetComponentID()
76 // see header file for class documentation
78 return "TRDClusterHisto";
81 void AliHLTTRDClusterHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
83 // see header file for class documentation
85 list.push_back( AliHLTTRDDefinitions::fgkClusterDataType );
88 AliHLTComponentDataType AliHLTTRDClusterHistoComponent::GetOutputDataType()
90 // see header file for class documentation
91 return kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD;
95 void AliHLTTRDClusterHistoComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
97 // see header file for class documentation
98 constBase = fOutputSize;
102 AliHLTComponent* AliHLTTRDClusterHistoComponent::Spawn()
104 // see header file for class documentation
105 return new AliHLTTRDClusterHistoComponent;
108 int AliHLTTRDClusterHistoComponent::DoInit(int argc, const char** argv)
110 // Initialize histograms
113 TString configuration="";
115 for (int i=0; i<argc && iResult>=0; i++) {
117 if (!configuration.IsNull()) configuration+=" ";
118 configuration+=argument;
121 if (!configuration.IsNull()) {
122 iResult=Configure(configuration.Data());
125 fClusterArray = new TClonesArray("AliTRDcluster");
127 fNClsDet = new TH1F("trdClsDet", ";detector", 540, -0.5, 539.5);
128 fClsAmp = new TH1F("trdClsAmp", ";amplitude", 200, -0.5, 1999.5);
129 fClsAmpDrift = new TH1F("trdClsAmpDrift", ";amplitude", 200, -0.5, 199.5) ;
130 fClsTB = new TH1F("trdClsTB", ";time bin", 35, -0.5, 34.5);
131 fClsAmpDist = new TH1F("trdClsAmpDist", "mean amplitude", 200, 0, 1000);
132 fSClsDist = new TH1F("sclsdist", "Super cluster spectrum", 200, 0, 8000);
133 fNScls = new TH1F("nscls", "No. of Kr clusters per event", 540, 0, 540);
135 for(int i=0; i<540; i++) {
136 fClsAmpDriftDet[i] = new TH1F(Form("trdClsDriftDet_%d",i), "", 200, -0.5, 199.5);
142 int AliHLTTRDClusterHistoComponent::DoDeinit()
144 // see header file for class documentation
146 fClusterArray->Delete();
147 delete fClusterArray;
150 if (fNClsDet) delete fNClsDet;
151 if (fClsAmp) delete fClsAmp;
152 if (fClsAmpDrift) delete fClsAmpDrift;
153 if (fClsTB) delete fClsTB;
154 if (fClsAmpDist) delete fClsAmpDist;
155 if (fSClsDist) delete fSClsDist;
156 if (fNScls) delete fNScls;
158 for(int i=0; i<540; i++){
159 if (fClsAmpDriftDet[i]) delete fClsAmpDriftDet[i];
165 int AliHLTTRDClusterHistoComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/,
166 AliHLTComponentTriggerData& /*trigData*/)
169 // if (GetFirstInputBlock(kAliHLTDataTypeSOR)) return 0;
170 // else if (GetFirstInputBlock(kAliHLTDataTypeEOR))
172 // TString fileName="/tmp/ClusterHistoDump_run";
173 // fileName+=AliCDBManager::Instance()->GetRun();
174 // fileName+=".root";
175 // HLTInfo("Dumping Histogram file to %s",fileName.Data());
176 // TFile* file = TFile::Open(fileName, "RECREATE");
177 // fNClsDet->Write();
179 // fClsAmpDrift->Write();
181 // fClsAmpDist->Write();
182 // fSClsDist->Write();
185 // HLTInfo("Histogram file dumped");
189 if (GetFirstInputBlock(kAliHLTDataTypeSOR) || GetFirstInputBlock(kAliHLTDataTypeEOR)) return 0;
191 const AliHLTComponentBlockData* iter = NULL;
192 Bool_t gotData = kFALSE;
193 AliHLTUInt32_t spec = 0;
195 for ( iter = GetFirstInputBlock(AliHLTTRDDefinitions::fgkClusterDataType);
196 iter != NULL; iter = GetNextInputBlock() ) {
198 HLTDebug("We get the right data type: Block Ptr: 0x%x; Block Size: %i",
199 iter->fPtr, iter->fSize);
201 AliHLTTRDUtils::ReadClusters(fClusterArray, iter->fPtr, iter->fSize);
203 spec = iter->fSpecification;
206 if(!gotData) return 0;
208 HLTDebug("TClonesArray of clusters: nbEntries = %i", fClusterArray->GetEntriesFast());
210 Float_t sClusterCharge[540] = { 0 };
213 // loop over clusters
214 for(int i=0;i<fClusterArray->GetEntriesFast();i++) {
216 cls=(AliTRDcluster*)fClusterArray->At(i);
218 fNClsDet->Fill(cls->GetDetector());
219 fClsAmp->Fill(cls->GetQ());
221 int tb = cls->GetPadTime();
223 if (tb > 5 && tb <25){
224 fClsAmpDrift->Fill(cls->GetQ());
227 //fClsAmpDriftDet[cls->GetDetector()]->Fill(cls->GetQ());
229 Int_t det = cls->GetDetector();
230 sClusterCharge[det] += cls->GetQ();
234 fClusterArray->Delete();
236 //fClsAmpDist->Reset();
237 //Int_t nSClusters = 0;
238 for(int det=0; det<540; det++) {
239 // if (fClsAmpDriftDet[det]->GetSum() > 0)
240 // fClsAmpDist->Fill(fClsAmpDriftDet[det]->GetMean());
241 if(sClusterCharge[det])
242 fSClsDist->Fill(sClusterCharge[det]);
245 //fNScls->Fill(nSClusters);
247 PushBack((TObject*)fNClsDet, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, spec);
248 PushBack((TObject*)fClsAmp, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, spec);
249 PushBack((TObject*)fClsAmpDrift, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, spec);
250 PushBack((TObject*)fClsTB, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, spec);
251 //PushBack((TObject*)fClsAmpDist, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, spec);
252 //PushBack((TObject*)fNScls, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, spec);
253 PushBack((TObject*)fSClsDist, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, spec);
258 int AliHLTTRDClusterHistoComponent::Configure(const char* arguments){
260 if (!arguments) return iResult;
262 TString allArgs=arguments;
266 TObjArray* pTokens=allArgs.Tokenize(" ");
268 for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
269 argument=((TObjString*)pTokens->At(i))->GetString();
270 if (argument.IsNull()) continue;
272 if (argument.CompareTo("output_size")==0) {
273 if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
274 HLTInfo("Setting output size to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
275 fOutputSize=((TObjString*)pTokens->At(i))->GetString().Atoi();
279 HLTError("unknown argument: %s", argument.Data());
287 HLTError("missing parameter for argument %s", argument.Data());