// RooAbsReal is the common abstract base class for objects that represent a
// real value and implements functionality common to all real-valued objects
// such as the ability to plot them, to construct integrals of them, the
// ability to advertise (partial) analytical integrals etc..
// Implementation of RooAbsReal may be derived, thus no interface
// is provided to modify the contents.
//
//
// END_HTML
#include <sys/types.h>
#include "RooFit.h"
#include "RooMsgService.h"
#include "RooAbsReal.h"
#include "RooAbsReal.h"
#include "RooArgSet.h"
#include "RooArgList.h"
#include "RooPlot.h"
#include "RooCurve.h"
#include "RooRealVar.h"
#include "RooArgProxy.h"
#include "RooFormulaVar.h"
#include "RooRealBinding.h"
#include "RooRealIntegral.h"
#include "RooAbsCategoryLValue.h"
#include "RooCustomizer.h"
#include "RooAbsData.h"
#include "RooScaledFunc.h"
#include "RooAddPdf.h"
#include "RooCmdConfig.h"
#include "RooCategory.h"
#include "RooNumIntConfig.h"
#include "RooAddition.h"
#include "RooDataSet.h"
#include "RooDataHist.h"
#include "RooDataWeightedAverage.h"
#include "RooNumRunningInt.h"
#include "RooGlobalFunc.h"
#include "RooParamBinning.h"
#include "RooProfileLL.h"
#include "Riostream.h"
#include "TObjString.h"
#include "TTree.h"
#include "TH1.h"
#include "TH2.h"
#include "TH3.h"
#include "TBranch.h"
#include "TLeaf.h"
#include "TAttLine.h"
#include <sstream>
using namespace std ;
ClassImp(RooAbsReal)
;
Bool_t RooAbsReal::_cacheCheck(kFALSE) ;
Bool_t RooAbsReal::_globalSelectComp = kFALSE ;
Bool_t RooAbsReal::_doLogEvalError ;
map<const RooAbsArg*,pair<string,list<RooAbsReal::EvalError> > > RooAbsReal::_evalErrorList ;
RooAbsReal::RooAbsReal() : _specIntegratorConfig(0), _treeVar(kFALSE), _selectComp(kTRUE)
{
}
RooAbsReal::RooAbsReal(const char *name, const char *title, const char *unit) :
RooAbsArg(name,title), _plotMin(0), _plotMax(0), _plotBins(100),
_value(0), _unit(unit), _forceNumInt(kFALSE), _specIntegratorConfig(0), _treeVar(kFALSE), _selectComp(kTRUE)
{
setValueDirty() ;
setShapeDirty() ;
}
RooAbsReal::RooAbsReal(const char *name, const char *title, Double_t inMinVal,
Double_t inMaxVal, const char *unit) :
RooAbsArg(name,title), _plotMin(inMinVal), _plotMax(inMaxVal), _plotBins(100),
_value(0), _unit(unit), _forceNumInt(kFALSE), _specIntegratorConfig(0), _treeVar(kFALSE), _selectComp(kTRUE)
{
setValueDirty() ;
setShapeDirty() ;
}
RooAbsReal::RooAbsReal(const RooAbsReal& other, const char* name) :
RooAbsArg(other,name), _plotMin(other._plotMin), _plotMax(other._plotMax),
_plotBins(other._plotBins), _value(other._value), _unit(other._unit), _forceNumInt(other._forceNumInt),
_treeVar(other._treeVar), _selectComp(other._selectComp)
{
if (other._specIntegratorConfig) {
_specIntegratorConfig = new RooNumIntConfig(*other._specIntegratorConfig) ;
} else {
_specIntegratorConfig = 0 ;
}
}
RooAbsReal::~RooAbsReal()
{
if (_specIntegratorConfig) delete _specIntegratorConfig ;
}
Bool_t RooAbsReal::operator==(Double_t value) const
{
return (getVal()==value) ;
}
Bool_t RooAbsReal::operator==(const RooAbsArg& other)
{
const RooAbsReal* otherReal = dynamic_cast<const RooAbsReal*>(&other) ;
return otherReal ? operator==(otherReal->getVal()) : kFALSE ;
}
TString RooAbsReal::getTitle(Bool_t appendUnit) const
{
TString title(GetTitle());
if(appendUnit && 0 != strlen(getUnit())) {
title.Append(" (");
title.Append(getUnit());
title.Append(")");
}
return title;
}
Double_t RooAbsReal::getVal(const RooArgSet* set) const
{
if (isValueDirty() || isShapeDirty()) {
_value = traceEval(set) ;
clearValueDirty() ;
clearShapeDirty() ;
} else if (_cacheCheck) {
Double_t checkValue = traceEval(set);
if (checkValue != _value) {
coutW(Eval) << "RooAbsReal::getVal(" << GetName() << ") WARNING: cache contains " << _value
<< " but evaluate() returns " << checkValue << endl ;
_value = checkValue ;
}
}
return _value ;
}
Double_t RooAbsReal::traceEval(const RooArgSet* ) const
{
Double_t value = evaluate() ;
cxcoutD(Tracing) << "RooAbsReal::getVal(" << GetName() << ") operMode = " << _operMode << " recalculated, new value = " << value << endl ;
if (!isValidReal(value)) {
coutW(Tracing) << "RooAbsReal::traceEval(" << GetName()
<< "): validation failed: " << value << endl ;
}
traceEvalHook(value) ;
return value ;
}
Int_t RooAbsReal::getAnalyticalIntegralWN(RooArgSet& allDeps, RooArgSet& analDeps,
const RooArgSet* , const char* rangeName) const
{
return _forceNumInt ? 0 : getAnalyticalIntegral(allDeps,analDeps,rangeName) ;
}
Int_t RooAbsReal::getAnalyticalIntegral(RooArgSet& , RooArgSet& , const char* ) const
{
return 0 ;
}
Double_t RooAbsReal::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const
{
if (code==0) return getVal(normSet) ;
return analyticalIntegral(code,rangeName) ;
}
Double_t RooAbsReal::analyticalIntegral(Int_t code, const char* ) const
{
coutF(Eval) << "RooAbsReal::analyticalIntegral(" << GetName() << ") code " << code << " not implemented" << endl ;
return 0 ;
}
const char *RooAbsReal::getPlotLabel() const
{
return _label.IsNull() ? fName.Data() : _label.Data();
}
void RooAbsReal::setPlotLabel(const char *label)
{
_label= label;
}
Bool_t RooAbsReal::readFromStream(istream& , Bool_t , Bool_t )
{
return kFALSE ;
}
void RooAbsReal::writeToStream(ostream& , Bool_t ) const
{
}
void RooAbsReal::printValue(ostream& os) const
{
os << getVal() ;
}
void RooAbsReal::printMultiline(ostream& os, Int_t contents, Bool_t verbose, TString indent) const
{
RooAbsArg::printMultiline(os,contents,verbose,indent) ;
os << indent << "--- RooAbsReal ---" << endl;
TString unit(_unit);
if(!unit.IsNull()) unit.Prepend(' ');
os << indent << " Value = " << getVal() << unit << endl;
os << endl << indent << " Plot label is \"" << getPlotLabel() << "\"" << endl;
}
Bool_t RooAbsReal::isValid() const
{
return isValidReal(_value) ;
}
Bool_t RooAbsReal::isValidReal(Double_t , Bool_t ) const
{
return kTRUE ;
}
RooAbsReal* RooAbsReal::createProfile(const RooArgSet& paramsOfInterest)
{
TString name(Form("%s_Profile[",GetName())) ;
TIterator* iter = paramsOfInterest.createIterator() ;
RooAbsArg* arg ;
Bool_t first(kTRUE) ;
while((arg=(RooAbsArg*)iter->Next())) {
if (first) {
first=kFALSE ;
} else {
name.Append(",") ;
}
name.Append(arg->GetName()) ;
}
delete iter ;
name.Append("]") ;
return new RooProfileLL(name.Data(),Form("Profile of %s",GetTitle()),*this,paramsOfInterest) ;
}
RooAbsReal* RooAbsReal::createIntegral(const RooArgSet& iset, const RooCmdArg arg1, const RooCmdArg arg2,
const RooCmdArg arg3, const RooCmdArg arg4, const RooCmdArg arg5,
const RooCmdArg arg6, const RooCmdArg arg7, const RooCmdArg arg8) const
{
RooCmdConfig pc(Form("RooAbsReal::createIntegral(%s)",GetName())) ;
pc.defineString("rangeName","RangeWithName",0,"",kTRUE) ;
pc.defineObject("normSet","NormSet",0,0) ;
pc.defineObject("numIntConfig","NumIntConfig",0,0) ;
pc.process(arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8) ;
if (!pc.ok(kTRUE)) {
return 0 ;
}
const char* rangeName = pc.getString("rangeName",0,kTRUE) ;
const RooArgSet* nset = static_cast<const RooArgSet*>(pc.getObject("normSet",0)) ;
const RooNumIntConfig* cfg = static_cast<const RooNumIntConfig*>(pc.getObject("numIntConfig",0)) ;
return createIntegral(iset,nset,cfg,rangeName) ;
}
RooAbsReal* RooAbsReal::createIntegral(const RooArgSet& iset, const RooArgSet* nset,
const RooNumIntConfig* cfg, const char* rangeName) const
{
if (!rangeName || strchr(rangeName,',')==0) {
return createIntObj(iset,nset,cfg,rangeName) ;
}
char* buf = new char[strlen(rangeName)+1] ;
strcpy(buf,rangeName) ;
char* range = strtok(buf,",") ;
RooArgSet components ;
while (range) {
RooAbsReal* compIntegral = createIntObj(iset,nset,cfg,range) ;
components.add(*compIntegral) ;
range = strtok(0,",") ;
}
delete[] buf ;
TString title(GetTitle()) ;
title.Prepend("Integral of ") ;
TString fullName(GetName()) ;
fullName.Append(integralNameSuffix(iset,nset,rangeName)) ;
return new RooAddition(fullName.Data(),title.Data(),components,kTRUE) ;
}
RooAbsReal* RooAbsReal::createIntObj(const RooArgSet& iset2, const RooArgSet* nset2,
const RooNumIntConfig* cfg, const char* rangeName) const
{
RooArgSet iset(iset2) ;
const RooArgSet* nset = nset2 ;
Bool_t error = kFALSE ;
const RooAbsReal* integrand = this ;
RooAbsReal* integral = 0 ;
if (iset.getSize()==0) {
TString title(GetTitle()) ;
title.Prepend("Integral of ") ;
TString name(GetName()) ;
name.Append(integralNameSuffix(iset,nset,rangeName)) ;
return new RooRealIntegral(name,title,*this,iset,nset,cfg,rangeName) ;
}
while(iset.getSize()>0) {
RooArgSet innerSet ;
findInnerMostIntegration(iset,innerSet,rangeName) ;
if (innerSet.getSize()==0) {
error = kTRUE ;
break ;
}
TString title(integrand->GetTitle()) ;
title.Prepend("Integral of ") ;
TString name(integrand->GetName()) ;
name.Append(integrand->integralNameSuffix(innerSet,nset,rangeName)) ;
integral = new RooRealIntegral(name,title,*integrand,innerSet,nset,cfg,rangeName) ;
if (integrand != this) {
integral->addOwnedComponents(*integrand) ;
}
iset.remove(innerSet) ;
if (integrand == this && iset.getSize()>0) {
coutI(Integration) << GetName() << " : multidimensional integration over observables with parameterized ranges in terms of other integrated observables detected, using recursive integration strategy to construct final integral" << endl ;
}
integrand = integral ;
nset = 0 ;
}
if (error) {
coutE(Integration) << GetName() << " : ERROR while defining recursive integral over observables with parameterized integration ranges, please check that integration rangs specify uniquely defined integral " << endl;
delete integral ;
integral = 0 ;
}
return integral ;
}
void RooAbsReal::findInnerMostIntegration(const RooArgSet& allObs, RooArgSet& innerObs, const char* rangeName) const
{
RooArgSet obsWithFixedRange(allObs) ;
RooArgSet obsWithParamRange ;
RooArgSet obsServingAsRangeParams ;
TIterator* oiter = allObs.createIterator() ;
RooAbsArg* aarg ;
while((aarg=(RooAbsArg*)oiter->Next())) {
RooAbsRealLValue* arglv = dynamic_cast<RooAbsRealLValue*>(aarg) ;
if (arglv) {
RooAbsBinning& binning = arglv->getBinning(rangeName,kTRUE,kTRUE) ;
if (binning.isParameterized()) {
RooArgSet* loBoundObs = binning.lowBoundFunc()->getObservables(allObs) ;
RooArgSet* hiBoundObs = binning.highBoundFunc()->getObservables(allObs) ;
if (loBoundObs->overlaps(allObs) || hiBoundObs->overlaps(allObs)) {
obsWithParamRange.add(*aarg) ;
obsWithFixedRange.remove(*aarg) ;
obsServingAsRangeParams.add(*loBoundObs,kFALSE) ;
obsServingAsRangeParams.add(*hiBoundObs,kFALSE) ;
}
delete loBoundObs ;
delete hiBoundObs ;
}
}
}
delete oiter ;
RooArgSet obsWithFixedRangeNP(obsWithFixedRange) ;
obsWithFixedRangeNP.remove(obsServingAsRangeParams) ;
RooArgSet obsWithParamRangeNP(obsWithParamRange) ;
obsWithParamRangeNP.remove(obsServingAsRangeParams) ;
innerObs.removeAll() ;
innerObs.add(obsWithFixedRangeNP) ;
innerObs.add(obsWithParamRangeNP) ;
}
TString RooAbsReal::integralNameSuffix(const RooArgSet& iset, const RooArgSet* nset, const char* rangeName, Bool_t omitEmpty) const
{
TString name ;
if (iset.getSize()>0) {
name.Append("_Int[") ;
TIterator* iter = iset.createIterator() ;
RooAbsArg* arg ;
Bool_t first(kTRUE) ;
while((arg=(RooAbsArg*)iter->Next())) {
if (first) {
first=kFALSE ;
} else {
name.Append(",") ;
}
name.Append(arg->GetName()) ;
}
delete iter ;
if (rangeName) {
name.Append("|") ;
name.Append(rangeName) ;
}
name.Append("]");
} else if (!omitEmpty) {
name.Append("_Int[]") ;
}
if (nset && nset->getSize()>0 ) {
name.Append("_Norm[") ;
Bool_t first(kTRUE);
TIterator* iter = nset->createIterator() ;
RooAbsArg* arg ;
while((arg=(RooAbsArg*)iter->Next())) {
if (first) {
first=kFALSE ;
} else {
name.Append(",") ;
}
name.Append(arg->GetName()) ;
}
delete iter ;
name.Append("]") ;
}
return name ;
}
const RooAbsReal* RooAbsReal::createPlotProjection(const RooArgSet& depVars, const RooArgSet& projVars,
RooArgSet*& cloneSet) const
{
return createPlotProjection(depVars,&projVars,cloneSet) ;
}
const RooAbsReal* RooAbsReal::createPlotProjection(const RooArgSet& depVars, const RooArgSet& projVars) const
{
RooArgSet* cloneSet = new RooArgSet() ;
return createPlotProjection(depVars,&projVars,cloneSet) ;
}
const RooAbsReal *RooAbsReal::createPlotProjection(const RooArgSet &dependentVars, const RooArgSet *projectedVars,
RooArgSet *&cloneSet, const char* rangeName) const
{
RooArgSet leafNodes;
RooArgSet treeNodes;
leafNodeServerList(&leafNodes,this);
treeNodeServerList(&treeNodes,this) ;
TIterator *dependentIterator= dependentVars.createIterator();
assert(0 != dependentIterator);
const RooAbsArg *arg = 0;
while((arg= (const RooAbsArg*)dependentIterator->Next())) {
if(!arg->isFundamental() && !dynamic_cast<const RooAbsLValue*>(arg)) {
coutE(Plotting) << ClassName() << "::" << GetName() << ":createPlotProjection: variable \"" << arg->GetName()
<< "\" of wrong type: " << arg->ClassName() << endl;
delete dependentIterator;
return 0;
}
RooAbsArg *found= treeNodes.find(arg->GetName());
if(!found) {
coutE(Plotting) << ClassName() << "::" << GetName() << ":createPlotProjection: \"" << arg->GetName()
<< "\" is not a dependent and will be ignored." << endl;
continue;
}
if(found != arg) {
if (leafNodes.find(found->GetName())) {
leafNodes.replace(*found,*arg);
} else {
leafNodes.add(*arg) ;
RooArgSet* lvDep = arg->getObservables(&leafNodes) ;
RooAbsArg* lvs ;
TIterator* iter = lvDep->createIterator() ;
while((lvs=(RooAbsArg*)iter->Next())) {
RooAbsArg* tmp = leafNodes.find(lvs->GetName()) ;
if (tmp) {
leafNodes.remove(*tmp) ;
leafNodes.add(*lvs) ;
}
}
delete iter ;
}
}
if(0 != projectedVars && projectedVars->find(arg->GetName())) {
coutE(Plotting) << ClassName() << "::" << GetName() << ":createPlotProjection: \"" << arg->GetName()
<< "\" cannot be both a dependent and a projected variable." << endl;
delete dependentIterator;
return 0;
}
}
if(0 != projectedVars) leafNodes.remove(*projectedVars,kTRUE);
cloneSet= (RooArgSet*)RooArgSet(*this).snapshot(kTRUE);
if (!cloneSet) {
coutE(Plotting) << "RooAbsPdf::createPlotProjection(" << GetName() << ") Couldn't deep-clone PDF, abort," << endl ;
return 0 ;
}
RooAbsReal *theClone= (RooAbsReal*)cloneSet->find(GetName());
RooArgSet* plotLeafNodes = (RooArgSet*) leafNodes.selectCommon(dependentVars) ;
theClone->recursiveRedirectServers(*plotLeafNodes,kFALSE,kFALSE,kFALSE);
delete plotLeafNodes ;
RooArgSet normSet(dependentVars);
if(0 != projectedVars) normSet.add(*projectedVars);
RooArgSet empty;
if(0 == projectedVars) projectedVars= ∅
TString name = GetName() ;
name += integralNameSuffix(*projectedVars,&normSet,rangeName,kTRUE) ;
TString title(GetTitle());
title.Prepend("Projection of ");
RooRealIntegral *projected= new RooRealIntegral(name.Data(),title.Data(),*theClone,*projectedVars,&normSet,0,rangeName);
if(0 == projected || !projected->isValid()) {
coutE(Plotting) << ClassName() << "::" << GetName() << ":createPlotProjection: cannot integrate out ";
projectedVars->printStream(cout,kName|kArgs,kSingleLine);
if(0 != projected) delete projected;
delete dependentIterator;
return 0;
}
cloneSet->addOwned(*projected);
delete dependentIterator;
return projected;
}
TH1 *RooAbsReal::fillHistogram(TH1 *hist, const RooArgList &plotVars,
Double_t scaleFactor, const RooArgSet *projectedVars, Bool_t scaleForDensity) const
{
if(0 == hist) {
coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillHistogram: no valid histogram to fill" << endl;
return 0;
}
Int_t hdim= hist->GetDimension();
if(hdim != plotVars.getSize()) {
coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillHistogram: plotVars has the wrong dimension" << endl;
return 0;
}
RooArgSet plotClones;
for(Int_t index= 0; index < plotVars.getSize(); index++) {
const RooAbsArg *var= plotVars.at(index);
const RooRealVar *realVar= dynamic_cast<const RooRealVar*>(var);
if(0 == realVar) {
coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillHistogram: cannot plot variable \"" << var->GetName()
<< "\" of type " << var->ClassName() << endl;
return 0;
}
if(!this->dependsOn(*realVar)) {
coutE(InputArguments) << ClassName() << "::" << GetName()
<< ":fillHistogram: WARNING: variable is not an explicit dependent: " << realVar->GetName() << endl;
}
plotClones.addClone(*realVar,kTRUE);
}
TIterator* pciter= plotClones.createIterator() ;
RooAbsArg* pc ;
while((pc=(RooAbsArg*)pciter->Next())) {
pc->recursiveRedirectServers(plotClones,kFALSE,kFALSE,kTRUE) ;
}
delete pciter ;
RooArgSet allDeps(plotClones) ;
if (projectedVars) allDeps.add(*projectedVars) ;
if (checkObservables(&allDeps)) {
coutE(InputArguments) << "RooAbsReal::fillHistogram(" << GetName() << ") error in checkObservables, abort" << endl ;
return hist ;
}
RooArgSet *cloneSet = 0;
const RooAbsReal *projected= createPlotProjection(plotClones,projectedVars,cloneSet);
cxcoutD(Plotting) << "RooAbsReal::fillHistogram(" << GetName() << ") plot projection object is " << projected->GetName() << endl ;
Int_t xbins(0),ybins(1),zbins(1);
RooRealVar *xvar = 0;
RooRealVar *yvar = 0;
RooRealVar *zvar = 0;
TAxis *xaxis = 0;
TAxis *yaxis = 0;
TAxis *zaxis = 0;
switch(hdim) {
case 3:
zbins= hist->GetNbinsZ();
zvar= dynamic_cast<RooRealVar*>(plotClones.find(plotVars.at(2)->GetName()));
zaxis= hist->GetZaxis();
assert(0 != zvar && 0 != zaxis);
if (scaleForDensity) {
scaleFactor*= (zaxis->GetXmax() - zaxis->GetXmin())/zbins;
}
case 2:
ybins= hist->GetNbinsY();
yvar= dynamic_cast<RooRealVar*>(plotClones.find(plotVars.at(1)->GetName()));
yaxis= hist->GetYaxis();
assert(0 != yvar && 0 != yaxis);
if (scaleForDensity) {
scaleFactor*= (yaxis->GetXmax() - yaxis->GetXmin())/ybins;
}
case 1:
xbins= hist->GetNbinsX();
xvar= dynamic_cast<RooRealVar*>(plotClones.find(plotVars.at(0)->GetName()));
xaxis= hist->GetXaxis();
assert(0 != xvar && 0 != xaxis);
if (scaleForDensity) {
scaleFactor*= (xaxis->GetXmax() - xaxis->GetXmin())/xbins;
}
break;
default:
coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillHistogram: cannot fill histogram with "
<< hdim << " dimensions" << endl;
break;
}
Int_t xbin(0),ybin(0),zbin(0);
Int_t bins= xbins*ybins*zbins;
for(Int_t bin= 0; bin < bins; bin++) {
switch(hdim) {
case 3:
if(bin % (xbins*ybins) == 0) {
zbin++;
zvar->setVal(zaxis->GetBinCenter(zbin));
}
case 2:
if(bin % xbins == 0) {
ybin= (ybin%ybins) + 1;
yvar->setVal(yaxis->GetBinCenter(ybin));
}
case 1:
xbin= (xbin%xbins) + 1;
xvar->setVal(xaxis->GetBinCenter(xbin));
break;
default:
coutE(InputArguments) << "RooAbsReal::fillHistogram: Internal Error!" << endl;
break;
}
Double_t result= scaleFactor*projected->getVal();
hist->SetBinContent(hist->GetBin(xbin,ybin,zbin),result);
}
delete cloneSet;
return hist;
}
RooDataHist* RooAbsReal::fillDataHist(RooDataHist *hist, const RooArgSet* normSet,
Double_t scaleFactor, Bool_t correctForBinSize, Bool_t showProgress) const
{
if(0 == hist) {
coutE(InputArguments) << ClassName() << "::" << GetName() << ":fillDataHist: no valid RooDataHist to fill" << endl;
return 0;
}
RooArgSet allDeps(*hist->get()) ;
if (checkObservables(&allDeps)) {
coutE(InputArguments) << "RooAbsReal::fillDataHist(" << GetName() << ") error in checkObservables, abort" << endl ;
return hist ;
}
RooArgSet* cloneSet = (RooArgSet*) RooArgSet(*this).snapshot(kTRUE) ;
RooAbsReal* theClone = (RooAbsReal*) cloneSet->find(GetName()) ;
theClone->attachDataSet(*hist) ;
Int_t onePct = hist->numEntries()/100 ;
if (onePct==0) {
onePct++ ;
}
for (Int_t i=0 ; i<hist->numEntries() ; i++) {
if (showProgress && (i%onePct==0)) {
ccoutP(Eval) << "." << flush ;
}
const RooArgSet* obs = hist->get(i) ;
Double_t binVal = theClone->getVal(normSet?normSet:obs)*scaleFactor ;
if (correctForBinSize) binVal*= hist->binVolume() ;
hist->set(binVal) ;
}
delete cloneSet ;
return hist;
}
TH1* RooAbsReal::createHistogram(const char* varNameList, Int_t xbins, Int_t ybins, Int_t zbins) const
{
char buf[1024] ;
strcpy(buf,varNameList) ;
char* varName = strtok(buf,",:") ;
RooArgSet* vars = getVariables() ;
RooRealVar* xvar = (RooRealVar*) vars->find(varName) ;
varName = strtok(0,",") ;
RooRealVar* yvar = varName ? (RooRealVar*) vars->find(varName) : 0 ;
varName = strtok(0,",") ;
RooRealVar* zvar = varName ? (RooRealVar*) vars->find(varName) : 0 ;
delete vars ;
RooLinkedList argList ;
if (xbins>0) {
argList.Add(RooFit::Binning(xbins).Clone()) ;
}
if (yvar) {
if (ybins>0) {
argList.Add(RooFit::YVar(*yvar,RooFit::Binning(ybins)).Clone()) ;
} else {
argList.Add(RooFit::YVar(*yvar).Clone()) ;
}
}
if (zvar) {
if (zbins>0) {
argList.Add(RooFit::ZVar(*zvar,RooFit::Binning(zbins)).Clone()) ;
} else {
argList.Add(RooFit::ZVar(*zvar).Clone()) ;
}
}
TH1* result = createHistogram(GetName(),*xvar,argList) ;
argList.Delete() ;
return result ;
}
TH1 *RooAbsReal::createHistogram(const char *name, const RooAbsRealLValue& xvar,
const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4,
const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8) const
{
RooLinkedList l ;
l.Add((TObject*)&arg1) ; l.Add((TObject*)&arg2) ;
l.Add((TObject*)&arg3) ; l.Add((TObject*)&arg4) ;
l.Add((TObject*)&arg5) ; l.Add((TObject*)&arg6) ;
l.Add((TObject*)&arg7) ; l.Add((TObject*)&arg8) ;
return createHistogram(name,xvar,l) ;
}
TH1* RooAbsReal::createHistogram(const char *name, const RooAbsRealLValue& xvar, RooLinkedList& argList) const
{
RooCmdConfig pc(Form("RooAbsReal::createHistogram(%s)",GetName())) ;
pc.defineInt("scaling","Scaling",0,1) ;
pc.defineObject("projObs","ProjectedObservables",0,0) ;
pc.defineObject("yvar","YVar",0,0) ;
pc.defineObject("zvar","ZVar",0,0) ;
pc.allowUndefined() ;
pc.process(argList) ;
if (!pc.ok(kTRUE)) {
return 0 ;
}
RooArgList vars(xvar) ;
RooAbsArg* yvar = static_cast<RooAbsArg*>(pc.getObject("yvar")) ;
if (yvar) {
vars.add(*yvar) ;
}
RooAbsArg* zvar = static_cast<RooAbsArg*>(pc.getObject("zvar")) ;
if (zvar) {
vars.add(*zvar) ;
}
RooArgSet* projObs = static_cast<RooArgSet*>(pc.getObject("projObs")) ;
Bool_t doScaling = pc.getInt("scaling") ;
RooLinkedList l2 ;
for (Int_t i=0 ; i<argList.GetSize() ; i++) {
if (TString(argList.At(i)->GetName()).CompareTo("Scaling")) {
l2.Add(argList.At(i)) ;
}
}
TH1* histo = xvar.createHistogram(name,l2) ;
fillHistogram(histo,vars,1.0,projObs,doScaling) ;
return histo ;
}
RooPlot* RooAbsReal::plotOn(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,
const RooCmdArg& arg9, const RooCmdArg& arg10) const
{
RooLinkedList l ;
l.Add((TObject*)&arg1) ; l.Add((TObject*)&arg2) ;
l.Add((TObject*)&arg3) ; l.Add((TObject*)&arg4) ;
l.Add((TObject*)&arg5) ; l.Add((TObject*)&arg6) ;
l.Add((TObject*)&arg7) ; l.Add((TObject*)&arg8) ;
l.Add((TObject*)&arg9) ; l.Add((TObject*)&arg10) ;
return plotOn(frame,l) ;
}
RooPlot* RooAbsReal::plotOn(RooPlot* frame, RooLinkedList& argList) const
{
RooCmdArg* rcmd = (RooCmdArg*) argList.FindObject("RangeWithName") ;
if (rcmd && TString(rcmd->getString(0)).Contains(",")) {
RooCmdArg rnorm = RooFit::NormRange(rcmd->getString(0)) ;
argList.Add(&rnorm) ;
list<string> rlist ;
char buf[1024] ;
strcpy(buf,rcmd->getString(0)) ;
char* oneRange = strtok(buf,",") ;
while(oneRange) {
rlist.push_back(oneRange) ;
oneRange = strtok(0,",") ;
}
for (list<string>::iterator riter=rlist.begin() ; riter!=rlist.end() ; ++riter) {
rcmd->setString(0,riter->c_str()) ;
RooAbsReal::plotOn(frame,argList) ;
}
return frame ;
}
RooCmdConfig pc(Form("RooAbsReal::plotOn(%s)",GetName())) ;
pc.defineString("drawOption","DrawOption",0,"L") ;
pc.defineString("projectionRangeName","ProjectionRange",0,"",kTRUE) ;
pc.defineString("curveNameSuffix","CurveNameSuffix",0,"") ;
pc.defineString("sliceCatState","SliceCat",0,"",kTRUE) ;
pc.defineDouble("scaleFactor","Normalization",0,1.0) ;
pc.defineObject("sliceSet","SliceVars",0) ;
pc.defineObject("sliceCatList","SliceCat",0,0,kTRUE) ;
pc.defineObject("projSet","Project",0) ;
pc.defineObject("asymCat","Asymmetry",0) ;
pc.defineDouble("precision","Precision",0,1e-3) ;
pc.defineDouble("evalErrorVal","EvalErrorValue",0,0) ;
pc.defineInt("doEvalError","EvalErrorValue",0,0) ;
pc.defineInt("shiftToZero","ShiftToZero",0,0) ;
pc.defineObject("projDataSet","ProjData",0) ;
pc.defineObject("projData","ProjData",1) ;
pc.defineInt("binProjData","ProjData",0,0) ;
pc.defineDouble("rangeLo","Range",0,-999.) ;
pc.defineDouble("rangeHi","Range",1,-999.) ;
pc.defineInt("numee","PrintEvalErrors",0,10) ;
pc.defineInt("rangeAdjustNorm","Range",0,0) ;
pc.defineInt("rangeWNAdjustNorm","RangeWithName",0,0) ;
pc.defineInt("VLines","VLines",0,2) ;
pc.defineString("rangeName","RangeWithName",0,"") ;
pc.defineString("normRangeName","NormRange",0,"") ;
pc.defineInt("lineColor","LineColor",0,-999) ;
pc.defineInt("lineStyle","LineStyle",0,-999) ;
pc.defineInt("lineWidth","LineWidth",0,-999) ;
pc.defineInt("fillColor","FillColor",0,-999) ;
pc.defineInt("fillStyle","FillStyle",0,-999) ;
pc.defineString("curveName","Name",0,"") ;
pc.defineInt("curveInvisible","Invisible",0,0) ;
pc.defineInt("numCPU","NumCPU",0,1) ;
pc.defineInt("interleave","NumCPU",1,0) ;
pc.defineString("addToCurveName","AddTo",0,"") ;
pc.defineDouble("addToWgtSelf","AddTo",0,1.) ;
pc.defineDouble("addToWgtOther","AddTo",1,1.) ;
pc.defineInt("moveToBack","MoveToBack",0,0) ;
pc.defineMutex("SliceVars","Project") ;
pc.defineMutex("AddTo","Asymmetry") ;
pc.defineMutex("Range","RangeWithName") ;
pc.process(argList) ;
if (!pc.ok(kTRUE)) {
return frame ;
}
PlotOpt o ;
o.numee = pc.getInt("numee") ;
o.drawOptions = pc.getString("drawOption") ;
o.curveNameSuffix = pc.getString("curveNameSuffix") ;
o.scaleFactor = pc.getDouble("scaleFactor") ;
o.projData = (const RooAbsData*) pc.getObject("projData") ;
o.binProjData = pc.getInt("binProjData") ;
o.projDataSet = (const RooArgSet*) pc.getObject("projDataSet") ;
o.numCPU = pc.getInt("numCPU") ;
o.interleave = pc.getInt("interleave") ;
o.eeval = pc.getDouble("evalErrorVal") ;
o.doeeval = pc.getInt("doEvalError") ;
const RooArgSet* sliceSetTmp = (const RooArgSet*) pc.getObject("sliceSet") ;
RooArgSet* sliceSet = sliceSetTmp ? ((RooArgSet*) sliceSetTmp->Clone()) : 0 ;
const RooArgSet* projSet = (const RooArgSet*) pc.getObject("projSet") ;
const RooAbsCategoryLValue* asymCat = (const RooAbsCategoryLValue*) pc.getObject("asymCat") ;
const char* sliceCatState = pc.getString("sliceCatState",0,kTRUE) ;
const RooLinkedList& sliceCatList = pc.getObjectList("sliceCatList") ;
if (sliceCatState) {
if (!sliceSet) {
sliceSet = new RooArgSet ;
}
char buf[1024] ;
strcpy(buf,sliceCatState) ;
const char* slabel = strtok(buf,",") ;
TIterator* iter = sliceCatList.MakeIterator() ;
RooCategory* scat ;
while((scat=(RooCategory*)iter->Next())) {
if (slabel) {
scat->setLabel(slabel) ;
sliceSet->add(*scat,kFALSE) ;
}
slabel = strtok(0,",") ;
}
delete iter ;
}
o.precision = pc.getDouble("precision") ;
o.shiftToZero = (pc.getInt("shiftToZero")!=0) ;
Int_t vlines = pc.getInt("VLines");
if (pc.hasProcessed("Range")) {
o.rangeLo = pc.getDouble("rangeLo") ;
o.rangeHi = pc.getDouble("rangeHi") ;
o.postRangeFracScale = pc.getInt("rangeAdjustNorm") ;
if (vlines==2) vlines=0 ;
} else if (pc.hasProcessed("RangeWithName")) {
o.normRangeName = pc.getString("rangeName",0,kTRUE) ;
o.rangeLo = frame->getPlotVar()->getMin(pc.getString("rangeName",0,kTRUE)) ;
o.rangeHi = frame->getPlotVar()->getMax(pc.getString("rangeName",0,kTRUE)) ;
o.postRangeFracScale = pc.getInt("rangeWNAdjustNorm") ;
if (vlines==2) vlines=0 ;
}
if (pc.hasProcessed("NormRange")) {
o.normRangeName = pc.getString("normRangeName") ;
o.postRangeFracScale = kTRUE ;
}
o.wmode = (vlines==2)?RooCurve::Extended:(vlines==1?RooCurve::Straight:RooCurve::NoWings) ;
o.projectionRangeName = pc.getString("projectionRangeName",0,kTRUE) ;
o.curveName = pc.getString("curveName",0,kTRUE) ;
o.curveInvisible = pc.getInt("curveInvisible") ;
o.addToCurveName = pc.getString("addToCurveName",0,kTRUE) ;
o.addToWgtSelf = pc.getDouble("addToWgtSelf") ;
o.addToWgtOther = pc.getDouble("addToWgtOther") ;
if (o.addToCurveName && !frame->findObject(o.addToCurveName,RooCurve::Class())) {
coutE(InputArguments) << "RooAbsReal::plotOn(" << GetName() << ") cannot find existing curve " << o.addToCurveName << " to add to in RooPlot" << endl ;
return frame ;
}
RooArgSet projectedVars ;
if (sliceSet) {
cxcoutD(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") Preprocessing: have slice " << *sliceSet << endl ;
makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,kTRUE) ;
TIterator* iter = sliceSet->createIterator() ;
RooAbsArg* sliceArg ;
while((sliceArg=(RooAbsArg*)iter->Next())) {
RooAbsArg* arg = projectedVars.find(sliceArg->GetName()) ;
if (arg) {
projectedVars.remove(*arg) ;
} else {
coutI(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") slice variable "
<< sliceArg->GetName() << " was not projected anyway" << endl ;
}
}
delete iter ;
} else if (projSet) {
cxcoutD(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") Preprocessing: have projSet " << *projSet << endl ;
makeProjectionSet(frame->getPlotVar(),projSet,projectedVars,kFALSE) ;
} else {
cxcoutD(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") Preprocessing: have neither sliceSet nor projSet " << endl ;
makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,kTRUE) ;
}
o.projSet = &projectedVars ;
cxcoutD(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") Preprocessing: projectedVars = " << projectedVars << endl ;
RooPlot* ret ;
if (!asymCat) {
ret = RooAbsReal::plotOn(frame,o) ;
} else {
ret = RooAbsReal::plotAsymOn(frame,*asymCat,o) ;
}
delete sliceSet ;
Int_t lineColor = pc.getInt("lineColor") ;
Int_t lineStyle = pc.getInt("lineStyle") ;
Int_t lineWidth = pc.getInt("lineWidth") ;
Int_t fillColor = pc.getInt("fillColor") ;
Int_t fillStyle = pc.getInt("fillStyle") ;
if (lineColor!=-999) ret->getAttLine()->SetLineColor(lineColor) ;
if (lineStyle!=-999) ret->getAttLine()->SetLineStyle(lineStyle) ;
if (lineWidth!=-999) ret->getAttLine()->SetLineWidth(lineWidth) ;
if (fillColor!=-999) ret->getAttFill()->SetFillColor(fillColor) ;
if (fillStyle!=-999) ret->getAttFill()->SetFillStyle(fillStyle) ;
if (pc.getInt("moveToBack") && frame->numItems()>1) {
frame->drawBefore(frame->getObject(0)->GetName(), frame->getCurve()->GetName());
}
return ret ;
}
RooPlot* RooAbsReal::plotOn(RooPlot *frame, PlotOpt o) const
{
if (plotSanityChecks(frame)) return frame ;
RooArgSet projDataVars ;
if (o.projData) {
cxcoutD(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") have ProjData with observables = " << *o.projData->get() << endl ;
if (o.projDataSet) {
RooArgSet* tmp = (RooArgSet*) o.projData->get()->selectCommon(*o.projDataSet) ;
projDataVars.add(*tmp) ;
cxcoutD(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") have ProjDataSet = " << *o.projDataSet << " will only use this subset of projData" << endl ;
delete tmp ;
} else {
cxcoutD(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") using full ProjData" << endl ;
projDataVars.add(*o.projData->get()) ;
}
}
cxcoutD(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") ProjDataVars = " << projDataVars << endl ;
RooArgSet projectedVars ;
RooArgSet sliceSet ;
if (o.projSet) {
cxcoutD(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") have input projSet = " << *o.projSet << endl ;
makeProjectionSet(frame->getPlotVar(),o.projSet,projectedVars,kFALSE) ;
cxcoutD(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") calculated projectedVars = " << *o.projSet << endl ;
if (frame->getNormVars()) {
RooArgSet *sliceSetTmp = getObservables(*frame->getNormVars()) ;
cxcoutD(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") frame->getNormVars() that are also observables = " << *sliceSetTmp << endl ;
sliceSetTmp->remove(projectedVars,kTRUE,kTRUE) ;
sliceSetTmp->remove(*frame->getPlotVar(),kTRUE,kTRUE) ;
if (o.projData) {
RooArgSet* tmp = (RooArgSet*) projDataVars.selectCommon(*o.projSet) ;
sliceSetTmp->remove(*tmp,kTRUE,kTRUE) ;
delete tmp ;
}
if (sliceSetTmp->getSize()) {
coutI(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") plot on "
<< frame->getPlotVar()->GetName() << " represents a slice in " << *sliceSetTmp << endl ;
}
sliceSet.add(*sliceSetTmp) ;
delete sliceSetTmp ;
}
} else {
makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,kTRUE) ;
}
cxcoutD(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") projectedVars = " << projectedVars << " sliceSet = " << sliceSet << endl ;
RooArgSet* projDataNeededVars = 0 ;
if (o.projData) {
projDataNeededVars = (RooArgSet*) projectedVars.selectCommon(projDataVars) ;
projectedVars.remove(projDataVars,kTRUE,kTRUE) ;
}
RooAbsReal* realVar = (RooRealVar*) frame->getPlotVar() ;
RooArgSet* plotCloneSet = (RooArgSet*) RooArgSet(*realVar).snapshot(kTRUE) ;
if (!plotCloneSet) {
coutE(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") Couldn't deep-clone self, abort," << endl ;
return frame ;
}
RooRealVar* plotVar = (RooRealVar*) plotCloneSet->find(realVar->GetName());
if (projectedVars.getSize()) {
coutI(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") plot on " << plotVar->GetName()
<< " integrates over variables " << projectedVars
<< (o.projectionRangeName?Form(" in range %s",o.projectionRangeName):"") << endl;
}
if (projDataNeededVars && projDataNeededVars->getSize()>0) {
coutI(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") plot on " << plotVar->GetName()
<< " averages using data variables " << *projDataNeededVars << endl ;
}
RooArgSet* projectionCompList ;
RooArgSet* deps = getObservables(frame->getNormVars()) ;
deps->remove(projectedVars,kTRUE,kTRUE) ;
if (projDataNeededVars) {
deps->remove(*projDataNeededVars,kTRUE,kTRUE) ;
}
deps->remove(*plotVar,kTRUE,kTRUE) ;
deps->add(*plotVar) ;
if (checkObservables(deps)) {
coutE(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") error in checkObservables, abort" << endl ;
delete deps ;
delete plotCloneSet ;
if (projDataNeededVars) delete projDataNeededVars ;
return frame ;
}
RooArgSet normSet(*deps) ;
RooAbsReal *projection = (RooAbsReal*) createPlotProjection(normSet, &projectedVars, projectionCompList, o.projectionRangeName) ;
cxcoutD(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") plot projection object is " << projection->GetName() << endl ;
if (dologD(Plotting)) {
projection->printStream(ccoutD(Plotting),0,kVerbose) ;
}
RooArgSet fullNormSet(*deps) ;
fullNormSet.add(projectedVars) ;
if (projDataNeededVars && projDataNeededVars->getSize()>0) {
fullNormSet.add(*projDataNeededVars) ;
}
RooArgSet* compSet = projection->getComponents() ;
TIterator* iter = compSet->createIterator() ;
RooAbsArg* arg ;
while((arg=(RooAbsArg*)iter->Next())) {
RooAbsPdf* pdf = dynamic_cast<RooAbsPdf*>(arg) ;
if (pdf) {
pdf->selectNormalization(&fullNormSet) ;
}
}
delete iter ;
delete compSet ;
if (o.projData && projDataNeededVars && projDataNeededVars->getSize()>0) {
RooAbsData* projDataSel = (RooAbsData*)o.projData;
if (projDataNeededVars->getSize()<o.projData->get()->getSize()) {
RooArgSet* sliceDataSet = (RooArgSet*) sliceSet.selectCommon(*o.projData->get()) ;
TString cutString ;
if (sliceDataSet->getSize()>0) {
TIterator* iter2 = sliceDataSet->createIterator() ;
RooAbsArg* sliceVar ;
Bool_t first(kTRUE) ;
while((sliceVar=(RooAbsArg*)iter2->Next())) {
if (!first) {
cutString.Append("&&") ;
} else {
first=kFALSE ;
}
RooAbsRealLValue* real ;
RooAbsCategoryLValue* cat ;
if ((real = dynamic_cast<RooAbsRealLValue*>(sliceVar))) {
cutString.Append(Form("%s==%f",real->GetName(),real->getVal())) ;
} else if ((cat = dynamic_cast<RooAbsCategoryLValue*>(sliceVar))) {
cutString.Append(Form("%s==%d",cat->GetName(),cat->getIndex())) ;
}
}
delete iter2 ;
}
delete sliceDataSet ;
if (!cutString.IsNull()) {
projDataSel = ((RooAbsData*)o.projData)->reduce(*projDataNeededVars,cutString) ;
coutI(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") reducing given projection dataset to entries with " << cutString << endl ;
} else {
projDataSel = ((RooAbsData*)o.projData)->reduce(*projDataNeededVars) ;
}
coutI(Plotting) << "RooAbsReal::plotOn(" << GetName()
<< ") only the following components of the projection data will be used: " << *projDataNeededVars << endl ;
}
if (!o.binProjData && dynamic_cast<RooDataSet*>(projDataSel)!=0) {
TIterator* iter2 = projDataSel->get()->createIterator() ;
Bool_t allCat(kTRUE) ;
RooAbsArg* arg2 ;
while((arg2=(RooAbsArg*)iter2->Next())) {
if (!dynamic_cast<RooCategory*>(arg2)) allCat = kFALSE ;
}
delete iter2 ;
if (allCat) {
o.binProjData = kTRUE ;
coutI(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") unbinned projection dataset consist only of discrete variables,"
<< " performing projection with binned copy for optimization." << endl ;
}
}
if (o.binProjData) {
RooAbsData* tmp = new RooDataHist(Form("%s_binned",projDataSel->GetName()),"Binned projection data",*projDataSel->get(),*projDataSel) ;
if (projDataSel!=o.projData) delete projDataSel ;
projDataSel = tmp ;
}
projection->getVal(projDataSel->get()) ;
projection->attachDataSet(*projDataSel) ;
RooDataWeightedAverage dwa(Form("%sDataWgtAvg",GetName()),"Data Weighted average",*projection,*projDataSel,o.numCPU,o.interleave,kTRUE) ;
dwa.constOptimizeTestStatistic(Activate) ;
RooRealBinding projBind(dwa,*plotVar) ;
RooScaledFunc scaleBind(projBind,o.scaleFactor);
if (o.rangeLo==0 && o.rangeHi==0) {
o.rangeLo = frame->GetXaxis()->GetXmin() ;
o.rangeHi = frame->GetXaxis()->GetXmax() ;
}
TString curveName(projection->GetName()) ;
curveName.Append(Form("_DataAvg[%s]",projDataSel->get()->contentsString().c_str())) ;
if (sliceSet.getSize()>0) {
curveName.Append(Form("_Slice[%s]",sliceSet.contentsString().c_str())) ;
}
if (o.curveNameSuffix) {
curveName.Append(o.curveNameSuffix) ;
}
RooAbsReal::enableEvalErrorLogging(kTRUE) ;
RooCurve *curve = new RooCurve(projection->GetName(),projection->GetTitle(),scaleBind,
o.rangeLo,o.rangeHi,frame->GetNbinsX(),o.precision,o.precision,o.shiftToZero,o.wmode,o.numee,o.doeeval,o.eeval) ;
RooAbsReal::enableEvalErrorLogging(kFALSE) ;
curve->SetName(curveName.Data()) ;
if (o.addToCurveName) {
RooCurve* otherCurve = static_cast<RooCurve*>(frame->findObject(o.addToCurveName,RooCurve::Class())) ;
RooCurve* sumCurve = new RooCurve(projection->GetName(),projection->GetTitle(),*curve,*otherCurve,o.addToWgtSelf,o.addToWgtOther) ;
sumCurve->SetName(Form("%s_PLUS_%s",curve->GetName(),otherCurve->GetName())) ;
delete curve ;
curve = sumCurve ;
}
if (o.curveName) {
curve->SetName(o.curveName) ;
}
frame->addPlotable(curve, o.drawOptions);
if (projDataSel!=o.projData) delete projDataSel ;
} else {
if (o.rangeLo==0 && o.rangeHi==0) {
o.rangeLo = frame->GetXaxis()->GetXmin() ;
o.rangeHi = frame->GetXaxis()->GetXmax() ;
}
if (o.postRangeFracScale) {
if (!o.normRangeName) {
o.normRangeName = "plotRange" ;
plotVar->setRange("plotRange",o.rangeLo,o.rangeHi) ;
}
Bool_t tmp = _globalSelectComp ;
globalSelectComp(kTRUE) ;
RooAbsReal* intFrac = projection->createIntegral(*plotVar,*plotVar,o.normRangeName) ;
globalSelectComp(kTRUE) ;
o.scaleFactor /= intFrac->getVal() ;
globalSelectComp(tmp) ;
delete intFrac ;
}
RooAbsReal::enableEvalErrorLogging(kTRUE) ;
RooCurve *curve = new RooCurve(*projection,*plotVar,o.rangeLo,o.rangeHi,frame->GetNbinsX(),
o.scaleFactor,0,o.precision,o.precision,o.shiftToZero,o.wmode,o.numee,o.doeeval,o.eeval);
RooAbsReal::enableEvalErrorLogging(kFALSE) ;
TString curveName(projection->GetName()) ;
if (sliceSet.getSize()>0) {
curveName.Append(Form("_Slice[%s]",sliceSet.contentsString().c_str())) ;
}
if (o.curveNameSuffix) {
curveName.Append(o.curveNameSuffix) ;
}
curve->SetName(curveName.Data()) ;
if (o.addToCurveName) {
RooCurve* otherCurve = static_cast<RooCurve*>(frame->findObject(o.addToCurveName,RooCurve::Class())) ;
RooCurve* sumCurve = new RooCurve(projection->GetName(),projection->GetTitle(),*curve,*otherCurve,o.addToWgtSelf,o.addToWgtOther) ;
sumCurve->SetName(Form("%s_PLUS_%s",curve->GetName(),otherCurve->GetName())) ;
delete curve ;
curve = sumCurve ;
}
if (o.curveName) {
curve->SetName(o.curveName) ;
}
frame->addPlotable(curve, o.drawOptions, o.curveInvisible);
}
if (projDataNeededVars) delete projDataNeededVars ;
delete deps ;
delete projectionCompList ;
delete plotCloneSet ;
return frame;
}
RooPlot* RooAbsReal::plotSliceOn(RooPlot *frame, const RooArgSet& sliceSet, Option_t* drawOptions,
Double_t scaleFactor, ScaleType stype, const RooAbsData* projData) const
{
RooArgSet projectedVars ;
makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,kTRUE) ;
TIterator* iter = sliceSet.createIterator() ;
RooAbsArg* sliceArg ;
while((sliceArg=(RooAbsArg*)iter->Next())) {
RooAbsArg* arg = projectedVars.find(sliceArg->GetName()) ;
if (arg) {
projectedVars.remove(*arg) ;
} else {
coutI(Plotting) << "RooAbsReal::plotSliceOn(" << GetName() << ") slice variable "
<< sliceArg->GetName() << " was not projected anyway" << endl ;
}
}
delete iter ;
PlotOpt o ;
o.drawOptions = drawOptions ;
o.scaleFactor = scaleFactor ;
o.stype = stype ;
o.projData = projData ;
o.projSet = &projectedVars ;
return plotOn(frame,o) ;
}
RooPlot* RooAbsReal::plotAsymOn(RooPlot *frame, const RooAbsCategoryLValue& asymCat, PlotOpt o) const
{
if (plotSanityChecks(frame)) return frame ;
RooArgSet projDataVars ;
if (o.projData) {
if (o.projDataSet) {
RooArgSet* tmp = (RooArgSet*) o.projData->get()->selectCommon(*o.projDataSet) ;
projDataVars.add(*tmp) ;
delete tmp ;
} else {
projDataVars.add(*o.projData->get()) ;
}
}
if (!dependsOn(asymCat)) {
coutE(Plotting) << "RooAbsReal::plotAsymOn(" << GetName()
<< ") function doesn't depend on asymmetry category " << asymCat.GetName() << endl ;
return frame ;
}
if (!asymCat.isSignType()) {
coutE(Plotting) << "RooAbsReal::plotAsymOn(" << GetName()
<< ") asymmetry category must have 2 or 3 states with index values -1,0,1" << endl ;
return frame ;
}
RooArgSet projectedVars ;
RooArgSet sliceSet ;
if (o.projSet) {
makeProjectionSet(frame->getPlotVar(),o.projSet,projectedVars,kFALSE) ;
if (frame->getNormVars()) {
RooArgSet *sliceSetTmp = getObservables(*frame->getNormVars()) ;
sliceSetTmp->remove(projectedVars,kTRUE,kTRUE) ;
sliceSetTmp->remove(*frame->getPlotVar(),kTRUE,kTRUE) ;
if (o.projData) {
RooArgSet* tmp = (RooArgSet*) projDataVars.selectCommon(*o.projSet) ;
sliceSetTmp->remove(*tmp,kTRUE,kTRUE) ;
delete tmp ;
}
if (sliceSetTmp->getSize()) {
coutI(Plotting) << "RooAbsReal::plotAsymOn(" << GetName() << ") plot on "
<< frame->getPlotVar()->GetName() << " represents a slice in " << *sliceSetTmp << endl ;
}
sliceSet.add(*sliceSetTmp) ;
delete sliceSetTmp ;
}
} else {
makeProjectionSet(frame->getPlotVar(),frame->getNormVars(),projectedVars,kTRUE) ;
}
RooArgSet* projDataNeededVars = 0 ;
if (o.projData) {
projDataNeededVars = (RooArgSet*) projectedVars.selectCommon(projDataVars) ;
projectedVars.remove(projDataVars,kTRUE,kTRUE) ;
}
if (projectedVars.find(asymCat.GetName())) {
projectedVars.remove(*projectedVars.find(asymCat.GetName())) ;
}
RooAbsReal* realVar = (RooRealVar*) frame->getPlotVar() ;
RooRealVar* plotVar = (RooRealVar*) realVar->Clone() ;
if (projectedVars.getSize()) {
coutI(Plotting) << "RooAbsReal::plotAsymOn(" << GetName() << ") plot on " << plotVar->GetName()
<< " projects variables " << projectedVars << endl ;
}
if (projDataNeededVars && projDataNeededVars->getSize()>0) {
coutI(Plotting) << "RooAbsReal::plotOn(" << GetName() << ") plot on " << plotVar->GetName()
<< " averages using data variables "<< *projDataNeededVars << endl ;
}
RooAbsCategoryLValue* asymPos = (RooAbsCategoryLValue*) asymCat.Clone("asym_pos") ;
RooAbsCategoryLValue* asymNeg = (RooAbsCategoryLValue*) asymCat.Clone("asym_neg") ;
asymPos->setIndex(1) ;
asymNeg->setIndex(-1) ;
RooCustomizer* custPos = new RooCustomizer(*this,"pos") ;
RooCustomizer* custNeg = new RooCustomizer(*this,"neg") ;
custPos->setOwning(kTRUE) ;
custNeg->setOwning(kTRUE) ;
custPos->replaceArg(asymCat,*asymPos) ;
custNeg->replaceArg(asymCat,*asymNeg) ;
RooAbsReal* funcPos = (RooAbsReal*) custPos->build() ;
RooAbsReal* funcNeg = (RooAbsReal*) custNeg->build() ;
RooArgSet *posProjCompList, *negProjCompList ;
RooArgSet depPos(*plotVar,*asymPos) ;
RooArgSet depNeg(*plotVar,*asymNeg) ;
depPos.add(projDataVars) ;
depNeg.add(projDataVars) ;
const RooAbsReal *posProj = funcPos->createPlotProjection(depPos, &projectedVars, posProjCompList) ;
const RooAbsReal *negProj = funcNeg->createPlotProjection(depNeg, &projectedVars, negProjCompList) ;
if (!posProj || !negProj) {
coutE(Plotting) << "RooAbsReal::plotAsymOn(" << GetName() << ") Unable to create projections, abort" << endl ;
return frame ;
}
TString asymName(GetName()) ;
asymName.Append("_Asym[") ;
asymName.Append(asymCat.GetName()) ;
asymName.Append("]") ;
TString asymTitle(asymCat.GetName()) ;
asymTitle.Append(" Asymmetry of ") ;
asymTitle.Append(GetTitle()) ;
RooFormulaVar* funcAsym = new RooFormulaVar(asymName,asymTitle,"(@0-@1)/(@0+@1)",RooArgSet(*posProj,*negProj)) ;
if (o.projData) {
RooAbsData* projDataSel = (RooAbsData*)o.projData;
if (projDataNeededVars->getSize()<o.projData->get()->getSize()) {
RooArgSet* sliceDataSet = (RooArgSet*) sliceSet.selectCommon(*o.projData->get()) ;
TString cutString ;
if (sliceDataSet->getSize()>0) {
TIterator* iter = sliceDataSet->createIterator() ;
RooAbsArg* sliceVar ;
Bool_t first(kTRUE) ;
while((sliceVar=(RooAbsArg*)iter->Next())) {
if (!first) {
cutString.Append("&&") ;
} else {
first=kFALSE ;
}
RooAbsRealLValue* real ;
RooAbsCategoryLValue* cat ;
if ((real = dynamic_cast<RooAbsRealLValue*>(sliceVar))) {
cutString.Append(Form("%s==%f",real->GetName(),real->getVal())) ;
} else if ((cat = dynamic_cast<RooAbsCategoryLValue*>(sliceVar))) {
cutString.Append(Form("%s==%d",cat->GetName(),cat->getIndex())) ;
}
}
delete iter ;
}
delete sliceDataSet ;
if (!cutString.IsNull()) {
projDataSel = ((RooAbsData*)o.projData)->reduce(*projDataNeededVars,cutString) ;
coutI(Plotting) << "RooAbsReal::plotAsymOn(" << GetName()
<< ") reducing given projection dataset to entries with " << cutString << endl ;
} else {
projDataSel = ((RooAbsData*)o.projData)->reduce(*projDataNeededVars) ;
}
coutI(Plotting) << "RooAbsReal::plotAsymOn(" << GetName()
<< ") only the following components of the projection data will be used: " << *projDataNeededVars << endl ;
}
RooDataWeightedAverage dwa(Form("%sDataWgtAvg",GetName()),"Data Weighted average",*funcAsym,*projDataSel,o.numCPU,o.interleave,kTRUE) ;
dwa.constOptimizeTestStatistic(Activate) ;
RooRealBinding projBind(dwa,*plotVar) ;
((RooAbsReal*)posProj)->attachDataSet(*projDataSel) ;
((RooAbsReal*)negProj)->attachDataSet(*projDataSel) ;
RooScaledFunc scaleBind(projBind,o.scaleFactor);
if (o.rangeLo==0 && o.rangeHi==0) {
o.rangeLo = frame->GetXaxis()->GetXmin() ;
o.rangeHi = frame->GetXaxis()->GetXmax() ;
}
TString curveName(funcAsym->GetName()) ;
curveName.Append(Form("_DataAvg[%s]",projDataSel->get()->contentsString().c_str())) ;
if (sliceSet.getSize()>0) {
curveName.Append(Form("_Slice[%s]",sliceSet.contentsString().c_str())) ;
}
if (o.curveNameSuffix) {
curveName.Append(o.curveNameSuffix) ;
}
RooAbsReal::enableEvalErrorLogging(kTRUE) ;
RooCurve *curve = new RooCurve(funcAsym->GetName(),funcAsym->GetTitle(),scaleBind,
o.rangeLo,o.rangeHi,frame->GetNbinsX(),o.precision,o.precision,kFALSE,o.wmode,o.numee,o.doeeval,o.eeval) ;
RooAbsReal::enableEvalErrorLogging(kFALSE) ;
dynamic_cast<TAttLine*>(curve)->SetLineColor(2) ;
frame->addPlotable(curve, o.drawOptions);
ccoutW(Eval) << endl ;
if (projDataSel!=o.projData) delete projDataSel ;
} else {
if (o.rangeLo==0 && o.rangeHi==0) {
o.rangeLo = frame->GetXaxis()->GetXmin() ;
o.rangeHi = frame->GetXaxis()->GetXmax() ;
}
RooAbsReal::enableEvalErrorLogging(kTRUE) ;
RooCurve* curve= new RooCurve(*funcAsym,*plotVar,o.rangeLo,o.rangeHi,frame->GetNbinsX(),
o.scaleFactor,0,o.precision,o.precision,kFALSE,o.wmode,o.numee,o.doeeval,o.eeval);
RooAbsReal::enableEvalErrorLogging(kFALSE) ;
dynamic_cast<TAttLine*>(curve)->SetLineColor(2) ;
TString curveName(funcAsym->GetName()) ;
if (sliceSet.getSize()>0) {
curveName.Append(Form("_Slice[%s]",sliceSet.contentsString().c_str())) ;
}
if (o.curveNameSuffix) {
curveName.Append(o.curveNameSuffix) ;
}
curve->SetName(curveName.Data()) ;
frame->addPlotable(curve, o.drawOptions);
}
delete custPos ;
delete custNeg ;
delete posProjCompList ;
delete negProjCompList ;
delete asymPos ;
delete asymNeg ;
delete funcAsym ;
delete plotVar ;
return frame;
}
Bool_t RooAbsReal::plotSanityChecks(RooPlot* frame) const
{
if(0 == frame) {
coutE(Plotting) << ClassName() << "::" << GetName() << ":plotOn: frame is null" << endl;
return kTRUE;
}
RooAbsReal* var = frame->getPlotVar() ;
if(!var) {
coutE(Plotting) << ClassName() << "::" << GetName()
<< ":plotOn: frame does not specify a plot variable" << endl;
return kTRUE;
}
if(!dynamic_cast<RooAbsRealLValue*>(var)) {
coutE(Plotting) << ClassName() << "::" << GetName() << ":plotOn: cannot plot variable \""
<< var->GetName() << "\" of type " << var->ClassName() << endl;
return kTRUE;
}
if(!this->dependsOn(*var)) {
coutE(Plotting) << ClassName() << "::" << GetName() << ":plotOn: WARNING: variable is not an explicit dependent: "
<< var->GetName() << endl;
}
return kFALSE ;
}
void RooAbsReal::makeProjectionSet(const RooAbsArg* plotVar, const RooArgSet* allVars,
RooArgSet& projectedVars, Bool_t silent) const
{
cxcoutD(Plotting) << "RooAbsReal::makeProjectionSet(" << GetName() << ") plotVar = " << plotVar->GetName()
<< " allVars = " << (allVars?(*allVars):RooArgSet()) << endl ;
projectedVars.removeAll() ;
if (!allVars) return ;
projectedVars.add(*allVars) ;
RooAbsArg *found= projectedVars.find(plotVar->GetName());
if(found) {
projectedVars.remove(*found);
RooArgSet* plotServers = plotVar->getObservables(&projectedVars) ;
TIterator* psIter = plotServers->createIterator() ;
RooAbsArg* ps ;
while((ps=(RooAbsArg*)psIter->Next())) {
RooAbsArg* tmp = projectedVars.find(ps->GetName()) ;
if (tmp) {
cxcoutD(Plotting) << "RooAbsReal::makeProjectionSet(" << GetName() << ") removing " << tmp->GetName()
<< " from projection set because it a server of " << plotVar->GetName() << endl ;
projectedVars.remove(*tmp) ;
}
}
delete psIter ;
delete plotServers ;
if (!silent) {
coutW(Plotting) << "RooAbsReal::plotOn(" << GetName()
<< ") WARNING: cannot project out frame variable ("
<< found->GetName() << "), ignoring" << endl ;
}
}
TIterator* iter = allVars->createIterator() ;
RooAbsArg* arg ;
while((arg=(RooAbsArg*)iter->Next())) {
if (!dependsOnValue(*arg)) {
projectedVars.remove(*arg,kTRUE) ;
cxcoutD(Plotting) << "RooAbsReal::plotOn(" << GetName()
<< ") function doesn't depend on projection variable "
<< arg->GetName() << ", ignoring" << endl ;
}
}
delete iter ;
}
Bool_t RooAbsReal::isSelectedComp() const
{
return _selectComp || _globalSelectComp ;
}
void RooAbsReal::globalSelectComp(Bool_t flag)
{
_globalSelectComp = flag ;
}
RooAbsFunc *RooAbsReal::bindVars(const RooArgSet &vars, const RooArgSet* nset, Bool_t clipInvalid) const
{
RooAbsFunc *binding= new RooRealBinding(*this,vars,nset,clipInvalid);
if(binding && !binding->isValid()) {
coutE(InputArguments) << ClassName() << "::" << GetName() << ":bindVars: cannot bind to " << vars << endl ;
delete binding;
binding= 0;
}
return binding;
}
void RooAbsReal::copyCache(const RooAbsArg* source)
{
RooAbsReal* other = static_cast<RooAbsReal*>(const_cast<RooAbsArg*>(source)) ;
if (!other->_treeVar) {
_value = other->_value ;
} else {
if (source->getAttribute("FLOAT_TREE_BRANCH")) {
_value = other->_floatValue ;
} else if (source->getAttribute("INTEGER_TREE_BRANCH")) {
_value = other->_intValue ;
} else if (source->getAttribute("BYTE_TREE_BRANCH")) {
_value = other->_byteValue ;
} else if (source->getAttribute("UNSIGNED_INTEGER_TREE_BRANCH")) {
_value = other->_uintValue ;
}
}
setValueDirty() ;
}
void RooAbsReal::attachToTree(TTree& t, Int_t bufSize)
{
TString cleanName(cleanBranchName()) ;
TBranch* branch = t.GetBranch(cleanName) ;
if (branch) {
TLeaf* leaf = (TLeaf*)branch->GetListOfLeaves()->At(0) ;
Int_t dummy ;
TLeaf* counterLeaf = leaf->GetLeafCounter(dummy) ;
if (counterLeaf) {
coutE(Eval) << "RooAbsReal::attachToTree(" << GetName() << ") ERROR: TTree branch " << GetName()
<< " is an array and cannot be attached to a RooAbsReal" << endl ;
return ;
}
TString typeName(leaf->GetTypeName()) ;
if (!typeName.CompareTo("Float_t")) {
coutI(Eval) << "RooAbsReal::attachToTree(" << GetName() << ") TTree Float_t branch " << GetName()
<< " will be converted to double precision" << endl ;
setAttribute("FLOAT_TREE_BRANCH",kTRUE) ;
_treeVar = kTRUE ;
t.SetBranchAddress(cleanName,&_floatValue) ;
} else if (!typeName.CompareTo("Int_t")) {
coutI(Eval) << "RooAbsReal::attachToTree(" << GetName() << ") TTree Int_t branch " << GetName()
<< " will be converted to double precision" << endl ;
setAttribute("INTEGER_TREE_BRANCH",kTRUE) ;
_treeVar = kTRUE ;
t.SetBranchAddress(cleanName,&_intValue) ;
} else if (!typeName.CompareTo("UChar_t")) {
coutI(Eval) << "RooAbsReal::attachToTree(" << GetName() << ") TTree UChar_t branch " << GetName()
<< " will be converted to double precision" << endl ;
setAttribute("BYTE_TREE_BRANCH",kTRUE) ;
_treeVar = kTRUE ;
t.SetBranchAddress(cleanName,&_byteValue) ;
} else if (!typeName.CompareTo("UInt_t")) {
coutI(Eval) << "RooAbsReal::attachToTree(" << GetName() << ") TTree UInt_t branch " << GetName()
<< " will be converted to double precision" << endl ;
setAttribute("UNSIGNED_INTEGER_TREE_BRANCH",kTRUE) ;
_treeVar = kTRUE ;
t.SetBranchAddress(cleanName,&_uintValue) ;
} else {
t.SetBranchAddress(cleanName,&_value) ;
}
if (branch->GetCompressionLevel()<0) {
branch->SetCompressionLevel(1) ;
}
} else {
TString format(cleanName);
format.Append("/D");
branch = t.Branch(cleanName, &_value, (const Text_t*)format, bufSize);
branch->SetCompressionLevel(1) ;
}
}
void RooAbsReal::fillTreeBranch(TTree& t)
{
TBranch* branch = t.GetBranch(cleanBranchName()) ;
if (!branch) {
coutE(Eval) << "RooAbsReal::fillTreeBranch(" << GetName() << ") ERROR: not attached to tree: " << cleanBranchName() << endl ;
assert(0) ;
}
branch->Fill() ;
}
void RooAbsReal::setTreeBranchStatus(TTree& t, Bool_t active)
{
TBranch* branch = t.GetBranch(cleanBranchName()) ;
if (branch) {
t.SetBranchStatus(cleanBranchName(),active?1:0) ;
}
}
RooAbsArg *RooAbsReal::createFundamental(const char* newname) const
{
RooRealVar *fund= new RooRealVar(newname?newname:GetName(),GetTitle(),_value,getUnit());
fund->removeRange();
fund->setPlotLabel(getPlotLabel());
fund->setAttribute("fundamentalCopy");
return fund;
}
Bool_t RooAbsReal::matchArgs(const RooArgSet& allDeps, RooArgSet& analDeps,
const RooArgProxy& a) const
{
TList nameList ;
nameList.Add(new TObjString(a.absArg()->GetName())) ;
Bool_t result = matchArgsByName(allDeps,analDeps,nameList) ;
nameList.Delete() ;
return result ;
}
Bool_t RooAbsReal::matchArgs(const RooArgSet& allDeps, RooArgSet& analDeps,
const RooArgProxy& a, const RooArgProxy& b) const
{
TList nameList ;
nameList.Add(new TObjString(a.absArg()->GetName())) ;
nameList.Add(new TObjString(b.absArg()->GetName())) ;
Bool_t result = matchArgsByName(allDeps,analDeps,nameList) ;
nameList.Delete() ;
return result ;
}
Bool_t RooAbsReal::matchArgs(const RooArgSet& allDeps, RooArgSet& analDeps,
const RooArgProxy& a, const RooArgProxy& b,
const RooArgProxy& c) const
{
TList nameList ;
nameList.Add(new TObjString(a.absArg()->GetName())) ;
nameList.Add(new TObjString(b.absArg()->GetName())) ;
nameList.Add(new TObjString(c.absArg()->GetName())) ;
Bool_t result = matchArgsByName(allDeps,analDeps,nameList) ;
nameList.Delete() ;
return result ;
}
Bool_t RooAbsReal::matchArgs(const RooArgSet& allDeps, RooArgSet& analDeps,
const RooArgProxy& a, const RooArgProxy& b,
const RooArgProxy& c, const RooArgProxy& d) const
{
TList nameList ;
nameList.Add(new TObjString(a.absArg()->GetName())) ;
nameList.Add(new TObjString(b.absArg()->GetName())) ;
nameList.Add(new TObjString(c.absArg()->GetName())) ;
nameList.Add(new TObjString(d.absArg()->GetName())) ;
Bool_t result = matchArgsByName(allDeps,analDeps,nameList) ;
nameList.Delete() ;
return result ;
}
Bool_t RooAbsReal::matchArgs(const RooArgSet& allDeps, RooArgSet& analDeps,
const RooArgSet& refset) const
{
TList nameList ;
TIterator* iter = refset.createIterator() ;
RooAbsArg* arg ;
while ((arg=(RooAbsArg*)iter->Next())) {
nameList.Add(new TObjString(arg->GetName())) ;
}
delete iter ;
Bool_t result = matchArgsByName(allDeps,analDeps,nameList) ;
nameList.Delete() ;
return result ;
}
Bool_t RooAbsReal::matchArgsByName(const RooArgSet &allArgs, RooArgSet &matchedArgs,
const TList &nameList) const
{
RooArgSet matched("matched");
TIterator *iterator= nameList.MakeIterator();
TObjString *name = 0;
Bool_t isMatched(kTRUE);
while((isMatched && (name= (TObjString*)iterator->Next()))) {
RooAbsArg *found= allArgs.find(name->String().Data());
if(found) {
matched.add(*found);
}
else {
isMatched= kFALSE;
}
}
if (isMatched && (matched.getSize()!=nameList.GetSize())) {
isMatched = kFALSE ;
}
delete iterator;
if(isMatched) matchedArgs.add(matched);
return isMatched;
}
RooNumIntConfig* RooAbsReal::defaultIntegratorConfig()
{
return &RooNumIntConfig::defaultConfig() ;
}
RooNumIntConfig* RooAbsReal::specialIntegratorConfig() const
{
return _specIntegratorConfig ;
}
const RooNumIntConfig* RooAbsReal::getIntegratorConfig() const
{
const RooNumIntConfig* config = specialIntegratorConfig() ;
if (config) return config ;
return defaultIntegratorConfig() ;
}
void RooAbsReal::setIntegratorConfig(const RooNumIntConfig& config)
{
if (_specIntegratorConfig) {
delete _specIntegratorConfig ;
}
_specIntegratorConfig = new RooNumIntConfig(config) ;
}
void RooAbsReal::setIntegratorConfig()
{
if (_specIntegratorConfig) {
delete _specIntegratorConfig ;
}
_specIntegratorConfig = 0 ;
}
void RooAbsReal::selectNormalization(const RooArgSet*, Bool_t)
{
}
void RooAbsReal::selectNormalizationRange(const char*, Bool_t)
{
}
void RooAbsReal::setCacheCheck(Bool_t flag)
{
_cacheCheck = flag ;
}
Int_t RooAbsReal::getMaxVal(const RooArgSet& ) const
{
return 0 ;
}
Double_t RooAbsReal::maxVal(Int_t )
{
assert(1) ;
return 0 ;
}
void RooAbsReal::logEvalError(const char* message, const char* serverValueString) const
{
EvalError ee ;
ee.setMessage(message) ;
if (serverValueString) {
ee.setServerValues(serverValueString) ;
} else {
string srvval ;
ostringstream oss ;
Bool_t first(kTRUE) ;
for (Int_t i=0 ; i<numProxies() ; i++) {
RooAbsProxy* p = getProxy(i) ;
if (p->name()[0]=='!') continue ;
if (first) {
first=kFALSE ;
} else {
oss << ", " ;
}
getProxy(i)->print(oss,kTRUE) ;
}
ee.setServerValues(oss.str().c_str()) ;
}
ostringstream oss2 ;
printStream(oss2,kName|kClassName|kArgs,kInline) ;
if (!_doLogEvalError) {
coutE(Eval) << "RooAbsReal::logEvalError(" << GetName() << ") evaluation error, " << endl
<< " origin : " << oss2.str() << endl
<< " message : " << ee._msg << endl
<< " server values: " << ee._srvval << endl ;
} else {
_evalErrorList[this].first = oss2.str().c_str() ;
_evalErrorList[this].second.push_back(ee) ;
}
}
void RooAbsReal::clearEvalErrorLog()
{
if (!_doLogEvalError) return ;
_evalErrorList.clear() ;
}
void RooAbsReal::printEvalErrors(ostream& os, Int_t maxPerNode)
{
if (maxPerNode<0) return ;
map<const RooAbsArg*,pair<string,list<EvalError> > >::iterator iter = _evalErrorList.begin() ;
for(;iter!=_evalErrorList.end() ; ++iter) {
if (maxPerNode==0) {
os << iter->second.first ;
os << " has " << iter->second.second.size() << " errors" << endl ;
} else {
os << iter->second.first << endl ;
Int_t i(0) ;
std::list<EvalError>::iterator iter2 = iter->second.second.begin() ;
for(;iter2!=iter->second.second.end() ; ++iter2, i++) {
os << " " << iter2->_msg << " @ " << iter2->_srvval << endl ;
if (i>maxPerNode) {
os << " ... (remaining " << iter->second.second.size() - maxPerNode << " messages suppressed)" << endl ;
break ;
}
}
}
}
}
Int_t RooAbsReal::numEvalErrors()
{
Int_t ntot(0) ;
map<const RooAbsArg*,pair<string,list<EvalError> > >::iterator iter = _evalErrorList.begin() ;
for(;iter!=_evalErrorList.end() ; ++iter) {
ntot += iter->second.second.size() ;
}
return ntot ;
}
void RooAbsReal::fixAddCoefNormalization(const RooArgSet& addNormSet, Bool_t force)
{
RooArgSet* compSet = getComponents() ;
TIterator* iter = compSet->createIterator() ;
RooAbsArg* arg ;
while((arg=(RooAbsArg*)iter->Next())) {
RooAbsPdf* pdf = dynamic_cast<RooAbsPdf*>(arg) ;
if (pdf) {
if (addNormSet.getSize()>0) {
pdf->selectNormalization(&addNormSet,force) ;
} else {
pdf->selectNormalization(0,force) ;
}
}
}
delete iter ;
delete compSet ;
}
void RooAbsReal::fixAddCoefRange(const char* rangeName, Bool_t force)
{
RooArgSet* compSet = getComponents() ;
TIterator* iter = compSet->createIterator() ;
RooAbsArg* arg ;
while((arg=(RooAbsArg*)iter->Next())) {
RooAbsPdf* pdf = dynamic_cast<RooAbsPdf*>(arg) ;
if (pdf) {
pdf->selectNormalizationRange(rangeName,force) ;
}
}
delete iter ;
delete compSet ;
}
void RooAbsReal::preferredObservableScanOrder(const RooArgSet& obs, RooArgSet& orderedObs) const
{
orderedObs.removeAll() ;
orderedObs.add(obs) ;
}
RooAbsReal* RooAbsReal::createRunningIntegral(const RooArgSet& iset, const RooArgSet& nset)
{
return createRunningIntegral(iset,RooFit::SupNormSet(nset)) ;
}
RooAbsReal* RooAbsReal::createRunningIntegral(const RooArgSet& iset, const RooCmdArg arg1, const RooCmdArg arg2,
const RooCmdArg arg3, const RooCmdArg arg4, const RooCmdArg arg5,
const RooCmdArg arg6, const RooCmdArg arg7, const RooCmdArg arg8)
{
RooCmdConfig pc(Form("RooAbsReal::createRunningIntegral(%s)",GetName())) ;
pc.defineObject("supNormSet","SupNormSet",0,0) ;
pc.defineInt("numScanBins","ScanParameters",0,1000) ;
pc.defineInt("intOrder","ScanParameters",1,2) ;
pc.defineInt("doScanNum","ScanNum",0,1) ;
pc.defineInt("doScanAll","ScanAll",0,0) ;
pc.defineInt("doScanNon","ScanNone",0,0) ;
pc.defineMutex("ScanNum","ScanAll","ScanNone") ;
pc.process(arg1,arg2,arg3,arg4,arg5,arg6,arg7,arg8) ;
if (!pc.ok(kTRUE)) {
return 0 ;
}
const RooArgSet* snset = static_cast<const RooArgSet*>(pc.getObject("supNormSet",0)) ;
RooArgSet nset ;
if (snset) {
nset.add(*snset) ;
}
Int_t numScanBins = pc.getInt("numScanBins") ;
Int_t intOrder = pc.getInt("intOrder") ;
Int_t doScanNum = pc.getInt("doScanNum") ;
Int_t doScanAll = pc.getInt("doScanAll") ;
Int_t doScanNon = pc.getInt("doScanNon") ;
if (doScanNon) {
return createIntRI(iset,nset) ;
}
if (doScanAll) {
return createScanRI(iset,nset,numScanBins,intOrder) ;
}
if (doScanNum) {
RooRealIntegral* tmp = (RooRealIntegral*) createIntegral(iset) ;
Int_t isNum= (tmp->numIntRealVars().getSize()==1) ;
delete tmp ;
if (isNum) {
coutI(NumIntegration) << "RooAbsPdf::createRunningIntegral(" << GetName() << ") integration over observable(s) " << iset << " involves numeric integration," << endl
<< " constructing cdf though numeric integration of sampled pdf in " << numScanBins << " bins and applying order "
<< intOrder << " interpolation on integrated histogram." << endl
<< " To override this choice of technique use argument ScanNone(), to change scan parameters use ScanParameters(nbins,order) argument" << endl ;
}
return isNum ? createScanRI(iset,nset,numScanBins,intOrder) : createIntRI(iset,nset) ;
}
return 0 ;
}
RooAbsReal* RooAbsReal::createScanRI(const RooArgSet& iset, const RooArgSet& nset, Int_t numScanBins, Int_t intOrder)
{
string name = string(GetName()) + "_NUMRUNINT_" + integralNameSuffix(iset,&nset).Data() ;
RooRealVar* ivar = (RooRealVar*) iset.first() ;
ivar->setBins(numScanBins,"numcdf") ;
RooNumRunningInt* ret = new RooNumRunningInt(name.c_str(),name.c_str(),*this,*ivar,"numrunint") ;
ret->setInterpolationOrder(intOrder) ;
return ret ;
}
RooAbsReal* RooAbsReal::createIntRI(const RooArgSet& iset, const RooArgSet& nset)
{
RooArgList ilist ;
TIterator* iter2 = iset.createIterator() ;
RooAbsArg* arg ;
while((arg=(RooAbsArg*)iter2->Next())) {
if (dynamic_cast<RooRealVar*>(arg)) {
ilist.add(*arg) ;
} else {
coutW(InputArguments) << "RooAbsPdf::createRunningIntegral(" << GetName() << ") WARNING ignoring non-RooRealVar input argument " << arg->GetName() << endl ;
}
}
delete iter2 ;
RooArgList cloneList ;
RooArgList loList ;
RooArgSet clonedBranchNodes ;
RooCustomizer cust(*this,"cdf") ;
cust.setCloneBranchSet(clonedBranchNodes) ;
cust.setOwning(kFALSE) ;
TIterator* iter = ilist.createIterator() ;
RooRealVar* rrv ;
while((rrv=(RooRealVar*)iter->Next())) {
RooRealVar* cloneArg = (RooRealVar*) rrv->clone(Form("%s_prime",rrv->GetName())) ;
cloneList.add(*cloneArg) ;
cust.replaceArg(*rrv,*cloneArg) ;
RooRealVar* cloneLo = (RooRealVar*) rrv->clone(Form("%s_lowbound",rrv->GetName())) ;
cloneLo->setVal(rrv->getMin()) ;
loList.add(*cloneLo) ;
RooParamBinning pb(*cloneLo,*rrv,100) ;
cloneArg->setBinning(pb,"CDF") ;
}
delete iter ;
RooAbsReal* tmp = (RooAbsReal*) cust.build() ;
RooArgSet finalNset(nset) ;
finalNset.add(cloneList,kTRUE) ;
RooAbsReal* cdf = tmp->createIntegral(cloneList,finalNset,"CDF") ;
cdf->addOwnedComponents(*tmp) ;
cdf->addOwnedComponents(cloneList) ;
cdf->addOwnedComponents(loList) ;
return cdf ;
}
Last change: Fri Dec 12 09:05:56 2008
Last generated: 2008-12-12 09:05
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.