library: libRooFit
#include "RooAbsPdf.h"

RooAbsPdf


class description - header file - source file
viewCVS header - viewCVS source

class RooAbsPdf: public RooAbsReal

Inheritance Inherited Members Includes Libraries
Class Charts

Function Members (Methods)

Display options:
Show inherited
Show non-public
 
    This is an abstract class, constructors will not be documented.
    Look at the header to check for available constructors.

public:
virtual~RooAbsPdf()
voidTObject::AbstractMethod(const char* method) const
virtual Double_tRooAbsReal::analyticalIntegral(Int_t code, const char* rangeName = "0") const
virtual Double_tanalyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName = "0") const
virtual voidTObject::AppendPad(Option_t* option = "")
voidRooAbsArg::attachDataSet(const RooAbsData& set)
TIterator*RooAbsArg::attribIterator() const
RooAbsFunc*RooAbsReal::bindVars(const RooArgSet& vars, const RooArgSet* nset = 0, Bool_t clipInvalid = kFALSE) const
voidRooAbsArg::branchNodeServerList(RooAbsCollection* list, const RooAbsArg* arg = 0) const
virtual voidTObject::Browse(TBrowser* b)
Bool_tcanBeExtended() const
Bool_tRooAbsArg::checkDependents(const RooArgSet* nset) const
virtual Bool_tRooAbsArg::checkObservables(const RooArgSet* nset) const
static TClass*Class()
virtual const char*TObject::ClassName() const
virtual voidTNamed::Clear(Option_t* option = "")
static voidclearEvalError()
TIterator*RooAbsArg::clientIterator() const
virtual TObject*RooAbsArg::clone(const char* newname) const
virtual TObject*RooAbsArg::Clone(const char* newname = "0") const
virtual Int_tRooAbsArg::Compare(const TObject* other) const
virtual voidRooAbsArg::constOptimize(RooAbsArg::ConstOpCode opcode)
virtual voidTNamed::Copy(TObject& named) const
static voidRooAbsArg::copyList(TList& dest, const TList& source)
virtual RooAbsArg*RooAbsReal::createFundamental(const char* newname = "0") const
TH1*RooAbsReal::createHistogram(const char* name, const RooAbsRealLValue& xvar, const RooCmdArg& arg1 = RooCmdArg::none, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none, const RooCmdArg& arg6 = RooCmdArg::none, const RooCmdArg& arg7 = RooCmdArg::none, const RooCmdArg& arg8 = RooCmdArg::none) const
RooAbsReal*RooAbsReal::createIntegral(const RooArgSet& iset, const char* rangeName) const
RooAbsReal*RooAbsReal::createIntegral(const RooArgSet& iset, const RooArgSet& nset, const char* rangeName = "0") const
RooAbsReal*RooAbsReal::createIntegral(const RooArgSet& iset, const RooNumIntConfig& cfg, const char* rangeName = "0") const
RooAbsReal*RooAbsReal::createIntegral(const RooArgSet& iset, const RooArgSet& nset, const RooNumIntConfig& cfg, const char* rangeName = "0") const
virtual RooAbsReal*RooAbsReal::createIntegral(const RooArgSet& iset, const RooArgSet* nset = 0, const RooNumIntConfig* cfg = 0, const char* rangeName = "0") const
RooAbsReal*RooAbsReal::createIntegral(const RooArgSet& iset, const RooCmdArg arg1, const RooCmdArg arg2 = RooCmdArg::none, const RooCmdArg arg3 = RooCmdArg::none, const RooCmdArg arg4 = RooCmdArg::none, const RooCmdArg arg5 = RooCmdArg::none, const RooCmdArg arg6 = RooCmdArg::none, const RooCmdArg arg7 = RooCmdArg::none, const RooCmdArg arg8 = RooCmdArg::none) const
const RooAbsReal*RooAbsReal::createProjection(const RooArgSet& depVars, const RooArgSet& projVars) const
const RooAbsReal*RooAbsReal::createProjection(const RooArgSet& depVars, const RooArgSet& projVars, RooArgSet*& cloneSet) const
virtual Double_tRooAbsReal::defaultErrorLevel() const
static RooNumIntConfig*RooAbsReal::defaultIntegratorConfig()
static ostream&RooPrintable::defaultStream(ostream* os = 0)
virtual voidTObject::Delete(Option_t* option = "")
Bool_tRooAbsArg::deleteWatch() const
Bool_tRooAbsArg::dependentOverlaps(const RooAbsData* dset, const RooAbsArg& testArg) const
Bool_tRooAbsArg::dependentOverlaps(const RooArgSet* depList, const RooAbsArg& testArg) const
Bool_tRooAbsArg::dependsOn(const RooAbsCollection& serverList, const RooAbsArg* ignoreArg = 0) const
Bool_tRooAbsArg::dependsOn(const RooAbsArg& server, const RooAbsArg* ignoreArg = 0) const
virtual Int_tTObject::DistancetoPrimitive(Int_t px, Int_t py)
virtual voidTObject::Draw(Option_t* option = "")
virtual voidTObject::DrawClass() const
virtual TObject*TObject::DrawClone(Option_t* option = "") const
virtual voidTObject::Dump() const
virtual voidTObject::Error(const char* method, const char* msgfmt) const
static Bool_tevalError()
virtual voidTObject::Execute(const char* method, const char* params, Int_t* error = 0)
virtual voidTObject::Execute(TMethod* method, TObjArray* params, Int_t* error = 0)
virtual voidTObject::ExecuteEvent(Int_t event, Int_t px, Int_t py)
virtual Double_texpectedEvents(const RooArgSet* nset) const
virtual Double_texpectedEvents(const RooArgSet& nset) const
virtual Double_textendedTerm(UInt_t observedEvents, const RooArgSet* nset = 0) const
virtual RooAbsPdf::ExtendModeextendMode() const
virtual voidTObject::Fatal(const char* method, const char* msgfmt) const
virtual voidTNamed::FillBuffer(char*& buffer)
TH1*RooAbsReal::fillHistogram(TH1* hist, const RooArgList& plotVars, Double_t scaleFactor = 1, const RooArgSet* projectedVars = 0) const
virtual TObject*TObject::FindObject(const char* name) const
virtual TObject*TObject::FindObject(const TObject* obj) const
RooAbsArg*RooAbsArg::findServer(const char* name) const
RooAbsArg*RooAbsArg::findServer(const RooAbsArg& arg) const
RooAbsArg*RooAbsArg::findServer(Int_t index) const
virtual RooFitResult*fitTo(RooAbsData& data, const RooLinkedList& cmdList)
virtual RooFitResult*fitTo(RooAbsData& data, Option_t* fitOpt = "", Option_t* optOpt = "c", const char* fitRange = "0")
virtual RooFitResult*fitTo(RooAbsData& data, const RooArgSet& projDeps, Option_t* fitOpt = "", Option_t* optOpt = "c", const char* fitRange = "0")
virtual RooFitResult*fitTo(RooAbsData& data, RooCmdArg arg1, RooCmdArg arg2 = RooCmdArg::none, RooCmdArg arg3 = RooCmdArg::none, RooCmdArg arg4 = RooCmdArg::none, RooCmdArg arg5 = RooCmdArg::none, RooCmdArg arg6 = RooCmdArg::none, RooCmdArg arg7 = RooCmdArg::none, RooCmdArg arg8 = RooCmdArg::none)
virtual voidfixAddCoefNormalization(const RooArgSet& addNormSet = RooArgSet())
virtual voidfixAddCoefRange(const char* rangeName = "0")
virtual Bool_tRooAbsReal::forceAnalyticalInt(const RooAbsArg&) const
virtual voidRooAbsReal::forceNumInt(Bool_t flag = kTRUE)
RooDataSet*generate(const RooArgSet& whatVars, Int_t nEvents = 0, Bool_t verbose = kFALSE) const
RooDataSet*generate(const RooArgSet& whatVars, const RooDataSet& prototype, Int_t nEvents = 0, Bool_t verbose = kFALSE, Bool_t randProtoOrder = kFALSE, Bool_t resampleProto = kFALSE) const
RooDataSet*generate(const RooArgSet& whatVars, Int_t nEvents, const RooCmdArg& arg1, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none)
RooDataSet*generate(const RooArgSet& whatVars, const RooCmdArg& arg1 = RooCmdArg::none, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none, const RooCmdArg& arg6 = RooCmdArg::none)
virtual voidgenerateEvent(Int_t code)
virtual Int_tRooAbsReal::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName = "0") const
virtual Int_tRooAbsReal::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet, const char* rangeName = "0") const
Bool_tRooAbsArg::getAttribute(const Text_t* name) const
RooLinkedListRooAbsArg::getCloningAncestors() const
RooArgSet*RooAbsArg::getComponents() const
RooArgSet*RooAbsArg::getDependents(const RooArgSet& set) const
RooArgSet*RooAbsArg::getDependents(const RooAbsData* set) const
RooArgSet*RooAbsArg::getDependents(const RooArgSet* depList) const
virtual Option_t*TObject::GetDrawOption() const
static Long_tTObject::GetDtorOnly()
virtual Int_tgetGenerator(const RooArgSet& directVars, RooArgSet& generateVars, Bool_t staticInitOK = kTRUE) const
virtual const char*TObject::GetIconName() const
const RooNumIntConfig*RooAbsReal::getIntegratorConfig() const
Double_tgetLogVal(const RooArgSet* set = 0) const
virtual Int_tRooAbsReal::getMaxVal(const RooArgSet& vars) const
virtual const char*TNamed::GetName() const
Double_tgetNorm(const RooArgSet& nset) const
virtual Double_tgetNorm(const RooArgSet* set = 0) const
virtual const RooAbsReal*getNormObj(const RooArgSet* set, const RooArgSet* iset, const TNamed* rangeName = 0) const
virtual char*TObject::GetObjectInfo(Int_t px, Int_t py) const
static Bool_tTObject::GetObjectStat()
RooArgSet*RooAbsArg::getObservables(const RooArgSet& set) const
RooArgSet*RooAbsArg::getObservables(const RooAbsData* data) const
RooArgSet*RooAbsArg::getObservables(const RooAbsData& data) const
virtual RooArgSet*RooAbsArg::getObservables(const RooArgSet* depList) const
virtual Option_t*TObject::GetOption() const
RooArgSet*RooAbsArg::getParameters(const RooAbsData* data) const
RooArgSet*RooAbsArg::getParameters(const RooAbsData& data) const
RooArgSet*RooAbsArg::getParameters(const RooArgSet& set) const
virtual RooArgSet*RooAbsArg::getParameters(const RooArgSet* depList) const
virtual Int_tRooAbsReal::getPlotBins() const
const char*RooAbsReal::getPlotLabel() const
Double_tRooAbsReal::getPlotMax() const
Double_tRooAbsReal::getPlotMin() const
virtual const char*TNamed::GetTitle() const
TStringRooAbsReal::getTitle(Bool_t appendUnit = kFALSE) const
virtual UInt_tTObject::GetUniqueID() const
const Text_t*RooAbsReal::getUnit() const
virtual Double_tgetVal(const RooArgSet* set = 0) const
RooArgSet*RooAbsArg::getVariables() const
virtual Bool_tTObject::HandleTimer(TTimer* timer)
virtual ULong_tTNamed::Hash() const
virtual Bool_tRooAbsArg::hasRange(const char*) const
virtual voidTObject::Info(const char* method, const char* msgfmt) const
virtual Bool_tTObject::InheritsFrom(const char* classname) const
virtual Bool_tTObject::InheritsFrom(const TClass* cl) const
virtual voidinitGenerator(Int_t code)
static voidRooPrintable::inLinePrint(ostream& os, const TNamed& named)
virtual Bool_tRooAbsReal::inPlotRange(Double_t value) const
virtual Bool_tRooAbsArg::inRange(const char*) const
virtual voidTObject::Inspect() const
voidTObject::InvertBit(UInt_t f)
virtual TClass*IsA() const
Bool_tRooAbsArg::isCloneOf(const RooAbsArg& other) const
Bool_tRooAbsArg::isConstant() const
virtual Bool_tRooAbsArg::isDerived() const
virtual Bool_tisDirectGenSafe(const RooAbsArg& arg) const
virtual Bool_tTObject::IsEqual(const TObject* obj) const
virtual Bool_tTObject::IsFolder() const
virtual Bool_tRooAbsArg::isFundamental() const
virtual Bool_tRooAbsArg::isLValue() const
Bool_tTObject::IsOnHeap() const
Bool_tisSelectedComp() const
Bool_tRooAbsArg::isShapeServer(const RooAbsArg& arg) const
Bool_tRooAbsArg::isShapeServer(const char* name) const
virtual Bool_tRooAbsArg::IsSortable() const
Bool_tRooAbsArg::isValueServer(const RooAbsArg& arg) const
Bool_tRooAbsArg::isValueServer(const char* name) const
Bool_tTObject::IsZombie() const
voidRooAbsArg::leafNodeServerList(RooAbsCollection* list, const RooAbsArg* arg = 0) const
RooPrintable::PrintOptionRooPrintable::lessVerbose(RooPrintable::PrintOption opt) const
virtual voidTNamed::ls(Option_t* option = "") const
virtual Double_tRooAbsReal::maxVal(Int_t code)
voidTObject::MayNotUse(const char* method) const
Bool_tmustBeExtended() const
static voidRooAbsArg::nameFieldLength(Int_t newLen)
virtual Bool_tTObject::Notify()
Bool_tRooAbsArg::observableOverlaps(const RooAbsData* dset, const RooAbsArg& testArg) const
Bool_tRooAbsArg::observableOverlaps(const RooArgSet* depList, const RooAbsArg& testArg) const
static voidRooPrintable::oneLinePrint(ostream& os, const TNamed& named)
static voidTObject::operator delete(void* ptr)
static voidTObject::operator delete(void* ptr, void* vp)
static voidTObject::operator delete[](void* ptr)
static voidTObject::operator delete[](void* ptr, void* vp)
void*TObject::operator new(size_t sz)
void*TObject::operator new(size_t sz, void* vp)
void*TObject::operator new[](size_t sz)
void*TObject::operator new[](size_t sz, void* vp)
TNamed&TNamed::operator=(const TNamed& rhs)
Bool_tRooAbsReal::operator==(Double_t value) const
virtual Bool_tRooAbsReal::operator==(const RooAbsArg& other)
Bool_tRooAbsArg::overlaps(const RooAbsArg& testArg) const
virtual voidTObject::Paint(Option_t* option = "")
virtual RooPlot*paramOn(RooPlot* frame, const RooAbsData* data, const char* label = "", Int_t sigDigits = 2, Option_t* options = "NELU", Double_t xmin = 0.65, Double_t xmax = 0.99, Double_t ymax = 0.95)
virtual RooPlot*paramOn(RooPlot* frame, const RooCmdArg& arg1 = RooCmdArg::none, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none, const RooCmdArg& arg6 = RooCmdArg::none, const RooCmdArg& arg7 = RooCmdArg::none, const RooCmdArg& arg8 = RooCmdArg::none)
RooPrintable::PrintOptionRooPrintable::parseOptions(Option_t* options) const
virtual RooPlot*plotCompOn(RooPlot* frame, const char* compNameList, Option_t* drawOptions = "L", Double_t scaleFactor = 1.0, RooAbsReal::ScaleType stype = Relative, const RooAbsData* projData = 0, const RooArgSet* projSet = 0) const
virtual RooPlot*plotCompOn(RooPlot* frame, const RooArgSet& compSet, Option_t* drawOptions = "L", Double_t scaleFactor = 1.0, RooAbsReal::ScaleType stype = Relative, const RooAbsData* projData = 0, const RooArgSet* projSet = 0) const
virtual RooPlot*plotCompSliceOn(RooPlot* frame, const char* compNameList, const RooArgSet& sliceSet, Option_t* drawOptions = "L", Double_t scaleFactor = 1.0, RooAbsReal::ScaleType stype = Relative, const RooAbsData* projData = 0) const
virtual RooPlot*plotCompSliceOn(RooPlot* frame, const RooArgSet& compSet, const RooArgSet& sliceSet, Option_t* drawOptions = "L", Double_t scaleFactor = 1.0, RooAbsReal::ScaleType stype = Relative, const RooAbsData* projData = 0) const
RooPlot*plotNLLOn(RooPlot* frame, RooDataSet* data, Option_t* drawOptions = "L", Double_t prec = 1e-2, Bool_t fixMinToZero = kTRUE)
virtual RooPlot*plotNLLOn(RooPlot* frame, RooDataSet* data, Bool_t extended, Option_t* drawOptions = "L", Double_t prec = 1e-2, Bool_t fixMinToZero = kTRUE)
virtual RooPlot*plotNLLOn(RooPlot* frame, RooDataSet* data, Bool_t extended, const RooArgSet& projDeps, Option_t* drawOptions = "L", Double_t prec = 1e-2, Bool_t fixMinToZero = kTRUE)
virtual RooPlot*plotOn(RooPlot* frame, const RooCmdArg& arg1 = RooCmdArg::none, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none, const RooCmdArg& arg6 = RooCmdArg::none, const RooCmdArg& arg7 = RooCmdArg::none, const RooCmdArg& arg8 = RooCmdArg::none, const RooCmdArg& arg9 = RooCmdArg::none, const RooCmdArg& arg10 = RooCmdArg::none) const
virtual RooPlot*RooAbsReal::plotSliceOn(RooPlot* frame, const RooArgSet& sliceSet, Option_t* drawOptions = "L", Double_t scaleFactor = 1.0, RooAbsReal::ScaleType stype = Relative, const RooAbsData* projData = 0) const
virtual voidTObject::Pop()
virtual voidRooAbsArg::Print(Option_t* options = "0") const
voidRooAbsArg::printCompactTree(const char* indent = "", const char* fileName = "0", const char* namePat = "0")
voidRooAbsArg::printCompactTree(ostream& os, const char* indent = "", const char* namePat = "0")
virtual voidRooAbsArg::printCompactTreeHook(ostream& os, const char* ind = "")
voidRooAbsArg::printDirty(Bool_t depth = kTRUE) const
virtual voidprintToStream(ostream& stream, RooPrintable::PrintOption opt = Standard, TString indent = ) const
virtual Int_tTObject::Read(const char* name)
virtual Bool_tRooAbsReal::readFromStream(istream& is, Bool_t compact, Bool_t verbose = kFALSE)
Bool_tRooAbsArg::recursiveCheckDependents(const RooArgSet* nset) const
Bool_tRooAbsArg::recursiveCheckObservables(const RooArgSet* nset) const
virtual voidTObject::RecursiveRemove(TObject* obj)
voidTObject::ResetBit(UInt_t f)
virtual voidresetErrorCounters(Int_t resetValue = 10)
virtual voidTObject::SaveAs(const char* filename = "", Option_t* option = "") const
virtual voidTObject::SavePrimitive(ostream& out, Option_t* option = "")
virtual Bool_tselfNormalized() const
TIterator*RooAbsArg::serverIterator() const
voidRooAbsArg::setAttribute(const Text_t* name, Bool_t value = kTRUE)
voidTObject::SetBit(UInt_t f)
voidTObject::SetBit(UInt_t f, Bool_t set)
static voidRooAbsReal::setCacheCheck(Bool_t flag)
voidRooAbsArg::setDeleteWatch(Bool_t flag = kTRUE)
static voidRooAbsArg::setDirtyInhibit(Bool_t flag)
virtual voidTObject::SetDrawOption(Option_t* option = "")
static voidTObject::SetDtorOnly(void* obj)
voidRooAbsReal::setIntegratorConfig()
voidRooAbsReal::setIntegratorConfig(const RooNumIntConfig& config)
virtual voidTNamed::SetName(const char* name)
virtual voidTNamed::SetNameTitle(const char* name, const char* title)
static voidTObject::SetObjectStat(Bool_t stat)
voidRooAbsReal::setPlotBins(Int_t value)
voidRooAbsReal::setPlotLabel(const char* label)
voidRooAbsReal::setPlotMax(Double_t value)
voidRooAbsReal::setPlotMin(Double_t value)
voidRooAbsReal::setPlotRange(Double_t min, Double_t max)
virtual voidTNamed::SetTitle(const char* title = "")
voidsetTraceCounter(Int_t value, Bool_t allNodes = kFALSE)
virtual voidTObject::SetUniqueID(UInt_t uid)
voidRooAbsReal::setUnit(const char* unit)
TIterator*RooAbsArg::shapeClientIterator() const
virtual voidShowMembers(TMemberInspector& insp, char* parent)
virtual Int_tTNamed::Sizeof() const
RooNumIntConfig*RooAbsReal::specialIntegratorConfig() const
virtual voidStreamer(TBuffer& b)
voidStreamerNVirtual(TBuffer& b)
virtual voidTObject::SysError(const char* method, const char* msgfmt) const
Bool_tTObject::TestBit(UInt_t f) const
Int_tTObject::TestBits(UInt_t f) const
virtual Bool_ttraceEvalHook(Double_t value) const
Bool_ttraceEvalPdf(Double_t value) const
voidRooAbsArg::treeNodeServerList(RooAbsCollection* list, const RooAbsArg* arg = 0, Bool_t doBranch = kTRUE, Bool_t doLeaf = kTRUE, Bool_t valueOnly = kFALSE) const
virtual voidTObject::UseCurrentStyle()
TIterator*RooAbsArg::valueClientIterator() const
static voidRooAbsArg::verboseDirty(Bool_t flag)
static voidverboseEval(Int_t stat)
virtual voidTObject::Warning(const char* method, const char* msgfmt) const
virtual Int_tTObject::Write(const char* name = "0", Int_t option = 0, Int_t bufsize = 0)
virtual Int_tTObject::Write(const char* name = "0", Int_t option = 0, Int_t bufsize = 0) const
virtual voidRooAbsReal::writeToStream(ostream& os, Bool_t compact) const
protected:
voidRooAbsArg::addServer(RooAbsArg& server, Bool_t valueProp = kTRUE, Bool_t shapeProp = kFALSE)
voidRooAbsArg::addServerList(RooAbsCollection& serverList, Bool_t valueProp = kTRUE, Bool_t shapeProp = kFALSE)
Bool_tRooAbsReal::allClientsCached(RooAbsArg* var, RooArgSet& cacheList)
virtual voidRooAbsReal::attachToTree(TTree& t, Int_t bufSize = 32000)
voidRooAbsArg::changeServer(RooAbsArg& server, Bool_t valueProp, Bool_t shapeProp)
TStringRooAbsArg::cleanBranchName() const
voidRooAbsArg::clearShapeDirty() const
voidRooAbsArg::clearValueDirty() const
virtual voidRooAbsReal::copyCache(const RooAbsArg* source)
UInt_tRooAbsArg::crc32(const char* data) const
const RooAbsReal*RooAbsReal::createProjection(const RooArgSet& dependentVars, const RooArgSet* projectedVars, RooArgSet*& cloneSet, const char* rangeName = "0") const
voidRooAbsReal::doConstOpt(RooAbsData& dataset, const RooArgSet* normSet, Bool_t verbose)
virtual voidTObject::DoError(int level, const char* location, const char* fmt, va_list va) const
virtual Double_tRooAbsReal::evaluate() const
virtual voidRooAbsReal::fillTreeBranch(TTree& t)
Bool_tRooAbsReal::findCacheableBranches(RooAbsArg* arg, RooAbsData* dset, RooArgSet& cacheList, const RooArgSet* normSet, Bool_t verbose)
RooAbsArg*RooAbsArg::findNewServer(const RooAbsCollection& newSet, Bool_t nameChange) const
voidRooAbsReal::findRedundantCacheServers(RooAbsData* dset, RooArgSet& cacheList, RooArgSet& pruneList, Bool_t verbose)
voidRooAbsReal::findUnusedDataVariables(RooAbsData* dset, RooArgSet& pruneList, Bool_t verbose)
virtual RooAbsGenContext*genContext(const RooArgSet& vars, const RooDataSet* prototype = 0, const RooArgSet* auxProto = 0, Bool_t verbose = kFALSE) const
virtual voidRooAbsArg::getObservablesHook(const RooArgSet*, RooArgSet*) const
virtual voidRooAbsArg::getParametersHook(const RooArgSet*, RooArgSet*) const
RooAbsProxy*RooAbsArg::getProxy(Int_t index) const
static voidglobalSelectComp(Bool_t flag)
TStringRooAbsReal::integralNameSuffix(const RooArgSet& iset, const RooArgSet* nset = 0, const char* rangeName = "0") const
Bool_tRooAbsArg::isShapeDirty() const
virtual Bool_tRooAbsReal::isValid() const
virtual Bool_tRooAbsReal::isValidReal(Double_t value, Bool_t printError = kFALSE) const
Bool_tRooAbsArg::isValueDirty() const
voidRooAbsReal::makeProjectionSet(const RooAbsArg* plotVar, const RooArgSet* allVars, RooArgSet& projectedVars, Bool_t silent) const
voidTObject::MakeZombie()
Bool_tRooAbsReal::matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, const RooArgProxy& a) const
Bool_tRooAbsReal::matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, const RooArgSet& set) const
Bool_tRooAbsReal::matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, const RooArgProxy& a, const RooArgProxy& b) const
Bool_tRooAbsReal::matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, const RooArgProxy& a, const RooArgProxy& b, const RooArgProxy& c) const
Bool_tRooAbsReal::matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, const RooArgProxy& a, const RooArgProxy& b, const RooArgProxy& c, const RooArgProxy& d) const
Int_tRooAbsArg::numProxies() const
RooAbsArg::OperModeRooAbsArg::operMode() const
virtual voidoperModeHook()
voidRooAbsReal::optimizeDirty(RooAbsData& dataset, const RooArgSet* normSet, Bool_t verbose)
virtual RooPlot*paramOn(RooPlot* frame, const RooArgSet& params, Bool_t showConstants = kFALSE, const char* label = "", Int_t sigDigits = 2, Option_t* options = "NELU", Double_t xmin = 0.65, Double_t xmax = 0.99, Double_t ymax = 0.95, const RooCmdArg* formatCmd = 0)
virtual RooPlot*RooAbsReal::plotAsymOn(RooPlot* frame, const RooAbsCategoryLValue& asymCat, RooAbsReal::PlotOpt o) const
virtual RooPlot*plotCompOnEngine(RooPlot* frame, RooArgSet* selNodes, Option_t* drawOptions = "L", Double_t scaleFactor = 1.0, RooAbsReal::ScaleType stype = Relative, const RooAbsData* projData = 0, const RooArgSet* projSet = 0) const
virtual RooPlot*plotOn(RooPlot* frame, RooLinkedList& cmdList) const
virtual RooPlot*plotOn(RooPlot* frame, RooAbsReal::PlotOpt o) const
voidplotOnCompSelect(RooArgSet* selNodes) const
Bool_tRooAbsReal::plotSanityChecks(RooPlot* frame) const
voidRooAbsArg::printAttribList(ostream& os) const
static voidraiseEvalError()
Int_t*randomizeProtoOrder(Int_t nProto, Int_t nGen, Bool_t resample = kFALSE) const
Bool_tRooAbsArg::recursiveRedirectServers(const RooAbsCollection& newServerList, Bool_t mustReplaceAll = kFALSE, Bool_t nameChange = kFALSE)
Bool_tRooAbsArg::redirectServers(const RooAbsCollection& newServerList, Bool_t mustReplaceAll = kFALSE, Bool_t nameChange = kFALSE, Bool_t isRecursionStep = kFALSE)
virtual Bool_tredirectServersHook(const RooAbsCollection& newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t isRecursive)
voidRooAbsArg::registerProxy(RooArgProxy& proxy)
voidRooAbsArg::registerProxy(RooSetProxy& proxy)
voidRooAbsArg::registerProxy(RooListProxy& proxy)
voidRooAbsArg::removeServer(RooAbsArg& server, Bool_t force = kFALSE)
voidRooAbsArg::replaceServer(RooAbsArg& oldServer, RooAbsArg& newServer, Bool_t valueProp, Bool_t shapeProp)
voidselectComp(Bool_t flag)
virtual voidRooAbsReal::selectNormalization(const RooArgSet* depSet = 0, Bool_t force = kFALSE)
virtual voidRooAbsReal::selectNormalizationRange(const char* rangeName = "0", Bool_t force = kFALSE)
virtual voidRooAbsArg::serverNameChangeHook(const RooAbsArg*, const RooAbsArg*)
voidRooAbsArg::setOperMode(RooAbsArg::OperMode mode, Bool_t recurseADirty = kTRUE)
voidRooAbsArg::setProxyNormSet(const RooArgSet* nset)
voidRooAbsArg::setShapeDirty() const
virtual voidRooAbsReal::setTreeBranchStatus(TTree& t, Bool_t active)
voidRooAbsArg::setValueDirty() const
virtual voidRooAbsReal::syncCache(const RooArgSet* set = 0)
virtual Bool_tsyncNormalization(const RooArgSet* dset, Bool_t adjustProxies = kTRUE) const
virtual voidsyncNormalizationPostHook(RooAbsReal* norm, const RooArgSet* dset) const
virtual Bool_tsyncNormalizationPreHook(RooAbsReal* norm, const RooArgSet* dset) const
Double_tRooAbsReal::traceEval(const RooArgSet* set) const
voidRooAbsReal::undoConstOpt(RooAbsData& dataset, const RooArgSet* normSet, Bool_t verbose)
voidRooAbsArg::unRegisterProxy(RooArgProxy& proxy)
voidRooAbsArg::unRegisterProxy(RooSetProxy& proxy)
voidRooAbsArg::unRegisterProxy(RooListProxy& proxy)

