00001 #include "THnSparse.h"
00002 #include "TH2.h"
00003 #include "TF2.h"
00004 #include "TF3.h"
00005 #include "TCanvas.h"
00006 #include "TApplication.h"
00007 #include "TMath.h"
00008
00009 #include "Fit/SparseData.h"
00010 #include "HFitInterface.h"
00011 #include "Fit/Fitter.h"
00012 #include "Math/WrappedMultiTF1.h"
00013
00014 #include "TRandom2.h"
00015
00016 #include <iostream>
00017 #include <iterator>
00018 #include <algorithm>
00019
00020 #include <list>
00021 #include <vector>
00022
00023 #include <cmath>
00024 #include <cassert>
00025
00026 using namespace std;
00027
00028 const bool __DRAW__ = 1;
00029
00030 TRandom2 r(17);
00031 Int_t bsize[] = { 10, 10, 10 };
00032 Double_t xmin[] = { 0., 0., 0. };
00033 Double_t xmax[] = { 10., 10., 10. };
00034
00035 double gaus1D(double *x, double *p)
00036 {
00037 return p[0]*TMath::Gaus(x[0],p[1],p[2]);
00038 }
00039
00040 double gaus2D(double *x, double *p)
00041 {
00042 return p[0]*TMath::Gaus(x[0],p[1],p[2]) * TMath::Gaus(x[1],p[3],p[4]);
00043 }
00044
00045 double gaus3D(double *x, double *p)
00046 {
00047 return p[0] * TMath::Gaus(x[0],p[1],p[2])
00048 * TMath::Gaus(x[1],p[3],p[4])
00049 * TMath::Gaus(x[2],p[5],p[6]);
00050 }
00051
00052 double pol2D(double *x, double *p)
00053 {
00054 return p[0]*x[0]+ p[1]*x[0]*x[0] + p[2]*x[1] + p[3]*x[1]*x[1] + p[4];
00055 }
00056
00057 ostream& operator << (ostream& out, ROOT::Fit::BinData& bd)
00058 {
00059 const unsigned int ndim( bd.NDim() );
00060 const unsigned int npoints( bd.NPoints() );
00061 for ( unsigned int i = 0; i < npoints; ++i )
00062 {
00063 double value, error;
00064 const double *x = bd.GetPoint(i, value, error);
00065 for ( unsigned int j = 0; j < ndim; ++j )
00066 {
00067 out << " x[" << j << "]: " << x[j];
00068 }
00069 out << " value: " << value;
00070 out << " error: " << error;
00071 out << endl;
00072 }
00073 return out;
00074 }
00075
00076 int findBin(ROOT::Fit::BinData& bd, const double *x)
00077 {
00078 const unsigned int ndim = bd.NDim();
00079 const unsigned int npoints = bd.NPoints();
00080
00081 for ( unsigned int i = 0; i < npoints; ++i )
00082 {
00083 double value1, error1;
00084 const double *x1 = bd.GetPoint(i, value1, error1);
00085 bool thisIsIt = true;
00086 for ( unsigned int j = 0; j < ndim; ++j )
00087 {
00088 thisIsIt &= fabs(x1[j] - x[j]) < 1E-15;
00089 }
00090 if ( thisIsIt ) return i;
00091 }
00092
00093 cout << "NO ENCONTRADO!";
00094 copy(x, x+ndim, ostream_iterator<double>(cout, " " ));
00095 cout << endl;
00096
00097 return -1;
00098 }
00099
00100 bool operator ==(ROOT::Fit::BinData& bd1, ROOT::Fit::BinData& bd2)
00101 {
00102 const unsigned int ndim = bd1.NDim();
00103 const unsigned int npoints = bd1.NPoints();
00104
00105 bool equals = true;
00106
00107 cout << "Equals" << endl;
00108
00109 for ( unsigned int i = 0; i < npoints && equals; ++i )
00110 {
00111 double value1, error1;
00112 const double *x1 = bd1.GetPoint(i, value1, error1);
00113
00114 int bin = findBin(bd2, x1);
00115
00116 double value2 = 0, error2;
00117 const double *x2 = bd2.GetPoint(bin, value2, error2);
00118
00119 equals &= ( value1 == value2 );
00120 cout << " v: " << equals;
00121 equals &= ( error1 == error2 );
00122 cout << " e: " << equals;
00123 for ( unsigned int j = 0; j < ndim; ++j )
00124 {
00125 equals &= fabs(x1[j] - x2[j]) < 1E-15;
00126 cout << " x[" << j << "]: " << equals;
00127 }
00128
00129 cout << " bd1: ";
00130 std::copy(x1, &x1[ndim], ostream_iterator<double>(cout, " "));
00131 cout << " value:" << value1;
00132 cout << " error:" << error1;
00133
00134 cout << " bd2: ";
00135 std::copy(x2, &x2[ndim], ostream_iterator<double>(cout, " "));
00136 cout << " value:" << value2;
00137 cout << " error:" << error2;
00138
00139 cout << " equals: " << equals;
00140
00141 cout << endl;
00142 }
00143
00144 return equals;
00145 }
00146
00147 void fillSparse(THnSparse* s, TF1* f, int nEvents = 5)
00148 {
00149 const unsigned int ndim = s->GetNdimensions();
00150
00151 for ( Int_t e = 0; e < nEvents; ++e ) {
00152 Double_t points[ndim];
00153 for ( UInt_t i = 0; i < ndim; ++ i )
00154 points[i] = r.Uniform( xmin[0] * .9 , xmax[0] * 1.1 );
00155 double value = gRandom->Poisson( f->EvalPar(points));
00156 s->Fill(points, value );
00157 cout << value << " " << s->GetNbins() << endl;
00158 }
00159 }
00160
00161 void DoFit(THnSparse* s, TF1* f, ROOT::Fit::BinData& bd)
00162 {
00163
00164 cout << "Retrieving the Sparse Data Structure" << endl;
00165
00166 ROOT::Fit::SparseData d(s->GetNdimensions(), xmin, xmax);
00167 ROOT::Fit::FillData(d, s, 0);
00168 d.GetBinData(bd);
00169
00170
00171
00172 bool ret;
00173 ROOT::Fit::Fitter fitter;
00174 ROOT::Math::WrappedMultiTF1 wf2(*f);
00175 ROOT::Math::IParametricFunctionMultiDim & if2 = wf2;
00176 fitter.Config().SetMinimizer("Minuit2");
00177
00178 ROOT::Fit::DataOptions opt;
00179 opt.fUseEmpty = true;
00180 opt.fIntegral = true;
00181
00182
00183 cout << "\n ******* Likelihood with BinData and NoCeros *******" << endl;
00184 ROOT::Fit::BinData bdNoCeros;
00185 d.GetBinDataNoZeros(bdNoCeros);
00186 ret = fitter.LikelihoodFit(bdNoCeros, if2);
00187 fitter.Result().Print(std::cout);
00188 if (!ret)
00189 std::cout << "Fit Failed " << std::endl;
00190
00191
00192 cout << "\n ******* Likelihood with BinData with Ceros *******" << endl;
00193 ROOT::Fit::BinData bdWithCeros(opt);
00194 d.GetBinDataIntegral(bdWithCeros);
00195 ret = fitter.LikelihoodFit(bdWithCeros, if2);
00196 fitter.Result().Print(std::cout);
00197 if (!ret)
00198 std::cout << "Fit Failed " << std::endl;
00199
00200
00201 }
00202
00203 void fitSparse1D()
00204 {
00205 const unsigned int ndim = 1;
00206
00207 THnSparseD* s1 = new THnSparseD("mfND-s1", "s1-Title", ndim, bsize, xmin, xmax);
00208
00209 TF1* f = new TF1("func1D", gaus1D, xmin[0] - 2, xmax[0] + 2, 3);
00210 f->SetParameters(10., 5., 2.);
00211
00212 fillSparse(s1,f,5);
00213
00214 cout << "Retrieving the Sparse Data Structure" << endl;
00215
00216 ROOT::Fit::SparseData d(ndim, xmin, xmax);
00217 ROOT::Fit::FillData(d, s1, 0);
00218 }
00219
00220 void fitSparse2D()
00221 {
00222 const unsigned int ndim = 2;
00223
00224 THnSparseD* s1 = new THnSparseD("mfND-s1", "s1-Title", ndim, bsize, xmin, xmax);
00225
00226 TF2* f = new TF2("func2D", gaus2D, xmin[0],xmax[0], xmin[1], xmax[1], 5);
00227 f->SetParameters(40,5,2,5,1);
00228
00229 for (int ix=1; ix <= bsize[0]; ++ix) {
00230 for (int iy=1; iy <= bsize[1]; ++iy) {
00231 double x = s1->GetAxis(0)->GetBinCenter(ix);
00232 double y = s1->GetAxis(1)->GetBinCenter(iy);
00233
00234 double value = gRandom->Poisson( f->Eval(x,y) );
00235 if ( value )
00236 {
00237 double points[] = {x,y};
00238 s1->Fill( points, value );
00239 }
00240 }
00241 }
00242
00243 ROOT::Fit::BinData bd;
00244 DoFit(s1, f, bd);
00245
00246 TH2D* h2 = new TH2D("2D Blanked Hist Fit", "h1-title",
00247 bsize[0], xmin[0], xmax[0],
00248 bsize[1], xmin[1], xmax[1]);
00249 cout << "Filling second histogram" << endl;
00250 for ( unsigned int i = 0; i < bd.NPoints(); ++i )
00251 {
00252 const double* x;
00253 double value, error;
00254 x = bd.GetPoint(i, value, error);
00255 value = (value)?value:-10;
00256 h2->Fill(x[0], x[1], value);
00257 }
00258
00259 h2->Draw("colZ");
00260 }
00261
00262 void fitSparse3D()
00263 {
00264 const unsigned int ndim = 3;
00265
00266 THnSparseD* s1 = new THnSparseD("mfND-s1", "s1-Title", ndim, bsize, xmin, xmax);
00267
00268 TF3* f = new TF3("func3D", gaus3D,
00269 xmin[0],xmax[0],
00270 xmin[1],xmax[1],
00271 xmin[2],xmax[2],
00272 7);
00273 f->SetParameters(100,5,2,5,1,5,2);
00274
00275 for (int ix=1; ix <= bsize[0]; ++ix) {
00276 for (int iy=1; iy <= bsize[1]; ++iy) {
00277 for (int iz=1; iz <= bsize[2]; ++iz) {
00278 double x = s1->GetAxis(0)->GetBinCenter(ix);
00279 double y = s1->GetAxis(1)->GetBinCenter(iy);
00280 double z = s1->GetAxis(2)->GetBinCenter(iz);
00281 double value = gRandom->Poisson( f->Eval(x,y,z) );
00282
00283 if ( value )
00284 {
00285 double points[] = {x,y,z};
00286 s1->Fill( points, value );
00287 }
00288 }
00289 }
00290 }
00291
00292 ROOT::Fit::BinData bd;
00293 DoFit(s1, f, bd);
00294 }
00295
00296
00297 int main(int argc, char** argv)
00298 {
00299 TApplication* theApp = 0;
00300
00301 if ( __DRAW__ )
00302 theApp = new TApplication("App",&argc,argv);
00303
00304 fitSparse1D();
00305 fitSparse2D();
00306 fitSparse3D();
00307
00308 cout << "C'est fini!" << endl;
00309
00310 if ( __DRAW__ ) {
00311 theApp->Run();
00312 delete theApp;
00313 theApp = 0;
00314 }
00315
00316 return 0;
00317 }
00318