rarSPlot Class Reference

Create an sPlot. More...

#include <rarSPlot.hh>

List of all members.

Public Member Functions

 rarSPlot ()
 Trivial ctor.
 rarSPlot (const char *name, const char *title, RooAbsReal *obs, RooDataSet &data, RooFitResult *fitRes, const RooArgList &pdfs, const RooArgList &yields, const RooArgList &pdf0s=RooArgList(), const RooArgList &yield0s=RooArgList(), const RooArgSet &projDeps=RooArgSet(), const Bool_t verbose=kTRUE)
 Default ctor.
virtual ~rarSPlot ()
 Trivial dtor.
virtual void setVerbose (Bool_t verbose=kTRUE)
 Set verbose mode.
RooDataSet * fill (RooAbsReal &yield, Int_t nbins, Double_t min, Double_t max, Bool_t doErrors=kTRUE)
 Fill the sPlot histogram and dataset for a yield.
virtual TH1F * getSPlotHist ()
 Get sPlot histogram.
virtual RooDataSet * getSPlotData ()
 Get sPlot dataset.

Protected Member Functions

virtual void init ()
 Initial function called by ctor.
virtual void fillsPn (Int_t compIdx)
 Fill sPn array.

Protected Attributes

RooAbsReal * _obs
 Obs to do sPlot.
RooDataSet _data
 Input dataset.
Int_t _nEvts
 Number of events in dataset.
RooFitResult * _fitRes
 Fit results.
RooArgList _fitPars
 Fit params.
RooArgList _pdfs
 Component pdfs.
RooArgList _yields
 Component yields.
RooArgList _pdf0s
 Component pdfs w/ fixed yields.
RooArgList _yield0s
 Fixed component yields.
RooArgSet _projDeps
RooArgSet _normVars
 Obs to be ignored for normalization.
Int_t _nComps
 Obs to be used for normalization Number of components.
Int_t _nComp0s
 Number of fixed components.
TArrayI _idxMap
 Array for index mapping between _yields and _fitPars.
TMatrixD _covM
 Covariance matrix from Minuit.
TMatrixD _iV
 Inverse of covM or from direct calculation.
TMatrixD _V
 covM or from inverse of iV
TArrayD _dens
 Denominators.
TArrayD _sPns
 Array for sPn.
TArrayI _sPnb
 Array for sPn fill bit.
TH1F * _sHist
 SPlot histogram.
RooDataSet * _sData
 SPlot dataset.
Bool_t _verbose
 Boolean to control debug output.

Private Member Functions

 ClassDef (rarSPlot, 0)


Detailed Description

Create an sPlot.

This creates an sPlot. See BAD 509 for more information on sPlots. An sPlot gives us the distribution of some variable, x in our data sample for a given species, for instance the distribution of signal mES. The result is similar to a likelihood projection plot, but no cuts are made, so every event contributes to the distribution.

To use this class, you first must perform your fit twice: The first time perform your nominal fit. For the second fit, fix your parameters at the minimum, float only your yields, and remove any PDFs correlated with the variable of interest. (In particular, if you want to make a signal mES sPlot, fix everything but yields, and perform the fit without the mES PDFs.) Be sure to save the RooFitResult. Now you can use rarSPlot.

Definition at line 45 of file rarSPlot.hh.


Constructor & Destructor Documentation

rarSPlot::rarSPlot (  ) 

Trivial ctor.

Usually the objects should be created using other ctors.

Definition at line 34 of file rarSPlot.cc.

00034                    :
00035   TNamed(),
00036   _obs(0), _fitRes(0),
00037   _nComps(0), _nComp0s(0), _idxMap(0), _covM(0,0), _iV(0,0), _V(0,0),
00038   _sHist(0), _sData(0),
00039   _verbose(kFALSE)
00040 {
00041   // Default constructor
00042 }

rarSPlot::rarSPlot ( const char *  name,
const char *  title,
RooAbsReal *  obs,
RooDataSet &  data,
RooFitResult *  fitRes,
const RooArgList &  pdfs,
const RooArgList &  yields,
const RooArgList &  pdf0s = RooArgList(),
const RooArgList &  yield0s = RooArgList(),
const RooArgSet &  projDeps = RooArgSet(),
const Bool_t  verbose = kTRUE 
)

Default ctor.

