82 #define alpha_Low "-5"
83 #define alpha_High "5"
84 #define NoHistConst_Low "0"
85 #define NoHistConst_High "2000"
88 using namespace RooFit ;
89 using namespace RooStats ;
92 ClassImp(RooStats::HistFactory::HistoToWorkspaceFactoryFast);
95 namespace HistFactory{
97 HistoToWorkspaceFactoryFast::HistoToWorkspaceFactoryFast() :
98 fNomLumi(1.0), fLumiError(0),
99 fLowBin(0), fHighBin(0)
102 HistoToWorkspaceFactoryFast::~HistoToWorkspaceFactoryFast(){
105 HistoToWorkspaceFactoryFast::HistoToWorkspaceFactoryFast(RooStats::HistFactory::Measurement& measurement ) :
106 fSystToFix( measurement.GetConstantParams() ),
107 fParamValues( measurement.GetParamValues() ),
108 fNomLumi( measurement.GetLumi() ),
109 fLumiError( measurement.GetLumi()*measurement.GetLumiRelErr() ),
110 fLowBin( measurement.GetBinLow() ),
111 fHighBin( measurement.GetBinHigh() ) {
114 SetFunctionsToPreprocess( measurement.GetPreprocessFunctions() );
118 void HistoToWorkspaceFactoryFast::ConfigureWorkspaceForMeasurement(
const std::string& ModelName, RooWorkspace* ws_single, Measurement& measurement ) {
125 ModelConfig * proto_config = (ModelConfig *) ws_single->obj(
"ModelConfig");
126 if( proto_config == NULL ) {
127 std::cout <<
"Error: Did not find 'ModelConfig' object in file: " << ws_single->GetName()
132 std::vector<std::string> poi_list = measurement.GetPOIList();
133 if( poi_list.size()==0 ) {
134 std::cout <<
"Warining: No Parametetrs of interest are set" << std::endl;
137 cout <<
"Setting Parameter(s) of Interest as: ";
138 for(
unsigned int i = 0; i < poi_list.size(); ++i) {
139 cout << poi_list.at(i) <<
" ";
144 for(
unsigned int i = 0; i < poi_list.size(); ++i ) {
145 std::string poi_name = poi_list.at(i);
146 RooRealVar* poi = (RooRealVar*) ws_single->var( poi_name.c_str() );
151 std::cout <<
"WARNING: Can't find parameter of interest: " << poi_name
152 <<
" in Workspace. Not setting in ModelConfig." << std::endl;
156 proto_config->SetParametersOfInterest(params);
159 std::string NewModelName =
"newSimPdf";
163 if( measurement.GetGammaSyst().size() > 0
164 || measurement.GetUniformSyst().size() > 0
165 || measurement.GetLogNormSyst().size() > 0
166 || measurement.GetNoSyst().size() > 0) {
167 HistoToWorkspaceFactoryFast::EditSyst( ws_single, (ModelName).c_str(),
168 measurement.GetGammaSyst(),
169 measurement.GetUniformSyst(),
170 measurement.GetLogNormSyst(),
171 measurement.GetNoSyst());
173 proto_config->SetPdf( *ws_single->pdf(
"newSimPdf" ) );
178 RooAbsData* expData = ws_single->data(
"asimovData");
180 std::cout <<
"Error: Failed to find dataset: " << expData
181 <<
" in workspace" << std::endl;
184 if(poi_list.size()!=0){
185 proto_config->GuessObsAndNuisance(*expData);
195 RooAbsPdf* pdf = ws_single->pdf(NewModelName.c_str());
196 if( !pdf ) pdf = ws_single->pdf( ModelName.c_str() );
197 const RooArgSet* observables = ws_single->set(
"observables");
200 std::string SnapShotName =
"NominalParamValues";
201 ws_single->saveSnapshot(SnapShotName.c_str(), ws_single->allVars());
203 for(
unsigned int i=0; i<measurement.GetAsimovDatasets().size(); ++i) {
206 RooStats::HistFactory::Asimov& asimov = measurement.GetAsimovDatasets().at(i);
207 std::string AsimovName = asimov.GetName();
209 std::cout <<
"Generating additional Asimov Dataset: " << AsimovName << std::endl;
210 asimov.ConfigureWorkspace(ws_single);
211 RooDataSet* asimov_dataset =
212 (RooDataSet*) AsymptoticCalculator::GenerateAsimovData(*pdf, *observables);
214 std::cout <<
"Importing Asimov dataset" << std::endl;
215 bool failure = ws_single->import(*asimov_dataset, Rename(AsimovName.c_str()));
217 std::cout <<
"Error: Failed to import Asimov dataset: " << AsimovName
219 delete asimov_dataset;
225 ws_single->loadSnapshot(SnapShotName.c_str());
228 delete asimov_dataset;
238 RooWorkspace* HistoToWorkspaceFactoryFast::MakeSingleChannelModel( Measurement& measurement, Channel& channel ) {
250 string ch_name = channel.GetName();
253 RooWorkspace* ws_single = this->MakeSingleChannelWorkspace(measurement, channel);
254 if( ws_single == NULL ) {
255 std::cout <<
"Error: Failed to make Single-Channel workspace for channel: " << ch_name
256 <<
" and measurement: " << measurement.GetName() << std::endl;
262 HistoToWorkspaceFactoryFast::ConfigureWorkspaceForMeasurement(
"model_"+ch_name, ws_single, measurement );
268 RooWorkspace* HistoToWorkspaceFactoryFast::MakeCombinedModel( Measurement& measurement ) {
281 HistoToWorkspaceFactoryFast factory( measurement );
284 vector<RooWorkspace*> channel_workspaces;
285 vector<string> channel_names;
287 for(
unsigned int chanItr = 0; chanItr < measurement.GetChannels().size(); ++chanItr ) {
289 HistFactory::Channel& channel = measurement.GetChannels().at( chanItr );
291 if( ! channel.CheckHistograms() ) {
292 std::cout <<
"MakeModelAndMeasurementsFast: Channel: " << channel.GetName()
293 <<
" has uninitialized histogram pointers" << std::endl;
297 string ch_name = channel.GetName();
298 channel_names.push_back(ch_name);
301 RooWorkspace* ws_single = factory.MakeSingleChannelModel( measurement, channel );
303 channel_workspaces.push_back(ws_single);
310 RooWorkspace* ws = factory.MakeCombinedModel( channel_names, channel_workspaces );
314 HistoToWorkspaceFactoryFast::ConfigureWorkspaceForMeasurement(
"simPdf", ws, measurement );
317 for (vector<RooWorkspace*>::iterator iter = channel_workspaces.begin() ; iter != channel_workspaces.end() ; ++iter) {
326 void HistoToWorkspaceFactoryFast::ProcessExpectedHisto(
const TH1* hist,RooWorkspace* proto,
327 string prefix,
string productPrefix,
330 cout <<
"processing hist " << hist->GetName() << endl;
332 cout <<
"hist is empty" << endl;
333 R__ASSERT(hist != 0);
338 if (fObsNameVec.empty() && !fObsName.empty()) { fObsNameVec.push_back(fObsName); }
339 R__ASSERT( fObsNameVec.size()>=1 && fObsNameVec.size()<=3 );
342 unsigned int histndim(1);
343 std::string classname = hist->ClassName();
344 if (classname.find(
"TH1")==0) { histndim=1; }
345 else if (classname.find(
"TH2")==0) { histndim=2; }
346 else if (classname.find(
"TH3")==0) { histndim=3; }
347 R__ASSERT( histndim==fObsNameVec.size() );
350 RooArgList observables;
351 std::vector<std::string>::iterator itr = fObsNameVec.begin();
352 for (
int idx=0; itr!=fObsNameVec.end(); ++itr, ++idx ) {
353 if ( !proto->var(itr->c_str()) ) {
354 const TAxis* axis(0);
355 if (idx==0) { axis = hist->GetXaxis(); }
356 if (idx==1) { axis = hist->GetYaxis(); }
357 if (idx==2) { axis = hist->GetZaxis(); }
358 Int_t nbins = axis->GetNbins();
359 Double_t xmin = axis->GetXmin();
360 Double_t xmax = axis->GetXmax();
362 proto->factory(Form(
"%s[%f,%f]",itr->c_str(),xmin,xmax));
363 proto->var(itr->c_str())->setBins(nbins);
365 observables.add( *proto->var(itr->c_str()) );
368 RooDataHist* histDHist =
new RooDataHist((prefix+
"nominalDHist").c_str(),
"",observables,hist);
369 RooHistFunc* histFunc =
new RooHistFunc((prefix+
"_nominal").c_str(),
"",observables,*histDHist,0) ;
371 proto->import(*histFunc);
374 proto->factory((
"prod:"+productPrefix+
"("+prefix+
"_nominal,"+systTerm+
")").c_str() );
381 void HistoToWorkspaceFactoryFast::AddMultiVarGaussConstraint(RooWorkspace* proto,
string prefix,
int lowBin,
int highBin, vector<string>& constraintTermNames){
385 TVectorD mean(highBin);
387 for(Int_t i=lowBin; i<highBin; ++i){
388 std::stringstream str;
390 RooRealVar* temp = proto->var((prefix+str.str()).c_str());
391 mean(i) = temp->getVal();
394 TMatrixDSym Cov(highBin-lowBin);
395 for(
int i=lowBin; i<highBin; ++i){
396 for(
int j=0; j<highBin-lowBin; ++j){
397 if(i==j) { Cov(i,j) = sqrt(mean(i)); }
398 else { Cov(i,j) = 0; }
403 RooArgList floating( *(proto->set(prefix.c_str() ) ) );
404 RooMultiVarGaussian constraint((prefix+
"Constraint").c_str(),
"",
405 floating, mean, Cov);
407 proto->import(constraint);
409 constraintTermNames.push_back(constraint.GetName());
412 void HistoToWorkspaceFactoryFast::LinInterpWithConstraint(RooWorkspace* proto,
const TH1* nominal,
413 std::vector<HistoSys> histoSysList,
414 string prefix,
string productPrefix,
416 vector<string>& constraintTermNames){
422 if (fObsNameVec.empty() && !fObsName.empty()) { fObsNameVec.push_back(fObsName); }
423 R__ASSERT( fObsNameVec.size()>=1 && fObsNameVec.size()<=3 );
426 unsigned int histndim(1);
427 std::string classname = nominal->ClassName();
428 if (classname.find(
"TH1")==0) { histndim=1; }
429 else if (classname.find(
"TH2")==0) { histndim=2; }
430 else if (classname.find(
"TH3")==0) { histndim=3; }
431 R__ASSERT( histndim==fObsNameVec.size() );
435 RooArgList observables;
436 std::vector<std::string>::iterator itr = fObsNameVec.begin();
437 for (
int idx=0; itr!=fObsNameVec.end(); ++itr, ++idx ) {
438 if ( !proto->var(itr->c_str()) ) {
439 const TAxis* axis(
nullptr);
440 if (idx==0) { axis = nominal->GetXaxis(); }
441 else if (idx==1) { axis = nominal->GetYaxis(); }
442 else if (idx==2) { axis = nominal->GetZaxis(); }
444 std::cout <<
"Error: Too many observables. "
445 <<
"HistFactory only accepts up to 3 observables (3d) "
449 Int_t nbins = axis->GetNbins();
450 Double_t xmin = axis->GetXmin();
451 Double_t xmax = axis->GetXmax();
453 proto->factory(Form(
"%s[%f,%f]",itr->c_str(),xmin,xmax));
454 proto->var(itr->c_str())->setBins(nbins);
456 observables.add( *proto->var(itr->c_str()) );
459 RooDataHist* nominalDHist =
new RooDataHist((prefix+
"nominalDHist").c_str(),
"",observables,nominal);
460 RooHistFunc* nominalFunc =
new RooHistFunc((prefix+
"nominal").c_str(),
"",observables,*nominalDHist,0) ;
463 RooArgList params( (
"alpha_Hist") );
465 string range=string(
"[")+alpha_Low+
","+alpha_High+
"]";
468 for(
unsigned int j=0; j<histoSysList.size(); ++j){
469 std::stringstream str;
472 HistoSys& histoSys = histoSysList.at(j);
473 string histoSysName = histoSys.GetName();
475 RooRealVar* temp = (RooRealVar*) proto->var((
"alpha_" + histoSysName).c_str());
478 temp = (RooRealVar*) proto->factory((
"alpha_" + histoSysName + range).c_str());
481 string command=(
"Gaussian::alpha_"+histoSysName+
"Constraint(alpha_"+histoSysName+
",nom_alpha_"+histoSysName+
"[0.,-10,10],1.)");
482 cout << command << endl;
483 constraintTermNames.push_back( proto->factory( command.c_str() )->GetName() );
484 proto->var((
"nom_alpha_"+histoSysName).c_str())->setConstant();
485 const_cast<RooArgSet*
>(proto->set(
"globalObservables"))->add(*proto->var((
"nom_alpha_"+histoSysName).c_str()));
492 vector<double> low, high;
493 RooArgSet lowSet, highSet;
495 for(
unsigned int j=0; j<histoSysList.size(); ++j){
496 std::stringstream str;
499 HistoSys& histoSys = histoSysList.at(j);
500 RooDataHist* lowDHist =
new RooDataHist((prefix+str.str()+
"lowDHist").c_str(),
"",observables, histoSys.GetHistoLow());
501 RooDataHist* highDHist =
new RooDataHist((prefix+str.str()+
"highDHist").c_str(),
"",observables, histoSys.GetHistoHigh());
502 RooHistFunc* lowFunc =
new RooHistFunc((prefix+str.str()+
"low").c_str(),
"",observables,*lowDHist,0) ;
503 RooHistFunc* highFunc =
new RooHistFunc((prefix+str.str()+
"high").c_str(),
"",observables,*highDHist,0) ;
504 lowSet.add(*lowFunc);
505 highSet.add(*highFunc);
509 PiecewiseInterpolation interp(prefix.c_str(),
"",*nominalFunc,lowSet,highSet,params);
510 interp.setPositiveDefinite();
511 interp.setAllInterpCodes(4);
513 RooArgSet observableSet(observables);
514 interp.setBinIntegrator(observableSet);
515 interp.forceNumInt();
517 proto->import(interp);
520 proto->factory((
"prod:"+productPrefix+
"("+prefix+
","+systTerm+
")").c_str() );
525 string HistoToWorkspaceFactoryFast::AddNormFactor(RooWorkspace* proto,
string& channel,
string& sigmaEpsilon, Sample& sample,
bool doRatio){
526 string overallNorm_times_sigmaEpsilon ;
529 vector<NormFactor> normList = sample.GetNormFactorList();
530 vector<string> normFactorNames, rangeNames;
532 if(normList.size() > 0){
534 for(vector<NormFactor>::iterator itr = normList.begin(); itr != normList.end(); ++itr){
536 NormFactor& norm = *itr;
539 if(!prodNames.empty()) prodNames +=
",";
541 varname = norm.GetName() +
"_" + channel;
544 varname=norm.GetName();
550 std::stringstream range;
551 range <<
"[" << norm.GetVal() <<
"," << norm.GetLow() <<
"," << norm.GetHigh() <<
"]";
553 if( proto->obj(varname.c_str()) == NULL) {
554 cout <<
"making normFactor: " << norm.GetName() << endl;
556 proto->factory((varname + range.str()).c_str());
559 if(norm.GetConst()) {
562 cout <<
"WARNING: Const attribute to <NormFactor> tag is deprecated, will ignore." <<
563 " Instead, add \n\t<ParamSetting Const=\"True\">" << varname <<
"</ParamSetting>\n" <<
564 " to your top-level XML's <Measurment> entry" << endl;
567 rangeNames.push_back(range.str());
568 normFactorNames.push_back(varname);
571 overallNorm_times_sigmaEpsilon = sample.GetName() +
"_" + channel +
"_overallNorm_x_sigma_epsilon";
572 proto->factory((
"prod::" + overallNorm_times_sigmaEpsilon +
"(" + prodNames +
"," + sigmaEpsilon +
")").c_str());
575 unsigned int rangeIndex=0;
576 for( vector<string>::iterator nit = normFactorNames.begin(); nit!=normFactorNames.end(); ++nit){
577 if( count (normFactorNames.begin(), normFactorNames.end(), *nit) > 1 ){
578 cout <<
"WARNING: <NormFactor Name =\""<<*nit<<
"\"> is duplicated for <Sample Name=\""
579 << sample.GetName() <<
"\">, but only one factor will be included. \n Instead, define something like"
580 <<
"\n\t<Function Name=\""<<*nit<<
"Squared\" Expresion=\""<<*nit<<
"*"<<*nit<<
"\" Var=\""<<*nit<<rangeNames.at(rangeIndex)
581 <<
"\"> \nin your top-level XML's <Measurment> entry and use <NormFactor Name=\""<<*nit<<
"Squared\" in your channel XML file."<< endl;
586 if(!overallNorm_times_sigmaEpsilon.empty())
587 return overallNorm_times_sigmaEpsilon;
592 void HistoToWorkspaceFactoryFast::AddConstraintTerms(RooWorkspace* proto, Measurement & meas,
string prefix,
594 std::vector<OverallSys>& systList,
595 vector<string>& constraintTermNames,
596 vector<string>& totSystTermNames) {
601 string range=string(
"[0,")+alpha_Low+
","+alpha_High+
"]";
602 totSystTermNames.push_back(prefix);
604 RooArgSet params(prefix.c_str());
605 vector<double> lowVec, highVec;
607 std::map<std::string, double>::iterator itconstr;
608 for(
unsigned int i = 0; i < systList.size(); ++i) {
610 OverallSys& sys = systList.at(i);
611 std::string strname = sys.GetName();
612 const char * name = strname.c_str();
615 if (meas.GetNoSyst().count(sys.GetName()) > 0 ) {
616 std::cout <<
"HistoToWorkspaceFast::AddConstraintTerm - skip systematic " << sys.GetName() << std::endl;
620 if (meas.GetGammaSyst().count(sys.GetName()) > 0 ) {
621 double relerr = meas.GetGammaSyst().find(sys.GetName() )->second;
623 std::cout <<
"HistoToWorkspaceFast::AddConstraintTerm - zero uncertainty assigned - skip systematic " << sys.GetName() << std::endl;
626 double tauVal = 1./(relerr*relerr);
627 double sqtau = 1./relerr;
628 RooAbsArg * beta = proto->factory(TString::Format(
"beta_%s[1,0,10]",name) );
630 RooAbsArg * yvar = proto->factory(TString::Format(
"nom_%s[%f,0,10]",beta->GetName(),tauVal)) ;
632 RooAbsArg * theta = proto->factory(TString::Format(
"theta_%s[%f]",name,1./tauVal));
634 RooAbsArg* alphaOfBeta = proto->factory(TString::Format(
"PolyVar::alphaOfBeta_%s(beta_%s,{%f,%f})",name,name,-sqtau,sqtau));
638 RooAbsArg * kappa = proto->factory(TString::Format(
"sum::k_%s(%s,1.)",name,yvar->GetName()) );
639 RooAbsArg * gamma = proto->factory(TString::Format(
"Gamma::%sConstraint(%s, %s, %s, 0.0)",beta->GetName(),beta->GetName(), kappa->GetName(), theta->GetName() ) );
640 alphaOfBeta->Print(
"t");
642 constraintTermNames.push_back(gamma->GetName());
644 RooRealVar * gobs =
dynamic_cast<RooRealVar*
>(yvar); assert(gobs);
645 gobs->setConstant(
true);
646 const_cast<RooArgSet*
>(proto->set(
"globalObservables"))->add(*yvar);
649 params.add(*alphaOfBeta);
650 std::cout <<
"Added a gamma constraint for " << name << std::endl;
659 double gaussSigma = 1;
660 if (meas.GetUniformSyst().count(sys.GetName()) > 0 ) {
662 std::cout <<
"Added a uniform constraint for " << name <<
" as a gaussian constraint with a very large sigma " << std::endl;
666 RooRealVar* alpha = (RooRealVar*) proto->var((prefix + sys.GetName()).c_str());
669 alpha = (RooRealVar*) proto->factory((prefix + sys.GetName() + range).c_str());
670 RooAbsArg * nomAlpha = proto->factory(TString::Format(
"nom_%s[0.,-10,10]",alpha->GetName() ) );
671 RooAbsArg * gausConstraint = proto->factory(TString::Format(
"Gaussian::%sConstraint(%s,%s,%f)",alpha->GetName(),alpha->GetName(), nomAlpha->GetName(), gaussSigma) );
673 constraintTermNames.push_back( gausConstraint->GetName() );
674 proto->var((
"nom_" + prefix + sys.GetName()).c_str())->setConstant();
675 const_cast<RooArgSet*
>(proto->set(
"globalObservables"))->add(*nomAlpha);
683 if (meas.GetLogNormSyst().count(sys.GetName()) == 0 && meas.GetGammaSyst().count(sys.GetName()) == 0 ) {
688 if (meas.GetLogNormSyst().count(sys.GetName()) > 0 ) {
690 double relerr = meas.GetLogNormSyst().find(sys.GetName() )->second;
691 double tauVal = 1./relerr;
692 std::string tauName =
"tau_" + sys.GetName();
693 proto->factory(TString::Format(
"%s[%f]",tauName.c_str(),tauVal ) );
694 double kappaVal = 1. + relerr;
695 std::string kappaName =
"kappa_" + sys.GetName();
696 proto->factory(TString::Format(
"%s[%f]",kappaName.c_str(),kappaVal ) );
697 const char * alphaName = alpha->GetName();
699 std::string alphaOfBetaName =
"alphaOfBeta_" + sys.GetName();
700 RooAbsArg * alphaOfBeta = proto->factory(TString::Format(
"expr::%s('%s*(pow(%s,%s)-1.)',%s,%s,%s)",alphaOfBetaName.c_str(),
701 tauName.c_str(),kappaName.c_str(),alphaName,
702 tauName.c_str(),kappaName.c_str(),alphaName ) );
704 std::cout <<
"Added a log-normal constraint for " << name << std::endl;
705 alphaOfBeta->Print(
"t");
706 params.add(*alphaOfBeta);
711 double low = sys.GetLow();
712 double high = sys.GetHigh();
713 lowVec.push_back(low);
714 highVec.push_back(high);
718 if(systList.size() > 0){
722 assert( params.getSize() > 0);
723 assert(
int(lowVec.size()) == params.getSize() );
725 FlexibleInterpVar interp( (interpName).c_str(),
"", params, 1., lowVec, highVec);
726 interp.setAllInterpCodes(4);
728 proto->import(interp);
732 RooConstVar interp( (interpName).c_str(),
"", 1.);
733 proto->import(interp);
742 void HistoToWorkspaceFactoryFast::MakeTotalExpected(RooWorkspace* proto,
string totName,
743 vector<string>& syst_x_expectedPrefixNames,
744 vector<string>& normByNames){
752 if (fObsNameVec.empty() && !fObsName.empty()) { fObsNameVec.push_back(fObsName); }
754 double binWidth(1.0);
755 std::string obsNameVecStr;
756 std::vector<std::string>::iterator itr = fObsNameVec.begin();
757 for (; itr!=fObsNameVec.end(); ++itr) {
758 std::string obsName = *itr;
759 binWidth *= proto->var(obsName.c_str())->numBins()/(proto->var(obsName.c_str())->getMax() - proto->var(obsName.c_str())->getMin()) ;
760 if (obsNameVecStr.size()>0) { obsNameVecStr +=
"_"; }
761 obsNameVecStr += obsName;
765 for(
unsigned int j=0; j<syst_x_expectedPrefixNames.size();++j){
766 std::stringstream str;
770 command=string(Form(
"binWidth_%s_%d[%e]",obsNameVecStr.c_str(),j,binWidth));
771 proto->factory(command.c_str());
772 proto->var(Form(
"binWidth_%s_%d",obsNameVecStr.c_str(),j))->setConstant();
773 coeffList+=prepend+
"binWidth_"+obsNameVecStr+str.str();
775 command=
"prod::L_x_"+syst_x_expectedPrefixNames.at(j)+
"("+normByNames.at(j)+
","+syst_x_expectedPrefixNames.at(j)+
")";
777 proto->factory(command.c_str());
778 shapeList+=prepend+
"L_x_"+syst_x_expectedPrefixNames.at(j);
787 proto->defineSet(
"coefList",coeffList.c_str());
788 proto->defineSet(
"shapeList",shapeList.c_str());
790 RooRealSumPdf tot(totName.c_str(),totName.c_str(),*proto->set(
"shapeList"),*proto->set(
"coefList"),kTRUE);
791 tot.specialIntegratorConfig(kTRUE)->method1D().setLabel(
"RooBinIntegrator") ;
792 tot.specialIntegratorConfig(kTRUE)->method2D().setLabel(
"RooBinIntegrator") ;
793 tot.specialIntegratorConfig(kTRUE)->methodND().setLabel(
"RooBinIntegrator") ;
797 tot.setAttribute(
"GenerateBinned");
821 void HistoToWorkspaceFactoryFast::AddPoissonTerms(RooWorkspace* proto,
string prefix,
string obsPrefix,
string expPrefix,
int lowBin,
int highBin,
822 vector<string>& likelihoodTermNames){
826 RooArgSet Pois(prefix.c_str());
827 for(Int_t i=lowBin; i<highBin; ++i){
828 std::stringstream str;
831 string command(
"Poisson::"+prefix+str.str()+
"("+obsPrefix+str.str()+
","+expPrefix+str.str()+
",1)");
832 RooAbsArg* temp = (proto->factory( command.c_str() ) );
835 cout <<
"Poisson Term " << command << endl;
836 ((RooAbsPdf*) temp)->setEvalErrorLoggingMode(RooAbsReal::PrintErrors);
839 likelihoodTermNames.push_back( temp->GetName() );
842 proto->defineSet(prefix.c_str(),Pois);
845 void HistoToWorkspaceFactoryFast::SetObsToExpected(RooWorkspace* proto,
string obsPrefix,
string expPrefix,
int lowBin,
int highBin){
848 TTree* tree =
new TTree();
849 Double_t* obsForTree =
new Double_t[highBin-lowBin];
850 RooArgList obsList(
"obsList");
852 for(Int_t i=lowBin; i<highBin; ++i){
853 std::stringstream str;
855 RooRealVar* obs = (RooRealVar*) proto->var((obsPrefix+str.str()).c_str());
856 cout <<
"expected number of events called: " << expPrefix << endl;
857 RooAbsReal* exp = proto->function((expPrefix+str.str()).c_str());
861 obs->setVal( exp->getVal() );
862 cout <<
"setting obs"+str.str()+
" to expected = " << exp->getVal() <<
" check: " << obs->getVal() << endl;
865 obsForTree[i] = exp->getVal();
866 tree->Branch((obsPrefix+str.str()).c_str(), obsForTree+i ,(obsPrefix+str.str()+
"/D").c_str());
869 cout <<
"problem retrieving obs or exp " << obsPrefix+str.str() << obs <<
" " << expPrefix+str.str() << exp << endl;
873 RooDataSet* data =
new RooDataSet(
"expData",
"", tree, obsList);
876 delete [] obsForTree;
878 proto->import(*data);
886 void HistoToWorkspaceFactoryFast::EditSyst(RooWorkspace* proto,
const char* pdfNameChar,
887 map<string,double> gammaSyst,
888 map<string,double> uniformSyst,
889 map<string,double> logNormSyst,
890 map<string,double> noSyst) {
891 string pdfName(pdfNameChar);
893 ModelConfig * combined_config = (ModelConfig *) proto->obj(
"ModelConfig");
894 if( combined_config==NULL ) {
895 std::cout <<
"Error: Failed to find object 'ModelConfig' in workspace: "
896 << proto->GetName() << std::endl;
901 string edit=
"EDIT::newSimPdf("+pdfName+
",";
903 string lastPdf=pdfName;
905 unsigned int numReplacements = 0;
906 unsigned int nskipped = 0;
907 map<string,double>::iterator it;
911 for(it=gammaSyst.begin(); it!=gammaSyst.end(); ++it) {
913 if(! proto->var((
"alpha_"+it->first).c_str())){
920 double relativeUncertainty = it->second;
921 double scale = 1/sqrt((1+1/pow(relativeUncertainty,2)));
924 proto->factory(Form(
"beta_%s[1,0,10]",it->first.c_str()));
925 proto->factory(Form(
"y_%s[%f]",it->first.c_str(),1./pow(relativeUncertainty,2))) ;
926 proto->factory(Form(
"theta_%s[%f]",it->first.c_str(),pow(relativeUncertainty,2))) ;
927 proto->factory(Form(
"Gamma::beta_%sConstraint(beta_%s,sum::k_%s(y_%s,one[1]),theta_%s,zero[0])",
932 it->first.c_str())) ;
948 proto->factory(Form(
"PolyVar::alphaOfBeta_%s(beta_%s,{%f,%f})",it->first.c_str(),it->first.c_str(),-1./scale,1./scale));
951 if(proto->var(Form(
"alpha_%s",it->first.c_str()))->isConstant()) {
952 proto->var(Form(
"beta_%s",it->first.c_str()))->setConstant(
true);
955 proto->var(Form(
"beta_%s",it->first.c_str()))->setConstant(
false);
961 editList+=precede +
"alpha_"+it->first+
"Constraint=beta_" + it->first+
"Constraint";
963 editList+=precede +
"alpha_"+it->first+
"=alphaOfBeta_"+ it->first;
973 if(numReplacements%10 == 0 && numReplacements+nskipped!=gammaSyst.size()){
974 edit=
"EDIT::"+lastPdf+
"_("+lastPdf+
","+editList+
")";
978 cout <<
"Going to issue this edit command\n" << edit<< endl;
979 proto->factory( edit.c_str() );
980 RooAbsPdf* newOne = proto->pdf(lastPdf.c_str());
982 cout <<
"\n\n ---------------------\n WARNING: failed to make EDIT\n\n" << endl;
988 for(it=uniformSyst.begin(); it!=uniformSyst.end(); ++it) {
989 cout <<
"edit for " << it->first <<
"with rel uncert = " << it->second << endl;
990 if(! proto->var((
"alpha_"+it->first).c_str())){
991 cout <<
"systematic not there" << endl;
998 proto->factory(Form(
"beta_%s[1,0,10]",it->first.c_str()));
999 proto->factory(Form(
"Uniform::beta_%sConstraint(beta_%s)",it->first.c_str(),it->first.c_str()));
1000 proto->factory(Form(
"PolyVar::alphaOfBeta_%s(beta_%s,{-1,1})",it->first.c_str(),it->first.c_str()));
1003 if(proto->var(Form(
"alpha_%s",it->first.c_str()))->isConstant())
1004 proto->var(Form(
"beta_%s",it->first.c_str()))->setConstant(
true);
1006 proto->var(Form(
"beta_%s",it->first.c_str()))->setConstant(
false);
1011 cout <<
"alpha_"+it->first+
"Constraint=beta_" + it->first+
"Constraint" << endl;
1012 editList+=precede +
"alpha_"+it->first+
"Constraint=beta_" + it->first+
"Constraint";
1014 cout <<
"alpha_"+it->first+
"=alphaOfBeta_"+ it->first << endl;
1015 editList+=precede +
"alpha_"+it->first+
"=alphaOfBeta_"+ it->first;
1017 if( proto->pdf((
"alpha_"+it->first+
"Constraint").c_str()) && proto->var((
"alpha_"+it->first).c_str()) )
1018 cout <<
" checked they are there" << proto->pdf((
"alpha_"+it->first+
"Constraint").c_str()) <<
" " << proto->var((
"alpha_"+it->first).c_str()) << endl;
1020 cout <<
"NOT THERE" << endl;
1023 if(numReplacements%10 == 0 && numReplacements+nskipped!=gammaSyst.size()){
1024 edit=
"EDIT::"+lastPdf+
"_("+lastPdf+
","+editList+
")";
1028 cout << edit<< endl;
1029 proto->factory( edit.c_str() );
1030 RooAbsPdf* newOne = proto->pdf(lastPdf.c_str());
1032 cout <<
"\n\n ---------------------\n WARNING: failed to make EDIT\n\n" << endl;
1042 for(it=logNormSyst.begin(); it!=logNormSyst.end(); ++it) {
1043 cout <<
"edit for " << it->first <<
"with rel uncert = " << it->second << endl;
1044 if(! proto->var((
"alpha_"+it->first).c_str())){
1045 cout <<
"systematic not there" << endl;
1051 double relativeUncertainty = it->second;
1052 double kappa = 1+relativeUncertainty;
1056 double scale = relativeUncertainty;
1059 const char * cname = it->first.c_str();
1062 proto->factory(TString::Format(
"beta_%s[1,0,10]",cname));
1063 proto->factory(TString::Format(
"nom_beta_%s[1]",cname));
1064 proto->factory(TString::Format(
"kappa_%s[%f]",cname,kappa));
1065 proto->factory(TString::Format(
"Lognormal::beta_%sConstraint(beta_%s,nom_beta_%s,kappa_%s)",
1066 cname, cname, cname, cname)) ;
1067 proto->factory(TString::Format(
"PolyVar::alphaOfBeta_%s(beta_%s,{%f,%f})",cname,cname,-1./scale,1./scale));
1071 if(proto->var(TString::Format(
"alpha_%s",cname))->isConstant())
1072 proto->var(TString::Format(
"beta_%s",cname))->setConstant(
true);
1074 proto->var(TString::Format(
"beta_%s",cname))->setConstant(
false);
1079 cout <<
"alpha_"+it->first+
"Constraint=beta_" + it->first+
"Constraint" << endl;
1080 editList+=precede +
"alpha_"+it->first+
"Constraint=beta_" + it->first+
"Constraint";
1082 cout <<
"alpha_"+it->first+
"=alphaOfBeta_"+ it->first << endl;
1083 editList+=precede +
"alpha_"+it->first+
"=alphaOfBeta_"+ it->first;
1085 if( proto->pdf((
"alpha_"+it->first+
"Constraint").c_str()) && proto->var((
"alpha_"+it->first).c_str()) )
1086 cout <<
" checked they are there" << proto->pdf((
"alpha_"+it->first+
"Constraint").c_str()) <<
" " << proto->var((
"alpha_"+it->first).c_str()) << endl;
1088 cout <<
"NOT THERE" << endl;
1091 if(numReplacements%10 == 0 && numReplacements+nskipped!=gammaSyst.size()){
1092 edit=
"EDIT::"+lastPdf+
"_("+lastPdf+
","+editList+
")";
1096 cout << edit<< endl;
1097 proto->factory( edit.c_str() );
1098 RooAbsPdf* newOne = proto->pdf(lastPdf.c_str());
1100 cout <<
"\n\n ---------------------\n WARNING: failed to make EDIT\n\n" << endl;
1104 const RooArgSet * gobs = proto->set(
"globalObservables");
1105 RooArgSet gobsNew(*gobs);
1106 gobsNew.add(*proto->var(TString::Format(
"nom_beta_%s",cname)) );
1107 proto->removeSet(
"globalObservables");
1108 proto->defineSet(
"globalObservables",gobsNew);
1116 for(it=noSyst.begin(); it!=noSyst.end(); ++it) {
1118 cout <<
"remove constraint for parameter" << it->first << endl;
1119 if(! proto->var((
"alpha_"+it->first).c_str()) || ! proto->pdf((
"alpha_"+it->first+
"Constraint").c_str()) ) {
1120 cout <<
"systematic not there" << endl;
1127 if ( !proto->var(
"one") ) { proto->factory(
"one[1.0]"); }
1128 proto->var(
"one")->setConstant();
1131 cout <<
"alpha_"+it->first+
"Constraint=one" << endl;
1132 editList+=precede +
"alpha_"+it->first+
"Constraint=one";
1136 if(numReplacements%10 == 0 && numReplacements+nskipped!=gammaSyst.size()){
1137 edit=
"EDIT::"+lastPdf+
"_("+lastPdf+
","+editList+
")";
1141 cout << edit << endl;
1142 proto->factory( edit.c_str() );
1143 RooAbsPdf* newOne = proto->pdf(lastPdf.c_str());
1144 if(!newOne) { cout <<
"\n\n ---------------------\n WARNING: failed to make EDIT\n\n" << endl; }
1151 edit=
"EDIT::newSimPdf("+lastPdf+
","+editList+
")";
1152 cout << edit<< endl;
1153 proto->factory( edit.c_str() );
1155 RooAbsPdf* newOne = proto->pdf(
"newSimPdf");
1158 combined_config->SetPdf(*newOne);
1161 cout <<
"\n\n ---------------------\n WARNING: failed to make EDIT\n\n" << endl;
1165 void HistoToWorkspaceFactoryFast::PrintCovarianceMatrix(RooFitResult* result, RooArgSet* params,
string filename){
1168 FILE* covFile = fopen ((filename).c_str(),
"w");
1170 TIter iti = params->createIterator();
1171 TIter itj = params->createIterator();
1172 RooRealVar *myargi, *myargj;
1173 fprintf(covFile,
" ") ;
1174 while ((myargi = (RooRealVar *)iti.Next())) {
1175 if(myargi->isConstant())
continue;
1176 fprintf(covFile,
" & %s", myargi->GetName());
1178 fprintf(covFile,
"\\\\ \\hline \n" );
1180 while ((myargi = (RooRealVar *)iti.Next())) {
1181 if(myargi->isConstant())
continue;
1182 fprintf(covFile,
"%s", myargi->GetName());
1184 while ((myargj = (RooRealVar *)itj.Next())) {
1185 if(myargj->isConstant())
continue;
1186 cout << myargi->GetName() <<
"," << myargj->GetName();
1187 fprintf(covFile,
" & %.2f", result->correlation(*myargi, *myargj));
1190 fprintf(covFile,
" \\\\\n");
1198 RooWorkspace* HistoToWorkspaceFactoryFast::MakeSingleChannelWorkspace(Measurement& measurement, Channel& channel) {
1202 if (channel.GetSamples().empty()) {
1203 Error(
"MakeSingleChannelWorkspace",
1204 "The input Channel does not contain any sample - return a nullptr");
1208 const TH1* channel_hist_template = channel.GetSamples().front().GetHisto();
1209 if (channel_hist_template ==
nullptr) {
1210 channel.CollectHistograms();
1211 channel_hist_template = channel.GetSamples().front().GetHisto();
1213 if (channel_hist_template ==
nullptr) {
1214 std::ostringstream stream;
1215 stream <<
"The sample " << channel.GetSamples().front().GetName()
1216 <<
" in channel " << channel.GetName() <<
" does not contain a histogram. This is the channel:\n";
1217 channel.Print(stream);
1218 Error(
"MakeSingleChannelWorkspace",
"%s", stream.str().c_str());
1222 if( ! channel.CheckHistograms() ) {
1223 std::cout <<
"MakeSingleChannelWorkspace: Channel: " << channel.GetName()
1224 <<
" has uninitialized histogram pointers" << std::endl;
1231 vector<string> systToFix = measurement.GetConstantParams();
1238 string channel_name = channel.GetName();
1241 fObsNameVec.clear();
1246 if (fObsNameVec.empty()) { GuessObsNameVec(channel_hist_template); }
1248 for (
unsigned int idx=0; idx<fObsNameVec.size(); ++idx ) {
1249 fObsNameVec[idx] =
"obs_" + fObsNameVec[idx] +
"_" + channel_name ;
1252 if (fObsNameVec.empty()) {
1253 fObsName=
"obs_" + channel_name;
1254 fObsNameVec.push_back( fObsName );
1257 R__ASSERT( fObsNameVec.size()>=1 && fObsNameVec.size()<=3 );
1259 cout <<
"\n\n-------------------\nStarting to process " << channel_name <<
" channel with " << fObsNameVec.size() <<
" observables" << endl;
1264 RooWorkspace* proto =
new RooWorkspace(channel_name.c_str(), (channel_name+
" workspace").c_str());
1265 auto proto_config = make_unique<ModelConfig>(
"ModelConfig", proto);
1266 proto_config->SetWorkspace(*proto);
1269 vector<string>::iterator funcIter = fPreprocessFunctions.begin();
1270 for(;funcIter!= fPreprocessFunctions.end(); ++funcIter){
1271 cout <<
"will preprocess this line: " << *funcIter <<endl;
1272 proto->factory(funcIter->c_str());
1276 RooArgSet likelihoodTerms(
"likelihoodTerms"), constraintTerms(
"constraintTerms");
1277 vector<string> likelihoodTermNames, constraintTermNames, totSystTermNames, syst_x_expectedPrefixNames, normalizationNames;
1279 vector< pair<string,string> > statNamePairs;
1280 vector< pair<const TH1*, const TH1*> > statHistPairs;
1281 std::string statFuncName;
1282 std::string statNodeName;
1286 string prefix, range;
1291 std::stringstream lumiStr;
1293 lumiStr<<
"["<<fNomLumi<<
",0,"<<10.*fNomLumi<<
"]";
1294 proto->factory((
"Lumi"+lumiStr.str()).c_str());
1295 cout <<
"lumi str = " << lumiStr.str() << endl;
1297 std::stringstream lumiErrorStr;
1298 lumiErrorStr <<
"nominalLumi["<<fNomLumi <<
",0,"<<fNomLumi+10*fLumiError<<
"]," << fLumiError ;
1299 proto->factory((
"Gaussian::lumiConstraint(Lumi,"+lumiErrorStr.str()+
")").c_str());
1300 proto->var(
"nominalLumi")->setConstant();
1301 proto->defineSet(
"globalObservables",
"nominalLumi");
1303 constraintTermNames.push_back(
"lumiConstraint");
1304 cout <<
"lumi Error str = " << lumiErrorStr.str() << endl;
1311 vector<Sample>::iterator it = channel.GetSamples().begin();
1312 for(; it!=channel.GetSamples().end(); ++it) {
1315 Sample& sample = (*it);
1316 string overallSystName = sample.GetName() +
"_" + channel_name +
"_epsilon";
1318 string systSourcePrefix =
"alpha_";
1322 AddConstraintTerms(proto,measurement, systSourcePrefix, overallSystName,
1323 sample.GetOverallSysList(), constraintTermNames , totSystTermNames);
1326 overallSystName = AddNormFactor(proto, channel_name, overallSystName, sample, doRatio);
1331 string syst_x_expectedPrefix =
"";
1335 const TH1* nominal = sample.GetHisto();
1348 if(sample.GetHistoSysList().size() == 0) {
1351 cout << sample.GetName() +
"_" + channel_name +
" has no variation histograms " << endl;
1352 string expPrefix = sample.GetName() +
"_" + channel_name;
1353 syst_x_expectedPrefix = sample.GetName() +
"_" + channel_name +
"_overallSyst_x_Exp";
1355 ProcessExpectedHisto(sample.GetHisto(), proto, expPrefix, syst_x_expectedPrefix,
1361 string constraintPrefix = sample.GetName() +
"_" + channel_name +
"_Hist_alpha";
1362 syst_x_expectedPrefix = sample.GetName() +
"_" + channel_name +
"_overallSyst_x_HistSyst";
1366 LinInterpWithConstraint(proto, nominal, sample.GetHistoSysList(),
1367 constraintPrefix, syst_x_expectedPrefix, overallSystName,
1368 constraintTermNames);
1375 if( sample.GetStatError().GetActivate() ) {
1377 if( fObsNameVec.size() > 3 ) {
1378 std::cout <<
"Cannot include Stat Error for histograms of more than 3 dimensions."
1387 std::cout <<
"Sample: " << sample.GetName() <<
" to be included in Stat Error "
1388 <<
"for channel " << channel_name
1418 string UncertName = syst_x_expectedPrefix +
"_StatAbsolUncert";
1419 TH1* statErrorHist = NULL;
1421 if( sample.GetStatError().GetErrorHist() == NULL ) {
1423 std::cout <<
"Making Statistical Uncertainty Hist for "
1424 <<
" Channel: " << channel_name
1425 <<
" Sample: " << sample.GetName()
1427 statErrorHist = MakeAbsolUncertaintyHist( UncertName, nominal );
1432 statErrorHist = (TH1*) sample.GetStatError().GetErrorHist()->Clone();
1436 std::cout <<
"Using external histogram for Stat Errors for "
1437 <<
" Channel: " << channel_name
1438 <<
" Sample: " << sample.GetName()
1440 std::cout <<
"Error Histogram: " << statErrorHist->GetName() << std::endl;
1443 statErrorHist->Multiply( nominal );
1444 statErrorHist->SetName( UncertName.c_str() );
1449 statHistPairs.push_back( std::make_pair(nominal, statErrorHist) );
1466 statFuncName =
"mc_stat_" + channel_name;
1467 ParamHistFunc* paramHist = (ParamHistFunc*) proto->function( statFuncName.c_str() );
1468 if( paramHist == NULL ) {
1472 RooArgList observables;
1473 std::vector<std::string>::iterator itr = fObsNameVec.begin();
1474 for (
int idx=0; itr!=fObsNameVec.end(); ++itr, ++idx ) {
1475 observables.add( *proto->var(itr->c_str()) );
1480 std::string ParamSetPrefix =
"gamma_stat_" + channel_name;
1481 Double_t gammaMin = 0.0;
1482 Double_t gammaMax = 10.0;
1483 RooArgList statFactorParams = ParamHistFunc::createParamSet(*proto,
1484 ParamSetPrefix.c_str(),
1486 gammaMin, gammaMax);
1488 ParamHistFunc statUncertFunc(statFuncName.c_str(), statFuncName.c_str(),
1489 observables, statFactorParams );
1491 proto->import( statUncertFunc, RecycleConflictNodes() );
1493 paramHist = (ParamHistFunc*) proto->function( statFuncName.c_str() );
1500 statNodeName = sample.GetName() +
"_" + channel_name +
"_overallSyst_x_StatUncert";
1502 RooAbsReal* expFunc = (RooAbsReal*) proto->function( syst_x_expectedPrefix.c_str() );
1503 RooProduct nodeWithMcStat(statNodeName.c_str(), statNodeName.c_str(),
1504 RooArgSet(*paramHist, *expFunc) );
1506 proto->import( nodeWithMcStat, RecycleConflictNodes() );
1511 syst_x_expectedPrefix = nodeWithMcStat.GetName();
1521 if( sample.GetShapeFactorList().size() > 0 ) {
1523 if( fObsNameVec.size() > 3 ) {
1524 std::cout <<
"Cannot include Stat Error for histograms of more than 3 dimensions."
1529 std::cout <<
"Sample: " << sample.GetName() <<
" in channel: " << channel_name
1530 <<
" to be include a ShapeFactor."
1533 std::vector<ParamHistFunc*> paramHistFuncList;
1534 std::vector<std::string> shapeFactorNameList;
1536 for(
unsigned int i=0; i < sample.GetShapeFactorList().size(); ++i) {
1538 ShapeFactor& shapeFactor = sample.GetShapeFactorList().at(i);
1540 std::string funcName = channel_name +
"_" + shapeFactor.GetName() +
"_shapeFactor";
1541 ParamHistFunc* paramHist = (ParamHistFunc*) proto->function( funcName.c_str() );
1542 if( paramHist == NULL ) {
1544 RooArgList observables;
1545 std::vector<std::string>::iterator itr = fObsNameVec.begin();
1546 for (
int idx=0; itr!=fObsNameVec.end(); ++itr, ++idx ) {
1547 observables.add( *proto->var(itr->c_str()) );
1551 std::string funcParams =
"gamma_" + shapeFactor.GetName();
1555 RooArgList shapeFactorParams = ParamHistFunc::createParamSet(*proto,
1557 observables, 0, 1000);
1560 ParamHistFunc shapeFactorFunc( funcName.c_str(), funcName.c_str(),
1561 observables, shapeFactorParams );
1564 if( shapeFactor.GetInitialShape() != NULL ) {
1565 TH1* initialShape =
static_cast<TH1*
>(shapeFactor.GetInitialShape()->Clone());
1566 std::cout <<
"Setting Shape Factor: " << shapeFactor.GetName()
1567 <<
" to have initial shape from hist: "
1568 << initialShape->GetName()
1570 shapeFactorFunc.setShape( initialShape );
1574 if( shapeFactor.IsConstant() ) {
1575 std::cout <<
"Setting Shape Factor: " << shapeFactor.GetName()
1576 <<
" to be constant" << std::endl;
1577 shapeFactorFunc.setConstant(
true);
1580 proto->import( shapeFactorFunc, RecycleConflictNodes() );
1581 paramHist = (ParamHistFunc*) proto->function( funcName.c_str() );
1585 paramHistFuncList.push_back(paramHist);
1586 shapeFactorNameList.push_back(funcName);
1595 std::string shapeFactorNodeName = syst_x_expectedPrefix;
1596 for(
unsigned int i=0; i < shapeFactorNameList.size(); ++i) {
1597 shapeFactorNodeName +=
"_x_" + shapeFactorNameList.at(i);
1600 RooAbsReal* expFunc = (RooAbsReal*) proto->function( syst_x_expectedPrefix.c_str() );
1601 RooArgSet nodesForProduct(*expFunc);
1602 for(
unsigned int i=0; i < paramHistFuncList.size(); ++i) {
1603 nodesForProduct.add( *paramHistFuncList.at(i) );
1608 RooProduct nodeWithShapeFactor(shapeFactorNodeName.c_str(),
1609 shapeFactorNodeName.c_str(),
1612 proto->import( nodeWithShapeFactor, RecycleConflictNodes() );
1617 syst_x_expectedPrefix = nodeWithShapeFactor.GetName();
1627 if( sample.GetShapeSysList().size() != 0 ) {
1629 if( fObsNameVec.size() > 3 ) {
1630 std::cout <<
"Cannot include Stat Error for histograms of more than 3 dimensions."
1636 std::vector<string> ShapeSysNames;
1638 for(
unsigned int i = 0; i < sample.GetShapeSysList().size(); ++i) {
1649 RooStats::HistFactory::ShapeSys& shapeSys = sample.GetShapeSysList().at(i);
1651 std::cout <<
"Sample: " << sample.GetName() <<
" in channel: " << channel_name
1652 <<
" to include a ShapeSys." << std::endl;
1654 std::string funcName = channel_name +
"_" + shapeSys.GetName() +
"_ShapeSys";
1655 ShapeSysNames.push_back( funcName );
1656 ParamHistFunc* paramHist = (ParamHistFunc*) proto->function( funcName.c_str() );
1657 if( paramHist == NULL ) {
1662 RooArgList observables;
1663 std::vector<std::string>::iterator itr = fObsNameVec.begin();
1664 for(; itr!=fObsNameVec.end(); ++itr ) {
1665 observables.add( *proto->var(itr->c_str()) );
1669 std::string funcParams =
"gamma_" + shapeSys.GetName();
1670 RooArgList shapeFactorParams = ParamHistFunc::createParamSet(*proto,
1672 observables, 0, 10);
1675 ParamHistFunc shapeFactorFunc( funcName.c_str(), funcName.c_str(),
1676 observables, shapeFactorParams );
1678 proto->import( shapeFactorFunc, RecycleConflictNodes() );
1679 paramHist = (ParamHistFunc*) proto->function( funcName.c_str() );
1688 const TH1* shapeErrorHist = shapeSys.GetErrorHist();
1691 Constraint::Type systype = shapeSys.GetConstraintType();
1692 if( systype == Constraint::Gaussian) {
1693 systype = Constraint::Gaussian;
1695 if( systype == Constraint::Poisson ) {
1696 systype = Constraint::Poisson;
1699 Double_t minShapeUncertainty = 0.0;
1700 RooArgList shapeConstraints = createStatConstraintTerms(proto, constraintTermNames,
1701 *paramHist, shapeErrorHist,
1703 minShapeUncertainty);
1711 std::string NodeName = syst_x_expectedPrefix;
1712 RooArgList ShapeSysForNode;
1713 RooAbsReal* expFunc = (RooAbsReal*) proto->function( syst_x_expectedPrefix.c_str() );
1714 ShapeSysForNode.add( *expFunc );
1715 for(
unsigned int i = 0; i < ShapeSysNames.size(); ++i ) {
1716 std::string ShapeSysName = ShapeSysNames.at(i);
1717 ShapeSysForNode.add( *proto->function(ShapeSysName.c_str()) );
1718 NodeName = NodeName +
"_x_" + ShapeSysName;
1722 RooProduct nodeWithShapeFactor(NodeName.c_str(), NodeName.c_str(), ShapeSysForNode );
1723 proto->import( nodeWithShapeFactor, RecycleConflictNodes() );
1728 syst_x_expectedPrefix = nodeWithShapeFactor.GetName();
1737 syst_x_expectedPrefixNames.push_back(syst_x_expectedPrefix);
1742 if( sample.GetNormalizeByTheory() ) {
1743 normalizationNames.push_back(
"Lumi" );
1746 TString lumiParamString;
1747 lumiParamString += measurement.GetLumi();
1748 lumiParamString.ReplaceAll(
' ', TString());
1749 normalizationNames.push_back(lumiParamString.Data());
1757 if( statHistPairs.size() > 0 ) {
1761 TH1* fracStatError = MakeScaledUncertaintyHist( statNodeName +
"_RelErr", statHistPairs );
1762 if( fracStatError == NULL ) {
1763 std::cout <<
"Error: Failed to make ScaledUncertaintyHist for: "
1764 << statNodeName << std::endl;
1770 ParamHistFunc* chanStatUncertFunc = (ParamHistFunc*) proto->function( statFuncName.c_str() );
1771 std::cout <<
"About to create Constraint Terms from: "
1772 << chanStatUncertFunc->GetName()
1773 <<
" params: " << chanStatUncertFunc->paramList()
1782 Constraint::Type statConstraintType = channel.GetStatErrorConfig().GetConstraintType();
1783 if( statConstraintType == Constraint::Gaussian) {
1784 std::cout <<
"Using Gaussian StatErrors in channel: " << channel.GetName() << std::endl;
1786 if( statConstraintType == Constraint::Poisson ) {
1787 std::cout <<
"Using Poisson StatErrors in channel: " << channel.GetName() << std::endl;
1790 double statRelErrorThreshold = channel.GetStatErrorConfig().GetRelErrorThreshold();
1791 RooArgList statConstraints = createStatConstraintTerms(proto, constraintTermNames,
1792 *chanStatUncertFunc, fracStatError,
1794 statRelErrorThreshold);
1798 for (
unsigned int i = 0; i < statHistPairs.size() ; ++i )
1799 delete statHistPairs[i].second;
1801 statHistPairs.clear();
1803 delete fracStatError;
1812 MakeTotalExpected(proto, channel_name+
"_model",
1813 syst_x_expectedPrefixNames, normalizationNames);
1814 likelihoodTermNames.push_back(channel_name+
"_model");
1818 for(
unsigned int i=0; i<systToFix.size(); ++i){
1819 RooRealVar* temp = proto->var((systToFix.at(i)).c_str());
1822 temp->setConstant();
1825 RooRealVar* auxMeas = NULL;
1826 if(systToFix.at(i)==
"Lumi"){
1827 auxMeas = proto->var(
"nominalLumi");
1829 auxMeas = proto->var(TString::Format(
"nom_%s",temp->GetName()));
1833 const_cast<RooArgSet*
>(proto->set(
"globalObservables"))->
remove(*auxMeas);
1835 cout <<
"could not corresponding auxiliary measurement "
1836 << TString::Format(
"nom_%s",temp->GetName()) << endl;
1839 cout <<
"could not find variable " << systToFix.at(i)
1840 <<
" could not set it to constant" << endl;
1846 for(
unsigned int i=0; i<constraintTermNames.size(); ++i){
1847 RooAbsArg* proto_arg = (proto->arg(constraintTermNames[i].c_str()));
1848 if( proto_arg==NULL ) {
1849 std::cout <<
"Error: Cannot find arg set: " << constraintTermNames.at(i)
1850 <<
" in workspace: " << proto->GetName() << std::endl;
1853 constraintTerms.add( *proto_arg );
1856 for(
unsigned int i=0; i<likelihoodTermNames.size(); ++i){
1857 RooAbsArg* proto_arg = (proto->arg(likelihoodTermNames[i].c_str()));
1858 if( proto_arg==NULL ) {
1859 std::cout <<
"Error: Cannot find arg set: " << likelihoodTermNames.at(i)
1860 <<
" in workspace: " << proto->GetName() << std::endl;
1863 likelihoodTerms.add( *proto_arg );
1865 proto->defineSet(
"constraintTerms",constraintTerms);
1866 proto->defineSet(
"likelihoodTerms",likelihoodTerms);
1870 RooArgList observables;
1871 std::string observablesStr;
1873 std::vector<std::string>::iterator itr = fObsNameVec.begin();
1874 for(; itr!=fObsNameVec.end(); ++itr ) {
1875 observables.add( *proto->var(itr->c_str()) );
1876 if (!observablesStr.empty()) { observablesStr +=
","; }
1877 observablesStr += *itr;
1883 proto->defineSet(
"observables", TString::Format(
"%s",observablesStr.c_str()));
1884 proto->defineSet(
"observablesSet", TString::Format(
"%s",observablesStr.c_str()));
1888 cout <<
"-----------------------------------------"<<endl;
1889 cout <<
"import model into workspace" << endl;
1891 auto model = make_unique<RooProdPdf>(
1892 (
"model_"+channel_name).c_str(),
1893 "product of Poissons accross bins for a single channel",
1894 constraintTerms, Conditional(likelihoodTerms,observables));
1895 proto->import(*model,RecycleConflictNodes());
1897 proto_config->SetPdf(*model);
1898 proto_config->SetObservables(observables);
1899 proto_config->SetGlobalObservables(*proto->set(
"globalObservables"));
1903 proto->import(*proto_config,proto_config->GetName());
1904 proto->importClassCode();
1909 const char* weightName=
"weightVar";
1910 proto->factory(TString::Format(
"%s[0,-1e10,1e10]",weightName));
1911 proto->defineSet(
"obsAndWeight",TString::Format(
"%s,%s",weightName,observablesStr.c_str()));
1919 RooDataSet* asimovDataUnbinned = new RooDataSet("asimovData","",*proto->set("obsAndWeight"),weightName);
1920 for(int i=0; i<asimov_data->numEntries(); ++i){
1921 asimov_data->get(i)->Print("v");
1922 //cout << "GREPME : " << i << " " << data->weight() <<endl;
1923 asimovDataUnbinned->add( *asimov_data->get(i), asimov_data->weight() );
1925 proto->import(*asimovDataUnbinned);
1930 unique_ptr<RooAbsData> asimov_dataset(AsymptoticCalculator::GenerateAsimovData(*model, observables));
1931 proto->import(dynamic_cast<RooDataSet&>(*asimov_dataset), Rename(
"asimovData"));
1934 if(channel.GetData().GetHisto() != NULL) {
1936 Data& data = channel.GetData();
1937 TH1* mnominal = data.GetHisto();
1939 std::cout <<
"Error: Data histogram for channel: " << channel.GetName()
1940 <<
" is NULL" << std::endl;
1945 auto obsDataUnbinned = make_unique<RooDataSet>(
"obsData",
"",*proto->set(
"obsAndWeight"),weightName);
1948 ConfigureHistFactoryDataset( obsDataUnbinned.get(), mnominal,
1949 proto, fObsNameVec );
1984 proto->import(*obsDataUnbinned);
1988 for(
unsigned int i=0; i < channel.GetAdditionalData().size(); ++i) {
1990 Data& data = channel.GetAdditionalData().at(i);
1991 std::string dataName = data.GetName();
1992 TH1* mnominal = data.GetHisto();
1994 std::cout <<
"Error: Additional Data histogram for channel: " << channel.GetName()
1995 <<
" with name: " << dataName <<
" is NULL" << std::endl;
2000 auto obsDataUnbinned = make_unique<RooDataSet>(dataName.c_str(), dataName.c_str(),
2001 *proto->set(
"obsAndWeight"), weightName);
2003 ConfigureHistFactoryDataset( obsDataUnbinned.get(), mnominal,
2004 proto, fObsNameVec );
2039 proto->import(*obsDataUnbinned);
2048 void HistoToWorkspaceFactoryFast::ConfigureHistFactoryDataset( RooDataSet* obsDataUnbinned,
2050 RooWorkspace* proto,
2051 std::vector<std::string> ObsNameVec) {
2057 if (ObsNameVec.empty() ) {
2058 Error(
"ConfigureHistFactoryDataset",
"Invalid input - return");
2064 TAxis* ax = mnominal->GetXaxis();
2065 TAxis* ay = mnominal->GetYaxis();
2066 TAxis* az = mnominal->GetZaxis();
2068 for (
int i=1; i<=ax->GetNbins(); ++i) {
2070 Double_t xval = ax->GetBinCenter(i);
2071 proto->var( ObsNameVec[0].c_str() )->setVal( xval );
2073 if(ObsNameVec.size()==1) {
2074 Double_t fval = mnominal->GetBinContent(i);
2075 obsDataUnbinned->add( *proto->set(
"obsAndWeight"), fval );
2078 for(
int j=1; j<=ay->GetNbins(); ++j) {
2079 Double_t yval = ay->GetBinCenter(j);
2080 proto->var( ObsNameVec[1].c_str() )->setVal( yval );
2082 if(ObsNameVec.size()==2) {
2083 Double_t fval = mnominal->GetBinContent(i,j);
2084 obsDataUnbinned->add( *proto->set(
"obsAndWeight"), fval );
2087 for(
int k=1; k<=az->GetNbins(); ++k) {
2088 Double_t zval = az->GetBinCenter(k);
2089 proto->var( ObsNameVec[2].c_str() )->setVal( zval );
2090 Double_t fval = mnominal->GetBinContent(i,j,k);
2091 obsDataUnbinned->add( *proto->set(
"obsAndWeight"), fval );
2099 void HistoToWorkspaceFactoryFast::GuessObsNameVec(
const TH1* hist)
2101 fObsNameVec.clear();
2104 unsigned int histndim(1);
2105 std::string classname = hist->ClassName();
2106 if (classname.find(
"TH1")==0) { histndim=1; }
2107 else if (classname.find(
"TH2")==0) { histndim=2; }
2108 else if (classname.find(
"TH3")==0) { histndim=3; }
2110 for (
unsigned int idx=0; idx<histndim; ++idx ) {
2111 if (idx==0) { fObsNameVec.push_back(
"x"); }
2112 if (idx==1) { fObsNameVec.push_back(
"y"); }
2113 if (idx==2) { fObsNameVec.push_back(
"z"); }
2118 RooWorkspace* HistoToWorkspaceFactoryFast::MakeCombinedModel(vector<string> ch_names, vector<RooWorkspace*> chs)
2123 if (ch_names.empty() || chs.empty() ) {
2124 Error(
"MakeCombinedModel",
"Input vectors are empty - return a nullptr");
2127 if (chs.size() < ch_names.size() ) {
2128 Error(
"MakeCombinedModel",
"Input vector of workspace has an invalid size - return a nullptr");
2136 map<string, RooAbsPdf*> pdfMap;
2137 vector<RooAbsPdf*> models;
2141 for(
unsigned int i = 0; i< ch_names.size(); ++i){
2142 ModelConfig * config = (ModelConfig *) chs[i]->obj(
"ModelConfig");
2143 obsList.add(*config->GetObservables());
2145 cout <<
"full list of observables:"<<endl;
2148 RooArgSet globalObs;
2149 for(
unsigned int i = 0; i< ch_names.size(); ++i){
2150 string channel_name=ch_names[i];
2152 if (ss.str().empty()) ss << channel_name ;
2153 else ss <<
',' << channel_name ;
2154 RooWorkspace * ch=chs[i];
2156 RooAbsPdf* model = ch->pdf((
"model_"+channel_name).c_str());
2157 if(!model) cout <<
"failed to find model for channel"<<endl;
2159 models.push_back(model);
2160 globalObs.add(*ch->set(
"globalObservables"));
2163 pdfMap[channel_name]=model;
2167 cout <<
"\n\n------------------\n Entering combination" << endl;
2168 RooWorkspace* combined =
new RooWorkspace(
"combined");
2172 RooCategory* channelCat = (RooCategory*) combined->factory((
"channelCat["+ss.str()+
"]").c_str());
2173 RooSimultaneous * simPdf=
new RooSimultaneous(
"simPdf",
"",pdfMap, *channelCat);
2174 ModelConfig * combined_config =
new ModelConfig(
"ModelConfig", combined);
2175 combined_config->SetWorkspace(*combined);
2178 combined->import(globalObs);
2179 combined->defineSet(
"globalObservables",globalObs);
2180 combined_config->SetGlobalObservables(*combined->set(
"globalObservables"));
2185 cout <<
"-----------------------------------------"<<endl;
2186 cout <<
"create toy data for " << ss.str() << endl;
2192 combined->factory(
"weightVar[0,-1e10,1e10]");
2193 obsList.add(*combined->var(
"weightVar"));
2212 RooDataSet* asimov_combined = (RooDataSet*) AsymptoticCalculator::GenerateAsimovData(*simPdf,
2214 if( asimov_combined ) {
2215 combined->import( *asimov_combined, Rename(
"asimovData"));
2218 std::cout <<
"Error: Failed to create combined asimov dataset" << std::endl;
2221 delete asimov_combined;
2224 if(chs[0]->data(
"obsData") != NULL) {
2225 MergeDataSets(combined, chs, ch_names,
"obsData", obsList, channelCat);
2273 obsList.add(*channelCat);
2274 combined->defineSet(
"observables",obsList);
2275 combined_config->SetObservables(*combined->set(
"observables"));
2279 cout <<
"\n\n----------------\n Importing combined model" << endl;
2280 combined->import(*simPdf,RecycleConflictNodes());
2285 std::map< std::string, double>::iterator param_itr = fParamValues.begin();
2286 for( ; param_itr != fParamValues.end(); ++param_itr ){
2288 std::string paramName = param_itr->first;
2289 double paramVal = param_itr->second;
2291 RooRealVar* temp = combined->var( paramName.c_str() );
2293 temp->setVal( paramVal );
2294 cout <<
"setting " << paramName <<
" to the value: " << paramVal << endl;
2296 cout <<
"could not find variable " << paramName <<
" could not set its value" << endl;
2300 for(
unsigned int i=0; i<fSystToFix.size(); ++i){
2302 RooRealVar* temp = combined->var((fSystToFix.at(i)).c_str());
2304 temp->setConstant();
2305 cout <<
"setting " << fSystToFix.at(i) <<
" constant" << endl;
2307 cout <<
"could not find variable " << fSystToFix.at(i) <<
" could not set it to constant" << endl;
2315 combined_config->SetPdf(*simPdf);
2318 combined->import(*combined_config,combined_config->GetName());
2319 combined->importClassCode();
2323 delete combined_config;
2330 RooDataSet* HistoToWorkspaceFactoryFast::MergeDataSets(RooWorkspace* combined,
2331 std::vector<RooWorkspace*> wspace_vec,
2332 std::vector<std::string> channel_names,
2333 std::string dataSetName,
2335 RooCategory* channelCat) {
2338 RooDataSet* simData=NULL;
2342 for(
unsigned int i = 0; i< channel_names.size(); ++i){
2345 std::cout <<
"Merging data for channel " << channel_names[i].c_str() << std::endl;
2346 RooDataSet* obsDataInChannel = (RooDataSet*) wspace_vec[i]->data(dataSetName.c_str());
2347 if( !obsDataInChannel ) {
2348 std::cout <<
"Error: Can't find DataSet: " << dataSetName
2349 <<
" in channel: " << channel_names.at(i)
2355 RooDataSet * tempData =
new RooDataSet(channel_names[i].c_str(),
"",
2356 obsList, Index(*channelCat),
2357 WeightVar(
"weightVar"),
2358 Import(channel_names[i].c_str(),*obsDataInChannel));
2360 simData->append(*tempData);
2371 combined->import(*simData, Rename(dataSetName.c_str()));
2373 simData = (RooDataSet*) combined->data(dataSetName.c_str() );
2376 std::cout <<
"Error: Unable to merge observable datasets" << std::endl;
2385 TH1* HistoToWorkspaceFactoryFast::MakeAbsolUncertaintyHist(
const std::string& Name,
const TH1* Nominal ) {
2391 TH1* ErrorHist = (TH1*) Nominal->Clone( Name.c_str() );
2394 Int_t numBins = Nominal->GetNbinsX()*Nominal->GetNbinsY()*Nominal->GetNbinsZ();
2395 Int_t binNumber = 0;
2398 for( Int_t i_bin = 0; i_bin < numBins; ++i_bin) {
2402 while( Nominal->IsBinUnderflow(binNumber) || Nominal->IsBinOverflow(binNumber) ){
2406 Double_t histError = Nominal->GetBinError( binNumber );
2409 if( histError != histError ) {
2410 std::cout <<
"Warning: In histogram " << Nominal->GetName()
2411 <<
" bin error for bin " << i_bin
2412 <<
" is NAN. Not using Error!!!"
2420 if( histError < 0 ) {
2421 std::cout <<
"Warning: In histogram " << Nominal->GetName()
2422 <<
" bin error for bin " << binNumber
2423 <<
" is < 0. Setting Error to 0"
2429 ErrorHist->SetBinContent( binNumber, histError );
2437 TH1* HistoToWorkspaceFactoryFast::MakeScaledUncertaintyHist(
const std::string& Name, std::vector< std::pair<const TH1*, const TH1*> > HistVec ) {
2449 unsigned int numHists = HistVec.size();
2451 if( numHists == 0 ) {
2452 std::cout <<
"Warning: Empty Hist Vector, cannot create total uncertainty" << std::endl;
2456 const TH1* HistTemplate = HistVec.at(0).first;
2457 Int_t numBins = HistTemplate->GetNbinsX()*HistTemplate->GetNbinsY()*HistTemplate->GetNbinsZ();
2461 for(
unsigned int i = 0; i < HistVec.size(); ++i ) {
2463 const TH1* nominal = HistVec.at(i).first;
2464 const TH1* error = HistVec.at(i).second;
2466 if( nominal->GetNbinsX()*nominal->GetNbinsY()*nominal->GetNbinsZ() != numBins ) {
2467 std::cout <<
"Error: Provided hists have unequal bins" << std::endl;
2470 if( error->GetNbinsX()*error->GetNbinsY()*error->GetNbinsZ() != numBins ) {
2471 std::cout <<
"Error: Provided hists have unequal bins" << std::endl;
2476 std::vector<double> TotalBinContent( numBins, 0.0);
2477 std::vector<double> HistErrorsSqr( numBins, 0.0);
2479 Int_t binNumber = 0;
2482 for( Int_t i_bins = 0; i_bins < numBins; ++i_bins) {
2485 while( HistTemplate->IsBinUnderflow(binNumber) || HistTemplate->IsBinOverflow(binNumber) ){
2489 for(
unsigned int i_hist = 0; i_hist < numHists; ++i_hist ) {
2491 const TH1* nominal = HistVec.at(i_hist).first;
2492 const TH1* error = HistVec.at(i_hist).second;
2496 Double_t histValue = nominal->GetBinContent( binNumber );
2497 Double_t histError = error->GetBinContent( binNumber );
2507 if( histError != histError ) {
2508 std::cout <<
"Warning: In histogram " << error->GetName()
2509 <<
" bin error for bin " << binNumber
2510 <<
" is NAN. Not using error!!"
2516 TotalBinContent.at(i_bins) += histValue;
2517 HistErrorsSqr.at(i_bins) += histError*histError;
2525 TH1* ErrorHist = (TH1*) HistTemplate->Clone( Name.c_str() );
2529 for( Int_t i = 0; i < numBins; ++i) {
2533 while( ErrorHist->IsBinUnderflow(binNumber) || ErrorHist->IsBinOverflow(binNumber) ){
2537 Double_t ErrorsSqr = HistErrorsSqr.at(i);
2538 Double_t TotalVal = TotalBinContent.at(i);
2540 if( TotalVal <= 0 ) {
2541 std::cout <<
"Warning: Sum of histograms for bin: " << binNumber
2542 <<
" is <= 0. Setting error to 0"
2545 ErrorHist->SetBinContent( binNumber, 0.0 );
2549 Double_t RelativeError = sqrt(ErrorsSqr) / TotalVal;
2553 if( RelativeError != RelativeError ) {
2554 std::cout <<
"Error: bin " << i <<
" error is NAN" << std::endl;
2555 std::cout <<
" HistErrorsSqr: " << ErrorsSqr
2556 <<
" TotalVal: " << TotalVal
2567 ErrorHist->SetBinError(binNumber, TotalVal);
2568 ErrorHist->SetBinContent(binNumber, RelativeError);
2570 std::cout <<
"Making Total Uncertainty for bin " << binNumber
2571 <<
" Error = " << sqrt(ErrorsSqr)
2572 <<
" Val = " << TotalVal
2573 <<
" RelativeError = " << RelativeError
2584 RooArgList HistoToWorkspaceFactoryFast::
2585 createStatConstraintTerms( RooWorkspace* proto, vector<string>& constraintTermNames,
2586 ParamHistFunc& paramHist,
const TH1* uncertHist,
2587 Constraint::Type type, Double_t minSigma ) {
2603 RooArgList ConstraintTerms;
2605 RooArgList paramSet = paramHist.paramList();
2609 Int_t numBins = uncertHist->GetNbinsX()*uncertHist->GetNbinsY()*uncertHist->GetNbinsZ();
2610 Int_t numParams = paramSet.getSize();
2615 if( numBins != numParams ) {
2616 std::cout <<
"Error: In createStatConstraintTerms, encountered bad number of bins" << std::endl;
2617 std::cout <<
"Given histogram with " << numBins <<
" bins,"
2618 <<
" but require exactly " << numParams << std::endl;
2622 Int_t TH1BinNumber = 0;
2623 for( Int_t i = 0; i < paramSet.getSize(); ++i) {
2627 while( uncertHist->IsBinUnderflow(TH1BinNumber) || uncertHist->IsBinOverflow(TH1BinNumber) ){
2631 RooRealVar& gamma = (RooRealVar&) (paramSet[i]);
2633 std::cout <<
"Creating constraint for: " << gamma.GetName()
2634 <<
". Type of constraint: " << type << std::endl;
2638 const double sigmaRel = uncertHist->GetBinContent(TH1BinNumber);
2642 if( sigmaRel <= 0 ){
2643 std::cout <<
"Not creating constraint term for "
2645 <<
" because sigma = " << sigmaRel
2647 <<
" (TH1 bin number = " << TH1BinNumber <<
")"
2649 gamma.setConstant(kTRUE);
2654 gamma.setMax( 1 + 5*sigmaRel );
2658 std::string constrName = string(gamma.GetName()) +
"_constraint";
2659 std::string nomName = string(
"nom_") + gamma.GetName();
2660 std::string sigmaName = string(gamma.GetName()) +
"_sigma";
2661 std::string poisMeanName = string(gamma.GetName()) +
"_poisMean";
2663 if( type == Constraint::Gaussian ) {
2669 RooConstVar constrSigma( sigmaName.c_str(), sigmaName.c_str(), sigmaRel );
2672 RooRealVar constrNom(nomName.c_str(), nomName.c_str(), 1.0,0,10);
2673 constrNom.setConstant(
true );
2676 RooGaussian gauss( constrName.c_str(), constrName.c_str(),
2677 constrNom, gamma, constrSigma );
2679 proto->import( gauss, RecycleConflictNodes() );
2683 gamma.setError(sigmaRel);
2684 }
else if( type == Constraint::Poisson ) {
2686 Double_t tau = 1/sigmaRel/sigmaRel;
2689 RooRealVar constrNom(nomName.c_str(), nomName.c_str(), tau);
2690 constrNom.setMin(0);
2691 constrNom.setConstant(
true );
2694 std::string scalingName = string(gamma.GetName()) +
"_tau";
2695 RooConstVar poissonScaling( scalingName.c_str(), scalingName.c_str(), tau);
2698 RooProduct constrMean( poisMeanName.c_str(), poisMeanName.c_str(), RooArgSet(gamma, poissonScaling) );
2703 RooPoisson pois(constrName.c_str(), constrName.c_str(), constrNom, constrMean);
2704 pois.setNoRounding(
true);
2705 proto->import( pois, RecycleConflictNodes() );
2709 gamma.setError(sigmaRel);
2713 std::cout <<
"Error: Did not recognize Stat Error constraint term type: "
2714 << type <<
" for : " << paramHist.GetName() << std::endl;
2721 if( sigmaRel < minSigma ) {
2722 std::cout <<
"Warning: Bin " << i <<
" = " << sigmaRel
2723 <<
" and is < " << minSigma
2724 <<
". Setting: " << gamma.GetName() <<
" to constant"
2726 gamma.setConstant(kTRUE);
2729 constraintTermNames.push_back( constrName );
2730 ConstraintTerms.add( *proto->pdf(constrName.c_str()) );
2734 RooArgSet* globalSet =
const_cast<RooArgSet*
>(proto->set(
"globalObservables"));
2736 RooRealVar* nomVarInWorkspace = proto->var(nomName.c_str());
2737 if( ! globalSet->contains(*nomVarInWorkspace) ) {
2738 globalSet->add( *nomVarInWorkspace );
2743 return ConstraintTerms;