]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/AliHLTTRDClusterHistoComponent.cxx
- minor fixes for the TRD histogram component (Theodor)
[u/mrichter/AliRoot.git] / HLT / TRD / AliHLTTRDClusterHistoComponent.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   AliHLTTRDClusterHistoComponent.cxx
18     @author Sylwester Radomski
19     @brief  Component for ploting charge in clusters
20 */
21
22 #if __GNUC__>= 3
23 using namespace std;
24 #endif
25
26 #include "AliHLTTRDClusterHistoComponent.h"
27 #include "AliHLTTRDDefinitions.h"
28 #include "AliHLTTRDCluster.h"
29 #include "AliTRDcluster.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 "AliHLTTRDUtils.h"
37
38 //#include "AliHLTTRD.h"
39 //#include <stdlib.h>
40 //#include <cerrno>
41
42 /** ROOT macro for the implementation of ROOT specific class methods */
43 ClassImp(AliHLTTRDClusterHistoComponent)
44
45 AliHLTTRDClusterHistoComponent::AliHLTTRDClusterHistoComponent()
46 : fClusterArray(NULL),
47   fNClsDet(NULL),
48   fClsAmp(NULL),
49   fClsAmpDrift(NULL),
50   fClsTB(NULL),
51   fClsAmpDist(NULL),
52   fSClsDist(NULL)
53 {
54   // see header file for class documentation
55   // or
56   // refer to README to build package
57   // or
58   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
59
60 }
61
62 AliHLTTRDClusterHistoComponent::~AliHLTTRDClusterHistoComponent()
63 {
64   // see header file for class documentation
65 }
66
67 // Public functions to implement AliHLTComponent's interface.
68 // These functions are required for the registration process
69
70 const char* AliHLTTRDClusterHistoComponent::GetComponentID()
71 {
72   // see header file for class documentation
73   
74   return "TRDClusterHisto";
75 }
76
77 void AliHLTTRDClusterHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
78 {
79   // see header file for class documentation
80   list.clear();
81   list.push_back( AliHLTTRDDefinitions::fgkClusterDataType );
82 }
83
84 AliHLTComponentDataType AliHLTTRDClusterHistoComponent::GetOutputDataType()
85 {
86   // see header file for class documentation
87   return kAliHLTDataTypeHistogram  | kAliHLTDataOriginTRD;
88
89 }
90
91 void AliHLTTRDClusterHistoComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
92 {
93   // see header file for class documentation
94   constBase = 5000;
95   inputMultiplier = 3;
96 }
97
98 AliHLTComponent* AliHLTTRDClusterHistoComponent::Spawn()
99 {
100   // see header file for class documentation
101   return new AliHLTTRDClusterHistoComponent;
102 }
103
104 int AliHLTTRDClusterHistoComponent::DoInit(int /*argc*/, const char** /*argv*/ )
105 {
106   // Initialize histograms
107
108   fClusterArray = new TClonesArray("AliTRDcluster");
109
110   fNClsDet = new TH1D("trdClsDet", ";detector", 540, -0.5, 539.5);
111   fClsAmp  = new TH1D("trdClsAmp", ";amplitude", 200, -0.5, 1999.5);
112   fClsAmpDrift = new TH1D("trdClsAmpDrift", ";amplitude", 200, -0.5, 199.5) ;
113   fClsTB = new TH1D("trdClsTB", ";time bin", 35, -0.5, 34.5);
114   fClsAmpDist = new TH1D("trdClsAmpDist", "mean amplitude", 200, 0, 1000);
115   fSClsDist = new TH1D("sclsdist", "Super cluster spectrum", 200, 0, 8000);
116
117   for(int i=0; i<540; i++)
118     fClsAmpDriftDet[i] = new TH1D(Form("trdClsDriftDet_%d",i), "", 200, -0.5, 199.5);
119     
120     return 0;
121 }
122   
123 int AliHLTTRDClusterHistoComponent::DoDeinit()
124 {
125   // see header file for class documentation
126
127   fClusterArray->Delete();
128   delete fClusterArray;
129
130   // delete histograms
131   if (fNClsDet) delete fNClsDet;
132   if (fClsAmp) delete fClsAmp;
133   if (fClsAmpDrift) delete fClsAmpDrift;
134   if (fClsTB) delete fClsTB;
135   if (fClsAmpDist) delete fClsAmpDist;
136   if (fSClsDist) delete fSClsDist;
137
138   for(int i=0; i<540; i++)
139     if (fClsAmpDriftDet[i]) delete fClsAmpDriftDet[i];
140
141   return 0;
142 }
143
144 int AliHLTTRDClusterHistoComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/, 
145                                             AliHLTComponentTriggerData& /*trigData*/)
146 {
147
148   if (GetFirstInputBlock(kAliHLTDataTypeSOR) || GetFirstInputBlock(kAliHLTDataTypeEOR)) return 0;
149   
150   const AliHLTComponentBlockData* iter = NULL;
151   
152   Float_t sClusterCharge[540] = { 0 };
153
154   for ( iter = GetFirstInputBlock(AliHLTTRDDefinitions::fgkClusterDataType); 
155         iter != NULL; iter = GetNextInputBlock() ) {
156
157     HLTDebug("We get the right data type: Block Ptr: 0x%x; Block Size: %i",
158              iter->fPtr, iter->fSize);
159
160     AliHLTTRDUtils::ReadClusters(fClusterArray, iter->fPtr, iter->fSize);
161     HLTDebug("TClonesArray of clusters: nbEntries = %i", fClusterArray->GetEntriesFast());
162
163     AliTRDcluster *cls;
164
165     // loop over clusters 
166     for(int i=0;i<fClusterArray->GetEntriesFast();i++) {
167
168       cls=(AliTRDcluster*)fClusterArray->At(i);
169       
170       fNClsDet->Fill(cls->GetDetector());
171       sClusterCharge[cls->GetDetector()] += cls->GetQ();
172       fClsAmp->Fill(cls->GetQ());
173       
174       int tb = cls->GetPadTime();
175       fClsTB->Fill(tb);
176       if (tb > 5 && tb <25)
177         fClsAmpDrift->Fill(cls->GetQ()); 
178       
179       fClsAmpDriftDet[cls->GetDetector()]->Fill(cls->GetQ());
180     }
181     
182     fClusterArray->Delete();
183     
184   }
185    
186   fClsAmpDist->Reset();
187   for(int det=0; det<540; det++) {
188     if (fClsAmpDriftDet[det]->GetSum() > 0) 
189       fClsAmpDist->Fill(fClsAmpDriftDet[det]->GetMean());
190     fSClsDist->Fill(sClusterCharge[det]);
191   }
192
193   PushBack((TObject*)fNClsDet, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);   
194   PushBack((TObject*)fClsAmp, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);  
195   PushBack((TObject*)fClsAmpDrift, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);   
196   PushBack((TObject*)fClsTB, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);  
197   PushBack((TObject*)fClsAmpDist, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);  
198   PushBack((TObject*)fSClsDist, kAliHLTDataTypeHistogram | kAliHLTDataOriginTRD, 0);  
199     
200   return 0;
201 }