Parameters:
name The name
title The title
obs The obs to get sPlot
data The dataset to get sPlot from
fitRes The fit result
pdfs The component pdfs
yields The component yields
pdf0s The component pdfs w/ fixed yields
yield0s The fixed component yields
projDeps Obs to be ignored for normalization
verbose Boolean to control debug output
Default ctor to set several common data members, and then call init.

Definition at line 59 of file rarSPlot.cc.

References init().

00063                                                                     :
00064   TNamed(name, title),
00065   _obs(obs), _data(data), _nEvts(_data.numEntries()),
00066   _fitRes(fitRes), _fitPars(fitRes->floatParsFinal()),
00067   _pdfs(pdfs), _yields(yields), _pdf0s(pdf0s), _yield0s(yield0s),
00068   _projDeps(projDeps), _normVars(*pdfs[0].getDependents(&data)),
00069   _nComps(0), _nComp0s(0), _idxMap(0), _covM(0,0), _iV(0,0), _V(0,0),
00070   _sHist(0), _sData(0),
00071   _verbose(verbose)
00072 {
00073   init();
00074 }

rarSPlot::~rarSPlot (  )  [virtual]

Trivial dtor.

Definition at line 77 of file rarSPlot.cc.

00078 {
00079 }


Member Function Documentation

rarSPlot::ClassDef ( rarSPlot  ,
 
) [private]

RooDataSet * rarSPlot::fill ( RooAbsReal &  yield,
Int_t  nbins,
Double_t  min,
Double_t  max,
Bool_t  doErrors = kTRUE 
)

Fill the sPlot histogram and dataset for a yield.

Parameters:
yield The yield to get sPlot for
nbins Number of bins for histogram
min Histogram min
max Histogram max
doErrors Use Sumw2() for histogram
Returns:
The sPlot dataset created
It fills the sPlot histogram and dataset

Definition at line 345 of file rarSPlot.cc.

References _data, _dens, _fitPars, _nComp0s, _nComps, _nEvts, _normVars, _obs, _pdfs, _sData, _sHist, _sPns, _V, _verbose, and fillsPn().

00347 {
00348   if (_nComps<=0) {
00349     _sHist=0;
00350     _sData=0;
00351     return 0;
00352   }
00353   
00354   // create sPlot histogram
00355   TString histName=Form("hist_%s_%s",GetName(), yield.GetName());
00356   _sHist=new TH1F(histName, histName, nbins, min, max);
00357   // Size of bins...
00358   Double_t binSize = (max-min)/nbins;
00359   if (_verbose) {
00360     cout <<" SPlot histogram bins: " << min << " to " << max
00361          <<" with " << nbins << " bins."
00362          <<" binSize = " << binSize << endl;
00363   }
00364   
00365   Int_t istar=_fitPars.index(_fitPars.find(yield.GetName()));
00366   // get sPns
00367   fillsPn(istar);
00368   if (_verbose) {
00369     Double_t sum = 0;
00370     for (Int_t j=0; j<_nComps; j++) {
00371       sum += _V(j,istar);
00372     }
00373     cout <<" Sum of star column in V: " << sum << endl
00374          <<" fitted yield "<<yield.GetName()<<" = "<<yield.getVal()<<endl;
00375   }
00376   // add a weight column to the dataset
00377   RooArgSet obs = *_data.get();
00378   RooRealVar evtWgt("evtWgt", "Weight for sPlot entry", 1.);
00379   obs.add(evtWgt); // add weight to list of observables
00380   _sData = new RooDataSet("sPlotData", "sPlotData", obs, "evtWgt");
00381   
00382   Double_t sum = 0;
00383   Double_t sumtmp = 0;
00384   // This forces the error in a bin to be calculated
00385   // as the sqrt of the sum of the squares of weights, as they should be.
00386   if (doErrors) _sHist->Sumw2();
00387   
00388   Double_t Nstar= yield.getVal();
00389   Double_t Vstar=0;
00390   for (Int_t j=0; j<_nComps; j++) {
00391     Vstar+=_V(istar, j);
00392   }
00393   Double_t Vsum=0;
00394   for(Int_t i=0; i<_nComps; i++) {
00395     for (Int_t j=0; j<_nComps; j++) {
00396       Vsum+=_V(i,j);
00397     }
00398   }
00399   Double_t Nratio=(Nstar-Vstar)/(_nEvts-Vsum);
00400   Double_t Nratio2=(Nstar-Vstar)*Nratio;
00401   if (_nComp0s>0) cout<<" Nratio = "<<Nratio<<" \tNratio2 = "<<Nratio2<<endl;
00402   cout<<" loop over dataset to fill sPlot ";
00403   //_data.write("/tmp/sd.txt");
00404   for (Int_t ievt=0; ievt<_nEvts; ievt++) {
00405     if (_verbose&&ievt%1000==0) {
00406       cout << ".";
00407       cout.flush();
00408     }
00409     // Read this event and find the value of x for this event.
00410     const RooArgSet *row=_data.get(ievt);
00411     Double_t obsVal = ((RooAbsReal*)row->find(_obs->GetName()))->getVal();
00412     Double_t p=1/Nstar*_sPns[istar*_nEvts+ievt];
00413     if (_nComp0s>0) {
00414       p=1/Nstar*(_sPns[istar*_nEvts+ievt]+Nratio*_sPns[_nComps*_nEvts+ievt]);
00415     }
00416     sumtmp += ((RooAbsPdf*)_pdfs.at(istar))->getVal(&_normVars)/_dens[ievt];
00417     
00418     //if (obsVal > min && obsVal < max)
00419     sum += p*Nstar;
00420     
00421     // Add weighted event to the sPlot histogram
00422     _sHist->Fill(obsVal, p*Nstar);
00423     // ... and to the dataSet
00424     _sData->add(*row, p*Nstar);
00425     
00426   }// end event loop
00427   cout<<endl;
00428   
00429   _sHist->SetEntries(sum*Double_t(binSize));
00430   
00431   if (_verbose)
00432     cout << " Entries should be: " << sum*Double_t(binSize)
00433          << " (" << sum << " events)" << endl
00434          << " Sum of likelihood ratios for nstar: " << sumtmp << endl;
00435   
00436   return _sData;
00437 }

