1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 /* $Id: AliTRDCalibViewer.cxx 40390 2010-04-14 09:43:23Z cblume $ */
18 ///////////////////////////////////////////////////////////////////////////////
20 // Class which implements AliBaseCalibViewer for the TRD //
21 // used for the calibration monitor //
23 // Authors: Marian Ivanov (Marian.Ivanov@cern.ch) //
24 // Jens Wiechula (Jens.Wiechula@cern.ch) //
25 // Ionut Arsene (iarsene@cern.ch) //
27 ///////////////////////////////////////////////////////////////////////////////
42 #include <THashTable.h>
43 #include <TObjString.h>
44 #include <TTimeStamp.h>
45 #include <TObjString.h>
46 #include <TTreeStream.h>
50 #include <TDirectory.h>
51 #include <TFriendElement.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"
71 #include "AliTRDCalibViewer.h"
76 ClassImp(AliTRDCalibViewer)
78 //_____________________________________________________________________________
79 AliTRDCalibViewer::AliTRDCalibViewer()
83 // Default constructor (just call base class constructor)
87 //_____________________________________________________________________________
88 AliTRDCalibViewer::AliTRDCalibViewer(const AliTRDCalibViewer &c)
89 :AliBaseCalibViewer(c)
92 // copy constructor (just call base class copy constructor)
96 //_____________________________________________________________________________
97 AliTRDCalibViewer::AliTRDCalibViewer(TTree* tree)
98 :AliBaseCalibViewer(tree)
101 // Constructor (just call the corresponding base constructor)
105 //_____________________________________________________________________________
106 AliTRDCalibViewer::AliTRDCalibViewer(const char* fileName, const char* treeName)
107 :AliBaseCalibViewer(fileName, treeName)
110 // Constructor (just call the corresponding base constructor)
114 //_____________________________________________________________________________
115 AliTRDCalibViewer& AliTRDCalibViewer::operator=(const AliTRDCalibViewer& param) {
117 // assignment operator
120 fTreeMustBeDeleted = param.fTreeMustBeDeleted;
121 fListOfObjectsToBeDeleted = param.fListOfObjectsToBeDeleted;
122 fAbbreviation = param.fAbbreviation;
123 fAppendString = param.fAppendString;
127 //_____________________________________________________________________________
128 AliTRDCalibViewer::~AliTRDCalibViewer()
131 // AliTRDCalibViewer destructor
132 // do nothing, the base class destructor will do the job
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) {
141 // Get time information from OCDB by calling the DumpOCDBtoTree.C macro
143 DumpOCDBtoTree(runList, outFile, firstRun, lastRun,
144 TESTBIT(infoFlags,1), TESTBIT(infoFlags,2),
145 TESTBIT(infoFlags,3), TESTBIT(infoFlags,4),
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~"
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"
164 // - Do actual replacing:
165 // - ReplaceAll <variable> with <variable><fAbbreviation>, e.g. "CEQmean" -> "CEQmean~"
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"
173 // Now all the missing "~" should be added.
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();
199 //Int_t *varLengths = new Int_t[nVariables];
200 int *varLengths = new int[nVariables];
201 for (Int_t i = 0; i < nVariables; i++) {
202 varLengths[i] = ((TObjString*)listOfVariables->At(i))->String().Length();
204 //Int_t *normLengths = new Int_t[nNorm];
205 int *normLengths = new int[nNorm];
206 for (Int_t i = 0; i < nNorm; i++) {
207 normLengths[i] = ((TObjString*)listOfNormalizationVariables->At(i))->String().Length();
209 //Int_t *varSort = new Int_t[nVariables];
210 int *varSort = new int[nVariables];
211 TMath::Sort(nVariables, varLengths, varSort, kTRUE);
212 //Int_t *normSort = new Int_t[nNorm];
213 int *normSort = new int[nNorm];
214 TMath::Sort(nNorm, normLengths, normSort, kTRUE);
216 for (Int_t ivar = 0; ivar < nVariables; ivar++) {
217 // ***** save correct tokens *****
218 // first get the next variable:
219 searchString = ((TObjString*)listOfVariables->At(varSort[ivar]))->String();
220 // form replaceString:
222 for (Int_t i = 0; i < searchString.Length(); i++) {
223 replaceString.Append(searchString[i]);
224 if (i == 0) replaceString.Append(removeString);
226 // go through normalization:
227 for (Int_t inorm = 0; inorm < nNorm; inorm++) {
228 normString = ((TObjString*)listOfNormalizationVariables->At(normSort[inorm]))->String();
229 str.ReplaceAll(searchString + normString, replaceString + normString);
230 // like: str.ReplaceAll("CEQmean_Mean", "C!EQmean_Mean");
232 str.ReplaceAll(searchString + fAbbreviation, replaceString + fAbbreviation);
233 // like: str.ReplaceAll("CEQmean~", "C!EQmean~");
234 str.ReplaceAll(searchString + fAppendString, replaceString + fAppendString);
235 // like: str.ReplaceAll("CEQmean.fElements", "C!EQmean.fElements");
237 // ***** add missing extensions *****
238 str.ReplaceAll(searchString, replaceString + fAbbreviation);
239 // like: str.ReplaceAll("CEQmean", "C!EQmean~");
242 // ***** undo saving *****
243 str.ReplaceAll(removeString, "");
245 if (printDrawCommand) std::cout << "The string looks now like: " << str.Data() << std::endl;
251 //_____________________________________________________________________________
252 TObjArray* AliTRDCalibViewer::GetListOfVariables(Bool_t printList) {
254 // scan the tree - produces a list of available variables in the tree
255 // printList: print the list to the screen, after the scan is done
257 TObjArray* arr = new TObjArray();
262 Int_t nentries = fTree->GetListOfBranches()->GetEntries();
263 for (Int_t i = 0; i < nentries; i++) {
264 str = new TObjString(fTree->GetListOfBranches()->At(i)->GetName());
265 str->String().ReplaceAll("_Mean", "");
266 str->String().ReplaceAll("_RMS", "");
267 str->String().ReplaceAll("_Median", "");
268 str->String().ReplaceAll(".", "");
269 str->String().ReplaceAll("_Run", "");
270 str->String().ReplaceAll("_SuperModule", "");
271 str->String().ReplaceAll("_Chamber", "");
272 // add all the variables in the tree to a list
273 // exception make variables which are used for mapping, specified in AddAbbreviations()
274 // These two functions should be kept synchronized with respect to the mapping variables
275 if (!arr->FindObject(str) &&
276 !(str->String() == "run" ||
277 str->String() == "SuperModule" || str->String() == "Layer" || str->String() == "Stack" ||
278 str->String() == "Chamber" ||
279 str->String() == "Channel" || str->String() == "Row" || str->String() == "Column" ||
280 str->String() == "PadSuperRow" || str->String() == "PadSuperColumn" || str->String() == "MCMSuperRow"
281 || str->String() == "MCMSuperColumn" || str->String() == "MCM" || str->String() == "ROB")) {
286 // loop over all friends (if there are some) and add them to the list
287 if (fTree->GetListOfFriends()) {
288 for (Int_t ifriend = 0; ifriend < fTree->GetListOfFriends()->GetEntries(); ifriend++){
289 // printf("iterating through friendlist, currently at %i\n", ifriend);
290 // printf("working with %s\n", fTree->GetListOfFriends()->At(ifriend)->ClassName());
291 if (TString(fTree->GetListOfFriends()->At(ifriend)->ClassName()) != "TFriendElement") continue; // no friendElement found
292 TFriendElement *friendElement = (TFriendElement*)fTree->GetListOfFriends()->At(ifriend);
293 if (friendElement->GetTree() == 0) continue; // no tree found in friendElement
294 // printf("friend found \n");
295 for (Int_t i = 0; i < friendElement->GetTree()->GetListOfBranches()->GetEntries(); i++) {
296 // printf("iterating through friendelement entries, currently at %i\n", i);
297 str = new TObjString(friendElement->GetTree()->GetListOfBranches()->At(i)->GetName());
298 str->String().ReplaceAll("_Mean", "");
299 str->String().ReplaceAll("_RMS", "");
300 str->String().ReplaceAll("_Median", "");
301 str->String().ReplaceAll(".", "");
302 str->String().ReplaceAll("_Run", "");
303 str->String().ReplaceAll("_SuperModule", "");
304 str->String().ReplaceAll("_Chamber", "");
305 if (!(str->String() == "run" ||
306 str->String() == "SuperModule" || str->String() == "Layer" || str->String() == "Stack" ||
307 str->String() == "Chamber" ||
308 str->String() == "Channel" || str->String() == "Row" || str->String() == "Column" ||
309 str->String() == "PadSuperRow" || str->String() == "PadSuperColumn" || str->String() == "MCMSuperRow"
310 || str->String() == "MCMSuperColumn" || str->String() == "MCM" || str->String() == "ROB")){
311 // insert "<friendName>." at the beginning: (<friendName> is per default "R")
312 str->String().Insert(0, ".");
313 str->String().Insert(0, friendElement->GetName());
314 if (!arr->FindObject(str)) {
317 // printf("added string %s \n", str->String().Data());
321 } // if (fTree->GetListOfFriends())
326 TIterator* iter = arr->MakeIterator();
328 TObjString* currentStr = 0;
329 while ( (currentStr = (TObjString*)(iter->Next())) ) {
330 std::cout << currentStr->GetString().Data() << std::endl;
337 TObjArray* AliTRDCalibViewer::GetListOfNormalizationVariables(Bool_t printList) const{
339 // produces a list of available variables for normalization in the tree
340 // printList: print the list to the screen, after the scan is done
342 TObjArray* arr = new TObjArray();
343 arr->Add(new TObjString("_Mean_Run"));
344 arr->Add(new TObjString("_Mean_SuperModule"));
345 arr->Add(new TObjString("_Mean_Chamber"));
346 arr->Add(new TObjString("_Median_Run"));
347 arr->Add(new TObjString("_Median_SuperModule"));
348 arr->Add(new TObjString("_Median_Chamber"));
351 TIterator* iter = arr->MakeIterator();
353 TObjString* currentStr = 0;
354 while ((currentStr = (TObjString*)(iter->Next()))) {
355 std::cout << currentStr->GetString().Data() << std::endl;
362 void AliTRDCalibViewer::GetLayerSectorStack(TString trdString, Int_t& layerNo, Int_t& sectorNo, Int_t& stackNo) const {
363 // Get the layer, sector and stack numbers out of a string
364 // encoded with the following format:
365 // Layer%dSector%dStack%d
367 sscanf(trdString.Data(), "Layer%1dSector%02dStack%1d", &layerNo, §orNo, &stackNo);
372 //_____________________________________________________________________________
373 Int_t AliTRDCalibViewer::EasyDraw(const char* drawCommand, const char* sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
375 // easy drawing of data, use '~' for abbreviation of '.fElements'
376 // example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
377 // sector: sector-number - only the specified sector will be drwawn
378 // 'A'/'C' or 'a'/'c' - side A/C will be drawn
379 // 'ALL' - whole TPC will be drawn, projected on one side
380 // cuts: specifies cuts
381 // drawOptions: draw options like 'same'
382 // writeDrawCommand: write the command, that is passed to TTree::Draw
385 TString drawStr(drawCommand);
387 TString sectorStr(sector);
391 GetLayerSectorStack(sectorStr, layerNo, sectorNo, stackNo);
393 Warning("EasyDraw", "The sector string must always contain the Layer number!");
396 if(layerNo<0 || layerNo>5) {
397 Warning("EasyDraw", "The Layer number must be in the range [0,5] !");
400 if(sectorNo!=-1 && (sectorNo<0 || sectorNo>17)) {
401 Warning("EasyDraw", "The SuperModule number must be in the range [0,17] !");
404 if(stackNo!=-1 && (stackNo<0 || stackNo>4)) {
405 Warning("EasyDraw", "The Stack number must be in the range [0,4] !");
411 Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
412 if (dangerousToDraw) {
413 Warning("EasyDraw", "The draw string must not contain ':' or '>>'. Using only first variable for drawing!");
414 drawStr.Resize(drawStr.First(":"));
417 TString drawOptionsStr("");
419 Int_t rndNumber = rnd.Integer(10000);
421 if (drawOptions && strcmp(drawOptions, "") != 0)
422 drawOptionsStr += drawOptions;
424 drawOptionsStr += "profcolz";
426 const Int_t gkNRows[ 5] = {16, 16, 12, 16, 16}; // number of pad rows in the chambers from each of the 5 stacks
428 // draw calibration stuff
429 if(drawStr.Contains("Status") || drawStr.Contains("Gain") || drawStr.Contains("Noise") ||
430 drawStr.Contains("Vdrift") || drawStr.Contains("T0") ||
431 drawStr.Contains("gain") || drawStr.Contains("chiSquare")) {
432 if(sectorNo==-1 && stackNo==-1) { // plot the entire layer
433 drawStr += Form(":PadSuperColumn%s:PadSuperRow%s>>prof", fAppendString.Data(), fAppendString.Data());
434 drawStr += rndNumber;
435 drawStr += "(76,-0.5,75.5,2592,-0.5,2591.5)";
436 cutStr += Form("Layer==%d", layerNo);
438 else if(sectorNo!=-1 && stackNo==-1) { // plot a sector from a layer
439 drawStr += Form(":Column%s:PadSuperRow%s>>prof", fAppendString.Data(), fAppendString.Data());
440 drawStr += rndNumber;
441 drawStr += "(76,-0.5,75.5,144,-0.5,143.5)";
442 cutStr += Form("Layer==%d && SuperModule==%d", layerNo, sectorNo);
444 else if(sectorNo==-1 && stackNo!=-1) { // plot a stack from a layer
445 drawStr += Form(":PadSuperColumn%s:Row%s>>prof", fAppendString.Data(), fAppendString.Data());
446 drawStr += rndNumber;
447 drawStr += Form("(%d,-0.5,%d-0.5,2592,-0.5,2591.5)", gkNRows[stackNo], gkNRows[stackNo]);
448 cutStr += Form("Layer==%d && Stack==%d", layerNo, stackNo);
450 else { // the layer, sector and stack are determined -> so plot a chamber
451 drawStr += Form(":Column%s:Row%s>>prof", fAppendString.Data(), fAppendString.Data());
452 drawStr += rndNumber;
453 drawStr += Form("(%d,-0.5,%d-0.5,144,-0.5,143.5)", gkNRows[stackNo], gkNRows[stackNo]);
454 cutStr += Form("Layer==%d && SuperModule==%d && Stack==%d", layerNo, sectorNo, stackNo);
458 else if(drawStr.Contains("SORandEOR") ||
459 drawStr.Contains("gsmSOR") || drawStr.Contains("gsmDelta") ||
460 drawStr.Contains("nimSOR") || drawStr.Contains("nimDelta") ||
461 drawStr.Contains("nevSOR") || drawStr.Contains("nevDelta") ||
462 drawStr.Contains("nptSOR") || drawStr.Contains("nptDelta")) {
463 if(sectorNo==-1 && stackNo==-1) { // plot the entire layer
464 drawStr += Form(":MCMSuperColumn%s:MCMSuperRow%s>>prof", fAppendString.Data(), fAppendString.Data());
465 drawStr += rndNumber;
466 drawStr += "(76,-0.5,75.5,144,-0.5,143.5)";
467 cutStr += Form("Layer==%d", layerNo);
469 else if(sectorNo!=-1 && stackNo==-1) { // plot a sector from a layer
470 drawStr += Form(":MCMSuperColumn%s:MCMSuperRow%s>>prof", fAppendString.Data(), fAppendString.Data());
471 drawStr += rndNumber;
472 drawStr += "(76,-0.5,75.5,144,-0.5,143.5)";
473 cutStr += Form("Layer==%d && SuperModule==%d", layerNo, sectorNo);
475 else if(sectorNo==-1 && stackNo!=-1) { // plot a stack from a layer
476 drawStr += Form(":MCMSuperColumn%s:MCMSuperRow%s>>prof", fAppendString.Data(), fAppendString.Data());
477 drawStr += rndNumber;
478 // drawStr += Form("(%d,-0.5,%d-0.5,2592,-0.5,2591.5)", gkNRows[stackNo], gkNRows[stackNo]);
479 drawStr += "(76,-0.5,75.5,144,-0.5,143.5)";
480 cutStr += Form("Layer==%d && Stack==%d", layerNo, stackNo);
482 else { // the layer, sector and stack are determined -> so plot a chamber
483 drawStr += Form(":ROB%s:MCM%s>>prof", fAppendString.Data(), fAppendString.Data());
484 drawStr += rndNumber;
485 drawStr += Form("(16,-0.5,15.5,%d,-0.5,%d-0.5)", gkNRows[stackNo]/2, gkNRows[stackNo]/2);
486 cutStr += Form("Layer==%d && SuperModule==%d && Stack==%d", layerNo, sectorNo, stackNo);
489 // draw alignment stuff
490 else if(drawStr.Contains("Align")) {
491 if(sectorNo==-1 && stackNo==-1) { // plot the entire layer
492 drawStr += ":SuperModule:Stack>>prof";
493 drawStr += rndNumber;
494 drawStr += "(5,-0.5,4.5,18,-0.5,17.5)";
495 cutStr += Form("Layer==%d", layerNo);
497 else if(sectorNo!=-1 && stackNo==-1) { // plot a sector from a layer
498 drawStr += ":SuperModule:Stack>>prof";
499 drawStr += rndNumber;
500 drawStr += Form("(5,-0.5,4.5,1,%f,%f)", sectorNo-0.5, sectorNo+0.5);
501 cutStr += Form("Layer==%d && SuperModule==%d", layerNo, sectorNo);
503 else if(sectorNo==-1 && stackNo!=-1) { // plot a stack from a layer
504 drawStr += ":SuperModule:Stack>>prof";
505 drawStr += rndNumber;
506 drawStr += Form("(1,%f,%f,18,-0.5,17.5)", stackNo-0.5, stackNo+0.5);
507 cutStr += Form("Layer==%d && Stack==%d", layerNo, stackNo);
509 else { // the layer, sector and stack are determined -> so plot a chamber
510 drawStr += ":SuperModule:Stack>>prof";
511 drawStr += rndNumber;
512 drawStr += Form("(1,%f,%f,1,%f,%f)", stackNo-0.5, stackNo+0.5, sectorNo-0.5, sectorNo+0.5);
513 cutStr += Form("Layer==%d && SuperModule==%d && Stack==%d", layerNo, sectorNo, stackNo);
518 if (cuts && cuts[0] != 0) {
519 if (cutStr.Length() != 0) cutStr += "&& ";
524 drawStr.ReplaceAll(fAbbreviation, fAppendString);
525 cutStr.ReplaceAll(fAbbreviation, fAppendString);
526 if (writeDrawCommand) std::cout << "fTree->Draw(\"" << drawStr << "\", \"" << cutStr << "\", \"" << drawOptionsStr << "\");" << std::endl;
527 Int_t returnValue = fTree->Draw(drawStr.Data(), cutStr.Data(), drawOptionsStr.Data());
528 TString profName("prof");
529 profName += rndNumber;
530 TObject *obj = gDirectory->Get(profName.Data());
531 // set the names of the axes
532 TH1 *histObj = (TH1*)obj;
533 if(drawStr.Contains("Status") || drawStr.Contains("Gain") || drawStr.Contains("Noise") ||
534 drawStr.Contains("Vdrift") || drawStr.Contains("T0") ||
535 drawStr.Contains("gain") || drawStr.Contains("chiSquare")) {
536 histObj->GetXaxis()->SetTitle("Row");
537 histObj->GetYaxis()->SetTitle("Column");
539 else if(drawStr.Contains("SORandEOR") ||
540 drawStr.Contains("gsmSOR") || drawStr.Contains("gsmDelta") ||
541 drawStr.Contains("nimSOR") || drawStr.Contains("nimDelta") ||
542 drawStr.Contains("nevSOR") || drawStr.Contains("nevDelta") ||
543 drawStr.Contains("nptSOR") || drawStr.Contains("nptDelta")) {
544 histObj->GetXaxis()->SetTitle("MCM Row");
545 histObj->GetYaxis()->SetTitle("MCM Column");
547 else if(drawStr.Contains("Align")) {
548 histObj->GetXaxis()->SetTitle("Stack");
549 histObj->GetYaxis()->SetTitle("Sector");
552 if (obj && obj->InheritsFrom("TH1")) FormatHistoLabels((TH1*)obj);
556 //_____________________________________________________________________________
557 Int_t AliTRDCalibViewer::EasyDraw1D(const char* drawCommand, const char* sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
559 // easy drawing of data, use '~' for abbreviation of '.fElements'
560 // example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
561 // sector: sector-number - the specified sector will be drwawn
562 // 'A'/'C' or 'a'/'c' - side A/C will be drawn
563 // 'ALL' - whole TPC will be drawn, projected on one side
564 // cuts: specifies cuts
565 // drawOptions: draw options like 'same'
566 // writeDrawCommand: write the command, that is passed to TTree::Draw
569 TString drawStr(drawCommand);
571 TString sectorStr(sector);
575 GetLayerSectorStack(sectorStr, layerNo, sectorNo, stackNo);
577 Warning("EasyDraw", "The sector string must always contain the Layer number!");
580 if(layerNo<0 || layerNo>5) {
581 Warning("EasyDraw", "The Layer number must be in the range [0,5] !");
584 if(sectorNo!=-1 && (sectorNo<0 || sectorNo>17)) {
585 Warning("EasyDraw", "The Sector number must be in the range [0,17] !");
588 if(stackNo!=-1 && (stackNo<0 || stackNo>4)) {
589 Warning("EasyDraw", "The Stack number must be in the range [0,4] !");
593 TString drawOptionsStr(drawOptions);
596 if(sectorNo==-1 && stackNo==-1) // plot the entire layer
597 cutStr += Form("Layer==%d", layerNo);
598 else if(sectorNo!=-1 && stackNo==-1) // plot a sector from a layer
599 cutStr += Form("Layer==%d && SuperModule==%d", layerNo, sectorNo);
600 else if(sectorNo==-1 && stackNo!=-1) // plot a stack from a layer
601 cutStr += Form("Layer==%d && Stack==%d", layerNo, stackNo);
602 else // the layer, sector and stack are determined -> so plot a chamber
603 cutStr += Form("Layer==%d && SuperModule==%d && Stack==%d", layerNo, sectorNo, stackNo);
605 if(cuts && cuts[0] != 0) {
606 if (cutStr.Length() != 0) cutStr += "&& ";
612 drawStr.ReplaceAll(fAbbreviation, fAppendString);
613 cutStr.ReplaceAll(fAbbreviation, fAppendString);
614 if (writeDrawCommand) std::cout << "fTree->Draw(\"" << drawStr << "\", \"" << cutStr << "\", \"" << drawOptionsStr << "\");" << std::endl;
615 Int_t returnValue = fTree->Draw(drawStr.Data(), cutStr.Data(), drawOptionsStr.Data());
616 if (returnValue == -1) return -1;
618 TObject *obj = (gPad) ? gPad->GetPrimitive("htemp") : 0;
619 if (!obj) obj = (TH1F*)gDirectory->Get("htemp");
620 if (!obj) obj = gPad->GetPrimitive("tempHist");
621 if (!obj) obj = (TH1F*)gDirectory->Get("tempHist");
622 if (!obj) obj = gPad->GetPrimitive("Graph");
623 if (!obj) obj = (TH1F*)gDirectory->Get("Graph");
624 if (obj && obj->InheritsFrom("TH1")) FormatHistoLabels((TH1*)obj);
628 //_____________________________________________________________________________
629 Int_t AliTRDCalibViewer::EasyDraw(const char* drawCommand, Int_t chamber, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
631 // easy drawing of data, use '~' for abbreviation of '.fElements'
632 // example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0")
633 // sector: sector-number - only the specified sector will be drwawn
634 // cuts: specifies cuts
635 // drawOptions: draw options like 'same'
636 // writeDrawCommand: write the command, that is passed to TTree::Draw
638 if(chamber >= 0 && chamber < 540) {
639 Int_t superModuleNo = chamber/30;
640 Int_t stackNo = (chamber%30)/6;
641 Int_t layerNo = (chamber%30)%6;
643 snprintf(sectorChr,22, "Layer%iSector%iStack%i", layerNo, superModuleNo, stackNo);
644 return EasyDraw(drawCommand, sectorChr, cuts, drawOptions, writeDrawCommand);
646 Error("EasyDraw","The TRD contains only chamber from 0 to 539");
650 //_____________________________________________________________________________
651 Int_t AliTRDCalibViewer::EasyDraw1D(const char* drawCommand, Int_t chamber, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
653 // easy drawing of data, use '~' for abbreviation of '.fElements'
654 // example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0")
655 // sector: sector-number - the specified sector will be drwawn
656 // cuts: specifies cuts
657 // drawOptions: draw options like 'same'
658 // writeDrawCommand: write the command, that is passed to TTree::Draw
661 if (chamber >= 0 && chamber < 539) {
662 Int_t superModuleNo = chamber/30;
663 Int_t stackNo = (chamber%30)/6;
664 Int_t layerNo = (chamber%30)%6;
666 snprintf(sectorChr,22, "Layer%iSector%iStack%i", layerNo, superModuleNo, stackNo);
667 return EasyDraw1D(drawCommand, sectorChr, cuts, drawOptions, writeDrawCommand);
669 Error("EasyDraw1D","The TRD contains only chambers from 0 to 539");
673 //_____________________________________________________________________________
674 Bool_t AliTRDCalibViewer::DumpOCDBtoTreeDetails(const Char_t* runListFilename,
675 const Char_t* outFilename,
676 Int_t firstRun, Int_t lastRun,
677 const Char_t* storage,
684 // Retrieve TRD OCDB information for a given run list/range
687 if(runListFilename[0]!='\0' && firstRun==-1 && lastRun==-1) {
688 cout << "AliTRDCalibViewer::DumpOCDBtoTreeDetails(): You must provide at least a run range or an ascii filename with run numbers"
692 // initialize the OCDB manager
693 TString storageString = storage;
694 if(storageString.Contains("alien://")) {
695 TGrid::Connect("alien://");
697 AliCDBManager *manager = AliCDBManager::Instance();
698 if(storage[0]!='\0') {
699 manager->SetDefaultStorage(storage);
702 if(!manager->IsDefaultStorageSet()) {
703 cout << "AliTRDCalibViewer::DumpOCDBtoTreeDetails(): Default OCDB storage not set!!" << endl;
709 // open the ascii file
711 if(runListFilename[0]!='\0')
712 in.open(runListFilename);
714 // initialize the tree streamer
715 if(outFilename[0]=='\0') outFilename = "trdDetails.root";
716 TString calibFilename = outFilename;
718 TTreeSRedirector *treeStreamer = new TTreeSRedirector(calibFilename.Data());
721 if(runListFilename[0]=='\0' && firstRun!=-1 && lastRun!=-1)
728 if(runListFilename[0]!='\0') {
729 if(!(in>>currRun)) continue;
730 if(currRun < (firstRun==-1 ? 0 : firstRun) ||
731 currRun > (lastRun==-1 ? 999999999 : lastRun))
735 if(currRun>lastRun) break;
737 cout << "run = " << currRun << endl;
738 manager->SetRun(currRun);
740 // Get GRP data. If there is no proper GRP object for this run than
741 // this run is aborted
742 AliCDBEntry *entry = manager->Get("GRP/GRP/Data");
743 AliGRPObject* grpObject = 0;
745 entry->SetOwner(kFALSE);
746 grpObject = dynamic_cast<AliGRPObject*>(entry->GetObject());
754 cout << "No GRP info available for this run " << endl;
756 time_t startTimeGRP = 0;
757 TObjString runType("");
759 startTimeGRP = grpObject->GetTimeStart();
760 TTimeStamp start(grpObject->GetTimeStart());
761 TTimeStamp end(grpObject->GetTimeEnd());
762 cout << "Start time: " << start.GetDate()/10000 << "/"
763 << (start.GetDate()/100)-(start.GetDate()/10000)*100 << "/"
764 << start.GetDate()%100 << " "
765 << start.GetTime()/10000 << ":"
766 << (start.GetTime()/100)-(start.GetTime()/10000)*100 << ":"
767 << start.GetTime()%100 << endl;
768 cout << "End time: " << end.GetDate()/10000 << "/"
769 << (end.GetDate()/100)-(end.GetDate()/10000)*100 << "/"
770 << end.GetDate()%100 << " "
771 << end.GetTime()/10000 << ":"
772 << (end.GetTime()/100)-(end.GetTime()/10000)*100 << ":"
773 << end.GetTime()%100 << endl;
774 cout << "Run type = " << grpObject->GetRunType().Data() << endl;
775 runType = grpObject->GetRunType().Data();
779 AliTRDCalDet *chamberGainFactor = 0;
781 entry = manager->Get("TRD/Calib/ChamberGainFactor", currRun, version, subVersion);
783 entry->SetOwner(kFALSE);
784 chamberGainFactor = (AliTRDCalDet*)entry->GetObject();
787 AliTRDCalPad *padGainFactor = 0;
789 entry = manager->Get("TRD/Calib/LocalGainFactor", currRun, version, subVersion);
791 entry->SetOwner(kFALSE);
792 padGainFactor = (AliTRDCalPad*)entry->GetObject();
795 Double_t runMeanGain, runRMSGain;
796 TVectorD chamberMeanGain(AliTRDcalibDB::kNdet);
797 TVectorD chamberRMSGain(AliTRDcalibDB::kNdet);
798 TVectorD smMeanGain(AliTRDcalibDB::kNsector);
799 TVectorD smRMSGain(AliTRDcalibDB::kNsector);
800 for(Int_t iNdet=0; iNdet<AliTRDcalibDB::kNdet; iNdet++) {chamberMeanGain[iNdet] = 0.0; chamberRMSGain[iNdet] = 0.0;}
801 for(Int_t iSm=0; iSm<AliTRDcalibDB::kNsector; iSm++) {smMeanGain[iSm] = 0.0; smRMSGain[iSm] = 0.0;}
802 TString parName("Gain");
804 ProcessTRDCalibArray(chamberGainFactor, padGainFactor,
806 runMeanGain, runRMSGain,
807 chamberMeanGain, chamberRMSGain,
808 smMeanGain, smRMSGain);
811 AliTRDCalDet *chamberNoise = 0;
813 entry = manager->Get("TRD/Calib/DetNoise", currRun, version, subVersion);
815 entry->SetOwner(kFALSE);
816 chamberNoise = (AliTRDCalDet*)entry->GetObject();
819 AliTRDCalPad *padNoise = 0;
821 entry = manager->Get("TRD/Calib/PadNoise", currRun, version, subVersion);
823 entry->SetOwner(kFALSE);
824 padNoise = (AliTRDCalPad*)entry->GetObject();
827 Double_t runMeanNoise, runRMSNoise;
828 TVectorD chamberMeanNoise(AliTRDcalibDB::kNdet);
829 TVectorD chamberRMSNoise(AliTRDcalibDB::kNdet);
830 TVectorD smMeanNoise(AliTRDcalibDB::kNsector);
831 TVectorD smRMSNoise(AliTRDcalibDB::kNsector);
832 for(Int_t iNdet=0; iNdet<AliTRDcalibDB::kNdet; iNdet++) {chamberMeanNoise[iNdet] = 0.0; chamberRMSNoise[iNdet] = 0.0;}
833 for(Int_t iSm=0; iSm<AliTRDcalibDB::kNsector; iSm++) {smMeanNoise[iSm] = 0.0; smRMSNoise[iSm] = 0.0;}
836 ProcessTRDCalibArray(chamberNoise, padNoise,
838 runMeanNoise, runRMSNoise,
839 chamberMeanNoise, chamberRMSNoise,
840 smMeanNoise, smRMSNoise);
843 AliTRDCalDet *chamberVdrift = 0;
845 entry = manager->Get("TRD/Calib/ChamberVdrift", currRun, version, subVersion);
847 entry->SetOwner(kFALSE);
848 chamberVdrift = (AliTRDCalDet*)entry->GetObject();
851 AliTRDCalPad *padVdrift = 0;
853 entry = manager->Get("TRD/Calib/LocalVdrift", currRun, version, subVersion);
855 entry->SetOwner(kFALSE);
856 padVdrift = (AliTRDCalPad*)entry->GetObject();
859 Double_t runMeanVdrift, runRMSVdrift;
860 TVectorD chamberMeanVdrift(AliTRDcalibDB::kNdet);
861 TVectorD chamberRMSVdrift(AliTRDcalibDB::kNdet);
862 TVectorD smMeanVdrift(AliTRDcalibDB::kNsector);
863 TVectorD smRMSVdrift(AliTRDcalibDB::kNsector);
864 for(Int_t iNdet=0; iNdet<AliTRDcalibDB::kNdet; iNdet++) {chamberMeanVdrift[iNdet] = 0.0; chamberRMSVdrift[iNdet] = 0.0;}
865 for(Int_t iSm=0; iSm<AliTRDcalibDB::kNsector; iSm++) {smMeanVdrift[iSm] = 0.0; smRMSVdrift[iSm] = 0.0;}
868 ProcessTRDCalibArray(chamberVdrift, padVdrift,
870 runMeanVdrift, runRMSVdrift,
871 chamberMeanVdrift, chamberRMSVdrift,
872 smMeanVdrift, smRMSVdrift);
875 AliTRDCalDet *chamberT0 = 0;
877 entry = manager->Get("TRD/Calib/ChamberT0", currRun, version, subVersion);
879 entry->SetOwner(kFALSE);
880 chamberT0 = (AliTRDCalDet*)entry->GetObject();
883 AliTRDCalPad *padT0 = 0;
885 entry = manager->Get("TRD/Calib/LocalT0", currRun, version, subVersion);
887 entry->SetOwner(kFALSE);
888 padT0 = (AliTRDCalPad*)entry->GetObject();
891 Double_t runMeanT0, runRMST0;
892 TVectorD chamberMeanT0(AliTRDcalibDB::kNdet);
893 TVectorD chamberRMST0(AliTRDcalibDB::kNdet);
894 TVectorD smMeanT0(AliTRDcalibDB::kNsector);
895 TVectorD smRMST0(AliTRDcalibDB::kNsector);
896 for(Int_t iNdet=0; iNdet<AliTRDcalibDB::kNdet; iNdet++) {chamberMeanT0[iNdet] = 0.0; chamberRMST0[iNdet] = 0.0;}
897 for(Int_t iSm=0; iSm<AliTRDcalibDB::kNsector; iSm++) {smMeanT0[iSm] = 0.0; smRMST0[iSm] = 0.0;}
900 ProcessTRDCalibArray(chamberT0, padT0,
903 chamberMeanT0, chamberRMST0,
907 AliTRDCalChamberStatus* chamberStatus = 0;
909 entry = manager->Get("TRD/Calib/ChamberStatus", currRun, version, subVersion);
911 entry->SetOwner(kFALSE);
912 chamberStatus = (AliTRDCalChamberStatus*)entry->GetObject();
915 AliTRDCalPadStatus *padStatus = 0;
917 entry = manager->Get("TRD/Calib/PadStatus", currRun, version, subVersion);
919 entry->SetOwner(kFALSE);
920 padStatus = (AliTRDCalPadStatus*)entry->GetObject();
924 // DCS FEE information
925 TObjArray *dcsArray = 0;
927 entry = manager->Get("TRD/Calib/DCS");
929 entry->SetOwner(kTRUE);
930 dcsArray = (TObjArray*)entry->GetObject();
933 AliTRDCalDCS *dcsSOR = 0;
934 AliTRDCalDCS *dcsEOR = 0;
935 if(getDCS && dcsArray) {
936 dcsSOR = (AliTRDCalDCS*)dcsArray->At(0);
937 dcsEOR = (AliTRDCalDCS*)dcsArray->At(1);
940 // Alignment information
941 // get the geometry from OCDB
942 TGeoManager *geoMan = 0x0;
944 entry=manager->Get("GRP/Geometry/Data");
946 geoMan=(TGeoManager*)entry->GetObject();
948 cout << "Cannot get an entry for the geometry storage" << endl;
950 // get the alignment from OCDB
951 AliTRDalignment *alignMan=0;
952 if(getAlign && geoMan) {
953 entry=manager->Get("TRD/Align/Data", currRun, version, subVersion);
955 alignMan = new AliTRDalignment();
956 cout << "storage for alignment = " << manager->GetDefaultStorage()->GetURI().Data() << endl;
957 alignMan->ReadDB(manager->GetDefaultStorage()->GetURI().Data(), "TRD/Align/Data", currRun, version, subVersion);
960 cout << "Cannot get an entry for the alignment info" << endl;
964 Int_t kSuperModuleStatus[18] = {1, 1, 0, 0, 0, 0, // super module status (1- installed, 0- not installed)
967 Int_t kNRows[ 5] = {16, 16, 12, 16, 16}; // number of pad rows in the chambers from each of the 5 stacks
968 Int_t kNCols = 144; // number of pad columns in the chambers from each of the 18 supermodules
969 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)
970 Int_t kMCM = 16; // number of MCMs per ROB
971 for(Short_t iLayer=0; iLayer<AliTRDgeometry::kNlayer; iLayer++) { // loop over layers
972 for(Short_t iSector=0; iSector<AliTRDgeometry::kNsector; iSector++) { // loop over supermodules
973 if(kSuperModuleStatus[iSector]==0)
975 Double_t alignSMPars[6];
976 for(Int_t ipar=0; ipar<6; ipar++) alignSMPars[ipar]=0.0;
977 if(getAlign && alignMan)
978 alignMan->GetSm(iSector, alignSMPars);
979 for(Short_t iStack=0; iStack<AliTRDgeometry::kNstack; iStack++) { // loop over stacks
980 Short_t chamberNo = AliTRDgeometry::GetDetector(iLayer, iStack, iSector);
981 AliTRDCalROC *gainROC = 0;
982 if(padGainFactor) gainROC = padGainFactor->GetCalROC(chamberNo);
983 AliTRDCalROC *noiseROC = 0;
984 if(padNoise) noiseROC = padNoise->GetCalROC(chamberNo);
985 AliTRDCalROC *vdriftROC = 0;
986 if(padVdrift) vdriftROC = padVdrift->GetCalROC(chamberNo);
987 AliTRDCalROC *t0ROC = 0;
988 if(t0ROC) t0ROC = padT0->GetCalROC(chamberNo);
989 AliTRDCalSingleChamberStatus *statusROC = 0;
990 if(padStatus) statusROC = padStatus->GetCalROC(chamberNo);
991 TVectorD channelVector(kNRows[iStack]*kNCols);
992 TVectorD rowVector(kNRows[iStack]*kNCols);
993 TVectorD colVector(kNRows[iStack]*kNCols);
994 TVectorD statusVector(kNRows[iStack]*kNCols);
995 TVectorD gainVector(kNRows[iStack]*kNCols);
996 TVectorD noiseVector(kNRows[iStack]*kNCols);
997 TVectorD vdriftVector(kNRows[iStack]*kNCols);
998 TVectorD t0Vector(kNRows[iStack]*kNCols);
999 TVectorD padSuperRowVector(kNRows[iStack]*kNCols);
1000 TVectorD padSuperColumnVector(kNRows[iStack]*kNCols);
1001 for(Int_t ipar=0; ipar<kNRows[iStack]*kNCols; ipar++) {
1002 channelVector[ipar] = 0; rowVector[ipar] = 0; colVector[ipar] = 0;
1003 statusVector[ipar] = 0; gainVector[ipar] = 0; noiseVector[ipar] = 0;
1004 vdriftVector[ipar] = 0; t0Vector[ipar] = 0; padSuperRowVector[ipar] = 0;
1005 padSuperColumnVector[ipar] = 0;
1009 for(Short_t iRow=0; iRow<kNRows[iStack]; iRow++) { // loop over pad rows
1010 for(Short_t iCol=0; iCol<kNCols; iCol++) { // loop over pad columns
1011 Short_t padSuperRow = iRow;
1012 for(Int_t i=0; i<iStack; i++) padSuperRow = padSuperRow + kNRows[i];
1013 padSuperRowVector[index] = padSuperRow;
1014 Short_t padSuperColumn = iCol;
1015 for(Int_t i=0; i<iSector; i++) padSuperColumn = padSuperColumn + kNCols;
1016 padSuperColumnVector[index] = padSuperColumn;
1017 Short_t channelNo = -1;
1018 Float_t gain = -99.;
1019 if(gainROC && chamberGainFactor) {
1020 channelNo = gainROC->GetChannel(iCol, iRow);
1021 gain = chamberGainFactor->GetValue(chamberNo) * gainROC->GetValue(iCol, iRow);
1023 Float_t noise = -99.;
1024 if(noiseROC && chamberNoise)
1025 noise = chamberNoise->GetValue(chamberNo) * noiseROC->GetValue(iCol, iRow);
1026 Float_t vdrift = -99.;
1027 if(vdriftROC && chamberVdrift)
1028 vdrift = chamberVdrift->GetValue(chamberNo) * vdriftROC->GetValue(iCol, iRow);
1030 if(t0ROC && chamberT0)
1031 t0 = chamberT0->GetValue(chamberNo) + t0ROC->GetValue(iCol, iRow);
1034 status = statusROC->GetStatus(iCol, iRow);
1035 channelVector[index] = channelNo;
1036 rowVector[index] = iRow;
1037 colVector[index] = iCol;
1038 statusVector[index] = status;
1039 gainVector[index] = gain;
1040 noiseVector[index] = noise;
1041 vdriftVector[index] = vdrift;
1042 t0Vector[index] = t0;
1044 } // end loop over pad columns
1045 } // end loop over pad rows
1046 } // end if(getCalibs)
1048 // get the dcs information
1049 AliTRDCalDCSFEE *dcsfeeSOR = 0;
1050 AliTRDCalDCSFEE *dcsfeeEOR = 0;
1052 if(dcsSOR) dcsfeeSOR = dcsSOR->GetCalDCSFEEObj(chamberNo);
1053 if(dcsEOR) dcsfeeEOR = dcsEOR->GetCalDCSFEEObj(chamberNo);
1056 Bool_t sorAndEor = kFALSE;
1057 if(getDCS && dcsfeeSOR && dcsfeeEOR) sorAndEor = kTRUE;
1058 if(getDCS && !dcsfeeSOR && dcsfeeEOR) dcsfeeSOR = dcsfeeEOR;
1059 TVectorD robVector(kROB[iStack]*kMCM);
1060 TVectorD mcmVector(kROB[iStack]*kMCM);
1061 TVectorD sorandeorVector(kROB[iStack]*kMCM);
1062 TVectorD gsmSorVector(kROB[iStack]*kMCM);
1063 TVectorD gsmDeltaVector(kROB[iStack]*kMCM);
1064 TVectorD nimSorVector(kROB[iStack]*kMCM);
1065 TVectorD nimDeltaVector(kROB[iStack]*kMCM);
1066 TVectorD nevSorVector(kROB[iStack]*kMCM);
1067 TVectorD nevDeltaVector(kROB[iStack]*kMCM);
1068 TVectorD nptSorVector(kROB[iStack]*kMCM);
1069 TVectorD nptDeltaVector(kROB[iStack]*kMCM);
1070 TVectorD mcmSuperRowVector(kROB[iStack]*kMCM);
1071 TVectorD mcmSuperColumnVector(kROB[iStack]*kMCM);
1072 for(Int_t ipar=0; ipar<kROB[iStack]*kMCM; ipar++) {
1073 robVector[ipar] = 0; mcmVector[ipar] = 0; sorandeorVector[ipar] = 0;
1074 gsmSorVector[ipar] = 0; gsmDeltaVector[ipar] = 0; nimSorVector[ipar] = 0;
1075 nimDeltaVector[ipar] = 0; nevSorVector[ipar] = 0; nevDeltaVector[ipar] = 0;
1076 nptSorVector[ipar] = 0; nptDeltaVector[ipar] = 0; mcmSuperRowVector[ipar] = 0;
1077 mcmSuperColumnVector[ipar] = 0;
1080 Int_t robsRowDirection = kNRows[iStack]/4; // 4 or 3 ROBs per chamber in row direction
1082 if(getDCS && (dcsfeeSOR || dcsfeeEOR) && dcsfeeSOR->GetStatusBit()==0) {
1083 for(Int_t iROB=0; iROB<kROB[iStack]; iROB++) { // loop over ROBs
1084 for(Int_t iMCM=0; iMCM<kMCM; iMCM++) { // loop over MCMs
1085 Short_t superRowMCM = iMCM%4; // 4 MCMs per ROB in row direction
1086 superRowMCM += 4*(iROB%robsRowDirection); // now we have the row of this MCM inside one chamber
1087 for(Int_t kk=0; kk<iStack; kk++) superRowMCM += kNRows[kk]; // add number of rows in previous stacks
1088 Short_t superColumnMCM = iMCM/4; // 4 MCMs per ROB in column direction
1089 superColumnMCM += 4*(iROB/robsRowDirection); // should yield 0 or 1 (2 ROBs per chamber in col direction)
1090 superColumnMCM += iSector*8;
1091 mcmSuperRowVector[index1] = superRowMCM;
1092 mcmSuperColumnVector[index1] = superColumnMCM;
1093 Int_t gsm = dcsfeeSOR->GetMCMGlobalState(iROB, iMCM);
1094 Int_t nim = dcsfeeSOR->GetMCMStateNI(iROB, iMCM);
1095 Int_t nev = dcsfeeSOR->GetMCMEventCnt(iROB, iMCM);
1096 Int_t npt = dcsfeeSOR->GetMCMPtCnt(iROB, iMCM);
1097 Int_t dgsm = -100000;
1098 Int_t dnim = -100000;
1099 Int_t dnev = -100000;
1100 Int_t dnpt = -100000;
1102 dgsm = gsm - dcsfeeEOR->GetMCMGlobalState(iROB, iMCM);
1103 dnim = nim - dcsfeeEOR->GetMCMStateNI(iROB, iMCM);
1104 dnev = nev - dcsfeeEOR->GetMCMEventCnt(iROB, iMCM);
1105 dnpt = npt - dcsfeeEOR->GetMCMPtCnt(iROB, iMCM);
1106 if(gsm==-1 && dgsm==0) dgsm = -100000;
1107 if(nim==-1 && dnim==0) dnim = -100000;
1108 if(nev==-1 && dnev==0) dnev = -100000;
1109 if(npt==-1 && dnpt==0) dnpt = -100000;
1111 robVector[index1] = iROB;
1112 mcmVector[index1] = iMCM;
1113 sorandeorVector[index1] = sorAndEor;
1114 gsmSorVector[index1] = gsm;
1115 gsmDeltaVector[index1] = dgsm;
1116 nimSorVector[index1] = nim;
1117 nimDeltaVector[index1] = dnim;
1118 nevSorVector[index1] = nev;
1119 nevDeltaVector[index1] = dnev;
1120 nptSorVector[index1] = npt;
1121 nptDeltaVector[index1] = dnpt;
1123 } // end loop over MCMs
1124 } // end loop over ROBs
1125 } // end if(getDCS ...)
1127 Double_t alignChamberPars[6];
1128 for(Int_t ipar=0; ipar<6; ipar++) alignChamberPars[ipar]=0;
1129 if(getAlign && alignMan)
1130 alignMan->GetCh(chamberNo, alignChamberPars);
1132 (*treeStreamer)<< "TRDcalibDetails"
1133 << "run=" << currRun
1134 << "SuperModule=" << iSector
1135 << "Stack=" << iStack
1136 << "Layer=" << iLayer
1137 << "Chamber=" << chamberNo;
1139 (*treeStreamer)<< "TRDcalibDetails"
1140 << "Align_SM_ShiftRphi=" << alignSMPars[0]
1141 << "Align_SM_ShiftZ=" << alignSMPars[1]
1142 << "Align_SM_ShiftR=" << alignSMPars[2]
1143 << "Align_SM_RotRphi=" << alignSMPars[3]
1144 << "Align_SM_RotZ=" << alignSMPars[4]
1145 << "Align_SM_RotR=" << alignSMPars[5]
1146 << "Align_Ch_ShiftRphi=" << alignChamberPars[0]
1147 << "Align_Ch_ShiftZ=" << alignChamberPars[1]
1148 << "Align_Ch_ShiftR=" << alignChamberPars[2]
1149 << "Align_Ch_RotRphi=" << alignChamberPars[3]
1150 << "Align_Ch_RotZ=" << alignChamberPars[4]
1151 << "Align_Ch_RotR=" << alignChamberPars[5];
1153 (*treeStreamer)<< "TRDcalibDetails"
1154 << "Gain_Mean_Run=" << runMeanGain
1155 << "Gain_RMS_Run=" << runRMSGain
1156 << "Gain_Mean_SuperModule=" << smMeanGain[iSector]
1157 << "Gain_RMS_SuperModule=" << smRMSGain[iSector]
1158 << "Gain_Mean_Chamber=" << chamberMeanGain[chamberNo]
1159 << "Gain_RMS_Chamber=" << chamberRMSGain[chamberNo]
1160 << "Noise_Mean_Run=" << runMeanNoise
1161 << "Noise_RMS_Run=" << runRMSNoise
1162 << "Noise_Mean_SuperModule=" << smMeanNoise[iSector]
1163 << "Noise_RMS_SuperModule=" << smRMSNoise[iSector]
1164 << "Noise_Mean_Chamber=" << chamberMeanNoise[chamberNo]
1165 << "Noise_RMS_Chamber=" << chamberRMSNoise[chamberNo]
1166 << "Vdrift_Mean_Run=" << runMeanVdrift
1167 << "Vdrift_RMS_Run=" << runRMSVdrift
1168 << "Vdrift_Mean_SuperModule=" << smMeanVdrift[iSector]
1169 << "Vdrift_RMS_SuperModule=" << smRMSVdrift[iSector]
1170 << "Vdrift_Mean_Chamber=" << chamberMeanVdrift[chamberNo]
1171 << "Vdrift_RMS_Chamber=" << chamberRMSVdrift[chamberNo]
1172 << "T0_Mean_Run=" << runMeanT0
1173 << "T0_RMS_Run=" << runRMST0
1174 << "T0_Mean_SuperModule=" << smMeanT0[iSector]
1175 << "T0_RMS_SuperModule=" << smRMST0[iSector]
1176 << "T0_Mean_Chamber=" << chamberMeanT0[chamberNo]
1177 << "T0_RMS_Chamber=" << chamberRMST0[chamberNo]
1178 << "Channel.=" << &channelVector
1179 << "Row.=" << &rowVector
1180 << "Column.=" << &colVector
1181 << "PadSuperRow.=" << &padSuperRowVector
1182 << "PadSuperColumn.=" << &padSuperColumnVector
1183 << "Status.=" << &statusVector
1184 << "Gain.=" << &gainVector
1185 << "Noise.=" << &noiseVector
1186 << "Vdrift.=" << &vdriftVector
1187 << "T0.=" << &t0Vector;
1189 (*treeStreamer)<< "TRDcalibDetails"
1190 << "ROB.=" << &robVector
1191 << "MCM.=" << &mcmVector
1192 << "MCMSuperRow.=" << &mcmSuperRowVector
1193 << "MCMSuperColumn.=" << &mcmSuperColumnVector
1194 << "SORandEOR.=" << &sorandeorVector
1195 << "gsmSOR.=" << &gsmSorVector
1196 << "gsmDelta.=" << &gsmDeltaVector
1197 << "nimSOR.=" << &nimSorVector
1198 << "nimDelta.=" << &nimDeltaVector
1199 << "nevSOR.=" << &nevSorVector
1200 << "nevDelta.=" << &nevDeltaVector
1201 << "nptSOR.=" << &nptSorVector
1202 << "nptDelta.=" << &nptDeltaVector;
1203 (*treeStreamer)<< "TRDcalibDetails"
1205 } // end loop over stacks
1206 } // end loop over supermodules
1207 } // end loop over layers
1209 // add the run number to the list of runs
1210 runs.ResizeTo(runs.GetNoElements()+1);
1211 runs[runs.GetNoElements()-1] = currRun;
1214 if(chamberGainFactor) delete chamberGainFactor;
1215 if(padGainFactor) delete padGainFactor;
1216 if(chamberNoise) delete chamberNoise;
1217 if(padNoise) delete padNoise;
1218 if(chamberVdrift) delete chamberVdrift;
1219 if(padVdrift) delete padVdrift;
1220 if(chamberT0) delete chamberT0;
1221 if(padT0) delete padT0;
1222 if(chamberStatus) delete chamberStatus;
1223 if(padStatus) delete padStatus;
1225 // check if we still have run numbers in the file or provided range
1226 if(runListFilename[0]=='\0' && firstRun!=-1 && lastRun!=-1) {
1228 if(currRun>lastRun) break;
1230 if(runListFilename[0]!='\0' && in.eof())
1232 } // end loop over runs
1234 treeStreamer->GetFile()->cd();
1236 delete treeStreamer;
1238 // delete treeStreamerDCS;
1241 //_________________________________________________________________________
1242 void AliTRDCalibViewer::DumpCalibToTree(const Char_t* inFilename, const Char_t* outFilename)
1245 // extract info from CalPad objects and dump them into a tree to be viewed
1247 TTreeSRedirector *treeStreamer = new TTreeSRedirector(outFilename);
1248 //open file and retrieve list of calPad objects
1249 TFile f(inFilename);
1250 TList *l=(TList*)f.GetListOfKeys();
1252 TObjArray arrCalPads;
1253 TObjArray arrSMmean;
1255 arrCalPads.SetOwner();
1256 arrSMmean.SetOwner();
1257 arrSMrms.SetOwner();
1261 while ( (k=(TKey*)next()) ){
1262 AliTRDCalPad *pad=dynamic_cast<AliTRDCalPad*>(k->ReadObj());
1264 arrCalPads.Add(pad);
1266 TVectorD *smMean=new TVectorD(AliTRDcalibDB::kNsector);
1267 TVectorD *smRMS=new TVectorD(AliTRDcalibDB::kNsector);
1269 arrSMmean.Add(smMean);
1270 arrSMrms.Add(smRMS);
1272 ProcessTRDCalibArray(pad, *smMean, *smRMS);
1275 Int_t kSuperModuleStatus[18] = {1, 1, 0, 0, 0, 0, // super module status (1- installed, 0- not installed)
1279 AliTRDgeometry trdGeom;
1280 Int_t kNRows[5] = {16, 16, 12, 16, 16}; // number of pad rows in the chambers from each of the 5 stacks
1282 for(Short_t iLayer=0; iLayer<AliTRDgeometry::kNlayer; iLayer++) { // loop over layers
1283 for(Short_t iSM=0; iSM<AliTRDgeometry::kNsector; iSM++) { // loop over supermodules
1284 if(kSuperModuleStatus[iSM]==0)
1287 for(Short_t iStack=0; iStack<AliTRDgeometry::kNstack; iStack++) { // loop over stacks
1288 AliTRDpadPlane &plane=*trdGeom.GetPadPlane(iLayer, iStack);
1290 Short_t chamberNo = AliTRDgeometry::GetDetector(iLayer, iStack, iSM);
1291 const Int_t nrows=plane.GetNrows();
1292 const Int_t ncols=plane.GetNcols();
1293 const Int_t nchannels=nrows*ncols;
1294 // printf("chamberNo: %d (%03d-%03d-%03d)\n", chamberNo,nrows,ncols,nchannels);
1296 TVectorD channelVector(nchannels);
1297 TVectorD rowVector(nchannels);
1298 TVectorD colVector(nchannels);
1300 TVectorD gxVector(nchannels);
1301 TVectorD gyVector(nchannels);
1302 TVectorD gzVector(nchannels);
1304 TVectorD padSuperRowVector(nchannels);
1305 TVectorD padSuperColumnVector(nchannels);
1308 for(Short_t iRow=0; iRow<nrows; iRow++) { // loop over pad rows
1309 for(Short_t iCol=0; iCol<ncols; iCol++) { // loop over pad columns
1310 Short_t padSuperRow = iRow;
1311 for(Int_t i=0; i<iStack; i++) padSuperRow = padSuperRow + kNRows[i];
1312 padSuperRowVector.GetMatrixArray()[index] = padSuperRow;
1314 Short_t padSuperColumn = iCol;
1315 for(Int_t i=0; i<iSM; i++) padSuperColumn = padSuperColumn + ncols;
1316 padSuperColumnVector.GetMatrixArray()[index] = padSuperColumn;
1318 rowVector.GetMatrixArray()[index] = iRow;
1319 colVector.GetMatrixArray()[index] = iCol;
1322 } // end loop over pad columns
1323 } // end loop over pad rows
1325 (*treeStreamer)<< "TRDcalibDetails"
1326 << "SuperModule=" << iSM
1327 << "Stack=" << iStack
1328 << "Layer=" << iLayer
1329 << "Chamber=" << chamberNo
1330 //geographical information
1331 << "Channel.=" << &channelVector
1332 << "Row.=" << &rowVector
1333 << "Column.=" << &colVector
1334 << "PadSuperRow.=" << &padSuperRowVector
1335 << "PadSuperColumn.=" << &padSuperColumnVector;
1336 // << "gx.=" << &gxVector
1337 // << "gy.=" << &gyVector
1338 // << "gz.=" << &gzVector;
1344 arrTrash.SetOwner();
1345 Int_t ncalib=arrCalPads.GetEntriesFast();
1346 for (Int_t iCalib=0; iCalib<ncalib; ++iCalib){
1347 AliTRDCalPad *pad=(AliTRDCalPad*)arrCalPads.UncheckedAt(iCalib);
1348 AliTRDCalROC *calROC=pad->GetCalROC(chamberNo);
1350 TVectorD &smMean=*((TVectorD*)arrSMmean.UncheckedAt(iCalib));
1351 TVectorD &smRMS=*((TVectorD*)arrSMrms.UncheckedAt(iCalib));
1353 TString calibName=pad->GetName();
1355 TVectorD *valueVector=new TVectorD(nchannels);
1356 arrTrash.Add(valueVector);
1360 Double_t rocMedian=0;
1364 for(Short_t iRow=0; iRow<nrows; iRow++) {
1365 for(Short_t iCol=0; iCol<ncols; iCol++) {
1366 valueVector->GetMatrixArray()[index2] = calROC->GetValue(iCol,iRow);
1370 rocMean = calROC->GetMean();
1371 rocRMS = calROC->GetRMS();
1372 rocMedian = calROC->GetMedian();
1374 if ( !(rocMean<1e30) ) rocMean=0;
1375 if ( !(rocRMS<1e30) ) rocRMS=0;
1377 // printf("mean: %f\n",rocMean);
1378 // printf("rms: %f\n",rocRMS);
1379 // printf("median: %f\n",rocMedian);
1382 (*treeStreamer)<< "TRDcalibDetails"
1383 // statistical information
1384 << (Char_t*)((calibName+"_Mean_SuperModule=").Data()) << smMean[iSM]
1385 << (Char_t*)((calibName+"_RMS_SuperModule=").Data()) << smRMS[iSM]
1386 << (Char_t*)((calibName+"_Mean_Chamber=").Data()) << rocMean
1387 << (Char_t*)((calibName+"_RMS_Chamber=").Data()) << rocRMS
1388 << (Char_t*)((calibName+"_Median_Chamber=").Data()) << rocMedian
1389 //pad by pad information
1390 << (Char_t*)((calibName+".=").Data()) << valueVector;
1392 } // end loop over calib objects
1394 (*treeStreamer)<< "TRDcalibDetails"
1397 } // end loop over stacks
1398 } // end loop over supermodules
1399 } // end loop over layers
1400 delete treeStreamer;
1403 //_____________________________________________________________________________
1404 void AliTRDCalibViewer::ProcessTRDCalibArray(AliTRDCalDet* chamberCalib, AliTRDCalPad *padCalib,
1406 Double_t &runValue, Double_t &runRMS,
1407 TVectorD &chamberValues, TVectorD &chamberValuesRMS,
1408 TVectorD &superModuleValues, TVectorD &superModuleValuesRMS) {
1409 // Process the calibrations for a given run.
1410 // Calculates the run and chamber wise averages
1412 if(!chamberCalib) return;
1413 if(!padCalib) return;
1414 Int_t kSuperModuleStatus[18] = {1, 1, 0, 0, 0, 0, // super module status (1- installed, 0- not installed)
1418 // initialize the histograms used for extracting the mean and RMS
1419 TH1F *runWiseHisto = new TH1F("runHisto", "runHisto", 200, -10.0, 10.0);
1420 TH1F *superModuleWiseHisto = new TH1F("smHisto", "smHisto", 200, -10.0, 10.0);
1421 TH1F *chamberWiseHisto = new TH1F("chamberHisto", "chamberHisto", 200, -10.0, 10.0);
1423 // check if the calibration parameter is multiplicative or additive
1424 Bool_t multiplicative = kTRUE;
1425 if(!parName.CompareTo("T0")) multiplicative = kFALSE;
1427 // first iteration (calculate all averages and RMS without discrimination on the SM average)
1428 runWiseHisto->Reset();
1429 for(Int_t iSM = 0; iSM<AliTRDcalibDB::kNsector; iSM++) { // loop over supermodules
1430 // reset the super module histogram
1431 superModuleWiseHisto->Reset();
1432 // check if SM is installed
1433 if(!kSuperModuleStatus[iSM]) continue;
1434 for(Int_t iChamber=iSM*AliTRDcalibDB::kNstack*AliTRDcalibDB::kNlayer;
1435 iChamber < (iSM+1)*AliTRDcalibDB::kNstack*AliTRDcalibDB::kNlayer;
1436 iChamber++) { // loop over chambers in this supermodule
1437 // get the chamber value
1438 Float_t chamberValue = chamberCalib->GetValue(iChamber);
1439 // get the ROC object
1440 AliTRDCalROC *chamberROC = padCalib->GetCalROC(iChamber);
1443 chamberWiseHisto->Reset();
1444 for(Int_t iChannel = 0; iChannel < chamberROC->GetNchannels(); iChannel++){ // loop over channels
1445 // calculate the calibration parameter for this pad
1448 padValue = chamberValue * chamberROC->GetValue(iChannel);
1450 padValue = chamberValue + chamberROC->GetValue(iChannel);
1451 // fill the run, SM and chamber wise histograms
1452 chamberWiseHisto->Fill(padValue);
1453 // if the parameter is Noise then check if the pad value is not a default one
1454 // Default value is now 1.2!!!! Check with Raphaelle for more informations
1455 if(parName.Contains("Noise") &&
1456 TMath::Abs(padValue-1.2)<0.00001) continue;
1457 superModuleWiseHisto->Fill(padValue);
1458 runWiseHisto->Fill(padValue);
1459 } // end loop over channels
1460 // get the chamber wise mean and RMS
1461 chamberValues[iChamber] = chamberWiseHisto->GetMean();
1462 chamberValuesRMS[iChamber] = chamberWiseHisto->GetRMS();
1463 } // end loop over chambers
1464 // SM wise mean and RMS
1465 superModuleValues[iSM] = superModuleWiseHisto->GetMean();
1466 superModuleValuesRMS[iSM] = superModuleWiseHisto->GetRMS();
1467 } // end loop over supermodules
1468 // run wise mean and RMS
1469 runValue = runWiseHisto->GetMean();
1470 runRMS = runWiseHisto->GetRMS();
1472 // Noise and Gain are finished processing
1473 if(parName.Contains("Noise") || parName.Contains("Gain"))
1475 // second iteration (calculate SM and run wise averages and RMS for Vdrift and T0)
1476 // The pads with calib parameter equal to the SM average are discarded (default value)
1477 runWiseHisto->Reset();
1478 for(Int_t iSM = 0; iSM<AliTRDcalibDB::kNsector; iSM++) { // loop over supermodules
1479 superModuleWiseHisto->Reset();
1480 // eliminate the uninstalled super modules
1481 if(!kSuperModuleStatus[iSM]) continue;
1482 for(Int_t iChamber=iSM*AliTRDcalibDB::kNstack*AliTRDcalibDB::kNlayer;
1483 iChamber < (iSM+1)*AliTRDcalibDB::kNstack*AliTRDcalibDB::kNlayer;
1484 iChamber++) { // loop over chambers
1485 // the chamber value
1486 Float_t chamberValue = chamberCalib->GetValue(iChamber);
1487 AliTRDCalROC *chamberROC = padCalib->GetCalROC(iChamber);
1491 for(Int_t iChannel = 0; iChannel < chamberROC->GetNchannels(); iChannel++){ // loop over channels in a chamber
1492 // get the pad value
1495 padValue = chamberValue * chamberROC->GetValue(iChannel);
1497 padValue = chamberValue + chamberROC->GetValue(iChannel);
1498 // eliminate from the average and RMS calculation all pads which
1499 // have the calib parameter equal with the SM average
1500 if((parName.Contains("Vdrift") || parName.Contains("T0")) &&
1501 TMath::Abs(padValue-superModuleValues[iSM])<0.00001) continue;
1502 superModuleWiseHisto->Fill(padValue);
1503 runWiseHisto->Fill(padValue);
1504 } // end loop over channels
1505 } // end loop over chambers
1507 superModuleValues[iSM] = superModuleWiseHisto->GetMean();
1508 superModuleValuesRMS[iSM] = superModuleWiseHisto->GetRMS();
1509 } // end loop over super modules
1510 runValue = runWiseHisto->GetMean();
1511 runRMS = runWiseHisto->GetRMS();
1513 delete chamberWiseHisto;
1514 delete superModuleWiseHisto;
1515 delete runWiseHisto;
1521 //_____________________________________________________________________________
1522 void AliTRDCalibViewer::ProcessTRDCalibArray(AliTRDCalPad *padCalib,
1523 TVectorD &superModuleValues, TVectorD &superModuleValuesRMS)
1525 // Process the calibrations for a given run.
1526 // Calculates the run and chamber wise averages
1528 if(!padCalib) return;
1529 Int_t kSuperModuleStatus[18] = {1, 1, 0, 0, 0, 0, // super module status (1- installed, 0- not installed)
1533 for(Int_t iSM = 0; iSM<AliTRDcalibDB::kNsector; iSM++) { // loop over supermodules
1537 if(!kSuperModuleStatus[iSM]) continue;
1539 // loop over chambers in this supermodule
1540 const Int_t nchambers=AliTRDcalibDB::kNstack*AliTRDcalibDB::kNlayer;
1541 for(Int_t iChamber=iSM*nchambers; iChamber < (iSM+1)*nchambers; iChamber++) {
1542 AliTRDCalROC *chamberROC = padCalib->GetCalROC(iChamber);
1546 // loop over channels
1547 for(Int_t iChannel = 0; iChannel < chamberROC->GetNchannels(); iChannel++){
1548 mean+=chamberROC->GetValue(iChannel);
1549 rms+=chamberROC->GetValue(iChannel)*chamberROC->GetValue(iChannel);
1551 } // end loop over channels
1552 } // end loop over chambers
1554 //calculate mean and rms
1558 rms=TMath::Sqrt(TMath::Abs(mean*mean-rms));
1560 // SM wise mean and RMS
1561 superModuleValues[iSM] = mean;
1562 superModuleValuesRMS[iSM] = rms;
1563 } // end loop over supermodules