#include "RooStats/SamplingDistPlot.h"
#include "RooRealVar.h"
#include "RooPlot.h"
#include <algorithm>
#include <iostream>
ClassImp(RooStats::SamplingDistPlot);
using namespace RooStats;
SamplingDistPlot::SamplingDistPlot() :
fhist(0) ,fItems()
{
fIterator = fItems.MakeIterator();
fbins = 100;
fMarkerType = 20;
fColor = 1;
}
SamplingDistPlot::SamplingDistPlot(const Int_t nbins) :
fhist(0) ,fItems()
{
fIterator = fItems.MakeIterator();
fbins = nbins;
fMarkerType = 20;
fColor = 1;
}
SamplingDistPlot::SamplingDistPlot(const char* name, const char* title, Int_t nbins, Double_t xmin, Double_t xmax) :
fhist(0) ,fItems()
{
fhist = new TH1F(name, title, nbins, xmin, xmax);
fbins = nbins;
fMarkerType = 20;
fColor = 1;
}
SamplingDistPlot::~SamplingDistPlot()
{
fSamplingDistr.clear();
fSampleWeights.clear();
fItems.Clear();
}
void SamplingDistPlot::AddSamplingDistribution(const SamplingDistribution *samplingDist, Option_t *drawOptions)
{
fSamplingDistr = samplingDist->GetSamplingDistribution();
SetSampleWeights(samplingDist);
TString options(drawOptions);
options.ToUpper();
if(!options.Contains("SAME")) options.Append("SAME");
if(!options.Contains("E1")) options.Append("E1");
const Double_t xlow = *(std::min_element(fSamplingDistr.begin(),fSamplingDistr.end()));
const Double_t xup = *(std::max_element(fSamplingDistr.begin(),fSamplingDistr.end()));
fhist = new TH1F(samplingDist->GetName(),samplingDist->GetTitle(),fbins,xlow,xup);
fhist->GetXaxis()->SetTitle(samplingDist->GetVarName().Data());
fVarName = samplingDist->GetVarName().Data();
std::vector<Double_t>::iterator valuesIt = fSamplingDistr.begin();
for(int w_idx = 0; valuesIt != fSamplingDistr.end(); ++valuesIt, ++w_idx)
{
if(fIsWeighted) fhist->Fill(*valuesIt,fSampleWeights[w_idx]);
else fhist->Fill(*valuesIt);
}
fhist->Sumw2() ;
fhist->SetMarkerStyle(fMarkerType);
fhist->SetMarkerColor(fColor);
fhist->SetLineColor(fColor);
fMarkerType++;
fColor++;
fhist->SetStats(kFALSE);
addObject(fhist,options.Data());
return;
}
void SamplingDistPlot::SetSampleWeights(const SamplingDistribution* samplingDist)
{
fIsWeighted = kFALSE;
if(samplingDist->GetSampleWeights().size() != 0){
fIsWeighted = kTRUE;
fSampleWeights = samplingDist->GetSampleWeights();
}
return;
}
void SamplingDistPlot::addObject(TObject *obj, Option_t *drawOptions)
{
if(0 == obj) {
std::cerr << fName << "::addObject: called with a null pointer" << std::endl;
return;
}
fItems.Add(obj,drawOptions);
return;
}
void SamplingDistPlot::Draw(const Option_t * )
{
Float_t theMin(0.), theMax(0.), theYMax(0.);
GetAbsoluteInterval(theMin,theMax,theYMax);
RooRealVar xaxis("xaxis",fVarName.Data(),theMin,theMax);
RooPlot* frame = xaxis.frame();
frame->SetTitle("");
frame->SetMaximum(theYMax);
fIterator->Reset();
TH1F *obj = 0;
while((obj= (TH1F*)fIterator->Next()))
frame->addTH1(obj,fIterator->GetOption());
frame->Draw();
return;
}
void SamplingDistPlot::GetAbsoluteInterval(Float_t &theMin, Float_t &theMax, Float_t &theYMax) const
{
Float_t tmpmin = 999.;
Float_t tmpmax = -999.;
Float_t tmpYmax = -999.;
fIterator->Reset();
TH1F *obj = 0;
while((obj = (TH1F*)fIterator->Next())) {
if(obj->GetXaxis()->GetXmin() < tmpmin) tmpmin = obj->GetXaxis()->GetXmin();
if(obj->GetXaxis()->GetXmax() > tmpmax) tmpmax = obj->GetXaxis()->GetXmax();
if(obj->GetMaximum() > tmpYmax) tmpYmax = obj->GetMaximum() + 0.1*obj->GetMaximum();
}
theMin = tmpmin;
theMax = tmpmax;
theYMax = tmpYmax;
return;
}
void SamplingDistPlot::SetLineColor(const Color_t color, const SamplingDistribution *samplDist)
{
if(samplDist == 0){
fhist->SetLineColor(color);
}
else{
fIterator->Reset();
TH1F *obj = 0;
while((obj = (TH1F*)fIterator->Next())) {
if(!strcmp(obj->GetName(),samplDist->GetName())){
obj->SetLineColor(color);
break;
}
}
}
return;
}
void SamplingDistPlot::SetLineWidth(const Width_t lwidth, const SamplingDistribution *samplDist)
{
if(samplDist == 0){
fhist->SetLineWidth(lwidth);
}
else{
fIterator->Reset();
TH1F *obj = 0;
while((obj = (TH1F*)fIterator->Next())) {
if(!strcmp(obj->GetName(),samplDist->GetName())){
obj->SetLineWidth(lwidth);
break;
}
}
}
return;
}
void SamplingDistPlot::SetLineStyle(const Style_t style, const SamplingDistribution *samplDist)
{
if(samplDist == 0){
fhist->SetLineStyle(style);
}
else{
fIterator->Reset();
TH1F *obj = 0;
while((obj = (TH1F*)fIterator->Next())) {
if(!strcmp(obj->GetName(),samplDist->GetName())){
obj->SetLineStyle(style);
break;
}
}
}
return;
}
void SamplingDistPlot::SetMarkerStyle(const Style_t style, const SamplingDistribution *samplDist)
{
if(samplDist == 0){
fhist->SetMarkerStyle(style);
}
else{
fIterator->Reset();
TH1F *obj = 0;
while((obj = (TH1F*)fIterator->Next())) {
if(!strcmp(obj->GetName(),samplDist->GetName())){
obj->SetMarkerStyle(style);
break;
}
}
}
return;
}
void SamplingDistPlot::SetMarkerColor(const Color_t color, const SamplingDistribution *samplDist)
{
if(samplDist == 0){
fhist->SetMarkerColor(color);
}
else{
fIterator->Reset();
TH1F *obj = 0;
while((obj = (TH1F*)fIterator->Next())) {
if(!strcmp(obj->GetName(),samplDist->GetName())){
obj->SetMarkerColor(color);
break;
}
}
}
return;
}
void SamplingDistPlot::SetMarkerSize(const Size_t size, const SamplingDistribution *samplDist)
{
if(samplDist == 0){
fhist->SetMarkerSize(size);
}
else{
fIterator->Reset();
TH1F *obj = 0;
while((obj = (TH1F*)fIterator->Next())) {
if(!strcmp(obj->GetName(),samplDist->GetName())){
obj->SetMarkerSize(size);
break;
}
}
}
return;
}
void SamplingDistPlot::RebinDistribution(const Int_t rebinFactor, const SamplingDistribution *samplDist)
{
if(samplDist == 0){
fhist->Rebin(rebinFactor);
}
else{
fIterator->Reset();
TH1F *obj = 0;
while((obj = (TH1F*)fIterator->Next())) {
if(!strcmp(obj->GetName(),samplDist->GetName())){
obj->Rebin(rebinFactor);
break;
}
}
}
return;
}