/* ======================================================================== *\ ! ! * ! * This file is part of MARS, the MAGIC Analysis and Reconstruction ! * Software. It is distributed to you in the hope that it can be a useful ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes. ! * It is distributed WITHOUT ANY WARRANTY. ! * ! * Permission to use, copy, modify and distribute this software and its ! * documentation for any purpose is hereby granted without fee, ! * provided that the above copyright notice appear in all copies and ! * that both that copyright notice and this permission notice appear ! * in supporting documentation. It is provided "as is" without express ! * or implied warranty. ! * ! ! ! Author(s): Thomas Bretz, 3/2004 ! ! Copyright: MAGIC Software Development, 2000-2004 ! ! \* ======================================================================== */ ////////////////////////////////////////////////////////////////////////////// // // MHFalseSource // // Create a false source plot. For the Binning in x,y the object MBinning // "BinningFalseSource" is taken from the paremeter list. // // The binning in alpha is currently fixed as 18bins from 0 to 90deg. // // If MTime, MObservatory and MPointingPos is available in the paremeter // list each image is derotated. // // MHFalseSource fills a 3D histogram with alpha distribution for // each (derotated) x and y position. // // Only a radius of 1.2deg around the camera center is taken into account. // // The displayed histogram is the projection of alpha90-fAlphaCut // // By using the context menu (find it in a small region between the single // pads) you can change the AlphaCut 'online' // // Each Z-Projection (Alpha-histogram) is scaled such, that the number // of entries fits the maximum number of entries in all Z-Projections. // This should correct for the different acceptance in the center // and at the vorder of the camera. This, however, produces more noise // at the border. // // Here is a slightly simplified version of the algorithm: // ------------------------------------------------------ // MHillas hil; // Taken as second argument in MFillH // // MSrcPosCam src; // MHillasSrc hsrc; // hsrc.SetSrcPos(&src); // // for (int ix=0; ix #include #include #include #include #include #include #include "MGeomCam.h" #include "MSrcPosCam.h" #include "MHillasSrc.h" #include "MTime.h" #include "MObservatory.h" #include "MPointingPos.h" #include "MAstroCatalog.h" #include "MAstroSky2Local.h" #include "MStatusDisplay.h" #include "MMath.h" #include "MBinning.h" #include "MParList.h" #include "MLog.h" #include "MLogManip.h" ClassImp(MHFalseSource); using namespace std; // -------------------------------------------------------------------------- // // Default Constructor // MHFalseSource::MHFalseSource(const char *name, const char *title) : fTime(0), fPointPos(0), fObservatory(0), fMm2Deg(-1), fAlphaCut(12.5), fBgMean(55), fMinDist(-1), fMaxDist(-1), fMinLD(-1), fMaxLD(-1) { // // set the name and title of this object // fName = name ? name : "MHFalseSource"; fTitle = title ? title : "3D-plot of Alpha vs x, y"; fHist.SetDirectory(NULL); fHist.SetName("Alpha"); fHist.SetTitle("3D-plot of Alpha vs x, y"); fHist.SetXTitle("x [\\circ]"); fHist.SetYTitle("y [\\circ]"); fHist.SetZTitle("\\alpha [\\circ]"); } void MHFalseSource::MakeSymmetric(TH1 *h) { const Float_t min = TMath::Abs(h->GetMinimum()); const Float_t max = TMath::Abs(h->GetMaximum()); const Float_t absmax = TMath::Max(min, max)*1.002; h->SetMaximum( absmax); h->SetMinimum(-absmax); } // -------------------------------------------------------------------------- // // Set the alpha cut (|alpha|FindObject("MGeomCam"); if (!geom) { *fLog << err << "MGeomCam not found... aborting." << endl; return kFALSE; } fMm2Deg = geom->GetConvMm2Deg(); MBinning binsa; binsa.SetEdges(18, 0, 90); const MBinning *bins = (MBinning*)plist->FindObject("BinningFalseSource"); if (!bins) { const Float_t r = (geom ? geom->GetMaxRadius()/3 : 200)*fMm2Deg; MBinning b; b.SetEdges(20, -r, r); SetBinning(&fHist, &b, &b, &binsa); } else SetBinning(&fHist, bins, bins, &binsa); fPointPos = (MPointingPos*)plist->FindObject(AddSerialNumber("MPointingPos")); if (!fPointPos) *fLog << warn << "MPointingPos not found... no derotation." << endl; fTime = (MTime*)plist->FindObject(AddSerialNumber("MTime")); if (!fTime) *fLog << warn << "MTime not found... no derotation." << endl; fObservatory = (MObservatory*)plist->FindObject(AddSerialNumber("MObservatory")); if (!fObservatory) *fLog << warn << "MObservatory not found... no derotation." << endl; // FIXME: Because the pointing position could change we must check // for the current pointing position and add a offset in the // Fill function! fRa = fPointPos->GetRa(); fDec = fPointPos->GetDec(); return kTRUE; } // -------------------------------------------------------------------------- // // Fill the histogram. For details see the code or the class description // Bool_t MHFalseSource::Fill(const MParContainer *par, const Stat_t w) { const MHillas *hil = dynamic_cast(par); if (!hil) { *fLog << err << "MHFalseSource::Fill: No container specified!" << endl; return kFALSE; } // Get max radius... const Double_t maxr = 0.98*TMath::Abs(fHist.GetBinCenter(1)); // Get camera rotation angle Double_t rho = 0; if (fTime && fObservatory && fPointPos) rho = fPointPos->RotationAngle(*fObservatory, *fTime); //if (fPointPos) // rho = fPointPos->RotationAngle(*fObservatory); // Create necessary containers for calculation MSrcPosCam src; MHillasSrc hsrc; hsrc.SetSrcPos(&src); // Get number of bins and bin-centers const Int_t nx = fHist.GetNbinsX(); const Int_t ny = fHist.GetNbinsY(); Axis_t cx[nx]; Axis_t cy[ny]; fHist.GetXaxis()->GetCenter(cx); fHist.GetYaxis()->GetCenter(cy); for (int ix=0; ixmaxr) continue; // rotate center of bin // precalculation of sin/cos doesn't accelerate TVector2 v(cx[ix], cy[iy]); if (rho!=0) v=v.Rotate(rho); // convert degrees to millimeters v *= 1./fMm2Deg; src.SetXY(v); // Source dependant hillas parameters if (!hsrc.Calc(hil)) { *fLog << warn << "Calculation of MHillasSrc failed for x=" << cx[ix] << " y=" << cy[iy] << endl; return kFALSE; } // FIXME: This should be replaced by an external MFilter // and/or MTaskList // Source dependant distance cut if (fMinDist>0 && hsrc.GetDist()*fMm2Deg0 && hsrc.GetDist()*fMm2Deg>fMaxDist) continue; if (fMaxLD>0 && hil->GetLength()>fMaxLD*hsrc.GetDist()) continue; if (fMinLD>0 && hil->GetLength()SetTitle(Form("Distribution of %.1f\\circ<|\\alpha|<%.1f\\circ in x,y", cut1, cut2)); // Get projection for range TH2D *p = (TH2D*)fHist.Project3D("yx_off"); // Reset range axe.SetRange(0,9999); // Move contents from projection to h2 h2->Reset(); h2->Add(p, all->GetMaximum()); h2->Divide(all); // Delete p delete p; // Set Minimum as minimum value Greater Than 0 h2->SetMinimum(GetMinimumGT(*h2)); } // -------------------------------------------------------------------------- // // Create projection for on data, taking sum of bin contents of // range (0, fAlphaCut) // void MHFalseSource::ProjectOn(TH2D *h3, TH2D *all) { TAxis &axe = *fHist.GetZaxis(); // Find range to cut axe.SetRangeUser(0, fAlphaCut); const Float_t cut = axe.GetBinUpEdge(axe.GetLast()); h3->SetTitle(Form("Distribution of |\\alpha|<%.1f\\circ in x,y", cut)); // Get projection for range TH2D *p = (TH2D*)fHist.Project3D("yx_on"); // Reset range axe.SetRange(0,9999); // Move contents from projection to h3 h3->Reset(); h3->Add(p, all->GetMaximum()); h3->Divide(all); // Delete p delete p; // Set Minimum as minimum value Greater Than 0 h3->SetMinimum(GetMinimumGT(*h3)); } // -------------------------------------------------------------------------- // // Create projection for all data, taking sum of bin contents of // range (0, 90) - corresponding to the number of entries in this slice. // void MHFalseSource::ProjectAll(TH2D *h3) { h3->SetTitle("Number of entries"); // Get projection for range TH2D *p = (TH2D*)fHist.Project3D("yx_all"); // Move contents from projection to h3 h3->Reset(); h3->Add(p); delete p; // Set Minimum as minimum value Greater Than 0 h3->SetMinimum(GetMinimumGT(*h3)); } // -------------------------------------------------------------------------- // // Update the projections and paint them // void MHFalseSource::Paint(Option_t *opt) { // Set pretty color palette gStyle->SetPalette(1, 0); TVirtualPad *padsave = gPad; TH1D* h1; TH2D* h0; TH2D* h2; TH2D* h3; TH2D* h4; TH2D* h5; // Update projection of all-events padsave->GetPad(2)->cd(3); if ((h0 = (TH2D*)gPad->FindObject("Alpha_yx_all"))) ProjectAll(h0); // Update projection of on-events padsave->GetPad(1)->cd(1); if ((h3 = (TH2D*)gPad->FindObject("Alpha_yx_on"))) ProjectOn(h3, h0); // Update projection of off-events padsave->GetPad(1)->cd(3); if ((h2 = (TH2D*)gPad->FindObject("Alpha_yx_off"))) ProjectOff(h2, h0); padsave->GetPad(2)->cd(2); if ((h5 = (TH2D*)gPad->FindObject("Alpha_yx_diff"))) { h5->Add(h3, h2, -1); MakeSymmetric(h5); } // Update projection of significance padsave->GetPad(1)->cd(2); if (h2 && h3 && (h4 = (TH2D*)gPad->FindObject("Alpha_yx_sig"))) { const Int_t nx = h4->GetXaxis()->GetNbins(); const Int_t ny = h4->GetYaxis()->GetNbins(); //const Int_t nr = nx*nx + ny*ny; Int_t maxx=nx/2; Int_t maxy=ny/2; Int_t max = h4->GetBin(nx, ny); for (int ix=1; ix<=nx; ix++) for (int iy=1; iy<=ny; iy++) { const Int_t n = h4->GetBin(ix, iy); const Double_t s = h3->GetBinContent(n); const Double_t b = h2->GetBinContent(n); const Double_t sig = SignificanceLiMa(s, b); h4->SetBinContent(n, sig); if (sig>h4->GetBinContent(max) && sig>0/* && (ix-nx/2)*(ix-nx/2)+(iy-ny/2)*(iy-ny/2)GetPad(2)->cd(1); if ((h1 = (TH1D*)gPad->FindObject("Alpha")) && max>0) { const Double_t x = h4->GetXaxis()->GetBinCenter(maxx); const Double_t y = h4->GetYaxis()->GetBinCenter(maxy); const Double_t s = h4->GetBinContent(max); // This is because a 3D histogram x and y are vice versa // Than for their projections TH1 *h = fHist.ProjectionZ("Alpha", maxx, maxx, maxy, maxy); h->SetTitle(Form("Distribution of \\alpha for x=%.2f y=%.2f (\\sigma_{max}=%.1f)", x, y, s)); } } gPad = padsave; } // -------------------------------------------------------------------------- // // Get the MAstroCatalog corresponding to fRa, fDec. The limiting magnitude // is set to 9, while the fov is equal to the current fov of the false // source plot. // TObject *MHFalseSource::GetCatalog() { const Double_t maxr = 0.98*TMath::Abs(fHist.GetBinCenter(1)); // Create catalog... MAstroCatalog stars; stars.SetLimMag(9); stars.SetGuiActive(kFALSE); stars.SetRadiusFOV(maxr); stars.SetRaDec(fRa*TMath::DegToRad()*15, fDec*TMath::DegToRad()); stars.ReadBSC("bsc5.dat"); TObject *o = (MAstroCatalog*)stars.Clone(); o->SetBit(kCanDelete); return o; } // -------------------------------------------------------------------------- // // Draw the histogram // void MHFalseSource::Draw(Option_t *opt) { TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this); pad->SetBorderMode(0); AppendPad(""); pad->Divide(1, 2, 0, 0.03); TObject *catalog = GetCatalog(); // Initialize upper part pad->cd(1); gPad->SetBorderMode(0); gPad->Divide(3, 1); // PAD #1 pad->GetPad(1)->cd(1); gPad->SetBorderMode(0); fHist.GetZaxis()->SetRangeUser(0,fAlphaCut); TH1 *h3 = fHist.Project3D("yx_on"); fHist.GetZaxis()->SetRange(0,9999); h3->SetDirectory(NULL); h3->SetXTitle(fHist.GetXaxis()->GetTitle()); h3->SetYTitle(fHist.GetYaxis()->GetTitle()); h3->Draw("colz"); h3->SetBit(kCanDelete); catalog->Draw("mirror same"); // PAD #2 pad->GetPad(1)->cd(2); gPad->SetBorderMode(0); fHist.GetZaxis()->SetRange(0,0); TH1 *h4 = fHist.Project3D("yx_sig"); // Do this to get the correct binning.... fHist.GetZaxis()->SetRange(0,9999); h4->SetTitle("Significance"); h4->SetDirectory(NULL); h4->SetXTitle(fHist.GetXaxis()->GetTitle()); h4->SetYTitle(fHist.GetYaxis()->GetTitle()); h4->Draw("colz"); h4->SetBit(kCanDelete); catalog->Draw("mirror same"); // PAD #3 pad->GetPad(1)->cd(3); gPad->SetBorderMode(0); fHist.GetZaxis()->SetRangeUser(fBgMean-fAlphaCut/2, fBgMean+fAlphaCut/2); TH1 *h2 = fHist.Project3D("yx_off"); h2->SetDirectory(NULL); h2->SetXTitle(fHist.GetXaxis()->GetTitle()); h2->SetYTitle(fHist.GetYaxis()->GetTitle()); h2->Draw("colz"); h2->SetBit(kCanDelete); catalog->Draw("mirror same"); // Initialize lower part pad->cd(2); gPad->SetBorderMode(0); gPad->Divide(3, 1); // PAD #4 pad->GetPad(2)->cd(1); gPad->SetBorderMode(0); TH1 *h1 = fHist.ProjectionZ("Alpha"); h1->SetDirectory(NULL); h1->SetTitle("Distribution of \\alpha"); h1->SetXTitle(fHist.GetZaxis()->GetTitle()); h1->SetYTitle("Counts"); h1->Draw(opt); h1->SetBit(kCanDelete); // PAD #5 pad->GetPad(2)->cd(2); gPad->SetBorderMode(0); TH1 *h5 = (TH1*)h3->Clone("Alpha_yx_diff"); h5->Add(h2, -1); h5->SetTitle("Difference of on- and off-distribution"); h5->SetDirectory(NULL); h5->Draw("colz"); h5->SetBit(kCanDelete); catalog->Draw("mirror same"); // PAD #6 pad->GetPad(2)->cd(3); gPad->SetBorderMode(0); TH1 *h0 = fHist.Project3D("yx_all"); h0->SetDirectory(NULL); h0->SetXTitle(fHist.GetXaxis()->GetTitle()); h0->SetYTitle(fHist.GetYaxis()->GetTitle()); h0->Draw("colz"); h0->SetBit(kCanDelete); catalog->Draw("mirror same"); } // -------------------------------------------------------------------------- // // Everything which is in the main pad belongs to this class! // Int_t MHFalseSource::DistancetoPrimitive(Int_t px, Int_t py) { return 0; } // -------------------------------------------------------------------------- // // Set all sub-pads to Modified() // void MHFalseSource::Modified() { if (!gPad) return; TVirtualPad *padsave = gPad; padsave->Modified(); padsave->GetPad(1)->cd(1); gPad->Modified(); padsave->GetPad(1)->cd(3); gPad->Modified(); padsave->GetPad(2)->cd(1); gPad->Modified(); padsave->GetPad(2)->cd(2); gPad->Modified(); padsave->GetPad(2)->cd(3); gPad->Modified(); gPad->cd(); } // -------------------------------------------------------------------------- // // This is a preliminary implementation of a alpha-fit procedure for // all possible source positions. It will be moved into its own // more powerfull class soon. // // The fit function is "gaus(0)+pol2(3)" which is equivalent to: // [0]*exp(-0.5*((x-[1])/[2])^2) + [3] + [4]*x + [5]*x^2 // or // A*exp(-0.5*((x-mu)/sigma)^2) + a + b*x + c*x^2 // // Parameter [1] is fixed to 0 while the alpha peak should be // symmetric around alpha=0. // // Parameter [4] is fixed to 0 because the first derivative at // alpha=0 should be 0, too. // // In a first step the background is fitted between bgmin and bgmax, // while the parameters [0]=0 and [2]=1 are fixed. // // In a second step the signal region (alphaSetDirectory(0); TH1 *hists = fHist.Project3D("yx_fit"); hists->SetDirectory(0); TH1 *histb = fHist.Project3D("yx_fit"); histb->SetDirectory(0); hist->Reset(); hists->Reset(); histb->Reset(); hist->SetNameTitle("Significance", Form("Fit Region: Signal<%.1f\\circ, %.1f\\circSetName("Excess"); histb->SetName("Background"); hist->SetXTitle(fHist.GetXaxis()->GetTitle()); hists->SetXTitle(fHist.GetXaxis()->GetTitle()); histb->SetXTitle(fHist.GetXaxis()->GetTitle()); hist->SetYTitle(fHist.GetYaxis()->GetTitle()); hists->SetYTitle(fHist.GetYaxis()->GetTitle()); histb->SetYTitle(fHist.GetYaxis()->GetTitle()); const Double_t w = fHist.GetZaxis()->GetBinWidth(1); // xmin, xmax, npar //TF1 func("MyFunc", fcn, 0, 90, 6); // Implementing the function yourself is only about 5% faster TF1 func("", Form("gaus(0) + pol%d(3)", polynom), 0, 90); //TF1 func("", Form("[0]*(TMath::Gaus(x, [1], [2])+TMath::Gaus(x, -[1], [2]))+pol%d(3)", polynom), 0, 90); TArrayD maxpar(func.GetNpar()); /* func.SetParName(0, "A"); * func.SetParName(1, "mu"); * func.SetParName(2, "sigma"); */ func.FixParameter(1, 0); func.FixParameter(4, 0); func.SetParLimits(2, 0, 90); func.SetParLimits(3, -1, 1); const Int_t nx = hist->GetXaxis()->GetNbins(); const Int_t ny = hist->GetYaxis()->GetNbins(); //const Int_t nr = nx*nx+ny*ny; Double_t maxalpha0=0; Double_t maxs=3; Int_t maxx=0; Int_t maxy=0; TStopwatch clk; clk.Start(); *fLog << inf; *fLog << "Signal fit: alpha < " << sigmax << endl; *fLog << "Integration: alpha < " << sigint << endl; *fLog << "Background fit: " << bgmin << " < alpha < " << bgmax << endl; *fLog << "Polynom order: " << (int)polynom << endl; *fLog << "Fitting False Source Plot..." << flush; TH1 *h0 = fHist.Project3D("yx_entries"); Float_t entries = h0->GetMaximum(); delete h0; TH1 *h=0; for (int ix=1; ix<=nx; ix++) for (int iy=1; iy<=ny; iy++) { // This is because a 3D histogram x and y are vice versa // Than for their projections h = fHist.ProjectionZ("AlphaFit", ix, ix, iy, iy); const Double_t alpha0 = h->GetBinContent(1); // Check for the regios which is not filled... if (alpha0==0) continue; h->Scale(entries/h->GetEntries()); if (alpha0>maxalpha0) maxalpha0=alpha0; // First fit a polynom in the off region func.FixParameter(0, 0); func.FixParameter(2, 1); func.ReleaseParameter(3); for (int i=5; iFit(&func, "N0Q", "", bgmin, bgmax); h4a.Fill(func.GetChisquare()); //h5a.Fill(func.GetProb()); //func.SetParLimits(0, 0.5*h->GetBinContent(1), 1.5*h->GetBinContent(1)); //func.SetParLimits(2, 5, 90); func.ReleaseParameter(0); //func.ReleaseParameter(1); func.ReleaseParameter(2); func.FixParameter(3, func.GetParameter(3)); for (int i=5; iFit(&func, "N0Q", "", 0, sigmax); TArrayD p(func.GetNpar(), func.GetParameters()); // Fill results into some histograms h0a.Fill(p[0]); h0b.Fill(p[3]); h1.Fill(p[1]); h2.Fill(p[2]); if (polynom>1) h3.Fill(p[5]); h4b.Fill(func.GetChisquare()); //h5b.Fill(func.GetProb()); // Implementing the integral as analytical function // gives the same result in the order of 10e-5 // and it is not faster at all... //const Bool_t ok = /*p[0]>=0 && /*p[0]1.75;// && p[2]<88.5; /* if (p[0]<-fabs(alpha0-p[3])*7 && p[2]GetBin(ix, iy); hists->SetBinContent(n, s-b); histb->SetBinContent(n, b); hist->SetBinContent(n, sig); if (sig!=0) h6.Fill(sig); if (sig>maxs/* && (ix-nx/2)*(ix-nx/2)+(iy-ny/2)*(iy-ny/2)SetRangeUser(0, maxalpha0*1.5); h0b.GetXaxis()->SetRangeUser(0, maxalpha0*1.5); hists->SetTitle(Form("Excess events for \\alpha<%.0f\\circ (N_{max}=%d)", sigint, (int)hists->GetMaximum())); histb->SetTitle(Form("Background events for \\alpha<%.0f\\circ", sigint)); //hists->SetMinimum(GetMinimumGT(*hists)); histb->SetMinimum(GetMinimumGT(*histb)); MakeSymmetric(hists); MakeSymmetric(hist); clk.Stop(); clk.Print("m"); TCanvas *c=new TCanvas; gStyle->SetPalette(1, 0); c->Divide(3,2, 0, 0); c->cd(1); gPad->SetBorderMode(0); hists->Draw("colz"); hists->SetBit(kCanDelete); catalog->Draw("mirror same"); c->cd(2); gPad->SetBorderMode(0); hist->Draw("colz"); hist->SetBit(kCanDelete); catalog->Draw("mirror same"); c->cd(3); gPad->SetBorderMode(0); histb->Draw("colz"); histb->SetBit(kCanDelete); catalog->Draw("mirror same"); c->cd(4); gPad->Divide(1,3, 0, 0); TVirtualPad *p=gPad; p->SetBorderMode(0); p->cd(1); gPad->SetBorderMode(0); h0b.DrawCopy(); h0a.DrawCopy("same"); p->cd(2); gPad->SetBorderMode(0); h3.DrawCopy(); p->cd(3); gPad->SetBorderMode(0); h2.DrawCopy(); c->cd(6); gPad->Divide(1,2, 0, 0); TVirtualPad *q=gPad; q->SetBorderMode(0); q->cd(1); gPad->SetBorderMode(0); gPad->SetBorderMode(0); h4b.DrawCopy(); h4a.DrawCopy("same"); h6.DrawCopy("same"); q->cd(2); gPad->SetBorderMode(0); //h5b.DrawCopy(); //h5a.DrawCopy("same"); c->cd(5); gPad->SetBorderMode(0); if (maxx>0 && maxy>0) { const char *title = Form(" \\alpha for x=%.2f y=%.2f (\\sigma_{max}=%.1f) ", hist->GetXaxis()->GetBinCenter(maxx), hist->GetYaxis()->GetBinCenter(maxy), maxs); TH1 *result = fHist.ProjectionZ("AlphaFit", maxx, maxx, maxy, maxy); result->Scale(entries/h->GetEntries()); result->SetDirectory(NULL); result->SetNameTitle("Result \\alpha", title); result->SetBit(kCanDelete); result->SetXTitle("\\alpha [\\circ]"); result->SetYTitle("Counts"); result->Draw(); TF1 f1("", func.GetExpFormula(), 0, 90); TF1 f2("", func.GetExpFormula(), 0, 90); f1.SetParameters(maxpar.GetArray()); f2.SetParameters(maxpar.GetArray()); f2.FixParameter(0, 0); f2.FixParameter(1, 0); f2.FixParameter(2, 1); f1.SetLineColor(kGreen); f2.SetLineColor(kRed); f2.DrawCopy("same"); f1.DrawCopy("same"); TPaveText *leg = new TPaveText(0.35, 0.10, 0.90, 0.35, "brNDC"); leg->SetBorderSize(1); leg->SetTextSize(0.04); leg->AddText(0.5, 0.82, Form("A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + pol%d", polynom))->SetTextAlign(22); //leg->AddText(0.5, 0.82, "A * exp(-(\\frac{x-\\mu}{\\sigma})^{2}/2) + b*x^{2} + a")->SetTextAlign(22); leg->AddLine(0, 0.65, 0, 0.65); leg->AddText(0.06, 0.54, Form("A=%.2f", maxpar[0]))->SetTextAlign(11); leg->AddText(0.06, 0.34, Form("\\sigma=%.2f", maxpar[2]))->SetTextAlign(11); leg->AddText(0.06, 0.14, Form("\\mu=%.2f (fix)", maxpar[1]))->SetTextAlign(11); leg->AddText(0.60, 0.54, Form("a=%.2f", maxpar[3]))->SetTextAlign(11); leg->AddText(0.60, 0.34, Form("b=%.2f (fix)", maxpar[4]))->SetTextAlign(11); if (polynom>1) leg->AddText(0.60, 0.14, Form("c=%.2f", maxpar[5]))->SetTextAlign(11); leg->SetBit(kCanDelete); leg->Draw(); q->cd(2); TGraph *g = new TGraph; g->SetPoint(0, 0, 0); Int_t max=0; Float_t maxsig=0; for (int i=1; i<89; i++) { const Double_t s = f1.Integral(0, (float)i)/w; const Double_t b = f2.Integral(0, (float)i)/w; const Double_t sig = SignificanceLiMa(s, b); g->SetPoint(g->GetN(), i, sig); if (sig>maxsig) { max = i; maxsig = sig; } } g->SetNameTitle("SigVs\\alpha", "Significance vs \\alpha"); g->GetHistogram()->SetXTitle("\\alpha_{0} [\\circ]"); g->GetHistogram()->SetYTitle("Significance"); g->SetBit(kCanDelete); g->Draw("AP"); leg = new TPaveText(0.35, 0.10, 0.90, 0.25, "brNDC"); leg->SetBorderSize(1); leg->SetTextSize(0.1); leg->AddText(Form("\\sigma_{max}=%.1f at \\alpha_{max}=%d\\circ", maxsig, max)); leg->SetBit(kCanDelete); leg->Draw(); } }