]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQASDDDataMakerRec.cxx
6be5fa10b0615c66f6a07e86bb5d8c7a90812bd9
[u/mrichter/AliRoot.git] / ITS / AliITSQASDDDataMakerRec.cxx
1 /**************************************************************************
2  * Copyright(c) 2007-2009, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //  *************************************************************
19 //  Checks the quality assurance 
20 //  by comparing with reference data
21 //  contained in a DB
22 //  -------------------------------------------------------------
23 //  W. Ferrarese + P. Cerello Feb 2008
24 //  M.Siciliano Aug 2008 QA RecPoints and HLT mode
25 //  INFN Torino
26
27 // --- ROOT system ---
28
29 #include <TProfile2D.h>
30 #include <TH2D.h>
31 #include <TBranch.h>
32 #include <TTree.h>
33 #include <TGaxis.h>
34 #include <TMath.h>
35 #include <TDirectory.h>
36 #include <TSystem.h>
37 // --- Standard library ---
38
39 // --- AliRoot header files ---
40 #include "AliITSQASDDDataMakerRec.h"
41 #include "AliLog.h"
42 #include "AliQAv1.h"
43 #include "AliQAChecker.h"
44 #include "AliRawReader.h"
45 #include "AliITSRawStream.h"
46 #include "AliITSRawStreamSDD.h"
47 #include "AliITSRawStreamSDDCompressed.h"
48 #include "AliITSDetTypeRec.h"
49 #include "AliITSRecPoint.h"
50 #include "AliITSgeomTGeo.h"
51 #include "AliITSHLTforSDD.h"
52 #include "AliCDBManager.h"
53 #include "AliCDBStorage.h"
54 #include "AliCDBEntry.h"
55 #include "Riostream.h"
56
57 ClassImp(AliITSQASDDDataMakerRec)
58
59 //____________________________________________________________________________ 
60 AliITSQASDDDataMakerRec::AliITSQASDDDataMakerRec(AliITSQADataMakerRec *aliITSQADataMakerRec, Bool_t kMode, Short_t ldc) :
61 TObject(),
62 fAliITSQADataMakerRec(aliITSQADataMakerRec),
63 fkOnline(kMode),
64 fLDC(ldc),
65 fSDDhRawsTask(0),
66 fSDDhRecPointsTask(0),
67 fGenRawsOffset(0),
68 fGenRecPointsOffset(0),
69 fTimeBinSize(1),
70 fDDLModuleMap(0),
71 fHLTMode(0),
72 fHLTSDD(0)
73 {
74   //ctor used to discriminate OnLine-Offline analysis
75   if(fLDC < 0 || fLDC > 4) {
76         AliError("Error: LDC number out of range; return\n");
77   }
78   if(!fkOnline){AliInfo("Offline mode: HLT set from AliITSDetTypeRec for SDD\n");}
79   else
80     if(fkOnline){
81       AliInfo("Online mode: HLT set from environment for SDD\n");
82       SetHLTModeFromEnvironment();
83     }
84   //fDDLModuleMap=NULL;
85 }
86
87 //____________________________________________________________________________ 
88 AliITSQASDDDataMakerRec::AliITSQASDDDataMakerRec(const AliITSQASDDDataMakerRec& qadm) :
89 TObject(),
90 fAliITSQADataMakerRec(qadm.fAliITSQADataMakerRec),
91 fkOnline(qadm.fkOnline),
92 fLDC(qadm.fLDC),
93 fSDDhRawsTask(qadm.fSDDhRawsTask),
94 fSDDhRecPointsTask(qadm.fSDDhRecPointsTask),
95 fGenRawsOffset(qadm.fGenRawsOffset),
96 fGenRecPointsOffset(qadm.fGenRecPointsOffset),
97 fTimeBinSize(1),
98 fDDLModuleMap(0),
99 fHLTMode(qadm.fHLTMode),
100 fHLTSDD( qadm.fHLTSDD)
101 {
102   //copy ctor 
103   fAliITSQADataMakerRec->SetName((const char*)qadm.fAliITSQADataMakerRec->GetName()) ; 
104   fAliITSQADataMakerRec->SetTitle((const char*)qadm.fAliITSQADataMakerRec->GetTitle());
105   fDDLModuleMap=NULL;
106 }
107
108 //____________________________________________________________________________ 
109 AliITSQASDDDataMakerRec::~AliITSQASDDDataMakerRec(){
110   // destructor
111   //    if(fDDLModuleMap) delete fDDLModuleMap;
112 }
113 //__________________________________________________________________
114 AliITSQASDDDataMakerRec& AliITSQASDDDataMakerRec::operator = (const AliITSQASDDDataMakerRec& qac )
115 {
116   // Equal operator.
117   this->~AliITSQASDDDataMakerRec();
118   new(this) AliITSQASDDDataMakerRec(qac);
119   return *this;
120 }
121
122 //____________________________________________________________________________ 
123 void AliITSQASDDDataMakerRec::StartOfDetectorCycle()
124 {
125   //Detector specific actions at start of cycle
126   AliDebug(1,"AliITSQADM::Start of SDD Cycle\n");
127 }
128
129 //____________________________________________________________________________ 
130 void AliITSQASDDDataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t /*task*/, TObjArray* /*list*/)
131 {
132   // launch the QA checking
133   AliDebug(1,"AliITSDM instantiates checker with Run(AliQAv1::kITS, task, list)\n"); 
134 }
135
136 //____________________________________________________________________________ 
137 void AliITSQASDDDataMakerRec::InitRaws()
138
139   // Initialization for RAW data - SDD -
140   const Bool_t expert   = kTRUE ; 
141   const Bool_t saveCorr = kTRUE ; 
142   const Bool_t image    = kTRUE ; 
143   
144   fGenRawsOffset = (fAliITSQADataMakerRec->fRawsQAList[AliRecoParam::kDefault])->GetEntries();
145   AliCDBEntry *ddlMapSDD = AliCDBManager::Instance()->Get("ITS/Calib/DDLMapSDD");
146   Bool_t cacheStatus = AliCDBManager::Instance()->GetCacheFlag();
147   if(!ddlMapSDD)
148     {
149       AliError("Calibration object retrieval failed! SDD will not be processed");
150       fDDLModuleMap = NULL;
151       return;
152     }
153   fDDLModuleMap = (AliITSDDLModuleMapSDD*)ddlMapSDD->GetObject();
154   if(!cacheStatus)ddlMapSDD->SetObject(NULL);
155   ddlMapSDD->SetOwner(kTRUE);
156   if(!cacheStatus)
157     {
158       delete ddlMapSDD;
159     }
160   
161   if(fkOnline==kFALSE){
162     AliInfo("Offline mode: HLTforSDDobject used \n");
163     AliCDBEntry *hltforSDD = AliCDBManager::Instance()->Get("ITS/Calib/HLTforSDD");
164     if(!hltforSDD){
165       AliError("Calibration object retrieval failed! SDD will not be processed");    
166       fHLTSDD=NULL;
167       return;
168     }  
169     fHLTSDD = (AliITSHLTforSDD*)hltforSDD->GetObject();
170     if(!cacheStatus)hltforSDD->SetObject(NULL);
171     hltforSDD->SetOwner(kTRUE);
172     if(!cacheStatus)
173       {
174         delete hltforSDD;
175       }
176   }
177   Int_t lay, lad, det;
178   Int_t indexlast = 0;
179   Int_t index1 = 0;
180
181   if(fkOnline) 
182     {
183       AliInfo("Book Online Histograms for SDD\n");
184     }
185   else 
186     {
187       AliInfo("Book Offline Histograms for SDD\n ");
188     }
189   TH1D *h0 = new TH1D("SDDModPattern","HW Modules pattern",fgknSDDmodules,239.5,499.5); //0
190   h0->GetXaxis()->SetTitle("Module Number");
191   h0->GetYaxis()->SetTitle("Counts");
192   fAliITSQADataMakerRec->Add2RawsList((new TH1D(*h0)),0+fGenRawsOffset, expert, !image, !saveCorr);
193   delete h0;
194   fSDDhRawsTask++;
195   
196   //zPhi distribution using ladder and modules numbers
197   TH2D *hphil3 = new TH2D("SDDphizL3","SDD #varphiz Layer3 ",6,0.5,6.5,14,0.5,14.5);
198   hphil3->GetXaxis()->SetTitle("z[#Module L3 ]");
199   hphil3->GetYaxis()->SetTitle("#varphi[#Ladder L3]");
200   fAliITSQADataMakerRec->Add2RawsList((new TH2D(*hphil3)),1+fGenRawsOffset, !expert, image, saveCorr); 
201   delete hphil3;
202   fSDDhRawsTask++;
203   
204   TH2D *hphil4 = new TH2D("SDDphizL4","SDD #varphiz Layer4 ",8,0.5,8.5,22,0.5,22.5); 
205   hphil4->GetXaxis()->SetTitle("z[#Module L4]");
206   hphil4->GetYaxis()->SetTitle("#varphi[#Ladder L4]");
207   fAliITSQADataMakerRec->Add2RawsList((new TH2D(*hphil4)),2+fGenRawsOffset, !expert, image, saveCorr); 
208   delete hphil4;
209   fSDDhRawsTask++;
210   
211
212   if(fkOnline) 
213     {
214
215       //DDL Pattern 
216       TH2D *hddl = new TH2D("SDDDDLPattern","SDD DDL Pattern ",24,-0.5,23.5,24,-0.5,23.5); 
217       hddl->GetXaxis()->SetTitle("Channel");
218       hddl->GetYaxis()->SetTitle("#DDL");
219       fAliITSQADataMakerRec->Add2RawsList((new TH2D(*hddl)),3+fGenRawsOffset, expert, !image, !saveCorr);
220       delete hddl;
221       fSDDhRawsTask++;
222       Int_t indexlast1 = 0;
223   
224       fTimeBinSize = 4;
225       indexlast = 0;
226       index1 = 0;
227       indexlast1 = fSDDhRawsTask;
228       char *hname[3];
229       for(Int_t i=0; i<3; i++) hname[i]= new char[50];
230       for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){
231         for(Int_t iside=0;iside<fgknSide;iside++){
232           AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);
233           sprintf(hname[0],"SDDchargeMapFSE_L%d_%d_%d_%d",lay,lad,det,iside);
234           sprintf(hname[1],"SDDChargeMapForSingleEvent_L%d_%d_%d_%d",lay,lad,det,iside);
235           sprintf(hname[2],"SDDhmonoDMap_L%d_%d_%d_%d",lay,lad,det,iside);
236           TProfile2D *fModuleChargeMapFSE = new TProfile2D(hname[0],hname[1],256/fTimeBinSize,-0.5,255.5,256,-0.5,255.5);
237           fModuleChargeMapFSE->GetXaxis()->SetTitle("Time Bin");
238           fModuleChargeMapFSE->GetYaxis()->SetTitle("Anode");
239           fAliITSQADataMakerRec->Add2RawsList((new TProfile2D(*fModuleChargeMapFSE)),indexlast1 + index1 + fGenRawsOffset, expert, !image, !saveCorr);
240           delete fModuleChargeMapFSE;
241           
242           fSDDhRawsTask++;
243           index1++;      
244         }
245       }
246       
247       for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){
248         for(Int_t iside=0;iside<fgknSide;iside++){
249           AliITSgeomTGeo::GetModuleId(moduleSDD+fgkmodoffset, lay, lad, det);
250           sprintf(hname[0],"SDDchargeMap_L%d_%d_%d_%d",lay,lad,det,iside);
251           sprintf(hname[1],"SDDChargeMap_L%d_%d_%d_%d",lay,lad,det,iside);
252           TProfile2D *fModuleChargeMap = new TProfile2D(hname[0],hname[1],256/fTimeBinSize,-0.5,255.5,256,-0.5,255.5);
253           fModuleChargeMap->GetXaxis()->SetTitle("Time Bin");
254           fModuleChargeMap->GetYaxis()->SetTitle("Anode");
255           fAliITSQADataMakerRec->Add2RawsList((new TProfile2D(*fModuleChargeMap)),indexlast1 + index1 + fGenRawsOffset, expert, !image, !saveCorr);
256           delete fModuleChargeMap;
257           
258           fSDDhRawsTask++;
259           index1++;      
260         }
261       }
262       
263     }  // kONLINE
264   
265   AliDebug(1,Form("%d SDD Raws histograms booked\n",fSDDhRawsTask));
266 }
267
268
269 //____________________________________________________________________________
270 void AliITSQASDDDataMakerRec::MakeRaws(AliRawReader* rawReader)
271
272   // Fill QA for RAW - SDD -
273   
274   if(!fDDLModuleMap){
275     AliError("SDD DDL module map not available - skipping SDD QA");
276     return;
277   }
278   if(rawReader->GetType() != 7) return;  // skips non physical triggers
279   AliDebug(1,"entering MakeRaws\n");                 
280   
281   rawReader->Reset();       
282   AliITSRawStream *stream;
283   
284   if(fkOnline==kTRUE)
285     {
286       if(GetHLTMode()==kTRUE)
287         {
288           //AliInfo("Online  mode: HLT C compressed mode used for SDD\n");
289           stream = new AliITSRawStreamSDDCompressed(rawReader); }
290       else{ 
291         //AliInfo("Online  mode: HLT A mode used for SDD\n");
292         stream = new AliITSRawStreamSDD(rawReader);}     
293     }
294   else 
295     {
296       if(fHLTSDD->IsHLTmodeC()==kTRUE){
297         //AliInfo("Offline  mode: HLT C compressed mode used for SDD\n");
298         stream = new AliITSRawStreamSDDCompressed(rawReader);
299       }else 
300         {
301           //AliInfo("Offline  mode: HLT A mode used for SDD\n");
302           stream = new AliITSRawStreamSDD(rawReader);
303         }
304     }
305   
306   //ckeck on HLT mode
307  
308   //  AliITSRawStreamSDD s(rawReader); 
309   stream->SetDDLModuleMap(fDDLModuleMap);
310   
311   Int_t lay, lad, det; 
312   
313   Int_t index=0;
314   if(fkOnline) {
315     for(Int_t moduleSDD =0; moduleSDD<fgknSDDmodules; moduleSDD++){
316       for(Int_t iside=0;iside<fgknSide;iside++) {
317         if(fSDDhRawsTask > 4 + index) fAliITSQADataMakerRec->GetRawsData(4 + index +fGenRawsOffset)->Reset();   
318         // 4  because the 2D histos for single events start after the fourth position
319         index++;
320       }
321     }
322   }
323   
324   Int_t cnt = 0;
325   Int_t ildcID = -1;
326   Int_t iddl = -1;
327   Int_t isddmod = -1;
328   Int_t coord1, coord2, signal, moduleSDD, activeModule, index1; 
329   
330   while(stream->Next()) {
331     ildcID = rawReader->GetLDCId();
332     iddl = rawReader->GetDDLID() - fgkDDLIDshift;
333     
334     isddmod = fDDLModuleMap->GetModuleNumber(iddl,stream->GetCarlosId());
335     if(isddmod==-1){
336       AliDebug(1,Form("Found module with iddl: %d, stream->GetCarlosId: %d \n",iddl,stream->GetCarlosId()));
337       continue;
338     }
339     if(stream->IsCompletedModule()) {
340       AliDebug(1,Form("IsCompletedModule == KTRUE\n"));
341       continue;
342     } 
343     if(stream->IsCompletedDDL()) {
344       AliDebug(1,Form("IsCompletedDDL == KTRUE\n"));
345       continue;
346     } 
347     
348     coord1 = stream->GetCoord1();
349     coord2 = stream->GetCoord2();
350     signal = stream->GetSignal();
351     
352     moduleSDD = isddmod - fgkmodoffset;
353     
354     if(isddmod <fgkmodoffset|| isddmod>fgknSDDmodules+fgkmodoffset-1) {
355       AliDebug(1,Form( "Module SDD = %d, resetting it to 1 \n",isddmod));
356       isddmod = 1;
357     }
358     
359     AliITSgeomTGeo::GetModuleId(isddmod, lay, lad, det);
360
361     
362     fAliITSQADataMakerRec->GetRawsData( 0 + fGenRawsOffset )->Fill(isddmod);   
363     
364     if(lay==3)    fAliITSQADataMakerRec->GetRawsData(1+fGenRawsOffset)->Fill(det,lad); 
365     if(lay==4) { 
366       fAliITSQADataMakerRec->GetRawsData(2+fGenRawsOffset)->Fill(det,lad);}  
367     
368     Short_t iside = stream->GetChannel();
369     
370     fAliITSQADataMakerRec->GetRawsData(3+fGenRawsOffset)->Fill(2*(stream->GetCarlosId())+iside,iddl);
371     
372
373     if(fkOnline) {
374
375       activeModule = moduleSDD;
376       index1 = activeModule * 2 + iside;
377       
378       if(index1<0){
379         AliDebug(1,Form("Wrong index number %d - patched to 0\n",index1));
380         index1 = 0;
381       }      
382       fAliITSQADataMakerRec->GetRawsData(3+fGenRawsOffset)->Fill(2*(stream->GetCarlosId())+iside,iddl);
383       if(fSDDhRawsTask > 4 + index1) {                                  
384         ((TProfile2D *)(fAliITSQADataMakerRec->GetRawsData(4 + index1 +fGenRawsOffset)))->Fill(coord2, coord1, signal);     
385         ((TProfile2D *)(fAliITSQADataMakerRec->GetRawsData(4 + index1 + 260*2 +fGenRawsOffset)))->Fill(coord2, coord1, signal); 
386       }
387     }
388     cnt++;
389     if(!(cnt%10000)) AliDebug(1,Form(" %d raw digits read",cnt));
390   }
391   AliDebug(1,Form("Event completed, %d raw digits read",cnt)); 
392   delete stream;
393   stream = NULL; 
394 }
395
396 //____________________________________________________________________________ 
397 void AliITSQASDDDataMakerRec::InitRecPoints()
398 {
399   // Initialization for RECPOINTS - SDD -
400   const Bool_t expert   = kTRUE ; 
401   const Bool_t image    = kTRUE ; 
402   
403   fGenRecPointsOffset = (fAliITSQADataMakerRec->fRecPointsQAList[AliRecoParam::kDefault])->GetEntries();
404
405   Int_t nOnline=1;
406   Int_t  nOnline2=1;
407   Int_t  nOnline3=1; 
408   Int_t  nOnline4=1;
409   if(fkOnline)
410     {
411       nOnline=4;
412       nOnline2=28;
413       nOnline3=64;
414       nOnline4=14;
415     }
416
417   
418   TH1F *h0 = new TH1F("SDDLay3TotCh","Layer 3 total charge",1000/nOnline,-0.5, 499.5); //position number 0
419   h0->GetXaxis()->SetTitle("ADC value");
420   h0->GetYaxis()->SetTitle("Entries");
421   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h0)), 0 +fGenRecPointsOffset, !expert, image);
422   delete h0;
423   fSDDhRecPointsTask++;
424  
425   TH1F *h1 = new TH1F("SDDLay4TotCh","Layer 4 total charge",1000/nOnline,-0.5, 499.5);//position number 1
426   h1->GetXaxis()->SetTitle("ADC value");
427   h1->GetYaxis()->SetTitle("Entries");
428   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h1)), 1 +fGenRecPointsOffset, !expert, image);
429   delete h1;
430   fSDDhRecPointsTask++;
431
432   char hisnam[50];
433   TH2F *h2 = new TH2F("SDDGlobalCoordDistribYX","YX Global Coord Distrib",5600/nOnline2,-28,28,5600/nOnline2,-28,28);//position number 2
434   h2->GetYaxis()->SetTitle("Y[cm]");
435   h2->GetXaxis()->SetTitle("X[cm]");
436   fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h2)),2+fGenRecPointsOffset, expert, !image);
437   delete h2;
438   fSDDhRecPointsTask++;
439
440   TH2F *h3 = new TH2F("SDDGlobalCoordDistribRZ","RZ Global Coord Distrib",6400/nOnline3,-32,32,1400/nOnline4,12,26);//position number 3
441   h3->GetYaxis()->SetTitle("R[cm]");
442   h3->GetXaxis()->SetTitle("Z[cm]");
443   fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h3)),3+fGenRecPointsOffset, expert, !image);
444   delete h3;
445   fSDDhRecPointsTask++;
446   
447   TH2F *h4 = new TH2F("SDDGlobalCoordDistribL3PHIZ","#varphi Z Global Coord Distrib L3",6400/nOnline3,-32,32,360/nOnline,-TMath::Pi(),TMath::Pi());//position number 4
448   h4->GetYaxis()->SetTitle("#phi[rad]");
449   h4->GetXaxis()->SetTitle("Z[cm]");
450   fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h4)),4+fGenRecPointsOffset, !expert, image);
451   delete h4;
452   fSDDhRecPointsTask++;
453
454   TH2F *h5 = new TH2F("SDDGlobalCoordDistribL4PHIZ","#varphi Z Global Coord Distrib L4",6400/nOnline3,-32,32,360/nOnline,-TMath::Pi(),TMath::Pi());//position number 5
455   h5->GetYaxis()->SetTitle("#phi[rad]");
456   h5->GetXaxis()->SetTitle("Z[cm]");
457   fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h5)),5+fGenRecPointsOffset, !expert, image);
458   delete h5;
459   fSDDhRecPointsTask++;
460   
461   TH1F *h6 = new TH1F("SDDModPatternRP","Modules pattern RP",fgknSDDmodules,239.5,499.5); //position number 6
462   h6->GetXaxis()->SetTitle("Module number"); //spd offset = 240
463   h6->GetYaxis()->SetTitle("Entries");
464   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h6)),6 +fGenRecPointsOffset, expert, !image);
465   delete h6;
466   fSDDhRecPointsTask++;
467   TH1F *h7 = new TH1F("SDDLadPatternL3RP","Ladder pattern L3 RP",14,0.5,14.5);  //position number 7
468   h7->GetXaxis()->SetTitle("Ladder #, Layer 3");
469   h7->GetYaxis()->SetTitle("Entries");
470   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h7)),7 +fGenRecPointsOffset, expert, !image);
471   delete h7;
472   fSDDhRecPointsTask++;
473   TH1F *h8 = new TH1F("SDDLadPatternL4RP","Ladder pattern L4 RP",22,0.5,22.5); //position number 8
474   h8->GetXaxis()->SetTitle("Ladder #, Layer 4");
475   h8->GetYaxis()->SetTitle("Entries");
476   fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h8)),8 +fGenRecPointsOffset, expert, !image);
477   delete h8;
478   fSDDhRecPointsTask++;
479   TH2F *h9 = new TH2F("SDDLocalCoordDistrib","Local Coord Distrib",1000/nOnline,-4,4,1000/nOnline,-4,4);//position number 9
480   h9->GetXaxis()->SetTitle("X local coord, drift, cm");
481   h9->GetYaxis()->SetTitle("Z local coord, anode, cm");
482   fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h9)),9 +fGenRecPointsOffset, expert, !image);
483   delete h9;
484   fSDDhRecPointsTask++;
485
486
487     TH1F *h10 = new TH1F("SDDrdistrib_Layer3" ,"SDD r distribution Layer3" ,100,14.,18.);//position number 10 (L3)
488     h10->GetXaxis()->SetTitle("r[cm]");
489     h10->GetXaxis()->CenterTitle();
490     h10->GetYaxis()->SetTitle("Entries");
491     fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h10)),10 +fGenRecPointsOffset, expert, !image);
492     delete h10;
493     fSDDhRecPointsTask++;
494
495     TH1F *h11 = new TH1F("SDDrdistrib_Layer4" ,"SDD r distribution Layer4" ,100,22.,26.);// and position number 11 (L4)
496     h11->GetXaxis()->SetTitle("r[cm]");
497     h11->GetXaxis()->CenterTitle();
498     h11->GetYaxis()->SetTitle("Entries");
499     fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h11)),11 +fGenRecPointsOffset, expert, !image);
500     delete h11;
501     fSDDhRecPointsTask++;
502
503   for(Int_t iLay=0; iLay<=1; iLay++){
504     sprintf(hisnam,"SDDphidistrib_Layer%d",iLay+3);
505     TH1F *h12 = new TH1F(hisnam,hisnam,180,-TMath::Pi(),TMath::Pi());//position number 12 (L3) and position number 13 (L4)
506     h12->GetXaxis()->SetTitle("#varphi[rad]");
507     h12->GetXaxis()->CenterTitle();
508     h12->GetYaxis()->SetTitle("Entries");
509     fAliITSQADataMakerRec->Add2RecPointsList((new TH1F(*h12)),iLay+12+fGenRecPointsOffset, expert, !image);
510     delete h12;
511     fSDDhRecPointsTask++;
512   }
513
514   if(fkOnline)
515     {
516       TH2F *h14 = new TH2F("SDDGlobalCoordDistribYXFSE","YX Global Coord Distrib FSE",5600/nOnline2,-28,28,5600/nOnline2,-28,28);//position number 14
517       h14->GetYaxis()->SetTitle("Y[cm]");
518       h14->GetXaxis()->SetTitle("X[cm]");
519       fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h14)),14+fGenRecPointsOffset, expert, !image);
520       delete h14;
521       fSDDhRecPointsTask++;
522       
523       TH2F *h15 = new TH2F("SDDGlobalCoordDistribRZFSE","RZ Global Coord Distrib FSE",Int_t(6400/nOnline3),-32,32,1400/nOnline4,12,26);//position number 15
524       h15->GetYaxis()->SetTitle("R[cm]");
525       h15->GetXaxis()->SetTitle("Z[cm]");
526       fAliITSQADataMakerRec->Add2RecPointsList((new TH2F(*h15)),15+fGenRecPointsOffset, expert, !image);
527       delete h15;
528       fSDDhRecPointsTask++;
529       
530     }//online
531
532   //printf("%d SDD Recs histograms booked\n",fSDDhRecPointsTask);
533
534
535   AliDebug(1,Form("%d SDD Recs histograms booked\n",fSDDhRecPointsTask));
536
537
538 }
539
540 //____________________________________________________________________________ 
541 void AliITSQASDDDataMakerRec::MakeRecPoints(TTree * clustersTree)
542 {
543
544
545  // Fill QA for RecPoints - SDD -
546   Int_t lay, lad, det; 
547   TBranch *branchRecP = clustersTree->GetBranch("ITSRecPoints");
548   if (!branchRecP) {
549     AliError("can't get the branch with the ITS clusters !");
550     return;
551   }
552
553   static TClonesArray statRecpoints("AliITSRecPoint") ;
554   TClonesArray *recpoints = &statRecpoints;
555   branchRecP->SetAddress(&recpoints);
556   Int_t npoints = 0;      
557   Float_t cluglo[3]={0.,0.,0.}; 
558   if(fkOnline)
559     {
560       for(Int_t i=14;i<16;i++)
561         {
562           fAliITSQADataMakerRec->GetRecPointsData(i+fGenRecPointsOffset)->Reset();
563         }
564     }
565   for(Int_t module=0; module<clustersTree->GetEntries();module++){
566     branchRecP->GetEvent(module);
567     npoints += recpoints->GetEntries();
568     AliITSgeomTGeo::GetModuleId(module, lay, lad, det);
569     //printf("modnumb %d, lay %d, lad %d, det %d \n",module, lay, lad, det);
570     for(Int_t j=0;j<recpoints->GetEntries();j++){
571       AliITSRecPoint *recp = (AliITSRecPoint*)recpoints->At(j);    
572       fAliITSQADataMakerRec->GetRecPointsData(6 +fGenRecPointsOffset)->Fill(module);//modpatternrp
573       recp->GetGlobalXYZ(cluglo);
574       Float_t rad=TMath::Sqrt(cluglo[0]*cluglo[0]+cluglo[1]*cluglo[1]); 
575       Float_t phi=TMath::ATan2(cluglo[1],cluglo[0]);
576       fAliITSQADataMakerRec->GetRecPointsData(9 +fGenRecPointsOffset)->Fill(recp->GetDetLocalX(),recp->GetDetLocalZ());//local distribution
577       fAliITSQADataMakerRec->GetRecPointsData(2 +fGenRecPointsOffset)->Fill(cluglo[0],cluglo[1]);//global distribution YX
578       fAliITSQADataMakerRec->GetRecPointsData(3 +fGenRecPointsOffset)->Fill(cluglo[2],rad);//global distribution rz
579       if(fkOnline)
580         {
581           fAliITSQADataMakerRec->GetRecPointsData(14 +fGenRecPointsOffset)->Fill(cluglo[0],cluglo[1]);//global distribution YX FSE
582           fAliITSQADataMakerRec->GetRecPointsData(15 +fGenRecPointsOffset)->Fill(cluglo[2],rad);//global distribution rz FSE
583         }
584       if(recp->GetLayer() ==2) {
585         fAliITSQADataMakerRec->GetRecPointsData(0  +fGenRecPointsOffset)->Fill(recp->GetQ()) ;//total charge of layer 3
586         fAliITSQADataMakerRec->GetRecPointsData(7  +fGenRecPointsOffset)->Fill(lad);//lad pattern layer 3
587         fAliITSQADataMakerRec->GetRecPointsData(10 +fGenRecPointsOffset)->Fill(rad);//r distribution layer 3
588         fAliITSQADataMakerRec->GetRecPointsData(12 +fGenRecPointsOffset)->Fill(phi);// phi distribution layer 3
589         fAliITSQADataMakerRec->GetRecPointsData(4  +fGenRecPointsOffset)->Fill(cluglo[2],phi);// phi distribution layer 3
590       }
591       else if(recp->GetLayer() ==3) {
592         fAliITSQADataMakerRec->GetRecPointsData(1  +fGenRecPointsOffset)->Fill(recp->GetQ()) ;//total charge layer 4
593         fAliITSQADataMakerRec->GetRecPointsData(8  +fGenRecPointsOffset)->Fill(lad);//ladpatternlayer4
594         fAliITSQADataMakerRec->GetRecPointsData(11 +fGenRecPointsOffset)->Fill(rad);//r distribution
595         fAliITSQADataMakerRec->GetRecPointsData(13 +fGenRecPointsOffset)->Fill(phi);//phi distribution
596         fAliITSQADataMakerRec->GetRecPointsData(5  +fGenRecPointsOffset)->Fill(cluglo[2],phi);// phi distribution layer 4
597       }
598     }
599   }
600   statRecpoints.Clear();
601
602 }
603
604 //_______________________________________________________________
605
606 void AliITSQASDDDataMakerRec::SetHLTModeFromEnvironment()
607 {
608
609    Int_t  hltmode= ::atoi(gSystem->Getenv("HLT_MODE"));
610
611    if(hltmode==1)
612      {
613        AliInfo("Online mode: HLT mode A selected from environment for SDD\n");
614        SetHLTMode(kFALSE);
615      }
616    else
617      if(hltmode==2)
618        {
619        AliInfo("Online mode: HLT mode C compressed selected from environment for SDD\n");
620        SetHLTMode(kTRUE);
621        }
622 }
623
624
625 //_______________________________________________________________
626
627 Int_t AliITSQASDDDataMakerRec::GetOffset(AliQAv1::TASKINDEX_t task)
628 {
629   Int_t offset=0;
630   if( task == AliQAv1::kRAWS )
631     {
632       offset=fGenRawsOffset;  
633     }
634   else
635     if( task == AliQAv1::kRECPOINTS )
636       {
637         offset=fGenRecPointsOffset;   
638       }
639     else AliInfo("No task has been selected. Offset set to zero.\n");
640   return offset;
641 }
642
643 //_______________________________________________________________
644
645 Int_t AliITSQASDDDataMakerRec::GetTaskHisto(AliQAv1::TASKINDEX_t task)
646 {
647
648   Int_t histotot=0;
649
650   if( task == AliQAv1::kRAWS )
651     {
652       histotot=fSDDhRawsTask ;  
653     }
654   else
655     if( task == AliQAv1::kRECPOINTS )
656       {
657         histotot=fSDDhRecPointsTask;   
658       }
659     else AliInfo("No task has been selected. TaskHisto set to zero.\n");
660   return histotot;
661 }