Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
HistoToWorkspaceFactoryFast.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id: cranmer $
2// Author: Kyle Cranmer, Akira Shibata
3/*************************************************************************
4 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11////////////////////////////////////////////////////////////////////////////////
12
13/** \class RooStats::HistFactory::HistoToWorkspaceFactoryFast
14 * \ingroup HistFactory
15 * This class provides helper functions for creating likelihood models from histograms.
16 * It is used by RooStats::HistFactory::MakeModelAndMeasurementFast.
17 *
18 * A tutorial showing how to create a HistFactory model is hf001_example.C
19 */
20
21
22#ifndef __CINT__
23#include "RooGlobalFunc.h"
24#endif
25
26#include "RooDataSet.h"
27#include "RooRealVar.h"
28#include "RooConstVar.h"
29#include "RooAddition.h"
30#include "RooProduct.h"
31#include "RooProdPdf.h"
32#include "RooAddPdf.h"
33#include "RooGaussian.h"
34#include "RooPoisson.h"
35#include "RooExponential.h"
36#include "RooRandom.h"
37#include "RooCategory.h"
38#include "RooSimultaneous.h"
39#include "RooMultiVarGaussian.h"
40#include "RooNumIntConfig.h"
41#include "RooNLLVar.h"
42#include "RooProfileLL.h"
43#include "RooFitResult.h"
44#include "RooDataHist.h"
45#include "RooHistFunc.h"
46#include "RooHistPdf.h"
47#include "RooRealSumPdf.h"
48#include "RooWorkspace.h"
49#include "RooCustomizer.h"
50#include "RooPlot.h"
51#include "RooHelpers.h"
52#include "RooBinning.h"
53#include "RooBinWidthFunction.h"
59#include "HFMsgService.h"
60
61#include "TH1.h"
62#include "TTree.h"
63#include "TStopwatch.h"
64#include "TVectorD.h"
65#include "TMatrixDSym.h"
66
67// specific to this package
72#include "Helper.h"
73
74#include <algorithm>
75#include <memory>
76#include <utility>
77
78#define VERBOSE
79
80#define alpha_Low "-5"
81#define alpha_High "5"
82#define NoHistConst_Low "0"
83#define NoHistConst_High "2000"
84
85// use this order for safety on library loading
86using namespace RooFit ;
87using namespace RooStats ;
88using namespace std ;
89
91
92namespace RooStats{
93namespace HistFactory{
94
96 fNomLumi(1.0), fLumiError(0),
97 fLowBin(0), fHighBin(0)
98 {}
99
101 }
102
104 fSystToFix( measurement.GetConstantParams() ),
105 fParamValues( measurement.GetParamValues() ),
106 fNomLumi( measurement.GetLumi() ),
107 fLumiError( measurement.GetLumi()*measurement.GetLumiRelErr() ),
108 fLowBin( measurement.GetBinLow() ),
109 fHighBin( measurement.GetBinHigh() ) {
110
111 // Set Preprocess functions
113
114 }
115
116 void HistoToWorkspaceFactoryFast::ConfigureWorkspaceForMeasurement( const std::string& ModelName, RooWorkspace* ws_single, Measurement& measurement ) {
117
118 // Configure a workspace by doing any
119 // necessary post-processing and by
120 // creating a ModelConfig
121
122 // Make a ModelConfig and configure it
123 ModelConfig * proto_config = (ModelConfig *) ws_single->obj("ModelConfig");
124 if( proto_config == NULL ) {
125 std::cout << "Error: Did not find 'ModelConfig' object in file: " << ws_single->GetName()
126 << std::endl;
127 throw hf_exc();
128 }
129
130 std::vector<std::string> poi_list = measurement.GetPOIList();
131 if( poi_list.size()==0 ) {
132 cxcoutWHF << "No Parametetrs of interest are set" << std::endl;
133 }
134
135
136 std::stringstream sstream;
137 sstream << "Setting Parameter(s) of Interest as: ";
138 for(unsigned int i = 0; i < poi_list.size(); ++i) {
139 sstream << poi_list.at(i) << " ";
140 }
141 cxcoutIHF << sstream.str() << endl;
142
143 RooArgSet params;
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() );
147 if(poi){
148 params.add(*poi);
149 }
150 else {
151 std::cout << "WARNING: Can't find parameter of interest: " << poi_name
152 << " in Workspace. Not setting in ModelConfig." << std::endl;
153 //throw hf_exc();
154 }
155 }
156 proto_config->SetParametersOfInterest(params);
157
158 // Name of an 'edited' model, if necessary
159 std::string NewModelName = "newSimPdf"; // <- This name is hard-coded in HistoToWorkspaceFactoryFast::EditSyt. Probably should be changed to : std::string("new") + ModelName;
160
161 // Set the ModelConfig's Params of Interest
162 RooAbsData* expData = ws_single->data("asimovData");
163 if( !expData ) {
164 std::cout << "Error: Failed to find dataset: " << expData
165 << " in workspace" << std::endl;
166 throw hf_exc();
167 }
168 if(poi_list.size()!=0){
169 proto_config->GuessObsAndNuisance(*expData, RooMsgService::instance().isActive(static_cast<TObject*>(nullptr), RooFit::HistFactory, RooFit::INFO));
170 }
171
172 // Now, let's loop over any additional asimov datasets
173 // that we need to make
174
175 // Get the pdf
176 // Notice that we get the "new" pdf, this is the one that is
177 // used in the creation of these asimov datasets since they
178 // are fitted (or may be, at least).
179 RooAbsPdf* pdf = ws_single->pdf(NewModelName.c_str());
180 if( !pdf ) pdf = ws_single->pdf( ModelName.c_str() );
181 const RooArgSet* observables = ws_single->set("observables");
182
183 // Create a SnapShot of the nominal values
184 std::string SnapShotName = "NominalParamValues";
185 ws_single->saveSnapshot(SnapShotName.c_str(), ws_single->allVars());
186
187 for( unsigned int i=0; i<measurement.GetAsimovDatasets().size(); ++i) {
188
189 // Set the variable values and "const" ness with the workspace
190 RooStats::HistFactory::Asimov& asimov = measurement.GetAsimovDatasets().at(i);
191 std::string AsimovName = asimov.GetName();
192
193 cxcoutPHF << "Generating additional Asimov Dataset: " << AsimovName << std::endl;
194 asimov.ConfigureWorkspace(ws_single);
195 std::unique_ptr<RooDataSet> asimov_dataset{static_cast<RooDataSet*>(AsymptoticCalculator::GenerateAsimovData(*pdf, *observables))};
196
197 cxcoutPHF << "Importing Asimov dataset" << std::endl;
198 bool failure = ws_single->import(*asimov_dataset, Rename(AsimovName.c_str()));
199 if( failure ) {
200 std::cout << "Error: Failed to import Asimov dataset: " << AsimovName
201 << std::endl;
202 throw hf_exc();
203 }
204
205 // Load the snapshot at the end of every loop iteration
206 // so we start each loop with a "clean" snapshot
207 ws_single->loadSnapshot(SnapShotName.c_str());
208 }
209
210 // Cool, we're done
211 return; // ws_single;
212 }
213
214
215 // We want to eliminate this interface and use the measurment directly
217
218 // This is a pretty light-weight wrapper function
219 //
220 // Take a fully configured measurement as well as
221 // one of its channels
222 //
223 // Return a workspace representing that channel
224 // Do this by first creating a vector of EstimateSummary's
225 // and this by configuring the workspace with any post-processing
226
227 // Get the channel's name
228 string ch_name = channel.GetName();
229
230 // Create a workspace for a SingleChannel from the Measurement Object
231 RooWorkspace* ws_single = this->MakeSingleChannelWorkspace(measurement, channel);
232 if( ws_single == NULL ) {
233 cxcoutF(HistFactory) << "Error: Failed to make Single-Channel workspace for channel: " << ch_name
234 << " and measurement: " << measurement.GetName() << std::endl;
235 throw hf_exc();
236 }
237
238 // Finally, configure that workspace based on
239 // properties of the measurement
240 HistoToWorkspaceFactoryFast::ConfigureWorkspaceForMeasurement( "model_"+ch_name, ws_single, measurement );
241
242 return ws_single;
243
244 }
245
247
248 // This function takes a fully configured measurement
249 // which may contain several channels and returns
250 // a workspace holding the combined model
251 //
252 // This can be used, for example, within a script to produce
253 // a combined workspace on-the-fly
254 //
255 // This is a static function (for now) to make
256 // it a one-liner
257
259
260 // First, we create an instance of a HistFactory
261 HistoToWorkspaceFactoryFast factory( measurement );
262
263 // Loop over the channels and create the individual workspaces
264 vector<std::unique_ptr<RooWorkspace>> channel_workspaces;
265 vector<string> channel_names;
266
267 for( unsigned int chanItr = 0; chanItr < measurement.GetChannels().size(); ++chanItr ) {
268
269 HistFactory::Channel& channel = measurement.GetChannels().at( chanItr );
270
271 if( ! channel.CheckHistograms() ) {
272 cxcoutFHF << "MakeModelAndMeasurementsFast: Channel: " << channel.GetName()
273 << " has uninitialized histogram pointers" << std::endl;
274 throw hf_exc();
275 }
276
277 string ch_name = channel.GetName();
278 channel_names.push_back(ch_name);
279
280 // GHL: Renaming to 'MakeSingleChannelWorkspace'
281 RooWorkspace* ws_single = factory.MakeSingleChannelModel( measurement, channel );
282
283 channel_workspaces.emplace_back(ws_single);
284
285 }
286
287
288 // Now, combine the individual channel workspaces to
289 // form the combined workspace
290 RooWorkspace* ws = factory.MakeCombinedModel( channel_names, channel_workspaces );
291
292
293 // Configure the workspace
295
296 // Done. Return the pointer
297 return ws;
298
299 }
300
301/// Create observables of type RooRealVar. Creates 1 to 3 observables, depending on the type of the histogram.
303 RooArgList observables;
304
305 for (unsigned int idx=0; idx < fObsNameVec.size(); ++idx) {
306 if (!proto->var(fObsNameVec[idx].c_str())) {
307 const TAxis *axis = (idx == 0) ? hist->GetXaxis() : (idx == 1 ? hist->GetYaxis() : hist->GetZaxis());
308 Int_t nbins = axis->GetNbins();
309 Double_t xmin = axis->GetXmin();
310 Double_t xmax = axis->GetXmax();
311 // create observable
312 auto obs = static_cast<RooRealVar*>(proto->factory(
313 Form("%s[%f,%f]", fObsNameVec[idx].c_str(), xmin, xmax)));
314 obs->setBins(nbins);
315 if (axis->IsVariableBinSize()) {
316 RooBinning binning(nbins, axis->GetXbins()->GetArray());
317 obs->setBinning(binning);
318 }
319 }
320
321 observables.add(*proto->var(fObsNameVec[idx].c_str()));
322 }
323
324 return observables;
325}
326
327 /// Create the nominal hist function from `hist`, and register it in the workspace.
329 const RooArgList& observables) const {
330 if(hist) {
331 cxcoutI(HistFactory) << "processing hist " << hist->GetName() << endl;
332 } else {
333 cxcoutF(HistFactory) << "hist is empty" << endl;
334 R__ASSERT(hist != 0);
335 return nullptr;
336 }
337
338 // determine histogram dimensionality
339 unsigned int histndim(1);
340 std::string classname = hist->ClassName();
341 if (classname.find("TH1")==0) { histndim=1; }
342 else if (classname.find("TH2")==0) { histndim=2; }
343 else if (classname.find("TH3")==0) { histndim=3; }
344 R__ASSERT( histndim==fObsNameVec.size() );
345
346 prefix += "_Hist_alphanominal";
347
348 RooDataHist histDHist((prefix + "DHist").c_str(),"",observables,hist);
349 RooHistFunc histFunc(prefix.c_str(),"",observables,histDHist,0);
350
351 proto->import(histFunc, RecycleConflictNodes());
352 auto histFuncInWS = static_cast<RooHistFunc*>(proto->arg(prefix.c_str()));
353
354 return histFuncInWS;
355 }
356
357 namespace {
358
359 void makeGaussianConstraint(RooAbsArg& param, RooWorkspace& proto, bool isUniform,
360 std::vector<std::string> & constraintTermNames) {
361 std::string paramName = param.GetName();
362 std::string constraintName = paramName + "Constraint";
363
364 // do nothing if the constraint term already exists
365 if(proto.pdf(constraintName.c_str())) return;
366
367 // case systematic is uniform (asssume they are like a Gaussian but with
368 // a large width (100 instead of 1)
369 const double gaussSigma = isUniform ? 100. : 1.0;
370 if (isUniform) {
371 cxcoutIHF << "Added a uniform constraint for " << paramName << " as a Gaussian constraint with a very large sigma " << std::endl;
372 }
373
374 std::stringstream command;
375 command << "Gaussian::" << constraintName << "(" << paramName << ",nom_" << paramName << "[0.,-10,10],"
376 << gaussSigma << ")";
377 constraintTermNames.emplace_back(proto.factory( command.str().c_str() )->GetName());
378 auto * normParam = proto.var((std::string("nom_") + paramName).c_str());
379 normParam->setConstant();
380 const_cast<RooArgSet*>(proto.set("globalObservables"))->add(*normParam);
381 }
382
383 /// Make list of abstract parameters that interpolate in space of variations.
384 RooArgList makeInterpolationParameters(std::vector<HistoSys> const& histoSysList, RooWorkspace& proto) {
385 RooArgList params( ("alpha_Hist") );
386
387 // range is set using defined macro (see top of the page)
388 string range=string("[")+alpha_Low+","+alpha_High+"]";
389
390 for(auto const& histoSys : histoSysList) {
391 const std::string histoSysName = histoSys.GetName();
392 RooRealVar* temp = proto.var(("alpha_" + histoSysName).c_str());
393 if(!temp){
394 temp = (RooRealVar*) proto.factory(("alpha_" + histoSysName + range).c_str());
395 }
396 params.add(* temp );
397 }
398
399 return params;
400 }
401
402 /// Create a linear interpolation object that holds nominal and systematics, import it into the workspace,
403 /// and return a pointer to it.
404 RooAbsArg* makeLinInterp(RooArgList const& interpolationParams,
405 RooHistFunc* nominalFunc,
406 RooWorkspace* proto, const std::vector<HistoSys>& histoSysList,
407 const string& prefix,
408 const RooArgList& observables) {
409
410 // now make function that linearly interpolates expectation between variations
411 // get low/high variations to interpolate between
412 vector<double> low, high;
413 RooArgSet lowSet, highSet;
414 //ES// for(unsigned int j=0; j<lowHist.size(); ++j){
415 for(unsigned int j=0; j<histoSysList.size(); ++j){
416 std::stringstream str;
417 str<<"_"<<j;
418
419 const HistoSys& histoSys = histoSysList.at(j);
420 RooDataHist* lowDHist = new RooDataHist((prefix+str.str()+"lowDHist").c_str(),"",observables, histoSys.GetHistoLow());
421 RooDataHist* highDHist = new RooDataHist((prefix+str.str()+"highDHist").c_str(),"",observables, histoSys.GetHistoHigh());
422 RooHistFunc* lowFunc = new RooHistFunc((prefix+str.str()+"low").c_str(),"",observables,*lowDHist,0) ;
423 RooHistFunc* highFunc = new RooHistFunc((prefix+str.str()+"high").c_str(),"",observables,*highDHist,0) ;
424 lowSet.add(*lowFunc);
425 highSet.add(*highFunc);
426 }
427
428 // this is sigma(params), a piece-wise linear interpolation
429 PiecewiseInterpolation interp(prefix.c_str(),"",*nominalFunc,lowSet,highSet,interpolationParams);
430 interp.setPositiveDefinite();
431 interp.setAllInterpCodes(4); // LM: change to 4 (piece-wise linear to 6th order polynomial interpolation + linear extrapolation )
432 // KC: interpo codes 1 etc. don't have proper analytic integral.
433 RooArgSet observableSet(observables);
434 interp.setBinIntegrator(observableSet);
435 interp.forceNumInt();
436
437 proto->import(interp, RecycleConflictNodes()); // individual params have already been imported in first loop of this function
438
439 auto interpInWS = proto->arg(prefix.c_str());
440 assert(interpInWS);
441
442 return interpInWS;
443 }
444
445 }
446
447 // GHL: Consider passing the NormFactor list instead of the entire sample
448 std::unique_ptr<RooProduct> HistoToWorkspaceFactoryFast::CreateNormFactor(RooWorkspace* proto, string& channel, string& sigmaEpsilon, Sample& sample, bool doRatio){
449
450 std::vector<string> prodNames;
451
452 vector<NormFactor> normList = sample.GetNormFactorList();
453 vector<string> normFactorNames, rangeNames;
454
455
456 string overallNorm_times_sigmaEpsilon = sample.GetName() + "_" + channel + "_scaleFactors";
457 auto sigEps = proto->arg(sigmaEpsilon.c_str());
458 assert(sigEps);
459 auto normFactor = std::make_unique<RooProduct>(overallNorm_times_sigmaEpsilon.c_str(), overallNorm_times_sigmaEpsilon.c_str(), RooArgList(*sigEps));
460
461 if(normList.size() > 0){
462
463 for(vector<NormFactor>::iterator itr = normList.begin(); itr != normList.end(); ++itr){
464
465 NormFactor& norm = *itr;
466
467 string varname;
468 if(doRatio) {
469 varname = norm.GetName() + "_" + channel;
470 }
471 else {
472 varname=norm.GetName();
473 }
474
475 // GHL: Check that the NormFactor doesn't already exist
476 // (it may have been created as a function expression
477 // during preprocessing)
478 std::stringstream range;
479 range << "[" << norm.GetVal() << "," << norm.GetLow() << "," << norm.GetHigh() << "]";
480
481 if( proto->obj(varname.c_str()) == NULL) {
482 cxcoutI(HistFactory) << "making normFactor: " << norm.GetName() << endl;
483 // remove "doRatio" and name can be changed when ws gets imported to the combined model.
484 proto->factory((varname + range.str()).c_str());
485 }
486
487 if(norm.GetConst()) {
488 // proto->var(varname.c_str())->setConstant();
489 // cout <<"setting " << varname << " constant"<<endl;
490 cxcoutW(HistFactory) << "Const attribute to <NormFactor> tag is deprecated, will ignore." <<
491 " Instead, add \n\t<ParamSetting Const=\"True\"> " << varname << " </ParamSetting>\n" <<
492 " to your top-level XML's <Measurement> entry" << endl;
493 }
494 prodNames.push_back(varname);
495 rangeNames.push_back(range.str());
496 normFactorNames.push_back(varname);
497 }
498
499
500 for (const auto& name : prodNames) {
501 auto arg = proto->arg(name.c_str());
502 assert(arg);
503 normFactor->addTerm(arg);
504 }
505
506 }
507
508 unsigned int rangeIndex=0;
509 for( vector<string>::iterator nit = normFactorNames.begin(); nit!=normFactorNames.end(); ++nit){
510 if( count (normFactorNames.begin(), normFactorNames.end(), *nit) > 1 ){
511 cxcoutI(HistFactory) <<"<NormFactor Name =\""<<*nit<<"\"> is duplicated for <Sample Name=\""
512 << sample.GetName() << "\">, but only one factor will be included. \n Instead, define something like"
513 << "\n\t<Function Name=\""<<*nit<<"Squared\" Expression=\""<<*nit<<"*"<<*nit<<"\" Var=\""<<*nit<<rangeNames.at(rangeIndex)
514 << "\"> \nin your top-level XML's <Measurment> entry and use <NormFactor Name=\""<<*nit<<"Squared\" in your channel XML file."<< endl;
515 }
516 ++rangeIndex;
517 }
518
519 return normFactor;
520 }
521
523 string interpName,
524 std::vector<OverallSys>& systList,
525 vector<string>& constraintTermNames,
526 vector<string>& totSystTermNames) {
527
528 // add variables for all the relative overall uncertainties we expect
529 // range is set using defined macro (see top of the page)
530
531 string range=string("[0,")+alpha_Low+","+alpha_High+"]";
532 totSystTermNames.push_back(prefix);
533
534 RooArgSet params(prefix.c_str());
535 vector<double> lowVec, highVec;
536
537 std::map<std::string, double>::iterator itconstr;
538 for(unsigned int i = 0; i < systList.size(); ++i) {
539
540 OverallSys& sys = systList.at(i);
541 std::string strname = sys.GetName();
542 const char * name = strname.c_str();
543
544 // case of no systematic (is it possible)
545 if (meas.GetNoSyst().count(sys.GetName()) > 0 ) {
546 cxcoutI(HistFactory) << "HistoToWorkspaceFast::AddConstraintTerm - skip systematic " << sys.GetName() << std::endl;
547 continue;
548 }
549 // case systematic is a gamma constraint
550 if (meas.GetGammaSyst().count(sys.GetName()) > 0 ) {
551 double relerr = meas.GetGammaSyst().find(sys.GetName() )->second;
552 if (relerr <= 0) {
553 cxcoutI(HistFactory) << "HistoToWorkspaceFast::AddConstraintTerm - zero uncertainty assigned - skip systematic " << sys.GetName() << std::endl;
554 continue;
555 }
556 double tauVal = 1./(relerr*relerr);
557 double sqtau = 1./relerr;
558 RooAbsArg * beta = proto->factory(TString::Format("beta_%s[1,0,10]",name) );
559 // the global observable (y_s)
560 RooAbsArg * yvar = proto->factory(TString::Format("nom_%s[%f,0,10]",beta->GetName(),tauVal)) ;
561 // the rate of the gamma distribution (theta)
562 RooAbsArg * theta = proto->factory(TString::Format("theta_%s[%f]",name,1./tauVal));
563 // find alpha as function of beta
564 RooAbsArg* alphaOfBeta = proto->factory(TString::Format("PolyVar::alphaOfBeta_%s(beta_%s,{%f,%f})",name,name,-sqtau,sqtau));
565
566 // add now the constraint itself Gamma_beta_constraint(beta, y+1, tau, 0 )
567 // build the gamma parameter k = as y_s + 1
568 RooAbsArg * kappa = proto->factory(TString::Format("sum::k_%s(%s,1.)",name,yvar->GetName()) );
569 RooAbsArg * gamma = proto->factory(TString::Format("Gamma::%sConstraint(%s, %s, %s, 0.0)",beta->GetName(),beta->GetName(), kappa->GetName(), theta->GetName() ) );
571 alphaOfBeta->Print("t");
572 gamma->Print("t");
573 }
574 constraintTermNames.push_back(gamma->GetName());
575 // set global observables
576 RooRealVar * gobs = dynamic_cast<RooRealVar*>(yvar); assert(gobs);
577 gobs->setConstant(true);
578 const_cast<RooArgSet*>(proto->set("globalObservables"))->add(*yvar);
579
580 // add alphaOfBeta in the list of params to interpolate
581 params.add(*alphaOfBeta);
582 cxcoutI(HistFactory) << "Added a gamma constraint for " << name << std::endl;
583
584 }
585 else {
586 RooRealVar* alpha = (RooRealVar*) proto->var((prefix + sys.GetName()).c_str());
587 if(!alpha) {
588 alpha = (RooRealVar*) proto->factory((prefix + sys.GetName() + range).c_str());
589 }
590 // add the Gaussian constraint part
591 const bool isUniform = meas.GetUniformSyst().count(sys.GetName()) > 0;
592 makeGaussianConstraint(*alpha, *proto, isUniform, constraintTermNames);
593
594 // check if exists a log-normal constraint
595 if (meas.GetLogNormSyst().count(sys.GetName()) == 0 && meas.GetGammaSyst().count(sys.GetName()) == 0 ) {
596 // just add the alpha for the parameters of the FlexibleInterpVar function
597 params.add(*alpha);
598 }
599 // case systematic is a log-normal constraint
600 if (meas.GetLogNormSyst().count(sys.GetName()) > 0 ) {
601 // log normal constraint for parameter
602 double relerr = meas.GetLogNormSyst().find(sys.GetName() )->second;
603 double tauVal = 1./relerr;
604 std::string tauName = "tau_" + sys.GetName();
605 proto->factory(TString::Format("%s[%f]",tauName.c_str(),tauVal ) );
606 double kappaVal = 1. + relerr;
607 std::string kappaName = "kappa_" + sys.GetName();
608 proto->factory(TString::Format("%s[%f]",kappaName.c_str(),kappaVal ) );
609 const char * alphaName = alpha->GetName();
610
611 std::string alphaOfBetaName = "alphaOfBeta_" + sys.GetName();
612 RooAbsArg * alphaOfBeta = proto->factory(TString::Format("expr::%s('%s*(pow(%s,%s)-1.)',%s,%s,%s)",alphaOfBetaName.c_str(),
613 tauName.c_str(),kappaName.c_str(),alphaName,
614 tauName.c_str(),kappaName.c_str(),alphaName ) );
615
616 cxcoutI(HistFactory) << "Added a log-normal constraint for " << name << std::endl;
617 if (RooMsgService::instance().isActive(static_cast<TObject*>(nullptr), RooFit::HistFactory, RooFit::DEBUG))
618 alphaOfBeta->Print("t");
619 params.add(*alphaOfBeta);
620 }
621
622 }
623 // add low/high vectors
624 double low = sys.GetLow();
625 double high = sys.GetHigh();
626 lowVec.push_back(low);
627 highVec.push_back(high);
628
629 } // end sys loop
630
631 if(systList.size() > 0){
632 // this is epsilon(alpha_j), a piece-wise linear interpolation
633 // LinInterpVar interp( (interpName).c_str(), "", params, 1., lowVec, highVec);
634
635 assert( params.getSize() > 0);
636 assert(int(lowVec.size()) == params.getSize() );
637
638 FlexibleInterpVar interp( (interpName).c_str(), "", params, 1., lowVec, highVec);
639 interp.setAllInterpCodes(4); // LM: change to 4 (piece-wise linear to 6th order polynomial interpolation + linear extrapolation )
640 //interp.setAllInterpCodes(0); // simple linear interpolation
641 proto->import(interp); // params have already been imported in first loop of this function
642 } else{
643 // some strange behavior if params,lowVec,highVec are empty.
644 //cout << "WARNING: No OverallSyst terms" << endl;
645 RooConstVar interp( (interpName).c_str(), "", 1.);
646 proto->import(interp); // params have already been imported in first loop of this function
647 }
648
649 // std::cout << "after creating FlexibleInterpVar " << std::endl;
650 // proto->Print();
651
652 }
653
654
656 const vector<RooProduct*>& sampleScaleFactors, std::vector<vector<RooAbsArg*>>& sampleHistFuncs) const {
657 assert(sampleScaleFactors.size() == sampleHistFuncs.size());
658
659 // for ith bin calculate totN_i = lumi * sum_j expected_j * syst_j
660
661 if (fObsNameVec.empty() && !fObsName.empty())
662 throw std::logic_error("HistFactory didn't process the observables correctly. Please file a bug report.");
663
664 auto firstHistFunc = dynamic_cast<const RooHistFunc*>(sampleHistFuncs.front().front());
665 if (!firstHistFunc) {
666 auto piecewiseInt = dynamic_cast<const PiecewiseInterpolation*>(sampleHistFuncs.front().front());
667 firstHistFunc = dynamic_cast<const RooHistFunc*>(piecewiseInt->nominalHist());
668 }
669 assert(firstHistFunc);
670
671 // Prepare a function to divide all bin contents by bin width to get a density:
672 const std::string binWidthFunctionName = totName + "_binWidth";
673 RooBinWidthFunction binWidth(binWidthFunctionName.c_str(), "Divide by bin width to obtain probability density", *firstHistFunc, true);
674 proto->import(binWidth);
675 auto binWidthWS = proto->function(binWidthFunctionName.c_str());
676 assert(binWidthWS);
677
678 // Loop through samples and create products of their functions:
679 RooArgSet coefList;
680 RooArgSet shapeList;
681 for (unsigned int i=0; i < sampleHistFuncs.size(); ++i) {
682 assert(!sampleHistFuncs[i].empty());
683 coefList.add(*sampleScaleFactors[i]);
684
685 std::vector<RooAbsArg*>& thisSampleHistFuncs = sampleHistFuncs[i];
686 thisSampleHistFuncs.push_back(binWidthWS);
687
688 if (thisSampleHistFuncs.size() == 1) {
689 // Just one function. Book it.
690 shapeList.add(*thisSampleHistFuncs.front());
691 } else {
692 // Have multiple functions. We need to multiply them.
693 std::string name = thisSampleHistFuncs.front()->GetName();
694 auto pos = name.find("Hist_alpha");
695 if (pos != std::string::npos) {
696 name = name.substr(0, pos) + "shapes";
697 } else if ( (pos = name.find("nominal")) != std::string::npos) {
698 name = name.substr(0, pos) + "shapes";
699 }
700
701 RooProduct shapeProduct(name.c_str(), name.c_str(), RooArgSet(thisSampleHistFuncs.begin(), thisSampleHistFuncs.end()));
702 proto->import(shapeProduct, RecycleConflictNodes());
703 shapeList.add(*proto->function(name.c_str()));
704 }
705 }
706
707 // Sum all samples
708 RooRealSumPdf tot(totName.c_str(), totName.c_str(), shapeList, coefList, true);
709 tot.specialIntegratorConfig(kTRUE)->method1D().setLabel("RooBinIntegrator") ;
710 tot.specialIntegratorConfig(kTRUE)->method2D().setLabel("RooBinIntegrator") ;
711 tot.specialIntegratorConfig(kTRUE)->methodND().setLabel("RooBinIntegrator") ;
712 tot.forceNumInt();
713
714 // for mixed generation in RooSimultaneous
715 tot.setAttribute("GenerateBinned"); // for use with RooSimultaneous::generate in mixed mode
716
717 proto->import(tot, RecycleConflictNodes());
718 }
719
720 void HistoToWorkspaceFactoryFast::AddPoissonTerms(RooWorkspace* proto, string prefix, string obsPrefix, string expPrefix, int lowBin, int highBin,
721 vector<string>& likelihoodTermNames){
722 /////////////////////////////////
723 // Relate observables to expected for each bin
724 // later modify variable named expPrefix_i to be product of terms
725 RooArgSet Pois(prefix.c_str());
726 for(Int_t i=lowBin; i<highBin; ++i){
727 std::stringstream str;
728 str<<"_"<<i;
729 //string command("Poisson::"+prefix+str.str()+"("+obsPrefix+str.str()+","+expPrefix+str.str()+")");
730 string command("Poisson::"+prefix+str.str()+"("+obsPrefix+str.str()+","+expPrefix+str.str()+",1)");//for no rounding
731 RooAbsArg* temp = (proto->factory( command.c_str() ) );
732
733 // output
734 cout << "Poisson Term " << command << endl;
736 //cout << temp << endl;
737
738 likelihoodTermNames.push_back( temp->GetName() );
739 Pois.add(* temp );
740 }
741 proto->defineSet(prefix.c_str(),Pois); // add argset to workspace
742 }
743
744 void HistoToWorkspaceFactoryFast::SetObsToExpected(RooWorkspace* proto, string obsPrefix, string expPrefix, int lowBin, int highBin){
745 /////////////////////////////////
746 // set observed to expected
747 TTree* tree = new TTree();
748 Double_t* obsForTree = new Double_t[highBin-lowBin];
749 RooArgList obsList("obsList");
750
751 for(Int_t i=lowBin; i<highBin; ++i){
752 std::stringstream str;
753 str<<"_"<<i;
754 RooRealVar* obs = (RooRealVar*) proto->var((obsPrefix+str.str()).c_str());
755 cout << "expected number of events called: " << expPrefix << endl;
756 RooAbsReal* exp = proto->function((expPrefix+str.str()).c_str());
757 if(obs && exp){
758
759 //proto->Print();
760 obs->setVal( exp->getVal() );
761 cout << "setting obs"+str.str()+" to expected = " << exp->getVal() << " check: " << obs->getVal() << endl;
762
763 // add entry to array and attach to tree
764 obsForTree[i] = exp->getVal();
765 tree->Branch((obsPrefix+str.str()).c_str(), obsForTree+i ,(obsPrefix+str.str()+"/D").c_str());
766 obsList.add(*obs);
767 }else{
768 cout << "problem retrieving obs or exp " << obsPrefix+str.str() << obs << " " << expPrefix+str.str() << exp << endl;
769 }
770 }
771 tree->Fill();
772 RooDataSet* data = new RooDataSet("expData","", tree, obsList); // one experiment
773
774 delete tree;
775 delete [] obsForTree;
776
777 proto->import(*data);
778
779 delete data;
780
781 }
782
783 //////////////////////////////////////////////////////////////////////////////
784
786 map<string,double> gammaSyst,
787 map<string,double> uniformSyst,
788 map<string,double> logNormSyst,
789 map<string,double> noSyst) {
790 string pdfName(pdfNameChar);
791
792 ModelConfig * combined_config = (ModelConfig *) proto->obj("ModelConfig");
793 if( combined_config==NULL ) {
794 std::cout << "Error: Failed to find object 'ModelConfig' in workspace: "
795 << proto->GetName() << std::endl;
796 throw hf_exc();
797 }
798 // const RooArgSet * constrainedParams=combined_config->GetNuisanceParameters();
799 // RooArgSet temp(*constrainedParams);
800 string edit="EDIT::newSimPdf("+pdfName+",";
801 string editList;
802 string lastPdf=pdfName;
803 string precede="";
804 unsigned int numReplacements = 0;
805 unsigned int nskipped = 0;
806 map<string,double>::iterator it;
807
808
809 // add gamma terms and their constraints
810 for(it=gammaSyst.begin(); it!=gammaSyst.end(); ++it) {
811 //cout << "edit for " << it->first << "with rel uncert = " << it->second << endl;
812 if(! proto->var(("alpha_"+it->first).c_str())){
813 //cout << "systematic not there" << endl;
814 nskipped++;
815 continue;
816 }
817 numReplacements++;
818
819 double relativeUncertainty = it->second;
820 double scale = 1/sqrt((1+1/pow(relativeUncertainty,2)));
821
822 // this is the Gamma PDF and in a form that doesn't have roundoff problems like the Poisson does
823 proto->factory(Form("beta_%s[1,0,10]",it->first.c_str()));
824 proto->factory(Form("y_%s[%f]",it->first.c_str(),1./pow(relativeUncertainty,2))) ;
825 proto->factory(Form("theta_%s[%f]",it->first.c_str(),pow(relativeUncertainty,2))) ;
826 proto->factory(Form("Gamma::beta_%sConstraint(beta_%s,sum::k_%s(y_%s,one[1]),theta_%s,zero[0])",
827 it->first.c_str(),
828 it->first.c_str(),
829 it->first.c_str(),
830 it->first.c_str(),
831 it->first.c_str())) ;
832
833 /*
834 // this has some problems because N in poisson is rounded to nearest integer
835 proto->factory(Form("Poisson::beta_%sConstraint(y_%s[%f],prod::taub_%s(taus_%s[%f],beta_%s[1,0,5]))",
836 it->first.c_str(),
837 it->first.c_str(),
838 1./pow(relativeUncertainty,2),
839 it->first.c_str(),
840 it->first.c_str(),
841 1./pow(relativeUncertainty,2),
842 it->first.c_str()
843 ) ) ;
844 */
845 // combined->factory(Form("expr::alphaOfBeta('(beta-1)/%f',beta)",scale));
846 // combined->factory(Form("expr::alphaOfBeta_%s('(beta_%s-1)/%f',beta_%s)",it->first.c_str(),it->first.c_str(),scale,it->first.c_str()));
847 proto->factory(Form("PolyVar::alphaOfBeta_%s(beta_%s,{%f,%f})",it->first.c_str(),it->first.c_str(),-1./scale,1./scale));
848
849 // set beta const status to be same as alpha
850 if(proto->var(Form("alpha_%s",it->first.c_str()))->isConstant()) {
851 proto->var(Form("beta_%s",it->first.c_str()))->setConstant(true);
852 }
853 else {
854 proto->var(Form("beta_%s",it->first.c_str()))->setConstant(false);
855 }
856 // set alpha const status to true
857 // proto->var(Form("alpha_%s",it->first.c_str()))->setConstant(true);
858
859 // replace alphas with alphaOfBeta and replace constraints
860 editList+=precede + "alpha_"+it->first+"Constraint=beta_" + it->first+ "Constraint";
861 precede=",";
862 editList+=precede + "alpha_"+it->first+"=alphaOfBeta_"+ it->first;
863
864 /*
865 if( proto->pdf(("alpha_"+it->first+"Constraint").c_str()) && proto->var(("alpha_"+it->first).c_str()) )
866 cout << " checked they are there" << proto->pdf(("alpha_"+it->first+"Constraint").c_str()) << " " << proto->var(("alpha_"+it->first).c_str()) << endl;
867 else
868 cout << "NOT THERE" << endl;
869 */
870
871 // EDIT seems to die if the list of edits is too long. So chunck them up.
872 if(numReplacements%10 == 0 && numReplacements+nskipped!=gammaSyst.size()){
873 edit="EDIT::"+lastPdf+"_("+lastPdf+","+editList+")";
874 lastPdf+="_"; // append an underscore for the edit
875 editList=""; // reset edit list
876 precede="";
877 cout << "Going to issue this edit command\n" << edit<< endl;
878 proto->factory( edit.c_str() );
879 RooAbsPdf* newOne = proto->pdf(lastPdf.c_str());
880 if(!newOne)
881 cxcoutWHF << "---------------------\n WARNING: failed to make EDIT\n\n" << endl;
882
883 }
884 }
885
886 // add uniform terms and their constraints
887 for(it=uniformSyst.begin(); it!=uniformSyst.end(); ++it) {
888 cout << "edit for " << it->first << "with rel uncert = " << it->second << endl;
889 if(! proto->var(("alpha_"+it->first).c_str())){
890 cout << "systematic not there" << endl;
891 nskipped++;
892 continue;
893 }
894 numReplacements++;
895
896 // this is the Uniform PDF
897 proto->factory(Form("beta_%s[1,0,10]",it->first.c_str()));
898 proto->factory(Form("Uniform::beta_%sConstraint(beta_%s)",it->first.c_str(),it->first.c_str()));
899 proto->factory(Form("PolyVar::alphaOfBeta_%s(beta_%s,{-1,1})",it->first.c_str(),it->first.c_str()));
900
901 // set beta const status to be same as alpha
902 if(proto->var(Form("alpha_%s",it->first.c_str()))->isConstant())
903 proto->var(Form("beta_%s",it->first.c_str()))->setConstant(true);
904 else
905 proto->var(Form("beta_%s",it->first.c_str()))->setConstant(false);
906 // set alpha const status to true
907 // proto->var(Form("alpha_%s",it->first.c_str()))->setConstant(true);
908
909 // replace alphas with alphaOfBeta and replace constraints
910 cout << "alpha_"+it->first+"Constraint=beta_" + it->first+ "Constraint" << endl;
911 editList+=precede + "alpha_"+it->first+"Constraint=beta_" + it->first+ "Constraint";
912 precede=",";
913 cout << "alpha_"+it->first+"=alphaOfBeta_"+ it->first << endl;
914 editList+=precede + "alpha_"+it->first+"=alphaOfBeta_"+ it->first;
915
916 if( proto->pdf(("alpha_"+it->first+"Constraint").c_str()) && proto->var(("alpha_"+it->first).c_str()) )
917 cout << " checked they are there" << proto->pdf(("alpha_"+it->first+"Constraint").c_str()) << " " << proto->var(("alpha_"+it->first).c_str()) << endl;
918 else
919 cout << "NOT THERE" << endl;
920
921 // EDIT seems to die if the list of edits is too long. So chunck them up.
922 if(numReplacements%10 == 0 && numReplacements+nskipped!=gammaSyst.size()){
923 edit="EDIT::"+lastPdf+"_("+lastPdf+","+editList+")";
924 lastPdf+="_"; // append an underscore for the edit
925 editList=""; // reset edit list
926 precede="";
927 cout << edit<< endl;
928 proto->factory( edit.c_str() );
929 RooAbsPdf* newOne = proto->pdf(lastPdf.c_str());
930 if(!newOne)
931 cxcoutWHF << "---------------------\n WARNING: failed to make EDIT\n\n" << endl;
932
933 }
934 }
935
936 /////////////////////////////////////////
937 ////////////////////////////////////
938
939
940 // add lognormal terms and their constraints
941 for(it=logNormSyst.begin(); it!=logNormSyst.end(); ++it) {
942 cout << "edit for " << it->first << "with rel uncert = " << it->second << endl;
943 if(! proto->var(("alpha_"+it->first).c_str())){
944 cout << "systematic not there" << endl;
945 nskipped++;
946 continue;
947 }
948 numReplacements++;
949
950 double relativeUncertainty = it->second;
951 double kappa = 1+relativeUncertainty;
952 // when transforming beta -> alpha, need alpha=1 to be +1sigma value.
953 // the P(beta>kappa*\hat(beta)) = 16%
954 // and \hat(beta) is 1, thus
955 double scale = relativeUncertainty;
956 //double scale = kappa;
957
958 const char * cname = it->first.c_str();
959
960 // this is the LogNormal
961 proto->factory(TString::Format("beta_%s[1,0,10]",cname));
962 proto->factory(TString::Format("nom_beta_%s[1]",cname));
963 proto->factory(TString::Format("kappa_%s[%f]",cname,kappa));
964 proto->factory(TString::Format("Lognormal::beta_%sConstraint(beta_%s,nom_beta_%s,kappa_%s)",
965 cname, cname, cname, cname)) ;
966 proto->factory(TString::Format("PolyVar::alphaOfBeta_%s(beta_%s,{%f,%f})",cname,cname,-1./scale,1./scale));
967
968
969 // set beta const status to be same as alpha
970 if(proto->var(TString::Format("alpha_%s",cname))->isConstant())
971 proto->var(TString::Format("beta_%s",cname))->setConstant(true);
972 else
973 proto->var(TString::Format("beta_%s",cname))->setConstant(false);
974 // set alpha const status to true
975 // proto->var(TString::Format("alpha_%s",cname))->setConstant(true);
976
977 // replace alphas with alphaOfBeta and replace constraints
978 cout << "alpha_"+it->first+"Constraint=beta_" + it->first+ "Constraint" << endl;
979 editList+=precede + "alpha_"+it->first+"Constraint=beta_" + it->first+ "Constraint";
980 precede=",";
981 cout << "alpha_"+it->first+"=alphaOfBeta_"+ it->first << endl;
982 editList+=precede + "alpha_"+it->first+"=alphaOfBeta_"+ it->first;
983
984 if( proto->pdf(("alpha_"+it->first+"Constraint").c_str()) && proto->var(("alpha_"+it->first).c_str()) )
985 cout << " checked they are there" << proto->pdf(("alpha_"+it->first+"Constraint").c_str()) << " " << proto->var(("alpha_"+it->first).c_str()) << endl;
986 else
987 cout << "NOT THERE" << endl;
988
989 // EDIT seems to die if the list of edits is too long. So chunck them up.
990 if(numReplacements%10 == 0 && numReplacements+nskipped!=gammaSyst.size()){
991 edit="EDIT::"+lastPdf+"_("+lastPdf+","+editList+")";
992 lastPdf+="_"; // append an underscore for the edit
993 editList=""; // reset edit list
994 precede="";
995 cout << edit<< endl;
996 proto->factory( edit.c_str() );
997 RooAbsPdf* newOne = proto->pdf(lastPdf.c_str());
998 if(!newOne)
999 cxcoutWHF << "\n\n ---------------------\n WARNING: failed to make EDIT\n\n" << endl;
1000
1001 }
1002 // add global observables
1003 const RooArgSet * gobs = proto->set("globalObservables");
1004 RooArgSet gobsNew(*gobs);
1005 gobsNew.add(*proto->var(TString::Format("nom_beta_%s",cname)) );
1006 proto->removeSet("globalObservables");
1007 proto->defineSet("globalObservables",gobsNew);
1008 gobsNew.Print();
1009
1010 }
1011
1012 /////////////////////////////////////////
1013
1014 // MB: remove a systematic constraint
1015 for(it=noSyst.begin(); it!=noSyst.end(); ++it) {
1016
1017 cout << "remove constraint for parameter" << it->first << endl;
1018 if(! proto->var(("alpha_"+it->first).c_str()) || ! proto->pdf(("alpha_"+it->first+"Constraint").c_str()) ) {
1019 cout << "systematic not there" << endl;
1020 nskipped++;
1021 continue;
1022 }
1023 numReplacements++;
1024
1025 // dummy replacement pdf
1026 if ( !proto->var("one") ) { proto->factory("one[1.0]"); }
1027 proto->var("one")->setConstant();
1028
1029 // replace constraints
1030 cout << "alpha_"+it->first+"Constraint=one" << endl;
1031 editList+=precede + "alpha_"+it->first+"Constraint=one";
1032 precede=",";
1033
1034 // EDIT seems to die if the list of edits is too long. So chunck them up.
1035 if(numReplacements%10 == 0 && numReplacements+nskipped!=gammaSyst.size()){
1036 edit="EDIT::"+lastPdf+"_("+lastPdf+","+editList+")";
1037 lastPdf+="_"; // append an underscore for the edit
1038 editList=""; // reset edit list
1039 precede="";
1040 cout << edit << endl;
1041 proto->factory( edit.c_str() );
1042 RooAbsPdf* newOne = proto->pdf(lastPdf.c_str());
1043 if(!newOne) {
1044 cxcoutWHF << "---------------------\n WARNING: failed to make EDIT\n\n" << endl;
1045 }
1046 }
1047 }
1048
1049 /////////////////////////////////////////
1050
1051 // commit last bunch of edits
1052 edit="EDIT::newSimPdf("+lastPdf+","+editList+")";
1053 cout << edit<< endl;
1054 proto->factory( edit.c_str() );
1055 // proto->writeToFile(("results/model_"+fRowTitle+"_edited.root").c_str());
1056 RooAbsPdf* newOne = proto->pdf("newSimPdf");
1057 if(newOne){
1058 // newOne->graphVizTree(("results/"+pdfName+"_"+fRowTitle+"newSimPdf.dot").c_str());
1059 combined_config->SetPdf(*newOne);
1060 }
1061 else{
1062 cxcoutWHF << "\n\n ---------------------\n WARNING: failed to make EDIT\n\n" << endl;
1063 }
1064 }
1065
1067 // Change-> Now a static utility
1068
1069 FILE* covFile = fopen ((filename).c_str(),"w");
1070
1071 TIter iti = params->createIterator();
1072 TIter itj = params->createIterator();
1073 RooRealVar *myargi, *myargj;
1074 fprintf(covFile," ") ;
1075 while ((myargi = (RooRealVar *)iti.Next())) {
1076 if(myargi->isConstant()) continue;
1077 fprintf(covFile," & %s", myargi->GetName());
1078 }
1079 fprintf(covFile,"\\\\ \\hline \n" );
1080 iti.Reset();
1081 while ((myargi = (RooRealVar *)iti.Next())) {
1082 if(myargi->isConstant()) continue;
1083 fprintf(covFile,"%s", myargi->GetName());
1084 itj.Reset();
1085 while ((myargj = (RooRealVar *)itj.Next())) {
1086 if(myargj->isConstant()) continue;
1087 cout << myargi->GetName() << "," << myargj->GetName();
1088 fprintf(covFile, " & %.2f", result->correlation(*myargi, *myargj));
1089 }
1090 cout << endl;
1091 fprintf(covFile, " \\\\\n");
1092 }
1093 fclose(covFile);
1094
1095 }
1096
1097
1098 ///////////////////////////////////////////////
1100
1101 // check inputs (see JIRA-6890 )
1102
1103 if (channel.GetSamples().empty()) {
1104 Error("MakeSingleChannelWorkspace",
1105 "The input Channel does not contain any sample - return a nullptr");
1106 return 0;
1107 }
1108
1109 const TH1* channel_hist_template = channel.GetSamples().front().GetHisto();
1110 if (channel_hist_template == nullptr) {
1111 channel.CollectHistograms();
1112 channel_hist_template = channel.GetSamples().front().GetHisto();
1113 }
1114 if (channel_hist_template == nullptr) {
1115 std::ostringstream stream;
1116 stream << "The sample " << channel.GetSamples().front().GetName()
1117 << " in channel " << channel.GetName() << " does not contain a histogram. This is the channel:\n";
1118 channel.Print(stream);
1119 Error("MakeSingleChannelWorkspace", "%s", stream.str().c_str());
1120 return 0;
1121 }
1122
1123 if( ! channel.CheckHistograms() ) {
1124 std::cout << "MakeSingleChannelWorkspace: Channel: " << channel.GetName()
1125 << " has uninitialized histogram pointers" << std::endl;
1126 throw hf_exc();
1127 }
1128
1129
1130
1131 // Set these by hand inside the function
1132 vector<string> systToFix = measurement.GetConstantParams();
1133 bool doRatio=false;
1134
1135 // to time the macro
1136 TStopwatch t;
1137 t.Start();
1138 //ES// string channel_name=summary[0].channel;
1139 string channel_name = channel.GetName();
1140
1141 /// MB: reset observable names for each new channel.
1142 fObsNameVec.clear();
1143
1144 /// MB: label observables x,y,z, depending on histogram dimensionality
1145 /// GHL: Give it the first sample's nominal histogram as a template
1146 /// since the data histogram may not be present
1147 if (fObsNameVec.empty()) { GuessObsNameVec(channel_hist_template); }
1148
1149 for ( unsigned int idx=0; idx<fObsNameVec.size(); ++idx ) {
1150 fObsNameVec[idx] = "obs_" + fObsNameVec[idx] + "_" + channel_name ;
1151 }
1152
1153 if (fObsNameVec.empty()) {
1154 fObsName= "obs_" + channel_name; // set name ov observable
1155 fObsNameVec.push_back( fObsName );
1156 }
1157
1158 if (fObsNameVec.empty() || fObsNameVec.size() > 3) {
1159 throw hf_exc("HistFactory is limited to 1- to 3-dimensional histograms.");
1160 }
1161
1162 cxcoutP(HistFactory) << "\n-----------------------------------------\n"
1163 << "\tStarting to process '"
1164 << channel_name << "' channel with " << fObsNameVec.size() << " observables"
1165 << "\n-----------------------------------------\n" << endl;
1166
1167 //
1168 // our main workspace that we are using to construct the model
1169 //
1170 RooWorkspace* proto = new RooWorkspace(channel_name.c_str(), (channel_name+" workspace").c_str());
1171 auto proto_config = make_unique<ModelConfig>("ModelConfig", proto);
1172 proto_config->SetWorkspace(*proto);
1173
1174 // preprocess functions
1175 vector<string>::iterator funcIter = fPreprocessFunctions.begin();
1176 for(;funcIter!= fPreprocessFunctions.end(); ++funcIter){
1177 cxcoutI(HistFactory) << "will preprocess this line: " << *funcIter <<endl;
1178 proto->factory(funcIter->c_str());
1179 proto->Print();
1180 }
1181
1182 RooArgSet likelihoodTerms("likelihoodTerms"), constraintTerms("constraintTerms");
1183 vector<string> likelihoodTermNames, constraintTermNames, totSystTermNames;
1184 // All histogram functions to be multiplied in each sample
1185 std::vector<std::vector<RooAbsArg*>> allSampleHistFuncs;
1186 std::vector<RooProduct*> sampleScaleFactors;
1187
1188 std::vector< pair<string,string> > statNamePairs;
1189 std::vector< pair<const TH1*, std::unique_ptr<TH1>> > statHistPairs; // <nominal, error>
1190 const std::string statFuncName = "mc_stat_" + channel_name;
1191
1192 string prefix, range;
1193
1194 /////////////////////////////
1195 // shared parameters
1196 // this is ratio of lumi to nominal lumi. We will include relative uncertainty in model
1197 std::stringstream lumiStr;
1198 // lumi range
1199 lumiStr << "Lumi[" << fNomLumi << ",0," << 10.*fNomLumi << "]";
1200 proto->factory(lumiStr.str().c_str());
1201 cxcoutI(HistFactory) << "lumi str = " << lumiStr.str() << endl;
1202
1203 std::stringstream lumiErrorStr;
1204 lumiErrorStr << "nominalLumi["<<fNomLumi << ",0,"<<fNomLumi+10*fLumiError<<"]," << fLumiError ;
1205 proto->factory(("Gaussian::lumiConstraint(Lumi,"+lumiErrorStr.str()+")").c_str());
1206 proto->var("nominalLumi")->setConstant();
1207 proto->defineSet("globalObservables","nominalLumi");
1208 //likelihoodTermNames.push_back("lumiConstraint");
1209 constraintTermNames.push_back("lumiConstraint");
1210 cxcoutI(HistFactory) << "lumi Error str = " << lumiErrorStr.str() << endl;
1211
1212 //proto->factory((string("SigXsecOverSM[1.,0.5,1..8]").c_str()));
1213 ///////////////////////////////////
1214 // loop through estimates, add expectation, floating bin predictions,
1215 // and terms that constrain floating to expectation via uncertainties
1216 // GHL: Loop over samples instead, which doesn't contain the data
1217 for (Sample& sample : channel.GetSamples()) {
1218 string overallSystName = sample.GetName() + "_" + channel_name + "_epsilon";
1219
1220 string systSourcePrefix = "alpha_";
1221
1222 // constraintTermNames and totSystTermNames are vectors that are passed
1223 // by reference and filled by this method
1224 AddConstraintTerms(proto,measurement, systSourcePrefix, overallSystName,
1225 sample.GetOverallSysList(), constraintTermNames , totSystTermNames);
1226
1227 allSampleHistFuncs.emplace_back();
1228 std::vector<RooAbsArg*>& sampleHistFuncs = allSampleHistFuncs.back();
1229
1230 // GHL: Consider passing the NormFactor list instead of the entire sample
1231 auto normFactors = CreateNormFactor(proto, channel_name, overallSystName, sample, doRatio);
1232 assert(normFactors);
1233
1234 // Create the string for the object
1235 // that is added to the RooRealSumPdf
1236 // for this channel
1237// string syst_x_expectedPrefix = "";
1238
1239 // get histogram
1240 //ES// TH1* nominal = it->nominal;
1241 const TH1* nominal = sample.GetHisto();
1242
1243 // MB : HACK no option to have both non-hist variations and hist variations ?
1244 // get histogram
1245 // GHL: Okay, this is going to be non-trivial.
1246 // We will loop over histosys's, which contain both
1247 // the low hist and the high hist together.
1248
1249 // Logic:
1250 // - If we have no HistoSys's, do part A
1251 // - else, if the histo syst's don't match, return (we ignore this case)
1252 // - finally, we take the syst's and apply the linear interpolation w/ constraint
1253 string expPrefix = sample.GetName() + "_" + channel_name;
1254 // create roorealvar observables
1255 RooArgList observables = createObservables(sample.GetHisto(), proto);
1256 RooHistFunc* nominalHistFunc = MakeExpectedHistFunc(sample.GetHisto(), proto, expPrefix, observables);
1257 assert(nominalHistFunc);
1258
1259 if(sample.GetHistoSysList().size() == 0) {
1260 // If no HistoSys
1261 cxcoutI(HistFactory) << sample.GetName() + "_" + channel_name + " has no variation histograms " << endl;
1262
1263 sampleHistFuncs.push_back(nominalHistFunc);
1264 } else {
1265 // If there ARE HistoSys(s)
1266 // name of source for variation
1267 string constraintPrefix = sample.GetName() + "_" + channel_name + "_Hist_alpha";
1268
1269 // make list of abstract parameters that interpolate in space of variations
1270 RooArgList interpParams = makeInterpolationParameters(sample.GetHistoSysList(), *proto);
1271
1272 // next, cerate the constraint terms
1273 for(std::size_t i = 0; i < interpParams.size(); ++i) {
1274 bool isUniform = measurement.GetUniformSyst().count(sample.GetHistoSysList()[i].GetName()) > 0;
1275 makeGaussianConstraint(interpParams[i], *proto, isUniform, constraintTermNames);
1276 }
1277
1278 // finally, create the interpolated function
1279 sampleHistFuncs.push_back( makeLinInterp(interpParams, nominalHistFunc, proto,
1280 sample.GetHistoSysList(), constraintPrefix, observables) );
1281 }
1282
1283 ////////////////////////////////////
1284 // Add StatErrors to this Channel //
1285 ////////////////////////////////////
1286
1287 if( sample.GetStatError().GetActivate() ) {
1288
1289 if( fObsNameVec.size() > 3 ) {
1290 cxcoutF(HistFactory) << "Cannot include Stat Error for histograms of more than 3 dimensions."
1291 << std::endl;
1292 throw hf_exc();
1293 } else {
1294
1295 // If we are using StatUncertainties, we multiply this object
1296 // by the ParamHistFunc and then pass that to the
1297 // RooRealSumPdf by appending it's name to the list
1298
1299 cxcoutI(HistFactory) << "Sample: " << sample.GetName() << " to be included in Stat Error "
1300 << "for channel " << channel_name
1301 << std::endl;
1302
1303 string UncertName = sample.GetName() + "_" + channel_name + "_StatAbsolUncert";
1304 std::unique_ptr<TH1> statErrorHist;
1305
1306 if( sample.GetStatError().GetErrorHist() == NULL ) {
1307 // Make the absolute stat error
1308 cxcoutI(HistFactory) << "Making Statistical Uncertainty Hist for "
1309 << " Channel: " << channel_name
1310 << " Sample: " << sample.GetName()
1311 << std::endl;
1312 statErrorHist.reset(MakeAbsolUncertaintyHist( UncertName, nominal));
1313 } else {
1314 // clone the error histograms because in case the sample has not error hist
1315 // it is created in MakeAbsolUncertainty
1316 // we need later to clean statErrorHist
1317 statErrorHist.reset(static_cast<TH1*>(sample.GetStatError().GetErrorHist()->Clone()));
1318 // We assume the (relative) error is provided.
1319 // We must turn it into an absolute error
1320 // using the nominal histogram
1321 cxcoutI(HistFactory) << "Using external histogram for Stat Errors for "
1322 << "\tChannel: " << channel_name
1323 << "\tSample: " << sample.GetName()
1324 << "\tError Histogram: " << statErrorHist->GetName() << std::endl;
1325 // Multiply the relative stat uncertainty by the
1326 // nominal to get the overall stat uncertainty
1327 statErrorHist->Multiply( nominal );
1328 statErrorHist->SetName( UncertName.c_str() );
1329 }
1330
1331 // Save the nominal and error hists
1332 // for the building of constraint terms
1333 statHistPairs.emplace_back(nominal, std::move(statErrorHist));
1334
1335 // To do the 'conservative' version, we would need to do some
1336 // intervention here. We would probably need to create a different
1337 // ParamHistFunc for each sample in the channel. The would nominally
1338 // use the same gamma's, so we haven't increased the number of parameters
1339 // However, if a bin in the 'nominal' histogram is 0, we simply need to
1340 // change the parameter in that bin in the ParamHistFunc for this sample.
1341 // We also need to add a constraint term.
1342 // Actually, we'd probably not use the ParamHistFunc...?
1343 // we could remove the dependence in this ParamHistFunc on the ith gamma
1344 // and then create the poisson term: Pois(tau | n_exp)Pois(data | n_exp)
1345
1346
1347 // Next, try to get the common ParamHistFunc (it may have been
1348 // created by another sample in this channel)
1349 // or create it if it doesn't yet exist:
1350 ParamHistFunc* paramHist = dynamic_cast<ParamHistFunc*>( proto->function(statFuncName.c_str()) );
1351 if( paramHist == nullptr ) {
1352
1353 // Get a RooArgSet of the observables:
1354 // Names in the list fObsNameVec:
1355 RooArgList theObservables;
1356 std::vector<std::string>::iterator itr = fObsNameVec.begin();
1357 for (int idx=0; itr!=fObsNameVec.end(); ++itr, ++idx ) {
1358 theObservables.add( *proto->var(itr->c_str()) );
1359 }
1360
1361 // Create the list of terms to
1362 // control the bin heights:
1363 std::string ParamSetPrefix = "gamma_stat_" + channel_name;
1364 Double_t gammaMin = 0.0;
1365 Double_t gammaMax = 10.0;
1367 ParamSetPrefix.c_str(),
1368 theObservables,
1369 gammaMin, gammaMax);
1370
1371 ParamHistFunc statUncertFunc(statFuncName.c_str(), statFuncName.c_str(),
1372 theObservables, statFactorParams );
1373
1374 proto->import( statUncertFunc, RecycleConflictNodes() );
1375
1376 paramHist = (ParamHistFunc*) proto->function( statFuncName.c_str() );
1377 }
1378
1379 // apply stat function to sample
1380 sampleHistFuncs.push_back(paramHist);
1381 }
1382 } // END: if DoMcStat
1383
1384
1385 ///////////////////////////////////////////
1386 // Create a ShapeFactor for this channel //
1387 ///////////////////////////////////////////
1388
1389 if( sample.GetShapeFactorList().size() > 0 ) {
1390
1391 if( fObsNameVec.size() > 3 ) {
1392 cxcoutF(HistFactory) << "Cannot include Stat Error for histograms of more than 3 dimensions."
1393 << std::endl;
1394 throw hf_exc();
1395 } else {
1396
1397 cxcoutI(HistFactory) << "Sample: " << sample.GetName() << " in channel: " << channel_name
1398 << " to be include a ShapeFactor."
1399 << std::endl;
1400
1401 for(unsigned int i=0; i < sample.GetShapeFactorList().size(); ++i) {
1402
1403 ShapeFactor& shapeFactor = sample.GetShapeFactorList().at(i);
1404
1405 std::string funcName = channel_name + "_" + shapeFactor.GetName() + "_shapeFactor";
1406 ParamHistFunc* paramHist = (ParamHistFunc*) proto->function( funcName.c_str() );
1407 if( paramHist == nullptr ) {
1408
1409 RooArgList theObservables;
1410 std::vector<std::string>::iterator itr = fObsNameVec.begin();
1411 for (int idx=0; itr!=fObsNameVec.end(); ++itr, ++idx ) {
1412 theObservables.add( *proto->var(itr->c_str()) );
1413 }
1414
1415 // Create the Parameters
1416 std::string funcParams = "gamma_" + shapeFactor.GetName();
1417
1418 // GHL: Again, we are putting hard ranges on the gamma's
1419 // We should change this to range from 0 to /inf
1420 RooArgList shapeFactorParams = ParamHistFunc::createParamSet(*proto,
1421 funcParams.c_str(),
1422 theObservables, 0, 1000);
1423
1424 // Create the Function
1425 ParamHistFunc shapeFactorFunc( funcName.c_str(), funcName.c_str(),
1426 theObservables, shapeFactorParams );
1427
1428 // Set an initial shape, if requested
1429 if( shapeFactor.GetInitialShape() != NULL ) {
1430 TH1* initialShape = static_cast<TH1*>(shapeFactor.GetInitialShape()->Clone());
1431 cxcoutI(HistFactory) << "Setting Shape Factor: " << shapeFactor.GetName()
1432 << " to have initial shape from hist: "
1433 << initialShape->GetName()
1434 << std::endl;
1435 shapeFactorFunc.setShape( initialShape );
1436 }
1437
1438 // Set the variables constant, if requested
1439 if( shapeFactor.IsConstant() ) {
1440 cxcoutI(HistFactory) << "Setting Shape Factor: " << shapeFactor.GetName()
1441 << " to be constant" << std::endl;
1442 shapeFactorFunc.setConstant(true);
1443 }
1444
1445 proto->import( shapeFactorFunc, RecycleConflictNodes() );
1446 paramHist = (ParamHistFunc*) proto->function( funcName.c_str() );
1447
1448 } // End: Create ShapeFactor ParamHistFunc
1449
1450 sampleHistFuncs.push_back(paramHist);
1451 } // End loop over ShapeFactor Systematics
1452 }
1453 } // End: if ShapeFactorName!=""
1454
1455
1456 ////////////////////////////////////////
1457 // Create a ShapeSys for this channel //
1458 ////////////////////////////////////////
1459
1460 if( sample.GetShapeSysList().size() != 0 ) {
1461
1462 if( fObsNameVec.size() > 3 ) {
1463 cxcoutF(HistFactory) << "Cannot include Stat Error for histograms of more than 3 dimensions."
1464 << std::endl;
1465 throw hf_exc();
1466 } else {
1467
1468 // List of ShapeSys ParamHistFuncs
1469 std::vector<string> ShapeSysNames;
1470
1471 for( unsigned int i = 0; i < sample.GetShapeSysList().size(); ++i) {
1472
1473 // Create the ParamHistFunc's
1474 // Create their constraint terms and add them
1475 // to the list of constraint terms
1476
1477 // Create a single RooProduct over all of these
1478 // paramHistFunc's
1479
1480 // Send the name of that product to the RooRealSumPdf
1481
1482 RooStats::HistFactory::ShapeSys& shapeSys = sample.GetShapeSysList().at(i);
1483
1484 cxcoutI(HistFactory) << "Sample: " << sample.GetName() << " in channel: " << channel_name
1485 << " to include a ShapeSys." << std::endl;
1486
1487 std::string funcName = channel_name + "_" + shapeSys.GetName() + "_ShapeSys";
1488 ShapeSysNames.push_back( funcName );
1489 ParamHistFunc* paramHist = (ParamHistFunc*) proto->function( funcName.c_str() );
1490 if( paramHist == NULL ) {
1491
1492 //std::string funcParams = "gamma_" + it->shapeFactorName;
1493 //paramHist = CreateParamHistFunc( proto, fObsNameVec, funcParams, funcName );
1494
1495 RooArgList theObservables;
1496 std::vector<std::string>::iterator itr = fObsNameVec.begin();
1497 for(; itr!=fObsNameVec.end(); ++itr ) {
1498 theObservables.add( *proto->var(itr->c_str()) );
1499 }
1500
1501 // Create the Parameters
1502 std::string funcParams = "gamma_" + shapeSys.GetName();
1503 RooArgList shapeFactorParams = ParamHistFunc::createParamSet(*proto,
1504 funcParams.c_str(),
1505 theObservables, 0, 10);
1506
1507 // Create the Function
1508 ParamHistFunc shapeFactorFunc( funcName.c_str(), funcName.c_str(),
1509 theObservables, shapeFactorParams );
1510
1511 proto->import( shapeFactorFunc, RecycleConflictNodes() );
1512 paramHist = (ParamHistFunc*) proto->function( funcName.c_str() );
1513
1514 } // End: Create ShapeFactor ParamHistFunc
1515
1516 // Create the constraint terms and add
1517 // them to the workspace (proto)
1518 // as well as the list of constraint terms (constraintTermNames)
1519
1520 // The syst should be a fractional error
1521 const TH1* shapeErrorHist = shapeSys.GetErrorHist();
1522
1523 // Constraint::Type shapeConstraintType = Constraint::Gaussian;
1524 Constraint::Type systype = shapeSys.GetConstraintType();
1525 if( systype == Constraint::Gaussian) {
1526 systype = Constraint::Gaussian;
1527 }
1528 if( systype == Constraint::Poisson ) {
1529 systype = Constraint::Poisson;
1530 }
1531
1532 Double_t minShapeUncertainty = 0.0;
1533 RooArgList shapeConstraints = createStatConstraintTerms(proto, constraintTermNames,
1534 *paramHist, shapeErrorHist,
1535 systype,
1536 minShapeUncertainty);
1537
1538 } // End: Loop over ShapeSys vector in this EstimateSummary
1539
1540 // Now that we have the list of ShapeSys ParamHistFunc names,
1541 // we create the total RooProduct
1542 // we multiply the expected functio
1543
1544
1545 for( unsigned int i = 0; i < ShapeSysNames.size(); ++i ) {
1546 auto func = proto->function(ShapeSysNames.at(i).c_str());
1547 assert(func);
1548 sampleHistFuncs.push_back(func);
1549 }
1550
1551 } // End: NumObsVar == 1
1552
1553 } // End: GetShapeSysList.size() != 0
1554
1555
1556 // GHL: This was pretty confusing before,
1557 // hopefully using the measurement directly
1558 // will improve it
1559 auto lumi = proto->arg("Lumi");
1560 if( !sample.GetNormalizeByTheory() ) {
1561 if (!lumi) {
1562 TString lumiParamString;
1563 lumiParamString += measurement.GetLumi();
1564 lumiParamString.ReplaceAll(' ', TString());
1565 lumi = proto->factory(("Lumi[" + lumiParamString + "]").Data());
1566 } else {
1567 static_cast<RooAbsRealLValue*>(lumi)->setVal(measurement.GetLumi());
1568 }
1569 }
1570 assert(lumi);
1571 normFactors->addTerm(lumi);
1572
1573 // Append the name of the "node"
1574 // that is to be summed with the
1575 // RooRealSumPdf
1576 proto->import(*normFactors, RecycleConflictNodes());
1577 auto normFactorsInWS = dynamic_cast<RooProduct*>(proto->arg(normFactors->GetName()));
1578 assert(normFactorsInWS);
1579
1580 sampleScaleFactors.push_back(normFactorsInWS);
1581 } // END: Loop over EstimateSummaries
1582
1583 // If a non-zero number of samples call for
1584 // Stat Uncertainties, create the statFactor functions
1585 if(!statHistPairs.empty()) {
1586
1587 // Create the histogram of (binwise)
1588 // stat uncertainties:
1589 unique_ptr<TH1> fracStatError( MakeScaledUncertaintyHist( channel_name + "_StatUncert" + "_RelErr", statHistPairs) );
1590 if( fracStatError == NULL ) {
1591 cxcoutE(HistFactory) << "Error: Failed to make ScaledUncertaintyHist for: "
1592 << channel_name + "_StatUncert" + "_RelErr" << std::endl;
1593 throw hf_exc();
1594 }
1595
1596 // Using this TH1* of fractinal stat errors,
1597 // create a set of constraint terms:
1598 ParamHistFunc* chanStatUncertFunc = (ParamHistFunc*) proto->function( statFuncName.c_str() );
1599 cxcoutI(HistFactory) << "About to create Constraint Terms from: "
1600 << chanStatUncertFunc->GetName()
1601 << " params: " << chanStatUncertFunc->paramList()
1602 << std::endl;
1603
1604 // Get the constraint type and the
1605 // rel error threshold from the (last)
1606 // EstimateSummary looped over (but all
1607 // should be the same)
1608
1609 // Get the type of StatError constraint from the channel
1610 Constraint::Type statConstraintType = channel.GetStatErrorConfig().GetConstraintType();
1611 if( statConstraintType == Constraint::Gaussian) {
1612 cxcoutI(HistFactory) << "Using Gaussian StatErrors in channel: " << channel.GetName() << std::endl;
1613 }
1614 if( statConstraintType == Constraint::Poisson ) {
1615 cxcoutI(HistFactory) << "Using Poisson StatErrors in channel: " << channel.GetName() << std::endl;
1616 }
1617
1618 double statRelErrorThreshold = channel.GetStatErrorConfig().GetRelErrorThreshold();
1619 RooArgList statConstraints = createStatConstraintTerms(proto, constraintTermNames,
1620 *chanStatUncertFunc, fracStatError.get(),
1621 statConstraintType,
1622 statRelErrorThreshold);
1623
1624 } // END: Loop over stat Hist Pairs
1625
1626
1627 ///////////////////////////////////
1628 // for ith bin calculate totN_i = lumi * sum_j expected_j * syst_j
1629 MakeTotalExpected(proto, channel_name+"_model",
1630 sampleScaleFactors, allSampleHistFuncs);
1631 likelihoodTermNames.push_back(channel_name+"_model");
1632
1633 //////////////////////////////////////
1634 // fix specified parameters
1635 for(unsigned int i=0; i<systToFix.size(); ++i){
1636 RooRealVar* temp = proto->var((systToFix.at(i)).c_str());
1637 if(temp) {
1638 // set the parameter constant
1639 temp->setConstant();
1640
1641 // remove the corresponding auxiliary observable from the global observables
1642 RooRealVar* auxMeas = NULL;
1643 if(systToFix.at(i)=="Lumi"){
1644 auxMeas = proto->var("nominalLumi");
1645 } else {
1646 auxMeas = proto->var(TString::Format("nom_%s",temp->GetName()));
1647 }
1648
1649 if(auxMeas){
1650 const_cast<RooArgSet*>(proto->set("globalObservables"))->remove(*auxMeas);
1651 } else{
1652 cxcoutE(HistFactory) << "could not corresponding auxiliary measurement "
1653 << TString::Format("nom_%s",temp->GetName()) << endl;
1654 }
1655 } else {
1656 cxcoutE(HistFactory) << "could not find variable " << systToFix.at(i)
1657 << " could not set it to constant" << endl;
1658 }
1659 }
1660
1661 //////////////////////////////////////
1662 // final proto model
1663 for(unsigned int i=0; i<constraintTermNames.size(); ++i){
1664 RooAbsArg* proto_arg = (proto->arg(constraintTermNames[i].c_str()));
1665 if( proto_arg==NULL ) {
1666 cxcoutF(HistFactory) << "Error: Cannot find arg set: " << constraintTermNames.at(i)
1667 << " in workspace: " << proto->GetName() << std::endl;
1668 throw hf_exc();
1669 }
1670 constraintTerms.add( *proto_arg );
1671 // constraintTerms.add(* proto_arg(proto->arg(constraintTermNames[i].c_str())) );
1672 }
1673 for(unsigned int i=0; i<likelihoodTermNames.size(); ++i){
1674 RooAbsArg* proto_arg = (proto->arg(likelihoodTermNames[i].c_str()));
1675 if( proto_arg==NULL ) {
1676 cxcoutF(HistFactory) << "Error: Cannot find arg set: " << likelihoodTermNames.at(i)
1677 << " in workspace: " << proto->GetName() << std::endl;
1678 throw hf_exc();
1679 }
1680 likelihoodTerms.add( *proto_arg );
1681 }
1682 proto->defineSet("constraintTerms",constraintTerms);
1683 proto->defineSet("likelihoodTerms",likelihoodTerms);
1684
1685 // list of observables
1686 RooArgList observables;
1687 std::string observablesStr;
1688
1689 std::vector<std::string>::iterator itr = fObsNameVec.begin();
1690 for(; itr!=fObsNameVec.end(); ++itr ) {
1691 observables.add( *proto->var(itr->c_str()) );
1692 if (!observablesStr.empty()) { observablesStr += ","; }
1693 observablesStr += *itr;
1694 }
1695
1696 // We create two sets, one for backwards compatability
1697 // The other to make a consistent naming convention
1698 // between individual channels and the combined workspace
1699 proto->defineSet("observables", TString::Format("%s",observablesStr.c_str()));
1700 proto->defineSet("observablesSet", TString::Format("%s",observablesStr.c_str()));
1701
1702 // Create the ParamHistFunc
1703 // after observables have been made
1704 cxcoutP(HistFactory) << "\n-----------------------------------------\n"
1705 << "\timport model into workspace"
1706 << "\n-----------------------------------------\n" << endl;
1707
1708 auto model = make_unique<RooProdPdf>(
1709 ("model_"+channel_name).c_str(), // MB : have changed this into conditional pdf. Much faster for toys!
1710 "product of Poissons accross bins for a single channel",
1711 constraintTerms, Conditional(likelihoodTerms,observables));
1712 proto->import(*model,RecycleConflictNodes());
1713
1714 proto_config->SetPdf(*model);
1715 proto_config->SetObservables(observables);
1716 proto_config->SetGlobalObservables(*proto->set("globalObservables"));
1717 // proto->writeToFile(("results/model_"+channel+".root").c_str());
1718 // fill out nuisance parameters in model config
1719 // proto_config->GuessObsAndNuisance(*proto->data("asimovData"));
1720 proto->import(*proto_config,proto_config->GetName());
1721 proto->importClassCode();
1722
1723 ///////////////////////////
1724 // make data sets
1725 // THis works and is natural, but the memory size of the simultaneous dataset grows exponentially with channels
1726 const char* weightName="weightVar";
1727 proto->factory(TString::Format("%s[0,-1e10,1e10]",weightName));
1728 proto->defineSet("obsAndWeight",TString::Format("%s,%s",weightName,observablesStr.c_str()));
1729
1730 // New Asimov Generation: Use the code in the Asymptotic calculator
1731 // Need to get the ModelConfig...
1732 int asymcalcPrintLevel = 0;
1733 if (RooMsgService::instance().isActive(static_cast<TObject*>(nullptr), RooFit::HistFactory, RooFit::INFO)) asymcalcPrintLevel = 1;
1734 if (RooMsgService::instance().isActive(static_cast<TObject*>(nullptr), RooFit::HistFactory, RooFit::DEBUG)) asymcalcPrintLevel = 2;
1735 AsymptoticCalculator::SetPrintLevel(asymcalcPrintLevel);
1736 unique_ptr<RooAbsData> asimov_dataset(AsymptoticCalculator::GenerateAsimovData(*model, observables));
1737 proto->import(dynamic_cast<RooDataSet&>(*asimov_dataset), Rename("asimovData"));
1738
1739 // GHL: Determine to use data if the hist isn't 'NULL'
1740 if(TH1 const* mnominal = channel.GetData().GetHisto()) {
1741 // This works and is natural, but the memory size of the simultaneous
1742 // dataset grows exponentially with channels.
1743 RooDataSet dataset{"obsData","",*proto->set("obsAndWeight"),weightName};
1744 ConfigureHistFactoryDataset( dataset, *mnominal, *proto, fObsNameVec );
1745 proto->import(dataset);
1746 } // End: Has non-null 'data' entry
1747
1748
1749 for(auto const& data : channel.GetAdditionalData()) {
1750 if(data.GetName().empty()) {
1751 cxcoutF(HistFactory) << "Error: Additional Data histogram for channel: " << channel.GetName()
1752 << " has no name! The name always needs to be set for additional datasets, "
1753 << "either via the \"Name\" tag in the XML or via RooStats::HistFactory::Data::SetName()." << std::endl;
1754 throw hf_exc();
1755 }
1756 std::string const& dataName = data.GetName();
1757 TH1 const* mnominal = data.GetHisto();
1758 if( !mnominal ) {
1759 cxcoutF(HistFactory) << "Error: Additional Data histogram for channel: " << channel.GetName()
1760 << " with name: " << dataName << " is NULL" << std::endl;
1761 throw hf_exc();
1762 }
1763
1764 // THis works and is natural, but the memory size of the simultaneous dataset grows exponentially with channels
1765 RooDataSet dataset{dataName.c_str(), "", *proto->set("obsAndWeight"), weightName};
1766 ConfigureHistFactoryDataset( dataset, *mnominal, *proto, fObsNameVec );
1767 proto->import(dataset);
1768
1769 }
1770
1771 if (RooMsgService::instance().isActive(static_cast<TObject*>(nullptr), RooFit::HistFactory, RooFit::INFO))
1772 proto->Print();
1773
1774 return proto;
1775 }
1776
1777
1779 TH1 const& mnominal,
1781 std::vector<std::string> const& obsNameVec) {
1782
1783 // Take a RooDataSet and fill it with the entries
1784 // from a TH1*, using the observable names to
1785 // determine the columns
1786
1787 if (obsNameVec.empty() ) {
1788 Error("ConfigureHistFactoryDataset","Invalid input - return");
1789 return;
1790 }
1791
1792 TAxis const* ax = mnominal.GetXaxis();
1793 TAxis const* ay = mnominal.GetYaxis();
1794 TAxis const* az = mnominal.GetZaxis();
1795
1796 for (int i=1; i<=ax->GetNbins(); ++i) { // 1 or more dimension
1797
1798 double xval = ax->GetBinCenter(i);
1799 proto.var( obsNameVec[0].c_str() )->setVal( xval );
1800
1801 if(obsNameVec.size()==1) {
1802 double fval = mnominal.GetBinContent(i);
1803 obsDataUnbinned.add( *proto.set("obsAndWeight"), fval );
1804 } else { // 2 or more dimensions
1805
1806 for(int j=1; j<=ay->GetNbins(); ++j) {
1807 double yval = ay->GetBinCenter(j);
1808 proto.var( obsNameVec[1].c_str() )->setVal( yval );
1809
1810 if(obsNameVec.size()==2) {
1811 double fval = mnominal.GetBinContent(i,j);
1812 obsDataUnbinned.add( *proto.set("obsAndWeight"), fval );
1813 } else { // 3 dimensions
1814
1815 for(int k=1; k<=az->GetNbins(); ++k) {
1816 double zval = az->GetBinCenter(k);
1817 proto.var( obsNameVec[2].c_str() )->setVal( zval );
1818 double fval = mnominal.GetBinContent(i,j,k);
1819 obsDataUnbinned.add( *proto.set("obsAndWeight"), fval );
1820 }
1821 }
1822 }
1823 }
1824 }
1825 }
1826
1828 {
1829 fObsNameVec.clear();
1830
1831 // determine histogram dimensionality
1832 unsigned int histndim(1);
1833 std::string classname = hist->ClassName();
1834 if (classname.find("TH1")==0) { histndim=1; }
1835 else if (classname.find("TH2")==0) { histndim=2; }
1836 else if (classname.find("TH3")==0) { histndim=3; }
1837
1838 for ( unsigned int idx=0; idx<histndim; ++idx ) {
1839 if (idx==0) { fObsNameVec.push_back("x"); }
1840 if (idx==1) { fObsNameVec.push_back("y"); }
1841 if (idx==2) { fObsNameVec.push_back("z"); }
1842 }
1843 }
1844
1845
1846 RooWorkspace* HistoToWorkspaceFactoryFast::MakeCombinedModel(vector<string> ch_names, vector<std::unique_ptr<RooWorkspace>>& chs)
1847 {
1849
1850 // check first the inputs (see JIRA-6890)
1851 if (ch_names.empty() || chs.empty() ) {
1852 Error("MakeCombinedModel","Input vectors are empty - return a nullptr");
1853 return 0;
1854 }
1855 if (chs.size() < ch_names.size() ) {
1856 Error("MakeCombinedModel","Input vector of workspace has an invalid size - return a nullptr");
1857 return 0;
1858 }
1859
1860 //
1861 /// These things were used for debugging. Maybe useful in the future
1862 //
1863
1864 map<string, RooAbsPdf*> pdfMap;
1865 vector<RooAbsPdf*> models;
1866
1867 RooArgList obsList;
1868 for(unsigned int i = 0; i< ch_names.size(); ++i){
1869 ModelConfig * config = (ModelConfig *) chs[i]->obj("ModelConfig");
1870 obsList.add(*config->GetObservables());
1871 }
1872 cxcoutI(HistFactory) <<"full list of observables:\n" << obsList << std::endl;
1873
1874 RooArgSet globalObs;
1875 stringstream channelString;
1876 channelString << "channelCat[";
1877 for(unsigned int i = 0; i< ch_names.size(); ++i){
1878 string channel_name=ch_names[i];
1879 if (i == 0 && isdigit(channel_name[0])) {
1880 throw std::invalid_argument("The first channel name for HistFactory cannot start with a digit. Got " + channel_name);
1881 }
1882 if (channel_name.find(',') != std::string::npos) {
1883 throw std::invalid_argument("Channel names for HistFactory cannot contain ','. Got " + channel_name);
1884 }
1885
1886 if (i == 0) channelString << channel_name ;
1887 else channelString << ',' << channel_name ;
1888 RooWorkspace * ch=chs[i].get();
1889
1890 RooAbsPdf* model = ch->pdf(("model_"+channel_name).c_str());
1891 if(!model) cout <<"failed to find model for channel"<<endl;
1892 // cout << "int = " << model->createIntegral(*obsN)->getVal() << endl;;
1893 models.push_back(model);
1894 globalObs.add(*ch->set("globalObservables"), /*silent=*/true); // silent because observables might exist in other channel.
1895
1896 // constrainedParams->add( * ch->set("constrainedParams") );
1897 pdfMap[channel_name]=model;
1898 }
1899 channelString << "]";
1900
1901 cxcoutP(HistFactory) << "\n-----------------------------------------\n"
1902 << "\tEntering combination"
1903 << "\n-----------------------------------------\n" << endl;
1904 RooWorkspace* combined = new RooWorkspace("combined");
1905 // RooWorkspace* combined = chs[0];
1906
1907
1908 RooCategory* channelCat = dynamic_cast<RooCategory*>( combined->factory(channelString.str().c_str()) );
1909 if (!channelCat) throw std::runtime_error("Unable to construct a category from string " + channelString.str());
1910
1911 auto simPdf= std::make_unique<RooSimultaneous>("simPdf","",pdfMap, *channelCat);
1912 auto combined_config = std::make_unique<ModelConfig>("ModelConfig", combined);
1913 combined_config->SetWorkspace(*combined);
1914 // combined_config->SetNuisanceParameters(*constrainedParams);
1915
1916 combined->import(globalObs);
1917 combined->defineSet("globalObservables",globalObs);
1918 combined_config->SetGlobalObservables(*combined->set("globalObservables"));
1919
1920
1921 ////////////////////////////////////////////
1922 // Make toy simultaneous dataset
1923 cxcoutP(HistFactory) << "\n-----------------------------------------\n"
1924 << "\tcreate toy data for " << channelString.str()
1925 << "\n-----------------------------------------\n" << endl;
1926
1927
1928 // now with weighted datasets
1929 // First Asimov
1930 //RooDataSet * simData=NULL;
1931 combined->factory("weightVar[0,-1e10,1e10]");
1932 obsList.add(*combined->var("weightVar"));
1933
1934 // Create Asimov data for the combined dataset
1935 std::unique_ptr<RooDataSet> asimov_combined{static_cast<RooDataSet*>(AsymptoticCalculator::GenerateAsimovData(*simPdf,
1936 obsList))};
1937 if( asimov_combined ) {
1938 combined->import( *asimov_combined, Rename("asimovData"));
1939 }
1940 else {
1941 std::cout << "Error: Failed to create combined asimov dataset" << std::endl;
1942 throw hf_exc();
1943 }
1944
1945 // Now merge the observable datasets across the channels
1946 for(RooAbsData * data : chs[0]->allData()) {
1947 // Only include RooDataSets where a data with the same name doesn't only
1948 // exist in the merged workspace (like the "asimovData").
1949 if(dynamic_cast<RooDataSet*>(data) && !combined->data(data->GetName())) {
1950 MergeDataSets(combined, chs, ch_names, data->GetName(), obsList, channelCat);
1951 }
1952 }
1953
1954
1955 obsList.add(*channelCat);
1956 combined->defineSet("observables",obsList);
1957 combined_config->SetObservables(*combined->set("observables"));
1958
1959 if (RooMsgService::instance().isActive(static_cast<TObject*>(nullptr), RooFit::HistFactory, RooFit::INFO))
1960 combined->Print();
1961
1962 cxcoutP(HistFactory) << "\n-----------------------------------------\n"
1963 << "\tImporting combined model"
1964 << "\n-----------------------------------------\n" << endl;
1965 combined->import(*simPdf,RecycleConflictNodes());
1966 //combined->import(*simPdf, RenameVariable("SigXsecOverSM","SigXsecOverSM_comb"));
1967 // cout << "check pointer " << simPdf << endl;
1968 // cout << "check val " << simPdf->getVal() << endl;
1969
1970 std::map< std::string, double>::iterator param_itr = fParamValues.begin();
1971 for( ; param_itr != fParamValues.end(); ++param_itr ){
1972 // make sure they are fixed
1973 std::string paramName = param_itr->first;
1974 double paramVal = param_itr->second;
1975
1976 if(RooRealVar* temp = combined->var( paramName.c_str() )) {
1977 temp->setVal( paramVal );
1978 cxcoutI(HistFactory) <<"setting " << paramName << " to the value: " << paramVal << endl;
1979 } else
1980 cxcoutE(HistFactory) << "could not find variable " << paramName << " could not set its value" << endl;
1981 }
1982
1983
1984 for(unsigned int i=0; i<fSystToFix.size(); ++i){
1985 // make sure they are fixed
1986 if(RooRealVar* temp = combined->var(fSystToFix[i].c_str())) {
1987 temp->setConstant();
1988 cxcoutI(HistFactory) <<"setting " << fSystToFix.at(i) << " constant" << endl;
1989 } else
1990 cxcoutE(HistFactory) << "could not find variable " << fSystToFix.at(i) << " could not set it to constant" << endl;
1991 }
1992
1993 ///
1994 /// writing out the model in graphViz
1995 ///
1996 // RooAbsPdf* customized=combined->pdf("simPdf");
1997 //combined_config->SetPdf(*customized);
1998 combined_config->SetPdf(*simPdf);
1999 // combined_config->GuessObsAndNuisance(*simData);
2000 // customized->graphVizTree(("results/"+fResultsPrefixStr.str()+"_simul.dot").c_str());
2001 combined->import(*combined_config,combined_config->GetName());
2002 combined->importClassCode();
2003 // combined->writeToFile("results/model_combined.root");
2004
2005 return combined;
2006 }
2007
2008
2010 std::vector<std::unique_ptr<RooWorkspace>>& wspace_vec,
2011 std::vector<std::string> const& channel_names,
2012 std::string const& dataSetName,
2013 RooArgList obsList,
2014 RooCategory* channelCat) {
2015
2016 // Create the total dataset
2017 std::unique_ptr<RooDataSet> simData;
2018
2019 // Loop through channels, get their individual datasets,
2020 // and add them to the combined dataset
2021 for(unsigned int i = 0; i< channel_names.size(); ++i){
2022
2023 // Grab the dataset for the existing channel
2024 cxcoutPHF << "Merging data for channel " << channel_names[i].c_str() << std::endl;
2025 RooDataSet* obsDataInChannel = (RooDataSet*) wspace_vec[i]->data(dataSetName.c_str());
2026 if( !obsDataInChannel ) {
2027 std::cout << "Error: Can't find DataSet: " << dataSetName
2028 << " in channel: " << channel_names.at(i)
2029 << std::endl;
2030 throw hf_exc();
2031 }
2032
2033 // Create the new Dataset
2034 auto tempData = std::make_unique<RooDataSet>(channel_names[i].c_str(),"",
2035 obsList, Index(*channelCat),
2036 WeightVar("weightVar"),
2037 Import(channel_names[i].c_str(),*obsDataInChannel));
2038 if(simData) {
2039 simData->append(*tempData);
2040 }
2041 else {
2042 simData = std::move(tempData);
2043 }
2044 } // End Loop Over Channels
2045
2046 // Check that we successfully created the dataset
2047 // and import it into the workspace
2048 if(simData) {
2049 combined->import(*simData, Rename(dataSetName.c_str()));
2050 return static_cast<RooDataSet*>(combined->data(dataSetName.c_str()));
2051 }
2052 else {
2053 std::cout << "Error: Unable to merge observable datasets" << std::endl;
2054 throw hf_exc();
2055 return nullptr;
2056 }
2057 }
2058
2059
2061
2062 // Take a nominal TH1* and create
2063 // a TH1 representing the binwise
2064 // errors (taken from the nominal TH1)
2065
2066 TH1* ErrorHist = (TH1*) Nominal->Clone( Name.c_str() );
2067 ErrorHist->Reset();
2068
2069 Int_t numBins = Nominal->GetNbinsX()*Nominal->GetNbinsY()*Nominal->GetNbinsZ();
2070 Int_t binNumber = 0;
2071
2072 // Loop over bins
2073 for( Int_t i_bin = 0; i_bin < numBins; ++i_bin) {
2074
2075 binNumber++;
2076 // Ignore underflow / overflow
2077 while( Nominal->IsBinUnderflow(binNumber) || Nominal->IsBinOverflow(binNumber) ){
2078 binNumber++;
2079 }
2080
2081 Double_t histError = Nominal->GetBinError( binNumber );
2082
2083 // Check that histError != NAN
2084 if( histError != histError ) {
2085 std::cout << "Warning: In histogram " << Nominal->GetName()
2086 << " bin error for bin " << i_bin
2087 << " is NAN. Not using Error!!!"
2088 << std::endl;
2089 throw hf_exc();
2090 //histError = sqrt( histContent );
2091 //histError = 0;
2092 }
2093
2094 // Check that histError ! < 0
2095 if( histError < 0 ) {
2096 std::cout << "Warning: In histogram " << Nominal->GetName()
2097 << " bin error for bin " << binNumber
2098 << " is < 0. Setting Error to 0"
2099 << std::endl;
2100 //histError = sqrt( histContent );
2101 histError = 0;
2102 }
2103
2104 ErrorHist->SetBinContent( binNumber, histError );
2105
2106 }
2107
2108 return ErrorHist;
2109
2110 }
2111
2112 // Take a list of < nominal, absolError > TH1* pairs
2113 // and construct a single histogram representing the
2114 // total fractional error as:
2115
2116 // UncertInQuad(bin i) = Sum: absolUncert*absolUncert
2117 // Total(bin i) = Sum: Value
2118 //
2119 // TotalFracError(bin i) = Sqrt( UncertInQuad(i) ) / TotalBin(i)
2120 std::unique_ptr<TH1> HistoToWorkspaceFactoryFast::MakeScaledUncertaintyHist( const std::string& Name, std::vector< std::pair<const TH1*, std::unique_ptr<TH1>> > const& HistVec ) const {
2121
2122
2123 unsigned int numHists = HistVec.size();
2124
2125 if( numHists == 0 ) {
2126 cxcoutE(HistFactory) << "Warning: Empty Hist Vector, cannot create total uncertainty" << std::endl;
2127 return NULL;
2128 }
2129
2130 const TH1* HistTemplate = HistVec.at(0).first;
2131 Int_t numBins = HistTemplate->GetNbinsX()*HistTemplate->GetNbinsY()*HistTemplate->GetNbinsZ();
2132
2133 // Check that all histograms
2134 // have the same bins
2135 for( unsigned int i = 0; i < HistVec.size(); ++i ) {
2136
2137 const TH1* nominal = HistVec.at(i).first;
2138 const TH1* error = HistVec.at(i).second.get();
2139
2140 if( nominal->GetNbinsX()*nominal->GetNbinsY()*nominal->GetNbinsZ() != numBins ) {
2141 cxcoutE(HistFactory) << "Error: Provided hists have unequal bins" << std::endl;
2142 return NULL;
2143 }
2144 if( error->GetNbinsX()*error->GetNbinsY()*error->GetNbinsZ() != numBins ) {
2145 cxcoutE(HistFactory) << "Error: Provided hists have unequal bins" << std::endl;
2146 return NULL;
2147 }
2148 }
2149
2150 std::vector<double> TotalBinContent( numBins, 0.0);
2151 std::vector<double> HistErrorsSqr( numBins, 0.0);
2152
2153 Int_t binNumber = 0;
2154
2155 // Loop over bins
2156 for( Int_t i_bins = 0; i_bins < numBins; ++i_bins) {
2157
2158 binNumber++;
2159 while( HistTemplate->IsBinUnderflow(binNumber) || HistTemplate->IsBinOverflow(binNumber) ){
2160 binNumber++;
2161 }
2162
2163 for( unsigned int i_hist = 0; i_hist < numHists; ++i_hist ) {
2164
2165 const TH1* nominal = HistVec.at(i_hist).first;
2166 const TH1* error = HistVec.at(i_hist).second.get();
2167
2168 //Int_t binNumber = i_bins + 1;
2169
2170 Double_t histValue = nominal->GetBinContent( binNumber );
2171 Double_t histError = error->GetBinContent( binNumber );
2172 /*
2173 std::cout << " Getting Bin content for Stat Uncertainty"
2174 << " Nom name: " << nominal->GetName()
2175 << " Err name: " << error->GetName()
2176 << " HistNumber: " << i_hist << " bin: " << binNumber
2177 << " Value: " << histValue << " Error: " << histError
2178 << std::endl;
2179 */
2180
2181 if( histError != histError ) {
2182 cxcoutE(HistFactory) << "In histogram " << error->GetName()
2183 << " bin error for bin " << binNumber
2184 << " is NAN. Not using error!!";
2185 throw hf_exc();
2186 }
2187
2188 TotalBinContent.at(i_bins) += histValue;
2189 HistErrorsSqr.at(i_bins) += histError*histError; // Add in quadrature
2190
2191 }
2192 }
2193
2194 binNumber = 0;
2195
2196 // Creat the output histogram
2197 TH1* ErrorHist = (TH1*) HistTemplate->Clone( Name.c_str() );
2198 ErrorHist->Reset();
2199
2200 // Fill the output histogram
2201 for( Int_t i = 0; i < numBins; ++i) {
2202
2203 // Int_t binNumber = i + 1;
2204 binNumber++;
2205 while( ErrorHist->IsBinUnderflow(binNumber) || ErrorHist->IsBinOverflow(binNumber) ){
2206 binNumber++;
2207 }
2208
2209 Double_t ErrorsSqr = HistErrorsSqr.at(i);
2210 Double_t TotalVal = TotalBinContent.at(i);
2211
2212 if( TotalVal <= 0 ) {
2213 cxcoutW(HistFactory) << "Warning: Sum of histograms for bin: " << binNumber
2214 << " is <= 0. Setting error to 0"
2215 << std::endl;
2216
2217 ErrorHist->SetBinContent( binNumber, 0.0 );
2218 continue;
2219 }
2220
2221 Double_t RelativeError = sqrt(ErrorsSqr) / TotalVal;
2222
2223 // If we otherwise get a NAN
2224 // it's an error
2225 if( RelativeError != RelativeError ) {
2226 cxcoutE(HistFactory) << "Error: bin " << i << " error is NAN\n"
2227 << " HistErrorsSqr: " << ErrorsSqr
2228 << " TotalVal: " << TotalVal;
2229 throw hf_exc();
2230 }
2231
2232 // 0th entry in vector is
2233 // the 1st bin in TH1
2234 // (we ignore underflow)
2235
2236 // Error and bin content are interchanged because for some reason, the other functions
2237 // use the bin content to convey the error ...
2238 ErrorHist->SetBinError(binNumber, TotalVal);
2239 ErrorHist->SetBinContent(binNumber, RelativeError);
2240
2241 cxcoutI(HistFactory) << "Making Total Uncertainty for bin " << binNumber
2242 << " Error = " << sqrt(ErrorsSqr)
2243 << " CentralVal = " << TotalVal
2244 << " RelativeError = " << RelativeError << "\n";
2245
2246 }
2247
2248 return std::unique_ptr<TH1>(ErrorHist);
2249}
2250
2251
2252
2254 createStatConstraintTerms( RooWorkspace* proto, vector<string>& constraintTermNames,
2255 ParamHistFunc& paramHist, const TH1* uncertHist,
2256 Constraint::Type type, Double_t minSigma ) {
2257
2258
2259 // Take a RooArgList of RooAbsReal's and
2260 // create N constraint terms (one for
2261 // each gamma) whose relative uncertainty
2262 // is the value of the ith RooAbsReal
2263 //
2264 // The integer "type" controls the type
2265 // of constraint term:
2266 //
2267 // type == 0 : NONE
2268 // type == 1 : Gaussian
2269 // type == 2 : Poisson
2270 // type == 3 : LogNormal
2271
2272 RooArgList ConstraintTerms;
2273
2274 RooArgList paramSet = paramHist.paramList();
2275
2276 // Must get the full size of the TH1
2277 // (No direct method to do this...)
2278 Int_t numBins = uncertHist->GetNbinsX()*uncertHist->GetNbinsY()*uncertHist->GetNbinsZ();
2279 Int_t numParams = paramSet.getSize();
2280 // Int_t numBins = uncertHist->GetNbinsX()*uncertHist->GetNbinsY()*uncertHist->GetNbinsZ();
2281
2282 // Check that there are N elements
2283 // in the RooArgList
2284 if( numBins != numParams ) {
2285 std::cout << "Error: In createStatConstraintTerms, encountered bad number of bins" << std::endl;
2286 std::cout << "Given histogram with " << numBins << " bins,"
2287 << " but require exactly " << numParams << std::endl;
2288 throw hf_exc();
2289 }
2290
2291 Int_t TH1BinNumber = 0;
2292 for( Int_t i = 0; i < paramSet.getSize(); ++i) {
2293
2294 TH1BinNumber++;
2295
2296 while( uncertHist->IsBinUnderflow(TH1BinNumber) || uncertHist->IsBinOverflow(TH1BinNumber) ){
2297 TH1BinNumber++;
2298 }
2299
2300 RooRealVar& gamma = (RooRealVar&) (paramSet[i]);
2301
2302 cxcoutI(HistFactory) << "Creating constraint for: " << gamma.GetName()
2303 << ". Type of constraint: " << type << std::endl;
2304
2305 // Get the sigma from the hist
2306 // (the relative uncertainty)
2307 const double sigmaRel = uncertHist->GetBinContent(TH1BinNumber);
2308
2309 // If the sigma is <= 0,
2310 // do cont create the term
2311 if( sigmaRel <= 0 ){
2312 cxcoutI(HistFactory) << "Not creating constraint term for "
2313 << gamma.GetName()
2314 << " because sigma = " << sigmaRel
2315 << " (sigma<=0)"
2316 << " (TH1 bin number = " << TH1BinNumber << ")"
2317 << std::endl;
2318 gamma.setConstant(kTRUE);
2319 continue;
2320 }
2321
2322 // set reasonable ranges for gamma parameters
2323 gamma.setMax( 1 + 5*sigmaRel );
2324 gamma.setMin( 0. );
2325
2326 // Make Constraint Term
2327 std::string constrName = string(gamma.GetName()) + "_constraint";
2328 std::string nomName = string("nom_") + gamma.GetName();
2329 std::string sigmaName = string(gamma.GetName()) + "_sigma";
2330 std::string poisMeanName = string(gamma.GetName()) + "_poisMean";
2331
2332 if( type == Constraint::Gaussian ) {
2333
2334 // Type 1 : RooGaussian
2335
2336 // Make sigma
2337
2338 RooConstVar constrSigma( sigmaName.c_str(), sigmaName.c_str(), sigmaRel );
2339
2340 // Make "observed" value
2341 RooRealVar constrNom(nomName.c_str(), nomName.c_str(), 1.0,0,10);
2342 constrNom.setConstant( true );
2343
2344 // Make the constraint:
2345 RooGaussian gauss( constrName.c_str(), constrName.c_str(),
2346 constrNom, gamma, constrSigma );
2347
2348 proto->import( gauss, RecycleConflictNodes() );
2349
2350 // Give reasonable starting point for pre-fit errors by setting it to the absolute sigma
2351 // Mostly useful for pre-fit plotting.
2352 gamma.setError(sigmaRel);
2353 } else if( type == Constraint::Poisson ) {
2354
2355 Double_t tau = 1/sigmaRel/sigmaRel; // this is correct Poisson equivalent to a Gaussian with mean 1 and stdev sigma
2356
2357 // Make nominal "observed" value
2358 RooRealVar constrNom(nomName.c_str(), nomName.c_str(), tau);
2359 constrNom.setMin(0);
2360 constrNom.setConstant( true );
2361
2362 // Make the scaling term
2363 std::string scalingName = string(gamma.GetName()) + "_tau";
2364 RooConstVar poissonScaling( scalingName.c_str(), scalingName.c_str(), tau);
2365
2366 // Make mean for scaled Poisson
2367 RooProduct constrMean( poisMeanName.c_str(), poisMeanName.c_str(), RooArgSet(gamma, poissonScaling) );
2368 //proto->import( constrSigma, RecycleConflictNodes() );
2369 //proto->import( constrSigma );
2370
2371 // Type 2 : RooPoisson
2372 RooPoisson pois(constrName.c_str(), constrName.c_str(), constrNom, constrMean);
2373 pois.setNoRounding(true);
2374 proto->import( pois, RecycleConflictNodes() );
2375
2376 if (std::string(gamma.GetName()).find("gamma_stat") != std::string::npos) {
2377 // Give reasonable starting point for pre-fit errors.
2378 // Mostly useful for pre-fit plotting.
2379 gamma.setError(sigmaRel);
2380 }
2381
2382 } else {
2383
2384 std::cout << "Error: Did not recognize Stat Error constraint term type: "
2385 << type << " for : " << paramHist.GetName() << std::endl;
2386 throw hf_exc();
2387 }
2388
2389 // If the sigma value is less
2390 // than a supplied threshold,
2391 // set the variable to constant
2392 if( sigmaRel < minSigma ) {
2393 cxcoutW(HistFactory) << "Warning: Bin " << i << " = " << sigmaRel
2394 << " and is < " << minSigma
2395 << ". Setting: " << gamma.GetName() << " to constant"
2396 << std::endl;
2397 gamma.setConstant(kTRUE);
2398 }
2399
2400 constraintTermNames.push_back( constrName );
2401 ConstraintTerms.add( *proto->pdf(constrName.c_str()) );
2402
2403 // Add the "observed" value to the
2404 // list of global observables:
2405 RooArgSet* globalSet = const_cast<RooArgSet*>(proto->set("globalObservables"));
2406
2407 RooRealVar* nomVarInWorkspace = proto->var(nomName.c_str());
2408 if( ! globalSet->contains(*nomVarInWorkspace) ) {
2409 globalSet->add( *nomVarInWorkspace );
2410 }
2411
2412 } // end loop over parameters
2413
2414 return ConstraintTerms;
2415
2416}
2417
2418} // namespace RooStats
2419} // namespace HistFactory
2420
#define cxcoutPHF
#define cxcoutFHF
#define cxcoutIHF
#define cxcoutWHF
#define alpha_Low
#define alpha_High
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define cxcoutI(a)
#define cxcoutW(a)
#define cxcoutF(a)
#define cxcoutE(a)
#define cxcoutP(a)
const Bool_t kTRUE
Definition RtypesCore.h:100
#define ClassImp(name)
Definition Rtypes.h:364
#define R__ASSERT(e)
Definition TError.h:118
char name[80]
Definition TGX11.cxx:110
int type
Definition TGX11.cxx:121
float xmin
float xmax
char * Form(const char *fmt,...)
const char * proto
Definition civetweb.c:16613
A class which maps the current values of a RooRealVar (or a set of RooRealVars) to one of a number of...
const RooArgSet * get(Int_t masterIdx) const
void setConstant(bool constant)
const RooArgList & paramList() const
static RooArgList createParamSet(RooWorkspace &w, const std::string &, const RooArgList &Vars)
Create the list of RooRealVar parameters which represent the height of the histogram bins.
void setShape(TH1 *shape)
The PiecewiseInterpolation is a class that can morph distributions into each other,...
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:69
virtual void Print(Option_t *options=0) const
Print the object to the defaultPrintStream().
Definition RooAbsArg.h:337
void setAttribute(const Text_t *name, Bool_t value=kTRUE)
Set (default) or clear a named boolean attribute of this object.
Bool_t isConstant() const
Check if the "Constant" attribute is set.
Definition RooAbsArg.h:377
Int_t getSize() const
Bool_t contains(const RooAbsArg &var) const
Check if collection contains an argument with the same name as var.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
Storage_t::size_type size() const
const_iterator begin() const
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
TIterator * createIterator(Bool_t dir=kIterForward) const
TIterator-style iteration over contained elements.
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:82
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
void setConstant(Bool_t value=kTRUE)
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:94
RooNumIntConfig * specialIntegratorConfig() const
Returns the specialized integrator configuration for this RooAbsReal.
virtual void forceNumInt(Bool_t flag=kTRUE)
Definition RooAbsReal.h:176
static void setEvalErrorLoggingMode(ErrorLoggingMode m)
Set evaluation error logging mode.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:35
RooBinWidthFunction is a class that returns the bin width (or volume) given a RooHistFunc.
Class RooBinning is an implements RooAbsBinning in terms of an array of boundary values,...
Definition RooBinning.h:28
RooCategory is an object to represent discrete states.
Definition RooCategory.h:27
virtual Bool_t setLabel(const char *label, bool printError=true) override
Set value by specifying the name of the desired state.
RooConstVar represent a constant real-valued object.
Definition RooConstVar.h:26
The RooDataHist is a container class to hold N-dimensional binned data.
Definition RooDataHist.h:45
RooDataSet is a container class to hold unbinned data.
Definition RooDataSet.h:36
virtual void add(const RooArgSet &row, Double_t weight=1.0, Double_t weightError=0) override
Add a data point, with its coordinates specified in the 'data' argset, to the data set.
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
Double_t correlation(const RooAbsArg &par1, const RooAbsArg &par2) const
Return correlation between par1 and par2.
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
Switches the message service to a different level while the instance is alive.
Definition RooHelpers.h:42
RooHistFunc implements a real-valued function sampled from a multidimensional histogram.
Definition RooHistFunc.h:30
static RooMsgService & instance()
Return reference to singleton instance.
Bool_t isActive(const RooAbsArg *self, RooFit::MsgTopic facility, RooFit::MsgLevel level)
Check if logging is active for given object/topic/RooFitMsgLevel combination.
RooCategory & method2D()
RooCategory & methodND()
RooCategory & method1D()
Poisson pdf.
Definition RooPoisson.h:19
void setNoRounding(bool flag=kTRUE)
Switch off/on rounding of x to the nearest integer.
Definition RooPoisson.h:34
A RooProduct represents the product of a given set of RooAbsReal objects.
Definition RooProduct.h:29
The class RooRealSumPdf implements a PDF constructed from a sum of functions:
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
void setMin(const char *name, Double_t value)
Set minimum of name range to given value.
void setBins(Int_t nBins, const char *name=0)
Create a uniform binning under name 'name' for this variable.
virtual void setVal(Double_t value)
Set value of variable to 'value'.
static void SetPrintLevel(int level)
set print level (static function)
static RooAbsData * GenerateAsimovData(const RooAbsPdf &pdf, const RooArgSet &observables)
generate the asimov data for the observables (not the global ones) need to deal with the case of a si...
TODO Here, we are missing some documentation.
Definition Asimov.h:22
void ConfigureWorkspace(RooWorkspace *)
Definition Asimov.cxx:22
This class encapsulates all information for the statistical interpretation of one experiment.
Definition Channel.h:30
std::vector< RooStats::HistFactory::Data > & GetAdditionalData()
retrieve vector of additional data objects
Definition Channel.h:65
void Print(std::ostream &=std::cout)
Definition Channel.cxx:75
HistFactory::StatErrorConfig & GetStatErrorConfig()
get information about threshold for statistical uncertainties and constraint term
Definition Channel.h:72
RooStats::HistFactory::Data & GetData()
get data object
Definition Channel.h:59
std::vector< RooStats::HistFactory::Sample > & GetSamples()
get vector of samples for this channel
Definition Channel.h:77
std::string GetName() const
get name of channel
Definition Channel.h:43
This class provides helper functions for creating likelihood models from histograms.
std::unique_ptr< RooProduct > CreateNormFactor(RooWorkspace *proto, std::string &channel, std::string &sigmaEpsilon, Sample &sample, bool doRatio)
static void EditSyst(RooWorkspace *proto, const char *pdfNameChar, std::map< std::string, double > gammaSyst, std::map< std::string, double > uniformSyst, std::map< std::string, double > logNormSyst, std::map< std::string, double > noSyst)
RooWorkspace * MakeSingleChannelModel(Measurement &measurement, Channel &channel)
RooWorkspace * MakeSingleChannelWorkspace(Measurement &measurement, Channel &channel)
void SetObsToExpected(RooWorkspace *proto, std::string obsPrefix, std::string expPrefix, int lowBin, int highBin)
std::unique_ptr< TH1 > MakeScaledUncertaintyHist(const std::string &Name, std::vector< std::pair< const TH1 *, std::unique_ptr< TH1 > > > const &HistVec) const
void SetFunctionsToPreprocess(std::vector< std::string > lines)
RooHistFunc * MakeExpectedHistFunc(const TH1 *hist, RooWorkspace *proto, std::string prefix, const RooArgList &observables) const
Create the nominal hist function from hist, and register it in the workspace.
TH1 * MakeAbsolUncertaintyHist(const std::string &Name, const TH1 *Hist)
static void ConfigureWorkspaceForMeasurement(const std::string &ModelName, RooWorkspace *ws_single, Measurement &measurement)
RooArgList createStatConstraintTerms(RooWorkspace *proto, std::vector< std::string > &constraintTerms, ParamHistFunc &paramHist, const TH1 *uncertHist, Constraint::Type type, Double_t minSigma)
void MakeTotalExpected(RooWorkspace *proto, const std::string &totName, const std::vector< RooProduct * > &sampleScaleFactors, std::vector< std::vector< RooAbsArg * > > &sampleHistFuncs) const
void AddPoissonTerms(RooWorkspace *proto, std::string prefix, std::string obsPrefix, std::string expPrefix, int lowBin, int highBin, std::vector< std::string > &likelihoodTermNames)
void ConfigureHistFactoryDataset(RooDataSet &obsData, TH1 const &nominal, RooWorkspace &proto, std::vector< std::string > const &obsNameVec)
static void PrintCovarianceMatrix(RooFitResult *result, RooArgSet *params, std::string filename)
RooDataSet * MergeDataSets(RooWorkspace *combined, std::vector< std::unique_ptr< RooWorkspace > > &wspace_vec, std::vector< std::string > const &channel_names, std::string const &dataSetName, RooArgList obsList, RooCategory *channelCat)
void AddConstraintTerms(RooWorkspace *proto, Measurement &measurement, std::string prefix, std::string interpName, std::vector< OverallSys > &systList, std::vector< std::string > &likelihoodTermNames, std::vector< std::string > &totSystTermNames)
RooArgList createObservables(const TH1 *hist, RooWorkspace *proto) const
Create observables of type RooRealVar. Creates 1 to 3 observables, depending on the type of the histo...
RooWorkspace * MakeCombinedModel(std::vector< std::string >, std::vector< std::unique_ptr< RooWorkspace > > &)
The RooStats::HistFactory::Measurement class can be used to construct a model by combining multiple R...
Definition Measurement.h:30
std::map< std::string, double > & GetGammaSyst()
std::map< std::string, double > & GetLogNormSyst()
std::map< std::string, double > & GetNoSyst()
std::vector< std::string > & GetPOIList()
get vector of PoI names
Definition Measurement.h:50
std::map< std::string, double > & GetUniformSyst()
std::vector< std::string > & GetConstantParams()
get vector of all constant parameters
Definition Measurement.h:59
std::vector< RooStats::HistFactory::Channel > & GetChannels()
std::vector< RooStats::HistFactory::Asimov > & GetAsimovDatasets()
get vector of defined Asimov Datasets
Definition Measurement.h:79
std::vector< std::string > GetPreprocessFunctions() const
Returns a list of defined preprocess function expressions.
double GetLumi()
retrieve integrated luminosity
Definition Measurement.h:88
Configuration for an un- constrained overall systematic to scale sample normalisations.
Definition Systematics.h:77
Configuration for a constrained overall systematic to scale sample normalisations.
Definition Systematics.h:49
const std::string & GetName() const
Definition Systematics.h:56
std::string GetName() const
get name of sample
Definition Sample.h:83
std::vector< RooStats::HistFactory::NormFactor > & GetNormFactorList()
Definition Sample.h:110
*Un*constrained bin-by-bin variation of affected histogram.
const TH1 * GetInitialShape() const
Constrained bin-by-bin variation of affected histogram.
Constraint::Type GetConstraintType() const
const TH1 * GetErrorHist() const
Constraint::Type GetConstraintType() const
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition ModelConfig.h:30
void GuessObsAndNuisance(const RooAbsData &data, bool printModelConfig=true)
Makes sensible guesses of observables, parameters of interest and nuisance parameters if one or multi...
virtual void SetParametersOfInterest(const RooArgSet &set)
Specify parameters of interest.
const RooArgSet * GetObservables() const
get RooArgSet for observables (return NULL if not existing)
virtual void SetPdf(const RooAbsPdf &pdf)
Set the Pdf, add to the the workspace if not already there.
Definition ModelConfig.h:81
The RooWorkspace is a persistable container for RooFit projects.
RooAbsData * data(const char *name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found.
RooArgSet allVars() const
Return set with all variable objects.
Bool_t import(const RooAbsArg &arg, const RooCmdArg &arg1=RooCmdArg(), const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg(), const RooCmdArg &arg9=RooCmdArg())
Import a RooAbsArg object, e.g.
Bool_t importClassCode(const char *pat="*", Bool_t doReplace=kFALSE)
Inport code of all classes in the workspace that have a class name that matches pattern 'pat' and whi...
void Print(Option_t *opts=0) const
Print contents of the workspace.
Bool_t defineSet(const char *name, const RooArgSet &aset, Bool_t importMissing=kFALSE)
Define a named RooArgSet with given constituents.
Bool_t saveSnapshot(const char *name, const char *paramNames)
Save snapshot of values and attributes (including "Constant") of given parameters.
Bool_t loadSnapshot(const char *name)
Load the values and attributes of the parameters in the snapshot saved with the given name.
RooRealVar * var(const char *name) const
Retrieve real-valued variable (RooRealVar) with given name. A null pointer is returned if not found.
RooFactoryWSTool & factory()
Return instance to factory tool.
const RooArgSet * set(const char *name)
Return pointer to previously defined named set with given nmame If no such set is found a null pointe...
TObject * obj(const char *name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
RooAbsPdf * pdf(const char *name) const
Retrieve p.d.f (RooAbsPdf) with given name. A null pointer is returned if not found.
const Double_t * GetArray() const
Definition TArrayD.h:43
Class to manage histogram axis.
Definition TAxis.h:30
Bool_t IsVariableBinSize() const
Definition TAxis.h:136
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition TAxis.cxx:478
const TArrayD * GetXbins() const
Definition TAxis.h:130
Double_t GetXmax() const
Definition TAxis.h:134
Double_t GetXmin() const
Definition TAxis.h:133
Int_t GetNbins() const
Definition TAxis.h:121
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
TAxis * GetZaxis()
Definition TH1.h:322
virtual Int_t GetNbinsY() const
Definition TH1.h:297
virtual Double_t GetBinError(Int_t bin) const
Return value of error associated to bin number bin.
Definition TH1.cxx:8893
virtual Int_t GetNbinsZ() const
Definition TH1.h:298
virtual void Reset(Option_t *option="")
Reset this histogram: contents, errors, etc.
Definition TH1.cxx:7058
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition TH1.h:320
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
Definition TH1.cxx:2741
virtual Int_t GetNbinsX() const
Definition TH1.h:296
virtual void SetBinError(Int_t bin, Double_t error)
Set the bin Error Note that this resets the bin eror option to be of Normal Type and for the non-empt...
Definition TH1.cxx:9036
TAxis * GetYaxis()
Definition TH1.h:321
Bool_t IsBinUnderflow(Int_t bin, Int_t axis=0) const
Return true if the bin is underflow.
Definition TH1.cxx:5147
Bool_t IsBinOverflow(Int_t bin, Int_t axis=0) const
Return true if the bin is overflow.
Definition TH1.cxx:5115
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition TH1.cxx:9052
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:4994
TObject * Next()
void Reset()
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Mother of all ROOT objects.
Definition TObject.h:41
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:200
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:963
Stopwatch class.
Definition TStopwatch.h:28
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Basic string class.
Definition TString.h:136
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:692
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2336
A TTree represents a columnar dataset.
Definition TTree.h:79
RooCmdArg RecycleConflictNodes(Bool_t flag=kTRUE)
RooCmdArg Rename(const char *suffix)
RooCmdArg Conditional(const RooArgSet &pdfSet, const RooArgSet &depSet, Bool_t depsAreCond=kFALSE)
RooCmdArg WeightVar(const char *name, Bool_t reinterpretAsWeight=kFALSE)
RooCmdArg Index(RooCategory &icat)
RooCmdArg Import(const char *state, TH1 &histo)
RooCmdArg Name(const char *name)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition Common.h:18
@ ObjectHandling
Namespace for the RooStats classes.
Definition Asimov.h:19
Definition tree.py:1
void ws()
Definition ws.C:66