rarTriGauss Class Reference

TripleGaussian PDF/Model/GexpShape builder. More...

#include <rarTriGauss.hh>

Inheritance diagram for rarTriGauss:

rarBasePdf rarConfig List of all members.

Public Member Functions

 rarTriGauss ()
 Trivial ctor.
 rarTriGauss (const char *configFile, const char *configSec, const char *configStr, rarDatasets *theDatasets, RooDataSet *theData, const char *name, const char *title)
 Default ctor.
virtual ~rarTriGauss ()

Protected Member Functions

void init ()
 Initial function called by ctor.

Protected Attributes

RooAbsReal * _x
 Default obs.
RooAbsReal * _tau
 B0 lifetime
RooAbsReal * _msSF
 Scale factor for mean and sigma.
RooAbsReal * _meanC
 Mean of the core Gaussian.
RooAbsReal * _sigmaC
 Sigma of the core Gaussian.
RooAbsReal * _meanT
 Mean of the tail Gaussian.
RooAbsReal * _sigmaT
 Sigma of the tail Gaussian.
RooAbsReal * _meanO
 Mean of the outliner Gaussian.
RooAbsReal * _sigmaO
 Sigma of the outliner Gaussian.
RooAbsReal * _fracC
 Fraction of the core Gaussian.
RooAbsReal * _fracO
 Fraction of the outliner Gaussian.

Private Member Functions

 rarTriGauss (const rarTriGauss &)
 ClassDef (rarTriGauss, 0)

Detailed Description

TripleGaussian PDF/Model/GexpShape builder.

Build TripleGaussian PDF/Model/GexpShape. This class is for convenience of three commonly used composite pdfs related to triple-Gaussian.

Config Directives:
See doc for TriGauss configs.

Definition at line 26 of file rarTriGauss.hh.


Constructor & Destructor Documentation

rarTriGauss::rarTriGauss (  ) 

Trivial ctor.

Usually the objects should be created using other ctors.

Definition at line 43 of file rarTriGauss.cc.

References init().

00044   : rarBasePdf(),
00045     _x(0), _tau(0), _msSF(0),
00046     _meanC(0), _sigmaC(0), _meanT(0), _sigmaT(0), _meanO(0), _sigmaO(0),
00047     _fracC(0), _fracO(0)
00048 {
00049   init();
00050 }

rarTriGauss::rarTriGauss ( const char *  configFile,
const char *  configSec,
const char *  configStr,
rarDatasets theDatasets,
RooDataSet *  theData,
const char *  name,
const char *  title 
)

Default ctor.

Parameters:
configFile The config file
configSec The config section
configStr The config string
theDatasets Available datasets
theData Default dataset for this PDF
name The name
title The title
The default ctor first initializes data members, and then calls init.

Definition at line 64 of file rarTriGauss.cc.

References init().

00067   : rarBasePdf(configFile, configSec, configStr,
00068                theDatasets, theData, name, title),
00069     _x(0), _tau(0), _msSF(0),
00070     _meanC(0), _sigmaC(0), _meanT(0), _sigmaT(0), _meanO(0), _sigmaO(0),
00071     _fracC(0), _fracO(0)
00072 {
00073   init();
00074 }

rarTriGauss::~rarTriGauss (  )  [virtual]

Definition at line 76 of file rarTriGauss.cc.

00077 {
00078 }

rarTriGauss::rarTriGauss ( const rarTriGauss  )  [private]


Member Function Documentation

rarTriGauss::ClassDef ( rarTriGauss  ,
 
) [private]

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

Initial function called by ctor.

init is called by the ctor. It first creates the parameters by calling createAbsReal, and finally it builds RooTripleGaussian PDF/Model based on _pdfType.

Reimplemented from rarBasePdf.

Definition at line 86 of file rarTriGauss.cc.

References _fracC, _fracO, _meanC, _meanO, _meanT, _msSF, rarBasePdf::_obsSet, rarBasePdf::_params, rarBasePdf::_pdfType, _sigmaC, _sigmaO, _sigmaT, rarBasePdf::_subPdfs, _tau, rarBasePdf::_thePdf, _x, rarConfig::createAbsReal(), rarBasePdf::getVarSec(), rarConfig::readConfStr(), and rarBasePdf::setControlBit().

Referenced by rarTriGauss().

