]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDCalibChamberStatus.cxx
Complete implementation of pools, see #88914. From rev. 53550,53557,53568 (Ruben)
[u/mrichter/AliRoot.git] / TRD / AliTRDCalibChamberStatus.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 ////////////////////////////////////////////////////////////////////////////
17 //                                                                        //
18 // AliTRDCalibChamberStatus: to determine which half chambers are off     //
19 // Produce a AliTRDCalChamberStatus calibration object                    //
20 // Check with the AliTRDCalDCSFEEv2 info                                  //
21 //                                                                        //
22 //                                                                        //
23 // Authors:                                                               //
24 //   J. Book (jbook@ikf.uni-frankfurt.de)                                 //
25 //   R. Bailhache (rbailhache@ikf.uni-frankfurt.de)                       //
26 //                                                                        //
27 ////////////////////////////////////////////////////////////////////////////
28
29
30 //Root includes
31 #include <THnSparse.h>
32
33 #include <TDirectory.h>
34 #include <TFile.h>
35 #include <TAxis.h>
36 #include <TH2F.h>
37 #include <TStyle.h>
38 #include <TCanvas.h>
39
40 //AliRoot includes
41 #include "AliRawReader.h"
42
43 //header file
44 #include "AliLog.h"
45 #include "AliTRDCalibChamberStatus.h"
46 #include "AliTRDgeometry.h"
47 #include "AliTRDfeeParam.h"
48 #include "AliTRDdigitsManager.h"
49 #include "AliTRDSignalIndex.h"
50 #include "AliTRDpadPlane.h"
51 #include "./Cal/AliTRDCalChamberStatus.h"
52 #include "./Cal/AliTRDCalDCSv2.h"
53 #include "./Cal/AliTRDCalDCSFEEv2.h"
54
55 #include "AliTRDrawStream.h"
56 #include "AliTRDseedV1.h"
57 #include "AliTRDcluster.h"
58
59 #ifdef ALI_DATE
60 #include "event.h"
61 #endif
62
63 ClassImp(AliTRDCalibChamberStatus) /*FOLD00*/
64
65 //_____________________________________________________________________
66 AliTRDCalibChamberStatus::AliTRDCalibChamberStatus() : /*FOLD00*/
67   TObject(),
68   fDetector(-1),
69   fNumberOfTimeBins(0),
70   fCounterEventNotEmpty(0),
71   fCalChamberStatus(0x0),
72   fHnSparseI(0x0),
73   fHnSparseHCM(0x0),
74   fHnSparseEvtDet(0x0),
75   fHnSparseDebug(0x0),
76   fHnSparseMCM(0x0),
77   fC1(0x0),
78   fDebugLevel(0)
79 {
80     //
81     // default constructor
82     //
83
84 }
85 //_____________________________________________________________________
86 AliTRDCalibChamberStatus::AliTRDCalibChamberStatus(const AliTRDCalibChamberStatus &ped) : /*FOLD00*/
87   TObject(ped),
88   fDetector(ped.fDetector),
89   fNumberOfTimeBins(ped.fNumberOfTimeBins),
90   fCounterEventNotEmpty(ped.fCounterEventNotEmpty),
91   fCalChamberStatus(ped.fCalChamberStatus),
92   fHnSparseI(ped.fHnSparseI),
93   fHnSparseHCM(ped.fHnSparseHCM),
94   fHnSparseEvtDet(ped.fHnSparseEvtDet),
95   fHnSparseDebug(ped.fHnSparseDebug),  
96   fHnSparseMCM(ped.fHnSparseMCM),
97   fC1(ped.fC1),
98   fDebugLevel(ped.fDebugLevel)
99 {
100     //
101     // copy constructor
102     //
103   
104 }
105 //_____________________________________________________________________
106 AliTRDCalibChamberStatus& AliTRDCalibChamberStatus::operator = (const  AliTRDCalibChamberStatus &source)
107 {
108   //
109   // assignment operator
110   //
111   if (&source == this) return *this;
112   new (this) AliTRDCalibChamberStatus(source);
113
114   return *this;
115 }
116 //_____________________________________________________________________
117 AliTRDCalibChamberStatus::~AliTRDCalibChamberStatus() /*FOLD00*/
118 {
119   //
120   // destructor
121   //
122   if(fCalChamberStatus){
123     delete fCalChamberStatus;
124   }
125   if(fHnSparseI) {
126     delete fHnSparseI;
127   }
128   if(fHnSparseHCM) {
129     delete fHnSparseHCM;
130   }
131   if(fHnSparseEvtDet) {
132     delete fHnSparseEvtDet;
133   }
134   if(fHnSparseDebug) {
135     delete fHnSparseDebug;
136   }
137   if(fHnSparseMCM) {
138    delete fHnSparseMCM;
139   }
140   if(fC1) {
141    delete fC1;
142   }
143 }
144
145 //_____________________________________________________________________
146 void AliTRDCalibChamberStatus::Init() 
147 {
148   //
149   // Init the different THnSparse
150   //
151   //
152
153   //
154   // Init the fHnSparseI
155   //
156
157   //create the map
158   Int_t thnDimEvt[4]; // sm, layer, stack, halfchamber
159   thnDimEvt[0] = 18;
160   thnDimEvt[1] = 6;
161   thnDimEvt[2] = 5;
162   thnDimEvt[3] = 2;
163   //arrays for lower bounds :
164   Double_t* binEdgesEvt[4];
165   for(Int_t ivar = 0; ivar < 4; ivar++)
166     binEdgesEvt[ivar] = new Double_t[thnDimEvt[ivar] + 1];
167   //values for bin lower bounds
168   for(Int_t i=0; i<=thnDimEvt[0]; i++) binEdgesEvt[0][i]= 0.0  + (18.0)/thnDimEvt[0]*(Double_t)i;
169   for(Int_t i=0; i<=thnDimEvt[1]; i++) binEdgesEvt[1][i]= 0.0  + (6.0)/thnDimEvt[1]*(Double_t)i;
170   for(Int_t i=0; i<=thnDimEvt[2]; i++) binEdgesEvt[2][i]= 0.0  + (5.0)/thnDimEvt[2]*(Double_t)i;
171   for(Int_t i=0; i<=thnDimEvt[3]; i++) binEdgesEvt[3][i]= 0.0  + (2.0)/thnDimEvt[3]*(Double_t)i;
172   
173   //create the THnSparse
174   fHnSparseI = new THnSparseI("NumberOfEntries","NumberOfEntries",4,thnDimEvt);
175   for (int k=0; k<4; k++) {
176     fHnSparseI->SetBinEdges(k,binEdgesEvt[k]);
177   }
178   fHnSparseI->Sumw2();
179
180   //
181   // Init the fHnSparseHCM (THnSparseI)
182   //
183
184   //create the THnSparse
185   fHnSparseHCM = new THnSparseI("HCMerrors","HCMerrors",4,thnDimEvt);
186   for (int k=0; k<4; k++) {
187     fHnSparseHCM->SetBinEdges(k,binEdgesEvt[k]);
188   }
189   fHnSparseHCM->Sumw2();
190
191
192   //---------//
193   //  Debug  //
194   if(fDebugLevel > 0) {
195   
196     //
197     // Init the fHnSparseEvtDet (THnSparseI)
198     //
199     
200     //create the map
201     Int_t thnDimEvts[3]; // event, detector, halfchamber
202     thnDimEvts[0] = 10000;
203     thnDimEvts[1] = 540;
204     thnDimEvts[2] = 2;
205     //arrays for lower bounds :
206     Double_t* binEdgesEvts[3];
207     for(Int_t ivar = 0; ivar < 3; ivar++)
208       binEdgesEvts[ivar] = new Double_t[thnDimEvts[ivar] + 1];
209     //values for bin lower bounds
210     for(Int_t i=0; i<=thnDimEvts[0]; i++) binEdgesEvts[0][i]= 0.0  + (10000.0)/thnDimEvts[0]*(Double_t)i;
211     for(Int_t i=0; i<=thnDimEvts[1]; i++) binEdgesEvts[1][i]= 0.0  + (540.0)/thnDimEvts[1]*(Double_t)i;
212     for(Int_t i=0; i<=thnDimEvts[2]; i++) binEdgesEvts[2][i]= 0.0  + (2.0)/thnDimEvts[2]*(Double_t)i;
213     
214     //create the THnSparse
215     fHnSparseEvtDet = new THnSparseI("NumberOfEntriesPerEvent","NumberOfEntriesPerEvent",3,thnDimEvts);
216     for (int k=0; k<3; k++) {
217       fHnSparseEvtDet->SetBinEdges(k,binEdgesEvts[k]);
218     }
219     fHnSparseEvtDet->Sumw2();
220
221     //
222     // Init the fHnSparseDebug (THnSparseI)
223     //
224     
225     //create the THnSparse
226     fHnSparseDebug = new THnSparseI("NumberOfDifferentDecisions","NumberOfDifferentDecisions",4,thnDimEvt);
227     for (int k=0; k<4; k++) {
228       fHnSparseDebug->SetBinEdges(k,binEdgesEvt[k]);
229     }
230     fHnSparseDebug->Sumw2();
231         
232     //
233     // Init the fHnSparseMCM (THnSparseI)
234     //
235     
236     //create the map
237     Int_t thnDimEvtt[6]; // sm, layer, stack, ROB, MCM
238     thnDimEvtt[0] = 18;
239     thnDimEvtt[1] = 6;
240     thnDimEvtt[2] = 5;
241     thnDimEvtt[3] = 8;
242     thnDimEvtt[4] = 18;
243     thnDimEvtt[5] = 18;
244     //arrays for lower bounds :
245     Double_t* binEdgesEvtt[6];
246     for(Int_t ivar = 0; ivar < 6; ivar++)
247       binEdgesEvtt[ivar] = new Double_t[thnDimEvtt[ivar] + 1];
248     //values for bin lower bounds
249     for(Int_t i=0; i<=thnDimEvtt[0]; i++) binEdgesEvtt[0][i]= 0.0  + (18.0)/thnDimEvtt[0]*(Double_t)i;
250     for(Int_t i=0; i<=thnDimEvtt[1]; i++) binEdgesEvtt[1][i]= 0.0  + (6.0)/thnDimEvtt[1]*(Double_t)i;
251     for(Int_t i=0; i<=thnDimEvtt[2]; i++) binEdgesEvtt[2][i]= 0.0  + (5.0)/thnDimEvtt[2]*(Double_t)i;
252     for(Int_t i=0; i<=thnDimEvtt[3]; i++) binEdgesEvtt[3][i]= 0.0  + (8.0)/thnDimEvtt[3]*(Double_t)i;
253     for(Int_t i=0; i<=thnDimEvtt[4]; i++) binEdgesEvtt[4][i]= 0.0  + (18.0)/thnDimEvtt[4]*(Double_t)i;
254     for(Int_t i=0; i<=thnDimEvtt[5]; i++) binEdgesEvtt[5][i]= 0.0  + (18.0)/thnDimEvtt[5]*(Double_t)i;
255     
256     //create the THnSparse
257     fHnSparseMCM = new THnSparseI("MCMerrorDCS","MCMerrorDCS",6,thnDimEvtt);
258     for (int k=0; k<6; k++) {
259       fHnSparseMCM->SetBinEdges(k,binEdgesEvtt[k]);
260     }
261     fHnSparseMCM->Sumw2();
262   
263   }
264   //  Debug  //
265   //---------//
266
267 }
268 //_____________________________________________________________________
269 void AliTRDCalibChamberStatus::ProcessTrack(AliTRDtrackV1 * trdTrack)
270 {
271   //
272   // Track Processing to get half chamber status
273   //
274   //
275   
276   
277   const AliTRDseedV1 *tracklet = 0x0;
278   AliTRDcluster *cluster;
279   //////////////////////////////////////
280   // Loop tracklets
281   ///////////////////////////////////// 
282   for(Int_t itr = 0; itr < 6; ++itr){
283         
284     if(!(tracklet = trdTrack->GetTracklet(itr))) continue;
285     if(!tracklet->IsOK()) continue;
286   
287     // Loop on clusters
288     for(int ic=0; ic<AliTRDseedV1::kNtb; ++ic){
289       if((cluster = tracklet->GetClusters(ic))) {
290         //printf("ic %d\n",ic);
291         break;
292       }
293     }
294     if(!cluster) continue;
295     
296     Int_t det     = cluster->GetDetector();       
297     Int_t layer   = AliTRDgeometry::GetLayer(det);
298     Int_t sm      = AliTRDgeometry::GetSector(det);
299     Int_t stac    = AliTRDgeometry::GetStack(det);
300     
301     Int_t col     = cluster->GetPadCol();
302     Int_t iMcm    = (Int_t)(col/18);   // current group of 18 col pads
303     Double_t rphi = 0.5;
304     if(iMcm > 3) rphi = 1.5;
305           
306     Double_t val[4] = {sm,layer,stac,rphi}; 
307     if(fHnSparseI->GetBinContent((const Int_t*)val)<2147483646) fHnSparseI->Fill(&val[0]); 
308   }
309   
310 }
311 //_____________________________________________________________________
312 void AliTRDCalibChamberStatus::ProcessEvent(AliRawReader * rawReader, Int_t nevents_physics)
313 {
314   //
315   // Event Processing loop with AliTRDrawStream
316   //
317   //
318   
319   Bool_t notEmpty = kFALSE;
320     
321   AliTRDdigitsManager *digitsManager = new AliTRDdigitsManager(kTRUE);
322   digitsManager->CreateArrays();
323   
324   AliTRDrawStream *rawStream = new AliTRDrawStream(rawReader);
325   rawStream->SetDigitsManager(digitsManager);
326   //  rawStream->SetNoErrorWarning();
327   //  rawStream->SetSharedPadReadout(kFALSE);
328
329   
330   Int_t det    = 0;
331   while ((det = rawStream->NextChamber(digitsManager, NULL, NULL)) >= 0) { 
332
333     //nextchamber loop
334     
335     // do the QA analysis
336     if (digitsManager->GetIndexes(det)->HasEntry()) {//QA
337       // printf("there is ADC data on this chamber!\n");
338       
339       AliTRDSignalIndex *indexes = digitsManager->GetIndexes(det);
340       if (indexes->IsAllocated() == kFALSE) {
341         // AliError("Indexes do not exist!");
342         break;
343       }
344       
345       Int_t iRow  = 0;
346       Int_t iCol  = 0;
347       indexes->ResetCounters();
348       
349       while (indexes->NextRCIndex(iRow, iCol)){
350         Int_t iMcm        = (Int_t)(iCol/18);   // current group of 18 col pads
351         
352         Int_t layer = AliTRDgeometry::GetLayer(det);
353         Int_t sm    = AliTRDgeometry::GetSector(det);
354         Int_t stac  = AliTRDgeometry::GetStack(det);
355         Double_t rphi = 0.5;
356         if(iMcm > 3) rphi = 1.5;
357
358         Double_t val[4] = {sm,layer,stac,rphi}; 
359         fHnSparseI->Fill(&val[0]); 
360         notEmpty = kTRUE;
361         
362         //---------//
363         //  Debug  //
364         if(fDebugLevel > 0) {
365           Int_t detector = AliTRDgeometry::GetDetector(layer,stac,sm);
366           Double_t valu[3] = {nevents_physics,detector,rphi};
367           fHnSparseEvtDet->Fill(&valu[0]); 
368         }
369         //  Debug  //
370         //---------//
371       }
372       
373     }
374     digitsManager->ClearArrays(det);
375   }
376
377   if(notEmpty) fCounterEventNotEmpty++;
378
379   delete digitsManager;
380   delete rawStream;
381    
382 }
383 //_____________________________________________________________________
384 Bool_t AliTRDCalibChamberStatus::TestEventHisto(Int_t nevent) /*FOLD00*/
385 {
386   //
387   //  Test event loop
388   // fill the fHnSparseI with entries
389   //
390   
391   AliTRDgeometry geo;
392
393
394   for(Int_t ievent=0; ievent<nevent; ievent++){
395     for (Int_t ism=0; ism<18; ism++){
396       for (Int_t istack=0; istack<5; istack++){
397         for (Int_t ipl=0; ipl<6; ipl++){
398           for (Int_t icol=0; icol<geo.GetColMax(ipl); icol++){
399             Int_t side = 0;
400             if(icol > 72) side = 1;
401             Double_t val[4] = {ism,ipl,istack,side}; 
402             fHnSparseI->Fill(&val[0]); 
403           }
404         }
405       }
406     }
407   }
408   
409   return kTRUE;
410
411 }
412 //_____________________________________________________________________
413 void AliTRDCalibChamberStatus::AnalyseHisto(Int_t limit) /*FOLD00*/
414 {
415   //
416   //  Create the AliTRDCalChamberStatus according to the fHnSparseI
417   //
418   
419   if(fCalChamberStatus) delete fCalChamberStatus;
420   fCalChamberStatus = new AliTRDCalChamberStatus();
421
422   //printf("test0\n");
423
424   // Check if enough events/tracklets per halfchamber to say something
425   Double_t mean=0.0; //number of tracklets per HCS
426   Int_t coord2[4];
427   for(Int_t bin = 0; bin < fHnSparseI->GetNbins(); bin++) {
428     //if(fHnSparseI->GetBinContent(bin,coord2)==0.0) printf(" bin shouldnt be empty!!\n");
429     mean+=fHnSparseI->GetBinContent(bin,coord2);
430   }
431   if(fHnSparseI->GetNbins() > 0.0) mean/=fHnSparseI->GetNbins();
432   //printf(" mean tracklets per halfchamber %f \n",mean);
433   if((fCounterEventNotEmpty < limit) && (mean < limit)) {
434     // Say all good
435     for (Int_t idet=0; idet<540; idet++) {
436       fCalChamberStatus->SetStatus(idet,AliTRDCalChamberStatus::kGood);
437     }
438     return;
439   }
440
441   //printf("test1\n");
442
443   // set all chambers to NoData
444   for (Int_t idet=0; idet<540; idet++) {
445     fCalChamberStatus->SetStatus(idet,AliTRDCalChamberStatus::kNoData);
446   }
447
448   // set status according to fHnSparseI 
449   Int_t coord[4];
450   for(Int_t bin = 0; bin < fHnSparseI->GetNbins(); bin++) {
451     
452     //Double_t content = fHnSparseI->GetBinContent(bin,coord);
453     fHnSparseI->GetBinContent(bin,coord);
454     // layer, stack, sector
455     Int_t detector = AliTRDgeometry::GetDetector(coord[1]-1,coord[2]-1,coord[0]-1);
456     //
457     //printf("Number of entries for detector %d: %f\n",detector,content);
458     //
459     // Check which halfchamber side corresponds to the bin number (0=A, 1=B)
460     // Change the status accordingly
461     //
462     if(coord[3]-1==0) { // HCS-A
463       fCalChamberStatus->SetStatus(detector,AliTRDCalChamberStatus::kGood);
464       fCalChamberStatus->UnsetStatusBit(detector,AliTRDCalChamberStatus::kNoDataHalfChamberSideA); // A has data
465       //fCalChamberStatus->UnsetStatusBit(detector,AliTRDCalChamberStatus::kNoData); 
466     }
467     else { //HCS-B
468       fCalChamberStatus->SetStatus(detector,AliTRDCalChamberStatus::kGood);
469       fCalChamberStatus->UnsetStatusBit(detector,AliTRDCalChamberStatus::kNoDataHalfChamberSideB); // B has data
470       //fCalChamberStatus->UnsetStatusBit(detector,AliTRDCalChamberStatus::kNoData); 
471     }
472   }
473
474   //printf("test2\n");
475
476   // printf
477   //for (Int_t idet=0; idet<540; idet++) {
478   //  if(fCalChamberStatus->IsNoData(idet)) printf("No Data: chamber %d\n",idet);
479   //}
480
481
482 }
483 //_____________________________________________________________________
484 void AliTRDCalibChamberStatus::CheckEORStatus(AliTRDCalDCSv2 *calDCS) /*FOLD00*/
485 {
486   //
487   //  Correct the AliTRDCalChamberStatus according to the AliTRDCalDCSv2
488   //  Using globale state of the HalfChamberMerger (HCM)
489   //
490   for(Int_t det = 0; det < 540; det++) {
491     AliTRDCalDCSFEEv2* calDCSFEEEOR = calDCS->GetCalDCSFEEObj(det);
492
493     if(!calDCSFEEEOR) continue;
494     
495     // MCM Global State Machine State Definitions
496     //  low_power =  0,
497     //  test      =  1,
498     //  wait_pre  =  3,
499     //  preproc   =  7,
500     //  zero_sp   =  8,
501     //  full_rd   =  9,
502     //  clear_st  = 11,
503     //  wait_L1   = 12,
504     //  tr_send   = 14,
505     //  tr_proc   = 15 
506
507     Int_t sm   = AliTRDgeometry::GetSector(det);
508     Int_t lay  = AliTRDgeometry::GetLayer(det);
509     Int_t stac = AliTRDgeometry::GetStack(det);
510     
511     Int_t stateA = 0;   // 0=bad, 1=good state
512     Int_t stateB = 0;
513
514     // loop over all mcm to define DCS-HCS
515     for(Int_t ii = 0; ii < 8; ii++) { //ROB loop
516       for(Int_t i = 0; i < 18; i++) { //MCM loop
517         
518         Int_t side = ii%2;  // 0=sideA, 1=sideB
519         Int_t cstate = calDCSFEEEOR->GetMCMGlobalState(ii,i); //current mcm state
520         
521         if(cstate==3) {
522           switch(side) {
523           case 0: stateA=1; break;
524           case 1: stateB=1; break; 
525           }
526         }
527       }
528     }
529
530     //---------//
531     //  Debug  //
532     if(fDebugLevel > 0) {
533       if( ((fCalChamberStatus->GetStatus(det) <= 1) && (stateA==0 || stateB==0)) || 
534           ((fCalChamberStatus->GetStatus(det) == 2) && (stateA==1 || stateB==1)) || 
535           ((fCalChamberStatus->GetStatus(det) == 3) && (stateA==1 || stateB==0)) ||
536           ((fCalChamberStatus->GetStatus(det) == 4) && (stateB==0 || stateB==1))  )
537         {
538           //printf(" Different half chamber status in DCS and DATA!!\n");
539           Double_t val[4] = {sm,lay,stac,1};
540           fHnSparseDebug->Fill(&val[0]); 
541           
542           // Fill MCM status map
543           for(Int_t ii = 0; ii < 8; ii++) { //ROB loop
544             for(Int_t i = 0; i < 18; i++) { //MCM loop
545               Double_t valss[6] = {sm,lay,stac,ii,i
546                                    ,calDCSFEEEOR->GetMCMGlobalState(ii,i)};
547               fHnSparseMCM->Fill(&valss[0]);
548               
549             } 
550           }
551         }
552     }
553     //---------//
554     //  Debug  //
555
556     //---------------------------------------
557     // Change the status according to DCS
558     //---------------------------------------
559     Int_t StatusData = fCalChamberStatus->GetStatus(det);
560     switch(StatusData) 
561       {
562       case 1: 
563         if(stateA==0 && stateB==0) fCalChamberStatus->SetStatus(det,2); // completely masked from DCS
564         if(stateA==1 && stateB==0) fCalChamberStatus->SetStatus(det,4); // Only B side masked from DCS
565         if(stateA==0 && stateB==1) fCalChamberStatus->SetStatus(det,3); // Only A side masked from DCS
566         if(stateA==1 && stateB==1) fCalChamberStatus->SetStatus(det,1);
567         break;
568       case 2: // completely masked from DATA
569         if(stateA==0 && stateB==0) fCalChamberStatus->SetStatus(det,2); // completely masked from DCS
570         break;
571       case 3: // Only A side masked from DATA
572         if(stateA==0 && stateB==0) fCalChamberStatus->SetStatus(det,2); // completely masked from DCS
573         if(stateA==1 && stateB==0) fCalChamberStatus->SetStatus(det,2); // Only B side masked from DCS
574         if(stateA==0 && stateB==1) fCalChamberStatus->SetStatus(det,3); // Only A side masked from DCS
575         if(stateA==1 && stateB==1) fCalChamberStatus->SetStatus(det,3);
576         break;
577       case 4: // Only B side masked from DATA
578         if(stateA==0 && stateB==0) fCalChamberStatus->SetStatus(det,2); // completely masked from DCS
579         if(stateA==1 && stateB==0) fCalChamberStatus->SetStatus(det,4); // Only B side masked from DCS
580         if(stateA==0 && stateB==1) fCalChamberStatus->SetStatus(det,2); // Only A side masked from DCS
581         if(stateA==1 && stateB==1) fCalChamberStatus->SetStatus(det,4);
582         break;
583       }
584     
585   }
586
587 }
588
589 //_____________________________________________________________________________________
590 void AliTRDCalibChamberStatus::Add(AliTRDCalibChamberStatus *calibChamberStatus) /*FOLD00*/
591 {
592     //
593     //  Add the THnSparseI of this calibChamberStatus
594     //
595
596   fCounterEventNotEmpty += calibChamberStatus->GetNumberEventNotEmpty();
597
598   THnSparseI *hnSparseI = calibChamberStatus->GetSparseI();
599   if(!hnSparseI) return;
600
601   if(!fHnSparseI) {
602     fHnSparseI = (THnSparseI *) hnSparseI->Clone();
603   }
604   else {
605     fHnSparseI->Add(hnSparseI);
606   }
607   
608
609 }
610 //_____________________________________________________________________
611 void AliTRDCalibChamberStatus::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
612 {
613     //
614     //  Write class to file
615     //
616
617     TString sDir(dir);
618     TString option;
619
620     if ( append )
621         option = "update";
622     else
623         option = "recreate";
624
625     TDirectory *backup = gDirectory;
626     TFile f(filename,option.Data());
627     f.cd();
628     if ( !sDir.IsNull() ){
629         f.mkdir(sDir.Data());
630         f.cd(sDir);
631     }
632     this->Write();
633     f.Close();
634
635     if ( backup ) backup->cd();
636 }
637 //_____________________________________________________________________________
638 TH2D* AliTRDCalibChamberStatus::PlotSparseI(Int_t sm,Int_t side) 
639 {
640   //
641   // Plot number of entries for supermodule sm 
642   // as a function of layer and stack
643   //
644
645   if(!fHnSparseI) return 0x0;
646   
647   fHnSparseI->GetAxis(0)->SetRange(sm+1,sm+1); 
648   fHnSparseI->GetAxis(3)->SetRange(side+1,side+1);  
649   TH2D *h2 = fHnSparseI->Projection(1,2);
650  
651
652   return h2;
653
654 }
655 //_____________________________________________________________________
656 TH2F *AliTRDCalibChamberStatus::MakeHisto2DSmPlEORStatus(AliTRDCalDCSv2 *calDCS, Int_t sm, Int_t pl) /*FOLD00*/
657 {
658   //
659   //  Plot globale state of the HalfChamberMerger (HCM)
660   //
661   AliTRDfeeParam *paramfee = AliTRDfeeParam::Instance();
662
663   AliTRDgeometry *trdGeo = new AliTRDgeometry();
664   AliTRDpadPlane *padPlane0 = trdGeo->GetPadPlane(pl,0);        // layer,stack
665   Double_t row0    = padPlane0->GetRow0();
666   Double_t col0    = padPlane0->GetCol0();
667
668   char  name[1000];
669   snprintf(name,1000,"%s DCS status sm %d pl %d",GetTitle(),sm,pl);
670   TH2F * his = new TH2F( name, name, 88,-TMath::Abs(row0),TMath::Abs(row0)
671                                    ,148,-TMath::Abs(col0),TMath::Abs(col0));
672
673
674   // Where we begin
675   Int_t offsetsmpl = 30*sm+pl;
676   Int_t nstack = 5;
677   Int_t ncols = 144;
678
679   for (Int_t k = 0; k < nstack; k++){
680     Int_t det = offsetsmpl+k*6;
681     Int_t stac = AliTRDgeometry::GetStack(det);
682     AliTRDCalDCSFEEv2* calDCSFEEEOR = calDCS->GetCalDCSFEEObj(det);
683     if(!calDCSFEEEOR) { continue;}
684     for (Int_t icol=0; icol<ncols; icol++){
685       Int_t nrows = 16;
686       if(stac==2) nrows = 12;
687       for (Int_t irow=0; irow<nrows; irow++){
688         Int_t binz     = 0;
689         Int_t kb       = 5-1-k;
690         Int_t krow     = nrows-1-irow;
691         Int_t kcol     = ncols-1-icol;
692         if(kb > 2) binz = 16*(kb-1)+12+krow+1+2*(kb+1);
693         else binz = 16*kb+krow+1+2*(kb+1); 
694         Int_t biny = kcol+1+2;
695         // Take the value
696         Int_t mcm = paramfee->GetMCMfromPad(irow,icol);
697         Int_t rob = paramfee->GetROBfromPad(irow,icol);
698         if(mcm < 0) AliWarning("Problem with mcm number");
699         Int_t state = calDCSFEEEOR->GetMCMGlobalState(rob,TMath::Abs(mcm)); 
700         his->SetBinContent(binz,biny,state);
701       }
702     }
703     for(Int_t icol = 1; icol < 147; icol++){
704       for(Int_t l = 0; l < 2; l++){
705         Int_t binz     = 0;
706         Int_t kb       = 5-1-k;
707         if(kb > 2) binz = 16*(kb-1)+12+1+2*(kb+1)-(l+1);
708         else binz = 16*kb+1+2*(kb+1)-(l+1); 
709         his->SetBinContent(binz,icol,16.0);
710       }
711     }
712   }
713   
714   for(Int_t icol = 1; icol < 147; icol++){
715     his->SetBinContent(88,icol,16.0);
716     his->SetBinContent(87,icol,16.0);
717   }
718   for(Int_t irow = 1; irow < 89; irow++){
719     his->SetBinContent(irow,1,16.0);
720     his->SetBinContent(irow,2,16.0);
721     his->SetBinContent(irow,147,16.0);
722     his->SetBinContent(irow,148,16.0);
723   }
724
725   his->SetXTitle("z (cm)");
726   his->SetYTitle("y (cm)");
727   his->SetMaximum(12);
728   his->SetMinimum(0.0);
729   his->SetStats(0);
730
731   return his;
732
733 }
734 //_____________________________________________________________________________
735 TCanvas* AliTRDCalibChamberStatus::PlotHistos2DSmEORStatus(AliTRDCalDCSv2 *calDCS, Int_t sm, const Char_t *name)
736 {
737   //
738   // Make 2D graph
739   //
740
741   gStyle->SetPalette(1);
742   fC1 = new TCanvas(name,name,50,50,600,800);
743   fC1->Divide(3,2);
744   fC1->cd(1);
745   MakeHisto2DSmPlEORStatus(calDCS,sm,0)->Draw("colz");
746   fC1->cd(2);
747   MakeHisto2DSmPlEORStatus(calDCS,sm,1)->Draw("colz");
748   fC1->cd(3);
749   MakeHisto2DSmPlEORStatus(calDCS,sm,2)->Draw("colz");
750   fC1->cd(4);
751   MakeHisto2DSmPlEORStatus(calDCS,sm,3)->Draw("colz");
752   fC1->cd(5);
753   MakeHisto2DSmPlEORStatus(calDCS,sm,4)->Draw("colz");
754   fC1->cd(6);
755   MakeHisto2DSmPlEORStatus(calDCS,sm,5)->Draw("colz");
756
757   return fC1;
758
759 }
760
761
762