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