peaks.C

Go to the documentation of this file.
00001 // Illustrates how to find peaks in histograms.
00002 // This script generates a random number of gaussian peaks
00003 // on top of a linear background.
00004 // The position of the peaks is found via TSpectrum and injected
00005 // as initial values of parameters to make a global fit.
00006 // The background is computed and drawn on top of the original histogram.
00007 //
00008 // To execute this example, do
00009 //  root > .x peaks.C  (generate 10 peaks by default)
00010 //  root > .x peaks.C++ (use the compiler)
00011 //  root > .x peaks.C++(30) (generates 30 peaks)
00012 //
00013 // To execute only the first part of the script (without fitting)
00014 // specify a negative value for the number of peaks, eg
00015 //  root > .x peaks.C(-20)
00016 //
00017 //Author: Rene Brun
00018 
00019 #include "TCanvas.h"
00020 #include "TMath.h"
00021 #include "TH1.h"
00022 #include "TF1.h"
00023 #include "TRandom.h"
00024 #include "TSpectrum.h"
00025 #include "TVirtualFitter.h"
00026    
00027 Int_t npeaks = 30;
00028 Double_t fpeaks(Double_t *x, Double_t *par) {
00029    Double_t result = par[0] + par[1]*x[0];
00030    for (Int_t p=0;p<npeaks;p++) {
00031       Double_t norm  = par[3*p+2];
00032       Double_t mean  = par[3*p+3];
00033       Double_t sigma = par[3*p+4];
00034       result += norm*TMath::Gaus(x[0],mean,sigma);
00035    }
00036    return result;
00037 }
00038 void peaks(Int_t np=10) {
00039    npeaks = TMath::Abs(np);
00040    TH1F *h = new TH1F("h","test",500,0,1000);
00041    //generate n peaks at random
00042    Double_t par[3000];
00043    par[0] = 0.8;
00044    par[1] = -0.6/1000;
00045    Int_t p;
00046    for (p=0;p<npeaks;p++) {
00047       par[3*p+2] = 1;
00048       par[3*p+3] = 10+gRandom->Rndm()*980;
00049       par[3*p+4] = 3+2*gRandom->Rndm();
00050    }
00051    TF1 *f = new TF1("f",fpeaks,0,1000,2+3*npeaks);
00052    f->SetNpx(1000);
00053    f->SetParameters(par);
00054    TCanvas *c1 = new TCanvas("c1","c1",10,10,1000,900);
00055    c1->Divide(1,2);
00056    c1->cd(1);
00057    h->FillRandom("f",200000);
00058    h->Draw();
00059    TH1F *h2 = (TH1F*)h->Clone("h2");
00060    //Use TSpectrum to find the peak candidates
00061    TSpectrum *s = new TSpectrum(2*npeaks);
00062    Int_t nfound = s->Search(h,2,"",0.10);
00063    printf("Found %d candidate peaks to fit\n",nfound);
00064    //Estimate background using TSpectrum::Background
00065    TH1 *hb = s->Background(h,20,"same");
00066    if (hb) c1->Update();
00067    if (np <0) return;
00068 
00069    //estimate linear background using a fitting method
00070    c1->cd(2);
00071    TF1 *fline = new TF1("fline","pol1",0,1000);
00072    h->Fit("fline","qn");
00073    //Loop on all found peaks. Eliminate peaks at the background level
00074    par[0] = fline->GetParameter(0);
00075    par[1] = fline->GetParameter(1);
00076    npeaks = 0;
00077    Float_t *xpeaks = s->GetPositionX();
00078    for (p=0;p<nfound;p++) {
00079       Float_t xp = xpeaks[p];
00080       Int_t bin = h->GetXaxis()->FindBin(xp);
00081       Float_t yp = h->GetBinContent(bin);
00082       if (yp-TMath::Sqrt(yp) < fline->Eval(xp)) continue;
00083       par[3*npeaks+2] = yp;
00084       par[3*npeaks+3] = xp;
00085       par[3*npeaks+4] = 3;
00086       npeaks++;
00087    }
00088    printf("Found %d useful peaks to fit\n",npeaks);
00089    printf("Now fitting: Be patient\n");
00090    TF1 *fit = new TF1("fit",fpeaks,0,1000,2+3*npeaks);
00091    //we may have more than the default 25 parameters
00092    TVirtualFitter::Fitter(h2,10+3*npeaks);
00093    fit->SetParameters(par);
00094    fit->SetNpx(1000);
00095    h2->Fit("fit");             
00096 }

Generated on Tue Jul 5 15:45:10 2011 for ROOT_528-00b_version by  doxygen 1.5.1