]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/MUON/utils/AliHLTMUONClusterHistoComponent.cxx
170de46afb637dd03ca76b02b20f33124c1efe77
[u/mrichter/AliRoot.git] / HLT / MUON / utils / AliHLTMUONClusterHistoComponent.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 Author: Arshad Ahmad Masoodi <Arshad.Ahmad@cern.ch>            *
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   AliHLTMUONClusterHistoComponent.cxx
18     @author Arshad Ahmad <Arshad.Ahmad@cern.ch>
19     @date  27 Nov 2009
20     @brief  Component for onlinehistograms
21 */
22
23 #if __GNUC__>= 3
24 using namespace std;
25 #endif
26
27 #include "AliHLTMUONClusterHistoComponent.h"
28 #include "AliCDBEntry.h"
29 #include "AliCDBManager.h"
30 #include "AliHLTDataTypes.h"
31 #include "AliHLTMUONConstants.h"
32 #include "AliHLTMUONClustersBlockStruct.h"
33 #include "AliHLTMUONDataBlockReader.h"
34 #include <TFile.h>
35 #include <TString.h>
36 #include "TObjString.h"
37 #include "TObjArray.h"
38
39
40 /** ROOT macro for the implementation of ROOT specific class methods */
41 ClassImp(AliHLTMUONClusterHistoComponent);
42
43 AliHLTMUONClusterHistoComponent::AliHLTMUONClusterHistoComponent() :
44   AliHLTMUONProcessor(),
45   fChargePerClusterBending(NULL),
46   fChargePerClusterNonBending(NULL),
47   fNumberOfClusters(NULL),
48   fPlotChargePerClusterBending(kTRUE),
49   fPlotChargePerClusterNonBending(kTRUE),
50   fPlotNClusters(kTRUE)
51 {
52   // see header file for class documentation
53 }
54
55 AliHLTMUONClusterHistoComponent::~AliHLTMUONClusterHistoComponent()
56 {
57   // see header file for class documentation
58   
59  if (fChargePerClusterBending)          delete fChargePerClusterBending;
60  if (fChargePerClusterNonBending)       delete fChargePerClusterNonBending;
61  if (fNumberOfClusters)                 delete fNumberOfClusters;
62 }
63
64 // Public functions to implement AliHLTComponent's interface.
65 // These functions are required for the registration process
66
67 const char* AliHLTMUONClusterHistoComponent::GetComponentID()
68 {
69   // see header file for class documentation
70   
71   return AliHLTMUONConstants::ClusterHistogrammerId();
72 //   return "MUONClusterHistogrammer";
73
74 }
75
76 void AliHLTMUONClusterHistoComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list)
77 {
78   // see header file for class documentation
79   list.clear();
80   list.push_back(AliHLTMUONConstants::ClusterBlockDataType() );
81   //list.push_back( AliHLTMUONConstants::RecHitsBlockDataType() );
82   //list.push_back( AliHLTMUONConstants::TriggerRecordsBlockDataType() );
83 }
84
85 AliHLTComponentDataType AliHLTMUONClusterHistoComponent::GetOutputDataType()
86 {
87   // see header file for class documentation
88   return kAliHLTDataTypeHistogram;
89 }
90
91
92 int AliHLTMUONClusterHistoComponent::GetOutputDataTypes(AliHLTComponentDataTypeList& tgtList)
93
94 {
95   // see header file for class documentation
96   tgtList.clear();
97   tgtList.push_back( AliHLTMUONConstants::HistogramDataType() );
98   return tgtList.size();
99 }
100
101 void AliHLTMUONClusterHistoComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
102 {
103   // see header file for class documentation
104   
105   // The total constant size is the size of the TH1F object, the space needed for
106   // the arrays of data points and finally we add 4k for any extra streamer meta data
107   // etc. added by ROOT.
108   constBase = sizeof(TH1F)*3 + sizeof(Double_t)*(400*2+101) + 4*1024*3;
109   
110   // The multiplier can be zero since the histograms generated are a constant size
111   // independent of the input data size.
112   inputMultiplier = 0;
113 }
114
115 AliHLTComponent* AliHLTMUONClusterHistoComponent::Spawn()
116 {
117   // see header file for class documentation
118   return new AliHLTMUONClusterHistoComponent;
119 }
120
121 int AliHLTMUONClusterHistoComponent::DoInit( int argc, const char** argv )
122 {
123   // see header file for class documentation
124   
125   fPlotChargePerClusterBending=kTRUE;
126   fPlotChargePerClusterNonBending=kTRUE;
127   fPlotNClusters=kTRUE;
128    
129   if(fPlotChargePerClusterBending)
130     fChargePerClusterBending = new TH1F("fChargePerClusterBending","Total Charge of clusters ",400,0,4000);
131   
132   if(fPlotChargePerClusterNonBending)
133     fChargePerClusterNonBending = new TH1F("fChargePerClusterNonBending","Total Charge of clusters ",400,0,4000);
134   
135   if(fPlotNClusters)
136     fNumberOfClusters = new TH1F("fNumberOfClusters","Total Number of Clusters",101,0,100);
137   
138   int iResult=0;
139   TString configuration="";
140   TString argument="";
141   for (int i=0; i<argc && iResult>=0; i++) {
142     argument=argv[i];
143     if (!configuration.IsNull()) configuration+=" ";
144     configuration+=argument;
145   }
146   
147   if (!configuration.IsNull()) {
148     iResult=Configure(configuration.Data());
149   }
150
151   return iResult;
152 }
153
154 int AliHLTMUONClusterHistoComponent::DoDeinit()
155 {
156   // see header file for class documentation
157   if(fChargePerClusterBending!=NULL) delete fChargePerClusterBending;
158   if(fChargePerClusterNonBending!=NULL) delete fChargePerClusterNonBending;
159   if(fNumberOfClusters!=NULL) delete fNumberOfClusters;
160   return 0;
161 }
162
163 int AliHLTMUONClusterHistoComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/,
164                                                 const AliHLTComponentBlockData* /*blocks*/,
165                                                 AliHLTComponentTriggerData& /*trigData*/,
166                                                 AliHLTUInt8_t* /*outputPtr*/,
167                                                 AliHLTUInt32_t& size,
168                                                 AliHLTComponentBlockDataList& /*outputBlocks*/)
169 {
170   AliHLTUInt32_t specification = 0x0;
171   
172   if ( GetFirstInputBlock( kAliHLTDataTypeSOR ) || GetFirstInputBlock( kAliHLTDataTypeEOR ) ) return 0;
173      
174   for (const AliHLTComponentBlockData* block = GetFirstInputBlock(AliHLTMUONConstants::ClusterBlockDataType());
175        block != NULL;block = GetNextInputBlock()){
176     
177     if (block->fDataType != AliHLTMUONConstants::ClusterBlockDataType()) continue;
178     specification |= block->fSpecification;  // The specification bit pattern should indicate all the DDLs that contributed.
179     
180     HLTDebug("Handling block: with fDataType = '%s', fPtr = %p and fSize = %u bytes.",
181              DataType2Text(block->fDataType).c_str(), block->fPtr, block->fSize
182              );
183     
184     AliHLTMUONClustersBlockReader clusterBlock(block->fPtr, block->fSize);
185     if (not BlockStructureOk(clusterBlock))
186     {
187       //FIXME: Need to inherit AliHLTMUONProcessor DoInit functionality properly for
188       // the following Dump feature to work properly. Like in AliHLTMUONRawDataHistoComponent.
189       //if (DumpDataOnError()) DumpEvent(evtData, blocks, trigData, outputPtr, size, outputBlocks);
190       continue;
191     }
192     
193     const AliHLTMUONClusterStruct *cluster=clusterBlock.GetArray();
194     for (AliHLTUInt32_t i = 0; i < clusterBlock.Nentries(); ++i){
195       
196       //commented for future use
197       //AliHLTInt32_t detelement=cluster->fDetElemId;  // Detector ID number from AliRoot geometry
198
199       AliHLTUInt16_t nchannelsB=cluster->fNchannelsB; // Number of channels/pads in the cluster in bending plane.
200       AliHLTUInt16_t nchannelsNB=cluster->fNchannelsNB; // Number of channels/pads in the cluster in non-bending plane.
201       
202
203       AliHLTFloat32_t BCharge= cluster->fChargeB; // Cluster charge in bending plane. Can be -1 if invalid or uncomputed.
204       AliHLTFloat32_t NBCharge= cluster->fChargeNB; // Cluster charge in bending plane. Can be -1 if invalid or uncomputed.
205       
206       if(fPlotChargePerClusterBending)
207         fChargePerClusterBending->Fill(BCharge);
208       if(fPlotChargePerClusterNonBending)       
209         fChargePerClusterNonBending->Fill(NBCharge);
210       cluster++;
211     }  // forloop clusterBlock.Nentries
212     if(fPlotNClusters and clusterBlock.Nentries()>0)
213             fNumberOfClusters->Fill(clusterBlock.Nentries());
214     
215   }//forloop clusterblock
216   
217   if( fPlotChargePerClusterBending ) PushBack( fChargePerClusterBending, AliHLTMUONConstants::HistogramDataType(), specification);
218   if( fPlotChargePerClusterNonBending ) PushBack( fChargePerClusterBending, AliHLTMUONConstants::HistogramDataType(), specification);
219   if( fPlotNClusters ) PushBack( fNumberOfClusters, AliHLTMUONConstants::HistogramDataType(), specification);
220
221   return 0;
222 }
223
224
225 int AliHLTMUONClusterHistoComponent::Reconfigure(const char* cdbEntry, const char* chainId)
226 {
227   // see header file for class documentation
228   int iResult=0;
229   const char* path="HLT/ConfigMUON/MUONHistoComponent";
230   const char* defaultNotify="";
231   if (cdbEntry) {
232     path=cdbEntry;
233     defaultNotify=" (default)";
234   }
235   if (path) {
236     HLTInfo("reconfigure from entry %s%s, chain id %s", path, defaultNotify,(chainId!=NULL && chainId[0]!=0)?chainId:"<none>");
237     AliCDBEntry *pEntry = AliCDBManager::Instance()->Get(path/*,GetRunNo()*/);
238         if (pEntry) {
239         TObjString* pString=dynamic_cast<TObjString*>(pEntry->GetObject());
240         if (pString) {
241         HLTInfo("received configuration object string: \'%s\'", pString->GetString().Data());
242         iResult=Configure(pString->GetString().Data());
243         } else {
244         HLTError("configuration object \"%s\" has wrong type, required TObjString", path);
245       }
246         } else {
247       HLTError("can not fetch object \"%s\" from CDB", path);
248     }
249   }
250
251   return iResult;
252 }
253
254
255 int AliHLTMUONClusterHistoComponent::Configure(const char* arguments)
256 {
257   
258   int iResult=0;
259   if (!arguments) return iResult;
260   
261   TString allArgs=arguments;
262   TString argument;
263   
264   TObjArray* pTokens=allArgs.Tokenize(" ");
265   
266   if (pTokens) {
267     for (int i=0; i<pTokens->GetEntries() && iResult>=0; i++) {
268       argument=((TObjString*)pTokens->At(i))->GetString();
269       if (argument.IsNull()) continue;
270       
271     }
272     delete pTokens;
273   }
274   
275   return iResult;
276 }