]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDCalibChamberStatus.cxx
1. Adding the new histograms TPC z vertex correlation in order
[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   if (!rawStream) return;
279   rawStream->SetNoErrorWarning();
280   rawStream->SetSharedPadReadout(kFALSE);
281
282   AliTRDdigitsManager *digitsManager = new AliTRDdigitsManager(kTRUE);
283   if (!digitsManager) return;
284   digitsManager->CreateArrays();
285   
286   Int_t det    = 0;
287   while ((det = rawStream->NextChamber(digitsManager, NULL, NULL)) >= 0) { 
288
289     //nextchamber loop
290     
291     // do the QA analysis
292     if (digitsManager->GetIndexes(det)->HasEntry()) {//QA
293       // printf("there is ADC data on this chamber!\n");
294       
295       AliTRDSignalIndex *indexes = digitsManager->GetIndexes(det);
296       if (indexes->IsAllocated() == kFALSE) {
297         // AliError("Indexes do not exist!");
298         break;
299       }
300       
301       Int_t iRow  = 0;
302       Int_t iCol  = 0;
303       indexes->ResetCounters();
304       
305       while (indexes->NextRCIndex(iRow, iCol)){
306         Int_t iMcm        = (Int_t)(iCol/18);   // current group of 18 col pads
307         
308         Int_t layer = AliTRDgeometry::GetLayer(det);
309         Int_t sm    = AliTRDgeometry::GetSector(det);
310         Int_t stac  = AliTRDgeometry::GetStack(det);
311         Double_t rphi = 0.5;
312         if(iMcm > 3) rphi = 1.5;
313
314         Double_t val[4] = {sm,layer,stac,rphi}; 
315         fHnSparseI->Fill(&val[0]); 
316         notEmpty = kTRUE;
317         
318         //---------//
319         //  Debug  //
320         if(fDebugLevel > 0) {
321           Int_t detector = AliTRDgeometry::GetDetector(layer,stac,sm);
322           Double_t valu[3] = {nevents_physics,detector,rphi};
323           fHnSparseEvtDet->Fill(&valu[0]); 
324         }
325         //  Debug  //
326         //---------//
327       }
328       
329     }
330     digitsManager->ClearArrays(det);
331   }
332
333   if(notEmpty) fCounterEventNotEmpty++;
334
335   delete digitsManager;
336   delete rawStream;
337    
338 }//_____________________________________________________________________
339 void AliTRDCalibChamberStatus::ProcessEvent3(AliRawReader * rawReader, Int_t nevents_physics)
340 {
341   //
342   // Event Processing loop with AliTRDrawStream
343   //
344   //
345   
346   Bool_t notEmpty = kFALSE;
347     
348   AliTRDdigitsManager *digitsManager = new AliTRDdigitsManager(kTRUE);
349   digitsManager->CreateArrays();
350   
351   AliTRDrawStream *rawStream = new AliTRDrawStream(rawReader);
352   rawStream->SetDigitsManager(digitsManager);
353   //  rawStream->SetNoErrorWarning();
354   //  rawStream->SetSharedPadReadout(kFALSE);
355
356   
357   Int_t det    = 0;
358   while ((det = rawStream->NextChamber(digitsManager, NULL, NULL)) >= 0) { 
359
360     //nextchamber loop
361     
362     // do the QA analysis
363     if (digitsManager->GetIndexes(det)->HasEntry()) {//QA
364       // printf("there is ADC data on this chamber!\n");
365       
366       AliTRDSignalIndex *indexes = digitsManager->GetIndexes(det);
367       if (indexes->IsAllocated() == kFALSE) {
368         // AliError("Indexes do not exist!");
369         break;
370       }
371       
372       Int_t iRow  = 0;
373       Int_t iCol  = 0;
374       indexes->ResetCounters();
375       
376       while (indexes->NextRCIndex(iRow, iCol)){
377         Int_t iMcm        = (Int_t)(iCol/18);   // current group of 18 col pads
378         
379         Int_t layer = AliTRDgeometry::GetLayer(det);
380         Int_t sm    = AliTRDgeometry::GetSector(det);
381         Int_t stac  = AliTRDgeometry::GetStack(det);
382         Double_t rphi = 0.5;
383         if(iMcm > 3) rphi = 1.5;
384
385         Double_t val[4] = {sm,layer,stac,rphi}; 
386         fHnSparseI->Fill(&val[0]); 
387         notEmpty = kTRUE;
388         
389         //---------//
390         //  Debug  //
391         if(fDebugLevel > 0) {
392           Int_t detector = AliTRDgeometry::GetDetector(layer,stac,sm);
393           Double_t valu[3] = {nevents_physics,detector,rphi};
394           fHnSparseEvtDet->Fill(&valu[0]); 
395         }
396         //  Debug  //
397         //---------//
398       }
399       
400     }
401     digitsManager->ClearArrays(det);
402   }
403
404   if(notEmpty) fCounterEventNotEmpty++;
405
406   delete digitsManager;
407   delete rawStream;
408    
409 }
410 //_____________________________________________________________________
411 Bool_t AliTRDCalibChamberStatus::TestEventHisto(Int_t nevent) /*FOLD00*/
412 {
413   //
414   //  Test event loop
415   // fill the fHnSparseI with entries
416   //
417   
418   AliTRDgeometry geo;
419
420
421   for(Int_t ievent=0; ievent<nevent; ievent++){
422     for (Int_t ism=0; ism<18; ism++){
423       for (Int_t istack=0; istack<5; istack++){
424         for (Int_t ipl=0; ipl<6; ipl++){
425           for (Int_t icol=0; icol<geo.GetColMax(ipl); icol++){
426             Int_t side = 0;
427             if(icol > 72) side = 1;
428             Double_t val[4] = {ism,ipl,istack,side}; 
429             fHnSparseI->Fill(&val[0]); 
430           }
431         }
432       }
433     }
434   }
435   
436   return kTRUE;
437
438 }
439 //_____________________________________________________________________
440 void AliTRDCalibChamberStatus::AnalyseHisto() /*FOLD00*/
441 {
442   //
443   //  Create the AliTRDCalChamberStatus according to the fHnSparseI
444   //
445   
446   if(fCalChamberStatus) delete fCalChamberStatus;
447   fCalChamberStatus = new AliTRDCalChamberStatus();
448
449   // Check if enough events to say something
450   if(fCounterEventNotEmpty < 200) {
451     // Say all installed
452     for (Int_t ism=0; ism<18; ism++) {
453       for (Int_t ipl=0; ipl<6; ipl++) {
454         for (Int_t istack=0; istack<5; istack++) {
455           // layer, stack, sector
456           Int_t det = AliTRDgeometry::GetDetector(ipl,istack,ism);
457           fCalChamberStatus->SetStatus(det,1);
458         }
459       }
460     }
461     return;
462   }
463
464   // Mask out all chambers
465   for (Int_t ism=0; ism<18; ism++) {
466     for (Int_t ipl=0; ipl<6; ipl++) {
467       for (Int_t istack=0; istack<5; istack++) {
468         // layer, stack, sector
469         Int_t det = AliTRDgeometry::GetDetector(ipl,istack,ism);
470         fCalChamberStatus->SetStatus(det,2);
471       }
472     }
473   }
474
475   // Unmask good chambers 
476   Int_t coord[4];
477   for(Int_t bin = 0; bin < fHnSparseI->GetNbins(); bin++) {
478     
479     fHnSparseI->GetBinContent(bin,coord);
480     // layer, stack, sector
481     Int_t detector = AliTRDgeometry::GetDetector(coord[1]-1,coord[2]-1,coord[0]-1);
482
483     //
484     // Check which halfchamber side corresponds to the bin number (0=A, 1=B)
485     // Change the status accordingly
486     //
487
488     switch(fCalChamberStatus->GetStatus(detector)) 
489       {    
490       case 1: break;  // no changes
491       case 2: 
492         if(coord[3]-1==0) {
493           fCalChamberStatus->SetStatus(detector,4); break;      // only SideB is masked
494         }
495         else {
496           fCalChamberStatus->SetStatus(detector,3); break;      // only SideA is masked
497         }
498       case 3:  fCalChamberStatus->SetStatus(detector,1); break;  // unmask SideA
499       case 4:  fCalChamberStatus->SetStatus(detector,1); break;  // unmask SideB
500       }
501   }
502
503
504 }
505 //_____________________________________________________________________
506 void AliTRDCalibChamberStatus::CheckEORStatus(AliTRDCalDCS *calDCS) /*FOLD00*/
507 {
508   //
509   //  Correct the AliTRDCalChamberStatus according to the AliTRDCalDCS
510   //  Using globale state of the HalfChamberMerger (HCM)
511   //
512   for(Int_t det = 0; det < 540; det++) {
513     AliTRDCalDCSFEE* calDCSFEEEOR = calDCS->GetCalDCSFEEObj(det);
514
515     if(!calDCSFEEEOR) continue;
516     
517     // MCM Global State Machine State Definitions
518     //  low_power =  0,
519     //  test      =  1,
520     //  wait_pre  =  3,
521     //  preproc   =  7,
522     //  zero_sp   =  8,
523     //  full_rd   =  9,
524     //  clear_st  = 11,
525     //  wait_L1   = 12,
526     //  tr_send   = 14,
527     //  tr_proc   = 15 
528
529     Int_t sm   = AliTRDgeometry::GetSector(det);
530     Int_t lay  = AliTRDgeometry::GetLayer(det);
531     Int_t stac = AliTRDgeometry::GetStack(det);
532     
533     Int_t stateA = 0;   // 0=bad, 1=good state
534     Int_t stateB = 0;
535
536     // loop over all mcm to define DCS-HCS
537     for(Int_t ii = 0; ii < 8; ii++) { //ROB loop
538       for(Int_t i = 0; i < 18; i++) { //MCM loop
539         
540         Int_t side = ii%2;  // 0=sideA, 1=sideB
541         Int_t cstate = calDCSFEEEOR->GetMCMGlobalState(ii,i); //current mcm state
542         
543         if(cstate==3) {
544           switch(side) {
545           case 0: stateA=1; break;
546           case 1: stateB=1; break; 
547           }
548         }
549       }
550     }
551
552     //---------//
553     //  Debug  //
554     if(fDebugLevel > 0) {
555       if( ((fCalChamberStatus->GetStatus(det) <= 1) && (stateA==0 || stateB==0)) || 
556           ((fCalChamberStatus->GetStatus(det) == 2) && (stateA==1 || stateB==1)) || 
557           ((fCalChamberStatus->GetStatus(det) == 3) && (stateA==1 || stateB==0)) ||
558           ((fCalChamberStatus->GetStatus(det) == 4) && (stateB==0 || stateB==1))  )
559         {
560           //printf(" Different half chamber status in DCS and DATA!!\n");
561           Double_t val[4] = {sm,lay,stac,1};
562           fHnSparseDebug->Fill(&val[0]); 
563           
564           // Fill MCM status map
565           for(Int_t ii = 0; ii < 8; ii++) { //ROB loop
566             for(Int_t i = 0; i < 18; i++) { //MCM loop
567               Double_t valss[6] = {sm,lay,stac,ii,i
568                                    ,calDCSFEEEOR->GetMCMGlobalState(ii,i)};
569               fHnSparseMCM->Fill(&valss[0]);
570               
571             } 
572           }
573         }
574     }
575     //---------//
576     //  Debug  //
577
578     //---------------------------------------
579     // Change the status according to DCS
580     //---------------------------------------
581     Int_t StatusData = fCalChamberStatus->GetStatus(det);
582     switch(StatusData) 
583       {
584       case 1: 
585         if(stateA==0 && stateB==0) fCalChamberStatus->SetStatus(det,2); // completely masked from DCS
586         if(stateA==1 && stateB==0) fCalChamberStatus->SetStatus(det,4); // Only B side masked from DCS
587         if(stateA==0 && stateB==1) fCalChamberStatus->SetStatus(det,3); // Only A side masked from DCS
588         if(stateA==1 && stateB==1) fCalChamberStatus->SetStatus(det,1);
589         break;
590       case 2: // completely masked from DATA
591         if(stateA==0 && stateB==0) fCalChamberStatus->SetStatus(det,2); // completely masked from DCS
592         break;
593       case 3: // Only A side masked from DATA
594         if(stateA==0 && stateB==0) fCalChamberStatus->SetStatus(det,2); // completely masked from DCS
595         if(stateA==1 && stateB==0) fCalChamberStatus->SetStatus(det,2); // Only B side masked from DCS
596         if(stateA==0 && stateB==1) fCalChamberStatus->SetStatus(det,3); // Only A side masked from DCS
597         if(stateA==1 && stateB==1) fCalChamberStatus->SetStatus(det,3);
598         break;
599       case 4: // Only B side masked from DATA
600         if(stateA==0 && stateB==0) fCalChamberStatus->SetStatus(det,2); // completely masked from DCS
601         if(stateA==1 && stateB==0) fCalChamberStatus->SetStatus(det,4); // Only B side masked from DCS
602         if(stateA==0 && stateB==1) fCalChamberStatus->SetStatus(det,2); // Only A side masked from DCS
603         if(stateA==1 && stateB==1) fCalChamberStatus->SetStatus(det,4);
604         break;
605       }
606     
607   }
608
609 }
610
611 //_____________________________________________________________________________________
612 void AliTRDCalibChamberStatus::Add(AliTRDCalibChamberStatus *calibChamberStatus) /*FOLD00*/
613 {
614     //
615     //  Add the THnSparseI of this calibChamberStatus
616     //
617
618   fCounterEventNotEmpty += calibChamberStatus->GetNumberEventNotEmpty();
619
620   THnSparseI *hnSparseI = calibChamberStatus->GetSparseI();
621   if(!hnSparseI) return;
622
623   if(!fHnSparseI) {
624     fHnSparseI = (THnSparseI *) hnSparseI->Clone();
625   }
626   else {
627     fHnSparseI->Add(hnSparseI);
628   }
629   
630
631 }
632 //_____________________________________________________________________
633 void AliTRDCalibChamberStatus::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
634 {
635     //
636     //  Write class to file
637     //
638
639     TString sDir(dir);
640     TString option;
641
642     if ( append )
643         option = "update";
644     else
645         option = "recreate";
646
647     TDirectory *backup = gDirectory;
648     TFile f(filename,option.Data());
649     f.cd();
650     if ( !sDir.IsNull() ){
651         f.mkdir(sDir.Data());
652         f.cd(sDir);
653     }
654     this->Write();
655     f.Close();
656
657     if ( backup ) backup->cd();
658 }
659 //_____________________________________________________________________________
660 TH2D* AliTRDCalibChamberStatus::PlotSparseI(Int_t sm,Int_t side) 
661 {
662   //
663   // Plot number of entries for supermodule sm 
664   // as a function of layer and stack
665   //
666
667   if(!fHnSparseI) return 0x0;
668   
669   fHnSparseI->GetAxis(0)->SetRange(sm+1,sm+1); 
670   fHnSparseI->GetAxis(3)->SetRange(side+1,side+1);  
671   TH2D *h2 = fHnSparseI->Projection(1,2);
672  
673
674   return h2;
675
676 }
677 //_____________________________________________________________________
678 TH2F *AliTRDCalibChamberStatus::MakeHisto2DSmPlEORStatus(AliTRDCalDCS *calDCS, Int_t sm, Int_t pl) /*FOLD00*/
679 {
680   //
681   //  Plot globale state of the HalfChamberMerger (HCM)
682   //
683   AliTRDfeeParam *paramfee = AliTRDfeeParam::Instance();
684
685   AliTRDgeometry *trdGeo = new AliTRDgeometry();
686   AliTRDpadPlane *padPlane0 = trdGeo->GetPadPlane(pl,0);        // layer,stack
687   Double_t row0    = padPlane0->GetRow0();
688   Double_t col0    = padPlane0->GetCol0();
689
690   char  name[1000];
691   snprintf(name,1000,"%s DCS status sm %d pl %d",GetTitle(),sm,pl);
692   TH2F * his = new TH2F( name, name, 88,-TMath::Abs(row0),TMath::Abs(row0)
693                                    ,148,-TMath::Abs(col0),TMath::Abs(col0));
694
695
696   // Where we begin
697   Int_t offsetsmpl = 30*sm+pl;
698   Int_t nstack = 5;
699   Int_t ncols = 144;
700
701   for (Int_t k = 0; k < nstack; k++){
702     Int_t det = offsetsmpl+k*6;
703     Int_t stac = AliTRDgeometry::GetStack(det);
704     AliTRDCalDCSFEE* calDCSFEEEOR = calDCS->GetCalDCSFEEObj(det);
705     if(!calDCSFEEEOR) { continue;}
706     for (Int_t icol=0; icol<ncols; icol++){
707       Int_t nrows = 16;
708       if(stac==2) nrows = 12;
709       for (Int_t irow=0; irow<nrows; irow++){
710         Int_t binz     = 0;
711         Int_t kb       = 5-1-k;
712         Int_t krow     = nrows-1-irow;
713         Int_t kcol     = ncols-1-icol;
714         if(kb > 2) binz = 16*(kb-1)+12+krow+1+2*(kb+1);
715         else binz = 16*kb+krow+1+2*(kb+1); 
716         Int_t biny = kcol+1+2;
717         // Take the value
718         Int_t mcm = paramfee->GetMCMfromPad(irow,icol);
719         Int_t rob = paramfee->GetROBfromPad(irow,icol);
720         if(mcm < 0) AliWarning("Problem with mcm number");
721         Int_t state = calDCSFEEEOR->GetMCMGlobalState(rob,TMath::Abs(mcm)); 
722         his->SetBinContent(binz,biny,state);
723       }
724     }
725     for(Int_t icol = 1; icol < 147; icol++){
726       for(Int_t l = 0; l < 2; l++){
727         Int_t binz     = 0;
728         Int_t kb       = 5-1-k;
729         if(kb > 2) binz = 16*(kb-1)+12+1+2*(kb+1)-(l+1);
730         else binz = 16*kb+1+2*(kb+1)-(l+1); 
731         his->SetBinContent(binz,icol,16.0);
732       }
733     }
734   }
735   
736   for(Int_t icol = 1; icol < 147; icol++){
737     his->SetBinContent(88,icol,16.0);
738     his->SetBinContent(87,icol,16.0);
739   }
740   for(Int_t irow = 1; irow < 89; irow++){
741     his->SetBinContent(irow,1,16.0);
742     his->SetBinContent(irow,2,16.0);
743     his->SetBinContent(irow,147,16.0);
744     his->SetBinContent(irow,148,16.0);
745   }
746
747   his->SetXTitle("z (cm)");
748   his->SetYTitle("y (cm)");
749   his->SetMaximum(12);
750   his->SetMinimum(0.0);
751   his->SetStats(0);
752
753   return his;
754
755 }
756 //_____________________________________________________________________________
757 TCanvas* AliTRDCalibChamberStatus::PlotHistos2DSmEORStatus(AliTRDCalDCS *calDCS, Int_t sm, const Char_t *name)
758 {
759   //
760   // Make 2D graph
761   //
762
763   gStyle->SetPalette(1);
764   fC1 = new TCanvas(name,name,50,50,600,800);
765   fC1->Divide(3,2);
766   fC1->cd(1);
767   MakeHisto2DSmPlEORStatus(calDCS,sm,0)->Draw("colz");
768   fC1->cd(2);
769   MakeHisto2DSmPlEORStatus(calDCS,sm,1)->Draw("colz");
770   fC1->cd(3);
771   MakeHisto2DSmPlEORStatus(calDCS,sm,2)->Draw("colz");
772   fC1->cd(4);
773   MakeHisto2DSmPlEORStatus(calDCS,sm,3)->Draw("colz");
774   fC1->cd(5);
775   MakeHisto2DSmPlEORStatus(calDCS,sm,4)->Draw("colz");
776   fC1->cd(6);
777   MakeHisto2DSmPlEORStatus(calDCS,sm,5)->Draw("colz");
778
779   return fC1;
780
781 }
782
783
784