TChi2FitData.cxx

Go to the documentation of this file.
00001 // @(#)root/minuit2:$Id: TChi2FitData.cxx 20880 2007-11-19 11:23:41Z rdm $
00002 // Author: L. Moneta    10/2005  
00003 
00004 /**********************************************************************
00005  *                                                                    *
00006  * Copyright (c) 2005 ROOT Foundation,  CERN/PH-SFT                   *
00007  *                                                                    *
00008  **********************************************************************/
00009 
00010 #include <cassert>
00011 
00012 #include "RConfig.h"
00013 #include "TChi2FitData.h"
00014 
00015 #include "TVirtualFitter.h" 
00016 
00017 //#define DEBUG
00018 #ifdef DEBUG
00019 #include <iostream>
00020 #endif
00021 
00022 #include "TList.h"
00023 #include "TF1.h"
00024 #include "TH1.h"
00025 #include "TGraph.h"
00026 #include "TGraph2D.h"
00027 #include "TMultiGraph.h"
00028 
00029 
00030 
00031 
00032 
00033 TChi2FitData::TChi2FitData( const TVirtualFitter & fitter, bool skipEmptyBins) : 
00034 fSize(0), fSkipEmptyBins(skipEmptyBins), fIntegral(false)
00035 {
00036    // constructor - create fit data from fitter ROOT data object (histogram, graph, etc...)
00037    // need a pointer to TF1 (model function) for avoid points in regions rejected by TF1 
00038    // option skipEmptyBins is used to skip the bin with zero content or to use them in the fit 
00039    // (in that case an error of 1 is set)
00040    
00041    TF1 * func = dynamic_cast<TF1 *> ( fitter.GetUserFunc() );  
00042    assert( func != 0);
00043    
00044    TObject * obj = fitter.GetObjectFit(); 
00045    // downcast to see type of object
00046    TH1 * hfit = dynamic_cast<TH1*> ( obj );
00047    if (hfit) { 
00048       GetFitData(hfit, func, &fitter);    
00049       return; 
00050    } 
00051    // case of TGraph
00052    TGraph * graph = dynamic_cast<TGraph*> ( obj );
00053    if (graph) { 
00054       GetFitData(graph, func, &fitter);    
00055       return; 
00056    } 
00057    // case of TGraph2D
00058    TGraph2D * graph2D = dynamic_cast<TGraph2D*> ( obj );
00059    if (graph2D) { 
00060       GetFitData(graph2D, func, &fitter);    
00061       return; 
00062    } 
00063    // case of TMultiGraph
00064    TMultiGraph * multigraph = dynamic_cast<TMultiGraph*> ( obj );
00065    if (multigraph) { 
00066       GetFitData(graph2D, func, &fitter);    
00067       return; 
00068    } 
00069    // else 
00070 #ifdef DEBUG
00071    std::cout << "other fit type are not yet supported- assert" << std::endl;
00072 #endif
00073    assert(hfit != 0); 
00074    
00075 }
00076 
00077 
00078 void TChi2FitData::GetFitData(const TH1 * hfit, const TF1 * func, const TVirtualFitter * hFitter) {
00079    // get Histogram Data
00080    
00081    assert(hfit != 0); 
00082    assert(hFitter != 0);
00083    assert(func != 0);
00084    
00085    //std::cout << "creating Fit Data from histogram " << hfit->GetName() << std::endl; 
00086    
00087    //  use TVirtual fitter to get fit range (should have a FitRange class ) 
00088    
00089    // first and last bin
00090    int hxfirst = hFitter->GetXfirst(); 
00091    int hxlast  = hFitter->GetXlast(); 
00092    int hyfirst = hFitter->GetYfirst(); 
00093    int hylast  = hFitter->GetYlast(); 
00094    int hzfirst = hFitter->GetZfirst(); 
00095    int hzlast  = hFitter->GetZlast(); 
00096    TAxis *xaxis  = hfit->GetXaxis();
00097    TAxis *yaxis  = hfit->GetYaxis();
00098    TAxis *zaxis  = hfit->GetZaxis();
00099    
00100    // get fit option 
00101    Foption_t fitOption = hFitter->GetFitOption();
00102    if (fitOption.Integral) fIntegral=true;
00103    
00104    int n = (hxlast-hxfirst+1)*(hylast-hyfirst+1)*(hzlast-hzfirst+1); 
00105    
00106 #ifdef DEBUG
00107    std::cout << "TChi2FitData: ifirst = " << hxfirst << " ilast =  " << hxlast 
00108              << "total bins  " << hxlast-hxfirst+1  
00109              << "skip empty bins "  << fSkipEmptyBins << std::endl; 
00110 #endif
00111    
00112    fInvErrors.reserve(n);
00113    fValues.reserve(n);
00114    fCoordinates.reserve(n);
00115    
00116    int ndim =  hfit->GetDimension();
00117    assert( ndim > 0 );
00118    CoordData x = CoordData( hfit->GetDimension() );
00119    int binx = 0; 
00120    int biny = 0; 
00121    int binz = 0; 
00122    
00123    for ( binx = hxfirst; binx <= hxlast; ++binx) {
00124       if (fIntegral) {
00125          x[0] = xaxis->GetBinLowEdge(binx);       
00126       }
00127       else
00128          x[0] = xaxis->GetBinCenter(binx);
00129       
00130       if ( ndim > 1 ) { 
00131          for ( biny = hyfirst; biny <= hylast; ++biny) {
00132             if (fIntegral) 
00133                x[1] = yaxis->GetBinLowEdge(biny);
00134             else
00135                x[1] = yaxis->GetBinCenter(biny);
00136             
00137             if ( ndim >  2 ) { 
00138                for ( binz = hzfirst; binz <= hzlast; ++binz) {
00139                   if (fIntegral) 
00140                      x[2] = zaxis->GetBinLowEdge(binz);
00141                   else
00142                      x[2] = zaxis->GetBinCenter(binz);
00143                   if (!func->IsInside(&x.front()) ) continue;
00144                   double error =  hfit->GetBinError(binx, biny, binz); 
00145                   if (fitOption.W1) error = 1;
00146                   SetDataPoint( x,  hfit->GetBinContent(binx, biny, binz), error );
00147                }  // end loop on z bins
00148             }
00149             else if (ndim == 2) { 
00150                // for dim == 2
00151                if (!func->IsInside(&x.front()) ) continue;
00152                double error =  hfit->GetBinError(binx, biny); 
00153                if (fitOption.W1) error = 1;
00154                SetDataPoint( x,  hfit->GetBinContent(binx, biny), error );
00155             }   
00156             
00157          }  // end loop on y bins
00158          
00159       }
00160       else if (ndim == 1) { 
00161          // for 1D 
00162          if (!func->IsInside(&x.front()) ) continue;
00163          double error =  hfit->GetBinError(binx); 
00164          if (fitOption.W1) error = 1;
00165          SetDataPoint( x,  hfit->GetBinContent(binx), error );
00166       }
00167       
00168    }   // end 1D loop 
00169    
00170    // in case of integral store additional point with upper x values
00171    if (fIntegral) { 
00172       x[0] = xaxis->GetBinLowEdge(hxlast) +  xaxis->GetBinWidth(hxlast); 
00173       if (ndim > 1) { 
00174          x[1] = yaxis->GetBinLowEdge(hylast) +  yaxis->GetBinWidth(hylast); 
00175       }
00176       if (ndim > 2) { 
00177          x[2] = zaxis->GetBinLowEdge(hzlast) +  zaxis->GetBinWidth(hzlast); 
00178       }
00179       fCoordinates.push_back(x);
00180    }
00181    
00182 #ifdef DEBUG
00183    std::cout << "TChi2FitData: Hist FitData size is " << fCoordinates.size() << std::endl;
00184 #endif
00185    
00186 }
00187 
00188 
00189 void TChi2FitData::GetFitData(const TGraph * gr, const TF1 * func, const TVirtualFitter * hFitter ) {
00190    // get TGraph data (neglect error in x, that is used in the extended method)
00191    
00192    assert(gr != 0); 
00193    assert(hFitter != 0);
00194    assert(func != 0);
00195    
00196    // fit options
00197    Foption_t fitOption = hFitter->GetFitOption();
00198    
00199    int  nPoints = gr->GetN();
00200    double *gx = gr->GetX();
00201    double *gy = gr->GetY();
00202    
00203    CoordData x = CoordData( 1 );  // 1D graph
00204    
00205    for ( int i = 0; i < nPoints; ++i) { 
00206       
00207       x[0] = gx[i];
00208       // neglect error in x (it is a different chi2) 
00209       if (!func->IsInside(&x.front()) ) continue;
00210       double errorY = gr->GetErrorY(i); 
00211       // consider error = 0 as 1 
00212       if (errorY <= 0) errorY = 1;
00213       if (fitOption.W1) errorY = 1;
00214       SetDataPoint( x, gy[i], errorY );
00215       
00216    }
00217 }
00218 
00219 
00220 void TChi2FitData::GetFitData(const TGraph2D * gr, const TF1 * func, const TVirtualFitter * hFitter ) {
00221    // fetch graph 2D data for CHI2 fit. 
00222    // neglect errors in x and y (one use the ExtendedChi2 method)
00223    
00224    assert(gr != 0); 
00225    assert(hFitter != 0);
00226    assert(func != 0);
00227    
00228    // fit options
00229    Foption_t fitOption = hFitter->GetFitOption();
00230    
00231    int  nPoints = gr->GetN();
00232    double *gx = gr->GetX();
00233    double *gy = gr->GetY();
00234    double *gz = gr->GetZ();
00235    
00236    CoordData x = CoordData( 2 );  // 2D graph
00237    
00238    for ( int i = 0; i < nPoints; ++i) { 
00239       
00240       x[0] = gx[i];
00241       x[1] = gy[i];
00242       if (!func->IsInside(&x.front()) ) continue;
00243       // neglect error in x (it is a different chi2) 
00244       double error = gr->GetErrorZ(i); 
00245       // consider error = 0 as 1 
00246       if (error <= 0) error = 1;
00247       if (fitOption.W1) error = 1;
00248       SetDataPoint( x, gz[i], error );
00249       
00250    }
00251 }
00252 
00253 
00254 
00255 void TChi2FitData::GetFitData(const TMultiGraph * mg, const TF1 * func, const TVirtualFitter * hFitter ) {
00256    // data from a multigraph
00257    
00258    assert(mg != 0); 
00259    assert(hFitter != 0);
00260    assert(func != 0);
00261    
00262    // fit options
00263    Foption_t fitOption = hFitter->GetFitOption();
00264    
00265    TGraph *gr;
00266    TIter next(mg->GetListOfGraphs());   
00267    
00268    int  nPoints;
00269    double *gx;
00270    double *gy;
00271    
00272    CoordData x = CoordData( 1 );  // 1D graph
00273    
00274    while ((gr = (TGraph*) next())) {
00275       nPoints = gr->GetN();
00276       gx      = gr->GetX();
00277       gy      = gr->GetY();
00278       for ( int i = 0; i < nPoints; ++i) { 
00279          
00280          x[0] = gx[i];
00281          // neglect error in x (it is a different chi2) 
00282          if (!func->IsInside(&x.front()) ) continue;
00283          double errorY = gr->GetErrorY(i); 
00284          // consider error = 0 as 1 
00285          if (errorY <= 0) errorY = 1;
00286          if (fitOption.W1) errorY = 1;
00287          SetDataPoint( x, gy[i], errorY );
00288          
00289       }
00290    }
00291    
00292 }
00293 
00294 void TChi2FitData::SetDataPoint( const CoordData & x, double y, double error) { 
00295    // set internally the data in the internal vectors and count them
00296    // if error is <=0 (like for zero content bin) skip or set to 1, according to the set option
00297    if (error <= 0) {  
00298       if (SkipEmptyBins() ) 
00299          return;
00300       else   
00301          // set errors to 1 
00302          error = 1;
00303    }
00304    
00305    fCoordinates.push_back(x);
00306    fValues.push_back(y);
00307    fInvErrors.push_back(1./error);
00308    fSize++;
00309 }

Generated on Tue Jul 5 14:37:11 2011 for ROOT_528-00b_version by  doxygen 1.5.1