Data Members

public:
enum ExtendMode { CanNotBeExtended
CanBeExtended
MustBeExtended
};
enum RooAbsReal::ScaleType { Raw
Relative
NumEvent
RelativeExpected
};
enum RooAbsArg::ConstOpCode { Activate
DeActivate
ConfigChange
ValueChange
};
enum RooAbsArg::OperMode { Auto
AClean
ADirty
};
enum TObject::EStatusBits { kCanDelete
kMustCleanup
kObjInCanvas
kIsReferenced
kHasUUID
kCannotPick
kNoContextMenu
kInvalidObject
};
enum TObject::[unnamed] { kIsOnHeap
kNotDeleted
kZombie
kBitMask
kSingleKey
kOverwrite
kWriteDelete
};
enum RooPrintable::PrintOption { InLine
OneLine
Standard
Shape
Verbose
};
protected:
static Int_t_verboseEval
Double_t_rawValue
RooAbsReal*_normNormalization integral (owned by _normMgr)
RooArgSet*_normSetNormalization set with for above integral
RooNormManager_normMgrNormalization manager
Int_t_errorCountNumber of errors remaining to print
Int_t_traceCountNumber of traces remaining to print
Int_t_negCountNumber of negative probablities remaining to print
Bool_t_selectCompComponent selection flag for RooAbsPdf::plotCompOn
static Bool_t_globalSelectCompGlobal activation switch for component selection
static Bool_t_evalError
Double_tRooAbsReal::_plotMinMinimum of plot range
Double_tRooAbsReal::_plotMaxMaximum of plot range
Int_tRooAbsReal::_plotBinsNumber of plot bins
Double_tRooAbsReal::_valueCache for current value of object
TStringRooAbsReal::_unitUnit for objects value
TStringRooAbsReal::_labelPlot label for objects value
Bool_tRooAbsReal::_forceNumIntForce numerical integration if flag set
RooNumIntConfig*RooAbsReal::_specIntegratorConfig! Numeric integrator configuration specific for this object
static Bool_tRooAbsReal::_cacheCheck
RooRefCountListRooAbsArg::_serverList! list of server objects
RooRefCountListRooAbsArg::_clientList! list of client objects
RooRefCountListRooAbsArg::_clientListShape! subset of clients that requested shape dirty flag propagation
RooRefCountListRooAbsArg::_clientListValue! subset of clients that requested value dirty flag propagation
TListRooAbsArg::_proxyList! list of proxies
TIterator*RooAbsArg::_clientShapeIter! Iterator over _clientListShape
TIterator*RooAbsArg::_clientValueIter! Iterator over _clientListValue
THashListRooAbsArg::_attribListList of string attributes
static Bool_tRooAbsArg::_verboseDirtyStatic flag controlling verbose messaging for dirty state changes
static Bool_tRooAbsArg::_inhibitDirtyStatic flag controlling global inhibit of dirty state propagation
Bool_tRooAbsArg::_deleteWatch! Delete watch flag
static Int_tRooAbsArg::_nameLength
TStringTNamed::fNameobject identifier
TStringTNamed::fTitleobject title

