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