//////////////////////////////////////////////////////////////////////////// // // // This program should be run under root : // // root unfold.C++ // // // // Author(s) : T. Bretz 02/2002 // // Author(s) : W. Wittek 09/2002 // // // //////////////////////////////////////////////////////////////////////////// #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include TH1 *DrawMatrixClone(const TMatrixD &m, Option_t *opt="") { const Int_t nrows = m.GetNrows(); const Int_t ncols = m.GetNcols(); TMatrix m2(nrows, ncols); for (int i=0; iSetBit(kCanDelete); hist->Draw(opt); hist->SetDirectory(NULL); return hist; } TH1 *DrawMatrixColClone(const TMatrixD &m, Option_t *opt="", Int_t col=0) { const Int_t nrows = m.GetNrows(); TVector vec(nrows); for (int i=0; iSetBinContent(i+1, vec(i)); } hist->SetBit(kCanDelete); hist->Draw(opt); hist->SetDirectory(NULL); return hist; } void PrintTH2Content(const TH2 &hist) { cout << hist.GetName() << ": " << hist.GetTitle() << endl; cout << "-----------------------------------------------------" << endl; for (Int_t i=1; i<=hist.GetNbinsX(); i++) for (Int_t j=1; j<=hist.GetNbinsY(); j++) cout << hist.GetBinContent(i,j) << " \t"; cout << endl << endl; } void PrintTH2Error(const TH2 &hist) { cout << hist.GetName() << ": " << hist.GetTitle() << " " << endl; cout << "-----------------------------------------------------" << endl; for (Int_t i=1; i<=hist.GetNbinsX(); i++) { for (Int_t j=1; j<=hist.GetNbinsY(); j++) cout << hist.GetBinError(i, j) << " \t"; cout << endl << endl; } } void PrintTH1Content(const TH1 &hist) { cout << hist.GetName() << ": " << hist.GetTitle() << endl; cout << "-----------------------------------------------------" << endl; for (Int_t i=1; i<=hist.GetNbinsX(); i++) cout << hist.GetBinContent(i) << " \t"; cout << endl << endl; } void PrintTH1Error(const TH1 &hist) { cout << hist.GetName() << ": " << hist.GetTitle() << " " << endl; cout << "-----------------------------------------------------" << endl; for (Int_t i=1; i<=hist.GetNbinsX(); i++) cout << hist.GetBinError(i) << " \t"; cout << endl << endl; } void CopyCol(TMatrixD &m, const TH1 &h, Int_t col=0) { const Int_t n = m.GetNrows(); for (Int_t i=0; iGetNbins(),1), fVEps0(fVEps, 0) { // ha is the distribution to be unfolded // hacov is the covariance matrix of ha // hmig is the migration matrix; // this matrix will be used in the unfolding // unless SmoothMigrationMatrix(*hmigrat) is called; // in the latter case hmigrat is smoothed // and the smoothed matrix is used in the unfolding // Eigen values of the matrix G, which are smaller than EpsLambda // will be considered as being zero EpsLambda = 1.e-10; fW = 0.0; fNa = hmig.GetXaxis()->GetNbins(); const Double_t alow = hmig.GetXaxis()->GetXmin(); const Double_t aup = hmig.GetXaxis()->GetXmax(); fNb = hmig.GetYaxis()->GetNbins(); const Double_t blow = hmig.GetYaxis()->GetXmin(); const Double_t bup = hmig.GetYaxis()->GetXmax(); UInt_t Na = ha.GetNbinsX(); if (fNa != Na) { cout << "MUnfold::MUnfold : dimensions do not match, fNa = "; cout << fNa << ", Na = " << Na << endl; } cout << "MUnfold::MUnfold :" << endl; cout << "==================" << endl; cout << " fNa = " << fNa << ", fNb = " << fNb << endl; // ------------------------ fVa.ResizeTo(fNa, 1); CopyCol(fVa, ha, 0); cout << " fVa = "; for (UInt_t i=0; iSumw2(); //------------------------------------- // smoothed migration matrix shmig = new TH2D("sMigrat", "Smoothed migration matrix", fNa, alow, aup, fNb, blow, bup); shmig->Sumw2(); //------------------------------------- // chi2 contributions for smoothing of migration matrix shmigChi2 = new TH2D("sMigratChi2", "Chi2 contr. for smoothing", fNa, alow, aup, fNb, blow, bup); //------------------------------------- // eigen values of matrix G = M * M(transposed) hEigen = new TH1D("Eigen", "Eigen values of M*MT", fNa, 0.5, fNa+0.5); //------------------------------------ // Ideal distribution fhb0 = new TH1D("fhb0", "Ideal distribution", fNb, blow, bup); fhb0->Sumw2(); //------------------------------------ // Distribution to be unfolded fha = new TH1D("fha", "Distribution to be unfolded", fNa, alow, aup); fha->Sumw2(); //------------------------------------ // Prior distribution hprior = new TH1D("Prior", "Prior distribution", fNb, blow, bup); //------------------------------------ // Unfolded distribution hb = new TH1D("DataSp", "Unfolded distribution", fNb, blow, bup); hb->Sumw2(); } // ----------------------------------------------------------------------- // // Define prior distribution to be a constant // void SetPriorConstant() { fVEps0 = 1./fNb; CopyCol(*hprior, fVEps); cout << "SetPriorConstant : Prior distribution fVEps = " << endl; cout << "==============================================" << endl; fVEps.Print(); } // ----------------------------------------------------------------------- // // Take prior distribution from the histogram 'ha' // which may have a different binning than 'hprior' // Bool_t SetPriorRebin(TH1D &ha) { // ------------------------------------------------------------------ // // fill the contents of histogram 'ha' into the histogram 'hrior'; // the histograms need not have the same binning; // if the binnings are different, the bin contents of histogram 'ha' // are distributed properly (linearly) over the bins of 'hprior' // const Int_t na = ha.GetNbinsX(); const Double_t alow = ha.GetBinLowEdge(1); const Double_t aup = ha.GetBinLowEdge(na+1); const Int_t nb = hprior->GetNbinsX(); const Double_t blow = hprior->GetBinLowEdge(1); const Double_t bup = hprior->GetBinLowEdge(nb+1); // check whether there is an overlap // between the x ranges of the 2 histograms if (alow>bup || aupGetBinLowEdge(j); const Double_t yh = hprior->GetBinLowEdge(j+1); // search bins of histogram ha which contribute // to bin j of histogram hb //---------------- Int_t il=0; Int_t ih=0; for (Int_t i=2; i<=na+1; i++) { const Double_t xl = ha.GetBinLowEdge(i); if (xl>yl) { il = i-1; //................................. ih = 0; for (Int_t k=(il+1); k<=(na+1); k++) { const Double_t xh = ha.GetBinLowEdge(k); if (xh >= yh) { ih = k-1; break; } } //................................. if (ih == 0) ih = na; break; } } //---------------- if (il == 0) { cout << "Something is wrong " << endl; cout << " na, alow, aup = " << na << ", " << alow << ", " << aup << endl; cout << " nb, blow, bup = " << nb << ", " << blow << ", " << bup << endl; return kFALSE; } Double_t content=0; // sum up the contribution to bin j for (Int_t i=il; i<=ih; i++) { const Double_t xl = ha.GetBinLowEdge(i); const Double_t xh = ha.GetBinLowEdge(i+1); const Double_t bina = xh-xl; if (xl=yh) content += ha.GetBinContent(i) * (yh-yl) / bina; else if (xl>=yl && xh=yl && xh>=yh) content += ha.GetBinContent(i) * (yh-xl) / bina; } hprior->SetBinContent(j, content); sum += content; } // normalize histogram hb if (sum==0) { cout << "histogram hb is empty; sum of weights in ha = "; cout << ha.GetSumOfWeights() << endl; return kFALSE; } for (Int_t j=1; j<=nb; j++) { const Double_t content = hprior->GetBinContent(j)/sum; hprior->SetBinContent(j, content); fVEps0(j-1) = content; } cout << "SetPriorRebin : Prior distribution fVEps = " << endl; cout << "===========================================" << endl; fVEps.Print(); return kTRUE; } // ----------------------------------------------------------------------- // // Set prior distribution to a given distribution 'hpr' // Bool_t SetPriorInput(TH1D &hpr) { CopyCol(fVEps, hpr); const Double_t sum = GetMatrixSumCol(fVEps, 0); if (sum<=0) { cout << "MUnfold::SetPriorInput: invalid prior distribution" << endl; return kFALSE; } // normalize prior distribution fVEps0 *= 1./sum; CopyCol(*hprior, fVEps); cout << "SetPriorInput : Prior distribution fVEps = " << endl; cout << "===========================================" << endl; fVEps.Print(); return kTRUE; } // ----------------------------------------------------------------------- // // Define prior distribution to be a power law // use input distribution 'hprior' only // for defining the histogram parameters // Bool_t SetPriorPower(Double_t gamma) { // generate distribution according to a power law // dN/dE = E^{-gamma} // or with y = lo10(E), E = 10^y : // dN/dy = ln10 * 10^{y*(1-gamma)} TH1D hpower(*hprior); const UInt_t nbin = hprior->GetNbinsX(); const Double_t xmin = hprior->GetBinLowEdge(1); const Double_t xmax = hprior->GetBinLowEdge(nbin+1); cout << "nbin, xmin, xmax = " << nbin << ", "; cout << xmin << ", " << xmax << endl; TF1* fpow = new TF1("fpow", "pow(10.0, x*(1.0-[0]))", xmin,xmax); fpow->SetParName (0,"gamma"); fpow->SetParameter(0, gamma ); hpower.FillRandom("fpow", 100000); // fill prior distribution CopyCol(fVEps, hpower); const Double_t sum = GetMatrixSumCol(fVEps, 0); if (sum <= 0) { cout << "MUnfold::SetPriorPower : invalid prior distribution" << endl; return kFALSE; } // normalize prior distribution fVEps0 *= 1./sum; CopyCol(*hprior, fVEps); cout << "SetPriorPower : Prior distribution fVEps = " << endl; cout << "===========================================" << endl; fVEps.Print(); return kTRUE; } // ----------------------------------------------------------------------- // // Set the initial weight // Bool_t SetInitialWeight(Double_t &weight) { if (weight == 0.0) { TMatrixD v1(fVa, TMatrixD::kTransposeMult, fVacovInv); TMatrixD v2(v1, TMatrixD::kMult, fVa); weight = 1./sqrt(v2(0,0)); } cout << "MUnfold::SetInitialWeight : Initial Weight = " << weight << endl; return kTRUE; } // ----------------------------------------------------------------------- // // Print the unfolded distribution // void PrintResults() { cout << "PrintResults : Unfolded distribution fResult " << endl; cout << "=============================================" << endl; cout << "val, eparab, eplus, eminus, gcc = " << endl; for (UInt_t i=0; iSetBinContent(i, hb0.GetBinContent(i)); fhb0->SetBinError (i, hb0.GetBinError(i)); } //----------------------------------------------------------------------- // Initialization // ============== Int_t numGiteration; Int_t MaxGiteration = 1000; TMatrixD alpha; alpha.ResizeTo(fNa, 1); //----------------------------------------------------------------------- // Newton iteration // ================ Double_t dga2; Double_t dga2old; Double_t EpsG = 1.e-12; TMatrixD wZdp_inv(fNa, fNa); TMatrixD d(fNb, 1); TMatrixD p(fNb, 1); TMatrixD gamma (fNa, 1); TMatrixD dgamma(fNa, 1); Double_t fWinitial; fWinitial = 0.0; SetInitialWeight(fWinitial); // for my example this fWinitial was not good; therefore : fWinitial = 1.0; Int_t ix; Double_t xiter; //-------- start scanning weights -------------------------- // if full == kFALSE only quantities necessary for the Gauss-Newton // iteration are calculated in SchmellCore // if full == kTRUE in addition the unfolded distribution, // its covariance matrix and quantities like // Chisq, SpurAR, etc. are computed in SchmellCore //Bool_t full; //full = kFALSE; Int_t full; dga2 = 1.e20; for (ix=0; ix EpsG) { gamma0 = 1; gamma.Zero(); } dga2 = 1.e20; while (1) { dga2old = dga2; full = 0; SchmellCore(full, xiter, gamma, dgamma, dga2); gamma += dgamma; //cout << "Schmelling : ix, numGiteration, dga2, dga2old = " // << ix << ", " << numGiteration << ", " // << dga2 << ", " << dga2old << endl; numGiteration += 1; // convergence if (dga2 < EpsG) break; // no convergence if (numGiteration > MaxGiteration) break; // gamma doesn't seem to change any more if (fabs(dga2-dga2old) < EpsG/100.) break; } //---------- end Gauss-Newton iteration ------------------------ if (dga2 EpsG) { cout << "Schmelling : Gauss-Newton iteration has not converged;" << " numGiteration = " << numGiteration << endl; cout << " ix, dga2, dga2old = " << ix << ", " << dga2 << ", " << dga2old << endl; continue; } //cout << "Schmelling : Gauss-Newton iteration has converged;" << endl; //cout << "==================================================" << endl; //cout << " numGiteration = " << numGiteration << endl; //cout << " ix, dga2 = " << ix << ", " << dga2 << endl; // calculated quantities which will be useful for determining // the best weight (Chisq, SpurAR, ...) //full = kTRUE; full = 1; SchmellCore(full, xiter, gamma, dgamma, dga2); // calculate difference between ideal and unfolded distribution Double_t D2bar = 0.0; for (UInt_t i = 0; iFindBin( xbin ); hBchisq->SetBinContent(bin,chisq(ix)); hBSpAR->SetBinContent(bin,SpAR(ix)); hBSpSig->SetBinContent(bin,SpSig(ix)/fSpurVacov); hBSecDeriv->SetBinContent(bin,SecDer(ix)); hBZerDeriv->SetBinContent(bin,ZerDer(ix)); hBEntropy->SetBinContent(bin,Entrop(ix)); hBDAR2->SetBinContent(bin,DAR2(ix)); hBD2bar->SetBinContent(bin,Dsqbar(ix)); if (ix > 0) { Double_t DSpAR = SpAR(ix) - SpAR(ix-1); hBDSpAR->SetBinContent(bin,DSpAR); Double_t diff = SpSig(ix) - SpSig(ix-1); Double_t DSpSig = diff; hBDSpSig->SetBinContent(bin, DSpSig/fSpurVacov); Double_t DEntrop = Entrop(ix) - Entrop(ix-1); hBDEntropy->SetBinContent(bin,DEntrop); Double_t DSecDer = SecDer(ix) - SecDer(ix-1); hBDSecDeriv->SetBinContent(bin,DSecDer); Double_t DZerDer = ZerDer(ix) - ZerDer(ix-1); hBDZerDeriv->SetBinContent(bin,DZerDer); } } // Select best weight SelectBestWeight(); if (ixbest < 0.0) { cout << "Schmelling : no solution found; " << endl; return kFALSE; } // do the unfolding using the best weight //full = kTRUE; xiter = pow(10.0,log10(xmin)+ixbest*dlogx) * fWinitial; //---------- start Gauss-Newton iteration ---------------------- numGiteration = 0; gamma.Zero(); dga2 = 1.e20; while (1) { full = 1; SchmellCore(full, xiter, gamma, dgamma, dga2); gamma += dgamma; //cout << "Schmelling : sum(dgamma^2) = " << dga2 << endl; numGiteration += 1; if (numGiteration > MaxGiteration) break; if (dga2 < EpsG) break; } //---------- end Gauss-Newton iteration ------------------------ //----------------------------------------------------------------------- // termination stage // ================= cout << "Schmelling : best solution found; " << endl; cout << "==================================" << endl; cout << " xiter, ixbest, numGiteration, Chisq = " << xiter << ", " << ixbest << ", " << numGiteration << ", " << Chisq << endl; //------------------------------------ //.............................................. // put unfolded distribution into fResult // fResult(i,0) value in bin i // fResult(i,1) error of value in bin i fNdf = SpurAR; fChisq = Chisq; for (UInt_t i=0; i0 ? TMath::Prob(fChisq, iNdf) : 0; fResult.ResizeTo(fNb, 5); for (UInt_t i=0; idmax) dmax = d(j,0); Double_t psum = 0.0; for (UInt_t j=0; j 0.0) Entropy += p(j,0) * log( p(j,0) ); //-------------------------------------------------------- } // ----------------------------------------------------------------------- // // Smooth migration matrix // by fitting a function to the migration matrix // Bool_t SmoothMigrationMatrix(TH2D &hmigorig) { // copy histograms into matrices; the matrices will be used in fcnSmooth // ------------------------ cout << "MUnfold::SmoothMigrationMatrix : fNa, fNb = " << fNa << ", " << fNb << endl; cout << "MUnfold::SmoothMigrationMatrix: fMigOrig = " << endl; cout << "========================================" << endl; for (UInt_t i=0; i = a0 + a1*log10(Etrue) + a2*SQR(log10(Etrue)) // + log10(Etrue) // b0RMS, b1RMS, b2RMS // RMS(log10(Eest)) = b0 + b1*log10(Etrue) + b2*SQR(log10(Etrue)) // UInt_t npar = 6; if (npar > 20) { cout << "MUnfold::SmoothMigrationMatrix : too many parameters, npar = " << npar << endl; return kFALSE; } //.............................................. // Find reasonable starting values for a0, a1 and b0, b1 Double_t xbar = 0.0; Double_t xxbar = 0.0; Double_t ybarm = 0.0; Double_t xybarm = 0.0; Double_t ybarR = 0.0; Double_t xybarR = 0.0; Double_t Sum = 0.0; for (UInt_t j=0; j 0.0) { meany = meany / sum; RMSy = RMSy / sum - meany*meany; RMSy = sqrt(RMSy); Sum += sum; xbar += x * sum; xxbar += x*x * sum; ybarm += meany * sum; xybarm += x*meany * sum; ybarR += RMSy * sum; xybarR += x*RMSy * sum; } } if (Sum > 0.0) { xbar /= Sum; xxbar /= Sum; ybarm /= Sum; xybarm /= Sum; ybarR /= Sum; xybarR /= Sum; } Double_t a1start = (xybarm - xbar*ybarm) / (xxbar - xbar*xbar); Double_t a0start = ybarm - a1start*xbar; a1start = a1start - 1.0; Double_t b1start = (xybarR - xbar*ybarR) / (xxbar - xbar*xbar); Double_t b0start = ybarR - b1start*xbar; cout << "MUnfold::SmoothMigrationMatrix : " << endl; cout << "============================" << endl; cout << "a0start, a1start = " << a0start << ", " << a1start << endl; cout << "b0start, b1start = " << b0start << ", " << b1start << endl; //.............................................. // Set starting values and step sizes for parameters char name[20][100]; Double_t vinit[20]; Double_t step[20]; Double_t limlo[20]; Double_t limup[20]; Int_t fix[20]; sprintf(&name[0][0], "a0mean"); vinit[0] = a0start; //vinit[0] = 1.0; step[0] = 0.1; limlo[0] = 0.0; limup[0] = 0.0; fix[0] = 0; sprintf(&name[1][0], "a1mean"); vinit[1] = a1start; //vinit[1] = 0.0; step[1] = 0.1; limlo[1] = 0.0; limup[1] = 0.0; fix[1] = 0; sprintf(&name[2][0], "a2mean"); vinit[2] = 0.0; step[2] = 0.1; limlo[2] = 0.0; limup[2] = 0.0; fix[2] = 0; sprintf(&name[3][0], "b0RMS"); vinit[3] = b0start; //vinit[3] = 0.8; step[3] = 0.1; limlo[3] = 1.e-20; limup[3] = 10.0; fix[3] = 0; sprintf(&name[4][0], "b1RMS"); vinit[4] = b1start; //vinit[4] = 0.0; step[4] = 0.1; limlo[4] = 0.0; limup[4] = 0.0; fix[4] = 0; sprintf(&name[5][0], "b2RMS"); vinit[5] = 0.0; step[5] = 0.1; limlo[5] = 0.0; limup[5] = 0.0; fix[5] = 0; if ( CallMinuit(fcnSmooth, npar, name, vinit, step, limlo, limup, fix) ) { // ------------------------ // fMigrat is the migration matrix to be used in the unfolding; // fMigrat, as set by the constructor, is overwritten // by the smoothed migration matrix for (UInt_t i=0; i 20) { cout << "MUnfold::Tikhonov2 : too many parameters, npar = " << npar << ", fNb = " << fNb << endl; return kFALSE; } // copy ideal distribution for (UInt_t i=1; i<=fNb; i++) { fhb0->SetBinContent(i, hb0.GetBinContent(i)); fhb0->SetBinError (i, hb0.GetBinError(i)); } //--- start w loop ----------------------------------- Int_t ix; Double_t xiter; for (ix=0; ixFindBin( log10(xiter) ); hBchisq->SetBinContent(bin,chisq(ix)); //hBSpAR->SetBinContent(bin,SpAR(ix)); hBSpAR->SetBinContent(bin,0.0); hBSpSig->SetBinContent(bin,SpSig(ix)/fSpurVacov); hBSecDeriv->SetBinContent(bin,SecDer(ix)); hBZerDeriv->SetBinContent(bin,ZerDer(ix)); hBEntropy->SetBinContent(bin,Entrop(ix)); //hBDAR2->SetBinContent(bin,DAR2(ix)); hBDAR2->SetBinContent(bin,0.0); hBD2bar->SetBinContent(bin,Dsqbar(ix)); if (ix > 0) { //Double_t DSpAR = SpAR(ix) - SpAR(ix-1); //hBDSpAR->SetBinContent(bin,DSpAR); Double_t diff = SpSig(ix) - SpSig(ix-1); Double_t DSpSig = diff; hBDSpSig->SetBinContent(bin, DSpSig/fSpurVacov); Double_t DEntrop = Entrop(ix) - Entrop(ix-1); hBDEntropy->SetBinContent(bin,DEntrop); Double_t DSecDer = SecDer(ix) - SecDer(ix-1); hBDSecDeriv->SetBinContent(bin,DSecDer); Double_t DZerDer = ZerDer(ix) - ZerDer(ix-1); hBDZerDeriv->SetBinContent(bin,DZerDer); } } } //--- end w loop ----------------------------------- // Select best weight SelectBestWeight(); cout << " Tikhonov2 : after SelectBestWeight" << endl; if (ixbest < 0.0) { cout << "Tikhonov2 : no result found; " << endl; return kFALSE; } cout << "Tikhonov2 : best result found; " << endl; cout << "===============================" << endl; cout << " ixbest = " << ixbest << endl; // do a final unfolding using the best weight fW = pow(10.0,log10(xmin)+ixbest*dlogx); //.............................................. // Set starting values and step sizes for parameters char name[20][100]; Double_t vinit[20]; Double_t step[20]; Double_t limlo[20]; Double_t limup[20]; Int_t fix[20]; for (UInt_t i=0; iSetBinContent(i, hb0.GetBinContent(i)); fhb0->SetBinError (i, hb0.GetBinError(i)); } TMatrixD bold(fNb, 1); bold.Zero(); //---------------------------------------------------------- Double_t db2 = 1.e20; TMatrixD aminusaest(fNa, 1); //------- scan number of iterations ----------------- Int_t ix; for (ix=0; ixFindBin( log10(xiter) ); hBchisq->SetBinContent(bin,chisq(ix)); hBSpAR->SetBinContent(bin,SpAR(ix)); hBSpSig->SetBinContent(bin,SpSig(ix)/fSpurVacov); hBSecDeriv->SetBinContent(bin,SecDer(ix)); hBZerDeriv->SetBinContent(bin,ZerDer(ix)); hBEntropy->SetBinContent(bin,Entrop(ix)); hBDAR2->SetBinContent(bin,DAR2(ix)); hBD2bar->SetBinContent(bin,Dsqbar(ix)); if (ix > 0) { Double_t DSpAR = SpAR(ix) - SpAR(ix-1); hBDSpAR->SetBinContent(bin,DSpAR); Double_t diff = SpSig(ix) - SpSig(ix-1); Double_t DSpSig = diff; hBDSpSig->SetBinContent(bin, DSpSig/fSpurVacov); Double_t DEntrop = Entrop(ix) - Entrop(ix-1); hBDEntropy->SetBinContent(bin,DEntrop); Double_t DSecDer = SecDer(ix) - SecDer(ix-1); hBDSecDeriv->SetBinContent(bin,DSecDer); Double_t DZerDer = ZerDer(ix) - ZerDer(ix-1); hBDZerDeriv->SetBinContent(bin,DZerDer); } } //------- end of scan of number of iterations ----------------- // Select best weight SelectBestWeight(); if (ixbest < 0.0) { cout << "Bertero : weight iteration has NOT converged; " << endl; return kFALSE; } cout << "Bertero : weight iteration has converged; " << endl; cout << " ixbest = " << ixbest << endl; // do a final unfolding using the best weight // calculate solution for the iteration number xiter Double_t xiter = pow(10.0,log10(xmin)+ixbest*dlogx); BertCore(xiter); cout << "Bertero : Values for best weight " << endl; cout << "=================================" << endl; cout << "fW, ixbest, Chisq, SpurAR, SpurSigma, SecDeriv, ZerDeriv, Entrop, DiffAR2, D2bar = " << endl; cout << " " << fW << ", " << ixbest << ", " << Chisq << ", " << SpurAR << ", " << SpurSigma << ", " << SecDeriv << ", " << ZerDeriv << ", " << Entropy << ", " << DiffAR2 << ", " << Dsqbar(ixbest) << endl; //---------------- fNdf = SpurAR; fChisq = Chisq; for (UInt_t i=0; i0 ? TMath::Prob(fChisq, iNdf) : 0; fResult.ResizeTo(fNb, 5); for (UInt_t i=0; i 0.0) p *= 1.0/sumb; Entropy = 0; for (UInt_t j=0; j 0.0) Entropy += p(j,0) * log( p(j,0) ); //----------------------------------------------------- TMatrixD Gb(fMigrat, TMatrixD::kMult, fVb); aminusaest = fVa; aminusaest -= Gb; TMatrixD v1(1,fNa); for (UInt_t i=0; i fVapoints ? chisqmin : fVapoints/2.0; // second loop over all weights : // consider only weights for which chisq(ix) < chisq0 ixbest = -1; for (ix=0; ix 0 && chisq(ix-1) != 0.0) { Double_t diff = SpSig(ix) - SpSig(ix-1); if (diff > DiffSpSigmax) { DiffSpSigmax = diff; ixDiffSpSigmax = ix; } } //---------------------------------- // search weight where Chisq is close // to the number of significant measurements Double_t DiffSigpoints = fabs(Chisq-fVapoints); if (DiffSigpoints < DiffSigpointsmin) { DiffSigpointsmin = DiffSigpoints; ixDiffSigpointsmin = ix; } //---------------------------------- // search weight where Chisq is close // to the rank of matrix G Double_t DiffRankG = fabs(Chisq-RankG); if (DiffRankG < DiffRankGmin) { DiffRankGmin = DiffRankG; ixDiffRankGmin = ix; } //---------------------------------- // search weight where SpurSigma is close to 1.0 Double_t DiffSpSig1 = fabs(SpurSigma/fSpurVacov-1.0); if (DiffSpSig1 < DiffSpSig1min) { DiffSpSig1min = DiffSpSig1; ixDiffSpSig1min = ix; } //---------------------------------- // search weight where D2bar is minimal if (D2bar < D2barmin) { D2barmin = D2bar; ixD2barmin = ix; } //---------------------------------- } } // choose solution where increase of SpurSigma is biggest //if ( DiffSpSigmax > 0.0) // ixbest = ixDiffSpSigmax; //else // ixbest = ixDiffSigpointsmin; // choose Least Squares Solution //$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ // ixbest = ixmax; //$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ // choose weight where chi2 is close to the number of significant // measurements // ixbest = ixDiffSigpointsmin; // choose weight where chi2 is close to the rank of matrix G // ixbest = ixDiffRankGmin; // choose weight where chi2 is close to the rank of matrix G //$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ ixbest = ixDiffSpSig1min; //$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ cout << "SelectBestWeight : ixDiffSpSigmax, DiffSpSigmax = " << ixDiffSpSigmax << ", " << DiffSpSigmax << endl; cout << "================== ixDiffSigpointsmin, DiffSigpointsmin = " << ixDiffSigpointsmin << ", " << DiffSigpointsmin << endl; cout << " ixDiffRankGmin, DiffRankGmin = " << ixDiffRankGmin << ", " << DiffRankGmin << endl; cout << " ixDiffSpSig1min, DiffSpSig1min = " << ixDiffSpSig1min << ", " << DiffSpSig1min << endl; cout << " ixD2barmin, D2barmin = " << ixD2barmin << ", " << D2barmin << endl; cout << " ixmax = " << ixmax << endl; cout << " ixbest = " << ixbest << endl; return kTRUE; } // ----------------------------------------------------------------------- // // Draw the plots // Bool_t DrawPlots() { // in the plots, mark the weight which has been selected Double_t xbin = log10(xmin)+ixbest*dlogx; TMarker *m = new TMarker(); m->SetMarkerSize(1); m->SetMarkerStyle(20); //------------------------------------- // draw the iteration plots TCanvas *c = new TCanvas("iter", "Plots versus weight", 900, 600); c->Divide(3,2); c->cd(1); hBchisq->Draw(); gPad->SetLogy(); hBchisq->SetXTitle("log10(iteration number)"); hBchisq->SetYTitle("chisq"); m->DrawMarker(xbin, log10(chisq(ixbest))); c->cd(2); hBD2bar->Draw(); gPad->SetLogy(); hBD2bar->SetXTitle("log10(iteration number)"); hBD2bar->SetYTitle("(b_unfolded-b_ideal)**2"); m->DrawMarker(xbin, log10(Dsqbar(ixbest))); /* c->cd(3); hBDAR2->Draw(); gPad->SetLogy(); strgx = "log10(iteration number)"; strgy = "norm(AR-AR+)"; hBDAR2->SetXTitle(strgx); hBDAR2->SetYTitle(strgy); m->DrawMarker(xbin, log10(DAR2(ixbest))); */ c->cd(3); hBSecDeriv->Draw(); hBSecDeriv->SetXTitle("log10(iteration number)"); hBSecDeriv->SetYTitle("Second Derivative squared"); m->DrawMarker(xbin, SecDer(ixbest)); /* c->cd(8); hBDSecDeriv->Draw(); strgx = "log10(iteration number)"; strgy = "Delta(Second Derivative squared)"; hBDSecDeriv->SetXTitle(strgx); hBDSecDeriv->SetYTitle(strgy); */ /* c->cd(4); hBZerDeriv->Draw(); strgx = "log10(iteration number)"; strgy = "Zero Derivative squared"; hBZerDeriv->SetXTitle(strgx); hBZerDeriv->SetYTitle(strgy); m->DrawMarker(xbin, ZerDer(ixbest)); */ /* c->cd(5); hBDZerDeriv->Draw(); strgx = "log10(iteration number)"; strgy = "Delta(Zero Derivative squared)"; hBDZerDeriv->SetXTitle(strgx); hBDZerDeriv->SetYTitle(strgy); */ c->cd(4); hBSpAR->Draw(); hBSpAR->SetXTitle("log10(iteration number)"); hBSpAR->SetYTitle("SpurAR"); m->DrawMarker(xbin, SpAR(ixbest)); /* c->cd(11); hBDSpAR->Draw(); strgx = "log10(iteration number)"; strgy = "Delta(SpurAR)"; hBDSpAR->SetXTitle(strgx); hBDSpAR->SetYTitle(strgy); */ c->cd(5); hBSpSig->Draw(); hBSpSig->SetXTitle("log10(iteration number)"); hBSpSig->SetYTitle("SpurSig/SpurC"); m->DrawMarker(xbin, SpSig(ixbest)/fSpurVacov); /* c->cd(14); hBDSpSig->Draw(); strgx = "log10(iteration number)"; strgy = "Delta(SpurSig/SpurC)"; hBDSpSig->SetXTitle(strgx); hBDSpSig->SetYTitle(strgy); */ c->cd(6); hBEntropy->Draw(); hBEntropy->SetXTitle("log10(iteration number)"); hBEntropy->SetYTitle("Entropy"); m->DrawMarker(xbin, Entrop(ixbest)); /* c->cd(17); hBDEntropy->Draw(); strgx = "log10(iteration number)"; strgy = "Delta(Entropy)"; hBDEntropy->SetXTitle(strgx); hBDEntropy->SetYTitle(strgy); */ //------------------------------------- for (UInt_t i=0; iSetBinContent(i+1, fVa(i, 0) ); fha->SetBinError (i+1, sqrt(fVacov(i, i))); for (UInt_t j=0; jSetBinContent(i+1, j+1, fMigOrig(i, j) ); fhmig->SetBinError (i+1, j+1, sqrt(fMigOrigerr2(i, j)) ); shmig->SetBinContent(i+1, j+1, fMigrat(i, j) ); shmig->SetBinError (i+1, j+1, sqrt(fMigraterr2(i, j)) ); shmigChi2->SetBinContent(i+1, j+1, fMigChi2(i, j) ); } } PrintTH2Content(*shmig); PrintTH2Content(*shmigChi2); //------------------------------------- CopyCol(*hprior, fVEps); CopyCol(*hb, fVb); for (UInt_t i=0; iSetBinError(i+1, sqrt(fVbcov(i, i))); PrintTH1Content(*hb); PrintTH1Error(*hb); //.............................................. for (UInt_t i=0; iSetBinContent(i+1, EigenValue(i)); //.............................................. // draw the plots TCanvas *cc = new TCanvas("input", "Unfolding input", 900, 600); cc->Divide(3, 2); // distribution to be unfolded cc->cd(1); fha->Draw(); gPad->SetLogy(); fha->SetXTitle("log10(E-est/GeV)"); fha->SetYTitle("Counts"); // superimpose unfolded distribution // hb->Draw("*HSAME"); // prior distribution cc->cd(2); hprior->Draw(); gPad->SetLogy(); hprior->SetXTitle("log10(E-true/GeV)"); hprior->SetYTitle("Counts"); // migration matrix cc->cd(3); fhmig->Draw("box"); fhmig->SetXTitle("log10(E-est/GeV)"); fhmig->SetYTitle("log10(E-true/GeV)"); // smoothed migration matrix cc->cd(4); shmig->Draw("box"); shmig->SetXTitle("log10(E-est/GeV)"); shmig->SetYTitle("log10(E-true/GeV)"); // chi2 contributions for smoothing cc->cd(5); shmigChi2->Draw("box"); shmigChi2->SetXTitle("log10(E-est/GeV)"); shmigChi2->SetYTitle("log10(E-true/GeV)"); // Eigenvalues of matrix M*M(transposed) cc->cd(6); hEigen->Draw(); hEigen->SetXTitle("l"); hEigen->SetYTitle("Eigen values Lambda_l of M*M(transposed)"); //.............................................. // draw the results TCanvas *cr = new TCanvas("results", "Unfolding results", 600, 600); cr->Divide(2, 2); // unfolded distribution cr->cd(1); hb->Draw(); gPad->SetLogy(); hb->SetXTitle("log10(E-true/GeV)"); hb->SetYTitle("Counts"); // covariance matrix of unfolded distribution cr->cd(2); TH1 *hbcov=DrawMatrixClone(fVbcov, "lego"); hbcov->SetBins(fNb, hb->GetBinLowEdge(1), hb->GetBinLowEdge(fNb+1), fNb, hb->GetBinLowEdge(1), hb->GetBinLowEdge(fNb+1)); hbcov->SetName("hbcov"); hbcov->SetTitle("Error matrix of distribution hb"); hbcov->SetXTitle("log10(E-true/GeV)"); hbcov->SetYTitle("log10(E-true/GeV)"); // chi2 contributions cr->cd(3); TH1 *hchi2=DrawMatrixColClone(fChi2); hchi2->SetBins(fNa, fha->GetBinLowEdge(1), fha->GetBinLowEdge(fNa+1)); hchi2->SetName("Chi2"); hchi2->SetTitle("chi2 contributions"); hchi2->SetXTitle("log10(E-est/GeV)"); hchi2->SetYTitle("Chisquared"); // ideal distribution cr->cd(4); fhb0->Draw(); gPad->SetLogy(); fhb0->SetXTitle("log10(E-true/GeV)"); fhb0->SetYTitle("Counts"); // superimpose unfolded distribution hb->Draw("*Hsame"); return kTRUE; } // ----------------------------------------------------------------------- // // Interface to MINUIT // // Bool_t CallMinuit( void (*fcnx)(Int_t &, Double_t *, Double_t &, Double_t *, Int_t), UInt_t npar, char name[20][100], Double_t vinit[20], Double_t step[20], Double_t limlo[20], Double_t limup[20], Int_t fix[20]) { // // Be carefull: This is not thread safe // UInt_t maxpar = 100; if (npar > maxpar) { cout << "MUnfold::CallMinuit : too many parameters, npar = " << fNb << ", maxpar = " << maxpar << endl; return kFALSE; } //.............................................. // Set the maximum number of parameters TMinuit minuit(maxpar); //.............................................. // Set the print level // -1 no output except SHOW comands // 0 minimum output // 1 normal output (default) // 2 additional ouput giving intermediate results // 3 maximum output, showing progress of minimizations // Int_t printLevel = -1; minuit.SetPrintLevel(printLevel); //.............................................. // Printout for warnings // SET WAR print warnings // SET NOW suppress warnings Int_t errWarn; Double_t tmpwar = 0; minuit.mnexcm("SET NOW", &tmpwar, 0, errWarn); //.............................................. // Set the address of the minimization function minuit.SetFCN(fcnx); //.............................................. // Set starting values and step sizes for parameters for (UInt_t i=0; i 0) { Int_t parNo = i; minuit.FixParameter(parNo); } } //.............................................. // Set maximum number of iterations (default = 500) Int_t maxiter = 100000; minuit.SetMaxIterations(maxiter); //.............................................. // minimization by the method of Migrad // Int_t errMigrad; // Double_t tmp = 0; // minuit.mnexcm("MIGRAD", &tmp, 0, errMigrad); //.............................................. // same minimization as by Migrad // but switches to the SIMPLEX method if MIGRAD fails to converge Int_t errMinimize; Double_t tmp = 0; minuit.mnexcm("MINIMIZE", &tmp, 0, errMinimize); //.............................................. // check quality of minimization // istat = 0 covariance matrix not calculated // 1 diagonal approximation only (not accurate) // 2 full matrix, but forced positive-definite // 3 full accurate covariance matrix // (indication of normal convergence) Double_t fmin, fedm, errdef; Int_t npari, nparx, istat; minuit.mnstat(fmin, fedm, errdef, npari, nparx, istat); if (errMinimize || istat < 3) { cout << "MUnfold::CallMinuit : Minimization failed" << endl; cout << " fmin = " << fmin << ", fedm = " << fedm << ", errdef = " << errdef << ", istat = " << istat << endl; return kFALSE; } //.............................................. // Minos error analysis // minuit.mnmnos(); //.............................................. // Print current status of minimization // if nkode = 0 only function value // 1 parameter values, errors, limits // 2 values, errors, step sizes, internal values // 3 values, errors, step sizes, 1st derivatives // 4 values, paraboloc errors, MINOS errors //Int_t nkode = 4; //minuit.mnprin(nkode, fmin); //.............................................. // call fcn with IFLAG = 3 (final calculation : calculate p(chi2)) // iflag = 1 initial calculations only // 2 calculate 1st derivatives and function // 3 calculate function only // 4 calculate function + final calculations const char *command = "CALL"; Double_t iflag = 3; Int_t errfcn3; minuit.mnexcm(command, &iflag, 1, errfcn3); return kTRUE; } // ----------------------------------------------------------------------- // // Return the unfolded distribution // TMatrixD &GetVb() { return fVb; } // ----------------------------------------------------------------------- // // Return the covariance matrix of the unfolded distribution // TMatrixD &GetVbcov() { return fVbcov; } // ----------------------------------------------------------------------- // // Return the unfolded distribution + various errors // TMatrixD &GetResult() { return fResult; } // ----------------------------------------------------------------------- // // Return the chisquared contributions // TMatrixD &GetChi2() { return fChi2; } // ----------------------------------------------------------------------- // // Return the total chisquared // Double_t &GetChisq() { return fChisq; } // ----------------------------------------------------------------------- // // Return the number of degrees of freedom // Double_t &GetNdf() { return fNdf; } // ----------------------------------------------------------------------- // // Return the chisquared probability // Double_t &GetProb() { return fProb; } // ----------------------------------------------------------------------- // // Return the smoothed migration matrix // TMatrixD &GetMigSmoo() { return fMigSmoo; } // ----------------------------------------------------------------------- // // Return the error2 of the smoothed migration matrix // TMatrixD &GetMigSmooerr2() { return fMigSmooerr2; } // ----------------------------------------------------------------------- // // Return the chi2 contributions for the smoothing // TMatrixD &GetMigChi2() { return fMigChi2; } }; // end of definition of class MUnfold /////////////////////////////////////////////////// // ----------------------------------------------------------------------- // // fcnSmooth (used by SmoothMigrationMatrix) // // is called by MINUIT // for given values of the parameters it calculates the function // to be minimized // void fcnSmooth(Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag) { MUnfold &gUnfold = *(MUnfold*)gMinuit->GetObjectFit(); Double_t a0 = par[0]; Double_t a1 = par[1]; Double_t a2 = par[2]; Double_t b0 = par[3]; Double_t b1 = par[4]; Double_t b2 = par[5]; // loop over bins of log10(E-true) Double_t chi2 = 0.0; Int_t npoints = 0; Double_t func[20]; for (UInt_t j=0; jGetObjectFit(); // (npar+1) is the number of bins of the unfolded distribuition (fNb) // npar is the number of free parameters (fNb-1) UInt_t npar1 = npar + 1; UInt_t fNa = gUnfold.fNa; UInt_t fNb = gUnfold.fNb; if (npar1 != fNb) { cout << "fcnTikhonov2 : inconsistency in number of parameters; npar, fNb = "; cout << npar << ", " << fNb << endl; //return; } npar1 = fNb; TMatrixD p(npar1, 1); TMatrixD &fVb = gUnfold.fVb; // p is the normalized unfolded distribution // sum(p(i,0)) from i=0 to npar is equal to 1 Double_t sum = 0.0; for (Int_t i=0; imnerrs(i, eplus, eminus, eparab, gcc); gUnfold.fVb(i, 0) = val * norm; gUnfold.fResult(i, 0) = val * norm; gUnfold.fResult(i, 1) = eparab * norm; gUnfold.fResult(i, 2) = eplus * norm; gUnfold.fResult(i, 3) = eminus * norm; gUnfold.fResult(i, 4) = gcc; sum += val; } gUnfold.fVb(gUnfold.fNb-1, 0) = (1.0-sum) * norm; gUnfold.fResult(gUnfold.fNb-1, 0) = (1.0-sum) * norm; gUnfold.fResult(gUnfold.fNb-1, 1) = sqrt(gUnfold.fVbcov(gUnfold.fNb-1,gUnfold.fNb-1)); gUnfold.fResult(gUnfold.fNb-1, 2) = 0; gUnfold.fResult(gUnfold.fNb-1, 3) = 0; gUnfold.fResult(gUnfold.fNb-1, 4) = 1; //.............................................. //----------------------------------------------------- // calculate 0th derivative squared gUnfold.ZerDeriv = 0; for (UInt_t j=0; j 0) gUnfold.Entropy += p(j,0) * log( p(j,0) ); //----------------------------------------------------- // calculate SpurSigma gUnfold.SpurSigma = 0.0; for (UInt_t m=0; m0 ? TMath::Prob(gUnfold.fChisq, iNdf) : 0; } } // ====================================================== // // SteerUnfold // void SteerUnfold(TH1D &ha, TH2D &hacov, TH2D &hmig, TH2D &hmigor, TH1D &hb0, TH1D *hpr=NULL) { // ha is the distribution to be unfolded // hacov is the covariance matrix of the distribution ha // hmig is the migration matrix; // it is used in the unfolding unless it is overwritten // by SmoothMigrationMatrix by the smoothed migration matrix // hmigor is the migration matrix to be smoothed; // the smoothed migration matrix will be used in the unfolding // hpr the prior distribution // it is only used if SetPriorInput(*hpr) is called //.............................................. // create an MUnfold object; // fill histograms into vectors and matrices MUnfold unfold(ha, hacov, hmig); //.............................................. // smooth the migration matrix; // the smoothed migration matrix will be used in the unfolding // hmig is the original (unsmoothed) migration matrix unfold.SmoothMigrationMatrix(hmigor); //.............................................. // define prior distribution (has always to be defined) // the alternatives are : // 1 SetPriorConstant() : isotropic distribution // 2 SetPriorPower(gamma) : dN/dE = E^{-gamma} // 3 SetPriorInput(*hpr): the distribution *hpr is used // 4 SetPriorRebin(*ha) : use rebinned histogram ha UInt_t flagprior = 4; cout << "SteerUnfold : flagprior = " << flagprior << endl; cout << "=========================="<< endl; Bool_t errorprior=kTRUE; switch (flagprior) { case 1: unfold.SetPriorConstant(); break; case 2: errorprior = unfold.SetPriorPower(1.5); break; case 3: if (!hpr) { cout << "Error: No hpr!" << endl; return; } errorprior = unfold.SetPriorInput(*hpr); break; case 4: errorprior = unfold.SetPriorRebin(ha); break; } if (!errorprior) { cout << "MUnfold::SetPrior... : failed. flagprior = " ; cout << flagprior << endl; return; } //.............................................. // calculate the matrix G = M * M(transposed) // M being the migration matrix unfold.CalculateG(); //.............................................. // call steering routine for the actual unfolding; // the alternatives are : // 1 Schmelling : minimize the function Z by Gauss-Newton iteration; // the parameters to be fitted are gamma(i) = lambda(i)/w; // 2 Tikhonov2 : regularization term is sum of (2nd deriv.)**2 ; // minimization by using MINUIT; // the parameters to be fitted are // the bin contents of the unfolded distribution // 3 Bertero: minimization by iteration // UInt_t flagunfold = 1; cout << "SteerUnfold : flagunfold = " << flagunfold << endl; cout << "===========================" << endl; switch (flagunfold) { case 1: // Schmelling cout << "" << endl; cout << "Unfolding algorithm : Schmelling" << endl; cout << "================================" << endl; if (!unfold.Schmelling(hb0)) cout << "MUnfold::Schmelling : failed." << endl; break; case 2: // Tikhonov2 cout << "" << endl; cout << "Unfolding algorithm : Tikhonov" << endl; cout << "================================" << endl; if (!unfold.Tikhonov2(hb0)) cout << "MUnfold::Tikhonov2 : failed." << endl; break; case 3: // Bertero cout << "" << endl; cout << "Unfolding algorithm : Bertero" << endl; cout << "================================" << endl; if (!unfold.Bertero(hb0)) cout << "MUnfold::Bertero : failed." << endl; break; } //.............................................. // Print fResult unfold.PrintResults(); //.............................................. // Draw the plots unfold.DrawPlots(); //.............................................. // get unfolded distribution //TMatrixD &Vb = unfold.GetVb(); //TMatrixD &Vbcov = unfold.GetVbcov(); } //__________________________________________________________________________ //////////////////////////////////////////////////////////////////////////// // // // Main program // // defines the ideal distribution (hb0) // // defines the migration matrix (hMigrat) // // defines the distribution to be unfolded (hVa) // // // // calls member functions of the class MUnfold // // to do the unfolding // // // //////////////////////////////////////////////////////////////////////////// void unfold() { // ----------------------------------------- // migration matrix : // x corresponds to measured quantity // y corresponds to true quantity //const Int_t na = 13; const Int_t na = 18; const Axis_t alow = 0.25; const Axis_t aup = 3.50; //const Int_t nb = 11; const Int_t nb = 22; const Axis_t blow = 0.50; const Axis_t bup = 3.25; TH2D hmig("Migrat", "Migration Matrix", na, alow, aup, nb, blow, bup); hmig.Sumw2(); // parametrize migration matrix as // = a0 + a1*log10(Etrue) + a2*log10(Etrue)**2 // + log10(Etrue) // RMS( log10(Eest) ) = b0 + b1*log10(Etrue) + b2*log10(Etrue)**2 Double_t a0 = 0.0; Double_t a1 = 0.0; Double_t a2 = 0.0; Double_t b0 = 0.26; Double_t b1 =-0.054; Double_t b2 = 0.0; TF1 f2("f2", "gaus(0)", alow, aup); f2.SetParName(0, "ampl"); f2.SetParName(1, "mean"); f2.SetParName(2, "sigma"); // loop over log10(Etrue) bins TAxis &yaxis = *hmig.GetYaxis(); for (Int_t j=1; j<=nb; j++) { Double_t yvalue = yaxis.GetBinCenter(j); const Double_t mean = a0 + a1*yvalue + a2*yvalue*yvalue + yvalue; const Double_t sigma = b0 + b1*yvalue + b2*yvalue*yvalue; const Double_t ampl = 1./ ( sigma*TMath::Sqrt(2.0*TMath::Pi())); // gaus(0) is a substitute for [0]*exp( -0.5*( (x-[1])/[2] )**2 ) f2.SetParameter(0, ampl); f2.SetParameter(1, mean); f2.SetParameter(2, sigma); // fill temporary 1-dim histogram with the function // fill the histogram using // - either FillRandom // - or using Freq TH1D htemp("temp", "temp", na, alow, aup); htemp.Sumw2(); for (Int_t k=0; k<1000000; k++) htemp.Fill(f2.GetRandom()); // copy it into the migration matrix Double_t sum = 0; for (Int_t i=1; i<=na; i++) { const Stat_t content = htemp.GetBinContent(i); hmig.SetBinContent(i, j, content); sum += content; } // normalize migration matrix if (sum==0) continue; for (Int_t i=1; i<=na; i++) { const Stat_t content = hmig.GetBinContent(i,j); hmig.SetBinContent(i,j, content/sum); hmig.SetBinError (i,j,sqrt(content)/sum); } } PrintTH2Content(hmig); PrintTH2Error(hmig); // ----------------------------------------- // ideal distribution TH1D hb0("hb0", "Ideal distribution", nb, blow, bup); hb0.Sumw2(); // fill histogram with random numbers according to // an exponential function dN/dE = E^{-gamma} // or with y = log10(E), E = 10^y : // dN/dy = ln10 * 10^{y*(1-gamma)} TF1 f1("f1", "pow(10.0, x*(1.0-[0]))", blow, bup); f1.SetParName(0,"gamma"); f1.SetParameter(0, 2.7); // ntimes is the number of entries for (Int_t k=0; k<10000; k++) hb0.Fill(f1.GetRandom()); // introduce energy threshold at 50 GeV const Double_t lgEth = 1.70; const Double_t dlgEth = 0.09; for (Int_t j=1; j<=nb; j++) { const Double_t lgE = hb0.GetBinCenter(j); const Double_t c = hb0.GetBinContent(j); const Double_t dc = hb0.GetBinError(j); const Double_t f = 1.0 / (1.0 + exp( -(lgE-lgEth)/dlgEth )); hb0.SetBinContent(j, f* c); hb0.SetBinError (j, f*dc); } PrintTH1Content(hb0); // ----------------------------------------- // here the prior distribution can be defined for the call // to SetPriorInput(*hpr) TH1D hpr("hpr", "Prior distribution" , nb, blow, bup); for (Int_t j=1; j<=nb; j++) hpr.SetBinContent(j, 1.0/nb); PrintTH1Content(hpr); // ----------------------------------------- // generate distribution to be unfolded (ha) // by smearing the ideal distribution (hb0) // TH1D ha("ha", "Distribution to be unfolded", na, alow, aup); ha.Sumw2(); for (Int_t i=1; i<=na; i++) { Double_t cont = 0; for (Int_t j=1; j<=nb; j++) cont += hmig.GetBinContent(i, j) * hb0.GetBinContent(j); ha.SetBinContent(i, cont); ha.SetBinError(i, sqrt(cont)); } PrintTH1Content(ha); PrintTH1Error(ha); // ----------------------------------------- // covariance matrix of the distribution ha TH2D hacov("hacov", "Error matrix of distribution ha", na, alow, aup, na, alow, aup); for (Int_t i=1; i<=na; i++) { const Double_t content = ha.GetBinContent(i); hacov.SetBinContent(i, i, content<3 ? 3.0 : content); } PrintTH2Content(hacov); SteerUnfold(ha, hacov, hmig, hmig, hb0, &hpr); } //========================================================================//