void rarSPlot::fillsPn ( Int_t  compIdx  )  [protected, virtual]

Fill sPn array.

Parameters:
compIdx Component index

Definition at line 287 of file rarSPlot.cc.

References _data, _dens, _fitPars, _idxMap, _nComps, _nEvts, _normVars, _pdfs, _sPnb, _sPns, _V, and _verbose.

Referenced by fill(), and init().

00288 {
00289   if ((compIdx>_nComps)||(compIdx<0)) {
00290     cout<<" Invalid comp index "<<compIdx<<endl;
00291     return;
00292   }
00293   // check if done already
00294   if (_sPnb[compIdx]) return;
00295   if (compIdx==_nComps) {
00296     // fill sPn's
00297     for (Int_t i=0; i<_nComps; i++) {
00298       fillsPn(i);
00299     }
00300     // fill sP0
00301     cout<<" fill sP0 ";
00302     for (Int_t i=0; i<_nEvts; i++) {
00303       if (_verbose&&i%1000==0) {
00304         cout << ".";
00305         cout.flush();
00306       }
00307       _sPns[compIdx*_nEvts+i]=1;
00308       for (Int_t j=0; j<_nComps; j++) {
00309         _sPns[compIdx*_nEvts+i]-=_sPns[j*_nEvts+i];
00310       }
00311     }
00312   } else {
00313     // fill sPn
00314     cout<<" fill sPn for fitPar #"<<compIdx
00315         <<" "<<_fitPars[compIdx].GetName()<<" ";
00316     for (Int_t i=0; i<_nEvts; i++) {
00317       if (_verbose&&i%1000==0) {
00318         cout << ".";
00319         cout.flush();
00320       }
00321       _data.get(i);
00322       Double_t numerator=0.;
00323       for (Int_t j=0; j<_nComps; j++) {
00324         numerator+=_V(compIdx,j)*
00325           ((RooAbsPdf*)_pdfs.at(_idxMap[j]))->getVal(&_normVars);
00326       }
00327       _sPns[compIdx*_nEvts+i]=numerator/_dens[i];
00328     }
00329   }
00330   cout<<endl;
00331   
00332   _sPnb[compIdx]=1;
00333   return;
00334 }

virtual RooDataSet* rarSPlot::getSPlotData (  )  [inline, virtual]

Get sPlot dataset.

Returns:
The sPlot histogram

Definition at line 68 of file rarSPlot.hh.

References _sData.

00068 {return _sData;}

virtual TH1F* rarSPlot::getSPlotHist (  )  [inline, virtual]

Get sPlot histogram.

Returns:
The sPlot histogram

Definition at line 65 of file rarSPlot.hh.

References _sHist.

00065 {return _sHist;}

void rarSPlot::init (  )  [protected, virtual]

Initial function called by ctor.