Class Description

~RooAbsPdf()
 Destructor
if (_norm) delete _norm ;
Double_t getVal(const RooArgSet* nset)
 Return current value, normalizated by integrating over
 the dependents in 'nset'. If 'nset' is 0, the unnormalized value. 
 is returned. All elements of 'nset' must be lvalues
Double_t analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName)
 Analytical integral with normalization (see RooAbsReal::analyticalIntegralWN() for further information)

 This function applies the normalization specified by 'normSet' to the integral returned
 by RooAbsReal::analyticalIntegral(). The passthrough scenario (code=0) is also changed
 to return a normalized answer
Bool_t traceEvalPdf(Double_t value)
 Check that passed value is positive and not 'not-a-number'.
 If not, print an error, until the error counter reaches
 its set maximum.
Double_t getNorm(const RooArgSet* nset)
 Return the integral of this PDF over all elements of 'nset'. 
const RooAbsReal* getNormObj(const RooArgSet* nset, const RooArgSet* iset, const TNamed* rangeName)
 Check normalization is already stored
Bool_t syncNormalizationPreHook(RooAbsReal*,const RooArgSet*)
void syncNormalizationPostHook(RooAbsReal*,const RooArgSet*)
Bool_t syncNormalization(const RooArgSet* nset, Bool_t adjustProxies)
 Verify that the normalization integral cached with this PDF
 is valid for given set of normalization dependents

 If not, the cached normalization integral (if any) is deleted
 and a new integral is constructed for use with 'nset'
 Elements in 'nset' can be discrete and real, but must be lvalues

 By default, only actual dependents of the PDF listed in 'nset'
 are integration. This behaviour can be modified in subclasses
 by overloading the syncNormalizationPreHook() function. 
 
 For functions that declare to be self-normalized by overloading the
 selfNormalized() function, a unit normalization is always constructed
