00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046 #include <stdlib.h>
00047 #include "TApplication.h"
00048 #include "TBenchmark.h"
00049 #include "TCanvas.h"
00050 #include "TH2.h"
00051 #include "TF2.h"
00052 #include "TRandom.h"
00053 #include "TSpectrum.h"
00054 #include "TSpectrum2.h"
00055 #include "TStyle.h"
00056 #include "Riostream.h"
00057 #include "TROOT.h"
00058 #include "TMath.h"
00059
00060 Int_t npeaks;
00061 Double_t fpeaks(Double_t *x, Double_t *par) {
00062 Double_t result = par[0] + par[1]*x[0];
00063 for (Int_t p=0;p<npeaks;p++) {
00064 Double_t norm = par[3*p+2];
00065 Double_t mean = par[3*p+3];
00066 Double_t sigma = par[3*p+4];
00067 result += norm*TMath::Gaus(x[0],mean,sigma);
00068 }
00069 return result;
00070 }
00071 Double_t fpeaks2(Double_t *x, Double_t *par) {
00072 Double_t result = 0.1;
00073 for (Int_t p=0;p<npeaks;p++) {
00074 Double_t norm = par[5*p+0];
00075 Double_t mean1 = par[5*p+1];
00076 Double_t sigma1 = par[5*p+2];
00077 Double_t mean2 = par[5*p+3];
00078 Double_t sigma2 = par[5*p+4];
00079 result += norm*TMath::Gaus(x[0],mean1,sigma1)*TMath::Gaus(x[1],mean2,sigma2);
00080 }
00081 return result;
00082 }
00083 void findPeaks(Int_t pmin, Int_t pmax, Int_t &nfound, Int_t &ngood, Int_t &nghost) {
00084 npeaks = (Int_t)gRandom->Uniform(pmin,pmax);
00085 Int_t nbins = 500;
00086 Double_t dxbins = 2;
00087 TH1F *h = new TH1F("h","test",nbins,0,nbins*dxbins);
00088
00089 Double_t par[3000];
00090 par[0] = 0.8;
00091 par[1] = -0.6/1000;
00092 Int_t p,pf;
00093 for (p=0;p<npeaks;p++) {
00094 par[3*p+2] = 1;
00095 par[3*p+3] = 10+gRandom->Rndm()*(nbins-20)*dxbins;
00096 par[3*p+4] = 3+2*gRandom->Rndm();
00097 }
00098 TF1 *f = new TF1("f",fpeaks,0,nbins*dxbins,2+3*npeaks);
00099 f->SetNpx(1000);
00100 f->SetParameters(par);
00101 h->FillRandom("f",200000);
00102 TSpectrum *s = new TSpectrum(4*npeaks);
00103 nfound = s->Search(h,2,"goff");
00104
00105 ngood = 0;
00106 Float_t *xpeaks = s->GetPositionX();
00107 for (p=0;p<npeaks;p++) {
00108 for (Int_t pf=0;pf<nfound;pf++) {
00109 Double_t dx = TMath::Abs(xpeaks[pf] - par[3*p+3]);
00110 if (dx <dxbins) ngood++;
00111 }
00112 }
00113
00114 nghost = 0;
00115 for (pf=0;pf<nfound;pf++) {
00116 Int_t nf=0;
00117 for (Int_t p=0;p<npeaks;p++) {
00118 Double_t dx = TMath::Abs(xpeaks[pf] - par[3*p+3]);
00119 if (dx <dxbins) nf++;
00120 }
00121 if (nf == 0) nghost++;
00122 }
00123 delete f;
00124 delete h;
00125 delete s;
00126 }
00127
00128 void stress1(Int_t ntimes) {
00129 Int_t pmin = 5;
00130 Int_t pmax = 55;
00131 TCanvas *c1 = new TCanvas("c1","Spectrum results",10,10,800,800);
00132 c1->Divide(2,2);
00133 gStyle->SetOptFit();
00134 TH1F *hpeaks = new TH1F("hpeaks","Number of peaks",pmax-pmin,pmin,pmax);
00135 TH1F *hfound = new TH1F("hfound","% peak founds",100,0,100);
00136 TH1F *hgood = new TH1F("hgood", "% good peaks",100,0,100);
00137 TH1F *hghost = new TH1F("hghost","% ghost peaks",100,0,100);
00138 Int_t nfound,ngood,nghost;
00139 for (Int_t i=0;i<ntimes;i++) {
00140 findPeaks(pmin,pmax,nfound,ngood,nghost);
00141 hpeaks->Fill(npeaks);
00142 hfound->Fill(100*Double_t(nfound)/Double_t(npeaks));
00143 hgood->Fill(100*Double_t(ngood)/Double_t(npeaks));
00144 hghost->Fill(100*Double_t(nghost)/Double_t(npeaks));
00145
00146 }
00147 c1->cd(1);
00148 hpeaks->Fit("pol1","lq");
00149 c1->cd(2);
00150 hfound->Fit("gaus","lq");
00151 c1->cd(3);
00152 hgood->Fit("gaus","lq");
00153 c1->cd(4);
00154 hghost->Fit("gaus","lq","",0,30);
00155 c1->cd();
00156 Double_t p1 = hfound->GetFunction("gaus")->GetParameter(1);
00157 Double_t ep1 = hfound->GetFunction("gaus")->GetParError(1);
00158 Double_t p2 = hgood->GetFunction("gaus")->GetParameter(1);
00159 Double_t ep2 = hgood->GetFunction("gaus")->GetParError(1);
00160 Double_t p3 = hghost->GetFunction("gaus")->GetParameter(1);
00161 Double_t ep3 = hghost->GetFunction("gaus")->GetParError(1);
00162 Double_t p1ref = 70.21;
00163 Double_t p2ref = 65.03;
00164 Double_t p3ref = 8.54;
00165
00166
00167
00168 char sok[20];
00169 if (TMath::Abs(p1ref-p1) < 2*ep1 && TMath::Abs(p2ref-p2) < 2*ep2 && TMath::Abs(p3ref-p3) < 2*ep3 ) {
00170 sprintf(sok,"OK");
00171 } else {
00172 sprintf(sok,"failed");
00173 }
00174 printf("Peak1 : found =%6.2f/%6.2f, good =%6.2f/%6.2f, ghost =%5.2f/%5.2f,--- %s\n",
00175 p1,p1ref,p2,p2ref,p3,p3ref,sok);
00176 }
00177 void stress2(Int_t np2) {
00178 npeaks = np2;
00179 TRandom r;
00180 Int_t nbinsx = 200;
00181 Int_t nbinsy = 200;
00182 Double_t xmin = 0;
00183 Double_t xmax = (Double_t)nbinsx;
00184 Double_t ymin = 0;
00185 Double_t ymax = (Double_t)nbinsy;
00186 Double_t dx = (xmax-xmin)/nbinsx;
00187 Double_t dy = (ymax-ymin)/nbinsy;
00188 TH2F *h2 = new TH2F("h2","test",nbinsx,xmin,xmax,nbinsy,ymin,ymax);
00189 h2->SetStats(0);
00190
00191 Double_t par[3000];
00192 Int_t p;
00193 for (p=0;p<npeaks;p++) {
00194 par[5*p+0] = r.Uniform(0.2,1);
00195 par[5*p+1] = r.Uniform(xmin,xmax);
00196 par[5*p+2] = r.Uniform(dx,5*dx);
00197 par[5*p+3] = r.Uniform(ymin,ymax);
00198 par[5*p+4] = r.Uniform(dy,5*dy);
00199 }
00200 TF2 *f2 = new TF2("f2",fpeaks2,xmin,xmax,ymin,ymax,5*npeaks);
00201 f2->SetNpx(100);
00202 f2->SetNpy(100);
00203 f2->SetParameters(par);
00204 h2->FillRandom("f2",500000);
00205
00206 TSpectrum2 *s = new TSpectrum2(2*npeaks);
00207 Int_t nfound = s->Search(h2,2,"goff noMarkov");
00208
00209
00210 Int_t pf,ngood = 0;
00211 Float_t *xpeaks = s->GetPositionX();
00212 Float_t *ypeaks = s->GetPositionY();
00213 for (p=0;p<npeaks;p++) {
00214 for (Int_t pf=0;pf<nfound;pf++) {
00215 Double_t diffx = TMath::Abs(xpeaks[pf] - par[5*p+1]);
00216 Double_t diffy = TMath::Abs(ypeaks[pf] - par[5*p+3]);
00217 if (diffx < 2*dx && diffy < 2*dy) ngood++;
00218 }
00219 }
00220 if (ngood > nfound) ngood = nfound;
00221
00222 Int_t nghost = 0;
00223 for (pf=0;pf<nfound;pf++) {
00224 Int_t nf=0;
00225 for (Int_t p=0;p<npeaks;p++) {
00226 Double_t diffx = TMath::Abs(xpeaks[pf] - par[5*p+1]);
00227 Double_t diffy = TMath::Abs(ypeaks[pf] - par[5*p+3]);
00228 if (diffx < 2*dx && diffy < 2*dy) nf++;
00229 }
00230 if (nf == 0) nghost++;
00231 }
00232
00233 delete s;
00234 delete f2;
00235 delete h2;
00236 Int_t nfoundRef = 163;
00237 Int_t ngoodRef = 163;
00238 Int_t nghostRef = 8;
00239 char sok[20];
00240 if ( TMath::Abs(nfound - nfoundRef) < 5
00241 && TMath::Abs(ngood - ngoodRef) < 5
00242 && TMath::Abs(nghost - nghostRef) < 5) {
00243 sprintf(sok,"OK");
00244 } else {
00245 sprintf(sok,"failed");
00246 }
00247 printf("Peak2 : found =%d/%d, good =%d, ghost =%2d,---------------------------- %s\n",
00248 nfound,npeaks,ngood,nghost,sok);
00249 }
00250
00251 #ifndef __CINT__
00252 void stressSpectrum(Int_t ntimes) {
00253 #else
00254 void stressSpectrum(Int_t ntimes=100) {
00255 #endif
00256 cout << "****************************************************************************" <<endl;
00257 cout << "* Starting stress S P E C T R U M *" <<endl;
00258 cout << "****************************************************************************" <<endl;
00259 gBenchmark->Start("stressSpectrum");
00260 stress1(ntimes);
00261 stress2(300);
00262 gBenchmark->Stop ("stressSpectrum");
00263 Double_t reftime100 = 19.04;
00264 Double_t ct = gBenchmark->GetCpuTime("stressSpectrum");
00265 const Double_t rootmarks = 800*reftime100*ntimes/(100*ct);
00266 printf("****************************************************************************\n");
00267
00268 gBenchmark->Print("stressSpectrum");
00269 printf("****************************************************************************\n");
00270 printf("* ROOTMARKS =%6.1f * Root%-8s %d/%d\n",rootmarks,gROOT->GetVersion(),
00271 gROOT->GetVersionDate(),gROOT->GetVersionTime());
00272 printf("****************************************************************************\n");
00273 }
00274
00275 #ifndef __CINT__
00276
00277 int main(int argc, char **argv)
00278 {
00279 TApplication theApp("App", &argc, argv);
00280 gROOT->SetBatch();
00281 gBenchmark = new TBenchmark();
00282 Int_t ntimes = 100;
00283 if (argc > 1) ntimes = atoi(argv[1]);
00284 stressSpectrum(ntimes);
00285 return 0;
00286 }
00287
00288 #endif