It will initialize sPlot Matrics, etc. based on input

Definition at line 84 of file rarSPlot.cc.

References _covM, _data, _dens, _fitPars, _fitRes, _idxMap, _iV, _nComp0s, _nComps, _nEvts, _normVars, _pdf0s, _pdfs, _projDeps, _sPnb, _sPns, _V, _verbose, _yield0s, _yields, and fillsPn().

Referenced by rarSPlot().

00085 {
00086   
00087   if (_verbose) {
00088     cout<<" rarSPlot:  yields:"<<endl;
00089     _yields.Print();
00090     cout<<" pdfs:"<<endl;
00091     _pdfs.Print();
00092     cout<<" data:"<<endl;
00093     _data.Print();
00094     //_data.Print("v");
00095     cout<<endl;
00096   }
00097   
00098   _nComps=_yields.getSize();
00099   if (_nComps<=0) {
00100     cout<<" No comp found"<<endl;
00101     _nComps=0;
00102     return;
00103   }
00104   if (_pdfs.getSize()!=_nComps) {
00105     cout<<" The size of yields and pdfs do not match"<<endl;
00106     _nComps=0;
00107     return;
00108   }
00109   if (_fitPars.getSize()!=_nComps) {
00110     cout<<" The size of yields and floating params do not match"<<endl;
00111     _yields.Print();
00112     _fitPars.Print();
00113     _nComps=0;
00114     return;
00115   }
00116   _nComp0s=_yield0s.getSize();
00117   if (_pdf0s.getSize()!=_nComp0s) {
00118     cout<<" The size of yield0s and pdf0s do not match"<<endl;
00119     _nComps=0;
00120     return;
00121   }
00122   if ((_nComp0s!=0)&&_verbose) {
00123     cout<<" yield0s:"<<endl;
00124     _yield0s.Print();
00125     cout<<" pdf0s:"<<endl;
00126     _pdf0s.Print();
00127     cout<<endl;
00128   }
00129 
00130   // Set the right size
00131   _idxMap.Set(_nComps);
00132   _covM.ResizeTo(_nComps, _nComps);
00133   _iV.ResizeTo(_nComps, _nComps);
00134   _V.ResizeTo(_nComps, _nComps);
00135   
00136   // Remember, the order of finalPars is not necessarily the same
00137   // as in #_yields, and we'd better map between them
00138   for (Int_t i=0; i<_nComps; i++) {
00139     RooAbsArg *theYield(0);
00140     if(theYield=_yields.find(_fitPars[i].GetName()))
00141       _idxMap[i]=_yields.index(theYield);
00142     else {
00143       cout<<" Can not find "<<_fitPars[i].GetName()
00144           <<" in yields:"<<endl;
00145       _yields.Print();
00146       _nComps=0;
00147       return;
00148     }
00149     if (_verbose)
00150       cout<<" fitPar #"<<i<<" --> yield #"<<_idxMap[i]
00151           <<" "<<_fitPars[i].GetName()<<endl;
00152   }
00153   if (_verbose) cout<<endl;
00154   
00155   // normalization ignored obs
00156   // The list of variables to normalize over when calculating PDF values.
00157   _normVars.remove(_projDeps, kTRUE, kTRUE);
00158   if (_verbose) {
00159     cout<<" Normalization variables" << endl;
00160     _normVars.Print();
00161     cout<<" Norm ignored variables"<<endl;
00162     _projDeps.Print();
00163     cout<<endl;
00164   }
00165 
00166   // create full snapshot
00167   RooArgList *pdfs=(RooArgList*)_pdfs.snapshot();
00168   RooArgList *yields=(RooArgList*)_yields.snapshot();
00169   RooArgList *pdf0s=(RooArgList*)_pdf0s.snapshot();
00170   RooArgList *yield0s=(RooArgList*)_yield0s.snapshot();
00171   _pdfs.removeAll();
00172   _pdfs.add(*pdfs);
00173   _yields.removeAll();
00174   _yields.add(*yields);
00175   _pdf0s.removeAll();
00176   _pdf0s.add(*pdf0s);
00177   _yield0s.removeAll();
00178   _yield0s.add(*yield0s);
00179   // attach dataset
00180   for (Int_t i = 0; i<_nComps; i++) {
00181     RooAbsPdf *pdf = (RooAbsPdf*)_pdfs.at(_idxMap[i]);
00182     pdf->attachDataSet(_data);
00183   }
00184   for (Int_t i = 0; i<_nComp0s; i++) {
00185     RooAbsPdf *pdf = (RooAbsPdf*)_pdf0s.at(i);
00186     pdf->attachDataSet(_data);
00187   }
00188   
00189   // Loop over all the parameters to make the covariance matrix.
00190   for (Int_t i=0; i<_nComps; i++) {
00191     for (Int_t j=0; j<_nComps; j++) {
00192       const RooRealVar *rowVar= (const RooRealVar*)_fitPars.at(i);
00193       const RooRealVar *colVar= (const RooRealVar*)_fitPars.at(j);
00194       assert(0 != rowVar && 0 != colVar);
00195       _covM(i,j) = rowVar->getError()*colVar->getError()
00196         *_fitRes->correlation(rowVar->GetName(),colVar->GetName());
00197       _V(i,j)=_covM(i,j);
00198     }
00199   }
00200   
00201   // Get the inverse of the covariance matrix
00202   // First check if it's singular
00203   if (_covM.Determinant() == 0) {
00204     cout << "rarSPlot Error: covariance matrix is singular;"
00205          << " I can't invert it!" << endl;
00206     _covM.Print();
00207     _nComps=0;
00208     return;
00209   }
00210   _iV=TMatrixD(TMatrixD::kInverted, _V);
00211   
00212   if (_verbose) {
00213     cout << "Covariance matrix:" << endl;
00214     _V.Print();
00215     cout << "Inverse of covariance matrix:" << endl;
00216     _iV.Print();
00217   }
00218   
00219   // calculate all denominators
00220   _dens.Set(_nEvts);
00221   cout<<" calculate denominators ";
00222   for (Int_t ievt=0; ievt<_nEvts; ievt++) {
00223     if (_verbose&&ievt%1000==0) {
00224       cout << ".";
00225       cout.flush();
00226     }
00227     _data.get(ievt);
00228     _dens[ievt]=0;
00229     for (Int_t k=0; k<_nComps; k++) {
00230       _dens[ievt]+=((RooAbsReal*)_yields.at(_idxMap[k]))->getVal()*
00231         ((RooAbsPdf*)_pdfs.at(_idxMap[k]))->getVal(&_normVars);
00232       if (0==((RooAbsPdf*)_pdfs.at(_idxMap[k]))->getVal(&_normVars)) {
00233         _data.get(ievt);
00234         _dens[ievt]+=((RooAbsReal*)_yields.at(_idxMap[k]))->getVal()*
00235           ((RooAbsPdf*)_pdfs.at(_idxMap[k]))->getVal(&_normVars);
00236       }
00237     }
00238     for (Int_t k=0; k<_nComp0s; k++) {
00239       _dens[ievt]+=((RooAbsReal*)_yield0s.at(k))->getVal()*
00240         ((RooAbsPdf*)_pdf0s.at(k))->getVal(&_normVars);
00241     }
00242   }
00243   cout<<endl;
00244   
00245   // init _sPns
00246   _sPns.Set((_nComps+1)*_nEvts);
00247   _sPnb.Set(_nComps+1);
00248   _sPnb.Reset();
00249   
00250   if (_nComp0s>0) {
00251     // now using eqn 15 in BAD 509 to calculate #_iV
00252     for (Int_t i=0; i<_nComps; i++) {
00253       RooAbsPdf *fi=(RooAbsPdf*)_pdfs.at(_idxMap[i]);
00254       for (Int_t j=i; j<_nComps; j++) {
00255         RooAbsPdf *fj=(RooAbsPdf*)_pdfs.at(_idxMap[j]);
00256         Double_t iVij(0);
00257         for(Int_t ievt=0; ievt<_nEvts; ievt++) {
00258           _data.get(ievt);
00259           iVij+=fi->getVal(&_normVars)*fj->getVal(&_normVars)
00260             /_dens[ievt]/_dens[ievt];
00261         }
00262         _iV(i,j)=iVij;
00263       }
00264     }
00265     for (Int_t i=0; i<_nComps; i++) {
00266       for (Int_t j=0; j<i; j++) {
00267         _iV(i,j)=_iV(j,i);
00268       }
00269     }
00270     
00271     _V=TMatrixD(TMatrixD::kInverted, _iV);
00272     
00273     if (_verbose) {
00274       cout << "Covariance matrix:" << endl;
00275       _V.Print();
00276       cout << "Inverse of covariance matrix:" << endl;
00277       _iV.Print();
00278     }
00279     
00280     // get sP0
00281     fillsPn(_nComps);
00282   }
00283 }

