Line data Source code
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 : /// \class AliTPCCalibViewer
18 : /// \brief Class for viewing/visualizing TPC calibration data
19 : ///
20 : /// base on TTree functionality for visualization
21 : ///
22 : /// Create a list of AliTPCCalPads, arrange them in an TObjArray.
23 : /// Pass this TObjArray to MakeTree and create the calibration Tree
24 : /// While craating this tree some statistical information are calculated
25 : /// Open the viewer with this Tree: AliTPCCalibViewer v("CalibTree.root")
26 : /// Have fun!
27 : /// EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
28 : ///
29 : /// If you like to click, we recommand you the
30 : /// AliTPCCalibViewerGUI
31 : ///
32 : /// THE DOCUMENTATION IS STILL NOT COMPLETED !!!!
33 : ///
34 : /// ROOT includes
35 :
36 : #include <fstream>
37 : #include <iostream>
38 :
39 : #include <TFile.h>
40 : #include <TFriendElement.h>
41 : #include <TGraph.h>
42 : #include <TKey.h>
43 : #include <TPad.h>
44 : //#include <TCanvas.h>
45 : #include <TProfile2D.h>
46 : #include <TH1.h>
47 : #include <TH1F.h>
48 : #include <TLegend.h>
49 : #include <TLine.h>
50 : #include <TMath.h>
51 : #include <TObjString.h>
52 : //#include <TROOT.h>
53 : #include <TRandom.h>
54 : #include <TString.h>
55 : #include <TStyle.h>
56 : #include <TTreeStream.h>
57 :
58 : #include "AliTPCCalibCE.h"
59 : #include "AliMathBase.h"
60 : #include "AliTPCCalPad.h"
61 : #include "AliTPCCalROC.h"
62 : #include "AliTPCCalibPedestal.h"
63 : #include "AliTPCCalibPulser.h"
64 :
65 : //
66 : // AliRoot includes
67 : //
68 : #include "AliTPCCalibViewer.h"
69 :
70 : using std::ifstream;
71 : /// \cond CLASSIMP
72 24 : ClassImp(AliTPCCalibViewer)
73 : /// \endcond
74 :
75 :
76 : AliTPCCalibViewer::AliTPCCalibViewer()
77 0 : :TObject(),
78 0 : fTree(0),
79 0 : fFile(0),
80 0 : fListOfObjectsToBeDeleted(0),
81 0 : fTreeMustBeDeleted(0),
82 0 : fAbbreviation(0),
83 0 : fAppendString(0)
84 0 : {
85 : //
86 : // Default constructor
87 : //
88 :
89 0 : }
90 :
91 : //_____________________________________________________________________________
92 : AliTPCCalibViewer::AliTPCCalibViewer(const AliTPCCalibViewer &c)
93 0 : :TObject(c),
94 0 : fTree(0),
95 0 : fFile(0),
96 0 : fListOfObjectsToBeDeleted(0),
97 0 : fTreeMustBeDeleted(0),
98 0 : fAbbreviation(0),
99 0 : fAppendString(0)
100 0 : {
101 : /// dummy AliTPCCalibViewer copy constructor
102 : /// not yet working!!!
103 :
104 0 : fTree = c.fTree;
105 0 : fTreeMustBeDeleted = c.fTreeMustBeDeleted;
106 : //fFile = new TFile(*(c.fFile));
107 0 : fListOfObjectsToBeDeleted = c.fListOfObjectsToBeDeleted;
108 0 : fAbbreviation = c.fAbbreviation;
109 0 : fAppendString = c.fAppendString;
110 0 : }
111 :
112 : //_____________________________________________________________________________
113 : AliTPCCalibViewer::AliTPCCalibViewer(TTree *const tree)
114 0 : :TObject(),
115 0 : fTree(0),
116 0 : fFile(0),
117 0 : fListOfObjectsToBeDeleted(0),
118 0 : fTreeMustBeDeleted(0),
119 0 : fAbbreviation(0),
120 0 : fAppendString(0)
121 0 : {
122 : /// Constructor that initializes the calibration viewer
123 :
124 0 : fTree = tree;
125 0 : fTreeMustBeDeleted = kFALSE;
126 0 : fListOfObjectsToBeDeleted = new TObjArray();
127 0 : fAbbreviation = "~";
128 0 : fAppendString = ".fElements";
129 0 : }
130 :
131 : //_____________________________________________________________________________
132 : AliTPCCalibViewer::AliTPCCalibViewer(const char* fileName, const char* treeName)
133 0 : :TObject(),
134 0 : fTree(0),
135 0 : fFile(0),
136 0 : fListOfObjectsToBeDeleted(0),
137 0 : fTreeMustBeDeleted(0),
138 0 : fAbbreviation(0),
139 0 : fAppendString(0)
140 :
141 0 : {
142 : /// Constructor to initialize the calibration viewer
143 : /// the file 'fileName' contains the tree 'treeName'
144 :
145 0 : fFile = new TFile(fileName, "read");
146 0 : fTree = (TTree*) fFile->Get(treeName);
147 0 : fTreeMustBeDeleted = kTRUE;
148 0 : fListOfObjectsToBeDeleted = new TObjArray();
149 0 : fAbbreviation = "~";
150 0 : fAppendString = ".fElements";
151 0 : }
152 :
153 : //____________________________________________________________________________
154 : AliTPCCalibViewer & AliTPCCalibViewer::operator =(const AliTPCCalibViewer & param)
155 : {
156 : /// assignment operator - dummy
157 : /// not yet working!!!
158 :
159 0 : if (this == ¶m) return (*this);
160 :
161 0 : fTree = param.fTree;
162 0 : fTreeMustBeDeleted = param.fTreeMustBeDeleted;
163 : //fFile = new TFile(*(param.fFile));
164 0 : fListOfObjectsToBeDeleted = param.fListOfObjectsToBeDeleted;
165 0 : fAbbreviation = param.fAbbreviation;
166 0 : fAppendString = param.fAppendString;
167 0 : return (*this);
168 0 : }
169 :
170 : //_____________________________________________________________________________
171 : AliTPCCalibViewer::~AliTPCCalibViewer()
172 0 : {
173 : /// AliTPCCalibViewer destructor
174 : /// all objects will be deleted, the file will be closed, the pictures will disappear
175 :
176 0 : if (fTree && fTreeMustBeDeleted) {
177 0 : fTree->SetCacheSize(0);
178 0 : fTree->Delete();
179 : //delete fTree;
180 : }
181 0 : if (fFile) {
182 0 : fFile->Close();
183 0 : fFile = 0;
184 0 : }
185 :
186 0 : for (Int_t i = fListOfObjectsToBeDeleted->GetEntriesFast()-1; i >= 0; i--) {
187 : //cout << "Index " << i << " trying to delete the following object: " << fListOfObjectsToBeDeleted->At(i)->GetName() << "..."<< endl;
188 0 : delete fListOfObjectsToBeDeleted->At(i);
189 : }
190 0 : delete fListOfObjectsToBeDeleted;
191 0 : }
192 :
193 : //_____________________________________________________________________________
194 : void AliTPCCalibViewer::Delete(Option_t* /*option*/) {
195 : /// Should be called from AliTPCCalibViewerGUI class only.
196 : /// If you use Delete() do not call the destructor.
197 : /// All objects (except those contained in fListOfObjectsToBeDeleted) will be deleted, the file will be closed.
198 :
199 0 : if (fTree && fTreeMustBeDeleted) {
200 0 : fTree->SetCacheSize(0);
201 0 : fTree->Delete();
202 0 : }
203 0 : delete fFile;
204 0 : delete fListOfObjectsToBeDeleted;
205 0 : }
206 :
207 :
208 : const char* AliTPCCalibViewer::AddAbbreviations(const Char_t *c, Bool_t printDrawCommand){
209 : /// Replace all "<variable>" with "<variable><fAbbreviation>" (Adds forgotten "~")
210 : /// but take care on the statistical information, like "CEQmean_Mean"
211 : /// and also take care on correct given variables, like "CEQmean~"
212 : ///
213 : /// For each variable out of "listOfVariables":
214 : /// - 'Save' correct items:
215 : /// - form <replaceString>, take <variable>'s first char, add <removeString>, add rest of <variable>, e.g. "C!#EQmean" (<removeString> = "!#")
216 : /// - For each statistical information in "listOfNormalizationVariables":
217 : /// - ReplaceAll <variable><statistical_Information> with <replaceString><statistical_Information>
218 : /// - ReplaceAll <variable><abbreviation> with <replaceString><abbreviation>, e.g. "CEQmean~" -> "C!#EQmean~"
219 : /// - ReplaceAll <variable><appendStr> with <replaceString><appendStr>, e.g. "CEQmean.fElements" -> "C!#EQmean.fElements"
220 : ///
221 : /// - Do actual replacing:
222 : /// - ReplaceAll <variable> with <variable><fAbbreviation>, e.g. "CEQmean" -> "CEQmean~"
223 : ///
224 : /// - Undo saving:
225 : /// - For each statistical information in "listOfNormalizationVariables":
226 : /// - ReplaceAll <replaceString><statistical_Information> with <variable><statistical_Information>
227 : /// - ReplaceAll <replaceString><abbreviation> with <variable><abbreviation>, e.g. "C!#EQmean~" -> "CEQmean~"
228 : /// - ReplaceAll <replaceString><appendStr> with <variable><appendStr>, e.g. "C!#EQmean.fElements" -> "CEQmean.fElements"
229 : ///
230 : /// Now all the missing "~" should be added.
231 :
232 0 : TString str(c);
233 0 : TString removeString = "!#"; // very unpropable combination of chars
234 0 : TString replaceString = "";
235 0 : TString searchString = "";
236 0 : TString normString = "";
237 0 : TObjArray *listOfVariables = GetListOfVariables();
238 0 : listOfVariables->Add(new TObjString("channel"));
239 0 : listOfVariables->Add(new TObjString("gx"));
240 0 : listOfVariables->Add(new TObjString("gy"));
241 0 : listOfVariables->Add(new TObjString("lx"));
242 0 : listOfVariables->Add(new TObjString("ly"));
243 0 : listOfVariables->Add(new TObjString("pad"));
244 0 : listOfVariables->Add(new TObjString("row"));
245 0 : listOfVariables->Add(new TObjString("rpad"));
246 0 : listOfVariables->Add(new TObjString("sector"));
247 0 : TObjArray *listOfNormalizationVariables = GetListOfNormalizationVariables();
248 0 : Int_t nVariables = listOfVariables->GetEntriesFast();
249 0 : Int_t nNorm = listOfNormalizationVariables->GetEntriesFast();
250 :
251 0 : Int_t *varLengths = new Int_t[nVariables];
252 0 : for (Int_t i = 0; i < nVariables; i++) {
253 0 : varLengths[i] = ((TObjString*)listOfVariables->At(i))->String().Length();
254 : }
255 0 : Int_t *normLengths = new Int_t[nNorm];
256 0 : for (Int_t i = 0; i < nNorm; i++) {
257 0 : normLengths[i] = ((TObjString*)listOfNormalizationVariables->At(i))->String().Length();
258 : // printf("normLengths[%i] (%s) = %i \n", i,((TObjString*)listOfNormalizationVariables->At(i))->String().Data(), normLengths[i]);
259 : }
260 0 : Int_t *varSort = new Int_t[nVariables];
261 0 : TMath::Sort(nVariables, varLengths, varSort, kTRUE);
262 0 : Int_t *normSort = new Int_t[nNorm];
263 0 : TMath::Sort(nNorm, normLengths, normSort, kTRUE);
264 : // for (Int_t i = 0; i<nNorm; i++) printf("normLengths: %i\n", normLengths[normSort[i]]);
265 : // for (Int_t i = 0; i<nVariables; i++) printf("varLengths: %i\n", varLengths[varSort[i]]);
266 :
267 0 : for (Int_t ivar = 0; ivar < nVariables; ivar++) {
268 : // ***** save correct tokens *****
269 : // first get the next variable:
270 0 : searchString = ((TObjString*)listOfVariables->At(varSort[ivar]))->String();
271 : // printf("searchString: %s ++++++++++++++\n", searchString.Data());
272 : // form replaceString:
273 0 : replaceString = "";
274 0 : for (Int_t i = 0; i < searchString.Length(); i++) {
275 0 : replaceString.Append(searchString[i]);
276 0 : if (i == 0) replaceString.Append(removeString);
277 : }
278 : // go through normalization:
279 : // printf("go through normalization\n");
280 0 : for (Int_t inorm = 0; inorm < nNorm; inorm++) {
281 : // printf(" inorm=%i, nNorm=%i, normSort[inorm]=%i \n", inorm, nNorm, normSort[inorm]);
282 0 : normString = ((TObjString*)listOfNormalizationVariables->At(normSort[inorm]))->String();
283 : // printf(" walking in normalization, i=%i, normString=%s \n", inorm, normString.Data());
284 0 : str.ReplaceAll(searchString + normString, replaceString + normString);
285 : // like: str.ReplaceAll("CEQmean_Mean", "C!EQmean_Mean");
286 : }
287 0 : str.ReplaceAll(searchString + fAbbreviation, replaceString + fAbbreviation);
288 : // like: str.ReplaceAll("CEQmean~", "C!EQmean~");
289 0 : str.ReplaceAll(searchString + fAppendString, replaceString + fAppendString);
290 : // like: str.ReplaceAll("CEQmean.fElements", "C!EQmean.fElements");
291 :
292 : // ***** add missing extensions *****
293 0 : str.ReplaceAll(searchString, replaceString + fAbbreviation);
294 : // like: str.ReplaceAll("CEQmean", "C!EQmean~");
295 : }
296 :
297 : // ***** undo saving *****
298 0 : str.ReplaceAll(removeString, "");
299 :
300 0 : if (printDrawCommand) std::cout << "The string looks now like: " << str.Data() << std::endl;
301 0 : delete [] varLengths;
302 0 : delete [] normLengths;
303 0 : delete [] varSort;
304 0 : delete [] normSort;
305 0 : return str.Data();
306 0 : }
307 :
308 :
309 :
310 :
311 : //_____________________________________________________________________________
312 : Int_t AliTPCCalibViewer::EasyDraw(const char* drawCommand, const char* sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
313 : /// easy drawing of data, use '~' for abbreviation of '.fElements'
314 : /// example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
315 : /// sector: sector-number - only the specified sector will be drwawn
316 : /// 'A'/'C' or 'a'/'c' - side A/C will be drawn
317 : /// 'ALL' - whole TPC will be drawn, projected on one side
318 : /// cuts: specifies cuts
319 : /// drawOptions: draw options like 'same'
320 : /// writeDrawCommand: write the command, that is passed to TTree::Draw
321 :
322 0 : TString drawStr(drawCommand);
323 0 : TString sectorStr(sector);
324 0 : sectorStr.ToUpper();
325 0 : TString cutStr("");
326 : //TString drawOptionsStr("profcolz ");
327 0 : Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
328 0 : if (dangerousToDraw) {
329 0 : Warning("EasyDraw", "The draw string must not contain ':' or '>>'. Using only first variable for drawing!");
330 : // return -1;
331 : // drawStr.Resize(drawStr.First(">"));
332 0 : drawStr.Resize(drawStr.First(":"));
333 : }
334 :
335 0 : TString drawOptionsStr("");
336 0 : TRandom rnd(0);
337 0 : Int_t rndNumber = rnd.Integer(10000);
338 :
339 0 : if (drawOptions && strcmp(drawOptions, "") != 0)
340 0 : drawOptionsStr += drawOptions;
341 : else
342 0 : drawOptionsStr += "profcolz";
343 :
344 : TProfile2D *prof2D=0x0;
345 0 : const TString profName=TString::Format("prof%d",rndNumber);
346 :
347 0 : TString secLabel = "";
348 :
349 0 : TString secStr=sectorStr;
350 0 : secStr.ReplaceAll("S","");
351 0 : if (secStr.IsDigit()) {
352 0 : const Int_t secDigit=secStr.Atoi();
353 0 : if ( secDigit < 36 )
354 0 : secLabel = "IROC";
355 : else
356 0 : secLabel = "OROC";
357 :
358 0 : if (sectorStr.Contains("S"))
359 0 : secLabel="Sector";
360 :
361 0 : if ( secDigit%36<18 ) //A-Side
362 0 : secLabel += ", A";
363 : else
364 0 : secLabel += ", C";
365 :
366 0 : secLabel += Form("%02d",secDigit%18);
367 0 : secLabel += ";row;pad";
368 0 : }
369 :
370 0 : if (sectorStr == "A") {
371 0 : drawStr += Form(":gy%s:gx%s>>prof", fAppendString.Data(), fAppendString.Data());
372 0 : drawStr += rndNumber;
373 : // drawStr += "(330,-250,250,330,-250,250)";
374 0 : cutStr += "(sector/18)%2==0 ";
375 0 : prof2D = new TProfile2D(profName,"A-Side;x (cm); y (cm)",330,-250,250,330,-250,250);
376 0 : }
377 0 : else if (sectorStr == "C") {
378 0 : drawStr += Form(":gy%s:gx%s>>prof", fAppendString.Data(), fAppendString.Data());
379 0 : drawStr += rndNumber;
380 : // drawStr += "(330,-250,250,330,-250,250)";
381 0 : cutStr += "(sector/18)%2==1 ";
382 0 : prof2D = new TProfile2D(profName,"C-Side;x (cm); y (cm)",330,-250,250,330,-250,250);
383 0 : }
384 0 : else if (sectorStr == "ALL") {
385 0 : drawStr += Form(":gy%s:gx%s>>prof", fAppendString.Data(), fAppendString.Data());
386 0 : drawStr += rndNumber;
387 : // drawStr += "(330,-250,250,330,-250,250)";
388 0 : prof2D = new TProfile2D(profName,"A&C-Side average;x (cm); y (cm)",330,-250,250,330,-250,250);
389 0 : }
390 0 : else if (sectorStr.Contains("S")) {
391 0 : drawStr += Form(":rpad%s:row%s+(sector>35)*63>>prof", fAppendString.Data(), fAppendString.Data());
392 0 : drawStr += rndNumber;
393 : // drawStr += "(159,0,159,140,-70,70)";
394 0 : TString sec=sectorStr;
395 0 : sec.Remove(0,1);
396 0 : Int_t isec = sec.Atoi();
397 0 : cutStr += "sector%36=="+sec+" ";
398 0 : prof2D = new TProfile2D(profName,secLabel,159,0,159,140,-70,70);
399 0 : }
400 0 : else if (sectorStr.IsDigit()) {
401 0 : Int_t isec = sectorStr.Atoi();
402 0 : drawStr += Form(":rpad%s:row%s>>prof", fAppendString.Data(), fAppendString.Data());
403 0 : drawStr += rndNumber;
404 0 : if (isec < 36 && isec >= 0){
405 : // drawStr += "(63,0,63,108,-54,54)";
406 0 : prof2D = new TProfile2D(profName,secLabel,63,0,63,108,-54,54);
407 0 : } else if (isec < 72 && isec >= 36) {
408 : // drawStr += "(96,0,96,140,-70,70)";
409 0 : prof2D = new TProfile2D(profName,secLabel,96,0,96,140,-70,70);
410 : } else {
411 0 : Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
412 0 : return -1;
413 : }
414 0 : cutStr += "(sector==";
415 0 : cutStr += isec;
416 0 : cutStr += ") ";
417 0 : }
418 :
419 0 : if (cuts && cuts[0] != 0) {
420 0 : if (cutStr.Length() != 0) cutStr += "&& ";
421 0 : cutStr += "(";
422 0 : cutStr += cuts;
423 0 : cutStr += ")";
424 : }
425 0 : drawStr.ReplaceAll(fAbbreviation, fAppendString);
426 0 : cutStr.ReplaceAll(fAbbreviation, fAppendString);
427 0 : if (writeDrawCommand) std::cout << "fTree->Draw(\"" << drawStr << "\", \"" << cutStr << "\", \"" << drawOptionsStr << "\");" << std::endl;
428 0 : Int_t returnValue = fTree->Draw(drawStr.Data(), cutStr.Data(), drawOptionsStr.Data());
429 : // TString profName("prof");
430 : // profName += rndNumber;
431 : // TObject *obj = gDirectory->Get(profName.Data());
432 0 : TObject *obj=prof2D;
433 0 : if (obj && obj->InheritsFrom("TH1")) FormatHistoLabels((TH1*)obj);
434 : return returnValue;
435 0 : }
436 :
437 :
438 : Int_t AliTPCCalibViewer::EasyDraw(const char* drawCommand, Int_t sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
439 : /// easy drawing of data, use '~' for abbreviation of '.fElements'
440 : /// example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0")
441 : /// sector: sector-number - only the specified sector will be drwawn
442 : /// cuts: specifies cuts
443 : /// drawOptions: draw options like 'same'
444 : /// writeDrawCommand: write the command, that is passed to TTree::Draw
445 :
446 0 : if (sector >= 0 && sector < 72) {
447 0 : return EasyDraw(drawCommand, Form("%i", sector), cuts, drawOptions, writeDrawCommand);
448 : }
449 0 : Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
450 0 : return -1;
451 0 : }
452 :
453 :
454 : //_____________________________________________________________________________
455 : Int_t AliTPCCalibViewer::EasyDraw1D(const char* drawCommand, const char* sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
456 : /// easy drawing of data, use '~' for abbreviation of '.fElements'
457 : /// example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
458 : /// sector: sector-number - the specified sector will be drwawn
459 : /// 'A'/'C' or 'a'/'c' - side A/C will be drawn
460 : /// 'ALL' - whole TPC will be drawn, projected on one side
461 : /// cuts: specifies cuts
462 : /// drawOptions: draw options like 'same'
463 : /// writeDrawCommand: write the command, that is passed to TTree::Draw
464 :
465 0 : TString drawStr(drawCommand);
466 0 : TString sectorStr(sector);
467 0 : TString drawOptionsStr(drawOptions);
468 0 : sectorStr.ToUpper();
469 0 : TString cutStr("");
470 :
471 0 : if (sectorStr == "A")
472 0 : cutStr += "(sector/18)%2==0 ";
473 0 : else if (sectorStr == "C")
474 0 : cutStr += "(sector/18)%2==1 ";
475 0 : else if (sectorStr.IsDigit()) {
476 0 : Int_t isec = sectorStr.Atoi();
477 0 : if (isec < 0 || isec > 71) {
478 0 : Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
479 0 : return -1;
480 : }
481 0 : cutStr += "(sector==";
482 0 : cutStr += isec;
483 0 : cutStr += ") ";
484 0 : }
485 0 : else if (sectorStr.Contains("S")) {
486 0 : TString sec=sectorStr;
487 0 : sec.Remove(0,1);
488 0 : cutStr += "sector%36=="+sec+" ";
489 0 : }
490 :
491 0 : if (cuts && cuts[0] != 0) {
492 0 : if (cutStr.Length() != 0) cutStr += "&& ";
493 0 : cutStr += "(";
494 0 : cutStr += cuts;
495 0 : cutStr += ")";
496 : }
497 :
498 0 : drawStr.ReplaceAll(fAbbreviation, fAppendString);
499 0 : cutStr.ReplaceAll(fAbbreviation, fAppendString);
500 0 : if (writeDrawCommand) std::cout << "fTree->Draw(\"" << drawStr << "\", \"" << cutStr << "\", \"" << drawOptionsStr << "\");" << std::endl;
501 0 : Int_t returnValue = fTree->Draw(drawStr.Data(), cutStr.Data(), drawOptionsStr.Data());
502 0 : if (returnValue == -1) return -1;
503 :
504 0 : TObject *obj = (gPad) ? gPad->GetPrimitive("htemp") : 0;
505 0 : if (!obj) obj = (TH1F*)gDirectory->Get("htemp");
506 0 : if (!obj) obj = gPad->GetPrimitive("tempHist");
507 0 : if (!obj) obj = (TH1F*)gDirectory->Get("tempHist");
508 0 : if (!obj) obj = gPad->GetPrimitive("Graph");
509 0 : if (!obj) obj = (TH1F*)gDirectory->Get("Graph");
510 0 : if (obj && obj->InheritsFrom("TH1")) FormatHistoLabels((TH1*)obj);
511 : return returnValue;
512 0 : }
513 :
514 :
515 : Int_t AliTPCCalibViewer::EasyDraw1D(const char* drawCommand, Int_t sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
516 : /// easy drawing of data, use '~' for abbreviation of '.fElements'
517 : /// example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0")
518 : /// sector: sector-number - the specified sector will be drwawn
519 : /// cuts: specifies cuts
520 : /// drawOptions: draw options like 'same'
521 : /// writeDrawCommand: write the command, that is passed to TTree::Draw
522 :
523 0 : if (sector >= 0 && sector < 72) {
524 0 : return EasyDraw1D(drawCommand, Form("%i",sector), cuts, drawOptions, writeDrawCommand);
525 : }
526 0 : Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
527 0 : return -1;
528 0 : }
529 :
530 :
531 : void AliTPCCalibViewer::FormatHistoLabels(TH1 *histo) const {
532 : /// formats title and axis labels of histo
533 : /// removes '.fElements'
534 :
535 0 : if (!histo) return;
536 0 : TString replaceString(fAppendString.Data());
537 0 : TString *str = new TString(histo->GetTitle());
538 0 : str->ReplaceAll(replaceString, "");
539 0 : histo->SetTitle(str->Data());
540 0 : delete str;
541 0 : if (histo->GetXaxis()) {
542 0 : str = new TString(histo->GetXaxis()->GetTitle());
543 0 : str->ReplaceAll(replaceString, "");
544 0 : histo->GetXaxis()->SetTitle(str->Data());
545 0 : delete str;
546 : }
547 0 : if (histo->GetYaxis()) {
548 0 : str = new TString(histo->GetYaxis()->GetTitle());
549 0 : str->ReplaceAll(replaceString, "");
550 0 : histo->GetYaxis()->SetTitle(str->Data());
551 0 : delete str;
552 : }
553 0 : if (histo->GetZaxis()) {
554 0 : str = new TString(histo->GetZaxis()->GetTitle());
555 0 : str->ReplaceAll(replaceString, "");
556 0 : histo->GetZaxis()->SetTitle(str->Data());
557 0 : delete str;
558 : }
559 0 : }
560 :
561 :
562 : Int_t AliTPCCalibViewer::DrawHisto1D(const char* drawCommand, Int_t sector, const char* cuts, const char *sigmas, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM) const {
563 : /// Easy drawing of data, in principle the same as EasyDraw1D
564 : /// Difference: A line for the mean / median / LTM is drawn
565 : /// in 'sigmas' you can specify in which distance to the mean/median/LTM you want to see a line in sigma-units, separated by ';'
566 : /// example: sigmas = "2; 4; 6;" at \f$2 \sigma\f$, \f$4 \sigma\f$ and \f$6 \sigma\f$ a line is drawn.
567 : /// "plotMean", "plotMedian" and "plotLTM": what kind of lines do you want to see?
568 :
569 0 : if (sector >= 0 && sector < 72) {
570 0 : return DrawHisto1D(drawCommand, Form("%i", sector), cuts, sigmas, plotMean, plotMedian, plotLTM);
571 : }
572 0 : Error("DrawHisto1D","The TPC contains only sectors between 0 and 71.");
573 0 : return -1;
574 0 : }
575 :
576 :
577 : Int_t AliTPCCalibViewer::DrawHisto1D(const char* drawCommand, const char* sector, const char* cuts, const char *sigmas, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM) const {
578 : /// Easy drawing of data, in principle the same as EasyDraw1D
579 : /// Difference: A line for the mean / median / LTM is drawn
580 : /// in 'sigmas' you can specify in which distance to the mean/median/LTM you want to see a line in sigma-units, separated by ';'
581 : /// example: sigmas = "2; 4; 6;" at \f$2 \sigma\f$, \f$4 \sigma\f$ and \f$6 \sigma\f$ a line is drawn.
582 : /// "plotMean", "plotMedian" and "plotLTM": what kind of lines do you want to see?
583 :
584 0 : Int_t oldOptStat = gStyle->GetOptStat();
585 0 : gStyle->SetOptStat(0000000);
586 : Double_t ltmFraction = 0.8;
587 :
588 0 : TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
589 0 : TVectorF nsigma(sigmasTokens->GetEntriesFast());
590 0 : for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
591 0 : TString str(((TObjString*)sigmasTokens->At(i))->GetString());
592 0 : Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
593 0 : nsigma[i] = sig;
594 0 : }
595 0 : delete sigmasTokens;
596 :
597 0 : TString drawStr(drawCommand);
598 0 : Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
599 0 : if (dangerousToDraw) {
600 0 : Warning("DrawHisto1D", "The draw string must not contain ':' or '>>'.");
601 0 : return -1;
602 : }
603 0 : drawStr += " >> tempHist";
604 0 : Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts);
605 0 : TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
606 : // FIXME is this histogram deleted automatically?
607 0 : Double_t *values = fTree->GetV1(); // value is the array containing 'entries' numbers
608 :
609 0 : Double_t mean = TMath::Mean(entries, values);
610 0 : Double_t median = TMath::Median(entries, values);
611 0 : Double_t sigma = TMath::RMS(entries, values);
612 0 : Double_t maxY = htemp->GetMaximum();
613 :
614 0 : TLegend * legend = new TLegend(.7,.7, .99, .99, "Statistical information");
615 : //fListOfObjectsToBeDeleted->Add(legend);
616 :
617 0 : if (plotMean) {
618 : // draw Mean
619 0 : TLine* line = new TLine(mean, 0, mean, maxY);
620 : //fListOfObjectsToBeDeleted->Add(line);
621 0 : line->SetLineColor(kRed);
622 0 : line->SetLineWidth(2);
623 0 : line->SetLineStyle(1);
624 0 : line->Draw();
625 0 : legend->AddEntry(line, Form("Mean: %f", mean), "l");
626 : // draw sigma lines
627 0 : for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
628 0 : TLine* linePlusSigma = new TLine(mean + nsigma[i] * sigma, 0, mean + nsigma[i] * sigma, maxY);
629 : //fListOfObjectsToBeDeleted->Add(linePlusSigma);
630 0 : linePlusSigma->SetLineColor(kRed);
631 0 : linePlusSigma->SetLineStyle(2 + i);
632 0 : linePlusSigma->Draw();
633 0 : TLine* lineMinusSigma = new TLine(mean - nsigma[i] * sigma, 0, mean - nsigma[i] * sigma, maxY);
634 : //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
635 0 : lineMinusSigma->SetLineColor(kRed);
636 0 : lineMinusSigma->SetLineStyle(2 + i);
637 0 : lineMinusSigma->Draw();
638 0 : legend->AddEntry(lineMinusSigma, Form("%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i] * sigma)), "l");
639 : }
640 0 : }
641 0 : if (plotMedian) {
642 : // draw median
643 0 : TLine* line = new TLine(median, 0, median, maxY);
644 : //fListOfObjectsToBeDeleted->Add(line);
645 0 : line->SetLineColor(kBlue);
646 0 : line->SetLineWidth(2);
647 0 : line->SetLineStyle(1);
648 0 : line->Draw();
649 0 : legend->AddEntry(line, Form("Median: %f", median), "l");
650 : // draw sigma lines
651 0 : for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
652 0 : TLine* linePlusSigma = new TLine(median + nsigma[i] * sigma, 0, median + nsigma[i]*sigma, maxY);
653 : //fListOfObjectsToBeDeleted->Add(linePlusSigma);
654 0 : linePlusSigma->SetLineColor(kBlue);
655 0 : linePlusSigma->SetLineStyle(2 + i);
656 0 : linePlusSigma->Draw();
657 0 : TLine* lineMinusSigma = new TLine(median - nsigma[i] * sigma, 0, median - nsigma[i]*sigma, maxY);
658 : //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
659 0 : lineMinusSigma->SetLineColor(kBlue);
660 0 : lineMinusSigma->SetLineStyle(2 + i);
661 0 : lineMinusSigma->Draw();
662 0 : legend->AddEntry(lineMinusSigma, Form("%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i] * sigma)), "l");
663 : }
664 0 : }
665 0 : if (plotLTM) {
666 : // draw LTM
667 0 : Double_t ltmRms = 0;
668 0 : Double_t ltm = GetLTM(entries, values, <mRms, ltmFraction);
669 0 : TLine* line = new TLine(ltm, 0, ltm, maxY);
670 : //fListOfObjectsToBeDeleted->Add(line);
671 0 : line->SetLineColor(kGreen+2);
672 0 : line->SetLineWidth(2);
673 0 : line->SetLineStyle(1);
674 0 : line->Draw();
675 0 : legend->AddEntry(line, Form("LTM: %f", ltm), "l");
676 : // draw sigma lines
677 0 : for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
678 0 : TLine* linePlusSigma = new TLine(ltm + nsigma[i] * ltmRms, 0, ltm + nsigma[i] * ltmRms, maxY);
679 : //fListOfObjectsToBeDeleted->Add(linePlusSigma);
680 0 : linePlusSigma->SetLineColor(kGreen+2);
681 0 : linePlusSigma->SetLineStyle(2+i);
682 0 : linePlusSigma->Draw();
683 :
684 0 : TLine* lineMinusSigma = new TLine(ltm - nsigma[i] * ltmRms, 0, ltm - nsigma[i] * ltmRms, maxY);
685 : //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
686 0 : lineMinusSigma->SetLineColor(kGreen+2);
687 0 : lineMinusSigma->SetLineStyle(2+i);
688 0 : lineMinusSigma->Draw();
689 0 : legend->AddEntry(lineMinusSigma, Form("%i #sigma = %f", (Int_t)(nsigma[i]), (Float_t)(nsigma[i] * ltmRms)), "l");
690 : }
691 0 : }
692 0 : if (!plotMean && !plotMedian && !plotLTM) return -1;
693 0 : legend->Draw();
694 0 : gStyle->SetOptStat(oldOptStat);
695 0 : return 1;
696 0 : }
697 :
698 :
699 : Int_t AliTPCCalibViewer::SigmaCut(const char* drawCommand, Int_t sector, const char* cuts, Float_t sigmaMax, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, Bool_t pm, const char *sigmas, Float_t sigmaStep) const {
700 : /// Creates a histogram \f$S(t, \mu, \sigma)\f$, where you can see, how much of the data are inside sigma-intervals around the mean value
701 : /// The data of the distribution \f$f(x, \mu, \sigma)\f$ are given in 'array', 'n' specifies the length of the array
702 : /// 'mean' and 'sigma' are \f$\mu\f$ and \f$\sigma\f$ of the distribution in 'array', to be specified by the user
703 : /// 'nbins': number of bins, 'binLow': first bin, 'binUp': last bin
704 : /// sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, \f$t \sigma\f$)
705 : /// sigmaStep: the binsize of the generated histogram
706 : /// \f[
707 : /// f(x, \mu, \sigma) \Rightarrow S(t, \mu, \sigma) = (\int_{\mu}^{\mu + t \sigma} f(x, \mu, \sigma) dx + \int_{\mu}^{\mu - t \sigma} f(x, \mu, \sigma) dx) / (\int_{-\infty}^{+\infty} f(x, \mu, \sigma) dx)
708 : /// \f]
709 : ///
710 : /// Creates a histogram, where you can see, how much of the data are inside sigma-intervals
711 : /// around the mean/median/LTM
712 : /// with drawCommand, sector and cuts you specify your input data, see EasyDraw
713 : /// sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma)
714 : /// sigmaStep: the binsize of the generated histogram
715 : /// plotMean/plotMedian/plotLTM: specifies where to put the center
716 :
717 0 : if (sector >= 0 && sector < 72) {
718 0 : return SigmaCut(drawCommand, Form("%i", sector), cuts, sigmaMax, plotMean, plotMedian, plotLTM, pm, sigmas, sigmaStep);
719 : }
720 0 : Error("SigmaCut","The TPC contains only sectors between 0 and 71.");
721 0 : return -1;
722 0 : }
723 :
724 :
725 : Int_t AliTPCCalibViewer::SigmaCut(const char* drawCommand, const char* sector, const char* cuts, Float_t sigmaMax, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, Bool_t pm, const char *sigmas, Float_t sigmaStep) const {
726 : /// Creates a histogram, where you can see, how much of the data are inside sigma-intervals
727 : /// around the mean/median/LTM
728 : /// with drawCommand, sector and cuts you specify your input data, see EasyDraw
729 : /// sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma)
730 : /// sigmaStep: the binsize of the generated histogram
731 : /// plotMean/plotMedian/plotLTM: specifies where to put the center
732 :
733 : Double_t ltmFraction = 0.8;
734 :
735 0 : TString drawStr(drawCommand);
736 0 : Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
737 0 : if (dangerousToDraw) {
738 0 : Warning("SigmaCut", "The draw string must not contain ':' or '>>'.");
739 0 : return -1;
740 : }
741 0 : drawStr += " >> tempHist";
742 :
743 0 : Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
744 0 : TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
745 : // FIXME is this histogram deleted automatically?
746 0 : Double_t *values = fTree->GetV1(); // value is the array containing 'entries' numbers
747 :
748 0 : Double_t mean = TMath::Mean(entries, values);
749 0 : Double_t median = TMath::Median(entries, values);
750 0 : Double_t sigma = TMath::RMS(entries, values);
751 :
752 0 : TLegend * legend = new TLegend(.7,.7, .99, .99, "Cumulative");
753 : //fListOfObjectsToBeDeleted->Add(legend);
754 : TH1F *cutHistoMean = 0;
755 : TH1F *cutHistoMedian = 0;
756 : TH1F *cutHistoLTM = 0;
757 :
758 0 : TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
759 0 : TVectorF nsigma(sigmasTokens->GetEntriesFast());
760 0 : for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
761 0 : TString str(((TObjString*)sigmasTokens->At(i))->GetString());
762 0 : Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
763 0 : nsigma[i] = sig;
764 0 : }
765 0 : delete sigmasTokens;
766 :
767 0 : if (plotMean) {
768 0 : cutHistoMean = AliTPCCalibViewer::SigmaCut(htemp, mean, sigma, sigmaMax, sigmaStep, pm);
769 0 : if (cutHistoMean) {
770 : //fListOfObjectsToBeDeleted->Add(cutHistoMean);
771 0 : cutHistoMean->SetLineColor(kRed);
772 0 : legend->AddEntry(cutHistoMean, "Mean", "l");
773 0 : cutHistoMean->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
774 0 : cutHistoMean->Draw();
775 0 : DrawLines(cutHistoMean, nsigma, legend, kRed, pm);
776 0 : } // if (cutHistoMean)
777 :
778 : }
779 0 : if (plotMedian) {
780 0 : cutHistoMedian = AliTPCCalibViewer::SigmaCut(htemp, median, sigma, sigmaMax, sigmaStep, pm);
781 0 : if (cutHistoMedian) {
782 : //fListOfObjectsToBeDeleted->Add(cutHistoMedian);
783 0 : cutHistoMedian->SetLineColor(kBlue);
784 0 : legend->AddEntry(cutHistoMedian, "Median", "l");
785 0 : cutHistoMedian->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
786 0 : if (plotMean && cutHistoMean) cutHistoMedian->Draw("same");
787 0 : else cutHistoMedian->Draw();
788 0 : DrawLines(cutHistoMedian, nsigma, legend, kBlue, pm);
789 0 : } // if (cutHistoMedian)
790 : }
791 0 : if (plotLTM) {
792 0 : Double_t ltmRms = 0;
793 0 : Double_t ltm = GetLTM(entries, values, <mRms, ltmFraction);
794 0 : cutHistoLTM = AliTPCCalibViewer::SigmaCut(htemp, ltm, ltmRms, sigmaMax, sigmaStep, pm);
795 0 : if (cutHistoLTM) {
796 : //fListOfObjectsToBeDeleted->Add(cutHistoLTM);
797 0 : cutHistoLTM->SetLineColor(kGreen+2);
798 0 : legend->AddEntry(cutHistoLTM, "LTM", "l");
799 0 : cutHistoLTM->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
800 0 : if ((plotMean && cutHistoMean) || (plotMedian && cutHistoMedian)) cutHistoLTM->Draw("same");
801 0 : else cutHistoLTM->Draw();
802 0 : DrawLines(cutHistoLTM, nsigma, legend, kGreen+2, pm);
803 0 : }
804 0 : }
805 0 : if (!plotMean && !plotMedian && !plotLTM) return -1;
806 0 : legend->Draw();
807 0 : return 1;
808 0 : }
809 :
810 :
811 : Int_t AliTPCCalibViewer::SigmaCutNew(const char* drawCommand, const char* sector, const char* cuts, Float_t /*sigmaMax*/, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, Bool_t /*pm*/, const char *sigmas, Float_t /*sigmaStep*/) const {
812 : /// Creates a histogram, where you can see, how much of the data are inside sigma-intervals
813 : /// around the mean/median/LTM
814 : /// with drawCommand, sector and cuts you specify your input data, see EasyDraw
815 : /// sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma)
816 : /// sigmaStep: the binsize of the generated histogram
817 : /// plotMean/plotMedian/plotLTM: specifies where to put the center
818 :
819 : // Double_t ltmFraction = 0.8; //unused
820 :
821 0 : TString drawStr(drawCommand);
822 0 : drawStr += " >> tempHist";
823 :
824 0 : Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
825 0 : TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
826 : TGraph *cutGraphMean = 0;
827 : // TGraph *cutGraphMedian = 0;
828 : // TGraph *cutGraphLTM = 0;
829 0 : Double_t *values = fTree->GetV1(); // value is the array containing 'entries' numbers
830 0 : Int_t *index = new Int_t[entries];
831 0 : Float_t *xarray = new Float_t[entries];
832 0 : Float_t *yarray = new Float_t[entries];
833 0 : TMath::Sort(entries, values, index, kFALSE);
834 :
835 0 : Double_t mean = TMath::Mean(entries, values);
836 : // Double_t median = TMath::Median(entries, values);
837 0 : Double_t sigma = TMath::RMS(entries, values);
838 :
839 0 : TLegend * legend = new TLegend(.7,.7, .99, .99, "Cumulative");
840 : //fListOfObjectsToBeDeleted->Add(legend);
841 :
842 : // parse sigmas string
843 0 : TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
844 0 : TVectorF nsigma(sigmasTokens->GetEntriesFast());
845 0 : for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
846 0 : TString str(((TObjString*)sigmasTokens->At(i))->GetString());
847 0 : Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
848 0 : nsigma[i] = sig;
849 0 : }
850 0 : delete sigmasTokens;
851 :
852 0 : if (plotMean) {
853 0 : for (Int_t i = 0; i < entries; i++) {
854 0 : xarray[i] = TMath::Abs(values[index[i]] - mean) / sigma;
855 0 : yarray[i] = float(i) / float(entries);
856 : }
857 0 : cutGraphMean = new TGraph(entries, xarray, yarray);
858 0 : if (cutGraphMean) {
859 : //fListOfObjectsToBeDeleted->Add(cutGraphMean);
860 0 : cutGraphMean->SetLineColor(kRed);
861 0 : legend->AddEntry(cutGraphMean, "Mean", "l");
862 0 : cutGraphMean->SetTitle(Form("%s, Cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
863 0 : cutGraphMean->Draw("alu");
864 0 : DrawLines(cutGraphMean, nsigma, legend, kRed, kTRUE);
865 0 : }
866 : }
867 0 : delete [] index;
868 0 : delete [] xarray;
869 0 : delete [] yarray;
870 : /*
871 : if (plotMedian) {
872 : cutHistoMedian = AliTPCCalibViewer::SigmaCut(htemp, median, sigma, sigmaMax, sigmaStep, pm);
873 : if (cutHistoMedian) {
874 : fListOfObjectsToBeDeleted->Add(cutHistoMedian);
875 : cutHistoMedian->SetLineColor(kBlue);
876 : legend->AddEntry(cutHistoMedian, "Median", "l");
877 : cutHistoMedian->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
878 : if (plotMean && cutHistoMean) cutHistoMedian->Draw("same");
879 : else cutHistoMedian->Draw();
880 : DrawLines(cutHistoMedian, nsigma, legend, kBlue, pm);
881 : } // if (cutHistoMedian)
882 : }
883 : if (plotLTM) {
884 : Double_t ltmRms = 0;
885 : Double_t ltm = GetLTM(entries, values, <mRms, ltmFraction);
886 : cutHistoLTM = AliTPCCalibViewer::SigmaCut(htemp, ltm, ltmRms, sigmaMax, sigmaStep, pm);
887 : if (cutHistoLTM) {
888 : fListOfObjectsToBeDeleted->Add(cutHistoLTM);
889 : cutHistoLTM->SetLineColor(kGreen+2);
890 : legend->AddEntry(cutHistoLTM, "LTM", "l");
891 : cutHistoLTM->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
892 : if (plotMean && cutHistoMean || plotMedian && cutHistoMedian) cutHistoLTM->Draw("same");
893 : else cutHistoLTM->Draw();
894 : DrawLines(cutHistoLTM, nsigma, legend, kGreen+2, pm);
895 : }
896 : }*/
897 0 : if (!plotMean && !plotMedian && !plotLTM) return -1;
898 0 : legend->Draw();
899 0 : return 1;
900 0 : }
901 :
902 :
903 : Int_t AliTPCCalibViewer::Integrate(const char* drawCommand, Int_t sector, const char* cuts, Float_t sigmaMax, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, const char *sigmas, Float_t sigmaStep) const {
904 : /// Creates an integrated histogram \f$S(t, \mu, \sigma)\f$, out of the input distribution distribution \f$f(x, \mu, \sigma)\f$, given in "histogram"
905 : /// "mean" and "sigma" are \f$\mu\f$ and \f$\sigma\f$ of the distribution in "histogram", to be specified by the user
906 : /// sigmaMax: up to which sigma around the mean/median/LTM you want to integrate
907 : /// if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
908 : /// "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
909 : /// The actual work is done on the array.
910 : ///
911 : /// \f[
912 : /// f(x, \mu, \sigma) \Rightarrow S(t, \mu, \sigma) = \int_{-\infty}^{\mu + t \sigma} f(x, \mu, \sigma) dx / \int_{-\infty}^{+\infty} f(x, \mu, \sigma) dx
913 : /// \f]
914 :
915 0 : if (sector >= 0 && sector < 72) {
916 0 : return Integrate(drawCommand, Form("%i", sector), cuts, sigmaMax, plotMean, plotMedian, plotLTM, sigmas, sigmaStep);
917 : }
918 0 : Error("Integrate","The TPC contains only sectors between 0 and 71.");
919 0 : return -1;
920 :
921 0 : }
922 :
923 :
924 : Int_t AliTPCCalibViewer::IntegrateOld(const char* drawCommand, const char* sector, const char* cuts, Float_t sigmaMax, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, const char *sigmas, Float_t sigmaStep) const {
925 : /// Creates an integrated histogram \f$S(t, \mu, \sigma)\f$, out of the input distribution distribution \f$f(x, \mu, \sigma)\f$, given in "histogram"
926 : /// "mean" and "sigma" are \f$\mu\f$ and \f$\sigma\f$ of the distribution in "histogram", to be specified by the user
927 : /// sigmaMax: up to which sigma around the mean/median/LTM you want to integrate
928 : /// if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
929 : /// "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
930 : /// The actual work is done on the array.
931 : ///
932 : /// \f[
933 : /// f(x, \mu, \sigma) \Rightarrow S(t, \mu, \sigma) = \int_{-\infty}^{\mu + t \sigma} f(x, \mu, \sigma) dx / \int_{-\infty}^{+\infty} f(x, \mu, \sigma) dx
934 : /// \f]
935 :
936 :
937 : Double_t ltmFraction = 0.8;
938 :
939 0 : TString drawStr(drawCommand);
940 0 : drawStr += " >> tempHist";
941 :
942 0 : Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
943 0 : TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
944 : // FIXME is this histogram deleted automatically?
945 0 : Double_t *values = fTree->GetV1(); // value is the array containing 'entries' numbers
946 :
947 0 : Double_t mean = TMath::Mean(entries, values);
948 0 : Double_t median = TMath::Median(entries, values);
949 0 : Double_t sigma = TMath::RMS(entries, values);
950 :
951 0 : TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
952 0 : TVectorF nsigma(sigmasTokens->GetEntriesFast());
953 0 : for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
954 0 : TString str(((TObjString*)sigmasTokens->At(i))->GetString());
955 0 : Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
956 0 : nsigma[i] = sig;
957 0 : }
958 0 : delete sigmasTokens;
959 :
960 0 : TLegend * legend = new TLegend(.7,.7, .99, .99, "Integrated histogram");
961 : //fListOfObjectsToBeDeleted->Add(legend);
962 : TH1F *integralHistoMean = 0;
963 : TH1F *integralHistoMedian = 0;
964 : TH1F *integralHistoLTM = 0;
965 :
966 0 : if (plotMean) {
967 0 : integralHistoMean = AliTPCCalibViewer::Integrate(htemp, mean, sigma, sigmaMax, sigmaStep);
968 0 : if (integralHistoMean) {
969 : //fListOfObjectsToBeDeleted->Add(integralHistoMean);
970 0 : integralHistoMean->SetLineColor(kRed);
971 0 : legend->AddEntry(integralHistoMean, "Mean", "l");
972 0 : integralHistoMean->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
973 0 : integralHistoMean->Draw();
974 0 : DrawLines(integralHistoMean, nsigma, legend, kRed, kTRUE);
975 0 : }
976 : }
977 0 : if (plotMedian) {
978 0 : integralHistoMedian = AliTPCCalibViewer::Integrate(htemp, median, sigma, sigmaMax, sigmaStep);
979 0 : if (integralHistoMedian) {
980 : //fListOfObjectsToBeDeleted->Add(integralHistoMedian);
981 0 : integralHistoMedian->SetLineColor(kBlue);
982 0 : legend->AddEntry(integralHistoMedian, "Median", "l");
983 0 : integralHistoMedian->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
984 0 : if (plotMean && integralHistoMean) integralHistoMedian->Draw("same");
985 0 : else integralHistoMedian->Draw();
986 0 : DrawLines(integralHistoMedian, nsigma, legend, kBlue, kTRUE);
987 0 : }
988 : }
989 0 : if (plotLTM) {
990 0 : Double_t ltmRms = 0;
991 0 : Double_t ltm = GetLTM(entries, values, <mRms, ltmFraction);
992 0 : integralHistoLTM = AliTPCCalibViewer::Integrate(htemp, ltm, ltmRms, sigmaMax, sigmaStep);
993 0 : if (integralHistoLTM) {
994 : //fListOfObjectsToBeDeleted->Add(integralHistoLTM);
995 0 : integralHistoLTM->SetLineColor(kGreen+2);
996 0 : legend->AddEntry(integralHistoLTM, "LTM", "l");
997 0 : integralHistoLTM->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
998 0 : if ((plotMean && integralHistoMean) || (plotMedian && integralHistoMedian)) integralHistoLTM->Draw("same");
999 0 : else integralHistoLTM->Draw();
1000 0 : DrawLines(integralHistoLTM, nsigma, legend, kGreen+2, kTRUE);
1001 0 : }
1002 0 : }
1003 0 : if (!plotMean && !plotMedian && !plotLTM) return -1;
1004 0 : legend->Draw();
1005 0 : return 1;
1006 0 : }
1007 :
1008 :
1009 : Int_t AliTPCCalibViewer::Integrate(const char* drawCommand, const char* sector, const char* cuts, Float_t /*sigmaMax*/, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, const char *sigmas, Float_t /*sigmaStep*/) const {
1010 : /// Creates an integrated histogram \f$S(t, \mu, \sigma)\f$, out of the input distribution distribution \f$f(x, \mu, \sigma)\f$, given in "histogram"
1011 : /// "mean" and "sigma" are \f$\mu\f$ and \f$\sigma\f$ of the distribution in "histogram", to be specified by the user
1012 : /// sigmaMax: up to which sigma around the mean/median/LTM you want to integrate
1013 : /// if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
1014 : /// "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
1015 : /// The actual work is done on the array.
1016 : ///
1017 : /// \f[
1018 : /// f(x, \mu, \sigma) \Rightarrow S(t, \mu, \sigma) = \int_{-\infty}^{\mu + t \sigma} f(x, \mu, \sigma) dx / \int_{-\infty}^{+\infty} f(x, \mu, \sigma) dx
1019 : /// \f]
1020 :
1021 : Double_t ltmFraction = 0.8;
1022 :
1023 0 : TString drawStr(drawCommand);
1024 0 : Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
1025 0 : if (dangerousToDraw) {
1026 0 : Warning("Integrate", "The draw string must not contain ':' or '>>'.");
1027 0 : return -1;
1028 : }
1029 0 : drawStr += " >> tempHist";
1030 :
1031 0 : Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
1032 0 : TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
1033 : TGraph *integralGraphMean = 0;
1034 : TGraph *integralGraphMedian = 0;
1035 : TGraph *integralGraphLTM = 0;
1036 0 : Double_t *values = fTree->GetV1(); // value is the array containing 'entries' numbers
1037 0 : Int_t *index = new Int_t[entries];
1038 0 : Float_t *xarray = new Float_t[entries];
1039 0 : Float_t *yarray = new Float_t[entries];
1040 0 : TMath::Sort(entries, values, index, kFALSE);
1041 :
1042 0 : Double_t mean = TMath::Mean(entries, values);
1043 0 : Double_t median = TMath::Median(entries, values);
1044 0 : Double_t sigma = TMath::RMS(entries, values);
1045 :
1046 : // parse sigmas string
1047 0 : TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
1048 0 : TVectorF nsigma(sigmasTokens->GetEntriesFast());
1049 0 : for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
1050 0 : TString str(((TObjString*)sigmasTokens->At(i))->GetString());
1051 0 : Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
1052 0 : nsigma[i] = sig;
1053 0 : }
1054 0 : delete sigmasTokens;
1055 :
1056 0 : TLegend * legend = new TLegend(.7,.7, .99, .99, "Integrated histogram");
1057 : //fListOfObjectsToBeDeleted->Add(legend);
1058 :
1059 0 : if (plotMean) {
1060 0 : for (Int_t i = 0; i < entries; i++) {
1061 0 : xarray[i] = (values[index[i]] - mean) / sigma;
1062 0 : yarray[i] = float(i) / float(entries);
1063 : }
1064 0 : integralGraphMean = new TGraph(entries, xarray, yarray);
1065 0 : if (integralGraphMean) {
1066 : //fListOfObjectsToBeDeleted->Add(integralGraphMean);
1067 0 : integralGraphMean->SetLineColor(kRed);
1068 0 : legend->AddEntry(integralGraphMean, "Mean", "l");
1069 0 : integralGraphMean->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
1070 0 : integralGraphMean->Draw("alu");
1071 0 : DrawLines(integralGraphMean, nsigma, legend, kRed, kTRUE);
1072 0 : }
1073 : }
1074 0 : if (plotMedian) {
1075 0 : for (Int_t i = 0; i < entries; i++) {
1076 0 : xarray[i] = (values[index[i]] - median) / sigma;
1077 0 : yarray[i] = float(i) / float(entries);
1078 : }
1079 0 : integralGraphMedian = new TGraph(entries, xarray, yarray);
1080 0 : if (integralGraphMedian) {
1081 : //fListOfObjectsToBeDeleted->Add(integralGraphMedian);
1082 0 : integralGraphMedian->SetLineColor(kBlue);
1083 0 : legend->AddEntry(integralGraphMedian, "Median", "l");
1084 0 : integralGraphMedian->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
1085 0 : if (plotMean && integralGraphMean) integralGraphMedian->Draw("samelu");
1086 0 : else integralGraphMedian->Draw("alu");
1087 0 : DrawLines(integralGraphMedian, nsigma, legend, kBlue, kTRUE);
1088 0 : }
1089 : }
1090 0 : if (plotLTM) {
1091 0 : Double_t ltmRms = 0;
1092 0 : Double_t ltm = GetLTM(entries, values, <mRms, ltmFraction);
1093 0 : for (Int_t i = 0; i < entries; i++) {
1094 0 : xarray[i] = (values[index[i]] - ltm) / ltmRms;
1095 0 : yarray[i] = float(i) / float(entries);
1096 : }
1097 0 : integralGraphLTM = new TGraph(entries, xarray, yarray);
1098 0 : if (integralGraphLTM) {
1099 : //fListOfObjectsToBeDeleted->Add(integralGraphLTM);
1100 0 : integralGraphLTM->SetLineColor(kGreen+2);
1101 0 : legend->AddEntry(integralGraphLTM, "LTM", "l");
1102 0 : integralGraphLTM->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
1103 0 : if ((plotMean && integralGraphMean) || (plotMedian && integralGraphMedian)) integralGraphLTM->Draw("samelu");
1104 0 : else integralGraphLTM->Draw("alu");
1105 0 : DrawLines(integralGraphLTM, nsigma, legend, kGreen+2, kTRUE);
1106 0 : }
1107 0 : }
1108 0 : delete [] index;
1109 0 : delete [] xarray;
1110 0 : delete [] yarray;
1111 :
1112 0 : if (!plotMean && !plotMedian && !plotLTM) return -1;
1113 0 : legend->Draw();
1114 0 : return entries;
1115 0 : }
1116 :
1117 :
1118 : void AliTPCCalibViewer::DrawLines(TH1F *histogram, TVectorF nsigma, TLegend *legend, Int_t color, Bool_t pm) const {
1119 : /// Private function for SigmaCut(...) and Integrate(...)
1120 : /// Draws lines into the given histogram, specified by "nsigma", the lines are addeed to the legend
1121 :
1122 : // start to draw the lines, loop over requested sigmas
1123 0 : for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
1124 0 : if (!pm) {
1125 0 : Int_t bin = histogram->GetXaxis()->FindBin(nsigma[i]);
1126 0 : TLine* lineUp = new TLine(nsigma[i], 0, nsigma[i], histogram->GetBinContent(bin));
1127 : //fListOfObjectsToBeDeleted->Add(lineUp);
1128 0 : lineUp->SetLineColor(color);
1129 0 : lineUp->SetLineStyle(2 + i);
1130 0 : lineUp->Draw();
1131 0 : TLine* lineLeft = new TLine(nsigma[i], histogram->GetBinContent(bin), 0, histogram->GetBinContent(bin));
1132 : //fListOfObjectsToBeDeleted->Add(lineLeft);
1133 0 : lineLeft->SetLineColor(color);
1134 0 : lineLeft->SetLineStyle(2 + i);
1135 0 : lineLeft->Draw();
1136 0 : legend->AddEntry(lineLeft, Form("Fraction(%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin)), "l");
1137 0 : }
1138 : else { // if (pm)
1139 0 : Int_t bin = histogram->GetXaxis()->FindBin(nsigma[i]);
1140 0 : TLine* lineUp1 = new TLine(nsigma[i], 0, nsigma[i], histogram->GetBinContent(bin));
1141 : //fListOfObjectsToBeDeleted->Add(lineUp1);
1142 0 : lineUp1->SetLineColor(color);
1143 0 : lineUp1->SetLineStyle(2 + i);
1144 0 : lineUp1->Draw();
1145 0 : TLine* lineLeft1 = new TLine(nsigma[i], histogram->GetBinContent(bin), histogram->GetBinLowEdge(0)+histogram->GetBinWidth(0), histogram->GetBinContent(bin));
1146 : //fListOfObjectsToBeDeleted->Add(lineLeft1);
1147 0 : lineLeft1->SetLineColor(color);
1148 0 : lineLeft1->SetLineStyle(2 + i);
1149 0 : lineLeft1->Draw();
1150 0 : legend->AddEntry(lineLeft1, Form("Fraction(+%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin)), "l");
1151 0 : bin = histogram->GetXaxis()->FindBin(-nsigma[i]);
1152 0 : TLine* lineUp2 = new TLine(-nsigma[i], 0, -nsigma[i], histogram->GetBinContent(bin));
1153 : //fListOfObjectsToBeDeleted->Add(lineUp2);
1154 0 : lineUp2->SetLineColor(color);
1155 0 : lineUp2->SetLineStyle(2 + i);
1156 0 : lineUp2->Draw();
1157 0 : TLine* lineLeft2 = new TLine(-nsigma[i], histogram->GetBinContent(bin), histogram->GetBinLowEdge(0)+histogram->GetBinWidth(0), histogram->GetBinContent(bin));
1158 : //fListOfObjectsToBeDeleted->Add(lineLeft2);
1159 0 : lineLeft2->SetLineColor(color);
1160 0 : lineLeft2->SetLineStyle(2 + i);
1161 0 : lineLeft2->Draw();
1162 0 : legend->AddEntry(lineLeft2, Form("Fraction(-%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin)), "l");
1163 : }
1164 : } // for (Int_t i = 0; i < nsigma.GetNoElements(); i++)
1165 0 : }
1166 :
1167 :
1168 : void AliTPCCalibViewer::DrawLines(TGraph *graph, TVectorF nsigma, TLegend *legend, Int_t color, Bool_t pm) const {
1169 : /// Private function for SigmaCut(...) and Integrate(...)
1170 : /// Draws lines into the given histogram, specified by "nsigma", the lines are addeed to the legend
1171 :
1172 : // start to draw the lines, loop over requested sigmas
1173 0 : for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
1174 0 : if (!pm) {
1175 0 : TLine* lineUp = new TLine(nsigma[i], 0, nsigma[i], graph->Eval(nsigma[i]));
1176 : //fListOfObjectsToBeDeleted->Add(lineUp);
1177 0 : lineUp->SetLineColor(color);
1178 0 : lineUp->SetLineStyle(2 + i);
1179 0 : lineUp->Draw();
1180 0 : TLine* lineLeft = new TLine(nsigma[i], graph->Eval(nsigma[i]), 0, graph->Eval(nsigma[i]));
1181 : //fListOfObjectsToBeDeleted->Add(lineLeft);
1182 0 : lineLeft->SetLineColor(color);
1183 0 : lineLeft->SetLineStyle(2 + i);
1184 0 : lineLeft->Draw();
1185 0 : legend->AddEntry(lineLeft, Form("Fraction(%f #sigma) = %f",nsigma[i], graph->Eval(nsigma[i])), "l");
1186 0 : }
1187 : else { // if (pm)
1188 0 : TLine* lineUp1 = new TLine(nsigma[i], 0, nsigma[i], graph->Eval(nsigma[i]));
1189 : //fListOfObjectsToBeDeleted->Add(lineUp1);
1190 0 : lineUp1->SetLineColor(color);
1191 0 : lineUp1->SetLineStyle(2 + i);
1192 0 : lineUp1->Draw();
1193 0 : TLine* lineLeft1 = new TLine(nsigma[i], graph->Eval(nsigma[i]), graph->GetHistogram()->GetXaxis()->GetBinLowEdge(0), graph->Eval(nsigma[i]));
1194 : //fListOfObjectsToBeDeleted->Add(lineLeft1);
1195 0 : lineLeft1->SetLineColor(color);
1196 0 : lineLeft1->SetLineStyle(2 + i);
1197 0 : lineLeft1->Draw();
1198 0 : legend->AddEntry(lineLeft1, Form("Fraction(+%f #sigma) = %f",nsigma[i], graph->Eval(nsigma[i])), "l");
1199 0 : TLine* lineUp2 = new TLine(-nsigma[i], 0, -nsigma[i], graph->Eval(-nsigma[i]));
1200 : //fListOfObjectsToBeDeleted->Add(lineUp2);
1201 0 : lineUp2->SetLineColor(color);
1202 0 : lineUp2->SetLineStyle(2 + i);
1203 0 : lineUp2->Draw();
1204 0 : TLine* lineLeft2 = new TLine(-nsigma[i], graph->Eval(-nsigma[i]), graph->GetHistogram()->GetXaxis()->GetBinLowEdge(0), graph->Eval(-nsigma[i]));
1205 : //fListOfObjectsToBeDeleted->Add(lineLeft2);
1206 0 : lineLeft2->SetLineColor(color);
1207 0 : lineLeft2->SetLineStyle(2 + i);
1208 0 : lineLeft2->Draw();
1209 0 : legend->AddEntry(lineLeft2, Form("Fraction(-%f #sigma) = %f",nsigma[i], graph->Eval(-nsigma[i])), "l");
1210 : }
1211 : } // for (Int_t i = 0; i < nsigma.GetNoElements(); i++)
1212 0 : }
1213 :
1214 :
1215 :
1216 :
1217 :
1218 : /////////////////
1219 : // Array tools //
1220 : /////////////////
1221 :
1222 :
1223 : Int_t AliTPCCalibViewer::GetBin(Float_t value, Int_t nbins, Double_t binLow, Double_t binUp){
1224 : /// Returns the 'bin' for 'value'
1225 : /// The interval between 'binLow' and 'binUp' is divided into 'nbins' equidistant bins
1226 : /// avoid index out of bounds error: 'if (bin < binLow) bin = binLow' and vice versa
1227 : ///
1228 : /// \f[
1229 : /// GetBin(value) = \frac{nbins - 1}{binUp - binLow} \upoint (value - binLow) +1
1230 : /// \f]
1231 :
1232 0 : Int_t bin = TMath::Nint( (Float_t)(value - binLow) / (Float_t)(binUp - binLow) * (nbins-1) ) + 1;
1233 : // avoid index out of bounds:
1234 0 : if (value < binLow) bin = 0;
1235 0 : if (value > binUp) bin = nbins + 1;
1236 0 : return bin;
1237 :
1238 : }
1239 :
1240 :
1241 : Double_t AliTPCCalibViewer::GetLTM(Int_t n, const Double_t *const array, Double_t *const sigma, Double_t fraction){
1242 : /// returns the LTM and sigma
1243 :
1244 0 : Double_t *ddata = new Double_t[n];
1245 0 : Double_t mean = 0, lsigma = 0;
1246 : UInt_t nPoints = 0;
1247 0 : for (UInt_t i = 0; i < (UInt_t)n; i++) {
1248 0 : ddata[nPoints]= array[nPoints];
1249 0 : nPoints++;
1250 : }
1251 0 : Int_t hh = TMath::Min(TMath::Nint(fraction * nPoints), Int_t(n));
1252 0 : AliMathBase::EvaluateUni(nPoints, ddata, mean, lsigma, hh);
1253 0 : if (sigma) *sigma = lsigma;
1254 0 : delete [] ddata;
1255 0 : return mean;
1256 0 : }
1257 :
1258 :
1259 : TH1F* AliTPCCalibViewer::SigmaCut(TH1F *const histogram, Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep, Bool_t pm) {
1260 : /// Creates a cumulative histogram \f$S(t, \mu, \sigma)\f$, where you can see, how much of the data are inside sigma-intervals around the mean value
1261 : /// The data of the distribution \f$f(x, \mu, \sigma)\f$ are given in 'histogram'
1262 : /// 'mean' and 'sigma' are \f$\mu\f$ and \f$\sigma\f$ of the distribution in 'histogram', to be specified by the user
1263 : /// sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, \f$t \sigma\f$)
1264 : /// sigmaStep: the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
1265 : /// pm: Decide weather \f$t > 0\f$ (first case) or \f$t\f$ arbitrary (secound case)
1266 : /// The actual work is done on the array.
1267 : ///
1268 : /// \f[
1269 : /// f(x, \mu, \sigma) \Rightarrow S(t, \mu, \sigma) = (\int_{\mu}^{\mu + t \sigma} f(x, \mu, \sigma) dx + \int_{\mu}^{\mu - t \sigma} f(x, \mu, \sigma) dx) / (\int_{-\infty}^{+\infty} f(x, \mu, \sigma) dx), for t > 0
1270 : /// or
1271 : /// f(x, \mu, \sigma) \Rightarrow S(t, \mu, \sigma) = \int_{\mu}^{\mu + t \sigma} f(x, \mu, \sigma) dx / \int_{-\infty}^{+\infty} f(x, \mu, \sigma) dx
1272 : /// \f]
1273 : /// 
1274 :
1275 0 : Float_t *array = histogram->GetArray();
1276 0 : Int_t nbins = histogram->GetXaxis()->GetNbins();
1277 0 : Float_t binLow = histogram->GetXaxis()->GetXmin();
1278 0 : Float_t binUp = histogram->GetXaxis()->GetXmax();
1279 0 : return AliTPCCalibViewer::SigmaCut(nbins, array, mean, sigma, nbins, binLow, binUp, sigmaMax, sigmaStep, pm);
1280 : }
1281 :
1282 :
1283 : TH1F* AliTPCCalibViewer::SigmaCut(Int_t n, const Float_t *array, Float_t mean, Float_t sigma, Int_t nbins, Float_t binLow, Float_t binUp, Float_t sigmaMax, Float_t sigmaStep, Bool_t pm){
1284 : /// Creates a histogram \f$S(t, \mu, \sigma)\f$, where you can see, how much of the data are inside sigma-intervals around the mean value
1285 : /// The data of the distribution \f$f(x, \mu, \sigma)\f$ are given in 'array', 'n' specifies the length of the array
1286 : /// 'mean' and 'sigma' are \f$\mu\f$ and \f$\sigma\f$ of the distribution in 'array', to be specified by the user
1287 : /// 'nbins': number of bins, 'binLow': first bin, 'binUp': last bin
1288 : /// sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, \f$t \sigma\f$)
1289 : /// sigmaStep: the binsize of the generated histogram
1290 : /// Here the actual work is done.
1291 :
1292 0 : if (sigma == 0) return 0;
1293 0 : Float_t binWidth = (binUp-binLow)/(nbins - 1);
1294 0 : if (sigmaStep <= 0) sigmaStep = binWidth;
1295 0 : Int_t kbins = (Int_t)(sigmaMax * sigma / sigmaStep) + 1; // + 1 due to overflow bin in histograms
1296 0 : if (pm) kbins = 2 * (Int_t)(sigmaMax * sigma / sigmaStep) + 1;
1297 0 : Float_t kbinLow = !pm ? 0 : -sigmaMax;
1298 : Float_t kbinUp = sigmaMax;
1299 0 : TH1F *hist = new TH1F("sigmaCutHisto","Cumulative; Multiples of #sigma; Fraction of included data", kbins, kbinLow, kbinUp);
1300 0 : hist->SetDirectory(0);
1301 0 : hist->Reset();
1302 :
1303 : // calculate normalization
1304 : Double_t normalization = 0;
1305 0 : for (Int_t i = 0; i <= n; i++) {
1306 0 : normalization += array[i];
1307 : }
1308 :
1309 : // given units: units from given histogram
1310 : // sigma units: in units of sigma
1311 : // iDelta: integrate in interval (mean +- iDelta), given units
1312 : // x: ofset from mean for integration, given units
1313 : // hist: needs
1314 :
1315 : // printf("nbins: %i, binLow: %f, binUp: %f \n", nbins, binLow, binUp);
1316 : // fill histogram
1317 0 : for (Float_t iDelta = 0; iDelta <= sigmaMax * sigma; iDelta += sigmaStep) {
1318 : // integrate array
1319 0 : Double_t valueP = array[GetBin(mean, nbins, binLow, binUp)];
1320 0 : Double_t valueM = array[GetBin(mean-binWidth, nbins, binLow, binUp)];
1321 : // add bin of mean value only once to the histogram
1322 : // printf("++ adding bins: ");
1323 0 : for (Float_t x = binWidth; x <= iDelta; x += binWidth) {
1324 0 : valueP += (mean + x <= binUp) ? array[GetBin(mean + x, nbins, binLow, binUp)] : 0;
1325 0 : valueM += (mean-binWidth - x >= binLow) ? array[GetBin(mean-binWidth - x, nbins, binLow, binUp)] : 0;
1326 : // printf("%i, ", GetBin(mean + x, nbins, binLow, binUp));
1327 : }
1328 : // printf("\n");
1329 0 : if (valueP / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail +++ \n", valueP, normalization);
1330 0 : if (valueP / normalization > 100) return hist;
1331 0 : if (valueM / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail +++ \n", valueM, normalization);
1332 0 : if (valueM / normalization > 100) return hist;
1333 : valueP = (valueP / normalization);
1334 : valueM = (valueM / normalization);
1335 0 : if (pm) {
1336 0 : Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
1337 0 : hist->SetBinContent(bin, valueP);
1338 0 : bin = GetBin(-iDelta/sigma, kbins, kbinLow, kbinUp);
1339 0 : hist->SetBinContent(bin, valueM);
1340 0 : }
1341 : else { // if (!pm)
1342 0 : Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
1343 0 : hist->SetBinContent(bin, valueP + valueM);
1344 : // printf(" first integration bin: %i, last integration bin in + direction: %i \n", GetBin(mean+binWidth, nbins, binLow, binUp), GetBin(iDelta, nbins, binLow, binUp));
1345 : // printf(" first integration bin: %i, last integration bin in - direction: %i \n", GetBin(mean+binWidth, nbins, binLow, binUp), GetBin(-iDelta, nbins, binLow, binUp));
1346 : // printf(" value: %f, normalization: %f, iDelta: %f, Bin: %i \n", valueP+valueM, normalization, iDelta, bin);
1347 : }
1348 0 : }
1349 : //hist->SetMaximum(0.7);
1350 0 : if (!pm) hist->SetMaximum(1.2);
1351 0 : return hist;
1352 0 : }
1353 :
1354 :
1355 : TH1F* AliTPCCalibViewer::SigmaCut(Int_t /*n*/, const Double_t */*array*/, Double_t /*mean*/, Double_t /*sigma*/, Int_t /*nbins*/, const Double_t */*xbins*/, Double_t /*sigmaMax*/){
1356 : /// SigmaCut for variable binsize
1357 : /// NOT YET IMPLEMENTED !!!
1358 :
1359 0 : printf("SigmaCut with variable binsize, Not yet implemented\n");
1360 :
1361 0 : return 0;
1362 : }
1363 :
1364 :
1365 : TH1F* AliTPCCalibViewer::Integrate(TH1F *const histogram, Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep){
1366 : /// Creates an integrated histogram \f$S(t, \mu, \sigma)\f$, out of the input distribution distribution \f$f(x, \mu, \sigma)\f$, given in "histogram"
1367 : /// "mean" and "sigma" are \f$\mu\f$ and \f$\sigma\f$ of the distribution in "histogram", to be specified by the user
1368 : /// sigmaMax: up to which sigma around the mean/median/LTM you want to integrate
1369 : /// if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
1370 : /// "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
1371 : /// The actual work is done on the array.
1372 : /// \f[
1373 : /// f(x, \mu, \sigma) \Rightarrow S(t, \mu, \sigma) = \int_{-\infty}^{\mu + t \sigma} f(x, \mu, \sigma) dx / \int_{-\infty}^{+\infty} f(x, \mu, \sigma) dx
1374 : /// \f]
1375 : /// 
1376 :
1377 :
1378 0 : Float_t *array = histogram->GetArray();
1379 0 : Int_t nbins = histogram->GetXaxis()->GetNbins();
1380 0 : Float_t binLow = histogram->GetXaxis()->GetXmin();
1381 0 : Float_t binUp = histogram->GetXaxis()->GetXmax();
1382 0 : return AliTPCCalibViewer::Integrate(nbins, array, nbins, binLow, binUp, mean, sigma, sigmaMax, sigmaStep);
1383 : }
1384 :
1385 :
1386 : TH1F* AliTPCCalibViewer::Integrate(Int_t n, const Float_t *const array, Int_t nbins, Float_t binLow, Float_t binUp, Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep){
1387 : /// Creates an integrated histogram \f$S(t, \mu, \sigma)\f$, out of the input distribution distribution \f$f(x, \mu, \sigma)\f$, given in "histogram"
1388 : /// "mean" and "sigma" are \f$\mu\f$ and \f$\sigma\f$ of the distribution in "histogram", to be specified by the user
1389 : /// sigmaMax: up to which sigma around the mean/median/LTM you want to integrate
1390 : /// if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
1391 : /// "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
1392 : /// Here the actual work is done.
1393 :
1394 : Bool_t givenUnits = kTRUE;
1395 0 : if (sigma != 0 && sigmaMax != 0) givenUnits = kFALSE;
1396 0 : if (givenUnits) {
1397 : sigma = 1;
1398 0 : sigmaMax = (binUp - binLow) / 2.;
1399 0 : }
1400 :
1401 0 : Float_t binWidth = (binUp-binLow)/(nbins - 1);
1402 0 : if (sigmaStep <= 0) sigmaStep = binWidth;
1403 0 : Int_t kbins = (Int_t)(sigmaMax * sigma / sigmaStep) + 1; // + 1 due to overflow bin in histograms
1404 0 : Float_t kbinLow = givenUnits ? binLow : -sigmaMax;
1405 0 : Float_t kbinUp = givenUnits ? binUp : sigmaMax;
1406 : TH1F *hist = 0;
1407 0 : if (givenUnits) hist = new TH1F("integratedHisto","Integrated Histogram; Given x; Fraction of included data", kbins, kbinLow, kbinUp);
1408 0 : if (!givenUnits) hist = new TH1F("integratedHisto","Integrated Histogram; Multiples of #sigma; Fraction of included data", kbins, kbinLow, kbinUp);
1409 0 : hist->SetDirectory(0);
1410 0 : hist->Reset();
1411 :
1412 : // calculate normalization
1413 : // printf("calculating normalization, integrating from bin 1 to %i \n", n);
1414 : Double_t normalization = 0;
1415 0 : for (Int_t i = 1; i <= n; i++) {
1416 0 : normalization += array[i];
1417 : }
1418 : // printf("normalization: %f \n", normalization);
1419 :
1420 : // given units: units from given histogram
1421 : // sigma units: in units of sigma
1422 : // iDelta: integrate in interval (mean +- iDelta), given units
1423 : // x: ofset from mean for integration, given units
1424 : // hist: needs
1425 :
1426 : // fill histogram
1427 0 : for (Float_t iDelta = mean - sigmaMax * sigma; iDelta <= mean + sigmaMax * sigma; iDelta += sigmaStep) {
1428 : // integrate array
1429 : Double_t value = 0;
1430 0 : for (Float_t x = mean - sigmaMax * sigma; x <= iDelta; x += binWidth) {
1431 0 : value += (x <= binUp && x >= binLow) ? array[GetBin(x, nbins, binLow, binUp)] : 0;
1432 : }
1433 0 : if (value / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail +++ \n", value, normalization);
1434 0 : if (value / normalization > 100) return hist;
1435 0 : Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
1436 : // printf("first integration bin: %i, last integration bin: %i \n", GetBin(mean - sigmaMax * sigma, nbins, binLow, binUp), GetBin(iDelta, nbins, binLow, binUp));
1437 : // printf("value: %f, normalization: %f, normalized value: %f, iDelta: %f, Bin: %i \n", value, normalization, value/normalization, iDelta, bin);
1438 : value = (value / normalization);
1439 0 : hist->SetBinContent(bin, value);
1440 0 : }
1441 0 : return hist;
1442 0 : }
1443 :
1444 :
1445 :
1446 :
1447 :
1448 : ////////////////////////
1449 : // end of Array tools //
1450 : ////////////////////////
1451 :
1452 :
1453 :
1454 : //_____________________________________________________________________________
1455 : AliTPCCalPad* AliTPCCalibViewer::GetCalPadOld(const char* desiredData, const char* cuts, const char* calPadName) const {
1456 : /// creates a AliTPCCalPad out of the 'desiredData'
1457 : /// the functionality of EasyDraw1D is used
1458 : /// calPadName specifies the name of the created AliTPCCalPad
1459 : /// - this takes a while -
1460 :
1461 0 : TString drawStr(desiredData);
1462 0 : drawStr.Append(":channel");
1463 0 : drawStr.Append(fAbbreviation);
1464 0 : AliTPCCalPad * createdCalPad = new AliTPCCalPad(calPadName, calPadName);
1465 : Int_t entries = 0;
1466 0 : for (Int_t sec = 0; sec < 72; sec++) {
1467 0 : AliTPCCalROC * roc = createdCalPad->GetCalROC(sec);
1468 0 : entries = EasyDraw1D(drawStr.Data(), (Int_t)sec, cuts, "goff");
1469 0 : if (entries == -1) return 0;
1470 0 : const Double_t *pchannel = fTree->GetV2();
1471 0 : const Double_t *pvalue = fTree->GetV1();
1472 0 : for (Int_t i = 0; i < entries; i++)
1473 0 : roc->SetValue((UInt_t)(pchannel[i]), (Float_t)(pvalue[i]));
1474 0 : }
1475 0 : return createdCalPad;
1476 0 : }
1477 :
1478 :
1479 : //_____________________________________________________________________________
1480 : AliTPCCalPad* AliTPCCalibViewer::GetCalPad(const char* desiredData, const char* cuts, const char* calPadName) const {
1481 : /// creates a AliTPCCalPad out of the 'desiredData'
1482 : /// the functionality of EasyDraw1D is used
1483 : /// calPadName specifies the name of the created AliTPCCalPad
1484 : /// - this takes a while -
1485 :
1486 0 : TString drawStr(desiredData);
1487 0 : drawStr.Append(":channel.fElements:sector");
1488 0 : AliTPCCalPad * createdCalPad = new AliTPCCalPad(calPadName, calPadName);
1489 : //
1490 0 : Int_t entries = fTree->Draw(drawStr, cuts,"goff");
1491 0 : const Double_t *pvalue = fTree->GetV1();
1492 0 : const Double_t *pchannel = fTree->GetV2();
1493 0 : const Double_t *psector = fTree->GetV3();
1494 :
1495 0 : for (Int_t ientry=0; ientry<entries; ientry++){
1496 0 : Int_t sector= TMath::Nint(psector[ientry]);
1497 0 : AliTPCCalROC * roc = createdCalPad->GetCalROC(sector);
1498 0 : if (roc) roc->SetValue((UInt_t)(pchannel[ientry]), (Float_t)(pvalue[ientry]));
1499 : }
1500 :
1501 : // for (Int_t sec = 0; sec < 72; sec++) {
1502 : // AliTPCCalROC * roc = createdCalPad->GetCalROC(sec);
1503 : // entries = EasyDraw1D(drawStr.Data(), (Int_t)sec, cuts, "goff");
1504 : // if (entries == -1) return 0;
1505 : // for (Int_t i = 0; i < entries; i++)
1506 : // roc->SetValue((UInt_t)(pchannel[i]), (Float_t)(pvalue[i]));
1507 : // }
1508 : return createdCalPad;
1509 0 : }
1510 :
1511 : //_____________________________________________________________________________
1512 : AliTPCCalROC* AliTPCCalibViewer::GetCalROC(const char* desiredData, UInt_t sector, const char* cuts) const {
1513 : /// creates a AliTPCCalROC out of the desiredData
1514 : /// the functionality of EasyDraw1D is used
1515 : /// sector specifies the sector of the created AliTPCCalROC
1516 :
1517 0 : TString drawStr(desiredData);
1518 0 : drawStr.Append(":channel");
1519 0 : drawStr.Append(fAbbreviation);
1520 0 : Int_t entries = EasyDraw1D(drawStr.Data(), (Int_t)sector, cuts, "goff");
1521 0 : if (entries == -1) return 0;
1522 0 : AliTPCCalROC * createdROC = new AliTPCCalROC(sector);
1523 0 : for (Int_t i = 0; i < entries; i++)
1524 0 : createdROC->SetValue((UInt_t)(fTree->GetV2()[i]), fTree->GetV1()[i]);
1525 : return createdROC;
1526 0 : }
1527 :
1528 :
1529 : TObjArray* AliTPCCalibViewer::GetListOfVariables(Bool_t printList) {
1530 : /// scan the tree - produces a list of available variables in the tree
1531 : /// printList: print the list to the screen, after the scan is done
1532 :
1533 0 : TObjArray* arr = new TObjArray();
1534 : TObjString* str = 0;
1535 0 : if (!fTree) return 0;
1536 0 : Int_t nentries = fTree->GetListOfBranches()->GetEntries();
1537 0 : for (Int_t i = 0; i < nentries; i++) {
1538 0 : str = new TObjString(fTree->GetListOfBranches()->At(i)->GetName());
1539 0 : str->String().ReplaceAll("_Median", "");
1540 0 : str->String().ReplaceAll("_Mean", "");
1541 0 : str->String().ReplaceAll("_RMS", "");
1542 0 : str->String().ReplaceAll("_LTM", "");
1543 0 : str->String().ReplaceAll("_OutlierCutted", "");
1544 0 : str->String().ReplaceAll(".", "");
1545 0 : if (!arr->FindObject(str) &&
1546 0 : !(str->String() == "channel" || str->String() == "gx" || str->String() == "gy" ||
1547 0 : str->String() == "lx" || str->String() == "ly" || str->String() == "pad" ||
1548 0 : str->String() == "row" || str->String() == "rpad" || str->String() == "sector" ))
1549 0 : arr->Add(str);
1550 : }
1551 :
1552 : // loop over all friends (if there are some) and add them to the list
1553 0 : if (fTree->GetListOfFriends()) {
1554 0 : for (Int_t ifriend = 0; ifriend < fTree->GetListOfFriends()->GetEntries(); ifriend++){
1555 : // printf("iterating through friendlist, currently at %i\n", ifriend);
1556 : // printf("working with %s\n", fTree->GetListOfFriends()->At(ifriend)->ClassName());
1557 0 : if (TString(fTree->GetListOfFriends()->At(ifriend)->ClassName()) != "TFriendElement") continue; // no friendElement found
1558 0 : TFriendElement *friendElement = (TFriendElement*)fTree->GetListOfFriends()->At(ifriend);
1559 0 : if (friendElement->GetTree() == 0) continue; // no tree found in friendElement
1560 : // printf("friend found \n");
1561 0 : for (Int_t i = 0; i < friendElement->GetTree()->GetListOfBranches()->GetEntries(); i++) {
1562 : // printf("iterating through friendelement entries, currently at %i\n", i);
1563 0 : str = new TObjString(friendElement->GetTree()->GetListOfBranches()->At(i)->GetName());
1564 0 : str->String().ReplaceAll("_Median", "");
1565 0 : str->String().ReplaceAll("_Mean", "");
1566 0 : str->String().ReplaceAll("_RMS", "");
1567 0 : str->String().ReplaceAll("_LTM", "");
1568 0 : str->String().ReplaceAll("_OutlierCutted", "");
1569 0 : str->String().ReplaceAll(".", "");
1570 0 : if (!(str->String() == "channel" || str->String() == "gx" || str->String() == "gy" ||
1571 0 : str->String() == "lx" || str->String() == "ly" || str->String() == "pad" ||
1572 0 : str->String() == "row" || str->String() == "rpad" || str->String() == "sector" )){
1573 : // insert "<friendName>." at the beginning: (<friendName> is per default "R")
1574 0 : str->String().Insert(0, ".");
1575 0 : str->String().Insert(0, friendElement->GetName());
1576 0 : if (!arr->FindObject(str)) arr->Add(str);
1577 : // printf("added string %s \n", str->String().Data());
1578 : }
1579 : }
1580 0 : }
1581 0 : } // if (fTree->GetListOfFriends())
1582 :
1583 0 : arr->Sort();
1584 : // ((TFriendElement*)gui->GetViewer()->GetTree()->GetListOfFriends()->At(0))->GetTree()->GetListOfBranches()->At(0)->GetName()
1585 : // ((TFriendElement*)gui->GetViewer()->GetTree()->GetListOfFriends()->At(0))->GetTree()->GetListOfBranches()
1586 :
1587 :
1588 0 : if (printList) {
1589 0 : TIterator* iter = arr->MakeIterator();
1590 0 : iter->Reset();
1591 : TObjString* currentStr = 0;
1592 0 : while ( (currentStr = (TObjString*)(iter->Next())) ) {
1593 0 : std::cout << currentStr->GetString().Data() << std::endl;
1594 : }
1595 0 : delete iter;
1596 0 : }
1597 : return arr;
1598 0 : }
1599 :
1600 :
1601 : TObjArray* AliTPCCalibViewer::GetListOfNormalizationVariables(Bool_t printList) const{
1602 : /// produces a list of available variables for normalization in the tree
1603 : /// printList: print the list to the screen, after the scan is done
1604 :
1605 0 : TObjArray* arr = new TObjArray();
1606 0 : arr->Add(new TObjString("_Mean"));
1607 0 : arr->Add(new TObjString("_Mean_OutlierCutted"));
1608 0 : arr->Add(new TObjString("_Median"));
1609 0 : arr->Add(new TObjString("_Median_OutlierCutted"));
1610 0 : arr->Add(new TObjString("_LTM"));
1611 0 : arr->Add(new TObjString("_LTM_OutlierCutted"));
1612 0 : arr->Add(new TObjString(Form("LFitIntern_4_8%s", fAppendString.Data())));
1613 0 : arr->Add(new TObjString(Form("GFitIntern_Lin%s", fAppendString.Data())));
1614 0 : arr->Add(new TObjString(Form("GFitIntern_Par%s", fAppendString.Data())));
1615 0 : arr->Add(new TObjString("FitLinLocal"));
1616 0 : arr->Add(new TObjString("FitLinGlobal"));
1617 0 : arr->Add(new TObjString("FitParLocal"));
1618 0 : arr->Add(new TObjString("FitParGlobal"));
1619 :
1620 0 : if (printList) {
1621 0 : TIterator* iter = arr->MakeIterator();
1622 0 : iter->Reset();
1623 : TObjString* currentStr = 0;
1624 0 : while ((currentStr = (TObjString*)(iter->Next()))) {
1625 0 : std::cout << currentStr->GetString().Data() << std::endl;
1626 : }
1627 0 : delete iter;
1628 0 : }
1629 0 : return arr;
1630 0 : }
1631 :
1632 :
1633 : TFriendElement* AliTPCCalibViewer::AddReferenceTree(const char* filename, const char* treename, const char* refname){
1634 : /// add a reference tree to the current tree
1635 : /// by default the treename is 'calPads' and the reference treename is 'R'
1636 :
1637 0 : TFile *file = new TFile(filename);
1638 0 : fListOfObjectsToBeDeleted->Add(file);
1639 0 : TTree * tree = (TTree*)file->Get(treename);
1640 0 : return AddFriend(tree, refname);
1641 0 : }
1642 :
1643 :
1644 : TObjArray* AliTPCCalibViewer::GetArrayOfCalPads(){
1645 : /// Returns a TObjArray with all AliTPCCalPads that are stored in the tree
1646 : /// - this takes a while -
1647 :
1648 0 : TObjArray *listOfCalPads = GetListOfVariables();
1649 0 : TObjArray *calPadsArray = new TObjArray();
1650 0 : Int_t numberOfCalPads = listOfCalPads->GetEntries();
1651 0 : for (Int_t i = 0; i < numberOfCalPads; i++) {
1652 0 : std::cout << "Creating calPad " << (i+1) << " of " << numberOfCalPads << "\r" << std::flush;
1653 0 : char* calPadName = (char*)((TObjString*)(listOfCalPads->At(i)))->GetString().Data();
1654 0 : TString drawCommand = ((TObjString*)(listOfCalPads->At(i)))->GetString();
1655 0 : drawCommand.Append(fAbbreviation.Data());
1656 0 : AliTPCCalPad* calPad = GetCalPad(drawCommand.Data(), "", calPadName);
1657 0 : calPadsArray->Add(calPad);
1658 0 : }
1659 0 : std::cout << std::endl;
1660 0 : listOfCalPads->Delete();
1661 0 : delete listOfCalPads;
1662 0 : return calPadsArray;
1663 0 : }
1664 :
1665 :
1666 : TString* AliTPCCalibViewer::Fit(const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, TVectorD &fitParam, TMatrixD &covMatrix){
1667 : /// fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
1668 : /// returns chi2, fitParam and covMatrix
1669 : /// returns TString with fitted formula
1670 :
1671 0 : TString formulaStr(formula);
1672 0 : TString drawStr(drawCommand);
1673 0 : TString cutStr(cuts);
1674 :
1675 : // abbreviations:
1676 0 : drawStr.ReplaceAll(fAbbreviation, fAppendString);
1677 0 : cutStr.ReplaceAll(fAbbreviation, fAppendString);
1678 0 : formulaStr.ReplaceAll(fAbbreviation, fAppendString);
1679 :
1680 0 : formulaStr.ReplaceAll("++", fAbbreviation);
1681 0 : TObjArray* formulaTokens = formulaStr.Tokenize(fAbbreviation.Data());
1682 0 : Int_t dim = formulaTokens->GetEntriesFast();
1683 :
1684 0 : fitParam.ResizeTo(dim);
1685 0 : covMatrix.ResizeTo(dim,dim);
1686 :
1687 0 : TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
1688 0 : fitter->StoreData(kTRUE);
1689 0 : fitter->ClearPoints();
1690 :
1691 0 : Int_t entries = Draw(drawStr.Data(), cutStr.Data(), "goff");
1692 0 : if (entries == -1) {
1693 0 : delete formulaTokens;
1694 0 : return new TString("An ERROR has occured during fitting!");
1695 : }
1696 0 : Double_t **values = new Double_t*[dim+1] ;
1697 :
1698 0 : for (Int_t i = 0; i < dim + 1; i++){
1699 : Int_t centries = 0;
1700 0 : if (i < dim) centries = fTree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff");
1701 0 : else centries = fTree->Draw(drawStr.Data(), cutStr.Data(), "goff");
1702 :
1703 0 : if (entries != centries) {
1704 0 : delete [] values;
1705 0 : delete formulaTokens;
1706 0 : return new TString("An ERROR has occured during fitting!");
1707 : }
1708 0 : values[i] = new Double_t[entries];
1709 0 : memcpy(values[i], fTree->GetV1(), entries*sizeof(Double_t));
1710 0 : }
1711 :
1712 : // add points to the fitter
1713 0 : for (Int_t i = 0; i < entries; i++){
1714 0 : Double_t x[1000];
1715 0 : for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
1716 0 : fitter->AddPoint(x, values[dim][i], 1);
1717 0 : }
1718 :
1719 0 : fitter->Eval();
1720 0 : fitter->GetParameters(fitParam);
1721 0 : fitter->GetCovarianceMatrix(covMatrix);
1722 0 : chi2 = fitter->GetChisquare();
1723 : // chi2 = chi2;
1724 :
1725 0 : TString *preturnFormula = new TString(Form("( %e+",fitParam[0])), &returnFormula = *preturnFormula;
1726 :
1727 0 : for (Int_t iparam = 0; iparam < dim; iparam++) {
1728 0 : returnFormula.Append(Form("%s*(%e)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]));
1729 0 : if (iparam < dim-1) returnFormula.Append("+");
1730 : }
1731 0 : returnFormula.Append(" )");
1732 0 : delete formulaTokens;
1733 0 : delete fitter;
1734 0 : for (Int_t i = 0; i < dim + 1; i++){
1735 0 : delete [] values[i];
1736 : }
1737 0 : delete [] values;
1738 : return preturnFormula;
1739 0 : }
1740 :
1741 :
1742 : void AliTPCCalibViewer::MakeTreeWithObjects(const char *fileName, const TObjArray *const array, const char * mapFileName) {
1743 : /// Write tree with all available information
1744 : /// im mapFileName is speciefied, the Map information are also written to the tree
1745 : /// AliTPCCalPad-Objects are written directly to the tree, so that they can be accessd later on
1746 : /// (does not work!!!)
1747 :
1748 0 : AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
1749 :
1750 : TObjArray* mapIROCs = 0;
1751 : TObjArray* mapOROCs = 0;
1752 : TVectorF *mapIROCArray = 0;
1753 : TVectorF *mapOROCArray = 0;
1754 : Int_t mapEntries = 0;
1755 : TString* mapNames = 0;
1756 :
1757 0 : if (mapFileName) {
1758 0 : TFile mapFile(mapFileName, "read");
1759 :
1760 0 : TList* listOfROCs = mapFile.GetListOfKeys();
1761 0 : mapEntries = listOfROCs->GetEntries()/2;
1762 0 : mapIROCs = new TObjArray(mapEntries*2);
1763 0 : mapOROCs = new TObjArray(mapEntries*2);
1764 0 : mapIROCArray = new TVectorF[mapEntries];
1765 0 : mapOROCArray = new TVectorF[mapEntries];
1766 :
1767 0 : mapNames = new TString[mapEntries];
1768 0 : for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
1769 0 : TString rocName(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
1770 0 : rocName.Remove(rocName.Length()-4, 4);
1771 0 : mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((rocName + "IROC").Data()), ivalue);
1772 0 : mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((rocName + "OROC").Data()), ivalue);
1773 0 : mapNames[ivalue].Append(rocName);
1774 0 : }
1775 :
1776 0 : for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
1777 0 : mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
1778 0 : mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
1779 :
1780 0 : for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
1781 0 : (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
1782 0 : for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
1783 0 : (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
1784 : }
1785 :
1786 0 : } // if (mapFileName)
1787 :
1788 0 : TTreeSRedirector cstream(fileName);
1789 0 : Int_t arrayEntries = array->GetEntries();
1790 :
1791 : // Read names of AliTPCCalPads and save them in names[]
1792 0 : TString* names = new TString[arrayEntries];
1793 0 : for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
1794 0 : names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
1795 :
1796 0 : for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
1797 :
1798 0 : TVectorF *vectorArray = new TVectorF[arrayEntries];
1799 0 : for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
1800 0 : vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
1801 :
1802 :
1803 : //
1804 : // fill vectors of variable per pad
1805 : //
1806 0 : TVectorF *posArray = new TVectorF[8];
1807 0 : for (Int_t ivalue = 0; ivalue < 8; ivalue++)
1808 0 : posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
1809 :
1810 0 : Float_t posG[3] = {0};
1811 0 : Float_t posL[3] = {0};
1812 : Int_t ichannel = 0;
1813 0 : for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
1814 0 : for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
1815 0 : tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
1816 0 : tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
1817 0 : posArray[0][ichannel] = irow;
1818 0 : posArray[1][ichannel] = ipad;
1819 0 : posArray[2][ichannel] = posL[0];
1820 0 : posArray[3][ichannel] = posL[1];
1821 0 : posArray[4][ichannel] = posG[0];
1822 0 : posArray[5][ichannel] = posG[1];
1823 0 : posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
1824 0 : posArray[7][ichannel] = ichannel;
1825 :
1826 : // loop over array containing AliTPCCalPads
1827 0 : for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
1828 0 : AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
1829 0 : AliTPCCalROC* calROC = calPad->GetCalROC(isector);
1830 0 : if (calROC)
1831 0 : (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
1832 : else
1833 0 : (vectorArray[ivalue])[ichannel] = 0;
1834 : }
1835 0 : ichannel++;
1836 : }
1837 : }
1838 0 : AliTPCCalROC dummyROC(0);
1839 0 : for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
1840 0 : AliTPCCalROC *roc = ((AliTPCCalPad*)array->At(ivalue))->GetCalROC(isector);
1841 0 : if (!roc) roc = &dummyROC;
1842 0 : cstream << "calPads" <<
1843 0 : (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
1844 0 : cstream << "calPads" <<
1845 0 : (Char_t*)((names[ivalue] + "Pad.=").Data()) << roc;
1846 : }
1847 :
1848 0 : if (mapFileName) {
1849 0 : for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
1850 0 : if (isector < 36)
1851 0 : cstream << "calPads" <<
1852 0 : (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
1853 : else
1854 0 : cstream << "calPads" <<
1855 0 : (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
1856 : }
1857 0 : }
1858 :
1859 0 : cstream << "calPads" <<
1860 0 : "sector=" << isector;
1861 :
1862 0 : cstream << "calPads" <<
1863 0 : "row.=" << &posArray[0] <<
1864 0 : "pad.=" << &posArray[1] <<
1865 0 : "lx.=" << &posArray[2] <<
1866 0 : "ly.=" << &posArray[3] <<
1867 0 : "gx.=" << &posArray[4] <<
1868 0 : "gy.=" << &posArray[5] <<
1869 0 : "rpad.=" << &posArray[6] <<
1870 0 : "channel.=" << &posArray[7];
1871 :
1872 0 : cstream << "calPads" <<
1873 : "\n";
1874 :
1875 0 : delete[] posArray;
1876 0 : delete[] vectorArray;
1877 0 : } //for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++)
1878 :
1879 0 : delete[] names;
1880 0 : if (mapFileName) {
1881 0 : delete mapIROCs;
1882 0 : delete mapOROCs;
1883 0 : delete[] mapIROCArray;
1884 0 : delete[] mapOROCArray;
1885 0 : delete[] mapNames;
1886 : }
1887 0 : }
1888 :
1889 :
1890 : void AliTPCCalibViewer::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad *const outlierPad, Float_t ltmFraction) {
1891 : /// Write a tree with all available information
1892 : /// if mapFileName is speciefied, the Map information are also written to the tree
1893 : /// pads specified in outlierPad are not used for calculating statistics
1894 : /// The following statistical information on the basis of a ROC are calculated:
1895 : /// "_Median", "_Mean", "_LTM", "_RMS_LTM"
1896 : /// "_Median_OutlierCutted", "_Mean_OutlierCutted", "_RMS_OutlierCutted", "_LTM_OutlierCutted", "_RMS_LTM_OutlierCutted"
1897 : /// The following position variables are available:
1898 : /// "row", "pad", "lx", "ly", "gx", "gy", "rpad", "channel"
1899 : ///
1900 : /// The tree out of this function is the basis for the AliTPCCalibViewer and the AliTPCCalibViewerGUI.
1901 :
1902 0 : AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
1903 :
1904 : TObjArray* mapIROCs = 0;
1905 : TObjArray* mapOROCs = 0;
1906 : TVectorF *mapIROCArray = 0;
1907 : TVectorF *mapOROCArray = 0;
1908 : Int_t mapEntries = 0;
1909 : TString* mapNames = 0;
1910 :
1911 0 : if (mapFileName) {
1912 0 : TFile mapFile(mapFileName, "read");
1913 :
1914 0 : TList* listOfROCs = mapFile.GetListOfKeys();
1915 0 : mapEntries = listOfROCs->GetEntries()/2;
1916 0 : mapIROCs = new TObjArray(mapEntries*2);
1917 0 : mapOROCs = new TObjArray(mapEntries*2);
1918 0 : mapIROCArray = new TVectorF[mapEntries];
1919 0 : mapOROCArray = new TVectorF[mapEntries];
1920 :
1921 0 : mapNames = new TString[mapEntries];
1922 0 : for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
1923 0 : TString rocName(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
1924 0 : rocName.Remove(rocName.Length()-4, 4);
1925 0 : mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((rocName + "IROC").Data()), ivalue);
1926 0 : mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((rocName + "OROC").Data()), ivalue);
1927 0 : mapNames[ivalue].Append(rocName);
1928 0 : }
1929 :
1930 0 : for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
1931 0 : mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
1932 0 : mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
1933 :
1934 0 : for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
1935 0 : (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
1936 0 : for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
1937 0 : (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
1938 : }
1939 :
1940 0 : } // if (mapFileName)
1941 :
1942 0 : TTreeSRedirector cstream(fileName,"recreate");
1943 : Int_t arrayEntries = 0;
1944 0 : if (array) arrayEntries = array->GetEntries();
1945 :
1946 0 : TString* names = new TString[arrayEntries];
1947 0 : for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
1948 0 : names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
1949 :
1950 0 : for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
1951 : //
1952 : // get statistic for given sector
1953 : //
1954 0 : TVectorF median(arrayEntries);
1955 0 : TVectorF mean(arrayEntries);
1956 0 : TVectorF rms(arrayEntries);
1957 0 : TVectorF ltm(arrayEntries);
1958 0 : TVectorF ltmrms(arrayEntries);
1959 0 : TVectorF medianWithOut(arrayEntries);
1960 0 : TVectorF meanWithOut(arrayEntries);
1961 0 : TVectorF rmsWithOut(arrayEntries);
1962 0 : TVectorF ltmWithOut(arrayEntries);
1963 0 : TVectorF ltmrmsWithOut(arrayEntries);
1964 :
1965 0 : TVectorF *vectorArray = new TVectorF[arrayEntries];
1966 0 : for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++){
1967 0 : vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
1968 0 : vectorArray[ivalue].SetUniqueID(array->UncheckedAt(ivalue)->GetUniqueID());
1969 : // printf("UniqueID: %d\n",vectorArray[ivalue].GetUniqueID());
1970 : }
1971 :
1972 0 : for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
1973 0 : AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
1974 0 : AliTPCCalROC* calROC = calPad->GetCalROC(isector);
1975 : AliTPCCalROC* outlierROC = 0;
1976 0 : if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
1977 0 : if (calROC) {
1978 0 : median[ivalue] = calROC->GetMedian();
1979 0 : mean[ivalue] = calROC->GetMean();
1980 0 : rms[ivalue] = calROC->GetRMS();
1981 0 : Double_t ltmrmsValue = 0;
1982 0 : ltm[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction);
1983 0 : ltmrms[ivalue] = ltmrmsValue;
1984 0 : if (outlierROC) {
1985 0 : medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
1986 0 : meanWithOut[ivalue] = calROC->GetMean(outlierROC);
1987 0 : rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
1988 0 : ltmrmsValue = 0;
1989 0 : ltmWithOut[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction, outlierROC);
1990 0 : ltmrmsWithOut[ivalue] = ltmrmsValue;
1991 0 : }
1992 0 : }
1993 : else {
1994 0 : median[ivalue] = 0.;
1995 0 : mean[ivalue] = 0.;
1996 0 : rms[ivalue] = 0.;
1997 0 : ltm[ivalue] = 0.;
1998 0 : ltmrms[ivalue] = 0.;
1999 0 : medianWithOut[ivalue] = 0.;
2000 0 : meanWithOut[ivalue] = 0.;
2001 0 : rmsWithOut[ivalue] = 0.;
2002 0 : ltmWithOut[ivalue] = 0.;
2003 0 : ltmrmsWithOut[ivalue] = 0.;
2004 : }
2005 : }
2006 :
2007 : //
2008 : // fill vectors of variable per pad
2009 : //
2010 0 : TVectorF *posArray = new TVectorF[8];
2011 0 : for (Int_t ivalue = 0; ivalue < 8; ivalue++)
2012 0 : posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
2013 :
2014 0 : Float_t posG[3] = {0};
2015 0 : Float_t posL[3] = {0};
2016 : Int_t ichannel = 0;
2017 0 : for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
2018 0 : for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
2019 0 : tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
2020 0 : tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
2021 0 : posArray[0][ichannel] = irow;
2022 0 : posArray[1][ichannel] = ipad;
2023 0 : posArray[2][ichannel] = posL[0];
2024 0 : posArray[3][ichannel] = posL[1];
2025 0 : posArray[4][ichannel] = posG[0];
2026 0 : posArray[5][ichannel] = posG[1];
2027 0 : posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
2028 0 : posArray[7][ichannel] = ichannel;
2029 :
2030 : // loop over array containing AliTPCCalPads
2031 0 : for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
2032 0 : AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
2033 0 : AliTPCCalROC* calROC = calPad->GetCalROC(isector);
2034 0 : if (calROC)
2035 0 : (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
2036 : else
2037 0 : (vectorArray[ivalue])[ichannel] = 0;
2038 : }
2039 0 : ichannel++;
2040 : }
2041 : }
2042 :
2043 0 : cstream << "calPads" <<
2044 0 : "sector=" << isector;
2045 :
2046 0 : for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
2047 0 : if (outlierPad==0) {
2048 0 : cstream << "calPads" <<
2049 0 : (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
2050 0 : (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
2051 0 : (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
2052 0 : (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
2053 0 : (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
2054 0 : }
2055 0 : if (outlierPad) {
2056 0 : cstream << "calPads" <<
2057 0 : (Char_t*)((names[ivalue] + "_Median=").Data()) << medianWithOut[ivalue] <<
2058 0 : (Char_t*)((names[ivalue] + "_Mean").Data()) << meanWithOut[ivalue] <<
2059 0 : (Char_t*)((names[ivalue] + "_RMS=").Data()) << rmsWithOut[ivalue] <<
2060 0 : (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltmWithOut[ivalue] <<
2061 0 : (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrmsWithOut[ivalue];
2062 0 : }
2063 : //timestamp and run, if given in title
2064 : /* TString title(((AliTPCCalPad*) array->At(ivalue))->GetTitle());
2065 : TObjArray *arrtitle=title.Tokenize(",");
2066 : Int_t run=-1;
2067 : UInt_t time=0;
2068 : TIter next(arrtitle);
2069 : TObject *o=0;
2070 : while ( (o=next()) ){
2071 : TString &entry=((TObjString*)o)->GetString();
2072 : entry.Remove(TString::kBoth,' ');
2073 : entry.Remove(TString::kBoth,'\t');
2074 : if (entry.BeginsWith("Run:")) {
2075 : run=entry(4,entry.Length()).Atoi();
2076 : } else if (entry.BeginsWith("Time:")) {
2077 : time=entry(6,entry.Length()).Atoi();
2078 : }
2079 : }
2080 : delete arrtitle;*/
2081 :
2082 : }
2083 :
2084 0 : for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
2085 0 : cstream << "calPads" <<
2086 0 : (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
2087 : }
2088 :
2089 0 : if (mapFileName) {
2090 0 : for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
2091 0 : if (isector < 36)
2092 0 : cstream << "calPads" <<
2093 0 : (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
2094 : else
2095 0 : cstream << "calPads" <<
2096 0 : (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
2097 : }
2098 0 : }
2099 :
2100 0 : cstream << "calPads" <<
2101 0 : "row.=" << &posArray[0] <<
2102 0 : "pad.=" << &posArray[1] <<
2103 0 : "lx.=" << &posArray[2] <<
2104 0 : "ly.=" << &posArray[3] <<
2105 0 : "gx.=" << &posArray[4] <<
2106 0 : "gy.=" << &posArray[5] <<
2107 0 : "rpad.=" << &posArray[6] <<
2108 0 : "channel.=" << &posArray[7];
2109 :
2110 0 : cstream << "calPads" <<
2111 : "\n";
2112 :
2113 0 : delete[] posArray;
2114 0 : delete[] vectorArray;
2115 0 : }
2116 0 : TTree * tree = (cstream << "calPads").GetTree();
2117 0 : MakeCalPadAliases(tree);
2118 0 : for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
2119 : //(Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
2120 0 : tree->SetAlias(TString::Format("%s_LTMRatio",names[ivalue].Data()).Data(), TString::Format("%s.fElements/%s_LTM",names[ivalue].Data(),names[ivalue].Data()).Data());
2121 0 : tree->SetAlias(TString::Format("%s_MeanRatio",names[ivalue].Data()).Data(), TString::Format("%s.fElements/%s_Mean",names[ivalue].Data(),names[ivalue].Data()).Data());
2122 0 : tree->SetAlias(TString::Format("%s_MedianRatio",names[ivalue].Data()).Data(), TString::Format("%s.fElements/%s_Median",names[ivalue].Data(),names[ivalue].Data()).Data());
2123 0 : tree->SetAlias(names[ivalue].Data(), TString::Format("%s.fElements",names[ivalue].Data()).Data());
2124 : }
2125 0 : delete[] names;
2126 0 : if (mapFileName) {
2127 0 : delete mapIROCs;
2128 0 : delete mapOROCs;
2129 0 : delete[] mapIROCArray;
2130 0 : delete[] mapOROCArray;
2131 0 : delete[] mapNames;
2132 : }
2133 0 : }
2134 : void AliTPCCalibViewer::MakeCalPadAliases(TTree * tree){
2135 : /// make standard aliases for standard calibration trees
2136 :
2137 0 : tree->SetAlias("Pad0","MapPadLength.fElements==7.5"); // pad types
2138 0 : tree->SetAlias("Pad1","MapPadLength.fElements==10.0");
2139 0 : tree->SetAlias("Pad2","MapPadLength.fElements==15.0");
2140 0 : tree->SetAlias("ly","(ly.fElements+((sector<36)*(2*((sector%36)<18)-1)*0.2)+((sector>=36)*(2*((sector%36)<18)-1)*0.3))"); // correction for bug in tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
2141 0 : tree->SetAlias("dRowNorm0","(row.fElements/32-1)"); // normalized distance to the center of the chamber in lx (-1,1)
2142 0 : tree->SetAlias("dRowNorm1","(row.fElements/32-1)"); //
2143 0 : tree->SetAlias("dRowNorm2","((row.fElements-63)/16-1)"); //
2144 0 : tree->SetAlias("dRowNorm","(row.fElements<63)*(row.fElements/32-1)+(row.fElements>63)*((row.fElements-63)/16-1)");
2145 0 : tree->SetAlias("dPhiNorm","(ly/lx.fElements)/(pi/18)"); // normalized distance to the center of the chamber in phi (-1,1)
2146 : //
2147 0 : tree->SetAlias("fsector","(sector%36)+(0.5+9*atan2(ly,lx.fElements)/pi)"); // float sector number
2148 0 : tree->SetAlias("posEdge","((pi/18.)-abs(atan(ly/lx.fElements)))*lx.fElements"); //distance to the edge
2149 0 : tree->SetAlias("lphi","atan(ly/lx.fElements)"); //local phi angle
2150 0 : }
2151 :
2152 : void AliTPCCalibViewer::MakeTree(const char *outPutFileName, const Char_t *inputFileName, AliTPCCalPad *outlierPad, Float_t ltmFraction, const char *mapFileName ){
2153 : /// Function to create a calibration Tree with all available information.
2154 : /// See also documentation to MakeTree
2155 : /// the file "inputFileName" must be in the following format (see also CreateObjectList):
2156 : /// (each colum separated by tabs, "dependingOnType" can have 2 or 3 colums)
2157 : ///
2158 : /// type path dependingOnType
2159 : ///
2160 : /// type == "CE":
2161 : /// dependingOnType = CETmean CEQmean CETrms
2162 : ///
2163 : /// type == "Pulser":
2164 : /// dependingOnType = PulserTmean PulsterQmean PulserTrms
2165 : ///
2166 : /// type == "Pedestals":
2167 : /// dependingOnType = Pedestals Noise
2168 : ///
2169 : /// type == "CalPad":
2170 : /// dependingOnType = NameInFile NameToWriteToFile
2171 :
2172 0 : TObjArray objArray;
2173 0 : CreateObjectList(inputFileName, &objArray);
2174 0 : MakeTree(outPutFileName, &objArray, mapFileName, outlierPad, ltmFraction);
2175 0 : }
2176 :
2177 :
2178 : void AliTPCCalibViewer::CreateObjectList(const Char_t *filename, TObjArray *calibObjects){
2179 : /// Function to create a TObjArray out of a given file
2180 : /// the file must be in the following format:
2181 : /// (each colum separated by tabs, "dependingOnType" can have 2 or 3 colums)
2182 : ///
2183 : /// type path dependingOnType
2184 : ///
2185 : /// type == "CE":
2186 : /// dependingOnType = CETmean CEQmean CETrms
2187 : ///
2188 : /// type == "Pulser":
2189 : /// dependingOnType = PulserTmean PulsterQmean PulserTrms
2190 : ///
2191 : /// type == "Pedestals":
2192 : /// dependingOnType = Pedestals Noise
2193 : ///
2194 : /// type == "CalPad":
2195 : /// dependingOnType = NameInFile NameToWriteToFile
2196 :
2197 0 : if ( calibObjects == 0x0 ) return;
2198 0 : ifstream in;
2199 0 : in.open(filename);
2200 0 : if ( !in.is_open() ){
2201 0 : fprintf(stderr,"Error: cannot open list file '%s'", filename);
2202 0 : return;
2203 : }
2204 :
2205 : AliTPCCalPad *calPad=0x0;
2206 :
2207 0 : TString sFile;
2208 0 : sFile.ReadFile(in);
2209 0 : in.close();
2210 :
2211 0 : TObjArray *arrFileLine = sFile.Tokenize("\n");
2212 0 : TIter nextLine(arrFileLine);
2213 :
2214 : TObjString *sObjLine = 0x0;
2215 0 : while ( (sObjLine = (TObjString*)nextLine()) ){
2216 0 : TString sLine(sObjLine->GetString());
2217 :
2218 0 : TObjArray *arrCol = sLine.Tokenize("\t");
2219 0 : Int_t nCols = arrCol->GetEntriesFast();
2220 :
2221 0 : TObjString *sObjType = (TObjString*)(arrCol->At(0));
2222 0 : TObjString *sObjFileName = (TObjString*)(arrCol->At(1));
2223 : TObjString *sObjName = 0x0;
2224 :
2225 0 : if ( !sObjType || !sObjFileName ) continue;
2226 0 : TString sType(sObjType->GetString());
2227 0 : TString sFileName(sObjFileName->GetString());
2228 0 : printf("Type %s, opening %s \n", sType.Data(), sFileName.Data());
2229 0 : TFile *fIn = TFile::Open(sFileName);
2230 0 : if ( !fIn ){
2231 0 : fprintf(stderr,"File not found: '%s'", sFileName.Data());
2232 0 : continue;
2233 : }
2234 :
2235 0 : if ( sType == "CE" ){ // next three colums are the names for CETmean, CEQmean and CETrms
2236 0 : AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
2237 0 : calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());
2238 0 : if (nCols > 2) { // check, if the name is provided
2239 0 : sObjName = (TObjString*)(arrCol->At(2));
2240 0 : calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2241 0 : }
2242 0 : else calPad->SetNameTitle("CETmean","CETmean");
2243 0 : calibObjects->Add(calPad);
2244 :
2245 0 : calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());
2246 0 : if (nCols > 3) { // check, if the name is provided
2247 0 : sObjName = (TObjString*)(arrCol->At(3));
2248 0 : calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2249 0 : }
2250 0 : else calPad->SetNameTitle("CEQmean","CEQmean");
2251 0 : calibObjects->Add(calPad);
2252 :
2253 0 : calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
2254 0 : if (nCols > 4) { // check, if the name is provided
2255 0 : sObjName = (TObjString*)(arrCol->At(4));
2256 0 : calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2257 0 : }
2258 0 : else calPad->SetNameTitle("CETrms","CETrms");
2259 0 : calibObjects->Add(calPad);
2260 :
2261 0 : } else if ( sType == "Pulser") {
2262 0 : AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
2263 :
2264 0 : calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());
2265 0 : if (nCols > 2) {
2266 0 : sObjName = (TObjString*)(arrCol->At(2));
2267 0 : calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2268 0 : }
2269 0 : else calPad->SetNameTitle("PulserTmean","PulserTmean");
2270 0 : calibObjects->Add(calPad);
2271 :
2272 0 : calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());
2273 0 : if (nCols > 3) {
2274 0 : sObjName = (TObjString*)(arrCol->At(3));
2275 0 : calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2276 0 : }
2277 0 : else calPad->SetNameTitle("PulserQmean","PulserQmean");
2278 0 : calibObjects->Add(calPad);
2279 :
2280 0 : calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
2281 0 : if (nCols > 4) {
2282 0 : sObjName = (TObjString*)(arrCol->At(4));
2283 0 : calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2284 0 : }
2285 0 : else calPad->SetNameTitle("PulserTrms","PulserTrms");
2286 0 : calibObjects->Add(calPad);
2287 :
2288 0 : } else if ( sType == "Pedestals") {
2289 0 : AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
2290 :
2291 0 : calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());
2292 0 : if (nCols > 2) {
2293 0 : sObjName = (TObjString*)(arrCol->At(2));
2294 0 : calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2295 0 : }
2296 0 : else calPad->SetNameTitle("Pedestals","Pedestals");
2297 0 : calibObjects->Add(calPad);
2298 :
2299 0 : calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());
2300 0 : if (nCols > 3) {
2301 0 : sObjName = (TObjString*)(arrCol->At(3));
2302 0 : calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2303 0 : }
2304 0 : else calPad->SetNameTitle("Noise","Noise");
2305 0 : calibObjects->Add(calPad);
2306 :
2307 0 : } else if ( sType == "CalPad") {
2308 0 : if (nCols > 2) sObjName = (TObjString*)(arrCol->At(2));
2309 0 : else continue;
2310 0 : calPad = new AliTPCCalPad(*(AliTPCCalPad*)fIn->Get(sObjName->GetString().Data()));
2311 0 : if (nCols > 3) {
2312 0 : sObjName = (TObjString*)(arrCol->At(3));
2313 0 : calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2314 0 : }
2315 0 : calibObjects->Add(calPad);
2316 : } else {
2317 0 : fprintf(stderr,"Undefined Type: '%s'",sType.Data());
2318 : }
2319 0 : delete fIn;
2320 0 : }
2321 0 : delete arrFileLine;
2322 0 : }
|