CMS 3D CMS Logo

Classes | Public Member Functions | Static Public Member Functions | Private Attributes | Static Private Attributes

ResidualFitter Class Reference

#include <ResidualFitter.h>

List of all members.

Classes

class  Point

Public Member Functions

void AddMeas (double x, double y, double z)
void AddTrue (double x, double y, double z)
void cd ()
void CreateCanvas ()
TH1D * CurBin ()
void ExecuteEvent (Int_t event, Int_t px, Int_t py, TObject *sel)
void Fill ()
void Fit (int x, int y, const char *opt="")
void FitSlicesZ (TF1 *f1=0)
TCanvas * GetCanvas ()
TF1 * GetFunction ()
 ResidualFitter (const char *name, const char *title, int nbinsx, double xlow, double xup, int nbinsy, double ylow, double yup, int nbinsz, double zlow, double zup)
void SetAutoRange (int n)
void SetFitFunction (TF1 *func)
void SetFitOptions (const char *opt)
void SetMinN (int n)
void ShowFit (TH2 *h2d, double x, double y)
 ~ResidualFitter ()

Static Public Member Functions

static void SetCanvas (int x, int y)

Private Attributes

int autoRangeN_
TCanvas * canvas_
TCanvas * canvasFit_
TH2D * chi2_
TH1D * curBin_
TF1 * fitFunction_
std::string fitOptions_
TH2D * mean_
std::vector
< ResidualFitter::Point
meas_
int minN_
int oldBinx_
int oldBiny_
TH2D * sigma_
std::vector
< ResidualFitter::Point
true_

Static Private Attributes

static int xCanvas_ = 600
static int yCanvas_ = 600

Detailed Description

Definition at line 11 of file ResidualFitter.h.


Constructor & Destructor Documentation

ResidualFitter::ResidualFitter ( const char *  name,
const char *  title,
int  nbinsx,
double  xlow,
double  xup,
int  nbinsy,
double  ylow,
double  yup,
int  nbinsz,
double  zlow,
double  zup 
)

Definition at line 21 of file ResidualFitter.cc.

References chi2_, gather_cfg::cout, CreateCanvas(), mean_, mergeVDriftHistosByStation::name, and sigma_.

  : TH3D( name, title, 
          nbinsx, xlow, xup, 
          nbinsy, ylow, yup, 
          nbinsz, zlow, zup ), 
    fitFunction_( new TF1("gaus", "gaus") ), 
    canvasFit_(0),
    curBin_(0), 
    autoRangeN_(0),
    minN_(5) {

  cout<<"creating residual fitter with name "<<name<<endl;
  
  string meanname = name; meanname += "_mean";
  mean_ = new TH2D(meanname.c_str(), meanname.c_str(),
                   nbinsx, xlow, xup, 
                   nbinsy, ylow, yup);
  mean_->SetStats(0);
  
  string sigmaname = name; sigmaname += "_sigma";
  sigma_ = new TH2D(sigmaname.c_str(), sigmaname.c_str(),
                    nbinsx, xlow, xup, 
                    nbinsy, ylow, yup);

  sigma_->SetStats(0);

  string chi2name = name; chi2name += "_chi2";
  chi2_ = new TH2D(chi2name.c_str(), chi2name.c_str(),
                   nbinsx, xlow, xup, 
                   nbinsy, ylow, yup);
  chi2_->SetStats(0);

  //   string nseenname = name; nseenname += "_nseen";
  //   nseen_ = new TH2D(nseenname.c_str(), nseenname.c_str(),
  //               nbinsx, xlow, xup, 
  //               nbinsy, ylow, yup);
  //   nseen_->SetStats(0);

  gDirectory->ls();

  CreateCanvas();
}
ResidualFitter::~ResidualFitter ( )

Definition at line 68 of file ResidualFitter.cc.

References canvas_, chi2_, curBin_, fitFunction_, mean_, and sigma_.

                                {
  delete fitFunction_;
  if(canvas_) delete canvas_;
  if(curBin_) delete curBin_;

  delete mean_;
  delete sigma_;
  delete chi2_;
  //   delete nseen_;
}

Member Function Documentation

void ResidualFitter::AddMeas ( double  x,
double  y,
double  z 
) [inline]

Definition at line 29 of file ResidualFitter.h.

References meas_.

    {meas_.push_back( Point(x,y,z) ); }
void ResidualFitter::AddTrue ( double  x,
double  y,
double  z 
) [inline]

Definition at line 27 of file ResidualFitter.h.

References true_.

    {true_.push_back( Point(x,y,z) ); }
void ResidualFitter::cd ( void  ) [inline]

Definition at line 49 of file ResidualFitter.h.

References canvas_.

{canvas_->cd(); }
void ResidualFitter::CreateCanvas ( )

Definition at line 79 of file ResidualFitter.cc.

References canvas_, xCanvas_, and yCanvas_.