virtual void rarSPlot::setVerbose ( Bool_t  verbose = kTRUE  )  [inline, virtual]

Set verbose mode.

Parameters:
verbose Verbose mode

Definition at line 58 of file rarSPlot.hh.

References _verbose.

00058 {_verbose=verbose;}


Member Data Documentation

TMatrixD rarSPlot::_covM [protected]

Covariance matrix from Minuit.

Definition at line 88 of file rarSPlot.hh.

Referenced by init().

RooDataSet rarSPlot::_data [protected]

Input dataset.

Definition at line 74 of file rarSPlot.hh.

Referenced by fill(), fillsPn(), and init().

TArrayD rarSPlot::_dens [protected]

Denominators.

Definition at line 91 of file rarSPlot.hh.

Referenced by fill(), fillsPn(), and init().

RooArgList rarSPlot::_fitPars [protected]

Fit params.

Definition at line 77 of file rarSPlot.hh.

Referenced by fill(), fillsPn(), and init().

RooFitResult* rarSPlot::_fitRes [protected]

Fit results.

Definition at line 76 of file rarSPlot.hh.

Referenced by init().

TArrayI rarSPlot::_idxMap [protected]

Array for index mapping between _yields and _fitPars.

Definition at line 87 of file rarSPlot.hh.

Referenced by fillsPn(), and init().

