]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDTrackHistoComponent.cxx
code cleanup (Theo)
[u/mrichter/AliRoot.git] / HLT / TRD / AliHLTTRDTrackHistoComponent.cxx
1 //**************************************************************************
2 //* This file is property of and copyright by the ALICE HLT Project        * 
3 //* ALICE Experiment at CERN, All rights reserved.                         *
4 //*                                                                        *
5 //* Primary Authors: Sylwester Radomski radomski@physi.uni-heidelberg.de    *
6 //*                  for The ALICE HLT Project.                            *
7 //*                                                                        *
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 //**************************************************************************
16
17 /** @file   AliHLTTRDTrackHistoComponent.cxx
18     @author Raphaelle and Theodor
19     @brief  Component for ploting charge in clusters
20 */
21
22 #if __GNUC__>= 3
23 using namespace std;
24 #endif
25
26 #include <time.h>
27
28 #include "AliHLTTRDTrackHistoComponent.h"
29 #include "AliHLTTRDDefinitions.h"
30 #include "AliCDBEntry.h"
31 #include "AliCDBManager.h"
32 #include <TFile.h>
33 #include <TString.h>
34 #include "TObjString.h"
35 #include "TClonesArray.h"
36 #include "TTimeStamp.h"
37 #include "AliHLTTRDUtils.h"
38 #include "TH1F.h"
39 #include "AliTRDcluster.h"
40 #include "AliTRDtrackV1.h"
41 #include "AliTRDseedV1.h"
42
43 //#include "AliHLTTRD.h"
44 //#include <stdlib.h>
45 //#include <cerrno>
46
47 /** ROOT macro for the implementation of ROOT specific class methods */
48 ClassImp(AliHLTTRDTrackHistoComponent)
49
50 AliHLTTRDTrackHistoComponent::AliHLTTRDTrackHistoComponent()
51 : AliHLTProcessor(),
52   fOutputSize(100000),
53   fTracksArray(NULL),
54   fClPerTrkl(NULL),
55   fTrklPerTrk(NULL)
56 {
57   // see header file for class documentation
58   // or
59   // refer to README to build package
60   // or
61   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
62
63 }
64
65 AliHLTTRDTrackHistoComponent::~AliHLTTRDTrackHistoComponent()
66 {
67   // see header file for class documentation
68 }
69
70 // Public functions to implement AliHLTComponent's interface.
71 // These functions are required for the registration process
72
73 const char* AliHLTTRDTrackHistoComponent::GetComponentID()
74 {
75   // see header file for class documentation
76   
77   return "TRDTrackHisto";
78 }
79
80 void AliHLTTRDTrackHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
81 {
82   // see header file for class documentation
83   list.clear();
84   list.push_back( AliHLTTRDDefinitions::fgkTracksDataType );
85 }
86
87 AliHLTComponentDataType AliHLTTRDTrackHistoComponent::GetOutputDataType()
88 {
89   // see header file for class documentation
90   return kAliHLTDataTypeHistogram  | kAliHLTDataOriginTRD;
91
92 }
93
94 void AliHLTTRDTrackHistoComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
95 {
96   // see header file for class documentation
97   constBase = fOutputSize;
98   inputMultiplier = 0;
99 }
100
101 AliHLTComponent* AliHLTTRDTrackHistoComponent::Spawn()
102 {
103   // see header file for class documentation
104   return new AliHLTTRDTrackHistoComponent;
105 }
106
107 int AliHLTTRDTrackHistoComponent::DoInit(int argc, const char** argv)
108 {
109   // Initialize histograms
110   int iResult=0;
111   
112   TString configuration="";
113   TString argument="";
114   for (int i=0; i<argc && iResult>=0; i++) {
115     argument=argv[i];
116     if (!configuration.IsNull()) configuration+=" ";
117     configuration+=argument;
118   }
119
120   if (!configuration.IsNull()) {
121     iResult=Configure(configuration.Data());
122   } 
123
124   fTracksArray = new TClonesArray("AliTRDtrackV1");
125
126   fClPerTrkl = new TH1F("fClPerTrkl","Clusters per Tracklet", AliTRDseedV1::kNtb, -0.5, AliTRDseedV1::kNtb - 0.5);
127   fTrklPerTrk = new TH1F("fTrklPerTrk","Tracklets per Track", 7, -0.5, 6.5);
128   
129   return 0;
130 }
131   
132 int AliHLTTRDTrackHistoComponent::DoDeinit()
133 {
134   // see header file for class documentation
135
136   fTracksArray->Delete();
137   delete fTracksArray;
138
139   // delete histograms
140   if (fClPerTrkl) delete fClPerTrkl;
141   if (fTrklPerTrk) delete fTrklPerTrk;
142   
143   return 0;
144 }
145
146 int AliHLTTRDTrackHistoComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/,
147                                             AliHLTComponentTriggerData& /*trigData*/)
148 {
149
150   // if (GetFirstInputBlock(kAliHLTDataTypeSOR)) return 0;
151   // else if (GetFirstInputBlock(kAliHLTDataTypeEOR))
152   //   {
153   //     TString fileName="/tmp/TracksHistoDump_run";
154   //     fileName+=AliCDBManager::Instance()->GetRun();
155   //     fileName+=".root";
156   //     HLTInfo("Dumping Histogram file to %s",fileName.Data());
157   //     TFile* file = TFile::Open(fileName, "RECREATE");
158   //     fClPerTrkl->Write();
159   //     fTrklPerTrk->Write();
160   //     file->Close();
161   //     HLTInfo("Histogram file dumped");
162   //     return 0;
163   //   }
164
165   if (GetFirstInputBlock(kAliHLTDataTypeSOR) || GetFirstInputBlock(kAliHLTDataTypeEOR)) return 0;
166
167   const AliHLTComponentBlockData* iter = NULL;
168   Bool_t gotData = kFALSE;
169   
170   for(iter = GetFirstInputBlock(AliHLTTRDDefinitions::fgkTracksDataType); 
171         iter != NULL; iter = GetNextInputBlock() ) {
172     
173     AliHLTTRDUtils::ReadTracks(fTracksArray, iter->fPtr, iter->fSize);
174     HLTDebug("TClonesArray of clusters: nbEntries = %i", fTracksArray->GetEntriesFast());
175     gotData=kTRUE;
176   }
177   
178   if(!gotData) return 0;
179   
180   AliTRDtrackV1 *trk;
181   
182   // loop over tracks
183   for(int i=0;i<fTracksArray->GetEntriesFast();i++) {
184     trk=(AliTRDtrackV1*)fTracksArray->At(i);
185     Int_t nrOfTrkls=0;
186     for(int seedNr=0; seedNr<6; seedNr++){
187       AliTRDseedV1* seed = trk->GetTracklet(seedNr);
188       if(!seed)continue;
189       nrOfTrkls++;
190       Int_t nrOfCls=0;
191       for(int clsNr=0; clsNr<AliTRDseedV1::kNtb; clsNr++)
192         if(seed->GetClusters(clsNr))nrOfCls++;
193       fClPerTrkl->Fill(nrOfCls);
194     }
195     fTrklPerTrk->Fill(nrOfTrkls);
196   }
197   
198   fTracksArray->Delete();
199   
200   PushBack((TObject*)fClPerTrkl, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);   
201   PushBack((TObject*)fTrklPerTrk, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);  
202   
203   return 0;
204 }
205
206 int AliHLTTRDTrackHistoComponent::Configure(const char* arguments){
207   int iResult=0;
208   if (!arguments) return iResult;
209   
210   TString allArgs=arguments;
211   TString argument;
212   int bMissingParam=0;
213
214   TObjArray* pTokens=allArgs.Tokenize(" ");
215   if (pTokens) {
216     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
217       argument=((TObjString*)pTokens->At(i))->GetString();
218       if (argument.IsNull()) continue;
219       
220       if (argument.CompareTo("output_size")==0) {
221         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
222         HLTInfo("Setting output size to: %s", ((TObjString*)pTokens->At(i))->GetString().Data());
223         fOutputSize=((TObjString*)pTokens->At(i))->GetString().Atoi();
224         continue;
225       } 
226       if (argument.CompareTo("-everyNevent")==0) {
227         if ((bMissingParam=(++i>=pTokens->GetEntries()))) break;
228         HLTInfo("Option -everyNevent depreceated");
229         continue;
230       } 
231       else {
232         HLTError("unknown argument: %s", argument.Data());
233         iResult=-EINVAL;
234         break;
235       }
236     }
237     delete pTokens;
238   }
239   if (bMissingParam) {
240     HLTError("missing parameter for argument %s", argument.Data());
241     iResult=-EINVAL;
242   }
243   return iResult;
244 }