Referenced by ResidualFitter().

                                  {
  string cname = "ResidualFitterCanvas_"; cname += GetName();
  canvas_ = new TCanvas(cname.c_str(), cname.c_str(),xCanvas_, yCanvas_);
  
  canvas_ ->Connect("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)", 
                    "ResidualFitter",
                    this, "ExecuteEvent(Int_t,Int_t,Int_t,TObject*)");

  
}
TH1D* ResidualFitter::CurBin ( ) [inline]

Definition at line 53 of file ResidualFitter.h.

References curBin_.

{return curBin_;}
void ResidualFitter::ExecuteEvent ( Int_t  event,
Int_t  px,
Int_t  py,
TObject *  sel 
)

Definition at line 124 of file ResidualFitter.cc.

References canvas_, EgammaValidation_Wenu_cff::sel, ShowFit(), x, and detailsBasic3DVector::y.

                                                                               {

  if( event != kButton1Down ) return;

  TH2* histo2d = dynamic_cast<TH2*>(sel);
  if(!histo2d) return;

  float x = canvas_->AbsPixeltoX(px);
  float y = canvas_->AbsPixeltoY(py);
  x = canvas_->PadtoX(x);
  y = canvas_->PadtoY(y);

  ShowFit(histo2d, x, y);
}
void ResidualFitter::Fill ( )

Definition at line 91 of file ResidualFitter.cc.

References meas_, HcalObjRepresent::Reset(), and true_.

                          {

  Reset();
  
  for(unsigned it = 0; it<true_.size(); it++) {
    for(unsigned im = 0; im<meas_.size(); im++) {
      TH3D::Fill(true_[it].x_, 
                 true_[it].y_, 
                 meas_[it].z_ - true_[it].z_ );
    }
  }
}
void ResidualFitter::Fit ( int  x,
int  y,
const char *  opt = "" 
)

Definition at line 166 of file ResidualFitter.cc.

References autoRangeN_, chi2_, gather_cfg::cout, curBin_, fitFunction_, fitOptions_, mean_, minN_, and sigma_.

Referenced by FitSlicesZ(), and ShowFit().

                                                            {
  TH1::AddDirectory(0);
  
  if(curBin_) delete curBin_;
  curBin_ = TH3::ProjectionZ("", binx, binx, biny, biny);  

  if(curBin_->GetEntries() < minN_ ) {
    TH1::AddDirectory(1);
    return;
  }

  string sopt = fitOptions_; sopt += opt;

  if( autoRangeN_ ) {
    double maxpos = curBin_->GetBinCenter( curBin_->GetMaximumBin() ); 
    
    double minrange = maxpos-curBin_->GetRMS()* autoRangeN_;
    double maxrange = maxpos+curBin_->GetRMS()* autoRangeN_;

    fitFunction_->SetRange( minrange, maxrange );
    cout<<"range : "<<minrange<<" "<<maxrange<<endl;
  }

  curBin_->Fit(fitFunction_, sopt.c_str() );


  double chi2overndf=0;
  if(fitFunction_->GetNDF() ) {
    chi2overndf = fitFunction_->GetChisquare()/ fitFunction_->GetNDF();
    mean_->SetBinContent(binx,biny, fitFunction_->GetParameter(1) );
    mean_->SetBinError(binx,biny, fitFunction_->GetParError(1) );
    sigma_->SetBinContent(binx,biny, fitFunction_->GetParameter(2) );
    sigma_->SetBinError(binx,biny, fitFunction_->GetParError(2) );
 
    chi2_->SetBinContent(binx,biny,chi2overndf);
  }
  //   nseen_->SetBinContent(binx, biny, 
  //                    fitFunction_->Integral( fitFunction_->GetXmin(), 
  //                                            fitFunction_->GetXmax())
  //                    /curBin_->GetBinWidth(1) );

  TH1::AddDirectory(1);
}
void ResidualFitter::FitSlicesZ ( TF1 *  f1 = 0)

Definition at line 104 of file ResidualFitter.cc.

References chi2_, gather_cfg::cout, Fit(), mean_, SetFitFunction(), and sigma_.

                                       {

  mean_->Reset();
  sigma_->Reset();
  chi2_->Reset();
  //   nseen_->Reset();

  cout<<"ResidualFitter::FitSlicesZ"<<endl;
  if(f1) SetFitFunction(f1);
  
  for(int ix=1; ix<=GetNbinsX(); ix++) {
    for(int iy=1; iy<=GetNbinsY(); iy++) {
      Fit(ix, iy, "Q0");
    }
  }

  //  TH3D::FitSlicesZ(f1);
}
TCanvas* ResidualFitter::GetCanvas ( ) [inline]

Definition at line 55 of file ResidualFitter.h.

References canvas_.

{return canvas_;}
TF1* ResidualFitter::GetFunction ( ) [inline]

Definition at line 51 of file ResidualFitter.h.

References fitFunction_.

