00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include <TMath.h>
00011 #include <TCanvas.h>
00012 #include <TRandom3.h>
00013 #include <TFitter.h>
00014 #include <TF1.h>
00015 #include <TStyle.h>
00016 #include <TVector.h>
00017 #include <TGraph.h>
00018
00019 #include <TUnfoldSys.h>
00020
00021 using namespace std;
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038 TRandom *rnd=0;
00039
00040 Int_t GenerateGenEvent(Int_t nmax,const Double_t *probability) {
00041
00042
00043
00044
00045 Double_t f=rnd->Rndm();
00046 Int_t r=0;
00047 while((r<nmax)&&(f>=probability[r])) {
00048 f -= probability[r];
00049 r++;
00050 }
00051 return r;
00052 }
00053
00054 Double_t GenerateRecEvent(const Double_t *shapeParm) {
00055
00056
00057
00058
00059
00060
00061
00062
00063 Double_t f=rnd->Rndm();
00064 Double_t r;
00065 if(f<shapeParm[0]) {
00066 r=rnd->Gaus(shapeParm[1],shapeParm[2]);
00067 } else {
00068 r=rnd->Rndm()*(shapeParm[4]-shapeParm[3])+shapeParm[3];
00069 }
00070 return r;
00071 }
00072
00073 void testUnfold4()
00074 {
00075
00076 TH1::SetDefaultSumw2();
00077
00078
00079 rnd=new TRandom3();
00080
00081
00082 Double_t const nData0= 1000.0;
00083 Double_t const nMC0 = 50000.0;
00084
00085
00086
00087 Int_t const nDet=15;
00088 Double_t const xminDet=0.0;
00089 Double_t const xmaxDet=15.0;
00090
00091
00092 Int_t const nGen=3;
00093 Double_t const xminGen=-0.5;
00094 Double_t const xmaxGen= 2.5;
00095
00096
00097
00098 static const Double_t genFrac[]={0.4,0.4,0.2};
00099
00100
00101 static const Double_t genShape[][5]=
00102 {{1.0,2.0,1.5,0.,15.},
00103 {1.0,7.0,2.5,0.,15.},
00104 {0.0,0.0,0.0,0.,15.}};
00105
00106
00107
00108 TH1D *histDetDATA=new TH1D("Yrec",";DATA(Yrec)",nDet,xminDet,xmaxDet);
00109
00110
00111
00112 TH2D *histGenDetMC=new TH2D("Yrec%Xgen","MC(Xgen,Yrec)",
00113 nGen,xminGen,xmaxGen,nDet,xminDet,xmaxDet);
00114
00115 TH1D *histUnfold=new TH1D("Xgen",";DATA(Xgen)",nGen,xminGen,xmaxGen);
00116
00117 TH1D **histPullNC=new TH1D* [nGen];
00118 TH1D **histPullArea=new TH1D* [nGen];
00119 for(int i=0;i<nGen;i++) {
00120 histPullNC[i]=new TH1D(TString::Format("PullNC%d",i),"pull",15,-3.,3.);
00121 histPullArea[i]=new TH1D(TString::Format("PullArea%d",i),"pull",15,-3.,3.);
00122 }
00123
00124
00125 cout<<"TUnfold version is "<<TUnfold::GetTUnfoldVersion()<<"\n";
00126
00127 for(int itoy=0;itoy<1000;itoy++) {
00128 histDetDATA->Reset();
00129 histGenDetMC->Reset();
00130
00131 Int_t nData=rnd->Poisson(nData0);
00132 for(Int_t i=0;i<nData;i++) {
00133 Int_t iGen=GenerateGenEvent(nGen,genFrac);
00134 Double_t yObs=GenerateRecEvent(genShape[iGen]);
00135 histDetDATA->Fill(yObs);
00136 }
00137
00138 Int_t nMC=rnd->Poisson(nMC0);
00139 for(Int_t i=0;i<nMC;i++) {
00140 Int_t iGen=GenerateGenEvent(nGen,genFrac);
00141 Double_t yObs=GenerateRecEvent(genShape[iGen]);
00142 histGenDetMC->Fill(iGen,yObs);
00143 }
00144
00145
00146
00147
00148
00149
00150
00151
00152 TUnfoldSys unfold(histGenDetMC,TUnfold::kHistMapOutputHoriz,
00153 TUnfold::kRegModeSize,TUnfold::kEConstraintNone);
00154
00155 unfold.SetInput(histDetDATA);
00156
00157
00158 unfold.ScanLcurve(50,0.,0.,0,0,0);
00159
00160
00161 unfold.GetOutput(histUnfold);
00162
00163 for(int i=0;i<nGen;i++) {
00164 histPullNC[i]->Fill((histUnfold->GetBinContent(i+1)-genFrac[i]*nData0)/
00165 histUnfold->GetBinError(i+1));
00166 }
00167
00168
00169 unfold.SetConstraint(TUnfold::kEConstraintArea);
00170
00171
00172 unfold.ScanLcurve(50,0.,0.,0,0,0);
00173
00174
00175 unfold.GetOutput(histUnfold);
00176
00177 for(int i=0;i<nGen;i++) {
00178 histPullArea[i]->Fill((histUnfold->GetBinContent(i+1)-genFrac[i]*nData0)/
00179 histUnfold->GetBinError(i+1));
00180 }
00181
00182 }
00183 TCanvas output;
00184 output.Divide(3,2);
00185
00186 for(int i=0;i<nGen;i++) {
00187 output.cd(i+1);
00188 histPullNC[i]->Fit("gaus");
00189 histPullNC[i]->Draw();
00190 }
00191 for(int i=0;i<nGen;i++) {
00192 output.cd(i+4);
00193 histPullArea[i]->Fit("gaus");
00194 histPullArea[i]->Draw();
00195 }
00196 output.SaveAs("testUnfold4.ps");
00197 }