TMatrixD rarSPlot::_iV [protected]

Inverse of covM or from direct calculation.

Definition at line 89 of file rarSPlot.hh.

Referenced by init().

Int_t rarSPlot::_nComp0s [protected]

Number of fixed components.

Definition at line 86 of file rarSPlot.hh.

Referenced by fill(), and init().

Int_t rarSPlot::_nComps [protected]

Obs to be used for normalization Number of components.

Definition at line 85 of file rarSPlot.hh.

Referenced by fill(), fillsPn(), and init().

Int_t rarSPlot::_nEvts [protected]

Number of events in dataset.

Definition at line 75 of file rarSPlot.hh.

Referenced by fill(), fillsPn(), and init().

RooArgSet rarSPlot::_normVars [protected]

Obs to be ignored for normalization.

Definition at line 83 of file rarSPlot.hh.

Referenced by fill(), fillsPn(), and init().

RooAbsReal* rarSPlot::_obs [protected]

Obs to do sPlot.

Definition at line 73 of file rarSPlot.hh.

Referenced by fill().

RooArgList rarSPlot::_pdf0s [protected]

Component pdfs w/ fixed yields.

Definition at line 80 of file rarSPlot.hh.

Referenced by init().

RooArgList rarSPlot::_pdfs [protected]

Component pdfs.

Definition at line 78 of file rarSPlot.hh.

Referenced by fill(), fillsPn(), and init().

RooArgSet rarSPlot::_projDeps [protected]

Definition at line 82 of file rarSPlot.hh.

Referenced by init().

RooDataSet* rarSPlot::_sData [protected]

SPlot dataset.

Definition at line 96 of file rarSPlot.hh.

Referenced by fill(), and getSPlotData().

TH1F* rarSPlot::_sHist [protected]

SPlot histogram.

Definition at line 95 of file rarSPlot.hh.

Referenced by fill(), and getSPlotHist().

TArrayI rarSPlot::_sPnb [protected]

Array for sPn fill bit.

Definition at line 93 of file rarSPlot.hh.

Referenced by fillsPn(), and init().

TArrayD rarSPlot::_sPns [protected]

Array for sPn.

Definition at line 92 of file rarSPlot.hh.

Referenced by fill(), fillsPn(), and init().

TMatrixD rarSPlot::_V [protected]

covM or from inverse of iV

Definition at line 90 of file rarSPlot.hh.

Referenced by fill(), fillsPn(), and init().

Bool_t rarSPlot::_verbose [protected]

Boolean to control debug output.

Definition at line 98 of file rarSPlot.hh.

Referenced by fill(), fillsPn(), init(), and setVerbose().

RooArgList rarSPlot::_yield0s [protected]

Fixed component yields.

Definition at line 81 of file rarSPlot.hh.

Referenced by init().

RooArgList rarSPlot::_yields [protected]

Component yields.

Definition at line 79 of file rarSPlot.hh.

Referenced by init().


The documentation for this class was generated from the following files:
Generated on 30 Oct 2013 for RooRarFit by  doxygen 1.4.7