{return fitFunction_;}
void ResidualFitter::SetAutoRange ( int  n) [inline]

Definition at line 42 of file ResidualFitter.h.

References autoRangeN_, fitOptions_, and n.

                           {
    autoRangeN_= n;
    fitOptions_ += "R";
  }
static void ResidualFitter::SetCanvas ( int  x,
int  y 
) [inline, static]

Definition at line 24 of file ResidualFitter.h.

References x, xCanvas_, detailsBasic3DVector::y, and yCanvas_.

void ResidualFitter::SetFitFunction ( TF1 *  func) [inline]

Definition at line 33 of file ResidualFitter.h.

References fitFunction_.

Referenced by FitSlicesZ().

{ fitFunction_ = func; }  
void ResidualFitter::SetFitOptions ( const char *  opt) [inline]

Definition at line 34 of file ResidualFitter.h.

References fitOptions_.

{fitOptions_ = opt; }
void ResidualFitter::SetMinN ( int  n) [inline]

Definition at line 47 of file ResidualFitter.h.

References minN_, and n.

{ minN_ = n;}
void ResidualFitter::ShowFit ( TH2 *  h2d,
double  x,
double  y 
)

Definition at line 140 of file ResidualFitter.cc.

References canvas_, canvasFit_, Fit(), oldBinx_, and oldBiny_.

Referenced by ExecuteEvent().

                                                             {
  
  if(!canvasFit_) {
    string cname = "ResidualFitterCanvasFit_"; cname += GetName();
    canvasFit_ = new TCanvas(cname.c_str(), cname.c_str(),300,300);
  }
  canvasFit_ ->cd();
  
  int binx = histo2d->GetXaxis()->FindBin(x); 
  int biny = histo2d->GetYaxis()->FindBin(y); 

  if(binx == oldBinx_ && biny == oldBiny_ ) return;
  oldBinx_ = binx;
  oldBiny_ = biny;

  Fit(binx, biny);
 
  canvasFit_->Modified();
  canvasFit_->Update();

  canvas_->Modified();
  canvas_->Update();
  canvas_->cd();
}

Member Data Documentation

Definition at line 94 of file ResidualFitter.h.

Referenced by Fit(), and SetAutoRange().

TCanvas* ResidualFitter::canvas_ [private]

Definition at line 77 of file ResidualFitter.h.

Referenced by cd(), CreateCanvas(), ExecuteEvent(), GetCanvas(), ShowFit(), and ~ResidualFitter().

TCanvas* ResidualFitter::canvasFit_ [private]

Definition at line 83 of file ResidualFitter.h.

Referenced by ShowFit().

TH2D* ResidualFitter::chi2_ [private]

Definition at line 88 of file ResidualFitter.h.

Referenced by Fit(), FitSlicesZ(), ResidualFitter(), and ~ResidualFitter().

TH1D* ResidualFitter::curBin_ [private]

Definition at line 85 of file ResidualFitter.h.

Referenced by CurBin(), Fit(), and ~ResidualFitter().

Definition at line 74 of file ResidualFitter.h.

Referenced by Fit(), GetFunction(), SetFitFunction(), and ~ResidualFitter().

std::string ResidualFitter::fitOptions_ [private]

Definition at line 75 of file ResidualFitter.h.

Referenced by Fit(), SetAutoRange(), and SetFitOptions().

TH2D* ResidualFitter::mean_ [private]

Definition at line 86 of file ResidualFitter.h.

Referenced by Fit(), FitSlicesZ(), ResidualFitter(), and ~ResidualFitter().

std::vector< ResidualFitter::Point > ResidualFitter::meas_ [private]

Definition at line 72 of file ResidualFitter.h.

Referenced by AddMeas(), and Fill().

int ResidualFitter::minN_ [private]

if the number of entries in a bin is lower than this, fit is not performed. set to 5 by default.

Definition at line 98 of file ResidualFitter.h.

Referenced by Fit(), and SetMinN().

int ResidualFitter::oldBinx_ [private]

Definition at line 91 of file ResidualFitter.h.

Referenced by ShowFit().

int ResidualFitter::oldBiny_ [private]

Definition at line 92 of file ResidualFitter.h.

Referenced by ShowFit().

TH2D* ResidualFitter::sigma_ [private]

Definition at line 87 of file ResidualFitter.h.

Referenced by Fit(), FitSlicesZ(), ResidualFitter(), and ~ResidualFitter().

std::vector< ResidualFitter::Point > ResidualFitter::true_ [private]

Definition at line 71 of file ResidualFitter.h.

Referenced by AddTrue(), and Fill().

int ResidualFitter::xCanvas_ = 600 [static, private]

Definition at line 79 of file ResidualFitter.h.

Referenced by CreateCanvas(), and SetCanvas().

int ResidualFitter::yCanvas_ = 600 [static, private]

Definition at line 80 of file ResidualFitter.h.

Referenced by CreateCanvas(), and SetCanvas().