]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDCalibViewer.cxx
Coverity stuff once more ...
[u/mrichter/AliRoot.git] / TRD / AliTRDCalibViewer.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 /* $Id: AliTRDCalibViewer.cxx 40390 2010-04-14 09:43:23Z cblume $ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 //  Class which implements AliBaseCalibViewer for the TRD                    //
21 //  used for the calibration monitor                                         //
22 //                                                                           //
23 //  Authors:     Marian Ivanov (Marian.Ivanov@cern.ch)                       //
24 //               Jens Wiechula (Jens.Wiechula@cern.ch)                       //
25 //               Ionut Arsene  (iarsene@cern.ch)                             //
26 //                                                                           //
27 ///////////////////////////////////////////////////////////////////////////////
28
29 #include <iostream>
30 #include <fstream>
31 #include <TString.h>
32 #include <TRandom.h>
33 #include <TLegend.h>
34 #include <TLine.h>
35 #include <TCanvas.h>
36 #include <TROOT.h>
37 #include <TStyle.h>
38 #include <TH1.h> 
39 #include <TH1F.h>
40 #include <TMath.h>
41 #include <TVectorD.h>
42 #include <THashTable.h>
43 #include <TObjString.h>
44 #include <TTimeStamp.h>
45 #include <TObjString.h>
46 #include <TTreeStream.h>
47 #include <TFile.h>
48 #include <TKey.h>
49 #include <TGraph.h>
50 #include <TDirectory.h>
51 #include <TFriendElement.h>
52 #include <TGrid.h>
53 #include <TGeoManager.h>
54 #include "AliTRDCalDet.h"
55 #include "AliTRDCalPad.h"
56 #include "AliTRDCalROC.h"
57 #include "AliTRDCalChamberStatus.h"
58 #include "AliTRDCalSingleChamberStatus.h"
59 #include "AliTRDCalPadStatus.h"
60 #include "AliTRDCalDCS.h"
61 #include "AliTRDCalDCSFEE.h"
62 #include "AliTRDcalibDB.h"
63 #include "AliCDBManager.h"
64 #include "AliCDBStorage.h"
65 #include "AliCDBEntry.h"
66 #include "AliGRPObject.h"
67 #include "AliTRDalignment.h"
68 #include "AliTRDgeometry.h"
69 #include "AliTRDpadPlane.h"
70
71 #include "AliTRDCalibViewer.h"
72
73
74 using namespace std;
75
76 ClassImp(AliTRDCalibViewer)
77
78 //_____________________________________________________________________________
79 AliTRDCalibViewer::AliTRDCalibViewer()
80                   :AliBaseCalibViewer()
81 {
82   //
83   // Default constructor (just call base class constructor)
84   //
85 }
86
87 //_____________________________________________________________________________
88 AliTRDCalibViewer::AliTRDCalibViewer(const AliTRDCalibViewer &c)
89                   :AliBaseCalibViewer(c)
90 {
91   //
92   // copy constructor (just call base class copy constructor)
93   // 
94 }
95
96 //_____________________________________________________________________________
97 AliTRDCalibViewer::AliTRDCalibViewer(TTree* tree)
98                   :AliBaseCalibViewer(tree)
99 {
100   //
101   // Constructor (just call the corresponding base constructor)
102   //
103 }
104
105 //_____________________________________________________________________________
106 AliTRDCalibViewer::AliTRDCalibViewer(const char* fileName, const char* treeName)
107                   :AliBaseCalibViewer(fileName, treeName)                   
108 {
109    //
110    // Constructor (just call the corresponding base constructor)
111    //
112 }
113
114 //_____________________________________________________________________________
115 AliTRDCalibViewer& AliTRDCalibViewer::operator=(const AliTRDCalibViewer& param) {
116   //
117   // assignment operator
118   //
119   fTree = param.fTree;
120   fTreeMustBeDeleted = param.fTreeMustBeDeleted;
121   fListOfObjectsToBeDeleted = param.fListOfObjectsToBeDeleted;
122   fAbbreviation = param.fAbbreviation;
123   fAppendString = param.fAppendString;
124   return(*this);
125 }
126
127 //_____________________________________________________________________________
128 AliTRDCalibViewer::~AliTRDCalibViewer()
129 {
130    //
131    // AliTRDCalibViewer destructor
132    // do nothing, the base class destructor will do the job
133 }
134
135 /*
136 //_____________________________________________________________________________
137 void AliTRDCalibViewer::GetTimeInfoOCDB(const Char_t* runList, const Char_t* outFile,
138                                         Int_t firstRun, Int_t lastRun, UInt_t infoFlags,
139                                         const Char_t* ocdbStorage) {
140 //
141 //  Get time information from OCDB by calling the DumpOCDBtoTree.C macro
142 //
143   DumpOCDBtoTree(runList, outFile, firstRun, lastRun,
144                  TESTBIT(infoFlags,1), TESTBIT(infoFlags,2), 
145                  TESTBIT(infoFlags,3), TESTBIT(infoFlags,4),
146                  ocdbStorage);
147 }
148 */
149
150 //_____________________________________________________________________________
151 const char* AliTRDCalibViewer::AddAbbreviations(char* c, Bool_t printDrawCommand){ 
152    // Replace all "<variable>" with "<variable><fAbbreviation>" (Adds forgotten "~")
153    // but take care on the statistical information, like "CEQmean_Mean"
154    // and also take care on correct given variables, like "CEQmean~"
155    // 
156    // For each variable out of "listOfVariables":
157    // - 'Save' correct items:
158    //   - form <replaceString>, take <variable>'s first char, add <removeString>, add rest of <variable>, e.g. "C!#EQmean" (<removeString> = "!#")
159    //   - For each statistical information in "listOfNormalizationVariables":
160    //     - ReplaceAll <variable><statistical_Information> with <replaceString><statistical_Information>
161    //   - ReplaceAll <variable><abbreviation> with <replaceString><abbreviation>, e.g. "CEQmean~" -> "C!#EQmean~"
162    //   - ReplaceAll <variable><appendStr> with <replaceString><appendStr>, e.g. "CEQmean.fElements" -> "C!#EQmean.fElements"
163    //
164    // - Do actual replacing:
165    //   - ReplaceAll <variable> with <variable><fAbbreviation>, e.g. "CEQmean" -> "CEQmean~"
166    //
167    // - Undo saving:
168    //   - For each statistical information in "listOfNormalizationVariables":
169    //     - ReplaceAll <replaceString><statistical_Information> with <variable><statistical_Information> 
170    //   - ReplaceAll <replaceString><abbreviation> with <variable><abbreviation>, e.g. "C!#EQmean~" -> "CEQmean~"
171    //   - ReplaceAll <replaceString><appendStr> with <variable><appendStr>, e.g. "C!#EQmean.fElements" -> "CEQmean.fElements"
172    // 
173    // Now all the missing "~" should be added.
174    
175    TString str(c);
176    TString removeString = "!#";  // very unprobable combination of chars
177    TString replaceString = "";
178    TString searchString = "";
179    TString normString = "";
180    TObjArray *listOfVariables = GetListOfVariables();
181    // variables used for mapping the pads, mcms, ...
182    listOfVariables->Add(new TObjString("SuperModule"));
183    listOfVariables->Add(new TObjString("Layer"));
184    listOfVariables->Add(new TObjString("Stack"));
185    listOfVariables->Add(new TObjString("Channel"));
186    listOfVariables->Add(new TObjString("Row"));
187    listOfVariables->Add(new TObjString("Column"));
188    listOfVariables->Add(new TObjString("Chamber"));
189    listOfVariables->Add(new TObjString("PadSuperRow"));
190    listOfVariables->Add(new TObjString("PadSuperColumn"));
191    listOfVariables->Add(new TObjString("MCMSuperRow"));
192    listOfVariables->Add(new TObjString("MCMSuperColumn"));
193    listOfVariables->Add(new TObjString("ROB"));
194    listOfVariables->Add(new TObjString("MCM"));
195    TObjArray *listOfNormalizationVariables = GetListOfNormalizationVariables();
196    Int_t nVariables = listOfVariables->GetEntriesFast();
197    Int_t nNorm = listOfNormalizationVariables->GetEntriesFast();
198    
199    Int_t *varLengths = new Int_t[nVariables];
200    for (Int_t i = 0; i < nVariables; i++) {
201       varLengths[i] = ((TObjString*)listOfVariables->At(i))->String().Length();
202    }
203    Int_t *normLengths = new Int_t[nNorm];
204    for (Int_t i = 0; i < nNorm; i++) {
205       normLengths[i] = ((TObjString*)listOfNormalizationVariables->At(i))->String().Length();
206    }
207    Int_t *varSort = new Int_t[nVariables];
208    TMath::Sort(nVariables, varLengths, varSort, kTRUE);
209    Int_t *normSort = new Int_t[nNorm];
210    TMath::Sort(nNorm, normLengths, normSort, kTRUE);
211       
212    for (Int_t ivar = 0; ivar < nVariables; ivar++) {
213       // ***** save correct tokens *****
214       // first get the next variable:
215       searchString = ((TObjString*)listOfVariables->At(varSort[ivar]))->String();
216       // form replaceString:
217       replaceString = "";
218       for (Int_t i = 0; i < searchString.Length(); i++) {
219          replaceString.Append(searchString[i]);
220          if (i == 0) replaceString.Append(removeString);
221       }
222       // go through normalization:
223       for (Int_t inorm = 0; inorm < nNorm; inorm++) {
224         normString = ((TObjString*)listOfNormalizationVariables->At(normSort[inorm]))->String();
225         str.ReplaceAll(searchString + normString, replaceString + normString);
226         // like: str.ReplaceAll("CEQmean_Mean", "C!EQmean_Mean");
227       }
228       str.ReplaceAll(searchString + fAbbreviation, replaceString + fAbbreviation);
229       // like: str.ReplaceAll("CEQmean~", "C!EQmean~");
230       str.ReplaceAll(searchString + fAppendString,    replaceString + fAppendString);
231       // like: str.ReplaceAll("CEQmean.fElements", "C!EQmean.fElements");
232       
233       // ***** add missing extensions *****
234       str.ReplaceAll(searchString, replaceString + fAbbreviation);
235       // like: str.ReplaceAll("CEQmean", "C!EQmean~");
236    }
237    
238    // ***** undo saving *****
239    str.ReplaceAll(removeString, "");
240   
241    if (printDrawCommand) std::cout << "The string looks now like: " << str.Data() << std::endl;
242    delete [] varSort;
243    delete [] normSort;
244    return str.Data();
245 }
246
247 //_____________________________________________________________________________
248 TObjArray* AliTRDCalibViewer::GetListOfVariables(Bool_t printList) {
249   //
250   // scan the tree  - produces a list of available variables in the tree
251   // printList: print the list to the screen, after the scan is done
252   //
253   TObjArray* arr = new TObjArray();
254   TObjString* str = 0;
255   if (!fTree) {
256     return 0;
257   }
258   Int_t nentries = fTree->GetListOfBranches()->GetEntries();
259   for (Int_t i = 0; i < nentries; i++) {
260     str = new TObjString(fTree->GetListOfBranches()->At(i)->GetName());
261     str->String().ReplaceAll("_Mean", "");
262     str->String().ReplaceAll("_RMS", "");
263     str->String().ReplaceAll("_Median", "");
264     str->String().ReplaceAll(".", "");
265     str->String().ReplaceAll("_Run", "");
266     str->String().ReplaceAll("_SuperModule", "");
267     str->String().ReplaceAll("_Chamber", "");
268     // add all the variables in the tree to a list
269     // exception make variables which are used for mapping, specified in AddAbbreviations()
270     // These two functions should be kept synchronized with respect to the mapping variables 
271     if (!arr->FindObject(str) && 
272         !(str->String() == "run" || 
273           str->String() == "SuperModule" || str->String() == "Layer" || str->String() == "Stack" || 
274           str->String() == "Chamber" ||
275           str->String() == "Channel" || str->String() == "Row" || str->String() == "Column" || 
276           str->String() == "PadSuperRow" || str->String() == "PadSuperColumn" || str->String() == "MCMSuperRow"
277           || str->String() == "MCMSuperColumn" || str->String() == "MCM" || str->String() == "ROB")) {
278       arr->Add(str);
279     }
280   }
281   
282   // loop over all friends (if there are some) and add them to the list
283   if (fTree->GetListOfFriends()) {
284     for (Int_t ifriend = 0; ifriend < fTree->GetListOfFriends()->GetEntries(); ifriend++){
285       // printf("iterating through friendlist, currently at %i\n", ifriend);
286       // printf("working with %s\n", fTree->GetListOfFriends()->At(ifriend)->ClassName());
287       if (TString(fTree->GetListOfFriends()->At(ifriend)->ClassName()) != "TFriendElement") continue; // no friendElement found
288       TFriendElement *friendElement = (TFriendElement*)fTree->GetListOfFriends()->At(ifriend);
289       if (friendElement->GetTree() == 0) continue; // no tree found in friendElement
290       // printf("friend found \n");
291       for (Int_t i = 0; i < friendElement->GetTree()->GetListOfBranches()->GetEntries(); i++) {
292         // printf("iterating through friendelement entries, currently at %i\n", i);
293         str = new TObjString(friendElement->GetTree()->GetListOfBranches()->At(i)->GetName());
294         str->String().ReplaceAll("_Mean", "");
295         str->String().ReplaceAll("_RMS", "");
296         str->String().ReplaceAll("_Median", "");
297         str->String().ReplaceAll(".", "");
298         str->String().ReplaceAll("_Run", "");
299         str->String().ReplaceAll("_SuperModule", "");
300         str->String().ReplaceAll("_Chamber", "");
301         if (!(str->String() == "run" || 
302               str->String() == "SuperModule" || str->String() == "Layer" || str->String() == "Stack" || 
303               str->String() == "Chamber" ||
304               str->String() == "Channel" || str->String() == "Row" || str->String() == "Column" || 
305               str->String() == "PadSuperRow" || str->String() == "PadSuperColumn" || str->String() == "MCMSuperRow"
306               || str->String() == "MCMSuperColumn" || str->String() == "MCM" || str->String() == "ROB")){
307           // insert "<friendName>." at the beginning: (<friendName> is per default "R")
308           str->String().Insert(0, ".");
309           str->String().Insert(0, friendElement->GetName());
310           if (!arr->FindObject(str)) {
311             arr->Add(str);
312           }
313           // printf("added string %s \n", str->String().Data());
314         }
315       }
316     }
317   } // if (fTree->GetListOfFriends())
318   
319   arr->Sort();
320     
321   if (printList) {
322     TIterator* iter = arr->MakeIterator();
323     iter->Reset();
324     TObjString* currentStr = 0;
325     while ( (currentStr = (TObjString*)(iter->Next())) ) {
326       std::cout << currentStr->GetString().Data() << std::endl;
327     }
328     delete iter;
329   }
330   return arr;
331 }
332
333 TObjArray* AliTRDCalibViewer::GetListOfNormalizationVariables(Bool_t printList) const{
334   //
335   // produces a list of available variables for normalization in the tree
336   // printList: print the list to the screen, after the scan is done
337   //
338    TObjArray* arr = new TObjArray();
339    arr->Add(new TObjString("_Mean_Run"));
340    arr->Add(new TObjString("_Mean_SuperModule"));
341    arr->Add(new TObjString("_Mean_Chamber"));
342    arr->Add(new TObjString("_Median_Run"));
343    arr->Add(new TObjString("_Median_SuperModule"));
344    arr->Add(new TObjString("_Median_Chamber"));
345    
346    if (printList) {
347      TIterator* iter = arr->MakeIterator();
348      iter->Reset();
349      TObjString* currentStr = 0;
350      while ((currentStr = (TObjString*)(iter->Next()))) {
351        std::cout << currentStr->GetString().Data() << std::endl;
352      }
353      delete iter;
354    }
355    return arr;
356 }
357
358 void AliTRDCalibViewer::GetLayerSectorStack(TString trdString, Int_t& layerNo, Int_t& sectorNo, Int_t& stackNo) const {
359   // Get the layer, sector and stack numbers out of a string
360   // encoded with the following format:
361   // Layer%dSector%dStack%d
362
363   sscanf(trdString.Data(), "Layer%1dSector%02dStack%1d", &layerNo, &sectorNo, &stackNo);
364
365   return;
366 }
367
368 //_____________________________________________________________________________
369 Int_t AliTRDCalibViewer::EasyDraw(const char* drawCommand, const char* sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
370   //
371   // easy drawing of data, use '~' for abbreviation of '.fElements'
372   // example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
373   // sector: sector-number - only the specified sector will be drwawn
374   //         'A'/'C' or 'a'/'c' - side A/C will be drawn
375   //         'ALL' - whole TPC will be drawn, projected on one side
376   // cuts: specifies cuts
377   // drawOptions: draw options like 'same'
378   // writeDrawCommand: write the command, that is passed to TTree::Draw
379   //
380
381   TString drawStr(drawCommand);
382
383   TString sectorStr(sector);
384   Int_t layerNo = -1; 
385   Int_t sectorNo = -1; 
386   Int_t stackNo = -1;
387   GetLayerSectorStack(sectorStr, layerNo, sectorNo, stackNo);
388   if(layerNo==-1) {
389      Warning("EasyDraw", "The sector string must always contain the Layer number!");
390      return -1;
391    }
392    if(layerNo<0 || layerNo>5) {
393      Warning("EasyDraw", "The Layer number must be in the range [0,5] !");
394      return -1;
395    }
396    if(sectorNo!=-1 && (sectorNo<0 || sectorNo>17)) {
397      Warning("EasyDraw", "The SuperModule number must be in the range [0,17] !");
398      return -1;
399    }
400    if(stackNo!=-1 && (stackNo<0 || stackNo>4)) {
401      Warning("EasyDraw", "The Stack number must be in the range [0,4] !");
402      return -1;
403    }
404
405    TString cutStr("");
406
407    Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
408    if (dangerousToDraw) {
409       Warning("EasyDraw", "The draw string must not contain ':' or '>>'. Using only first variable for drawing!");
410       drawStr.Resize(drawStr.First(":"));
411    }
412
413    TString drawOptionsStr("");
414    TRandom rnd(0);
415    Int_t rndNumber = rnd.Integer(10000);
416
417    if (drawOptions && strcmp(drawOptions, "") != 0)
418       drawOptionsStr += drawOptions;
419    else
420       drawOptionsStr += "profcolz";
421
422    const Int_t gkNRows[ 5] = {16, 16, 12, 16, 16};  // number of pad rows in the chambers from each of the 5 stacks
423
424    // draw calibration stuff
425    if(drawStr.Contains("Status") || drawStr.Contains("Gain") || drawStr.Contains("Noise") ||
426       drawStr.Contains("Vdrift") || drawStr.Contains("T0") || 
427       drawStr.Contains("gain") || drawStr.Contains("chiSquare")) {
428      if(sectorNo==-1 && stackNo==-1) {    // plot the entire layer
429        drawStr += Form(":PadSuperColumn%s:PadSuperRow%s>>prof", fAppendString.Data(), fAppendString.Data());
430        drawStr += rndNumber;
431        drawStr += "(76,-0.5,75.5,2592,-0.5,2591.5)";
432        cutStr += Form("Layer==%d", layerNo);
433      }
434      else if(sectorNo!=-1 && stackNo==-1) {     // plot a sector from a layer
435        drawStr += Form(":Column%s:PadSuperRow%s>>prof", fAppendString.Data(), fAppendString.Data());
436        drawStr += rndNumber;
437        drawStr += "(76,-0.5,75.5,144,-0.5,143.5)";
438        cutStr += Form("Layer==%d && SuperModule==%d", layerNo, sectorNo);
439      }
440      else if(sectorNo==-1 && stackNo!=-1) {       // plot a stack from a layer
441        drawStr += Form(":PadSuperColumn%s:Row%s>>prof", fAppendString.Data(), fAppendString.Data());
442        drawStr += rndNumber;
443        drawStr += Form("(%d,-0.5,%d-0.5,2592,-0.5,2591.5)", gkNRows[stackNo], gkNRows[stackNo]);
444        cutStr += Form("Layer==%d && Stack==%d", layerNo, stackNo);
445      }
446      else {                // the layer, sector and stack are determined -> so plot a chamber
447        drawStr += Form(":Column%s:Row%s>>prof", fAppendString.Data(), fAppendString.Data());
448        drawStr += rndNumber;
449        drawStr += Form("(%d,-0.5,%d-0.5,144,-0.5,143.5)", gkNRows[stackNo], gkNRows[stackNo]);
450        cutStr += Form("Layer==%d && SuperModule==%d && Stack==%d", layerNo, sectorNo, stackNo);
451      }
452    }
453    // draw FEE stuff
454    else if(drawStr.Contains("SORandEOR") || 
455            drawStr.Contains("gsmSOR") || drawStr.Contains("gsmDelta") ||
456            drawStr.Contains("nimSOR") || drawStr.Contains("nimDelta") ||
457            drawStr.Contains("nevSOR") || drawStr.Contains("nevDelta") ||
458            drawStr.Contains("nptSOR") || drawStr.Contains("nptDelta")) {
459      if(sectorNo==-1 && stackNo==-1) {    // plot the entire layer
460        drawStr += Form(":MCMSuperColumn%s:MCMSuperRow%s>>prof", fAppendString.Data(), fAppendString.Data());
461        drawStr += rndNumber;
462        drawStr += "(76,-0.5,75.5,144,-0.5,143.5)";
463        cutStr += Form("Layer==%d", layerNo);
464      }
465      else if(sectorNo!=-1 && stackNo==-1) {     // plot a sector from a layer
466        drawStr += Form(":MCMSuperColumn%s:MCMSuperRow%s>>prof", fAppendString.Data(), fAppendString.Data());
467        drawStr += rndNumber;
468        drawStr += "(76,-0.5,75.5,144,-0.5,143.5)";
469        cutStr += Form("Layer==%d && SuperModule==%d", layerNo, sectorNo);
470      }
471      else if(sectorNo==-1 && stackNo!=-1) {       // plot a stack from a layer
472        drawStr += Form(":MCMSuperColumn%s:MCMSuperRow%s>>prof", fAppendString.Data(), fAppendString.Data());
473        drawStr += rndNumber;
474        //       drawStr += Form("(%d,-0.5,%d-0.5,2592,-0.5,2591.5)", gkNRows[stackNo], gkNRows[stackNo]);
475        drawStr += "(76,-0.5,75.5,144,-0.5,143.5)";
476        cutStr += Form("Layer==%d && Stack==%d", layerNo, stackNo);
477      }
478      else {                // the layer, sector and stack are determined -> so plot a chamber
479        drawStr += Form(":ROB%s:MCM%s>>prof", fAppendString.Data(), fAppendString.Data());
480        drawStr += rndNumber;
481        drawStr += Form("(16,-0.5,15.5,%d,-0.5,%d-0.5)", gkNRows[stackNo]/2, gkNRows[stackNo]/2);
482        cutStr += Form("Layer==%d && SuperModule==%d && Stack==%d", layerNo, sectorNo, stackNo);
483      }
484    }
485    // draw alignment stuff
486    else if(drawStr.Contains("Align")) {
487      if(sectorNo==-1 && stackNo==-1) {  // plot the entire layer
488        drawStr += ":SuperModule:Stack>>prof";
489        drawStr += rndNumber;
490        drawStr += "(5,-0.5,4.5,18,-0.5,17.5)";
491        cutStr += Form("Layer==%d", layerNo);
492      }
493      else if(sectorNo!=-1 && stackNo==-1) {     // plot a sector from a layer
494        drawStr += ":SuperModule:Stack>>prof";
495        drawStr += rndNumber;
496        drawStr += Form("(5,-0.5,4.5,1,%f,%f)", sectorNo-0.5, sectorNo+0.5);
497        cutStr += Form("Layer==%d && SuperModule==%d", layerNo, sectorNo);
498      }
499      else if(sectorNo==-1 && stackNo!=-1) {       // plot a stack from a layer
500        drawStr += ":SuperModule:Stack>>prof";
501        drawStr += rndNumber;
502        drawStr += Form("(1,%f,%f,18,-0.5,17.5)", stackNo-0.5, stackNo+0.5);
503        cutStr += Form("Layer==%d && Stack==%d", layerNo, stackNo);
504      }
505      else {                // the layer, sector and stack are determined -> so plot a chamber
506        drawStr += ":SuperModule:Stack>>prof";
507        drawStr += rndNumber;
508        drawStr += Form("(1,%f,%f,1,%f,%f)", stackNo-0.5, stackNo+0.5, sectorNo-0.5, sectorNo+0.5);
509        cutStr += Form("Layer==%d && SuperModule==%d && Stack==%d", layerNo, sectorNo, stackNo);
510      }
511    }
512
513
514    if (cuts && cuts[0] != 0) {
515       if (cutStr.Length() != 0) cutStr += "&& ";
516       cutStr += "(";
517       cutStr += cuts;
518       cutStr += ")";
519    }
520    drawStr.ReplaceAll(fAbbreviation, fAppendString);
521    cutStr.ReplaceAll(fAbbreviation, fAppendString);
522    if (writeDrawCommand) std::cout << "fTree->Draw(\"" << drawStr << "\", \"" <<  cutStr << "\", \"" << drawOptionsStr << "\");" << std::endl;
523    Int_t returnValue = fTree->Draw(drawStr.Data(), cutStr.Data(), drawOptionsStr.Data());
524    TString profName("prof");
525    profName += rndNumber;
526    TObject *obj = gDirectory->Get(profName.Data());
527    // set the names of the axes 
528    TH1 *histObj = (TH1*)obj;
529    if(drawStr.Contains("Status") || drawStr.Contains("Gain") || drawStr.Contains("Noise") ||
530       drawStr.Contains("Vdrift") || drawStr.Contains("T0") ||
531       drawStr.Contains("gain") || drawStr.Contains("chiSquare")) {
532      histObj->GetXaxis()->SetTitle("Row");
533      histObj->GetYaxis()->SetTitle("Column");
534    }
535    else if(drawStr.Contains("SORandEOR") || 
536            drawStr.Contains("gsmSOR") || drawStr.Contains("gsmDelta") ||
537            drawStr.Contains("nimSOR") || drawStr.Contains("nimDelta") ||
538            drawStr.Contains("nevSOR") || drawStr.Contains("nevDelta") ||
539            drawStr.Contains("nptSOR") || drawStr.Contains("nptDelta")) {
540      histObj->GetXaxis()->SetTitle("MCM Row");
541      histObj->GetYaxis()->SetTitle("MCM Column");
542    }
543    else if(drawStr.Contains("Align")) {
544      histObj->GetXaxis()->SetTitle("Stack");
545      histObj->GetYaxis()->SetTitle("Sector");
546    }
547
548    if (obj && obj->InheritsFrom("TH1")) FormatHistoLabels((TH1*)obj);
549    return returnValue;
550 }
551
552 //_____________________________________________________________________________
553 Int_t AliTRDCalibViewer::EasyDraw1D(const char* drawCommand, const char* sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
554   //
555   // easy drawing of data, use '~' for abbreviation of '.fElements'
556   // example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
557   // sector: sector-number - the specified sector will be drwawn
558   //         'A'/'C' or 'a'/'c' - side A/C will be drawn
559   //         'ALL' - whole TPC will be drawn, projected on one side
560   // cuts: specifies cuts
561   // drawOptions: draw options like 'same'
562   // writeDrawCommand: write the command, that is passed to TTree::Draw
563   //
564
565    TString drawStr(drawCommand);
566
567    TString sectorStr(sector);
568    Int_t layerNo = -1; 
569    Int_t sectorNo = -1; 
570    Int_t stackNo = -1;
571    GetLayerSectorStack(sectorStr, layerNo, sectorNo, stackNo);
572    if(layerNo==-1) {
573      Warning("EasyDraw", "The sector string must always contain the Layer number!");
574      return -1;
575    }
576    if(layerNo<0 || layerNo>5) {
577      Warning("EasyDraw", "The Layer number must be in the range [0,5] !");
578      return -1;
579    }
580    if(sectorNo!=-1 && (sectorNo<0 || sectorNo>17)) {
581      Warning("EasyDraw", "The Sector number must be in the range [0,17] !");
582      return -1;
583    }
584    if(stackNo!=-1 && (stackNo<0 || stackNo>4)) {
585      Warning("EasyDraw", "The Stack number must be in the range [0,4] !");
586      return -1;
587    }
588
589    TString drawOptionsStr(drawOptions);
590    TString cutStr("");
591
592    if(sectorNo==-1 && stackNo==-1)     // plot the entire layer
593      cutStr += Form("Layer==%d", layerNo);
594    else if(sectorNo!=-1 && stackNo==-1)      // plot a sector from a layer
595      cutStr += Form("Layer==%d && SuperModule==%d", layerNo, sectorNo);
596    else if(sectorNo==-1 && stackNo!=-1)        // plot a stack from a layer
597      cutStr += Form("Layer==%d && Stack==%d", layerNo, stackNo);
598    else                 // the layer, sector and stack are determined -> so plot a chamber
599      cutStr += Form("Layer==%d && SuperModule==%d && Stack==%d", layerNo, sectorNo, stackNo);
600    
601    if(cuts && cuts[0] != 0) {
602       if (cutStr.Length() != 0) cutStr += "&& ";
603       cutStr += "(";
604       cutStr += cuts;
605       cutStr += ")";
606    }
607
608    drawStr.ReplaceAll(fAbbreviation, fAppendString);
609    cutStr.ReplaceAll(fAbbreviation, fAppendString);
610    if (writeDrawCommand) std::cout << "fTree->Draw(\"" << drawStr << "\", \"" <<  cutStr << "\", \"" << drawOptionsStr << "\");" << std::endl;
611    Int_t returnValue = fTree->Draw(drawStr.Data(), cutStr.Data(), drawOptionsStr.Data());
612    if (returnValue == -1) return -1;
613    
614    TObject *obj = (gPad) ? gPad->GetPrimitive("htemp") : 0; 
615    if (!obj) obj = (TH1F*)gDirectory->Get("htemp");
616    if (!obj) obj = gPad->GetPrimitive("tempHist");
617    if (!obj) obj = (TH1F*)gDirectory->Get("tempHist");
618    if (!obj) obj = gPad->GetPrimitive("Graph");
619    if (!obj) obj = (TH1F*)gDirectory->Get("Graph");
620    if (obj && obj->InheritsFrom("TH1")) FormatHistoLabels((TH1*)obj);
621    return returnValue;
622 }
623
624 //_____________________________________________________________________________
625 Int_t AliTRDCalibViewer::EasyDraw(const char* drawCommand, Int_t chamber, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
626   //
627   // easy drawing of data, use '~' for abbreviation of '.fElements'
628   // example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0")
629   // sector: sector-number - only the specified sector will be drwawn
630   // cuts: specifies cuts
631   // drawOptions: draw options like 'same'
632   // writeDrawCommand: write the command, that is passed to TTree::Draw
633   //
634   if(chamber >= 0 && chamber < 540) {
635     Int_t superModuleNo = chamber/30;
636     Int_t stackNo = (chamber%30)/6;
637     Int_t layerNo = (chamber%30)%6;
638     char sectorChr[22];
639     snprintf(sectorChr,22, "Layer%iSector%iStack%i", layerNo, superModuleNo, stackNo);
640     return EasyDraw(drawCommand, sectorChr, cuts, drawOptions, writeDrawCommand);
641   }
642   Error("EasyDraw","The TRD contains only chamber from 0 to 539");
643   return -1;
644 }
645
646 //_____________________________________________________________________________
647 Int_t AliTRDCalibViewer::EasyDraw1D(const char* drawCommand, Int_t chamber, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
648   //
649   // easy drawing of data, use '~' for abbreviation of '.fElements'
650   // example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0")
651   // sector: sector-number - the specified sector will be drwawn
652   // cuts: specifies cuts
653   // drawOptions: draw options like 'same'
654   // writeDrawCommand: write the command, that is passed to TTree::Draw
655   //
656
657   if (chamber >= 0 && chamber < 539) {
658     Int_t superModuleNo = chamber/30;
659     Int_t stackNo = (chamber%30)/6;
660     Int_t layerNo = (chamber%30)%6;
661     char sectorChr[22];
662     snprintf(sectorChr,22, "Layer%iSector%iStack%i", layerNo, superModuleNo, stackNo);
663     return EasyDraw1D(drawCommand, sectorChr, cuts, drawOptions, writeDrawCommand);
664   }
665   Error("EasyDraw1D","The TRD contains only chambers from 0 to 539");
666   return -1;
667 }
668
669 //_____________________________________________________________________________
670 Bool_t AliTRDCalibViewer::DumpOCDBtoTreeDetails(const Char_t* runListFilename,
671                                                 const Char_t* outFilename,
672                                                 Int_t firstRun, Int_t lastRun,
673                                                 const Char_t* storage,
674                                                 Int_t version,
675                                                 Int_t subVersion,
676                                                 Bool_t getCalibs,
677                                                 Bool_t getDCS,
678                                                 Bool_t getAlign) {
679   //
680   // Retrieve TRD OCDB information for a given run list/range
681   //
682
683   if(runListFilename[0]!='\0' && firstRun==-1 && lastRun==-1) {
684     cout << "AliTRDCalibViewer::DumpOCDBtoTreeDetails(): You must provide at least a run range or an ascii filename with run numbers" 
685          << endl;
686     return kFALSE;
687   }
688   // initialize the OCDB manager
689   TString storageString = storage;
690   if(storageString.Contains("alien://")) {
691     TGrid::Connect("alien://");
692   }
693   AliCDBManager *manager = AliCDBManager::Instance();
694   if(storage[0]!='\0') {
695     manager->SetDefaultStorage(storage);
696   }
697   else {
698     if(!manager->IsDefaultStorageSet()) {
699       cout << "AliTRDCalibViewer::DumpOCDBtoTreeDetails(): Default OCDB storage not set!!" << endl;
700       return kFALSE;
701     }
702   }
703   manager->SetRun(1);
704
705   // open the ascii file
706   ifstream in;
707   if(runListFilename[0]!='\0')
708     in.open(runListFilename);
709
710   // initialize the tree streamer
711   if(outFilename[0]=='\0') outFilename = "trdDetails.root";
712   TString calibFilename = outFilename;
713   
714   TTreeSRedirector *treeStreamer = new TTreeSRedirector(calibFilename.Data());
715
716   Int_t currRun;
717   if(runListFilename[0]=='\0' && firstRun!=-1 && lastRun!=-1)
718     currRun = firstRun;
719
720   TVectorD runs;
721
722   // loop over runs
723   while(1) {
724     if(runListFilename[0]!='\0') {
725       if(!(in>>currRun)) continue;
726       if(currRun < (firstRun==-1 ? 0 : firstRun) ||
727          currRun > (lastRun==-1 ? 999999999 : lastRun))
728         continue;
729     }
730     else {
731       if(currRun>lastRun) break;
732     }
733     cout << "run = " << currRun << endl;
734     manager->SetRun(currRun);
735
736     // Get GRP data. If there is no proper GRP object for this run than
737     // this run is aborted
738     AliCDBEntry *entry = manager->Get("GRP/GRP/Data");
739     AliGRPObject* grpObject = 0;
740     if(entry) {
741       entry->SetOwner(kFALSE);
742       grpObject = dynamic_cast<AliGRPObject*>(entry->GetObject());
743     }
744     else {
745       currRun++;
746       //      continue;
747       //      return kFALSE;
748     }
749     if(!grpObject)
750       cout << "No GRP info available for this run " << endl;
751
752     time_t startTimeGRP = 0;
753     TObjString runType("");
754     if(grpObject) {
755       startTimeGRP = grpObject->GetTimeStart();
756       TTimeStamp start(grpObject->GetTimeStart());
757       TTimeStamp end(grpObject->GetTimeEnd());
758       cout << "Start time: " << start.GetDate()/10000 << "/" 
759            << (start.GetDate()/100)-(start.GetDate()/10000)*100 << "/" 
760            << start.GetDate()%100 << "   "
761            << start.GetTime()/10000 << ":"
762            << (start.GetTime()/100)-(start.GetTime()/10000)*100 << ":" 
763            << start.GetTime()%100 << endl;
764       cout << "End time: " << end.GetDate()/10000 << "/" 
765            << (end.GetDate()/100)-(end.GetDate()/10000)*100 << "/" 
766            << end.GetDate()%100 << "   "
767            << end.GetTime()/10000 << ":"
768            << (end.GetTime()/100)-(end.GetTime()/10000)*100 << ":"
769            << end.GetTime()%100 << endl;
770       cout << "Run type = " << grpObject->GetRunType().Data() << endl;
771       runType = grpObject->GetRunType().Data();
772     }
773
774     // gain
775     AliTRDCalDet *chamberGainFactor = 0;
776     if(getCalibs) {
777       entry = manager->Get("TRD/Calib/ChamberGainFactor", currRun, version, subVersion);
778       if(entry) {
779         entry->SetOwner(kFALSE);
780         chamberGainFactor = (AliTRDCalDet*)entry->GetObject();
781       }
782     }
783     AliTRDCalPad *padGainFactor = 0;
784     if(getCalibs) {
785       entry = manager->Get("TRD/Calib/LocalGainFactor", currRun, version, subVersion);
786       if(entry) {
787         entry->SetOwner(kFALSE);
788         padGainFactor = (AliTRDCalPad*)entry->GetObject();
789       }
790     }
791     Double_t runMeanGain, runRMSGain;
792     TVectorD chamberMeanGain(AliTRDcalibDB::kNdet);
793     TVectorD chamberRMSGain(AliTRDcalibDB::kNdet);
794     TVectorD smMeanGain(AliTRDcalibDB::kNsector);
795     TVectorD smRMSGain(AliTRDcalibDB::kNsector);
796     for(Int_t iNdet=0; iNdet<AliTRDcalibDB::kNdet; iNdet++) {chamberMeanGain[iNdet] = 0.0; chamberRMSGain[iNdet] = 0.0;}
797     for(Int_t iSm=0; iSm<AliTRDcalibDB::kNsector; iSm++) {smMeanGain[iSm] = 0.0; smRMSGain[iSm] = 0.0;}
798     TString parName("Gain");
799     if(getCalibs)
800       ProcessTRDCalibArray(chamberGainFactor, padGainFactor, 
801                            parName,
802                            runMeanGain, runRMSGain,
803                            chamberMeanGain, chamberRMSGain,
804                            smMeanGain, smRMSGain);
805
806     // noise/pedestals
807     AliTRDCalDet *chamberNoise = 0;
808     if(getCalibs) {
809       entry = manager->Get("TRD/Calib/DetNoise", currRun, version, subVersion);
810       if(entry) {
811         entry->SetOwner(kFALSE);
812         chamberNoise = (AliTRDCalDet*)entry->GetObject();
813       }
814     }
815     AliTRDCalPad *padNoise = 0;
816     if(getCalibs) {
817       entry = manager->Get("TRD/Calib/PadNoise", currRun, version, subVersion);
818       if(entry) {
819         entry->SetOwner(kFALSE);
820         padNoise = (AliTRDCalPad*)entry->GetObject();
821       }
822     }
823     Double_t runMeanNoise, runRMSNoise;
824     TVectorD chamberMeanNoise(AliTRDcalibDB::kNdet);
825     TVectorD chamberRMSNoise(AliTRDcalibDB::kNdet);
826     TVectorD smMeanNoise(AliTRDcalibDB::kNsector);
827     TVectorD smRMSNoise(AliTRDcalibDB::kNsector);
828     for(Int_t iNdet=0; iNdet<AliTRDcalibDB::kNdet; iNdet++) {chamberMeanNoise[iNdet] = 0.0; chamberRMSNoise[iNdet] = 0.0;}
829     for(Int_t iSm=0; iSm<AliTRDcalibDB::kNsector; iSm++) {smMeanNoise[iSm] = 0.0; smRMSNoise[iSm] = 0.0;}
830     parName = "Noise";
831     if(getCalibs)
832       ProcessTRDCalibArray(chamberNoise, padNoise, 
833                            parName,
834                            runMeanNoise, runRMSNoise,
835                            chamberMeanNoise, chamberRMSNoise,
836                            smMeanNoise, smRMSNoise);
837
838     // vdrift
839     AliTRDCalDet *chamberVdrift = 0;
840     if(getCalibs) {
841       entry = manager->Get("TRD/Calib/ChamberVdrift", currRun, version, subVersion);
842       if(entry) {
843         entry->SetOwner(kFALSE);
844         chamberVdrift = (AliTRDCalDet*)entry->GetObject();
845       }
846     }
847     AliTRDCalPad *padVdrift = 0;
848     if(getCalibs) {
849       entry = manager->Get("TRD/Calib/LocalVdrift", currRun, version, subVersion);
850       if(entry) {
851         entry->SetOwner(kFALSE);
852         padVdrift = (AliTRDCalPad*)entry->GetObject();
853       }
854     }
855     Double_t runMeanVdrift, runRMSVdrift;
856     TVectorD chamberMeanVdrift(AliTRDcalibDB::kNdet);
857     TVectorD chamberRMSVdrift(AliTRDcalibDB::kNdet);
858     TVectorD smMeanVdrift(AliTRDcalibDB::kNsector);
859     TVectorD smRMSVdrift(AliTRDcalibDB::kNsector);
860     for(Int_t iNdet=0; iNdet<AliTRDcalibDB::kNdet; iNdet++) {chamberMeanVdrift[iNdet] = 0.0; chamberRMSVdrift[iNdet] = 0.0;}
861     for(Int_t iSm=0; iSm<AliTRDcalibDB::kNsector; iSm++) {smMeanVdrift[iSm] = 0.0; smRMSVdrift[iSm] = 0.0;}
862     parName = "Vdrift";
863     if(getCalibs)
864       ProcessTRDCalibArray(chamberVdrift, padVdrift, 
865                            parName,
866                            runMeanVdrift, runRMSVdrift,
867                            chamberMeanVdrift, chamberRMSVdrift,
868                            smMeanVdrift, smRMSVdrift);
869
870     // T0
871     AliTRDCalDet *chamberT0 = 0;
872     if(getCalibs) {
873       entry = manager->Get("TRD/Calib/ChamberT0", currRun, version, subVersion);
874       if(entry) {
875         entry->SetOwner(kFALSE);
876         chamberT0 = (AliTRDCalDet*)entry->GetObject();
877       }
878     }
879     AliTRDCalPad *padT0 = 0;
880     if(getCalibs) {
881       entry = manager->Get("TRD/Calib/LocalT0", currRun, version, subVersion);
882       if(entry) {
883         entry->SetOwner(kFALSE);
884         padT0 = (AliTRDCalPad*)entry->GetObject();
885       }
886     }
887     Double_t runMeanT0, runRMST0;
888     TVectorD chamberMeanT0(AliTRDcalibDB::kNdet);
889     TVectorD chamberRMST0(AliTRDcalibDB::kNdet);
890     TVectorD smMeanT0(AliTRDcalibDB::kNsector);
891     TVectorD smRMST0(AliTRDcalibDB::kNsector);
892     for(Int_t iNdet=0; iNdet<AliTRDcalibDB::kNdet; iNdet++) {chamberMeanT0[iNdet] = 0.0; chamberRMST0[iNdet] = 0.0;}
893     for(Int_t iSm=0; iSm<AliTRDcalibDB::kNsector; iSm++) {smMeanT0[iSm] = 0.0; smRMST0[iSm] = 0.0;}
894     parName = "T0";
895     if(getCalibs)
896       ProcessTRDCalibArray(chamberT0, padT0, 
897                            parName,
898                            runMeanT0, runRMST0,
899                            chamberMeanT0, chamberRMST0,
900                            smMeanT0, smRMST0);
901
902     // status
903     AliTRDCalChamberStatus* chamberStatus = 0;
904     if(getCalibs) {
905       entry = manager->Get("TRD/Calib/ChamberStatus", currRun, version, subVersion);
906       if(entry) {
907         entry->SetOwner(kFALSE);
908         chamberStatus = (AliTRDCalChamberStatus*)entry->GetObject();
909       }
910     }
911     AliTRDCalPadStatus *padStatus = 0;
912     if(getCalibs) {
913       entry = manager->Get("TRD/Calib/PadStatus", currRun, version, subVersion);
914       if(entry) {
915         entry->SetOwner(kFALSE);
916         padStatus = (AliTRDCalPadStatus*)entry->GetObject();
917       }
918     }
919
920     // DCS FEE information
921     TObjArray *dcsArray = 0;
922     if(getDCS) {
923       entry = manager->Get("TRD/Calib/DCS");
924       if(entry) {
925         entry->SetOwner(kTRUE);
926         dcsArray = (TObjArray*)entry->GetObject();
927       }
928     }
929     AliTRDCalDCS *dcsSOR = 0;
930     AliTRDCalDCS *dcsEOR = 0;
931     if(getDCS && dcsArray) {
932       dcsSOR = (AliTRDCalDCS*)dcsArray->At(0);
933       dcsEOR = (AliTRDCalDCS*)dcsArray->At(1);
934     }
935
936     // Alignment information
937     // get the geometry from OCDB
938     TGeoManager *geoMan = 0x0;
939     if(getAlign) {
940       entry=manager->Get("GRP/Geometry/Data");
941       if(entry)
942         geoMan=(TGeoManager*)entry->GetObject();
943       else
944         cout << "Cannot get an entry for the geometry storage" << endl;
945     }
946     // get the alignment from OCDB
947     AliTRDalignment *alignMan=0;
948     if(getAlign && geoMan) {
949       entry=manager->Get("TRD/Align/Data", currRun, version, subVersion);
950       if(entry) {
951         alignMan = new AliTRDalignment();
952         cout << "storage for alignment = " << manager->GetDefaultStorage()->GetURI().Data() << endl;
953         alignMan->ReadDB(manager->GetDefaultStorage()->GetURI().Data(), "TRD/Align/Data", currRun, version, subVersion);
954       }
955       else {
956         cout << "Cannot get an entry for the alignment info" << endl;
957       }
958     }
959
960     Int_t kSuperModuleStatus[18] = {1, 1, 0, 0, 0, 0,     // super module status (1- installed, 0- not installed)
961                                     0, 1, 1, 1, 1, 0, 
962                                     0, 0, 0, 0, 0, 1};
963     Int_t kNRows[ 5] = {16, 16, 12, 16, 16};  // number of pad rows in the chambers from each of the 5 stacks
964     Int_t kNCols = 144;          // number of pad columns in the chambers from each of the 18 supermodules
965     Int_t kROB[5] = {8, 8, 6, 8, 8};   // number of read out boards(ROB) per chamber (6 in stack 2 and 8 in the rest)
966     Int_t kMCM = 16;                   // number of MCMs per ROB
967     for(Short_t iLayer=0; iLayer<AliTRDgeometry::kNlayer; iLayer++) {   // loop over layers
968       for(Short_t iSector=0; iSector<AliTRDgeometry::kNsector; iSector++) {  // loop over supermodules
969         if(kSuperModuleStatus[iSector]==0) 
970           continue;
971         Double_t alignSMPars[6];
972         for(Int_t ipar=0; ipar<6; ipar++) alignSMPars[ipar]=0.0;
973         if(getAlign && alignMan)
974           alignMan->GetSm(iSector, alignSMPars);
975         for(Short_t iStack=0; iStack<AliTRDgeometry::kNstack; iStack++) {    // loop over stacks
976           Short_t chamberNo = AliTRDgeometry::GetDetector(iLayer, iStack, iSector);
977           AliTRDCalROC *gainROC = 0;
978           if(padGainFactor) gainROC = padGainFactor->GetCalROC(chamberNo);
979           AliTRDCalROC *noiseROC = 0;
980           if(padNoise) noiseROC = padNoise->GetCalROC(chamberNo);
981           AliTRDCalROC *vdriftROC = 0;
982           if(padVdrift) vdriftROC = padVdrift->GetCalROC(chamberNo);
983           AliTRDCalROC *t0ROC = 0;
984           if(t0ROC) t0ROC = padT0->GetCalROC(chamberNo);
985           AliTRDCalSingleChamberStatus *statusROC = 0;
986           if(padStatus) statusROC = padStatus->GetCalROC(chamberNo);
987           TVectorD channelVector(kNRows[iStack]*kNCols);
988           TVectorD rowVector(kNRows[iStack]*kNCols);
989           TVectorD colVector(kNRows[iStack]*kNCols);
990           TVectorD statusVector(kNRows[iStack]*kNCols);
991           TVectorD gainVector(kNRows[iStack]*kNCols);
992           TVectorD noiseVector(kNRows[iStack]*kNCols);
993           TVectorD vdriftVector(kNRows[iStack]*kNCols);
994           TVectorD t0Vector(kNRows[iStack]*kNCols);
995           TVectorD padSuperRowVector(kNRows[iStack]*kNCols);
996           TVectorD padSuperColumnVector(kNRows[iStack]*kNCols);
997           for(Int_t ipar=0; ipar<kNRows[iStack]*kNCols; ipar++) {
998             channelVector[ipar] = 0; rowVector[ipar] = 0; colVector[ipar] = 0; 
999             statusVector[ipar] = 0; gainVector[ipar] = 0; noiseVector[ipar] = 0; 
1000             vdriftVector[ipar] = 0; t0Vector[ipar] = 0; padSuperRowVector[ipar] = 0; 
1001             padSuperColumnVector[ipar] = 0;
1002           }
1003           Int_t index = 0;
1004           if(getCalibs) {
1005             for(Short_t iRow=0; iRow<kNRows[iStack]; iRow++) {   // loop over pad rows
1006               for(Short_t iCol=0; iCol<kNCols; iCol++) {    // loop over pad columns
1007                 Short_t padSuperRow = iRow;
1008                 for(Int_t i=0; i<iStack; i++) padSuperRow = padSuperRow + kNRows[i]; 
1009                 padSuperRowVector[index] = padSuperRow;
1010                 Short_t padSuperColumn = iCol;
1011                 for(Int_t i=0; i<iSector; i++) padSuperColumn = padSuperColumn + kNCols;
1012                 padSuperColumnVector[index] = padSuperColumn;
1013                 Short_t channelNo = -1;
1014                 Float_t gain = -99.;
1015                 if(gainROC && chamberGainFactor) {
1016                   channelNo = gainROC->GetChannel(iCol, iRow);
1017                   gain = chamberGainFactor->GetValue(chamberNo) * gainROC->GetValue(iCol, iRow);
1018                 }
1019                 Float_t noise = -99.;
1020                 if(noiseROC && chamberNoise)
1021                   noise = chamberNoise->GetValue(chamberNo) * noiseROC->GetValue(iCol, iRow);
1022                 Float_t vdrift = -99.;
1023                 if(vdriftROC && chamberVdrift)
1024                   vdrift = chamberVdrift->GetValue(chamberNo) * vdriftROC->GetValue(iCol, iRow);
1025                 Float_t t0 = -99.;
1026                 if(t0ROC && chamberT0)
1027                   t0 = chamberT0->GetValue(chamberNo) + t0ROC->GetValue(iCol, iRow);
1028                 Int_t status = -99;
1029                 if(statusROC)
1030                   status = statusROC->GetStatus(iCol, iRow);
1031                 channelVector[index] = channelNo;
1032                 rowVector[index] = iRow;
1033                 colVector[index] = iCol;
1034                 statusVector[index] = status;
1035                 gainVector[index] = gain;
1036                 noiseVector[index] = noise;
1037                 vdriftVector[index] = vdrift;
1038                 t0Vector[index] = t0;
1039                 index++;
1040               }  // end loop over pad columns
1041             }  // end loop over pad rows
1042           }   // end if(getCalibs)
1043
1044           // get the dcs information
1045           AliTRDCalDCSFEE *dcsfeeSOR = 0;
1046           AliTRDCalDCSFEE *dcsfeeEOR = 0;
1047           if(getDCS) {
1048             if(dcsSOR) dcsfeeSOR = dcsSOR->GetCalDCSFEEObj(chamberNo);
1049             if(dcsEOR) dcsfeeEOR = dcsEOR->GetCalDCSFEEObj(chamberNo);
1050           }
1051           
1052           Bool_t sorAndEor = kFALSE;
1053           if(getDCS && dcsfeeSOR && dcsfeeEOR) sorAndEor = kTRUE;
1054           if(getDCS && !dcsfeeSOR && dcsfeeEOR) dcsfeeSOR = dcsfeeEOR;
1055           TVectorD robVector(kROB[iStack]*kMCM);
1056           TVectorD mcmVector(kROB[iStack]*kMCM);
1057           TVectorD sorandeorVector(kROB[iStack]*kMCM);
1058           TVectorD gsmSorVector(kROB[iStack]*kMCM);
1059           TVectorD gsmDeltaVector(kROB[iStack]*kMCM);
1060           TVectorD nimSorVector(kROB[iStack]*kMCM);
1061           TVectorD nimDeltaVector(kROB[iStack]*kMCM);
1062           TVectorD nevSorVector(kROB[iStack]*kMCM);
1063           TVectorD nevDeltaVector(kROB[iStack]*kMCM);
1064           TVectorD nptSorVector(kROB[iStack]*kMCM);
1065           TVectorD nptDeltaVector(kROB[iStack]*kMCM);
1066           TVectorD mcmSuperRowVector(kROB[iStack]*kMCM);
1067           TVectorD mcmSuperColumnVector(kROB[iStack]*kMCM);
1068           for(Int_t ipar=0; ipar<kROB[iStack]*kMCM; ipar++) {
1069             robVector[ipar] = 0; mcmVector[ipar] = 0; sorandeorVector[ipar] = 0; 
1070             gsmSorVector[ipar] = 0; gsmDeltaVector[ipar] = 0; nimSorVector[ipar] = 0; 
1071             nimDeltaVector[ipar] = 0; nevSorVector[ipar] = 0; nevDeltaVector[ipar] = 0; 
1072             nptSorVector[ipar] = 0; nptDeltaVector[ipar] = 0; mcmSuperRowVector[ipar] = 0; 
1073             mcmSuperColumnVector[ipar] = 0;
1074           }
1075
1076           Int_t robsRowDirection = kNRows[iStack]/4;    // 4 or 3 ROBs per chamber in row direction
1077           Int_t index1 = 0; 
1078           if(getDCS && (dcsfeeSOR || dcsfeeEOR) && dcsfeeSOR->GetStatusBit()==0) {
1079             for(Int_t iROB=0; iROB<kROB[iStack]; iROB++) { // loop over ROBs
1080               for(Int_t iMCM=0; iMCM<kMCM; iMCM++) {  // loop over MCMs
1081                 Short_t superRowMCM = iMCM%4;            // 4 MCMs per ROB in row direction
1082                 superRowMCM += 4*(iROB%robsRowDirection);   // now we have the row of this MCM inside one chamber
1083                 for(Int_t kk=0; kk<iStack; kk++) superRowMCM += kNRows[kk];   // add number of rows in previous stacks
1084                 Short_t superColumnMCM = iMCM/4;        // 4 MCMs per ROB in column direction
1085                 superColumnMCM += 4*(iROB/robsRowDirection);    // should yield 0 or 1 (2 ROBs per chamber in col direction)
1086                 superColumnMCM += iSector*8;
1087                 mcmSuperRowVector[index1] = superRowMCM;
1088                 mcmSuperColumnVector[index1] = superColumnMCM;
1089                 Int_t gsm = dcsfeeSOR->GetMCMGlobalState(iROB, iMCM);
1090                 Int_t nim = dcsfeeSOR->GetMCMStateNI(iROB, iMCM);
1091                 Int_t nev = dcsfeeSOR->GetMCMEventCnt(iROB, iMCM);
1092                 Int_t npt = dcsfeeSOR->GetMCMPtCnt(iROB, iMCM);
1093                 Int_t dgsm = -100000;
1094                 Int_t dnim = -100000;
1095                 Int_t dnev = -100000;
1096                 Int_t dnpt = -100000;
1097                 if(sorAndEor) {
1098                   dgsm = gsm - dcsfeeEOR->GetMCMGlobalState(iROB, iMCM);
1099                   dnim = nim - dcsfeeEOR->GetMCMStateNI(iROB, iMCM);
1100                   dnev = nev - dcsfeeEOR->GetMCMEventCnt(iROB, iMCM);
1101                   dnpt = npt - dcsfeeEOR->GetMCMPtCnt(iROB, iMCM);
1102                   if(gsm==-1 && dgsm==0) dgsm = -100000;
1103                   if(nim==-1 && dnim==0) dnim = -100000;
1104                   if(nev==-1 && dnev==0) dnev = -100000;
1105                   if(npt==-1 && dnpt==0) dnpt = -100000;
1106                 }
1107                 robVector[index1] = iROB;
1108                 mcmVector[index1] = iMCM;
1109                 sorandeorVector[index1] = sorAndEor;
1110                 gsmSorVector[index1] = gsm;
1111                 gsmDeltaVector[index1] = dgsm;
1112                 nimSorVector[index1] = nim;
1113                 nimDeltaVector[index1] = dnim;
1114                 nevSorVector[index1] = nev;
1115                 nevDeltaVector[index1] = dnev;
1116                 nptSorVector[index1] = npt;
1117                 nptDeltaVector[index1] = dnpt;
1118                 index1++;
1119               }  // end loop over MCMs
1120             }  // end loop over ROBs
1121           }  // end if(getDCS ...)
1122
1123           Double_t alignChamberPars[6];
1124           for(Int_t ipar=0; ipar<6; ipar++) alignChamberPars[ipar]=0;
1125           if(getAlign && alignMan)
1126             alignMan->GetCh(chamberNo, alignChamberPars);
1127
1128           (*treeStreamer)<< "TRDcalibDetails"
1129                          << "run=" << currRun
1130                          << "SuperModule=" << iSector
1131                          << "Stack=" << iStack
1132                          << "Layer=" << iLayer
1133                          << "Chamber=" << chamberNo;
1134           if(getAlign)
1135             (*treeStreamer)<< "TRDcalibDetails"
1136                            << "Align_SM_ShiftRphi=" << alignSMPars[0]
1137                            << "Align_SM_ShiftZ=" << alignSMPars[1]
1138                            << "Align_SM_ShiftR=" << alignSMPars[2]
1139                            << "Align_SM_RotRphi=" << alignSMPars[3]
1140                            << "Align_SM_RotZ=" << alignSMPars[4]
1141                            << "Align_SM_RotR=" << alignSMPars[5]
1142                            << "Align_Ch_ShiftRphi=" << alignChamberPars[0]
1143                            << "Align_Ch_ShiftZ=" << alignChamberPars[1]
1144                            << "Align_Ch_ShiftR=" << alignChamberPars[2]
1145                            << "Align_Ch_RotRphi=" << alignChamberPars[3]
1146                            << "Align_Ch_RotZ=" << alignChamberPars[4]
1147                            << "Align_Ch_RotR=" << alignChamberPars[5];
1148           if(getCalibs)
1149             (*treeStreamer)<< "TRDcalibDetails"
1150                            << "Gain_Mean_Run=" << runMeanGain
1151                            << "Gain_RMS_Run=" << runRMSGain
1152                            << "Gain_Mean_SuperModule=" << smMeanGain[iSector]
1153                            << "Gain_RMS_SuperModule=" << smRMSGain[iSector]
1154                            << "Gain_Mean_Chamber=" << chamberMeanGain[chamberNo]
1155                            << "Gain_RMS_Chamber=" << chamberRMSGain[chamberNo]
1156                            << "Noise_Mean_Run=" << runMeanNoise
1157                            << "Noise_RMS_Run=" << runRMSNoise
1158                            << "Noise_Mean_SuperModule=" << smMeanNoise[iSector]
1159                            << "Noise_RMS_SuperModule=" << smRMSNoise[iSector]
1160                            << "Noise_Mean_Chamber=" << chamberMeanNoise[chamberNo]
1161                            << "Noise_RMS_Chamber=" << chamberRMSNoise[chamberNo]
1162                            << "Vdrift_Mean_Run=" << runMeanVdrift
1163                            << "Vdrift_RMS_Run=" << runRMSVdrift
1164                            << "Vdrift_Mean_SuperModule=" << smMeanVdrift[iSector]
1165                            << "Vdrift_RMS_SuperModule=" << smRMSVdrift[iSector]
1166                            << "Vdrift_Mean_Chamber=" << chamberMeanVdrift[chamberNo]
1167                            << "Vdrift_RMS_Chamber=" << chamberRMSVdrift[chamberNo]
1168                            << "T0_Mean_Run=" << runMeanT0
1169                            << "T0_RMS_Run=" << runRMST0
1170                            << "T0_Mean_SuperModule=" << smMeanT0[iSector]
1171                            << "T0_RMS_SuperModule=" << smRMST0[iSector]
1172                            << "T0_Mean_Chamber=" << chamberMeanT0[chamberNo]
1173                            << "T0_RMS_Chamber=" << chamberRMST0[chamberNo]
1174                            << "Channel.=" << &channelVector
1175                            << "Row.=" << &rowVector
1176                            << "Column.=" << &colVector
1177                            << "PadSuperRow.=" << &padSuperRowVector
1178                            << "PadSuperColumn.=" << &padSuperColumnVector
1179                            << "Status.=" << &statusVector
1180                            << "Gain.=" << &gainVector
1181                            << "Noise.=" << &noiseVector
1182                            << "Vdrift.=" << &vdriftVector
1183                            << "T0.=" << &t0Vector;
1184           if(getDCS)
1185             (*treeStreamer)<< "TRDcalibDetails"
1186                            << "ROB.=" << &robVector
1187                            << "MCM.=" << &mcmVector
1188                            << "MCMSuperRow.=" << &mcmSuperRowVector
1189                            << "MCMSuperColumn.=" << &mcmSuperColumnVector
1190                            << "SORandEOR.=" << &sorandeorVector
1191                            << "gsmSOR.=" << &gsmSorVector
1192                            << "gsmDelta.=" << &gsmDeltaVector
1193                            << "nimSOR.=" << &nimSorVector
1194                            << "nimDelta.=" << &nimDeltaVector
1195                            << "nevSOR.=" << &nevSorVector
1196                            << "nevDelta.=" << &nevDeltaVector
1197                            << "nptSOR.=" << &nptSorVector
1198                            << "nptDelta.=" << &nptDeltaVector;
1199           (*treeStreamer)<< "TRDcalibDetails"
1200                          << "\n";
1201         }  // end loop over stacks
1202       }  // end loop over supermodules
1203     }  // end loop over layers
1204     
1205     // add the run number to the list of runs
1206     runs.ResizeTo(runs.GetNoElements()+1);
1207     runs[runs.GetNoElements()-1] = currRun;
1208     
1209     // do cleaning
1210     if(chamberGainFactor) delete chamberGainFactor;
1211     if(padGainFactor) delete padGainFactor;
1212     if(chamberNoise) delete chamberNoise;
1213     if(padNoise) delete padNoise;
1214     if(chamberVdrift) delete chamberVdrift;
1215     if(padVdrift) delete padVdrift;
1216     if(chamberT0) delete chamberT0;
1217     if(padT0) delete padT0;
1218     if(chamberStatus) delete chamberStatus;
1219     if(padStatus) delete padStatus;
1220
1221     // check if we still have run numbers in the file or provided range
1222     if(runListFilename[0]=='\0' && firstRun!=-1 && lastRun!=-1) {
1223       currRun++;
1224       if(currRun>lastRun) break;
1225     }
1226     if(runListFilename[0]!='\0' && in.eof())
1227       break;
1228   }   // end loop over runs
1229
1230   treeStreamer->GetFile()->cd();
1231   runs.Write("runs");
1232   delete treeStreamer;
1233   return kTRUE;
1234   //  delete treeStreamerDCS;
1235 }
1236
1237 //_________________________________________________________________________
1238 void AliTRDCalibViewer::DumpCalibToTree(const Char_t* inFilename, const Char_t* outFilename)
1239 {
1240   //
1241   //  extract info from CalPad objects and dump them into a tree to be viewed
1242   //
1243   TTreeSRedirector *treeStreamer = new TTreeSRedirector(outFilename);
1244   //open file and retrieve list of calPad objects
1245   TFile f(inFilename);
1246   TList *l=(TList*)f.GetListOfKeys();
1247
1248   TObjArray arrCalPads;
1249   TObjArray arrSMmean;
1250   TObjArray arrSMrms;
1251   arrCalPads.SetOwner();
1252   arrSMmean.SetOwner();
1253   arrSMrms.SetOwner();
1254
1255   TIter next(l);
1256   TKey *k=0x0;
1257   while ( (k=(TKey*)next()) ){
1258     AliTRDCalPad *pad=dynamic_cast<AliTRDCalPad*>(k->ReadObj());
1259     if (!pad) continue;
1260     arrCalPads.Add(pad);
1261
1262     TVectorD *smMean=new TVectorD(AliTRDcalibDB::kNsector);
1263     TVectorD *smRMS=new TVectorD(AliTRDcalibDB::kNsector);
1264
1265     arrSMmean.Add(smMean);
1266     arrSMrms.Add(smRMS);
1267
1268     ProcessTRDCalibArray(pad, *smMean, *smRMS);
1269   }
1270
1271   Int_t kSuperModuleStatus[18] = {1, 1, 0, 0, 0, 0,     // super module status (1- installed, 0- not installed)
1272       0, 1, 1, 1, 1, 0,
1273       0, 0, 0, 0, 0, 1};
1274
1275   AliTRDgeometry trdGeom;
1276   Int_t kNRows[5] = {16, 16, 12, 16, 16};  // number of pad rows in the chambers from each of the 5 stacks
1277   
1278   for(Short_t iLayer=0; iLayer<AliTRDgeometry::kNlayer; iLayer++) {   // loop over layers
1279     for(Short_t iSM=0; iSM<AliTRDgeometry::kNsector; iSM++) {  // loop over supermodules
1280       if(kSuperModuleStatus[iSM]==0)
1281         continue;
1282       
1283       for(Short_t iStack=0; iStack<AliTRDgeometry::kNstack; iStack++) {    // loop over stacks
1284         AliTRDpadPlane &plane=*trdGeom.GetPadPlane(iLayer, iStack);
1285
1286         Short_t chamberNo = AliTRDgeometry::GetDetector(iLayer, iStack, iSM);
1287         const Int_t nrows=plane.GetNrows();
1288         const Int_t ncols=plane.GetNcols();
1289         const Int_t nchannels=nrows*ncols;
1290 //         printf("chamberNo: %d (%03d-%03d-%03d)\n", chamberNo,nrows,ncols,nchannels);
1291         
1292         TVectorD channelVector(nchannels);
1293         TVectorD rowVector(nchannels);
1294         TVectorD colVector(nchannels);
1295         
1296         TVectorD gxVector(nchannels);
1297         TVectorD gyVector(nchannels);
1298         TVectorD gzVector(nchannels);
1299         
1300         TVectorD padSuperRowVector(nchannels);
1301         TVectorD padSuperColumnVector(nchannels);
1302         
1303         Int_t index = 0;
1304         for(Short_t iRow=0; iRow<nrows; iRow++) {   // loop over pad rows
1305           for(Short_t iCol=0; iCol<ncols; iCol++) {    // loop over pad columns
1306             Short_t padSuperRow = iRow;
1307             for(Int_t i=0; i<iStack; i++) padSuperRow = padSuperRow + kNRows[i];
1308             padSuperRowVector.GetMatrixArray()[index] = padSuperRow;
1309             
1310             Short_t padSuperColumn = iCol;
1311             for(Int_t i=0; i<iSM; i++) padSuperColumn = padSuperColumn + ncols;
1312             padSuperColumnVector.GetMatrixArray()[index] = padSuperColumn;
1313             
1314             rowVector.GetMatrixArray()[index] = iRow;
1315             colVector.GetMatrixArray()[index] = iCol;
1316             
1317             index++;
1318           }  // end loop over pad columns
1319         }  // end loop over pad rows
1320         
1321         (*treeStreamer)<< "TRDcalibDetails"
1322           << "SuperModule=" << iSM
1323           << "Stack=" << iStack
1324           << "Layer=" << iLayer
1325           << "Chamber=" << chamberNo
1326             //geographical information
1327           << "Channel.=" << &channelVector
1328           << "Row.=" << &rowVector
1329           << "Column.=" << &colVector
1330           << "PadSuperRow.=" << &padSuperRowVector
1331           << "PadSuperColumn.=" << &padSuperColumnVector;
1332 //          << "gx.=" << &gxVector
1333 //          << "gy.=" << &gyVector
1334 //          << "gz.=" << &gzVector;
1335           
1336         //
1337         // pad calibrations
1338         //
1339         TObjArray arrTrash;
1340         arrTrash.SetOwner();
1341         Int_t ncalib=arrCalPads.GetEntriesFast();
1342         for (Int_t iCalib=0; iCalib<ncalib; ++iCalib){
1343           AliTRDCalPad *pad=(AliTRDCalPad*)arrCalPads.UncheckedAt(iCalib);
1344           AliTRDCalROC *calROC=pad->GetCalROC(chamberNo);
1345           
1346           TVectorD &smMean=*((TVectorD*)arrSMmean.UncheckedAt(iCalib));
1347           TVectorD &smRMS=*((TVectorD*)arrSMrms.UncheckedAt(iCalib));
1348           
1349           TString calibName=pad->GetName();
1350           
1351           TVectorD *valueVector=new TVectorD(nchannels);
1352           arrTrash.Add(valueVector);
1353           
1354           Double_t rocMean=0;
1355           Double_t rocRMS=0;
1356           Double_t rocMedian=0;
1357           
1358           if (calROC){
1359             Int_t index2 = 0;
1360             for(Short_t iRow=0; iRow<nrows; iRow++) {
1361               for(Short_t iCol=0; iCol<ncols; iCol++) {
1362                 valueVector->GetMatrixArray()[index2] = calROC->GetValue(iCol,iRow);
1363                 index2++;
1364               }
1365             }
1366             rocMean   = calROC->GetMean();
1367             rocRMS    = calROC->GetRMS();
1368             rocMedian = calROC->GetMedian();
1369             //check for NaN
1370             if ( !(rocMean<1e30) ) rocMean=0;
1371             if ( !(rocRMS<1e30) ) rocRMS=0;
1372             
1373 //             printf("mean:   %f\n",rocMean);
1374 //             printf("rms:    %f\n",rocRMS);
1375 //             printf("median: %f\n",rocMedian);
1376           }
1377           
1378           (*treeStreamer)<< "TRDcalibDetails"
1379             // statistical information
1380             << (Char_t*)((calibName+"_Mean_SuperModule=").Data()) << smMean[iSM]
1381             << (Char_t*)((calibName+"_RMS_SuperModule=").Data())  << smRMS[iSM]
1382             << (Char_t*)((calibName+"_Mean_Chamber=").Data())     << rocMean
1383             << (Char_t*)((calibName+"_RMS_Chamber=").Data())      << rocRMS
1384             << (Char_t*)((calibName+"_Median_Chamber=").Data())   << rocMedian
1385             //pad by pad information
1386             << (Char_t*)((calibName+".=").Data()) << valueVector;
1387           
1388         }   // end loop over calib objects
1389         
1390         (*treeStreamer)<< "TRDcalibDetails"
1391           << "\n";
1392         arrTrash.Delete();
1393       }  // end loop over stacks
1394     }  // end loop over supermodules
1395   }  // end loop over layers
1396   delete treeStreamer;
1397 }
1398
1399 //_____________________________________________________________________________
1400 void AliTRDCalibViewer::ProcessTRDCalibArray(AliTRDCalDet* chamberCalib, AliTRDCalPad *padCalib,
1401                                              TString parName,
1402                                              Double_t &runValue, Double_t &runRMS,
1403                                              TVectorD &chamberValues, TVectorD &chamberValuesRMS,
1404                                              TVectorD &superModuleValues, TVectorD &superModuleValuesRMS) {
1405   // Process the calibrations for a given run.
1406   // Calculates the run and chamber wise averages
1407   //
1408   if(!chamberCalib) return;
1409   if(!padCalib) return;
1410   Int_t kSuperModuleStatus[18] = {1, 1, 0, 0, 0, 0,     // super module status (1- installed, 0- not installed)
1411                                   0, 1, 1, 1, 1, 0, 
1412                                   0, 0, 0, 0, 0, 1};
1413
1414   // initialize the histograms used for extracting the mean and RMS
1415   TH1F *runWiseHisto = new TH1F("runHisto", "runHisto", 200, -10.0, 10.0);
1416   TH1F *superModuleWiseHisto = new TH1F("smHisto", "smHisto", 200, -10.0, 10.0);
1417   TH1F *chamberWiseHisto = new TH1F("chamberHisto", "chamberHisto", 200, -10.0, 10.0);
1418
1419   // check if the calibration parameter is multiplicative or additive
1420   Bool_t multiplicative = kTRUE;
1421   if(!parName.CompareTo("T0")) multiplicative = kFALSE;
1422
1423   // first iteration (calculate all averages and RMS without discrimination on the SM average)
1424   runWiseHisto->Reset();
1425   for(Int_t iSM = 0; iSM<AliTRDcalibDB::kNsector; iSM++) {   // loop over supermodules
1426     // reset the super module histogram
1427     superModuleWiseHisto->Reset();
1428     // check if SM is installed
1429     if(!kSuperModuleStatus[iSM]) continue;
1430     for(Int_t iChamber=iSM*AliTRDcalibDB::kNstack*AliTRDcalibDB::kNlayer; 
1431         iChamber < (iSM+1)*AliTRDcalibDB::kNstack*AliTRDcalibDB::kNlayer; 
1432         iChamber++) {  // loop over chambers in this supermodule
1433       // get the chamber value
1434       Float_t chamberValue = chamberCalib->GetValue(iChamber);
1435       // get the ROC object
1436       AliTRDCalROC *chamberROC = padCalib->GetCalROC(iChamber);
1437       if(!chamberROC) 
1438         continue;
1439       chamberWiseHisto->Reset();
1440       for(Int_t iChannel = 0; iChannel < chamberROC->GetNchannels(); iChannel++){ // loop over channels
1441         // calculate the calibration parameter for this pad
1442         Float_t padValue;
1443         if(multiplicative)
1444           padValue = chamberValue * chamberROC->GetValue(iChannel);
1445         else
1446           padValue = chamberValue + chamberROC->GetValue(iChannel);
1447         // fill the run, SM and chamber wise histograms
1448         chamberWiseHisto->Fill(padValue);
1449         // if the parameter is Noise then check if the pad value is not a default one
1450         // Default value is now 1.2!!!! Check with Raphaelle for more informations
1451         if(parName.Contains("Noise") &&
1452            TMath::Abs(padValue-1.2)<0.00001) continue;
1453         superModuleWiseHisto->Fill(padValue);
1454         runWiseHisto->Fill(padValue);
1455       }  // end loop over channels
1456       // get the chamber wise mean and RMS
1457       chamberValues[iChamber] = chamberWiseHisto->GetMean();
1458       chamberValuesRMS[iChamber] = chamberWiseHisto->GetRMS();
1459     }  // end loop over chambers
1460     // SM wise mean and RMS
1461     superModuleValues[iSM] = superModuleWiseHisto->GetMean();
1462     superModuleValuesRMS[iSM] = superModuleWiseHisto->GetRMS();
1463   }  // end loop over supermodules
1464   // run wise mean and RMS
1465   runValue = runWiseHisto->GetMean();
1466   runRMS = runWiseHisto->GetRMS();
1467
1468   // Noise and Gain are finished processing
1469   if(parName.Contains("Noise") || parName.Contains("Gain"))
1470     return;
1471   // second iteration (calculate SM and run wise averages and RMS for Vdrift and T0)
1472   // The pads with calib parameter equal to the SM average are discarded (default value)
1473   runWiseHisto->Reset();
1474   for(Int_t iSM = 0; iSM<AliTRDcalibDB::kNsector; iSM++) {   // loop over supermodules
1475     superModuleWiseHisto->Reset();
1476     // eliminate the uninstalled super modules
1477     if(!kSuperModuleStatus[iSM]) continue;
1478     for(Int_t iChamber=iSM*AliTRDcalibDB::kNstack*AliTRDcalibDB::kNlayer; 
1479         iChamber < (iSM+1)*AliTRDcalibDB::kNstack*AliTRDcalibDB::kNlayer; 
1480         iChamber++) {  // loop over chambers
1481       // the chamber value
1482       Float_t chamberValue = chamberCalib->GetValue(iChamber);
1483       AliTRDCalROC *chamberROC = padCalib->GetCalROC(iChamber);
1484       if(!chamberROC) 
1485         continue;
1486       
1487       for(Int_t iChannel = 0; iChannel < chamberROC->GetNchannels(); iChannel++){ // loop over channels in a chamber
1488         // get the pad value
1489         Float_t padValue;
1490         if(multiplicative)
1491           padValue = chamberValue * chamberROC->GetValue(iChannel);
1492         else
1493           padValue = chamberValue + chamberROC->GetValue(iChannel);
1494         // eliminate from the average and RMS calculation all pads which
1495         // have the calib parameter equal with the SM average
1496         if((parName.Contains("Vdrift") || parName.Contains("T0")) && 
1497            TMath::Abs(padValue-superModuleValues[iSM])<0.00001) continue;
1498         superModuleWiseHisto->Fill(padValue);
1499         runWiseHisto->Fill(padValue);
1500       }   // end loop over channels
1501     }   // end loop over chambers 
1502     
1503     superModuleValues[iSM] = superModuleWiseHisto->GetMean();
1504     superModuleValuesRMS[iSM] = superModuleWiseHisto->GetRMS();
1505   }   // end loop over super modules
1506   runValue = runWiseHisto->GetMean();
1507   runRMS = runWiseHisto->GetRMS();
1508
1509   delete chamberWiseHisto;
1510   delete superModuleWiseHisto;
1511   delete runWiseHisto;
1512
1513   return;
1514 }
1515
1516
1517 //_____________________________________________________________________________
1518 void AliTRDCalibViewer::ProcessTRDCalibArray(AliTRDCalPad *padCalib,
1519                                              TVectorD &superModuleValues, TVectorD &superModuleValuesRMS)
1520 {
1521 // Process the calibrations for a given run.
1522 // Calculates the run and chamber wise averages
1523   
1524   if(!padCalib) return;
1525   Int_t kSuperModuleStatus[18] = {1, 1, 0, 0, 0, 0,     // super module status (1- installed, 0- not installed)
1526       0, 1, 1, 1, 1, 0,
1527       0, 0, 0, 0, 0, 1};
1528
1529   for(Int_t iSM = 0; iSM<AliTRDcalibDB::kNsector; iSM++) {   // loop over supermodules
1530     Double_t mean=0;
1531     Double_t rms=0;
1532     Double_t count=0;
1533     if(!kSuperModuleStatus[iSM]) continue;
1534
1535     // loop over chambers in this supermodule
1536     const Int_t nchambers=AliTRDcalibDB::kNstack*AliTRDcalibDB::kNlayer;
1537     for(Int_t iChamber=iSM*nchambers; iChamber < (iSM+1)*nchambers; iChamber++) {
1538       AliTRDCalROC *chamberROC = padCalib->GetCalROC(iChamber);
1539       if(!chamberROC)
1540         continue;
1541
1542       // loop over channels
1543       for(Int_t iChannel = 0; iChannel < chamberROC->GetNchannels(); iChannel++){ 
1544         mean+=chamberROC->GetValue(iChannel);
1545         rms+=chamberROC->GetValue(iChannel)*chamberROC->GetValue(iChannel);
1546         ++count;
1547       }  // end loop over channels
1548     }  // end loop over chambers
1549
1550     //calculate mean and rms
1551     if (count>0){
1552       mean/=count;
1553       rms/=count;
1554       rms=TMath::Sqrt(TMath::Abs(mean*mean-rms));
1555     }
1556     // SM wise mean and RMS
1557     superModuleValues[iSM]    = mean;
1558     superModuleValuesRMS[iSM] = rms;
1559   }  // end loop over supermodules
1560 }