00087 {
00088   cout<<"init of rarTriGauss for "<<GetName()<<":"<<endl;
00089   
00090   // first get its dependent/observable
00091   _x=createAbsReal("x", "observable"); assert(_x);
00092   RooRealVar *x=(RooRealVar *)RooArgList(_obsSet).at(0); assert(x);
00093   TString msSFStr=readConfStr("msSF", "notSet", getVarSec());
00094   if (_x!=x) {
00095     // make sure theDep is not derived for most of pdfTypes
00096     if ("TriGauss"!=_pdfType) {
00097       cout <<"No derived dependent allowed for "<<_pdfType<<endl;
00098       exit(-1);
00099     } else if ("notSet"!=msSFStr) {
00100       cout <<"No derived dependent allowed for TriGauss with scaling"<<endl;
00101       exit(-1);
00102     }
00103   }
00104   // Config pdf params
00105   if (("notSet"!=msSFStr)||("TriGauss"!=_pdfType))
00106     _msSF=createAbsReal("msSF", "SF_{#mu#sigma}", 1.);
00107   if ("GexpShape"==_pdfType) _tau=createAbsReal("tau","#tau", 1.5, 0,10, "ps");
00108   _meanC=createAbsReal("meanC", "#mu_{C}",
00109                        (x->getMin()+x->getMax())/2,
00110                        x->getMin(), x->getMax(), _x->getUnit());
00111   _sigmaC=createAbsReal("sigmaC", "#sigma_{C}", .05, 0.,1., _x->getUnit());
00112   _meanT=createAbsReal("meanT", "#mu_{T}",
00113                        (x->getMin()+x->getMax())/2,
00114                        x->getMin(), x->getMax(), _x->getUnit());
00115   _sigmaT=createAbsReal("sigmaT", "#sigma_{T}", .09, 0.,1., _x->getUnit());
00116   _meanO=createAbsReal("meanO", "#mu_{O}",
00117                        (x->getMin()+x->getMax())/2,
00118                        x->getMin(), x->getMax(), _x->getUnit());
00119   _sigmaO=createAbsReal("sigmaO", "#sigma_{O}", .20, 0.,1., _x->getUnit());
00120   _fracC=createAbsReal("fracC", "f_{C}", .80, 0., 1.);
00121   _fracO=createAbsReal("fracO", "f_{O}", .05, 0., 1.);
00122   
00123   _params.Print("v");
00124   
00125   // create pdf
00126   RooAbsPdf *corePdf(0), *tailPdf(0), *outlPdf(0);
00127   if (("TriGauss"==_pdfType)&&("notSet"==msSFStr)) {
00128     // core, tail, and outliner
00129     corePdf=new RooGaussian(Form("core_%s",GetName()),
00130                             Form("Core Gaussian %s", GetTitle()),
00131                             *_x, *_meanC, *_sigmaC);
00132     tailPdf=new RooGaussian(Form("tail_%s",GetName()),
00133                             Form("Tail Gaussian %s", GetTitle()),
00134                             *_x, *_meanT, *_sigmaT);
00135     outlPdf=new RooGaussian(Form("outl_%s",GetName()),
00136                             Form("Outl Gaussian %s", GetTitle()),
00137                             *_x, *_meanO, *_sigmaO);
00138     _thePdf=new
00139       RooAddPdf(Form("the_%s", GetName()), _pdfType+" "+GetTitle(),
00140                 RooArgList(*corePdf, *outlPdf, *tailPdf),
00141                 RooArgList(*_fracC, *_fracO));
00142   } else {
00143     // FlatScaleFactorIntegral?
00144     Bool_t FSFIntegral(kTRUE);
00145     if ("no"==readConfStr("FlatScaleFactorIntegral", "yes", getVarSec()))
00146       FSFIntegral=kFALSE;
00147     // core, tail, and outliner
00148     corePdf=new RooGaussModel(Form("core_%s",GetName()),
00149                               Form("Core Gaussian %s", GetTitle()),
00150                               *x, *_meanC, *_sigmaC, *_msSF);
00151     ((RooGaussModel*)corePdf)->advertiseFlatScaleFactorIntegral(FSFIntegral);
00152     tailPdf=new RooGaussModel(Form("tail_%s",GetName()),
00153                               Form("Tail Gaussian %s", GetTitle()),
00154                               *x, *_meanT, *_sigmaT, *_msSF);
00155     ((RooGaussModel*)tailPdf)->advertiseFlatScaleFactorIntegral(FSFIntegral);
00156     outlPdf=new RooGaussModel(Form("outl_%s",GetName()),
00157                               Form("Outl Gaussian %s", GetTitle()),
00158                               *x, *_meanO, *_sigmaO);
00159     ((RooGaussModel*)outlPdf)->advertiseFlatScaleFactorIntegral(FSFIntegral);
00160     
00161     if ("TriGauss"==_pdfType) {
00162       _thePdf=new
00163         RooAddPdf(Form("the_%s", GetName()), _pdfType+" "+GetTitle(),
00164                   RooArgList(*corePdf, *outlPdf, *tailPdf),
00165                   RooArgList(*_fracC, *_fracO));
00166     } else if ("TriGaussModel"==_pdfType) {
00167       _thePdf=new
00168         RooAddModel(Form("the_%s", GetName()), _pdfType+" "+GetTitle(),
00169                     RooArgList(*corePdf, *outlPdf, *tailPdf),
00170                     RooArgList(*_fracC, *_fracO));
00171     } else if ("GexpShape"==_pdfType) {
00172       // parse the decay type
00173       RooDecay::DecayType decayType=RooDecay::DoubleSided;
00174       TString decayTypeStr=readConfStr("decayType","DoubleSided",getVarSec());
00175       if ("SingleSided"==decayTypeStr) {
00176         decayType=RooDecay::SingleSided;
00177       } else if ("DoubleSided"==decayTypeStr) {
00178         decayType=RooDecay::DoubleSided;
00179       } else if ("Flipped"==decayTypeStr) {
00180         decayType=RooDecay::Flipped;
00181       }
00182       RooAbsPdf *expDecay=new RooDecay(Form("decay_%s",GetName()),
00183                                        Form("decay %s", GetTitle()),
00184                                        *x, *_tau,
00185                                        *((RooGaussModel*)corePdf), decayType);
00186       _subPdfs.add(*expDecay);
00187       _thePdf=new
00188         RooAddPdf(Form("the_%s", GetName()), _pdfType+" "+GetTitle(),
00189                   RooArgList(*expDecay, *outlPdf, *tailPdf),
00190                   RooArgList(*_fracC, *_fracO));
00191     }
00192   }
00193   if ("GexpShape"!=_pdfType) _subPdfs.add(*corePdf);
00194   _subPdfs.add(*tailPdf);
00195   _subPdfs.add(*outlPdf);
00196   // by default do comp plot
00197   setControlBit("CompsOnPlot", "compsOnPlot");
00198 }


