]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCCalibViewer.cxx
Using AliGeomManager in the macros (Raffaele)
[u/mrichter/AliRoot.git] / TPC / AliTPCCalibViewer.cxx
CommitLineData
39bcd65d 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16
17///////////////////////////////////////////////////////////////////////////////
18// //
19// Class for viewing/visualizing TPC calibration data //
20// base on TTree functionality for visualization //
21///////////////////////////////////////////////////////////////////////////////
22
23//
24// ROOT includes
25//
26#include <iostream>
27#include <TString.h>
28#include <TRandom.h>
29#include <TLegend.h>
30#include <TLine.h>
31#include <TCanvas.h>
32#include <TROOT.h>
33#include <TStyle.h>
34#include <TH1F.h>
35#include <THashTable.h>
36#include <TObjString.h>
37#include "TTreeStream.h"
38#include "TFile.h"
39#include "TKey.h"
40
41
42//
43// AliRoot includes
44//
45#include "AliTPCCalibViewer.h"
46
47ClassImp(AliTPCCalibViewer)
48
49AliTPCCalibViewer::AliTPCCalibViewer()
50 :TObject(),
51 fTree(0),
52 fFile(0),
a6d2bd0c 53 fListOfObjectsToBeDeleted(0),
54 fTreeMustBeDeleted(0)
39bcd65d 55{
56 //
57 // Default constructor
58 //
59
60}
61
62//_____________________________________________________________________________
63AliTPCCalibViewer::AliTPCCalibViewer(const AliTPCCalibViewer &c)
64 :TObject(c),
65 fTree(0),
66 fFile(0),
a6d2bd0c 67 fListOfObjectsToBeDeleted(0),
68 fTreeMustBeDeleted(0)
39bcd65d 69{
70 //
71 // dummy AliTPCCalibViewer copy constructor
72 // not yet working!!!
73 //
74 fTree = c.fTree;
a6d2bd0c 75 fTreeMustBeDeleted = c.fTreeMustBeDeleted;
39bcd65d 76 //fFile = new TFile(*(c.fFile));
77 fListOfObjectsToBeDeleted = c.fListOfObjectsToBeDeleted;
78}
79
80//_____________________________________________________________________________
81AliTPCCalibViewer::AliTPCCalibViewer(TTree* tree)
82 :TObject(),
83 fTree(0),
84 fFile(0),
a6d2bd0c 85 fListOfObjectsToBeDeleted(0),
86 fTreeMustBeDeleted(0)
39bcd65d 87{
88 //
89 // Constructor that initializes the calibration viewer
90 //
91 fTree = tree;
a6d2bd0c 92 fTreeMustBeDeleted = kFALSE;
39bcd65d 93 fListOfObjectsToBeDeleted = new TObjArray();
94}
95
96//_____________________________________________________________________________
97AliTPCCalibViewer::AliTPCCalibViewer(char* fileName, char* treeName)
98 :TObject(),
99 fTree(0),
100 fFile(0),
a6d2bd0c 101 fListOfObjectsToBeDeleted(0),
102 fTreeMustBeDeleted(0)
39bcd65d 103{
104 //
105 // Constructor to initialize the calibration viewer
106 // the file 'fileName' contains the tree 'treeName'
107 //
108 fFile = new TFile(fileName, "read");
109 fTree = (TTree*) fFile->Get(treeName);
a6d2bd0c 110 fTreeMustBeDeleted = kTRUE;
39bcd65d 111 fListOfObjectsToBeDeleted = new TObjArray();
112}
113
114//____________________________________________________________________________
115AliTPCCalibViewer & AliTPCCalibViewer::operator =(const AliTPCCalibViewer & param)
116{
117 //
118 // assignment operator - dummy
119 // not yet working!!!
120 //
121 fTree = param.fTree;
a6d2bd0c 122 fTreeMustBeDeleted = param.fTreeMustBeDeleted;
39bcd65d 123 //fFile = new TFile(*(param.fFile));
124 fListOfObjectsToBeDeleted = param.fListOfObjectsToBeDeleted;
125 return (*this);
126}
127
128//_____________________________________________________________________________
129AliTPCCalibViewer::~AliTPCCalibViewer()
130{
131 //
132 // AliTPCCalibViewer destructor
a6d2bd0c 133 // all objects will be deleted, the file will be closed, the pictures will disappear
39bcd65d 134 //
a6d2bd0c 135 if (fTree && fTreeMustBeDeleted) {
136 fTree->SetCacheSize(0);
137 fTree->Delete();
138 //delete fTree;
139 }
39bcd65d 140 if (fFile) {
141 fFile->Close();
142 fFile = 0;
143 }
144
145 for (Int_t i = fListOfObjectsToBeDeleted->GetEntriesFast()-1; i >= 0; i--) {
146 //cout << "Index " << i << " trying to delete the following object: " << fListOfObjectsToBeDeleted->At(i)->GetName() << "..."<< endl;
147 delete fListOfObjectsToBeDeleted->At(i);
148 }
149 delete fListOfObjectsToBeDeleted;
150}
151
a6d2bd0c 152//_____________________________________________________________________________
153void AliTPCCalibViewer::Delete(Option_t* option) {
154 //
155 // Should be called from AliTPCCalibViewerGUI class only.
156 // If you use Delete() do not call the destructor.
157 // All objects (except those contained in fListOfObjectsToBeDeleted) will be deleted, the file will be closed.
158 //
159
160 option = option; // to avoid warnings on compiling
161 if (fTree && fTreeMustBeDeleted) {
162 fTree->SetCacheSize(0);
163 fTree->Delete();
164 }
165 if (fFile)
166 delete fFile;
167 delete fListOfObjectsToBeDeleted;
168}
169
39bcd65d 170//_____________________________________________________________________________
171Int_t AliTPCCalibViewer::EasyDraw(const char* drawCommand, const char* sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
172 //
173 // easy drawing of data, use '~' for abbreviation of '.fElements'
174 // example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
175 // sector: sector-number - only the specified sector will be drwawn
176 // 'A'/'C' or 'a'/'c' - side A/C will be drawn
177 // 'ALL' - whole TPC will be drawn, projected on one side
178 // cuts: specifies cuts
179 // drawOptions: draw options like 'same'
180 // writeDrawCommand: write the command, that is passed to TTree::Draw
181 //
182 TString drawStr(drawCommand);
183 TString sectorStr(sector);
184 sectorStr.ToUpper();
185 TString cutStr("");
186 TString drawOptionsStr("profcolz ");
187 TRandom rnd(0);
188 Int_t rndNumber = rnd.Integer(10000);
189 if (drawOptions && drawOptions != "")
190 drawOptionsStr += drawOptions;
191
192 if (sectorStr == "A") {
193 drawStr += ":gy.fElements:gx.fElements>>prof";
194 drawStr += rndNumber;
195 drawStr += "(330,-250,250,330,-250,250)";
196 cutStr += "(sector/18)%2==0 ";
197 }
198 else if (sectorStr == "C") {
199 drawStr += ":gy.fElements:gx.fElements>>prof";
200 drawStr += rndNumber;
201 drawStr += "(330,-250,250,330,-250,250)";
202 cutStr += "(sector/18)%2==1 ";
203 }
204 else if (sectorStr == "ALL") {
205 drawStr += ":gy.fElements:gx.fElements>>prof";
206 drawStr += rndNumber;
207 drawStr += "(330,-250,250,330,-250,250)";
208 }
209 else if (sectorStr.IsDigit()) {
210 Int_t isec = sectorStr.Atoi();
211 drawStr += ":rpad.fElements:row.fElements>>prof";
212 drawStr += rndNumber;
213 if (isec < 36 && isec >= 0)
214 drawStr += "(63,0,63,108,-54,54)";
215 else if (isec < 72 && isec >= 36)
216 drawStr += "(96,0,96,140,-70,70)";
217 else {
218 Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
219 return -1;
220 }
221 cutStr += "(sector==";
222 cutStr += isec;
223 cutStr += ") ";
224 }
225
226 if (cuts && cuts[0] != 0) {
227 if (cutStr.Length() != 0) cutStr += "&& ";
228 cutStr += "(";
229 cutStr += cuts;
230 cutStr += ")";
231 }
232 drawStr.ReplaceAll("~", ".fElements");
233 cutStr.ReplaceAll("~", ".fElements");
234 if (writeDrawCommand) std::cout << "fTree->Draw(\"" << drawStr << "\", \"" << cutStr << "\", \"" << drawOptionsStr << "\");" << std::endl;
235 return fTree->Draw(drawStr.Data(), cutStr.Data(), drawOptionsStr.Data());
236}
237
238Int_t AliTPCCalibViewer::EasyDraw(const char* drawCommand, Int_t sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
239 //
240 // easy drawing of data, use '~' for abbreviation of '.fElements'
241 // example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0")
242 // sector: sector-number - only the specified sector will be drwawn
243 // cuts: specifies cuts
244 // drawOptions: draw options like 'same'
245 // writeDrawCommand: write the command, that is passed to TTree::Draw
246 //
247 if (sector >= 0 && sector < 72) {
248 char sectorChr[3];
249 sprintf(sectorChr, "%i", sector);
250 return EasyDraw(drawCommand, sectorChr, cuts, drawOptions, writeDrawCommand);
251 }
252 Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
253 return -1;
254}
255
256//_____________________________________________________________________________
257Int_t AliTPCCalibViewer::EasyDraw1D(const char* drawCommand, const char* sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
258 //
259 // easy drawing of data, use '~' for abbreviation of '.fElements'
260 // example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
261 // sector: sector-number - the specified sector will be drwawn
262 // 'A'/'C' or 'a'/'c' - side A/C will be drawn
263 // 'ALL' - whole TPC will be drawn, projected on one side
264 // cuts: specifies cuts
265 // drawOptions: draw options like 'same'
266 // writeDrawCommand: write the command, that is passed to TTree::Draw
267 //
268
269 TString drawStr(drawCommand);
270 TString sectorStr(sector);
271 TString drawOptionsStr(drawOptions);
272 sectorStr.ToUpper();
273 TString cutStr("");
274
275 if (sectorStr == "A")
276 cutStr += "(sector/18)%2==0 ";
277 else if (sectorStr == "C")
278 cutStr += "(sector/18)%2==1 ";
279 else if (sectorStr.IsDigit()) {
280 Int_t isec = sectorStr.Atoi();
281 if (isec < 0 || isec > 71) {
282 Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
283 return -1;
284 }
285 cutStr += "(sector==";
286 cutStr += isec;
287 cutStr += ") ";
288 }
289
290 if (cuts && cuts[0] != 0) {
291 if (cutStr.Length() != 0) cutStr += "&& ";
292 cutStr += "(";
293 cutStr += cuts;
294 cutStr += ")";
295 }
296
297 drawStr.ReplaceAll("~", ".fElements");
298 cutStr.ReplaceAll("~", ".fElements");
299 if (writeDrawCommand) std::cout << "fTree->Draw(\"" << drawStr << "\", \"" << cutStr << "\", \"" << drawOptionsStr << "\");" << std::endl;
300 return fTree->Draw(drawStr.Data(), cutStr.Data(), drawOptionsStr.Data());
301}
302
303Int_t AliTPCCalibViewer::EasyDraw1D(const char* drawCommand, Int_t sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
304 //
305 // easy drawing of data, use '~' for abbreviation of '.fElements'
306 // example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0")
307 // sector: sector-number - the specified sector will be drwawn
308 // cuts: specifies cuts
309 // drawOptions: draw options like 'same'
310 // writeDrawCommand: write the command, that is passed to TTree::Draw
311 //
312
313 if (sector >= 0 && sector < 72) {
314 char sectorChr[3];
315 sprintf(sectorChr, "%i", sector);
316 return EasyDraw1D(drawCommand, sectorChr, cuts, drawOptions, writeDrawCommand);
317 }
318 Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
319 return -1;
320}
321
322//_____________________________________________________________________________
323Int_t AliTPCCalibViewer::DrawHisto1D(const char* type, Int_t sector, TVectorF& nsigma, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM) const {
324 //
325 // draws a 1-dimensional histogram of 'type' for sector 'sector'
326 // TVectorF nsigma: Specifies, for which distances from the mean/median/LTM lines should be drawn, in units of sigma
327 // example: nsigma={2, 4, 6}: Three lines will be drawn, distance to mean/median/LTM: 2, 3 and 6 sigma
328 // plotMean, plotMedian, plotLTM: specifies, if mean, median and LTM should be drawn as lines into the histogram
329 //
330
331 TString typeStr(type);
332 TString sectorStr("sector==");
333 sectorStr += sector;
334
335 TCanvas* canvas = ((TCanvas*)gROOT->GetListOfCanvases()->Last());
336 Int_t oldOptStat = gStyle->GetOptStat();
337 gStyle->SetOptStat(0000000);
338
339 if (!canvas) {
340 canvas = new TCanvas();
341 fListOfObjectsToBeDeleted->Add(canvas);
342 }
343
344 char c[500];
345 sprintf(c, "%s, sector: %i", type, sector);
346 TLegend * legend = new TLegend(.8,.6, .99, .99, c);
347 fListOfObjectsToBeDeleted->Add(legend);
348
349 Int_t nentries = fTree->Draw((typeStr+".fElements").Data(), sectorStr.Data(), "");
350 ((TH1F*)canvas->GetPrimitive("htemp"))->SetTitle("");
351
352 //****************************************************************
353 //!!!!!!!!!!!!!!!!! Needs further investigaton !!!!!!!!!!!!!!!!!!!
354 //****************************************************************
355 //fListOfObjectsToBeDeleted->Add(canvas->GetPrimitive("htemp"));
356/*
357 By default the temporary histogram created is called "htemp", but only in
358 the one dimensional Draw("e1") it contains the TTree's data points. For
359 a two dimensional Draw, the data is filled into a TGraph which is named
360 "Graph". They can be retrieved by calling
361 TH1F *htemp = (TH1F*)gPad->GetPrimitive("htemp"); // 1D
362 TGraph *graph = (TGraph*)gPad->GetPrimitive("Graph"); // 2D
363*/
364
365 canvas->Update();
366 Double_t sigma = 0;
367
368 if (plotMean) {
369 fTree->Draw((typeStr+"_Mean").Data(), sectorStr.Data(), "goff");
370 Double_t lineX = fTree->GetV1()[0];
371 fTree->Draw((typeStr+"_RMS").Data(), sectorStr.Data(), "goff");
372 sigma = fTree->GetV1()[0];
373 TLine* line = new TLine(lineX, 0, lineX, canvas->GetUymax());
374 fListOfObjectsToBeDeleted->Add(line);
375 line->SetLineColor(kRed);
376 line->SetLineWidth(2);
377 line->SetLineStyle(1);
378 line->Draw();
379 sprintf(c, "Mean: %f", lineX);
380 legend->AddEntry(line, c, "l");
381
382 for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
383 TLine* linePlusSigma = new TLine(lineX+nsigma[i]*sigma, 0, lineX+nsigma[i]*sigma, canvas->GetUymax());
384 fListOfObjectsToBeDeleted->Add(linePlusSigma);
385 linePlusSigma->SetLineColor(kRed);
386 linePlusSigma->SetLineStyle(2+i);
387 linePlusSigma->Draw();
388
389 TLine* lineMinusSigma = new TLine(lineX-nsigma[i]*sigma, 0, lineX-nsigma[i]*sigma, canvas->GetUymax());
390 fListOfObjectsToBeDeleted->Add(lineMinusSigma);
391 lineMinusSigma->SetLineColor(kRed);
392 lineMinusSigma->SetLineStyle(2+i);
393 lineMinusSigma->Draw();
394 sprintf(c, "%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i]*sigma));
395 std::cout << "nsigma-char*: " << c << std::endl;
396 legend->AddEntry(lineMinusSigma, c, "l");
397 }
398 }
399
400 if (plotMedian) {
401 fTree->Draw((typeStr+"_Median").Data(), sectorStr.Data(), "goff");
402 Double_t lineX = fTree->GetV1()[0];
403 fTree->Draw((typeStr+"_RMS").Data(), sectorStr.Data(), "goff");
404 sigma = fTree->GetV1()[0];
405 TLine* line = new TLine(lineX, 0, lineX, canvas->GetUymax());
406 fListOfObjectsToBeDeleted->Add(line);
407 line->SetLineColor(kBlue);
408 line->SetLineWidth(2);
409 line->SetLineStyle(1);
410 line->Draw();
411 sprintf(c, "Median: %f", lineX);
412 legend->AddEntry(line, c, "l");
413
414 for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
415 TLine* linePlusSigma = new TLine(lineX+nsigma[i]*sigma, 0, lineX+nsigma[i]*sigma, canvas->GetUymax());
416 fListOfObjectsToBeDeleted->Add(linePlusSigma);
417 linePlusSigma->SetLineColor(kBlue);
418 linePlusSigma->SetLineStyle(2+i);
419 linePlusSigma->Draw();
420
421 TLine* lineMinusSigma = new TLine(lineX-nsigma[i]*sigma, 0, lineX-nsigma[i]*sigma, canvas->GetUymax());
422 fListOfObjectsToBeDeleted->Add(lineMinusSigma);
423 lineMinusSigma->SetLineColor(kBlue);
424 lineMinusSigma->SetLineStyle(2+i);
425 lineMinusSigma->Draw();
426 sprintf(c, "%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i]*sigma));
427 legend->AddEntry(lineMinusSigma, c, "l");
428 }
429 }
430
431 if (plotLTM) {
432 fTree->Draw((typeStr+"_LTM").Data(), sectorStr.Data(), "goff");
433 Double_t lineX = fTree->GetV1()[0];
434 fTree->Draw((typeStr+"_RMS_LTM").Data(), sectorStr.Data(), "goff");
435 sigma = fTree->GetV1()[0];
436 TLine* line = new TLine(lineX, 0, lineX, canvas->GetUymax());
437 fListOfObjectsToBeDeleted->Add(line);
438 line->SetLineColor(kGreen+2);
439 line->SetLineWidth(2);
440 line->SetLineStyle(1);
441 line->Draw();
442 sprintf(c, "LTM: %f", lineX);
443 legend->AddEntry(line, c, "l");
444
445 for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
446 TLine* linePlusSigma = new TLine(lineX+nsigma[i]*sigma, 0, lineX+nsigma[i]*sigma, canvas->GetUymax());
447 fListOfObjectsToBeDeleted->Add(linePlusSigma);
448 linePlusSigma->SetLineColor(kGreen+2);
449 linePlusSigma->SetLineStyle(2+i);
450 linePlusSigma->Draw();
451
452 TLine* lineMinusSigma = new TLine(lineX-nsigma[i]*sigma, 0, lineX-nsigma[i]*sigma, canvas->GetUymax());
453 fListOfObjectsToBeDeleted->Add(lineMinusSigma);
454 lineMinusSigma->SetLineColor(kGreen+2);
455 lineMinusSigma->SetLineStyle(2+i);
456 lineMinusSigma->Draw();
457 sprintf(c, "%i #sigma = %f", (Int_t)(nsigma[i]), (Float_t)(nsigma[i]*sigma));
458 legend->AddEntry(lineMinusSigma, c, "l");
459 }
460 }
461 legend->Draw();
462 gStyle->SetOptStat(oldOptStat);
463 return nentries;
464}
465
466//_____________________________________________________________________________
467void AliTPCCalibViewer::SigmaCut(const char* type, Int_t sector, Float_t sigmaMax, Float_t sigmaStep, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM) const {
468 //
469 // Creates a histogram, where you can see, how much of the data are inside sigma-intervals around the mean/median/LTM
470 // type: For which type of data the histogram is generated, e.g. 'CEQmean'
471 // sector: For which sector the histogram is generated
472 // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated
473 // sigmaStep: the binsize of the generated histogram
474 // plotMean/plotMedian/plotLTM: specifies where to put the center
475 //
476 Int_t oldOptStat = gStyle->GetOptStat();
477 gStyle->SetOptStat(0000000);
478
479 TString typeStr(type);
480 TString sectorStr("sector==");
481 sectorStr += sector;
482 Int_t entries = fTree->Draw((typeStr+".fElements").Data(), sectorStr.Data(), "goff");
483 char headline[500];
484 sprintf(headline, "%s in sector %i; Multiples of #sigma; Fraction of used pads", type, sector);
485 TH1F *histMean = new TH1F("histMean",headline, (Int_t)(sigmaMax/sigmaStep+1), 0, sigmaMax+sigmaStep);
486 sprintf(headline, "%s in sector %i; Multiples of #sigma; Fraction of used pads", type, sector);
487 TH1F *histMedian = new TH1F("histMedian",headline, (Int_t)(sigmaMax/sigmaStep+1), 0, sigmaMax+sigmaStep);
488 sprintf(headline, "%s in sector %i; Multiples of #sigma; Fraction of used pads", type, sector);
489 TH1F *histLTM = new TH1F("histLTM",headline, (Int_t)(sigmaMax/sigmaStep+1), 0, sigmaMax+sigmaStep);
490 histMean->SetDirectory(0);
491 histMedian->SetDirectory(0);
492 histLTM->SetDirectory(0);
493 fListOfObjectsToBeDeleted->Add(histMean);
494 fListOfObjectsToBeDeleted->Add(histMedian);
495 fListOfObjectsToBeDeleted->Add(histLTM);
496
497
498 // example-cut: sector==34 && TMath::Abs(CEQmean.fElements - CEQmean_Mean) < nsigma * CEQmean_RMS
499 for (Float_t nsigma = 0; nsigma <= sigmaMax; nsigma += sigmaStep) {
500 std::cout << "Calculating histograms, step: " << (Int_t)(nsigma/sigmaStep) << " of: " << (Int_t)(sigmaMax/sigmaStep) << "\r" << std::flush;
501 char cuts[5000];
502
503 if (plotMean) {
504 sprintf(cuts, "sector==%i && ( %s.fElements - %s_Median) < %f * %s_RMS", sector, type, type, nsigma, type );
505 sprintf(cuts, "%s && (-%s.fElements + %s_Median) < %f * %s_RMS", cuts, type, type, nsigma, type );
506 Float_t value = (Float_t)fTree->Draw((typeStr+".fElements").Data(), cuts, "goff")/entries;
507 histMean->Fill(nsigma, value);
508 }
509 if (plotMedian) {
510 sprintf(cuts, "sector==%i && ( %s.fElements - %s_Mean) < %f * %s_RMS", sector, type, type, nsigma, type );
511 sprintf(cuts, "%s && (-%s.fElements + %s_Mean) < %f * %s_RMS", cuts, type, type, nsigma, type );
512 Float_t value = (Float_t)fTree->Draw((typeStr+".fElements").Data(), cuts, "goff")/entries;
513 histMedian->Fill(nsigma, value);
514 }
515 if (plotLTM) {
516 sprintf(cuts, "sector==%i && ( %s.fElements - %s_LTM) < %f * %s_RMS_LTM", sector, type, type, nsigma, type );
517 sprintf(cuts, "%s && (-%s.fElements + %s_LTM) < %f * %s_RMS_LTM", cuts, type, type, nsigma, type );
518 Float_t value = (Float_t)fTree->Draw((typeStr+".fElements").Data(), cuts, "goff")/entries;
519 histLTM->Fill(nsigma, value);
520 }
521 }
522
523 char c[500];
524 sprintf(c, "Sigma Cut");
525 TLegend * legend = new TLegend(.85,.8, .99, .99, c);
526 fListOfObjectsToBeDeleted->Add(legend);
527
528 if (plotMean){
529 histMean->SetLineColor(kBlack);
530 sprintf(c, "Mean");
531 legend->AddEntry(histMean, c, "l");
532 histMean->Draw();
533 }
534 if (plotMedian){
535 histMedian->SetLineColor(kRed);
536 sprintf(c, "Median");
537 legend->AddEntry(histMedian, c, "l");
538 histMedian->Draw("same");
539 }
540 if (plotLTM){
541 histLTM->SetLineColor(kBlue);
542 sprintf(c, "LTM");
543 legend->AddEntry(histLTM, c, "l");
544 histLTM->Draw("same");
545 }
546
547 legend->Draw();
548 gStyle->SetOptStat(oldOptStat);
549}
550
551
552//_____________________________________________________________________________
553AliTPCCalPad* AliTPCCalibViewer::GetCalPad(const char* desiredData, char* cuts, char* calPadName) const {
554 //
555 // creates a AliTPCCalPad out of the 'desiredData'
556 // the functionality of EasyDraw1D is used
557 // calPadName specifies the name of the created AliTPCCalPad
558 // - this takes a while -
559 //
560 TString drawStr(desiredData);
561 drawStr.Append(":channel~");
562 AliTPCCalPad * createdCalPad = new AliTPCCalPad(calPadName, calPadName);
563 Int_t entries = 0;
564 for (Int_t sec = 0; sec < 72; sec++) {
565 entries = EasyDraw1D(drawStr.Data(), (Int_t)sec, cuts, "goff");
a6d2bd0c 566 if (entries == -1) return 0;
39bcd65d 567 for (Int_t i = 0; i < entries; i++)
568 createdCalPad->GetCalROC(sec)->SetValue((UInt_t)(fTree->GetV2()[i]), (Float_t)(fTree->GetV1()[i]));
569 }
570 return createdCalPad;
571}
572
573//_____________________________________________________________________________
574AliTPCCalROC* AliTPCCalibViewer::GetCalROC(const char* desiredData, UInt_t sector, char* cuts) const {
575 //
576 // creates a AliTPCCalROC out of the desiredData
577 // the functionality of EasyDraw1D is used
578 // sector specifies the sector of the created AliTPCCalROC
579 //
580 TString drawStr(desiredData);
581 drawStr.Append(":channel~");
582 Int_t entries = EasyDraw1D(drawStr.Data(), (Int_t)sector, cuts, "goff");
a6d2bd0c 583 if (entries == -1) return 0;
39bcd65d 584 AliTPCCalROC * createdROC = new AliTPCCalROC(sector);
585 for (Int_t i = 0; i < entries; i++)
586 createdROC->SetValue((UInt_t)(fTree->GetV2()[i]), fTree->GetV1()[i]);
587 return createdROC;
588}
589
590
591TObjArray* AliTPCCalibViewer::GetListOfVariables(Bool_t printList) {
592 //
593 // scan the tree - produces a list of available variables in the tree
594 // printList: print the list to the screen, after the scan is done
595 //
596 TObjArray* arr = new TObjArray();
597 TObjString* str = 0;
598 Int_t nentries = fTree->GetListOfBranches()->GetEntries();
599 for (Int_t i = 0; i < nentries; i++) {
600 str = new TObjString(fTree->GetListOfBranches()->At(i)->GetName());
601 str->String().ReplaceAll("_Median", "");
602 str->String().ReplaceAll("_Mean", "");
603 str->String().ReplaceAll("_RMS", "");
604 str->String().ReplaceAll("_LTM", "");
605 str->String().ReplaceAll("_OutlierCutted", "");
606 str->String().ReplaceAll(".", "");
607 if (!arr->FindObject(str) &&
608 !(str->String() == "channel" || str->String() == "gx" || str->String() == "gy" ||
609 str->String() == "lx" || str->String() == "ly" || str->String() == "pad" ||
610 str->String() == "row" || str->String() == "rpad" || str->String() == "sector" ))
611 arr->Add(str);
612 }
613 arr->Sort();
614
615 if (printList) {
616 TIterator* iter = arr->MakeIterator();
617 iter->Reset();
618 TObjString* currentStr = 0;
619 while ( (currentStr = (TObjString*)(iter->Next())) ) {
620 std::cout << currentStr->GetString().Data() << std::endl;
621 }
622 delete iter;
623 }
624 return arr;
625}
626
627
628TObjArray* AliTPCCalibViewer::GetListOfNormalizationVariables(Bool_t printList) {
629 //
630 // produces a list of available variables for normalization in the tree
631 // printList: print the list to the screen, after the scan is done
632 //
633 TObjArray* arr = new TObjArray();
634 arr->Add(new TObjString("_Mean"));
635 arr->Add(new TObjString("_Mean_OutlierCutted"));
636 arr->Add(new TObjString("_Median"));
637 arr->Add(new TObjString("_Median_OutlierCutted"));
638 arr->Add(new TObjString("_LTM"));
639 arr->Add(new TObjString("_LTM_OutlierCutted"));
640 arr->Add(new TObjString("LFitIntern_4_8.fElements"));
641 arr->Add(new TObjString("GFitIntern_Lin.fElements"));
642 arr->Add(new TObjString("GFitIntern_Par.fElements"));
a6d2bd0c 643 arr->Add(new TObjString("FitLinLocal"));
644 arr->Add(new TObjString("FitLinGlobal"));
39bcd65d 645
646 if (printList) {
647 TIterator* iter = arr->MakeIterator();
648 iter->Reset();
649 TObjString* currentStr = 0;
650 while ((currentStr = (TObjString*)(iter->Next()))) {
651 std::cout << currentStr->GetString().Data() << std::endl;
652 }
653 delete iter;
654 }
655 return arr;
656}
657
658
659TFriendElement* AliTPCCalibViewer::AddReferenceTree(const char* filename, const char* treename, const char* refname){
660 //
661 // add a reference tree to the current tree
662 // by default the treename is 'calPads' and the reference treename is 'R'
663 //
664 TFile *file = new TFile(filename);
665 fListOfObjectsToBeDeleted->Add(file);
666 TTree * tree = (TTree*)file->Get(treename);
667 return AddFriend(tree, refname);
668}
669
670
671TObjArray* AliTPCCalibViewer::GetArrayOfCalPads(){
672 //
673 // Returns a TObjArray with all AliTPCCalPads that are stored in the tree
674 // - this takes a while -
675 //
676 TObjArray *listOfCalPads = GetListOfVariables();
677 TObjArray *calPadsArray = new TObjArray();
678 Int_t numberOfCalPads = listOfCalPads->GetEntries();
679 for (Int_t i = 0; i < numberOfCalPads; i++) {
680 std::cout << "Creating calPad " << (i+1) << " of " << numberOfCalPads << "\r" << std::flush;
681 char* calPadName = (char*)((TObjString*)(listOfCalPads->At(i)))->GetString().Data();
682 TString drawCommand = ((TObjString*)(listOfCalPads->At(i)))->GetString();
683 drawCommand.Append("~");
684 AliTPCCalPad* calPad = GetCalPad(drawCommand.Data(), "", calPadName);
685 calPadsArray->Add(calPad);
686 }
687 std::cout << std::endl;
688 listOfCalPads->Delete();
689 delete listOfCalPads;
690 return calPadsArray;
691}
692
693
a6d2bd0c 694TString* AliTPCCalibViewer::Fit(const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, TVectorD &fitParam, TMatrixD &covMatrix){
695 //
696 // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
697 // returns chi2, fitParam and covMatrix
698 // returns TString with fitted formula
699 //
700
701 TString formulaStr(formula);
702 TString drawStr(drawCommand);
703 TString cutStr(cuts);
704
705 // abbreviations:
706 drawStr.ReplaceAll("~",".fElements");
707 cutStr.ReplaceAll("~",".fElements");
708 formulaStr.ReplaceAll("~", ".fElements");
709
710 formulaStr.ReplaceAll("++", "~");
711 TObjArray* formulaTokens = formulaStr.Tokenize("~");
712 Int_t dim = formulaTokens->GetEntriesFast();
713
714 fitParam.ResizeTo(dim);
715 covMatrix.ResizeTo(dim,dim);
716
717 TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
718 fitter->StoreData(kTRUE);
719 fitter->ClearPoints();
720
721 Int_t entries = Draw(drawStr.Data(), cutStr.Data(), "goff");
722 if (entries == -1) return new TString("An ERROR has occured during fitting!");
723 Double_t **values = new Double_t*[dim+1] ;
724
725 for (Int_t i = 0; i < dim + 1; i++){
726 Int_t centries = 0;
727 if (i < dim) centries = fTree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff");
728 else centries = fTree->Draw(drawStr.Data(), cutStr.Data(), "goff");
729
730 if (entries != centries) return new TString("An ERROR has occured during fitting!");
731 values[i] = new Double_t[entries];
732 memcpy(values[i], fTree->GetV1(), entries*sizeof(Double_t));
733 }
734
735 // add points to the fitter
736 for (Int_t i = 0; i < entries; i++){
737 Double_t x[1000];
738 for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
739 fitter->AddPoint(x, values[dim][i], 1);
740 }
741
742 fitter->Eval();
743 fitter->GetParameters(fitParam);
744 fitter->GetCovarianceMatrix(covMatrix);
745 chi2 = fitter->GetChisquare();
746 chi2 = chi2;
747
748 TString *preturnFormula = new TString(Form("( %f+",fitParam[0])), &returnFormula = *preturnFormula;
749
750 for (Int_t iparam = 0; iparam < dim; iparam++) {
751 returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]));
752 if (iparam < dim-1) returnFormula.Append("+");
753 }
754 returnFormula.Append(" )");
755 delete formulaTokens;
756 delete fitter;
757 delete[] values;
758 return preturnFormula;
759}
760
761
762
39bcd65d 763void AliTPCCalibViewer::MakeTreeWithObjects(const char * fileName, TObjArray * array, const char * mapFileName) {
764 //
765 // Write tree with all available information
766 // im mapFileName is speciefied, the Map information are also written to the tree
767 // AliTPCCalPad-Objects are written directly to the tree, so that they can be accessd later on
768 // (does not work!!!)
769 //
770 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
771
772 TObjArray* mapIROCs = 0;
773 TObjArray* mapOROCs = 0;
774 TVectorF *mapIROCArray = 0;
775 TVectorF *mapOROCArray = 0;
776 Int_t mapEntries = 0;
777 TString* mapNames = 0;
778
779 if (mapFileName) {
780 TFile mapFile(mapFileName, "read");
781
782 TList* listOfROCs = mapFile.GetListOfKeys();
783 mapEntries = listOfROCs->GetEntries()/2;
784 mapIROCs = new TObjArray(mapEntries*2);
785 mapOROCs = new TObjArray(mapEntries*2);
786 mapIROCArray = new TVectorF[mapEntries];
787 mapOROCArray = new TVectorF[mapEntries];
788
789 mapNames = new TString[mapEntries];
790 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
791 TString ROCname(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
792 ROCname.Remove(ROCname.Length()-4, 4);
793 mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((ROCname + "IROC").Data()), ivalue);
794 mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((ROCname + "OROC").Data()), ivalue);
795 mapNames[ivalue].Append(ROCname);
796 }
797
798 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
799 mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
800 mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
801
802 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
803 (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
804 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
805 (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
806 }
807
808 } // if (mapFileName)
809
810 TTreeSRedirector cstream(fileName);
811 Int_t arrayEntries = array->GetEntries();
812
813 // Read names of AliTPCCalPads and save them in names[]
814 TString* names = new TString[arrayEntries];
815 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
816 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
817
818 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
819
820 TVectorF *vectorArray = new TVectorF[arrayEntries];
821 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
822 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
823
824
825 //
826 // fill vectors of variable per pad
827 //
828 TVectorF *posArray = new TVectorF[8];
829 for (Int_t ivalue = 0; ivalue < 8; ivalue++)
830 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
831
832 Float_t posG[3] = {0};
833 Float_t posL[3] = {0};
834 Int_t ichannel = 0;
835 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
836 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
837 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
838 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
839 posArray[0][ichannel] = irow;
840 posArray[1][ichannel] = ipad;
841 posArray[2][ichannel] = posL[0];
842 posArray[3][ichannel] = posL[1];
843 posArray[4][ichannel] = posG[0];
844 posArray[5][ichannel] = posG[1];
845 posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
846 posArray[7][ichannel] = ichannel;
847
848 // loop over array containing AliTPCCalPads
849 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
850 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
851 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
852 if (calROC)
853 (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
854 else
855 (vectorArray[ivalue])[ichannel] = 0;
856 }
857 ichannel++;
858 }
859 }
860 AliTPCCalROC dummyROC(0);
861 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
862 AliTPCCalROC *roc = ((AliTPCCalPad*)array->At(ivalue))->GetCalROC(isector);
863 if (!roc) roc = &dummyROC;
864 cstream << "calPads" <<
865 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
866 cstream << "calPads" <<
867 (Char_t*)((names[ivalue] + "Pad.=").Data()) << roc;
868 }
869
870 if (mapFileName) {
871 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
872 if (isector < 36)
873 cstream << "calPads" <<
874 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
875 else
876 cstream << "calPads" <<
877 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
878 }
879 }
880
881 cstream << "calPads" <<
882 "sector=" << isector;
883
884 cstream << "calPads" <<
885 "row.=" << &posArray[0] <<
886 "pad.=" << &posArray[1] <<
887 "lx.=" << &posArray[2] <<
888 "ly.=" << &posArray[3] <<
889 "gx.=" << &posArray[4] <<
890 "gy.=" << &posArray[5] <<
891 "rpad.=" << &posArray[6] <<
892 "channel.=" << &posArray[7];
893
894 cstream << "calPads" <<
895 "\n";
896
897 delete[] posArray;
898 delete[] vectorArray;
899 } //for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++)
900
901 delete[] names;
902 if (mapFileName) {
903 delete mapIROCs;
904 delete mapOROCs;
905 delete[] mapIROCArray;
906 delete[] mapOROCArray;
907 delete[] mapNames;
908 }
909}
910
911void AliTPCCalibViewer::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
912 //
913 // Write a tree with all available information
914 // im mapFileName is speciefied, the Map information are also written to the tree
915 // pads specified in outlierPad are not used for calculating statistics
916 // - the same function as AliTPCCalPad::MakeTree -
917 //
918 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
919
920 TObjArray* mapIROCs = 0;
921 TObjArray* mapOROCs = 0;
922 TVectorF *mapIROCArray = 0;
923 TVectorF *mapOROCArray = 0;
924 Int_t mapEntries = 0;
925 TString* mapNames = 0;
926
927 if (mapFileName) {
928 TFile mapFile(mapFileName, "read");
929
930 TList* listOfROCs = mapFile.GetListOfKeys();
931 mapEntries = listOfROCs->GetEntries()/2;
932 mapIROCs = new TObjArray(mapEntries*2);
933 mapOROCs = new TObjArray(mapEntries*2);
934 mapIROCArray = new TVectorF[mapEntries];
935 mapOROCArray = new TVectorF[mapEntries];
936
937 mapNames = new TString[mapEntries];
938 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
939 TString ROCname(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
940 ROCname.Remove(ROCname.Length()-4, 4);
941 mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((ROCname + "IROC").Data()), ivalue);
942 mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((ROCname + "OROC").Data()), ivalue);
943 mapNames[ivalue].Append(ROCname);
944 }
945
946 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
947 mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
948 mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
949
950 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
951 (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
952 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
953 (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
954 }
955
956 } // if (mapFileName)
957
958 TTreeSRedirector cstream(fileName);
959 Int_t arrayEntries = array->GetEntries();
960
961 TString* names = new TString[arrayEntries];
962 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
963 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
964
965 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
966 //
967 // get statistic for given sector
968 //
969 TVectorF median(arrayEntries);
970 TVectorF mean(arrayEntries);
971 TVectorF rms(arrayEntries);
972 TVectorF ltm(arrayEntries);
973 TVectorF ltmrms(arrayEntries);
974 TVectorF medianWithOut(arrayEntries);
975 TVectorF meanWithOut(arrayEntries);
976 TVectorF rmsWithOut(arrayEntries);
977 TVectorF ltmWithOut(arrayEntries);
978 TVectorF ltmrmsWithOut(arrayEntries);
979
980 TVectorF *vectorArray = new TVectorF[arrayEntries];
981 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
982 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
983
984 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
985 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
986 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
987 AliTPCCalROC* outlierROC = 0;
988 if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
989 if (calROC) {
990 median[ivalue] = calROC->GetMedian();
991 mean[ivalue] = calROC->GetMean();
992 rms[ivalue] = calROC->GetRMS();
993 Double_t ltmrmsValue = 0;
994 ltm[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction);
995 ltmrms[ivalue] = ltmrmsValue;
996 if (outlierROC) {
997 medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
998 meanWithOut[ivalue] = calROC->GetMean(outlierROC);
999 rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
1000 ltmrmsValue = 0;
1001 ltmWithOut[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction, outlierROC);
1002 ltmrmsWithOut[ivalue] = ltmrmsValue;
1003 }
1004 }
1005 else {
1006 median[ivalue] = 0.;
1007 mean[ivalue] = 0.;
1008 rms[ivalue] = 0.;
1009 ltm[ivalue] = 0.;
1010 ltmrms[ivalue] = 0.;
1011 medianWithOut[ivalue] = 0.;
1012 meanWithOut[ivalue] = 0.;
1013 rmsWithOut[ivalue] = 0.;
1014 ltmWithOut[ivalue] = 0.;
1015 ltmrmsWithOut[ivalue] = 0.;
1016 }
1017 }
1018
1019 //
1020 // fill vectors of variable per pad
1021 //
1022 TVectorF *posArray = new TVectorF[8];
1023 for (Int_t ivalue = 0; ivalue < 8; ivalue++)
1024 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
1025
1026 Float_t posG[3] = {0};
1027 Float_t posL[3] = {0};
1028 Int_t ichannel = 0;
1029 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
1030 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
1031 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
1032 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
1033 posArray[0][ichannel] = irow;
1034 posArray[1][ichannel] = ipad;
1035 posArray[2][ichannel] = posL[0];
1036 posArray[3][ichannel] = posL[1];
1037 posArray[4][ichannel] = posG[0];
1038 posArray[5][ichannel] = posG[1];
1039 posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
1040 posArray[7][ichannel] = ichannel;
1041
1042 // loop over array containing AliTPCCalPads
1043 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
1044 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
1045 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
1046 if (calROC)
1047 (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
1048 else
1049 (vectorArray[ivalue])[ichannel] = 0;
1050 }
1051 ichannel++;
1052 }
1053 }
1054
1055 cstream << "calPads" <<
1056 "sector=" << isector;
1057
1058 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
1059 cstream << "calPads" <<
1060 (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
1061 (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
1062 (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
1063 (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
1064 (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
1065 if (outlierPad) {
1066 cstream << "calPads" <<
1067 (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
1068 (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
1069 (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
1070 (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
1071 (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
1072 }
1073 }
1074
1075 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
1076 cstream << "calPads" <<
1077 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
1078 }
1079
1080 if (mapFileName) {
1081 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
1082 if (isector < 36)
1083 cstream << "calPads" <<
1084 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
1085 else
1086 cstream << "calPads" <<
1087 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
1088 }
1089 }
1090
1091 cstream << "calPads" <<
1092 "row.=" << &posArray[0] <<
1093 "pad.=" << &posArray[1] <<
1094 "lx.=" << &posArray[2] <<
1095 "ly.=" << &posArray[3] <<
1096 "gx.=" << &posArray[4] <<
1097 "gy.=" << &posArray[5] <<
1098 "rpad.=" << &posArray[6] <<
1099 "channel.=" << &posArray[7];
1100
1101 cstream << "calPads" <<
1102 "\n";
1103
1104 delete[] posArray;
1105 delete[] vectorArray;
1106 }
1107
1108
1109 delete[] names;
1110 if (mapFileName) {
1111 delete mapIROCs;
1112 delete mapOROCs;
1113 delete[] mapIROCArray;
1114 delete[] mapOROCArray;
1115 delete[] mapNames;
1116 }
1117}
1118