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 #include <TMath.h>
00027 #include <TCanvas.h>
00028 #include <TRandom3.h>
00029 #include <TFitter.h>
00030 #include <TF1.h>
00031 #include <TStyle.h>
00032 #include <TVector.h>
00033 #include <TGraph.h>
00034 #include <TUnfold.h>
00035
00036 using namespace std;
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061 TRandom *rnd=0;
00062
00063
00064
00065
00066
00067 Double_t GenerateEvent(Double_t bgr,
00068 Double_t mass,
00069 Double_t gamma)
00070 {
00071 Double_t t;
00072 if(rnd->Rndm()>bgr) {
00073
00074
00075 do {
00076 do {
00077 t=rnd->Rndm();
00078 } while(t>=1.0);
00079 t=TMath::Tan((t-0.5)*TMath::Pi())*gamma+mass;
00080 } while(t<=0.0);
00081 return t;
00082 } else {
00083
00084
00085
00086 static Double_t const E0=2.4;
00087 static Double_t const N0=2.9;
00088 do {
00089 do {
00090 t=rnd->Rndm();
00091 } while(t>=1.0);
00092
00093
00094 t= -(TMath::Power(1.-t,1./(1.-N0))-1.0)*E0;
00095 } while(t>=0.0);
00096 return t;
00097 }
00098 }
00099
00100
00101
00102
00103
00104
00105 Double_t DetectorEvent(Double_t mTrue) {
00106
00107 static Double_t frac=0.1;
00108 static Double_t wideBias=0.03;
00109 static Double_t wideSigma=0.5;
00110 static Double_t smallBias=0.0;
00111 static Double_t smallSigma=0.1;
00112 if(rnd->Rndm()>frac) {
00113 return rnd->Gaus(mTrue+smallBias,smallSigma);
00114 } else {
00115 return rnd->Gaus(mTrue+wideBias,wideSigma);
00116 }
00117 }
00118
00119 int testUnfold2()
00120 {
00121
00122 TH1::SetDefaultSumw2();
00123
00124
00125 rnd=new TRandom3();
00126
00127
00128 Double_t const luminosityData=100000;
00129 Double_t const luminosityMC=1000000;
00130 Double_t const crossSection=1.0;
00131
00132 Int_t const nDet=250;
00133 Int_t const nGen=100;
00134 Double_t const xminDet=0.0;
00135 Double_t const xmaxDet=10.0;
00136 Double_t const xminGen=0.0;
00137 Double_t const xmaxGen=10.0;
00138
00139
00140
00141
00142 TH1D *histMgenMC=new TH1D("MgenMC",";mass(gen)",nGen,xminGen,xmaxGen);
00143 TH1D *histMdetMC=new TH1D("MdetMC",";mass(det)",nDet,xminDet,xmaxDet);
00144 TH2D *histMdetGenMC=new TH2D("MdetgenMC",";mass(det);mass(gen)",nDet,xminDet,xmaxDet,
00145 nGen,xminGen,xmaxGen);
00146 Int_t neventMC=rnd->Poisson(luminosityMC*crossSection);
00147 for(Int_t i=0;i<neventMC;i++) {
00148 Double_t mGen=GenerateEvent(0.3,
00149 4.0,
00150 0.2);
00151 Double_t mDet=DetectorEvent(TMath::Abs(mGen));
00152
00153
00154
00155
00156
00157
00158
00159
00160 histMgenMC->Fill(mGen,luminosityData/luminosityMC);
00161
00162 histMdetMC->Fill(mDet,luminosityData/luminosityMC);
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172
00173
00174
00175
00176
00177 histMdetGenMC->Fill(mDet,mGen,luminosityData/luminosityMC);
00178 }
00179
00180
00181
00182
00183 TH1D *histMgenData=new TH1D("MgenData",";mass(gen)",nGen,xminGen,xmaxGen);
00184 TH1D *histMdetData=new TH1D("MdetData",";mass(det)",nDet,xminDet,xmaxDet);
00185 Int_t neventData=rnd->Poisson(luminosityData*crossSection);
00186 for(Int_t i=0;i<neventData;i++) {
00187 Double_t mGen=GenerateEvent(0.4,
00188 3.8,
00189 0.15);
00190 Double_t mDet=DetectorEvent(TMath::Abs(mGen));
00191
00192
00193 histMgenData->Fill(mGen);
00194
00195
00196 histMdetData->Fill(mDet);
00197 }
00198
00199
00200
00201 TUnfold unfold(histMdetGenMC,TUnfold::kHistMapOutputVert,
00202 TUnfold::kRegModeNone);
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217 Double_t estimatedPeakPosition=3.8;
00218 Int_t nPeek=3;
00219 TUnfold::ERegMode regMode=TUnfold::kRegModeCurvature;
00220
00221 Int_t iPeek=(Int_t)(nGen*(estimatedPeakPosition-xminGen)/(xmaxGen-xminGen)
00222
00223
00224
00225 +1.5);
00226
00227 unfold.RegularizeBins(1,1,iPeek-nPeek,regMode);
00228
00229 unfold.RegularizeBins(iPeek+nPeek,1,nGen-(iPeek+nPeek),regMode);
00230
00231
00232
00233
00234
00235 if(unfold.SetInput(histMdetData,0.0)>=10000) {
00236 std::cout<<"Unfolding result may be wrong\n";
00237 }
00238
00239
00240 Double_t tauMin=0.0;
00241 Double_t tauMax=0.0;
00242 Int_t nScan=30;
00243 Int_t iBest;
00244 TSpline *logTauX,*logTauY;
00245 TGraph *lCurve;
00246
00247
00248 iBest=unfold.ScanLcurve(nScan,tauMin,tauMax,&lCurve,&logTauX,&logTauY);
00249 std::cout<<"tau="<<unfold.GetTau()<<"\n";
00250 std::cout<<"chi**2="<<unfold.GetChi2A()<<"+"<<unfold.GetChi2L()
00251 <<" / "<<unfold.GetNdf()<<"\n";
00252
00253
00254 Double_t t[1],x[1],y[1];
00255 logTauX->GetKnot(iBest,t[0],x[0]);
00256 logTauY->GetKnot(iBest,t[0],y[0]);
00257 TGraph *bestLcurve=new TGraph(1,x,y);
00258 TGraph *bestLogTauX=new TGraph(1,t,x);
00259
00260
00261
00262
00263
00264
00265
00266 Int_t *binMap=new Int_t[nGen+2];
00267 for(Int_t i=1;i<=nGen;i++) binMap[i]=i;
00268 binMap[0]=-1;
00269 binMap[nGen+1]=-1;
00270
00271 TH1D *histMunfold=new TH1D("Unfolded",";mass(gen)",nGen,xminGen,xmaxGen);
00272 unfold.GetOutput(histMunfold,binMap);
00273 TH1D *histMdetFold=unfold.GetFoldedOutput("FoldedBack","mass(det)",
00274 xminDet,xmaxDet);
00275
00276
00277 TH1D *histRhoi=new TH1D("rho_I","mass",nGen,xminGen,xmaxGen);
00278 unfold.GetRhoI(histRhoi,0,binMap);
00279
00280 delete[] binMap;
00281 binMap=0;
00282
00283
00284
00285 TCanvas output;
00286
00287
00288 output.Divide(3,2);
00289
00290
00291
00292
00293
00294
00295
00296 output.cd(1);
00297 histMdetGenMC->Draw("BOX");
00298
00299
00300
00301
00302
00303 output.cd(2);
00304 histMunfold->SetLineColor(kBlue);
00305 histMunfold->Draw();
00306 histMgenData->SetLineColor(kRed);
00307 histMgenData->Draw("SAME");
00308 histMgenMC->Draw("SAME HIST");
00309
00310
00311
00312
00313
00314 output.cd(3);
00315 histMdetFold->SetLineColor(kBlue);
00316 histMdetFold->Draw();
00317 histMdetData->SetLineColor(kRed);
00318 histMdetData->Draw("SAME");
00319 histMdetMC->Draw("SAME HIST");
00320
00321
00322
00323
00324
00325
00326 output.cd(4);
00327 histRhoi->Draw();
00328
00329
00330 output.cd(5);
00331 logTauX->Draw();
00332 bestLogTauX->SetMarkerColor(kRed);
00333 bestLogTauX->Draw("*");
00334
00335 output.cd(6);
00336 lCurve->Draw("AL");
00337 bestLcurve->SetMarkerColor(kRed);
00338 bestLcurve->Draw("*");
00339
00340 output.SaveAs("testUnfold2.ps");
00341 return 0;
00342 }