Member Data Documentation

RooAbsReal* rarTriGauss::_fracC [protected]

Fraction of the core Gaussian.

Definition at line 47 of file rarTriGauss.hh.

Referenced by init().

RooAbsReal* rarTriGauss::_fracO [protected]

Fraction of the outliner Gaussian.

Definition at line 48 of file rarTriGauss.hh.

Referenced by init().

RooAbsReal* rarTriGauss::_meanC [protected]

Mean of the core Gaussian.

Definition at line 41 of file rarTriGauss.hh.

Referenced by init().

RooAbsReal* rarTriGauss::_meanO [protected]

Mean of the outliner Gaussian.

Definition at line 45 of file rarTriGauss.hh.

Referenced by init().

RooAbsReal* rarTriGauss::_meanT [protected]

Mean of the tail Gaussian.

Definition at line 43 of file rarTriGauss.hh.

Referenced by init().

RooAbsReal* rarTriGauss::_msSF [protected]

Scale factor for mean and sigma.

Definition at line 40 of file rarTriGauss.hh.

Referenced by init().

RooAbsReal* rarTriGauss::_sigmaC [protected]

Sigma of the core Gaussian.

Definition at line 42 of file rarTriGauss.hh.

Referenced by init().

RooAbsReal* rarTriGauss::_sigmaO [protected]

Sigma of the outliner Gaussian.

Definition at line 46 of file rarTriGauss.hh.

Referenced by init().

RooAbsReal* rarTriGauss::_sigmaT [protected]

Sigma of the tail Gaussian.

Definition at line 44 of file rarTriGauss.hh.

Referenced by init().

RooAbsReal* rarTriGauss::_tau [protected]

B0 lifetime

Definition at line 39 of file rarTriGauss.hh.

Referenced by init().

RooAbsReal* rarTriGauss::_x [protected]

Default obs.

Definition at line 38 of file rarTriGauss.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