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