#include "TROOT.h"
#include "TMultiGraph.h"
#include "TGraph.h"
#include "TH1.h"
#include "TVirtualPad.h"
#include "Riostream.h"
#include "TVirtualFitter.h"
#include "TPluginManager.h"
#include "TClass.h"
#include "TMath.h"
#include "TSystem.h"
#include <stdlib.h>
#include <ctype.h>
extern void H1LeastSquareSeqnd(Int_t n, Double_t *a, Int_t idim, Int_t &ifail, Int_t k, Double_t *b);
ClassImp(TMultiGraph)
TMultiGraph::TMultiGraph(): TNamed()
{
fGraphs = 0;
fFunctions = 0;
fHistogram = 0;
fMaximum = -1111;
fMinimum = -1111;
}
TMultiGraph::TMultiGraph(const char *name, const char *title)
: TNamed(name,title)
{
fGraphs = 0;
fFunctions = 0;
fHistogram = 0;
fMaximum = -1111;
fMinimum = -1111;
}
TMultiGraph::TMultiGraph(const TMultiGraph& mg) :
TNamed (mg),
fGraphs(mg.fGraphs),
fFunctions(mg.fFunctions),
fHistogram(mg.fHistogram),
fMaximum(mg.fMaximum),
fMinimum(mg.fMinimum)
{
}
TMultiGraph& TMultiGraph::operator=(const TMultiGraph& mg)
{
if(this!=&mg) {
TNamed::operator=(mg);
fGraphs=mg.fGraphs;
fFunctions=mg.fFunctions;
fHistogram=mg.fHistogram;
fMaximum=mg.fMaximum;
fMinimum=mg.fMinimum;
}
return *this;
}
TMultiGraph::~TMultiGraph()
{
if (!fGraphs) return;
TGraph *g;
TIter next(fGraphs);
while ((g = (TGraph*) next())) {
g->ResetBit(kMustCleanup);
}
fGraphs->Delete();
delete fGraphs;
fGraphs = 0;
delete fHistogram;
fHistogram = 0;
if (fFunctions) {
fFunctions->SetBit(kInvalidObject);
TObject *obj;
while ((obj = fFunctions->First())) {
while(fFunctions->Remove(obj)) { }
delete obj;
}
delete fFunctions;
}
}
void TMultiGraph::Add(TGraph *graph, Option_t *chopt)
{
if (!fGraphs) fGraphs = new TList();
graph->SetBit(kMustCleanup);
fGraphs->Add(graph,chopt);
}
void TMultiGraph::Add(TMultiGraph *multigraph, Option_t *chopt)
{
TList *graphlist = multigraph->GetListOfGraphs();
if (!graphlist) return;
if (!fGraphs) fGraphs = new TList();
TGraph *gr;
gr = (TGraph*)graphlist->First();
fGraphs->Add(gr,chopt);
for(Int_t i = 1; i < graphlist->GetSize(); i++){
gr = (TGraph*)graphlist->After(gr);
fGraphs->Add(gr,chopt);
}
}
void TMultiGraph::Browse(TBrowser *)
{
Draw("alp");
gPad->Update();
}
Int_t TMultiGraph::DistancetoPrimitive(Int_t px, Int_t py)
{
const Int_t kMaxDiff = 10;
Int_t distance = 9999;
if (fHistogram) {
distance = fHistogram->DistancetoPrimitive(px,py);
if (distance <= 0) return distance;
}
if (!fGraphs) return distance;
TGraph *g;
TIter next(fGraphs);
while ((g = (TGraph*) next())) {
Int_t dist = g->DistancetoPrimitive(px,py);
if (dist <= 0) return 0;
if (dist < kMaxDiff) {gPad->SetSelected(g); return dist;}
}
return distance;
}
void TMultiGraph::Draw(Option_t *option)
{
if (!fGraphs) {
Error("Draw", "Cannot draw an empty TMultiGraph");
return;
}
AppendPad(option);
}
Int_t TMultiGraph::Fit(const char *fname, Option_t *option, Option_t *, Axis_t xmin, Axis_t xmax)
{
char *linear;
linear= (char*)strstr(fname, "++");
TF1 *f1=0;
if (linear)
f1=new TF1(fname, fname, xmin, xmax);
else {
f1 = (TF1*)gROOT->GetFunction(fname);
if (!f1) { Printf("Unknown function: %s",fname); return -1; }
}
return Fit(f1,option,"",xmin,xmax);
}
Int_t TMultiGraph::Fit(TF1 *f1, Option_t *option, Option_t *goption, Axis_t rxmin, Axis_t rxmax)
{
return DoFit(f1,option,goption,rxmin,rxmax);
}
void TMultiGraph::FitPanel()
{
if (!gPad)
gROOT->MakeDefCanvas();
if (!gPad) {
Error("FitPanel", "Unable to create a default canvas");
return;
}
TPluginHandler *handler = gROOT->GetPluginManager()->FindHandler("TFitEditor");
if (handler && handler->LoadPlugin() != -1) {
if (handler->ExecPlugin(2, gPad, this) == 0)
Error("FitPanel", "Unable to crate the FitPanel");
}
else
Error("FitPanel", "Unable to find the FitPanel plug-in");
}
Option_t *TMultiGraph::GetGraphDrawOption(const TGraph *gr) const
{
if (!fGraphs || !gr) return "";
TListIter next(fGraphs);
TObject *obj;
while ((obj = next())) {
if (obj == (TObject*)gr) return next.GetOption();
}
return "";
}
void TMultiGraph::InitGaus(Double_t xmin, Double_t xmax)
{
Double_t allcha, sumx, sumx2, x, val, rms, mean;
Int_t bin;
const Double_t sqrtpi = 2.506628;
Int_t np = 0;
allcha = sumx = sumx2 = 0;
TGraph *g;
TIter next(fGraphs);
Double_t *px, *py;
Int_t npp;
while ((g = (TGraph*) next())) {
px=g->GetX();
py=g->GetY();
npp=g->GetN();
for (bin=0; bin<npp; bin++){
x=px[bin];
if (x<xmin || x>xmax) continue;
np++;
val=py[bin];
sumx+=val*x;
sumx2+=val*x*x;
allcha+=val;
}
}
if (np == 0 || allcha == 0) return;
mean = sumx/allcha;
rms = TMath::Sqrt(sumx2/allcha - mean*mean);
Double_t binwidx = TMath::Abs((xmax-xmin)/np);
if (rms == 0) rms = 1;
TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
TF1 *f1 = (TF1*)grFitter->GetUserFunc();
f1->SetParameter(0,binwidx*allcha/(sqrtpi*rms));
f1->SetParameter(1,mean);
f1->SetParameter(2,rms);
f1->SetParLimits(2,0,10*rms);
}
void TMultiGraph::InitExpo(Double_t xmin, Double_t xmax)
{
Double_t constant, slope;
Int_t ifail;
LeastSquareLinearFit(-1, constant, slope, ifail, xmin, xmax);
TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
TF1 *f1 = (TF1*)grFitter->GetUserFunc();
f1->SetParameter(0,constant);
f1->SetParameter(1,slope);
}
void TMultiGraph::InitPolynom(Double_t xmin, Double_t xmax)
{
Double_t fitpar[25];
TVirtualFitter *grFitter = TVirtualFitter::GetFitter();
TF1 *f1 = (TF1*)grFitter->GetUserFunc();
Int_t npar = f1->GetNpar();
LeastSquareFit(npar, fitpar, xmin, xmax);
for (Int_t i=0;i<npar;i++) f1->SetParameter(i, fitpar[i]);
}
void TMultiGraph::LeastSquareFit(Int_t m, Double_t *a, Double_t xmin, Double_t xmax)
{
const Double_t zero = 0.;
const Double_t one = 1.;
const Int_t idim = 20;
Double_t b[400] ;
Int_t i, k, l, ifail, bin;
Double_t power;
Double_t da[20], xk, yk;
TGraph *g;
TIter next(fGraphs);
Double_t *px, *py;
Int_t n=0;
Int_t npp;
while ((g = (TGraph*) next())) {
px=g->GetX();
py=g->GetY();
npp=g->GetN();
for (bin=0; bin<npp; bin++){
xk=px[bin];
if (xk < xmin || xk > xmax) continue;
n++;
}
}
if (m <= 2) {
LeastSquareLinearFit(n, a[0], a[1], ifail, xmin, xmax);
return;
}
if (m > idim || m > n) return;
da[0] = zero;
for (l = 2; l <= m; ++l) {
b[l-1] = zero;
b[m + l*20 - 21] = zero;
da[l-1] = zero;
}
Int_t np = 0;
next.Reset();
while ((g = (TGraph*) next())) {
px=g->GetX();
py=g->GetY();
npp=g->GetN();
for (k = 0; k <= npp; ++k) {
xk = px[k];
if (xk < xmin || xk > xmax) continue;
np++;
yk = py[k];
power = one;
da[0] += yk;
for (l = 2; l <= m; ++l) {
power *= xk;
b[l-1] += power;
da[l-1] += power*yk;
}
for (l = 2; l <= m; ++l) {
power *= xk;
b[m + l*20 - 21] += power;
}
}
}
b[0] = Double_t(np);
for (i = 3; i <= m; ++i) {
for (k = i; k <= m; ++k) {
b[k - 1 + (i-1)*20 - 21] = b[k + (i-2)*20 - 21];
}
}
H1LeastSquareSeqnd(m, b, idim, ifail, 1, da);
if (ifail < 0) {
py=((TGraph *)fGraphs->First())->GetY();
a[0]=py[0];
for (i=1; i<m; ++i) a[i] = 0;
return;
}
for (i=0; i<m; ++i) a[i] = da[i];
}
void TMultiGraph::LeastSquareLinearFit(Int_t ndata, Double_t &a0, Double_t &a1, Int_t &ifail, Double_t xmin, Double_t xmax)
{
Double_t xbar, ybar, x2bar;
Int_t i;
Double_t xybar;
Double_t fn, xk, yk;
Double_t det;
ifail = -2;
xbar = ybar = x2bar = xybar = 0;
Int_t np = 0;
TGraph *g;
TIter next(fGraphs);
Double_t *px, *py;
Int_t npp;
while ((g = (TGraph*) next())) {
px=g->GetX();
py=g->GetY();
npp=g->GetN();
for (i = 0; i < npp; ++i) {
xk = px[i];
if (xk < xmin || xk > xmax) continue;
np++;
yk = py[i];
if (ndata < 0) {
if (yk <= 0) yk = 1e-9;
yk = TMath::Log(yk);
}
xbar += xk;
ybar += yk;
x2bar += xk*xk;
xybar += xk*yk;
}
}
fn = Double_t(np);
det = fn*x2bar - xbar*xbar;
ifail = -1;
if (det <= 0) {
if (fn > 0) a0 = ybar/fn;
else a0 = 0;
a1 = 0;
return;
}
ifail = 0;
a0 = (x2bar*ybar - xbar*xybar) / det;
a1 = (fn*xybar - xbar*ybar) / det;
}
TH1F *TMultiGraph::GetHistogram() const
{
if (fHistogram) return fHistogram;
if (!gPad) return 0;
gPad->Modified();
gPad->Update();
if (fHistogram) return fHistogram;
TH1F *h1 = (TH1F*)gPad->FindObject("hframe");
return h1;
}
TF1 *TMultiGraph::GetFunction(const char *name) const
{
if (!fFunctions) return 0;
return (TF1*)fFunctions->FindObject(name);
}
TList *TMultiGraph::GetListOfFunctions()
{
if (!fFunctions) fFunctions = new TList();
return fFunctions;
}
TAxis *TMultiGraph::GetXaxis() const
{
if (!gPad) return 0;
TH1 *h = GetHistogram();
if (!h) return 0;
return h->GetXaxis();
}
TAxis *TMultiGraph::GetYaxis() const
{
if (!gPad) return 0;
TH1 *h = GetHistogram();
if (!h) return 0;
return h->GetYaxis();
}
void TMultiGraph::Paint(Option_t *option)
{
if (fGraphs->GetSize() == 0) return;
char *l;
static char chopt[33];
Int_t nch = strlen(option);
Int_t i;
for (i=0;i<nch;i++) chopt[i] = toupper(option[i]);
chopt[nch] = 0;
TGraph *g;
l = strstr(chopt,"A");
if (l) {
*l = ' ';
TIter next(fGraphs);
Int_t npt = 100;
Double_t maximum, minimum, rwxmin, rwxmax, rwymin, rwymax, uxmin, uxmax, dx, dy;
rwxmin = gPad->GetUxmin();
rwxmax = gPad->GetUxmax();
rwymin = gPad->GetUymin();
rwymax = gPad->GetUymax();
char *xtitle = 0;
char *ytitle = 0;
Int_t firstx = 0;
Int_t lastx = 0;
if (fHistogram) {
if (fHistogram->GetMinimum() >= fHistogram->GetMaximum()) {
nch = strlen(fHistogram->GetXaxis()->GetTitle());
firstx = fHistogram->GetXaxis()->GetFirst();
lastx = fHistogram->GetXaxis()->GetLast();
if (nch) {
xtitle = new char[nch+1];
strcpy(xtitle,fHistogram->GetXaxis()->GetTitle());
}
nch = strlen(fHistogram->GetYaxis()->GetTitle());
if (nch) {
ytitle = new char[nch+1];
strcpy(ytitle,fHistogram->GetYaxis()->GetTitle());
}
delete fHistogram;
fHistogram = 0;
}
}
if (fHistogram) {
minimum = fHistogram->GetYaxis()->GetXmin();
maximum = fHistogram->GetYaxis()->GetXmax();
uxmin = gPad->PadtoX(rwxmin);
uxmax = gPad->PadtoX(rwxmax);
} else {
g = (TGraph*) next();
g->ComputeRange(rwxmin, rwymin, rwxmax, rwymax);
while ((g = (TGraph*) next())) {
Double_t rx1,ry1,rx2,ry2;
g->ComputeRange(rx1, ry1, rx2, ry2);
if (rx1 < rwxmin) rwxmin = rx1;
if (ry1 < rwymin) rwymin = ry1;
if (rx2 > rwxmax) rwxmax = rx2;
if (ry2 > rwymax) rwymax = ry2;
if (g->GetN() > npt) npt = g->GetN();
}
if (rwxmin == rwxmax) rwxmax += 1.;
if (rwymin == rwymax) rwymax += 1.;
dx = 0.05*(rwxmax-rwxmin);
dy = 0.05*(rwymax-rwymin);
uxmin = rwxmin - dx;
uxmax = rwxmax + dx;
if (gPad->GetLogy()) {
if (rwymin <= 0) rwymin = 0.001*rwymax;
minimum = rwymin/(1+0.5*TMath::Log10(rwymax/rwymin));
maximum = rwymax*(1+0.2*TMath::Log10(rwymax/rwymin));
} else {
minimum = rwymin - dy;
maximum = rwymax + dy;
}
if (minimum < 0 && rwymin >= 0) minimum = 0;
if (maximum > 0 && rwymax <= 0) maximum = 0;
}
if (fMinimum != -1111) rwymin = minimum = fMinimum;
if (fMaximum != -1111) rwymax = maximum = fMaximum;
if (uxmin < 0 && rwxmin >= 0) {
if (gPad->GetLogx()) uxmin = 0.9*rwxmin;
}
if (uxmax > 0 && rwxmax <= 0) {
if (gPad->GetLogx()) uxmax = 1.1*rwxmax;
}
if (minimum < 0 && rwymin >= 0) {
if(gPad->GetLogy()) minimum = 0.9*rwymin;
}
if (maximum > 0 && rwymax <= 0) {
if(gPad->GetLogy()) maximum = 1.1*rwymax;
}
if (minimum <= 0 && gPad->GetLogy()) minimum = 0.001*maximum;
if (uxmin <= 0 && gPad->GetLogx()) {
if (uxmax > 1000) uxmin = 1;
else uxmin = 0.001*uxmax;
}
rwymin = minimum;
rwymax = maximum;
if (fHistogram) {
fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
}
if (!fHistogram) {
rwxmin = uxmin;
rwxmax = uxmax;
fHistogram = new TH1F(GetName(),GetTitle(),npt,rwxmin,rwxmax);
if (!fHistogram) return;
fHistogram->SetMinimum(rwymin);
fHistogram->SetBit(TH1::kNoStats);
fHistogram->SetMaximum(rwymax);
fHistogram->GetYaxis()->SetLimits(rwymin,rwymax);
fHistogram->SetDirectory(0);
if (xtitle) {fHistogram->GetXaxis()->SetTitle(xtitle); delete [] xtitle;}
if (ytitle) {fHistogram->GetYaxis()->SetTitle(ytitle); delete [] ytitle;}
if (firstx != lastx) fHistogram->GetXaxis()->SetRange(firstx,lastx);
}
fHistogram->Paint("0");
}
TGraph *gfit = 0;
if (fGraphs) {
TObjOptLink *lnk = (TObjOptLink*)fGraphs->FirstLink();
TObject *obj = 0;
while (lnk) {
obj = lnk->GetObject();
if (strlen(lnk->GetOption())) obj->Paint(lnk->GetOption());
else obj->Paint(chopt);
lnk = (TObjOptLink*)lnk->Next();
}
gfit = (TGraph*)obj;
}
TObject *f;
TF1 *fit = 0;
if (fFunctions) {
TIter next(fFunctions);
while ((f = (TObject*) next())) {
if (f->InheritsFrom(TF1::Class())) {
if (f->TestBit(TF1::kNotDraw) == 0) f->Paint("lsame");
fit = (TF1*)f;
} else {
f->Paint();
}
}
}
if (fit) gfit->PaintStats(fit);
}
void TMultiGraph::Print(Option_t *option) const
{
TGraph *g;
if (fGraphs) {
TIter next(fGraphs);
while ((g = (TGraph*) next())) {
g->Print(option);
}
}
}
void TMultiGraph::RecursiveRemove(TObject *obj)
{
if (!fGraphs) return;
TObject *objr = fGraphs->Remove(obj);
if (!objr) return;
delete fHistogram; fHistogram = 0;
if (gPad) gPad->Modified();
}
void TMultiGraph::SavePrimitive(ostream &out, Option_t *option )
{
char quote = '"';
out<<" "<<endl;
if (gROOT->ClassSaved(TMultiGraph::Class())) {
out<<" ";
} else {
out<<" TMultiGraph *";
}
out<<"multigraph = new TMultiGraph();"<<endl;
out<<" multigraph->SetName("<<quote<<GetName()<<quote<<");"<<endl;
out<<" multigraph->SetTitle("<<quote<<GetTitle()<<quote<<");"<<endl;
if (fGraphs) {
TObjOptLink *lnk = (TObjOptLink*)fGraphs->FirstLink();
TObject *g;
while (lnk) {
g = lnk->GetObject();
g->SavePrimitive(out, Form("multigraph%s",lnk->GetOption()));
lnk = (TObjOptLink*)lnk->Next();
}
}
out<<" multigraph->Draw(" <<quote<<option<<quote<<");"<<endl;
TAxis *xaxis = GetXaxis();
TAxis *yaxis = GetYaxis();
if (xaxis) xaxis->SaveAttributes(out, "multigraph","->GetXaxis()");
if (yaxis) yaxis->SaveAttributes(out, "multigraph","->GetYaxis()");
}
void TMultiGraph::SetMaximum(Double_t maximum)
{
fMaximum = maximum;
if (fHistogram) fHistogram->SetMaximum(maximum);
}
void TMultiGraph::SetMinimum(Double_t minimum)
{
fMinimum = minimum;
if (fHistogram) fHistogram->SetMinimum(minimum);
}
Last change: Fri Oct 31 14:49:42 2008
Last generated: 2008-10-31 14:49
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.