]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/ITS/AliHLTITSSSDQARecPointsComponent.cxx
dynamic adjustment of the output buffer size estimator
[u/mrichter/AliRoot.git] / HLT / ITS / AliHLTITSSSDQARecPointsComponent.cxx
1 // $Id$
2 //**************************************************************************
3 //* This file is property of and copyright by the ALICE HLT Project        * 
4 //* ALICE Experiment at CERN, All rights reserved.                         *
5 //*                                                                        *
6 //* Primary Authors: Ingrid Kielen                                         *
7 //*                  for The ALICE HLT Project.                            *
8 //*                                                                        *
9 //* Permission to use, copy, modify and distribute this software and its   *
10 //* documentation strictly for non-commercial purposes is hereby granted   *
11 //* without fee, provided that the above copyright notice appears in all   *
12 //* copies and that both the copyright notice and this permission notice   *
13 //* appear in the supporting documentation. The authors make no claims     *
14 //* about the suitability of this software for any purpose. It is          *
15 //* provided "as is" without express or implied warranty.                  *
16 //**************************************************************************
17
18 /** @file   AliHLTITSSSDQARecPointsComponent.cxx
19     @author Ingrid Kielen - Panos Christakoglou (Panos.Christakoglou@cern.ch)
20     @brief  Component for the SSD clusters QA
21 */
22
23 #if __GNUC__>= 3
24 using namespace std;
25 #endif
26
27 #include "AliHLTITSSSDQARecPointsComponent.h"
28 #include "AliHLTITSClusterDataFormat.h"
29 #include "AliCDBEntry.h"
30 #include "AliCDBManager.h"
31 #include "AliITSRecPoint.h"
32 #include <TFile.h>
33 #include <TMath.h>
34 #include <TString.h>
35 #include "TObjString.h"
36 #include "TObjArray.h"
37 #include "TH1.h"
38 #include "TH2.h"
39 #include "AliITSgeomTGeo.h"
40
41 //#include <stdlib.h>
42 //#include <cerrno>
43
44 /** ROOT macro for the implementation of ROOT specific class methods */
45 ClassImp(AliHLTITSSSDQARecPointsComponent)
46
47 AliHLTITSSSDQARecPointsComponent::AliHLTITSSSDQARecPointsComponent()
48   :AliHLTProcessor(),
49    fHistSSDArray(NULL) {
50   // see header file for class documentation
51   // or
52   // refer to README to build package
53   // or
54   // visit http://web.ift.uib.no/~kjeks/doc/alice-hlt
55
56 }
57
58 AliHLTITSSSDQARecPointsComponent::~AliHLTITSSSDQARecPointsComponent()
59 {
60   // see header file for class documentation
61 }
62
63 // Public functions to implement AliHLTComponent's interface.
64 // These functions are required for the registration process
65
66 const char* AliHLTITSSSDQARecPointsComponent::GetComponentID()
67 {
68   // see header file for class documentation
69   
70   return "ITSSSDQARecPoints";
71 }
72
73 void AliHLTITSSSDQARecPointsComponent::GetInputDataTypes(AliHLTComponentDataTypeList& list) {
74   // see header file for class documentation
75   list.clear();
76   list.push_back( kAliHLTDataTypeClusters|kAliHLTDataOriginITSSSD );
77 }
78
79 AliHLTComponentDataType AliHLTITSSSDQARecPointsComponent::GetOutputDataType() {
80   // see header file for class documentation
81   return kAliHLTDataTypeTObjArray;
82 }
83
84 void AliHLTITSSSDQARecPointsComponent::GetOutputDataSize( unsigned long& constBase, double& inputMultiplier )
85 {
86   // see header file for class documentation
87   // XXX TODO: Find more realistic values.
88   constBase = 80000;
89   inputMultiplier = 10;
90 }
91
92 AliHLTComponent* AliHLTITSSSDQARecPointsComponent::Spawn()
93 {
94   // see header file for class documentation
95   return new AliHLTITSSSDQARecPointsComponent;
96 }
97
98 int AliHLTITSSSDQARecPointsComponent::DoInit(int /*argc*/, const char** /*argv*/) {
99   //Inititalization of histograms for the SSD QA component
100   fHistSSDArray = new TObjArray();  
101   fHistSSDArray->SetName("ssdArray");
102
103   Int_t nModuleOffset = 500;
104   Int_t nITSTotalModules = AliITSgeomTGeo::GetNModules();  
105   TH1F *gHistSSDModuleIdLayer5 = new TH1F("fHistSSDModuleIdLayer5",
106                                           "Module Id - Layer 5;Module Id;Entries",
107                                           fgkSSDMODULESLAYER5,
108                                           nModuleOffset - 0.5,
109                                           nITSTotalModules-fgkSSDMODULESLAYER6+0.5);
110   fHistSSDArray->Add(gHistSSDModuleIdLayer5);//entry 0
111
112   TH1F *gHistSSDModuleIdLayer6 = new TH1F("fHistSSDModuleIdLayer6",
113                                           "Module Id - Layer 6;Module Id;Entries",
114                                           fgkSSDMODULESLAYER6,
115                                           nModuleOffset+fgkSSDMODULESLAYER5 - 0.5,
116                                           nITSTotalModules + 0.5);
117   fHistSSDArray->Add(gHistSSDModuleIdLayer6);//entry 1
118   
119   TH1F *gHistSSDClusterPerEventLayer5 = new TH1F("fHistSSDClusterPerEventLayer5",
120                                                  "N_{clusters} - Layer 5;N_{clusters};Entries;",
121                                                  100,0.1,5000);
122   fHistSSDArray->Add(gHistSSDClusterPerEventLayer5);//entry 2
123   
124   TH1F *gHistSSDClusterPerEventLayer6 = new TH1F("fHistSSDClusterPerEventLayer6",
125                                                  "N_{clusters} - Layer 6;N_{clusters};Entries;",
126                                                  100,0.1,5000);
127   fHistSSDArray->Add(gHistSSDClusterPerEventLayer6);//entry 3
128   
129   TH1F *gHistSSDLocalXLayer5 = new TH1F("fHistSSDLocalXLayer5",
130                                         "Local x coord.- Layer 5;x [cm];Entries;",
131                                         100,-4.,4.);
132   fHistSSDArray->Add(gHistSSDLocalXLayer5);//entry 4
133   
134   TH1F *gHistSSDLocalXLayer6 = new TH1F("fHistSSDLocalXLayer6",
135                                         "Local x coord.- Layer 6;x [cm];Entries;",
136                                         100,-4.,4.);
137   fHistSSDArray->Add(gHistSSDLocalXLayer6);//entry 5
138   
139   TH1F *gHistSSDLocalZLayer5 = new TH1F("fHistSSDLocalZLayer5",
140                                         "Local z coord.- Layer 5;z [cm];Entries;",
141                                         100,-4.,4.);
142   fHistSSDArray->Add(gHistSSDLocalZLayer5);//entry 6
143   
144   TH1F *gHistSSDLocalZLayer6 = new TH1F("fHistSSDLocalZLayer6",
145                                         "Local z coord.- Layer 6;z [cm];Entries;",
146                                         100,-4.,4.);
147   fHistSSDArray->Add(gHistSSDLocalZLayer6);//entry 7
148   
149   TH1F *gHistSSDGlobalXLayer5 = new TH1F("fHistSSDGlobalXLayer5",
150                                          "Global x - Layer 5;x [cm];Entries;",
151                                          100,-40.,40.);
152   fHistSSDArray->Add(gHistSSDGlobalXLayer5);//entry 8
153   
154   TH1F *gHistSSDGlobalXLayer6 = new TH1F("fHistSSDGlobalXLayer6",
155                                          "Global x - Layer 6;x [cm];Entries;",
156                                          100,-45.,45.);
157   fHistSSDArray->Add(gHistSSDGlobalXLayer6);//entry 9
158   
159   TH1F *gHistSSDGlobalYLayer5 = new TH1F("fHistSSDGlobalYLayer5",
160                                          "Global y - Layer 5;y [cm];Entries;",
161                                          100,-40.,40);
162   fHistSSDArray->Add(gHistSSDGlobalYLayer5);//entry 10
163   
164   TH1F *gHistSSDGlobalYLayer6 = new TH1F("fHistSSDGlobalYLayer6",
165                                          "Global y - Layer 6;y [cm];Entries;",
166                                          100,-45.,45.);
167   fHistSSDArray->Add(gHistSSDGlobalYLayer6);//entry 11
168   
169   TH1F *gHistSSDGlobalZLayer5 = new TH1F("fHistSSDGlobalZLayer5",
170                                          "Global z - Layer 5;z [cm];Entries;",
171                                          100,-45.,45);
172   fHistSSDArray->Add(gHistSSDGlobalZLayer5);//entry 12
173   
174   TH1F *gHistSSDGlobalZLayer6 = new TH1F("fHistSSDGlobalZLayer6",
175                                          "Global z - Layer 6;z [cm];Entries;",
176                                          100,-55.,55.);
177   fHistSSDArray->Add(gHistSSDGlobalZLayer6);//entry 13
178   
179   TH1F *gHistSSDPhiLayer5 = new TH1F("fHistSSDPhiLayer5",
180                                      "#phi - Layer 5;#phi [rad];Entries;",
181                                      100,-TMath::Pi(),TMath::Pi());
182   fHistSSDArray->Add(gHistSSDPhiLayer5);//entry 14
183   
184   TH1F *gHistSSDPhiLayer6 = new TH1F("fHistSSDPhiLayer6",
185                                      "#phi - Layer 6;#phi [rad];Entries;",
186                                      100,-TMath::Pi(),TMath::Pi());
187   fHistSSDArray->Add(gHistSSDPhiLayer6);//entry 15
188   
189   TH1F *gHistSSDThetaLayer5 = new TH1F("fHistSSDThetaLayer5",
190                                        "#theta - Layer 5;#theta [rad];Entries;",
191                                        100,-TMath::Pi(),TMath::Pi());
192   fHistSSDArray->Add(gHistSSDThetaLayer5);//entry 16
193   
194   TH1F *gHistSSDThetaLayer6 = new TH1F("fHistSSDThetaLayer6",
195                                        "#theta - Layer 6;#theta [rad];Entries;",
196                                        100,-TMath::Pi(),TMath::Pi());
197   fHistSSDArray->Add(gHistSSDThetaLayer6);//entry 17
198   
199   TH1F *gHistSSDRadiusLayer5 = new TH1F("fHistSSDRadiusLayer5",
200                                         "r - Layer 5;r [cm];Entries;",
201                                         100,35.,50.);
202   fHistSSDArray->Add(gHistSSDRadiusLayer5);//entry 18
203   
204   TH1F *gHistSSDRadiusLayer6 = new TH1F("fHistSSDRadiusLayer6",
205                                         "r - Layer 6;r [cm];Entries;",
206                                         100,35.,50.);
207   fHistSSDArray->Add(gHistSSDRadiusLayer6);//entry 19
208   
209   TH1F *gHistSSDClusterTypeLayer5 = new TH1F("Layer5/fHistSSDClusterTypeLayer5",
210                                              "CL type - Layer 5;Cluster type;Entries;",
211                                              150,0,150);
212   fHistSSDArray->Add(gHistSSDClusterTypeLayer5);//entry 20
213   
214   TH1F *gHistSSDClusterTypeLayer6 = new TH1F("fHistSSDClusterTypeLayer6",
215                                              "CL type - Layer 6;Cluster type;Entries;",
216                                              150,0,150);
217   fHistSSDArray->Add(gHistSSDClusterTypeLayer6);//entry 21
218   
219   TH1F *gHistSSDChargeRatioLayer5 = new TH1F("fHistSSDChargeRatioLayer5",
220                                              "Charge ratio - Layer 5;q_{ratio};Entries;",
221                                              100,-2.0,2.0);
222   fHistSSDArray->Add(gHistSSDChargeRatioLayer5);//entry 22
223   
224   TH1F *gHistSSDChargeRatioLayer6 = new TH1F("fHistSSDChargeRatioLayer6",
225                                              "Charge ratio - Layer 6;q_{ratio};Entries;",
226                                              100,-2.0,2.0);
227   fHistSSDArray->Add(gHistSSDChargeRatioLayer6);//entry 23
228   
229   TH1F *gHistSSDChargekeVLayer5 = new TH1F("fHistSSDChargekeVLayer5",
230                                            "Charge - Layer 5;q [keV];Entries;",
231                                            100,0.,300.);
232   fHistSSDArray->Add(gHistSSDChargekeVLayer5);//entry 24
233   
234   TH1F *gHistSSDChargekeVLayer6 = new TH1F("fHistSSDChargekeVLayer6",
235                                            "Charge - Layer 6;q [keV];Entries;",
236                                            100,0.,300.);
237   fHistSSDArray->Add(gHistSSDChargekeVLayer6);//entry 25
238   
239   TH1F *gHistSSDChargePSideLayer5 = new TH1F("fHistSSDChargePSideLayer5",
240                                              "Charge P- Layer 5;q_{P} [keV];Entries;",
241                                              100,0.,300.);
242   fHistSSDArray->Add(gHistSSDChargePSideLayer5);//entry 26
243   
244   TH1F *gHistSSDChargePSideLayer6 = new TH1F("fHistSSDChargePSideLayer6",
245                                              "Charge P- Layer 6;q_{P} [keV];Entries;",
246                                              100,0.,300.);
247   fHistSSDArray->Add(gHistSSDChargePSideLayer6);//entry 27
248   
249   TH1F *gHistSSDChargeNSideLayer5 = new TH1F("fHistSSDChargeNSideLayer5",
250                                              "Charge N- Layer 5;q_{N} [keV];Entries;",
251                                              100,0.,300.);
252   fHistSSDArray->Add(gHistSSDChargeNSideLayer5);//entry 28
253   
254   TH1F *gHistSSDChargeNSideLayer6 = new TH1F("fHistSSDChargeNSideLayer6",
255                                              "Charge N- Layer 6;q_{N} [keV];Entries;",
256                                              100,0.,300.);
257   fHistSSDArray->Add(gHistSSDChargeNSideLayer6);//entry 29
258   
259   TH1F *gHistSSDChargeRatio2Layer5 = new TH1F("fHistSSDChargeRatio2Layer5",
260                                               "Charge Ratio qN/qP - Layer 5;q_{N}/q_{P};Entries;",
261                                               100,0,2);
262   fHistSSDArray->Add(gHistSSDChargeRatio2Layer5);//entry 30
263   
264   TH1F *gHistSSDChargeRatio2Layer6 = new TH1F("fHistSSDChargeRatio2Layer6",
265                                               "Charge Ratio qN/qP - Layer 6;q_{N}/q_{P};Entries;",
266                                               100,0,2);
267   fHistSSDArray->Add(gHistSSDChargeRatio2Layer6);//entry 31
268   
269   TH2F *gHistSSDChargePNSideLayer5 = new TH2F("fHistSSDChargePNSideLayer5",
270                                               "Charge correlation - Layer 5;q_{P} [keV];q_{N} [keV]",
271                                               100,0.,300.,
272                                               100,0.,300.);
273   fHistSSDArray->Add(gHistSSDChargePNSideLayer5);//entry 32
274   
275   TH2F *gHistSSDChargePNSideLayer6 = new TH2F("fHistSSDChargePNSideLayer6",
276                                               "Charge correlation - Layer 6;q_{P} [keV];q_{N} [keV]",
277                                               100,0.,300.,
278                                               100,0.,300.);
279   fHistSSDArray->Add(gHistSSDChargePNSideLayer6);//entry 33
280   
281   TH2F *gHistSSDChargeMapLayer5 = new TH2F("fHistSSDChargeMapLayer5",
282                                            "Charge map;N_{modules};N_{Ladders}",
283                                            fgkSSDMODULESPERLADDERLAYER5,
284                                            -0.5,fgkSSDMODULESPERLADDERLAYER5+0.5,
285                                            3*fgkSSDLADDERSLAYER5,
286                                            -0.5,fgkSSDLADDERSLAYER5+0.5);
287   fHistSSDArray->Add(gHistSSDChargeMapLayer5);//entry 34
288   
289   TH2F *gHistSSDChargeMapLayer6 = new TH2F("fHistSSDChargeMapLayer6",
290                                            "Charge map;N_{modules};N_{Ladders}",
291                                            fgkSSDMODULESPERLADDERLAYER6,
292                                            -0.5,fgkSSDMODULESPERLADDERLAYER6+0.5,
293                                            3*fgkSSDLADDERSLAYER6,
294                                            -0.5,fgkSSDLADDERSLAYER6+0.5);
295   fHistSSDArray->Add(gHistSSDChargeMapLayer6);//entry 35
296   
297   TH2F *gHistSSDClusterMapLayer5 = new TH2F("fHistSSDClusterMapLayer5",
298                                             "Layer 5;N_{module};N_{ladder}",
299                                             22,1,23,
300                                             34,500,534);
301   gHistSSDClusterMapLayer5->GetXaxis()->SetTitleColor(1);
302   gHistSSDClusterMapLayer5->SetStats(kFALSE);
303   gHistSSDClusterMapLayer5->GetYaxis()->SetTitleOffset(1.8);
304   gHistSSDClusterMapLayer5->GetXaxis()->SetNdivisions(22);
305   gHistSSDClusterMapLayer5->GetYaxis()->SetNdivisions(34);
306   gHistSSDClusterMapLayer5->GetXaxis()->SetLabelSize(0.03);
307   gHistSSDClusterMapLayer5->GetYaxis()->SetLabelSize(0.03);
308   gHistSSDClusterMapLayer5->GetZaxis()->SetTitleOffset(1.4);
309   gHistSSDClusterMapLayer5->GetZaxis()->SetTitle("N_{clusters}");
310   fHistSSDArray->Add(gHistSSDClusterMapLayer5); //entry 36
311   
312   TH2F *gHistSSDClusterMapLayer6 = new TH2F("fHistSSDClusterMapLayer6",
313                                             "Layer 6;N_{module};N_{ladder}",
314                                             25,1,26,
315                                             38,600,638);
316   gHistSSDClusterMapLayer6->GetXaxis()->SetTitleColor(1);
317   gHistSSDClusterMapLayer6->SetStats(kFALSE);
318   gHistSSDClusterMapLayer6->GetYaxis()->SetTitleOffset(1.8);
319   gHistSSDClusterMapLayer6->GetXaxis()->SetNdivisions(25);
320   gHistSSDClusterMapLayer6->GetYaxis()->SetNdivisions(38);
321   gHistSSDClusterMapLayer6->GetXaxis()->SetLabelSize(0.03);
322   gHistSSDClusterMapLayer6->GetYaxis()->SetLabelSize(0.03);
323   gHistSSDClusterMapLayer6->GetZaxis()->SetTitleOffset(1.4);
324   gHistSSDClusterMapLayer6->GetZaxis()->SetTitle("N_{clusters}");
325   fHistSSDArray->Add(gHistSSDClusterMapLayer6);//entry 37
326
327   HLTInfo("Finished initialization of SSD histograms");  
328
329   return 0;  
330 }
331   
332 int AliHLTITSSSDQARecPointsComponent::DoDeinit() {
333   // see header file for class documentation
334   if(fHistSSDArray) delete fHistSSDArray;
335
336   return 0;
337 }
338
339 int AliHLTITSSSDQARecPointsComponent::DoEvent(const AliHLTComponentEventData& /*evtData*/, AliHLTComponentTriggerData& /*trigData*/)
340 {
341   
342   Int_t nClustersLayer5 = 0, nClustersLayer6 = 0;   
343   
344   const AliHLTComponentBlockData* iter = NULL;
345   
346   if(!IsDataEvent())
347     return 0;
348   
349   for ( iter = GetFirstInputBlock(kAliHLTDataTypeClusters); iter != NULL; iter = GetNextInputBlock() ) {
350     
351     if(iter->fDataType!=(kAliHLTAnyDataType|kAliHLTDataOriginITSSSD))
352       continue;
353     
354     const AliHLTITSClusterData* clusterData = (const AliHLTITSClusterData*) iter->fPtr;
355     Int_t nSpacepoint = (Int_t) clusterData->fSpacePointCnt;
356     AliHLTITSSpacePointData *clusters = (AliHLTITSSpacePointData*) clusterData->fSpacePoints;
357     
358     for(int i = 0; i < nSpacepoint; i++) { //cluster loop
359       Int_t lab[4]={0,0,0,0};
360       Float_t hit[6]={0,0,0,0,0,0};
361       Int_t info[3]={0,0,0};
362       
363       lab[0]=clusters[i].fTracks[0];
364       lab[1]=clusters[i].fTracks[1];
365       lab[2]=clusters[i].fTracks[2];
366       lab[3]=clusters[i].fIndex;
367       hit[0]=clusters[i].fY;
368       hit[1]=clusters[i].fZ;
369       hit[2]=clusters[i].fSigmaY2;
370       hit[3]=clusters[i].fSigmaZ2;
371       hit[4]=clusters[i].fQ;
372       hit[5]=clusters[i].fSigmaYZ;
373       info[0]=clusters[i].fNy;
374       info[1]=clusters[i].fNz;
375       info[2]=clusters[i].fLayer;
376       
377       AliITSRecPoint recp(lab,hit,info);
378       
379       Int_t module = 0;
380       Int_t gLayer = 0, gLadder = 0, gModule = 0;
381       Int_t lLadderLocationY = 0;
382       
383       Float_t cluglo[3]={0.,0.,0.}; 
384       
385       Int_t layer = recp.GetLayer();
386       if (layer == 4) module = recp.GetDetectorIndex() + 500;
387       if (layer == 5) module = recp.GetDetectorIndex() + 1248;
388       
389       AliITSgeomTGeo::GetModuleId(module,gLayer,gLadder,gModule);
390       lLadderLocationY = 3*gLadder; 
391       
392       recp.GetGlobalXYZ(cluglo);
393       Float_t radius = TMath::Sqrt(cluglo[0]*cluglo[0]+cluglo[1]*cluglo[1]); 
394       Float_t phi = TMath::ATan2(cluglo[1],cluglo[0]);
395       Float_t theta = TMath::ATan2(radius,cluglo[2]);
396       Double_t chargeRatio = recp.GetChargeRatio();
397       Double_t clusterCharge = recp.GetQ();
398       Double_t chargePSide = clusterCharge*(1. + chargeRatio);
399       Double_t chargeNSide = clusterCharge*(1. - chargeRatio);
400       
401       if(layer == 4) {
402         dynamic_cast<TH1F *>(fHistSSDArray->At(0))->Fill(module);
403         dynamic_cast<TH1F *>(fHistSSDArray->At(4))->Fill(recp.GetDetLocalX());
404         dynamic_cast<TH1F *>(fHistSSDArray->At(6))->Fill(recp.GetDetLocalZ());
405         dynamic_cast<TH1F *>(fHistSSDArray->At(8))->Fill(cluglo[0]);
406         dynamic_cast<TH1F *>(fHistSSDArray->At(10))->Fill(cluglo[1]);
407         dynamic_cast<TH1F *>(fHistSSDArray->At(12))->Fill(cluglo[2]);
408         dynamic_cast<TH1F *>(fHistSSDArray->At(14))->Fill(phi);
409         dynamic_cast<TH1F *>(fHistSSDArray->At(16))->Fill(theta);
410         dynamic_cast<TH1F *>(fHistSSDArray->At(18))->Fill(radius);
411         dynamic_cast<TH1F *>(fHistSSDArray->At(20))->Fill(recp.GetType());
412         dynamic_cast<TH1F *>(fHistSSDArray->At(22))->Fill(recp.GetChargeRatio());
413         dynamic_cast<TH1F *>(fHistSSDArray->At(24))->Fill(recp.GetQ());
414         dynamic_cast<TH1F *>(fHistSSDArray->At(26))->Fill(chargePSide);
415         dynamic_cast<TH1F *>(fHistSSDArray->At(28))->Fill(chargeNSide);
416         if(chargePSide != 0.)
417           dynamic_cast<TH1F *>(fHistSSDArray->At(30))->Fill(chargeNSide/chargePSide);
418         dynamic_cast<TH2F *>(fHistSSDArray->At(32))->Fill(chargePSide,
419                                                           chargeNSide);
420         dynamic_cast<TH2F *>(fHistSSDArray->At(34))->SetBinContent(gModule,lLadderLocationY,recp.GetQ());
421         dynamic_cast<TH2F *>(fHistSSDArray->At(36))->Fill(gModule,499+gLadder,1);
422         nClustersLayer5 += 1;
423       }//layer 5 histograms
424       if(layer == 5) {
425         dynamic_cast<TH1F *>(fHistSSDArray->At(1))->Fill(module);
426         dynamic_cast<TH1F *>(fHistSSDArray->At(5))->Fill(recp.GetDetLocalX());
427         dynamic_cast<TH1F *>(fHistSSDArray->At(6))->Fill(recp.GetDetLocalZ());
428         dynamic_cast<TH1F *>(fHistSSDArray->At(9))->Fill(cluglo[0]);
429         dynamic_cast<TH1F *>(fHistSSDArray->At(11))->Fill(cluglo[1]);
430         dynamic_cast<TH1F *>(fHistSSDArray->At(13))->Fill(cluglo[2]);
431         dynamic_cast<TH1F *>(fHistSSDArray->At(15))->Fill(phi);
432         dynamic_cast<TH1F *>(fHistSSDArray->At(17))->Fill(theta);
433         dynamic_cast<TH1F *>(fHistSSDArray->At(19))->Fill(radius);
434         dynamic_cast<TH1F *>(fHistSSDArray->At(21))->Fill(recp.GetType());
435         dynamic_cast<TH1F *>(fHistSSDArray->At(23))->Fill(recp.GetChargeRatio());
436         dynamic_cast<TH1F *>(fHistSSDArray->At(25))->Fill(recp.GetQ());
437         dynamic_cast<TH1F *>(fHistSSDArray->At(27))->Fill(chargePSide);
438         dynamic_cast<TH1F *>(fHistSSDArray->At(29))->Fill(chargeNSide);
439         if(chargePSide != 0.)
440           dynamic_cast<TH1F *>(fHistSSDArray->At(31))->Fill(chargeNSide/chargePSide);
441         dynamic_cast<TH2F *>(fHistSSDArray->At(33))->Fill(chargePSide,
442                                                           chargeNSide);
443         dynamic_cast<TH2F *>(fHistSSDArray->At(35))->SetBinContent(gModule,lLadderLocationY,recp.GetQ());
444         dynamic_cast<TH2F *>(fHistSSDArray->At(37))->Fill(gModule,599+gLadder,1);
445         nClustersLayer6 += 1;
446       }//layer 6 histograms
447     }//cluster loop
448     dynamic_cast<TH1F *>(fHistSSDArray->At(2))->Fill(nClustersLayer5);
449     dynamic_cast<TH1F *>(fHistSSDArray->At(3))->Fill(nClustersLayer6);
450  }//input block
451  
452  //Publish array
453  AliHLTUInt32_t fSpecification = 0x0;
454  PushBack( (TObjArray*) fHistSSDArray,kAliHLTDataTypeTObjArray|kAliHLTDataOriginITSSSD,fSpecification);
455  
456  HLTInfo("ITSSSDQARecPoints found %d SSD clusters", nClustersLayer5+nClustersLayer6);
457
458  return 0;
459 }