Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGraphBentErrors.cxx
Go to the documentation of this file.
1// @(#)root/hist:$Id$
2// Author: Dave Morrison 30/06/2003
3
4/*************************************************************************
5 * Copyright (C) 1995-2003, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#include <cstring>
13#include <iostream>
14
15#include "TROOT.h"
16#include "TGraphBentErrors.h"
17#include "TMath.h"
18#include "TVirtualPad.h"
19#include "TH1.h"
20#include "TF1.h"
21
23
24
25////////////////////////////////////////////////////////////////////////////////
26
27/** \class TGraphBentErrors
28 \ingroup Graphs
29A TGraphBentErrors is a TGraph with bent, asymmetric error bars.
30
31The TGraphBentErrors painting is performed thanks to the TGraphPainter
32class. All details about the various painting options are given in this class.
33
34The picture below gives an example:
35Begin_Macro(source)
36{
37 auto c1 = new TCanvas("c1","A Simple Graph with bent error bars",200,10,700,500);
38 const Int_t n = 10;
39 Double_t x[n] = {-0.22, 0.05, 0.25, 0.35, 0.5, 0.61,0.7,0.85,0.89,0.95};
40 Double_t y[n] = {1,2.9,5.6,7.4,9,9.6,8.7,6.3,4.5,1};
41 Double_t exl[n] = {.05,.1,.07,.07,.04,.05,.06,.07,.08,.05};
42 Double_t eyl[n] = {.8,.7,.6,.5,.4,.4,.5,.6,.7,.8};
43 Double_t exh[n] = {.02,.08,.05,.05,.03,.03,.04,.05,.06,.03};
44 Double_t eyh[n] = {.6,.5,.4,.3,.2,.2,.3,.4,.5,.6};
45 Double_t exld[n] = {.0,.0,.0,.0,.0,.0,.0,.0,.0,.0};
46 Double_t eyld[n] = {.0,.0,.05,.0,.0,.0,.0,.0,.0,.0};
47 Double_t exhd[n] = {.0,.0,.0,.0,.0,.0,.0,.0,.0,.0};
48 Double_t eyhd[n] = {.0,.0,.0,.0,.0,.0,.0,.0,.05,.0};
49 auto gr = new TGraphBentErrors(n,x,y,exl,exh,eyl,eyh,exld,exhd,eyld,eyhd);
50 gr->SetTitle("TGraphBentErrors Example");
51 gr->SetMarkerColor(4);
52 gr->SetMarkerStyle(21);
53 gr->Draw("ALP");
54}
55End_Macro
56*/
57
58
59////////////////////////////////////////////////////////////////////////////////
60/// TGraphBentErrors default constructor.
61
63{
64 if (!CtorAllocate()) return;
65}
66
67
68////////////////////////////////////////////////////////////////////////////////
69/// TGraphBentErrors copy constructor
70
72 : TGraph(gr)
73{
74 if (!CtorAllocate()) return;
75 Int_t n = fNpoints*sizeof(Double_t);
76 memcpy(fEXlow, gr.fEXlow, n);
77 memcpy(fEYlow, gr.fEYlow, n);
78 memcpy(fEXhigh, gr.fEXhigh, n);
79 memcpy(fEYhigh, gr.fEYhigh, n);
80 memcpy(fEXlowd, gr.fEXlowd, n);
81 memcpy(fEYlowd, gr.fEYlowd, n);
82 memcpy(fEXhighd, gr.fEXhighd, n);
83 memcpy(fEYhighd, gr.fEYhighd, n);
84}
85
86
87////////////////////////////////////////////////////////////////////////////////
88/// TGraphBentErrors normal constructor.
89///
90/// the arrays are preset to zero
91
93 : TGraph(n)
94{
95 if (!CtorAllocate()) return;
97}
98
99
100////////////////////////////////////////////////////////////////////////////////
101/// TGraphBentErrors normal constructor.
102///
103/// if exl,h or eyl,h are null, the corresponding arrays are preset to zero
104
106 const Float_t *x, const Float_t *y,
107 const Float_t *exl, const Float_t *exh,
108 const Float_t *eyl, const Float_t *eyh,
109 const Float_t *exld, const Float_t *exhd,
110 const Float_t *eyld, const Float_t *eyhd)
111 : TGraph(n,x,y)
112{
113 if (!CtorAllocate()) return;
114
115 for (Int_t i=0;i<n;i++) {
116 if (exl) fEXlow[i] = exl[i];
117 else fEXlow[i] = 0;
118 if (exh) fEXhigh[i] = exh[i];
119 else fEXhigh[i] = 0;
120 if (eyl) fEYlow[i] = eyl[i];
121 else fEYlow[i] = 0;
122 if (eyh) fEYhigh[i] = eyh[i];
123 else fEYhigh[i] = 0;
124
125 if (exld) fEXlowd[i] = exld[i];
126 else fEXlowd[i] = 0;
127 if (exhd) fEXhighd[i] = exhd[i];
128 else fEXhighd[i] = 0;
129 if (eyld) fEYlowd[i] = eyld[i];
130 else fEYlowd[i] = 0;
131 if (eyhd) fEYhighd[i] = eyhd[i];
132 else fEYhighd[i] = 0;
133 }
134}
135
136
137////////////////////////////////////////////////////////////////////////////////
138/// TGraphBentErrors normal constructor.
139///
140/// if exl,h or eyl,h are null, the corresponding arrays are preset to zero
141
143 const Double_t *x, const Double_t *y,
144 const Double_t *exl, const Double_t *exh,
145 const Double_t *eyl, const Double_t *eyh,
146 const Double_t *exld, const Double_t *exhd,
147 const Double_t *eyld, const Double_t *eyhd)
148 : TGraph(n,x,y)
149{
150 if (!CtorAllocate()) return;
151 auto memsz = sizeof(Double_t)*fNpoints;
152
153 if (exl) memcpy(fEXlow, exl, memsz);
154 else memset(fEXlow, 0, memsz);
155 if (exh) memcpy(fEXhigh, exh, memsz);
156 else memset(fEXhigh, 0, memsz);
157 if (eyl) memcpy(fEYlow, eyl, memsz);
158 else memset(fEYlow, 0, memsz);
159 if (eyh) memcpy(fEYhigh, eyh, memsz);
160 else memset(fEYhigh, 0, memsz);
161
162 if (exld) memcpy(fEXlowd, exld, memsz);
163 else memset(fEXlowd, 0, memsz);
164 if (exhd) memcpy(fEXhighd, exhd, memsz);
165 else memset(fEXhighd, 0, memsz);
166 if (eyld) memcpy(fEYlowd, eyld, memsz);
167 else memset(fEYlowd, 0, memsz);
168 if (eyhd) memcpy(fEYhighd, eyhd, memsz);
169 else memset(fEYhighd, 0, memsz);
170}
171
172
173////////////////////////////////////////////////////////////////////////////////
174/// TGraphBentErrors default destructor.
175
177{
178 delete [] fEXlow;
179 delete [] fEXhigh;
180 delete [] fEYlow;
181 delete [] fEYhigh;
182
183 delete [] fEXlowd;
184 delete [] fEXhighd;
185 delete [] fEYlowd;
186 delete [] fEYhighd;
187}
188
189////////////////////////////////////////////////////////////////////////////////
190/// Add a point with bent errors to the graph.
191
193 Double_t exld, Double_t exhd, Double_t eyld, Double_t eyhd)
194{
195 AddPoint(x, y);
196 SetPointError(fNpoints - 1, exl, exh, eyl, eyh, exld, exhd, eyld, eyhd);
197}
198
199////////////////////////////////////////////////////////////////////////////////
200/// Apply a function to all data points \f$ y = f(x,y) \f$.
201///
202/// Errors are calculated as \f$ eyh = f(x,y+eyh)-f(x,y) \f$ and
203/// \f$ eyl = f(x,y)-f(x,y-eyl) \f$.
204///
205/// Special treatment has to be applied for the functions where the
206/// role of "up" and "down" is reversed.
207///
208/// Function suggested/implemented by Miroslav Helbich <helbich@mail.desy.de>
209
211{
212 Double_t x,y,exl,exh,eyl,eyh,eyl_new,eyh_new,fxy;
213
214 if (fHistogram) {
215 delete fHistogram;
216 fHistogram = nullptr;
217 }
218 for (Int_t i = 0; i < GetN(); i++) {
219 GetPoint(i, x, y);
220 exl = GetErrorXlow(i);
221 exh = GetErrorXhigh(i);
222 eyl = GetErrorYlow(i);
223 eyh = GetErrorYhigh(i);
224
225 fxy = f->Eval(x, y);
226 SetPoint(i, x, fxy);
227
228 // in the case of the functions like y-> -1*y the roles of the
229 // upper and lower error bars is reversed
230 if (f->Eval(x,y-eyl) < f->Eval(x,y+eyh)) {
231 eyl_new = TMath::Abs(fxy - f->Eval(x,y-eyl));
232 eyh_new = TMath::Abs(f->Eval(x,y+eyh) - fxy);
233 } else {
234 eyh_new = TMath::Abs(fxy - f->Eval(x,y-eyl));
235 eyl_new = TMath::Abs(f->Eval(x,y+eyh) - fxy);
236 }
237
238 //error on x doesn't change
239 SetPointError(i,exl,exh,eyl_new,eyh_new);
240 }
241 if (gPad) gPad->Modified();
242}
243
244
245////////////////////////////////////////////////////////////////////////////////
246/// Compute range.
247
249{
251
252 for (Int_t i=0;i<fNpoints;i++) {
253 if (fX[i] -fEXlow[i] < xmin) {
254 if (gPad && gPad->GetLogx()) {
255 if (fEXlow[i] < fX[i]) xmin = fX[i]-fEXlow[i];
256 else xmin = TMath::Min(xmin,fX[i]/3);
257 } else {
258 xmin = fX[i]-fEXlow[i];
259 }
260 }
261 if (fX[i] +fEXhigh[i] > xmax) xmax = fX[i]+fEXhigh[i];
262 if (fY[i] -fEYlow[i] < ymin) {
263 if (gPad && gPad->GetLogy()) {
264 if (fEYlow[i] < fY[i]) ymin = fY[i]-fEYlow[i];
265 else ymin = TMath::Min(ymin,fY[i]/3);
266 } else {
267 ymin = fY[i]-fEYlow[i];
268 }
269 }
270 if (fY[i] +fEYhigh[i] > ymax) ymax = fY[i]+fEYhigh[i];
271 }
272}
273
274
275////////////////////////////////////////////////////////////////////////////////
276/// Copy and release.
277
279 Int_t ibegin, Int_t iend, Int_t obegin)
280{
281 CopyPoints(newarrays, ibegin, iend, obegin);
282 if (newarrays) {
283 delete[] fEXlow;
284 fEXlow = newarrays[0];
285 delete[] fEXhigh;
286 fEXhigh = newarrays[1];
287 delete[] fEYlow;
288 fEYlow = newarrays[2];
289 delete[] fEYhigh;
290 fEYhigh = newarrays[3];
291 delete[] fEXlowd;
292 fEXlowd = newarrays[4];
293 delete[] fEXhighd;
294 fEXhighd = newarrays[5];
295 delete[] fEYlowd;
296 fEYlowd = newarrays[6];
297 delete[] fEYhighd;
298 fEYhighd = newarrays[7];
299 delete[] fX;
300 fX = newarrays[8];
301 delete[] fY;
302 fY = newarrays[9];
303 delete[] newarrays;
304 }
305}
306
307
308////////////////////////////////////////////////////////////////////////////////
309/// Copy errors from `fE*** `to `arrays[***]`
310/// or to `f***` Copy points.
311
313 Int_t ibegin, Int_t iend, Int_t obegin)
314{
315 if (TGraph::CopyPoints(arrays ? arrays+8 : nullptr, ibegin, iend, obegin)) {
316 Int_t n = (iend - ibegin)*sizeof(Double_t);
317 if (arrays) {
318 memmove(&arrays[0][obegin], &fEXlow[ibegin], n);
319 memmove(&arrays[1][obegin], &fEXhigh[ibegin], n);
320 memmove(&arrays[2][obegin], &fEYlow[ibegin], n);
321 memmove(&arrays[3][obegin], &fEYhigh[ibegin], n);
322 memmove(&arrays[4][obegin], &fEXlowd[ibegin], n);
323 memmove(&arrays[5][obegin], &fEXhighd[ibegin], n);
324 memmove(&arrays[6][obegin], &fEYlowd[ibegin], n);
325 memmove(&arrays[7][obegin], &fEYhighd[ibegin], n);
326 } else {
327 memmove(&fEXlow[obegin], &fEXlow[ibegin], n);
328 memmove(&fEXhigh[obegin], &fEXhigh[ibegin], n);
329 memmove(&fEYlow[obegin], &fEYlow[ibegin], n);
330 memmove(&fEYhigh[obegin], &fEYhigh[ibegin], n);
331 memmove(&fEXlowd[obegin], &fEXlowd[ibegin], n);
332 memmove(&fEXhighd[obegin], &fEXhighd[ibegin], n);
333 memmove(&fEYlowd[obegin], &fEYlowd[ibegin], n);
334 memmove(&fEYhighd[obegin], &fEYhighd[ibegin], n);
335 }
336 return kTRUE;
337 } else {
338 return kFALSE;
339 }
340}
341
342
343////////////////////////////////////////////////////////////////////////////////
344/// Should be called from ctors after `fNpoints` has been set.
345
347{
348 if (!fNpoints) {
349 fEXlow = fEYlow = fEXhigh = fEYhigh = nullptr;
350 fEXlowd = fEYlowd = fEXhighd = fEYhighd = nullptr;
351 return kFALSE;
352 }
353 fEXlow = new Double_t[fMaxSize];
354 fEYlow = new Double_t[fMaxSize];
361 return kTRUE;
362}
363
364////////////////////////////////////////////////////////////////////////////////
365/// Protected function to perform the merge operation of a graph with asymmetric errors.
366
368{
369 if (g->GetN() == 0) return kFALSE;
370
371 Double_t *exl = g->GetEXlow();
372 Double_t *exh = g->GetEXhigh();
373 Double_t *eyl = g->GetEYlow();
374 Double_t *eyh = g->GetEYhigh();
375
376 Double_t *exld = g->GetEXlowd();
377 Double_t *exhd = g->GetEXhighd();
378 Double_t *eyld = g->GetEYlowd();
379 Double_t *eyhd = g->GetEYhighd();
380
381 if (!exl || !exh || !eyl || !eyh ||
382 !exld || !exhd || !eyld || !eyhd) {
383 if (g->IsA() != TGraph::Class() )
384 Warning("DoMerge", "Merging a %s is not compatible with a TGraphBentErrors - errors will be ignored", g->IsA()->GetName());
385 return TGraph::DoMerge(g);
386 }
387 for (Int_t i = 0 ; i < g->GetN(); i++) {
388 Int_t ipoint = GetN();
389 Double_t x = g->GetX()[i];
390 Double_t y = g->GetY()[i];
391 SetPoint(ipoint, x, y);
392 SetPointError(ipoint, exl[i], exh[i], eyl[i], eyh[i],
393 exld[i], exhd[i], eyld[i], eyhd[i]);
394 }
395
396 return kTRUE;
397
398}
399////////////////////////////////////////////////////////////////////////////////
400/// It returns the error along X at point `i`.
401
403{
404 if (i < 0 || i >= fNpoints) return -1;
405 if (!fEXlow && !fEXhigh) return -1;
406 Double_t elow = 0, ehigh = 0;
407 if (fEXlow) elow = fEXlow[i];
408 if (fEXhigh) ehigh = fEXhigh[i];
409 return TMath::Sqrt(0.5*(elow*elow + ehigh*ehigh));
410}
411
412
413////////////////////////////////////////////////////////////////////////////////
414/// It returns the error along Y at point `i`.
415
417{
418 if (i < 0 || i >= fNpoints) return -1;
419 if (!fEYlow && !fEYhigh) return -1;
420 Double_t elow=0, ehigh=0;
421 if (fEYlow) elow = fEYlow[i];
422 if (fEYhigh) ehigh = fEYhigh[i];
423 return TMath::Sqrt(0.5*(elow*elow + ehigh*ehigh));
424}
425
426
427////////////////////////////////////////////////////////////////////////////////
428/// Get high error on X[i].
429
431{
432 if (i<0 || i>fNpoints) return -1;
433 if (fEXhigh) return fEXhigh[i];
434 return -1;
435}
436
437
438////////////////////////////////////////////////////////////////////////////////
439/// Get low error on X[i].
440
442{
443 if (i<0 || i>fNpoints) return -1;
444 if (fEXlow) return fEXlow[i];
445 return -1;
446}
447
448
449////////////////////////////////////////////////////////////////////////////////
450/// Get high error on Y[i].
451
453{
454 if (i<0 || i>fNpoints) return -1;
455 if (fEYhigh) return fEYhigh[i];
456 return -1;
457}
458
459
460////////////////////////////////////////////////////////////////////////////////
461/// Get low error on Y[i].
462
464{
465 if (i<0 || i>fNpoints) return -1;
466 if (fEYlow) return fEYlow[i];
467 return -1;
468}
469
470
471////////////////////////////////////////////////////////////////////////////////
472/// Set zero values for point arrays in the range `[begin, end]`
473
475 Bool_t from_ctor)
476{
477 if (!from_ctor) {
478 TGraph::FillZero(begin, end, from_ctor);
479 }
480 Int_t n = (end - begin)*sizeof(Double_t);
481 memset(fEXlow + begin, 0, n);
482 memset(fEXhigh + begin, 0, n);
483 memset(fEYlow + begin, 0, n);
484 memset(fEYhigh + begin, 0, n);
485 memset(fEXlowd + begin, 0, n);
486 memset(fEXhighd + begin, 0, n);
487 memset(fEYlowd + begin, 0, n);
488 memset(fEYhighd + begin, 0, n);
489}
490
491
492////////////////////////////////////////////////////////////////////////////////
493/// Print graph and errors values.
494
496{
497 for (Int_t i=0;i<fNpoints;i++) {
498 printf("x[%d]=%g, y[%d]=%g, exl[%d]=%g, exh[%d]=%g, eyl[%d]=%g, eyh[%d]=%g\n"
499 ,i,fX[i],i,fY[i],i,fEXlow[i],i,fEXhigh[i],i,fEYlow[i],i,fEYhigh[i]);
500 }
501}
502
503////////////////////////////////////////////////////////////////////////////////
504/// Multiply the values and errors of a TGraphBentErrors by a constant c1.
505///
506/// If option contains "x" the x values and errors are scaled
507/// If option contains "y" the y values and errors are scaled
508/// If option contains "xy" both x and y values and errors are scaled
509
511{
513 TString opt = option; opt.ToLower();
514 if (opt.Contains("x") && GetEXlow()) {
515 for (Int_t i=0; i<GetN(); i++)
516 GetEXlow()[i] *= c1;
517 }
518 if (opt.Contains("x") && GetEXhigh()) {
519 for (Int_t i=0; i<GetN(); i++)
520 GetEXhigh()[i] *= c1;
521 }
522 if (opt.Contains("y") && GetEYlow()) {
523 for (Int_t i=0; i<GetN(); i++)
524 GetEYlow()[i] *= c1;
525 }
526 if (opt.Contains("y") && GetEYhigh()) {
527 for (Int_t i=0; i<GetN(); i++)
528 GetEYhigh()[i] *= c1;
529 }
530 if (opt.Contains("x") && GetEXlowd()) {
531 for (Int_t i=0; i<GetN(); i++)
532 GetEXlowd()[i] *= c1;
533 }
534 if (opt.Contains("x") && GetEXhighd()) {
535 for (Int_t i=0; i<GetN(); i++)
536 GetEXhighd()[i] *= c1;
537 }
538 if (opt.Contains("y") && GetEYlowd()) {
539 for (Int_t i=0; i<GetN(); i++)
540 GetEYlowd()[i] *= c1;
541 }
542 if (opt.Contains("y") && GetEYhighd()) {
543 for (Int_t i=0; i<GetN(); i++)
544 GetEYhighd()[i] *= c1;
545 }
546}
547
548////////////////////////////////////////////////////////////////////////////////
549/// Save primitive as a C++ statement(s) on output stream out.
550
551void TGraphBentErrors::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
552{
553 out << " " << std::endl;
554 static Int_t frameNumber = 2000;
555 frameNumber++;
556
557 auto fXName = SaveArray(out, "fx", frameNumber, fX);
558 auto fYName = SaveArray(out, "fy", frameNumber, fY);
559 auto fElXName = SaveArray(out, "felx", frameNumber, fEXlow);
560 auto fElYName = SaveArray(out, "fely", frameNumber, fEYlow);
561 auto fEhXName = SaveArray(out, "fehx", frameNumber, fEXhigh);
562 auto fEhYName = SaveArray(out, "fehy", frameNumber, fEYhigh);
563 auto fEldXName = SaveArray(out, "feldx", frameNumber, fEXlowd);
564 auto fEldYName = SaveArray(out, "feldy", frameNumber, fEYlowd);
565 auto fEhdXName = SaveArray(out, "fehdx", frameNumber, fEXhighd);
566 auto fEhdYName = SaveArray(out, "fehdy", frameNumber, fEYhighd);
567
568 if (gROOT->ClassSaved(TGraphBentErrors::Class()))
569 out << " ";
570 else
571 out << " TGraphBentErrors *";
572 out << "grbe = new TGraphBentErrors("<< fNpoints << ","
573 << fXName << "," << fYName << ","
574 << fElXName << "," << fEhXName << ","
575 << fElYName << "," << fEhYName << ","
576 << fEldXName << "," << fEhdXName << ","
577 << fEldYName << "," << fEhdYName << ");"
578 << std::endl;
579
580 SaveHistogramAndFunctions(out, "grbe", frameNumber, option);
581}
582
583
584////////////////////////////////////////////////////////////////////////////////
585/// Set ex and ey values for point pointed by the mouse.
586
588 Double_t exld, Double_t exhd, Double_t eyld, Double_t eyhd)
589{
590 if (!gPad) {
591 Error("SetPointError", "Cannot be used without gPad, requires last mouse position");
592 return;
593 }
594
595 Int_t px = gPad->GetEventX();
596 Int_t py = gPad->GetEventY();
597
598 //localize point to be deleted
599 Int_t ipoint = -2;
600 // start with a small window (in case the mouse is very close to one point)
601 for (Int_t i = 0; i < fNpoints; i++) {
602 Int_t dpx = px - gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
603 Int_t dpy = py - gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
604 if (dpx*dpx+dpy*dpy < 25) {ipoint = i; break;}
605 }
606 if (ipoint == -2) return;
607
608 fEXlow[ipoint] = exl;
609 fEYlow[ipoint] = eyl;
610 fEXhigh[ipoint] = exh;
611 fEYhigh[ipoint] = eyh;
612 fEXlowd[ipoint] = exld;
613 fEXhighd[ipoint] = exhd;
614 fEYlowd[ipoint] = eyld;
615 fEYhighd[ipoint] = eyhd;
616
617 gPad->Modified();
618}
619
620
621////////////////////////////////////////////////////////////////////////////////
622/// Set ex and ey values for point number `i`.
623
625 Double_t exld, Double_t exhd, Double_t eyld, Double_t eyhd)
626{
627 if (i < 0) return;
628 if (i >= fNpoints) {
629 // re-allocate the object
631 }
632 fEXlow[i] = exl;
633 fEYlow[i] = eyl;
634 fEXhigh[i] = exh;
635 fEYhigh[i] = eyh;
636 fEXlowd[i] = exld;
637 fEXhighd[i] = exhd;
638 fEYlowd[i] = eyld;
639 fEYhighd[i] = eyhd;
640}
641
642
643////////////////////////////////////////////////////////////////////////////////
644/// Swap points.
645
647{
648 SwapValues(fEXlow, pos1, pos2);
649 SwapValues(fEXhigh, pos1, pos2);
650 SwapValues(fEYlow, pos1, pos2);
651 SwapValues(fEYhigh, pos1, pos2);
652
653 SwapValues(fEXlowd, pos1, pos2);
654 SwapValues(fEXhighd, pos1, pos2);
655 SwapValues(fEYlowd, pos1, pos2);
656 SwapValues(fEYhighd, pos1, pos2);
657
658 TGraph::SwapPoints(pos1, pos2);
659}
660
661////////////////////////////////////////////////////////////////////////////////
662/// Update the fX, fY, fEXlow, fEXhigh, fEXlowd, fEXhighd, fEYlow, fEYhigh, fEYlowd,
663/// and fEYhighd arrays with the sorted values.
664
665void TGraphBentErrors::UpdateArrays(const std::vector<Int_t> &sorting_indices, Int_t numSortedPoints, Int_t low)
666{
667 std::vector<Double_t> fEXlowSorted(numSortedPoints);
668 std::vector<Double_t> fEXhighSorted(numSortedPoints);
669 std::vector<Double_t> fEXlowdSorted(numSortedPoints);
670 std::vector<Double_t> fEXhighdSorted(numSortedPoints);
671
672 std::vector<Double_t> fEYlowSorted(numSortedPoints);
673 std::vector<Double_t> fEYhighSorted(numSortedPoints);
674 std::vector<Double_t> fEYlowdSorted(numSortedPoints);
675 std::vector<Double_t> fEYhighdSorted(numSortedPoints);
676
677 // Fill the sorted X and Y error values based on the sorted indices
678 std::generate(fEXlowSorted.begin(), fEXlowSorted.end(),
679 [begin = low, &sorting_indices, this]() mutable { return fEXlow[sorting_indices[begin++]]; });
680 std::generate(fEXhighSorted.begin(), fEXhighSorted.end(),
681 [begin = low, &sorting_indices, this]() mutable { return fEXhigh[sorting_indices[begin++]]; });
682 std::generate(fEXlowdSorted.begin(), fEXlowdSorted.end(),
683 [begin = low, &sorting_indices, this]() mutable { return fEXlowd[sorting_indices[begin++]]; });
684 std::generate(fEXhighdSorted.begin(), fEXhighdSorted.end(),
685 [begin = low, &sorting_indices, this]() mutable { return fEXhighd[sorting_indices[begin++]]; });
686
687 std::generate(fEYlowSorted.begin(), fEYlowSorted.end(),
688 [begin = low, &sorting_indices, this]() mutable { return fEYlow[sorting_indices[begin++]]; });
689 std::generate(fEYhighSorted.begin(), fEYhighSorted.end(),
690 [begin = low, &sorting_indices, this]() mutable { return fEYhigh[sorting_indices[begin++]]; });
691 std::generate(fEYlowdSorted.begin(), fEYlowdSorted.end(),
692 [begin = low, &sorting_indices, this]() mutable { return fEYlowd[sorting_indices[begin++]]; });
693 std::generate(fEYhighdSorted.begin(), fEYhighdSorted.end(),
694 [begin = low, &sorting_indices, this]() mutable { return fEYhighd[sorting_indices[begin++]]; });
695
696 // Copy the sorted X and Y error values back to the original arrays
697 std::copy(fEXlowSorted.begin(), fEXlowSorted.end(), fEXlow + low);
698 std::copy(fEXhighSorted.begin(), fEXhighSorted.end(), fEXhigh + low);
699 std::copy(fEXlowdSorted.begin(), fEXlowdSorted.end(), fEXlowd + low);
700 std::copy(fEXhighdSorted.begin(), fEXhighdSorted.end(), fEXhighd + low);
701
702 std::copy(fEYlowSorted.begin(), fEYlowSorted.end(), fEYlow + low);
703 std::copy(fEYhighSorted.begin(), fEYhighSorted.end(), fEYhigh + low);
704 std::copy(fEYlowdSorted.begin(), fEYlowdSorted.end(), fEYlowd + low);
705 std::copy(fEYhighdSorted.begin(), fEYhighdSorted.end(), fEYhighd + low);
706
707 TGraph::UpdateArrays(sorting_indices, numSortedPoints, low);
708}
#define f(i)
Definition RSha256.hxx:104
#define g(i)
Definition RSha256.hxx:105
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:382
Option_t Option_t option
float xmin
float ymin
float xmax
float ymax
#define gROOT
Definition TROOT.h:406
#define gPad
1-Dim function class
Definition TF1.h:233
A TGraphBentErrors is a TGraph with bent, asymmetric error bars.
void UpdateArrays(const std::vector< Int_t > &sorting_indices, Int_t numSortedPoints, Int_t low) override
Update the fX, fY, fEXlow, fEXhigh, fEXlowd, fEXhighd, fEYlow, fEYhigh, fEYlowd, and fEYhighd array...
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitive as a C++ statement(s) on output stream out.
TGraphBentErrors()
TGraphBentErrors default constructor.
Double_t * GetEXlow() const override
Double_t * fEYhigh
[fNpoints] array of Y high errors
Double_t * fEXlowd
[fNpoints] array of X low displacements
Double_t * fEYlowd
[fNpoints] array of Y low displacements
Bool_t CtorAllocate()
Should be called from ctors after fNpoints has been set.
Bool_t DoMerge(const TGraph *g) override
Protected function to perform the merge operation of a graph with asymmetric errors.
Double_t * GetEYhigh() const override
Double_t GetErrorXlow(Int_t bin) const override
Get low error on X[i].
Double_t * GetEYlow() const override
Double_t * GetEYhighd() const override
void ComputeRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const override
Compute range.
Double_t * GetEXhigh() const override
static TClass * Class()
Double_t * GetEYlowd() const override
Double_t * fEXhighd
[fNpoints] array of X high displacements
Double_t * GetEXhighd() const override
Double_t GetErrorXhigh(Int_t bin) const override
Get high error on X[i].
virtual void AddPointError(Double_t x, Double_t y, Double_t exl, Double_t exh, Double_t eyl, Double_t eyh, Double_t exld=0, Double_t exhd=0, Double_t eyld=0, Double_t eyhd=0)
Add a point with bent errors to the graph.
virtual void SetPointError(Double_t exl, Double_t exh, Double_t eyl, Double_t eyh, Double_t exld=0, Double_t exhd=0, Double_t eyld=0, Double_t eyhd=0)
Set ex and ey values for point pointed by the mouse.
Double_t * fEXlow
[fNpoints] array of X low errors
Double_t * fEXhigh
[fNpoints] array of X high errors
~TGraphBentErrors() override
TGraphBentErrors default destructor.
void Print(Option_t *chopt="") const override
Print graph and errors values.
void Scale(Double_t c1=1., Option_t *option="y") override
Multiply the values and errors of a TGraphBentErrors by a constant c1.
void FillZero(Int_t begin, Int_t end, Bool_t from_ctor=kTRUE) override
Set zero values for point arrays in the range [begin, end]
Double_t * GetEXlowd() const override
Double_t GetErrorY(Int_t bin) const override
It returns the error along Y at point i.
Double_t * fEYlow
[fNpoints] array of Y low errors
Double_t * fEYhighd
[fNpoints] array of Y high displacements
void CopyAndRelease(Double_t **newarrays, Int_t ibegin, Int_t iend, Int_t obegin) override
Copy and release.
void SwapPoints(Int_t pos1, Int_t pos2) override
Swap points.
Bool_t CopyPoints(Double_t **arrays, Int_t ibegin, Int_t iend, Int_t obegin) override
Copy errors from fE***to arrays[***] or to f*** Copy points.
Double_t GetErrorYlow(Int_t bin) const override
Get low error on Y[i].
Double_t GetErrorX(Int_t bin) const override
It returns the error along X at point i.
Double_t GetErrorYhigh(Int_t bin) const override
Get high error on Y[i].
void Apply(TF1 *f) override
Apply a function to all data points .
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
static TClass * Class()
virtual void AddPoint(Double_t x, Double_t y)
Append a new point to the graph.
Definition TGraph.h:98
Int_t fNpoints
Number of points <= fMaxSize.
Definition TGraph.h:46
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition TGraph.cxx:2347
Int_t fMaxSize
!Current dimension of arrays fX and fY
Definition TGraph.h:45
void SaveHistogramAndFunctions(std::ostream &out, const char *varname, Int_t &frameNumber, Option_t *option)
Save histogram and list of functions of TGraph as C++ statement Used in all TGraph-derived classes.
Definition TGraph.cxx:2204
TH1F * fHistogram
Pointer to histogram used for drawing axis.
Definition TGraph.h:50
virtual void UpdateArrays(const std::vector< Int_t > &sorting_indices, Int_t numSortedPoints, Int_t low)
Update the fX and fY arrays with the sorted values.
Definition TGraph.cxx:2597
Int_t GetN() const
Definition TGraph.h:132
Double_t * fY
[fNpoints] array of Y points
Definition TGraph.h:48
TString SaveArray(std::ostream &out, const char *suffix, Int_t frameNumber, Double_t *arr)
Save array as C++ code Returns name of created array.
Definition TGraph.cxx:2177
virtual void ComputeRange(Double_t &xmin, Double_t &ymin, Double_t &xmax, Double_t &ymax) const
Compute the x/y range of the points in this graph.
Definition TGraph.cxx:733
virtual void Scale(Double_t c1=1., Option_t *option="y")
Multiply the values of a TGraph by a constant c1.
Definition TGraph.cxx:2264
static void SwapValues(Double_t *arr, Int_t pos1, Int_t pos2)
Swap values.
Definition TGraph.cxx:2616
virtual Bool_t DoMerge(const TGraph *g)
protected function to perform the merge operation of a graph
Definition TGraph.cxx:2681
virtual void SwapPoints(Int_t pos1, Int_t pos2)
Swap points.
Definition TGraph.cxx:2588
virtual void FillZero(Int_t begin, Int_t end, Bool_t from_ctor=kTRUE)
Set zero values for point arrays in the range [begin, end) Should be redefined in descendant classes.
Definition TGraph.cxx:1104
Double_t * fX
[fNpoints] array of X points
Definition TGraph.h:47
virtual Int_t GetPoint(Int_t i, Double_t &x, Double_t &y) const
Get x and y values for point number i.
Definition TGraph.cxx:1535
virtual Bool_t CopyPoints(Double_t **newarrays, Int_t ibegin, Int_t iend, Int_t obegin)
Copy points from fX and fY to arrays[0] and arrays[1] or to fX and fY if arrays == 0 and ibegin !...
Definition TGraph.cxx:781
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:991
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1005
Basic string class.
Definition TString.h:139
void ToLower()
Change string to lower-case.
Definition TString.cxx:1182
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:632
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
TGraphErrors * gr
Definition legend1.C:25
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:666
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123