Bool_t traceEvalHook(Double_t value)
 WVE 08/21/01 Probably obsolete now.
void resetErrorCounters(Int_t resetValue)
 Reset error counter to given value, limiting the number
 of future error messages for this pdf to 'resetValue'
void setTraceCounter(Int_t value, Bool_t allNodes)
 Reset trace counter to given value, limiting the
 number of future trace messages for this pdf to 'value'
void operModeHook()
 WVE 08/21/01 Probably obsolete now
Double_t getLogVal(const RooArgSet* nset)
 Return the log of the current value with given normalization
 An error message is printed if the argument of the log is negative.
Double_t extendedTerm(UInt_t observed, const RooArgSet* nset)
 Returned the extended likelihood term (Nexpect - Nobserved*log(NExpected)
 of this PDF for the given number of observed events

 For successfull operation the PDF implementation must indicate
 it is extendable by overloading canBeExtended() and must
 implemented the expectedEvents() function.
RooFitResult* fitTo(RooAbsData& data, RooCmdArg arg1, RooCmdArg arg2, RooCmdArg arg3, RooCmdArg arg4, RooCmdArg arg5, RooCmdArg arg6, RooCmdArg arg7, RooCmdArg arg8)
 Fit PDF to given dataset. If dataset is unbinned, an unbinned maximum likelihood is performed. If the dataset
 is binned, a binned maximum likelihood is performed. By default the fit is executed through the MINUIT
 commands MIGRAD, HESSE and MINOS in succession.

 The following named arguments are supported

 Options to control construction of -log(L)

 ConditionalObservables(const RooArgSet& set) -- Do not normalize PDF over listed observables
 Extended(Bool_t flag)           -- Add extended likelihood term, off by default
 Range(const char* name)         -- Fit only data inside range with given name
 Range(Double_t lo, Double_t hi) -- Fit only data inside given range. A range named "fit" is created on the fly on all observables.
 NumCPU(int num)                 -- Parallelize NLL calculation on num CPUs
 Optimize(Bool_t flag)           -- Activate constant term optimization (on by default)
 SplitRange(Bool_t flag)         -- Use separate fit ranges in a simultaneous fit. Actual range name for each
                                    subsample is assumed to by rangeName_{indexState} where indexState
                                    is the state of the master index category of the simultaneous fit

 Options to control flow of fit procedure

 InitialHesse(Bool_t flag)      -- Flag controls if HESSE before MIGRAD as well, off by default
 Hesse(Bool_t flag)             -- Flag controls if HESSE is run after MIGRAD, on by default
 Minos(Bool_t flag)             -- Flag controls if MINOS is run after HESSE, on by default
 Minos(const RooArgSet& set)    -- Only run MINOS on given subset of arguments
 Save(Bool_t flag)              -- Flac controls if RooFitResult object is produced and returned, off by default
 Strategy(Int_t flag)           -- Set Minuit strategy (0 through 2, default is 1)
 FitOptions(const char* optStr) -- Steer fit with classic options string (for backward compatibility). Use of this option
                                   excludes use of any of the new style steering options.

 Options to control informational output

 Verbose(Bool_t flag)           -- Flag controls if verbose output is printed (NLL, parameter changes during fit
 Timer(Bool_t flag)             -- Time CPU and wall clock consumption of fit steps, off by default
 PrintLevel(Int_t level)        -- Set Minuit print level (-1 through 3, default is 1). At -1 all RooFit informational 
                                   messages are suppressed as well
 
 
RooFitResult* fitTo(RooAbsData& data, const RooLinkedList& cmdList)
 Fit PDF to given dataset. If dataset is unbinned, an unbinned maximum likelihood is performed. If the dataset
 is binned, a binned maximum likelihood is performed. By default the fit is executed through the MINUIT
 commands MIGRAD, HESSE and MINOS in succession.

 See RooAbsPdf::fitTo(RooAbsData& data, RooCmdArg arg1, RooCmdArg arg2, RooCmdArg arg3, RooCmdArg arg4, 
                                         RooCmdArg arg5, RooCmdArg arg6, RooCmdArg arg7, RooCmdArg arg8) 

 for documentation of options
RooFitResult* fitTo(RooAbsData& data, Option_t *fitOpt, Option_t *optOpt, const char* fitRange)
RooFitResult* fitTo(RooAbsData& data, const RooArgSet& projDeps, Option_t *fitOpt, Option_t *optOpt, const char* fitRange)
 Fit this PDF to given data set

 OLD STYLE INTERFACE, PLEASE USE NEW INTERFACE fitTo(RooAbsData& data, RooCmdArg arg1,...,RooCmdArg arg8) 

 The dataset can be either binned, in which case a binned maximum likelihood fit
 is performed, or unbinned, in which case an unbinned maximum likelihood fit is performed

 Available fit options:

  "m" = MIGRAD only, i.e. no MINOS 
  "s" = estimate step size with HESSE before starting MIGRAD
  "h" = run HESSE after MIGRAD
  "e" = Perform extended MLL fit
  "0" = Run MIGRAD with strategy MINUIT 0 (no correlation matrix calculation at end)
        Does not apply to HESSE or MINOS, if run afterwards.
 
  "q" = Switch off verbose mode
  "l" = Save log file with parameter values at each MINUIT step
  "v" = Show changed parameters at each MINUIT step
  "t" = Time fit 
  "r" = Save fit output in RooFitResult object (return value is object RFR pointer)

 Available optimizer options

  "c" = Cache and precalculate components of PDF that exclusively depend on constant parameters
  "2" = Do NLL calculation in multi-processor mode on 2 processors
  "3" = Do NLL calculation in multi-processor mode on 3 processors
  "4" = Do NLL calculation in multi-processor mode on 4 processors

 The actual fit is performed to a temporary copy of both PDF and data set. Several optimization
 algorithm are run to increase the efficiency of the likelihood calculation and may increase
 the speed of complex fits up to an order of magnitude. All optimizations are exact, i.e the fit result
 of any fit should _exactly_ the same with and without optimization. We strongly encourage
 to stick to the default optimizer setting (all on). If for any reason you see a difference in the result
 with and without optimizer, please file a bug report.

 The function always return null unless the "r" fit option is specified. In that case a pointer to a RooFitResult
 is returned. The RooFitResult object contains the full fit output, including the correlation matrix.
void printToStream(ostream& os, PrintOption opt, TString indent)
 Print info about this object to the specified stream. In addition to the info
 from RooAbsArg::printToStream() we add:

     Shape : value, units, plot range
   Verbose : default binning and print label
RooAbsGenContext* genContext(const RooArgSet &vars, const RooDataSet *prototype, const RooArgSet* auxProto, Bool_t verbose)
RooDataSet * generate(const RooArgSet& whatVars, Int_t nEvents, const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3,const RooCmdArg& arg4, const RooCmdArg& arg5)
 Generate a new dataset containing the specified variables with events sampled from our distribution. 
 Generate the specified number of events or expectedEvents() if not specified.

 Any variables of this PDF that are not in whatVars will use their
 current values and be treated as fixed parameters. Returns zero
 in case of an error. The caller takes ownership of the returned
 dataset.

 The following named arguments are supported

 Verbose(Bool_t flag)               -- Print informational messages during event generation
 Extended()                         -- The actual number of events generated will be sampled from a Poisson distribution
                                       with mu=nevt. For use with extended maximum likelihood fits
 ProtoData(const RooDataSet& data,  -- Use specified dataset as prototype dataset. If randOrder is set to true
                 Bool_t randOrder)     the order of the events in the dataset will be read in a random order
                                       if the requested number of events to be generated does not match the
                                       number of events in the prototype dataset

 If ProtoData() is used, the specified existing dataset as a prototype: the new dataset will contain 
 the same number of events as the prototype (unless otherwise specified), and any prototype variables not in
 whatVars will be copied into the new dataset for each generated event and also used to set our PDF parameters. 
 The user can specify a  number of events to generate that will override the default. The result is a
 copy of the prototype dataset with only variables in whatVars randomized. Variables in whatVars that 
 are not in the prototype will be added as new columns to the generated dataset.  
RooDataSet * generate(const RooArgSet& whatVars, const RooCmdArg& arg1,const RooCmdArg& arg2, const RooCmdArg& arg3,const RooCmdArg& arg4, const RooCmdArg& arg5,const RooCmdArg& arg6)
 Generate a new dataset containing the specified variables with events sampled from our distribution. 
 Generate the specified number of events or expectedEvents() if not specified.

 Any variables of this PDF that are not in whatVars will use their
 current values and be treated as fixed parameters. Returns zero
 in case of an error. The caller takes ownership of the returned
 dataset.

 The following named arguments are supported

 Verbose(Bool_t flag)               -- Print informational messages during event generation
 NumEvent(int nevt)                 -- Generate specified number of events
 Extended()                         -- The actual number of events generated will be sampled from a Poisson distribution
                                       with mu=nevt. For use with extended maximum likelihood fits
 ProtoData(const RooDataSet& data,  -- Use specified dataset as prototype dataset. If randOrder is set to true
                 Bool_t randOrder,     the order of the events in the dataset will be read in a random order
                 Bool_t resample)      if the requested number of events to be generated does not match the
                                       number of events in the prototype dataset. If resample is also set to 
                                       true, the prototype dataset will be resampled rather than be strictly
                                       reshuffled. In this mode events of the protodata may be used more than
                                       once.

 If ProtoData() is used, the specified existing dataset as a prototype: the new dataset will contain 
 the same number of events as the prototype (unless otherwise specified), and any prototype variables not in
 whatVars will be copied into the new dataset for each generated event and also used to set our PDF parameters. 
 The user can specify a  number of events to generate that will override the default. The result is a
 copy of the prototype dataset with only variables in whatVars randomized. Variables in whatVars that 
 are not in the prototype will be added as new columns to the generated dataset.  
RooDataSet * generate(const RooArgSet &whatVars, Int_t nEvents, Bool_t verbose)
 Generate a new dataset containing the specified variables with
 events sampled from our distribution. Generate the specified
 number of events or else try to use expectedEvents() if nEvents <= 0.
 Any variables of this PDF that are not in whatVars will use their
 current values and be treated as fixed parameters. Returns zero
 in case of an error. The caller takes ownership of the returned
 dataset.
RooDataSet * generate(const RooArgSet &whatVars, const RooDataSet &prototype, Int_t nEvents, Bool_t verbose, Bool_t randProtoOrder, Bool_t resampleProto)
 Generate a new dataset with values of the whatVars variables
 sampled from our distribution. Use the specified existing dataset
 as a prototype: the new dataset will contain the same number of
 events as the prototype (by default), and any prototype variables not in
 whatVars will be copied into the new dataset for each generated
 event and also used to set our PDF parameters. The user can specify a
 number of events to generate that will override the default. The result is a
 copy of the prototype dataset with only variables in whatVars
 randomized. Variables in whatVars that are not in the prototype
 will be added as new columns to the generated dataset.  Returns
 zero in case of an error. The caller takes ownership of the
 returned dataset.
Int_t* randomizeProtoOrder(Int_t nProto, Int_t, Bool_t resampleProto)
 Return lookup table with randomized access order for prototype events,
 given nProto prototype data events and nGen events that will actually
 be accessed
Int_t getGenerator(const RooArgSet &/*directVars*/, RooArgSet &/*generatedVars*/, Bool_t /*staticInitOK*/)
 Load generatedVars with the subset of directVars that we can generate events for,
 and return a code that specifies the generator algorithm we will use. A code of
 zero indicates that we cannot generate any of the directVars (in this case, nothing
 should be added to generatedVars). Any non-zero codes will be passed to our generateEvent()
 implementation, but otherwise its value is arbitrary. The default implemetation of
 this method returns zero. Subclasses will usually implement this method using the
 matchArgs() methods to advertise the algorithms they provide.
void initGenerator(Int_t /*code*/)
 One-time initialization to setup the generator for the specified code.
void generateEvent(Int_t /*code*/)
 Generate an event using the algorithm corresponding to the specified code. The
 meaning of each code is defined by the getGenerator() implementation. The default
 implementation does nothing.
Bool_t isDirectGenSafe(const RooAbsArg& arg)
 Check if PDF depends via more than route on given arg
RooPlot* plotOn(RooPlot* frame, RooLinkedList& cmdList)
 Plot (project) PDF on specified frame. If a PDF is plotted in an empty frame, it
 will show a unit normalized curve in the frame variable, taken at the present value 
 of other observables defined for this PDF

 If a PDF is plotted in a frame in which a dataset has already been plotted, it will
 show a projected curve integrated over all variables that were present in the shown
 dataset except for the one on the x-axis. The normalization of the curve will also
 be adjusted to the event count of the plotted dataset. An informational message
 will be printed for each projection step that is performed

 This function takes the following named arguments

 Projection control

 Slice(const RooArgSet& set)     -- Override default projection behaviour by omittting observables listed 
                                    in set from the projection, resulting a 'slice' plot. Slicing is usually
                                    only sensible in discrete observables
 Project(const RooArgSet& set)   -- Override default projection behaviour by projecting over observables
                                    given in set and complete ignoring the default projection behavior. Advanced use only.
 ProjWData(const RooAbsData& d)  -- Override default projection _technique_ (integration). For observables present in given dataset
                                    projection of PDF is achieved by constructing an average over all observable values in given set.
                                    Consult RooFit plotting tutorial for further explanation of meaning & use of this technique
 ProjWData(const RooArgSet& s,   -- As above but only consider subset 's' of observables in dataset 'd' for projection through data averaging
           const RooAbsData& d)
 ProjectionRange(const char* rn) -- Override default range of projection integrals to a different range speficied by given range name.
                                    This technique allows you to project a finite width slice in a real-valued observable
 
 Misc content control

 Normalization(Double_t scale,   -- Adjust normalization by given scale factor. Interpretation of number depends on code: Relative:
                ScaleType code)     relative adjustment factor, NumEvent: scale to match given number of events.
 Name(const chat* name)          -- Give curve specified name in frame. Useful if curve is to be referenced later
 Asymmetry(const RooCategory& c) -- Show the asymmetry of the PDF in given two-state category [F(+)-F(-)] / [F(+)+F(-)] rather than
                                    the PDF projection. Category must have two states with indices -1 and +1 or three states with
                                    indeces -1,0 and +1.
 ShiftToZero(Bool_t flag)        -- Shift entire curve such that lowest visible point is at exactly zero. Mostly useful when
                                    plotting -log(L) or chi^2 distributions
 AddTo(const char* name,         -- Add constructed projection to already existing curve with given name and relative weight factors
       double_t wgtSelf, double_t wgtOther)

 Plotting control 

 LineStyle(Int_t style)          -- Select line style by ROOT line style code, default is solid
 LineColor(Int_t color)          -- Select line color by ROOT color code, default is blue
 LineWidth(Int_t width)          -- Select line with in pixels, default is 3
 FillStyle(Int_t style)          -- Select fill style, default is not filled. If a filled style is selected, also use VLines()
                                    to add vertical downward lines at end of curve to ensure proper closure
 FillColor(Int_t color)          -- Select fill color by ROOT color code
 Range(const char* name)         -- Only draw curve in range defined by given name
 Range(double lo, double hi)     -- Only draw curve in specified range
 VLines()                        -- Add vertical lines to y=0 at end points of curve
 Precision(Double_t eps)         -- Control precision of drawn curve w.r.t to scale of plot, default is 1e-3. Higher precision
                                    will result in more and more densely spaced curve points
 Invisble(Bool_t flag)           -- Add curve to frame, but do not display. Useful in combination AddTo()
cout << " plotOn(" << GetName() << ")
void plotOnCompSelect(RooArgSet* selNodes)
 Get complete set of tree branch nodes
RooPlot* plotOn(RooPlot *frame, PlotOpt o)
 Plot oneself on 'frame'. In addition to features detailed in  RooAbsReal::plotOn(),
 the scale factor for a PDF can be interpreted in three different ways. The interpretation
 is controlled by ScaleType

  Relative  -  Scale factor is applied on top of PDF normalization scale factor 
  NumEvent  -  Scale factor is interpreted as a number of events. The surface area
               under the PDF curve will match that of a histogram containing the specified
               number of event
  Raw       -  Scale factor is applied to the raw (projected) probability density.
               Not too useful, option provided for completeness.
RooPlot* plotCompOn(RooPlot *frame, const RooArgSet& compSet, Option_t* drawOptions, Double_t scaleFactor, ScaleType stype, const RooAbsData* projData, const RooArgSet* projSet)
 THIS FUNCTION IS OBSOLETE AND ONLY RETAINED FOR BACKWARD COMPATIBILITY. 
 PLEASE USE plotOn(frame,Componenents(...),...)

 Plot only the PDF components listed in 'compSet' of this PDF on 'frame'. 
 See RooAbsReal::plotOn() for a description of the remaining arguments and other features
RooPlot* plotCompOn(RooPlot *frame, const char* compNameList, Option_t* drawOptions, Double_t scaleFactor, ScaleType stype, const RooAbsData* projData, const RooArgSet* projSet)
 THIS FUNCTION IS OBSOLETE AND ONLY RETAINED FOR BACKWARD COMPATIBILITY. 
 PLEASE USE plotOn(frame,Componenents(...),...)

 Plot only the PDF components listed in 'compSet' of this PDF on 'frame'. 
 See RooAbsReal::plotOn() for a description of the remaining arguments and other features
RooPlot* plotCompOnEngine(RooPlot *frame, RooArgSet* selNodes, Option_t* drawOptions, Double_t scaleFactor, ScaleType stype, const RooAbsData* projData, const RooArgSet* projSet)
 Get complete set of tree branch nodes
RooPlot* plotCompSliceOn(RooPlot *frame, const char* compNameList, const RooArgSet& sliceSet, Option_t* drawOptions, Double_t scaleFactor, ScaleType stype, const RooAbsData* projData)
 THIS FUNCTION IS OBSOLETE AND ONLY RETAINED FOR BACKWARD COMPATIBILITY. 
 PLEASE USE plotOn(frame,Componenents(...),Slice(...),...)

 Plot ourselves on given frame, as done in plotOn(), except that the variables 
 listed in 'sliceSet' are taken out from the default list of projected dimensions created
 by plotOn().
RooPlot* plotCompSliceOn(RooPlot *frame, const RooArgSet& compSet, const RooArgSet& sliceSet, Option_t* drawOptions, Double_t scaleFactor, ScaleType stype, const RooAbsData* projData)
 THIS FUNCTION IS OBSOLETE AND ONLY RETAINED FOR BACKWARD COMPATIBILITY. 
 PLEASE USE plotOn(frame,Componenents(...),Slice(...),...)

 Plot ourselves on given frame, as done in plotOn(), except that the variables 
 listed in 'sliceSet' are taken out from the default list of projected dimensions created
 by plotOn().
RooPlot* paramOn(RooPlot* frame, const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4, const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8)
 Add a box with parameter values (and errors) to the specified frame

 The following named arguments are supported

   Parameters(const RooArgSet& param) -- Only the specified subset of parameters will be shown. 
                                         By default all non-contant parameters are shown
   ShowConstant(Bool_t flag)          -- Also display constant parameters
   Format(const char* optStr)         -- Classing [arameter formatting options, provided for backward compatibility
   Format(const char* what,...)       -- Parameter formatting options, details given below
   Label(const chat* label)           -- Add header label to parameter box
   Layout(Double_t xmin,              -- Specify relative position of left,right side of box and top of box. Position of 
       Double_t xmax, Double_t ymax)     bottom of box is calculated automatically from number lines in box


 The Format(const char* what,...) has the following structure

   const char* what      -- Controls what is shown. "N" adds name, "E" adds error, 
                            "A" shows asymmetric error, "U" shows unit, "H" hides the value
   FixedPrecision(int n) -- Controls precision, set fixed number of digits
   AutoPrecision(int n)  -- Controls precision. Number of shown digits is calculated from error 
                            + n specified additional digits (1 is sensible default)

 Example use: pdf.paramOn(frame, Label("fit result"), Format("NEU",AutoPrecision(1)) ) ;

RooPlot* paramOn(RooPlot* frame, const RooAbsData* data, const char *label, Int_t sigDigits, Option_t *options, Double_t xmin, Double_t xmax ,Double_t ymax)
 OBSOLETE FUNCTION PROVIDED FOR BACKWARD COMPATIBILITY
RooPlot* paramOn(RooPlot* frame, const RooArgSet& params, Bool_t showConstants, const char *label, Int_t sigDigits, Option_t *options, Double_t xmin, Double_t xmax ,Double_t ymax, const RooCmdArg* formatCmd)
 Add a text box with the current parameter values and their errors to the frame.
 Dependents of this PDF appearing in the 'data' dataset will be omitted.

 Optional label will be inserted as first line of the text box. Use 'sigDigits'
 to modify the default number of significant digits printed. The 'xmin,xmax,ymax'
 values specify the inital relative position of the text box in the plot frame  
void fixAddCoefNormalization(const RooArgSet& addNormSet)
void fixAddCoefRange(const char* rangeName)
RooPlot* plotNLLOn(RooPlot* frame, RooDataSet* data, Bool_t extended, const RooArgSet& /*projDeps*/, Option_t* /*drawOptions*/, Double_t prec, Bool_t fixMinToZero)
Bool_t redirectServersHook(const RooAbsCollection& newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t /*isRecursive*/)
Double_t expectedEvents(const RooArgSet*)
RooPlot * plotNLLOn(RooPlot* frame, RooDataSet* data, Option_t* drawOptions="L", Double_t prec=1e-2, Bool_t fixMinToZero=kTRUE)
RooPlot * plotNLLOn(RooPlot* frame, RooDataSet* data, Bool_t extended, Option_t* drawOptions="L", Double_t prec=1e-2, Bool_t fixMinToZero=kTRUE)
Double_t getNorm(const RooArgSet& nset)
{ return getNorm(&nset) ; }
Bool_t selfNormalized()
{ return kFALSE ; }
ExtendMode extendMode()
{ return CanNotBeExtended ; }
Bool_t canBeExtended()
{ return (extendMode() != CanNotBeExtended) ; }
Bool_t mustBeExtended()
{ return (extendMode() == MustBeExtended) ; }
Double_t expectedEvents(const RooArgSet* nset)
void verboseEval(Int_t stat)
{ _verboseEval = stat ; }
Bool_t isSelectedComp()
{ return _selectComp || _globalSelectComp ; }
void clearEvalError()
{ _evalError = kFALSE ; }
Bool_t evalError()
{ return _evalError ; }
void selectComp(Bool_t flag)
{ _selectComp = flag ; }
void globalSelectComp(Bool_t flag)
{ _globalSelectComp = flag ; }
void raiseEvalError()
{ _evalError = kTRUE ; }

Last update: Sat Dec 9 09:56:19 2006
Copyright (c) 2000-2005, Regents of the University of California *


ROOT page - Class index - Class Hierarchy - Top of the page

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.