Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
HistFactoryNavigation.cxx
Go to the documentation of this file.
1
2/** \class RooStats::HistFactory::HistFactoryNavigation
3 * \ingroup HistFactory
4 */
5
6#include <iomanip>
7#include <sstream>
8
9#include "TFile.h"
10#include "TRegexp.h"
11#include "TMath.h"
12
13#include "RooRealSumPdf.h"
14#include "RooProduct.h"
15#include "RooMsgService.h"
16#include "RooCategory.h"
17#include "RooSimultaneous.h"
18#include "RooWorkspace.h"
19#include "RooFit/ModelConfig.h"
20
23
24
26
27
28namespace RooStats {
29 namespace HistFactory {
30
31
32 // CONSTRUCTOR
34 : _minBinToPrint(-1), _maxBinToPrint(-1),
35 _label_print_width(20), _bin_print_width(12) {
36
37 if( !mc ) {
38 std::cout << "Error: The supplied ModelConfig is nullptr " << std::endl;
39 throw hf_exc();
40 }
41
42 // Save the model pointer
43 RooAbsPdf* pdf_in_mc = mc->GetPdf();
44 if( !pdf_in_mc ) {
45 std::cout << "Error: The pdf found in the ModelConfig: " << mc->GetName()
46 << " is nullptr" << std::endl;
47 throw hf_exc();
48 }
49
50 // Set the PDF member
51 fModel = mc->GetPdf();
52
53 // Get the observables
54 RooArgSet* observables_in_mc = const_cast<RooArgSet*>(mc->GetObservables());
55 if( !observables_in_mc ) {
56 std::cout << "Error: Observable set in the ModelConfig: " << mc->GetName()
57 << " is nullptr" << std::endl;
58 throw hf_exc();
59 }
60 if( observables_in_mc->empty() ) {
61 std::cout << "Error: Observable list: " << observables_in_mc->GetName()
62 << " found in ModelConfig: " << mc->GetName()
63 << " has no entries." << std::endl;
64 throw hf_exc();
65 }
66
67 // Set the observables member
68 fObservables = observables_in_mc;
69
70 // Initialize the rest of the members
72
73 }
74
75
76 // CONSTRUCTOR
78 const std::string& WorkspaceName,
79 const std::string& ModelConfigName) :
80 _minBinToPrint(-1), _maxBinToPrint(-1),
81 _label_print_width(20), _bin_print_width(12) {
82
83 // Open the File
84 auto file = std::make_unique<TFile>(FileName.c_str());
85 if( !file ) {
86 std::cout << "Error: Failed to open file: " << FileName << std::endl;
87 throw hf_exc();
88 }
89
90 // Get the workspace
91 RooWorkspace* wspace = static_cast<RooWorkspace*>(file->Get(WorkspaceName.c_str()));
92 if( !wspace ) {
93 std::cout << "Error: Failed to get workspace: " << WorkspaceName
94 << " from file: " << FileName << std::endl;
95 throw hf_exc();
96 }
97
98 // Get the ModelConfig
99 ModelConfig* mc = static_cast<ModelConfig*>(wspace->obj(ModelConfigName));
100 if( !mc ) {
101 std::cout << "Error: Failed to find ModelConfig: " << ModelConfigName
102 << " from workspace: " << WorkspaceName
103 << " in file: " << FileName << std::endl;
104 throw hf_exc();
105 }
106
107 // Save the model pointer
108 RooAbsPdf* pdf_in_mc = mc->GetPdf();
109 if( !pdf_in_mc ) {
110 std::cout << "Error: The pdf found in the ModelConfig: " << ModelConfigName
111 << " is nullptr" << std::endl;
112 throw hf_exc();
113 }
114
115 // Set the PDF member
116 fModel = pdf_in_mc;
117
118 // Get the observables
119 RooArgSet* observables_in_mc = const_cast<RooArgSet*>(mc->GetObservables());
120 if( !observables_in_mc ) {
121 std::cout << "Error: Observable set in the ModelConfig: " << ModelConfigName
122 << " is nullptr" << std::endl;
123 throw hf_exc();
124 }
125 if( observables_in_mc->empty() ) {
126 std::cout << "Error: Observable list: " << observables_in_mc->GetName()
127 << " found in ModelConfig: " << ModelConfigName
128 << " in file: " << FileName
129 << " has no entries." << std::endl;
130 throw hf_exc();
131 }
132
133 // Set the observables member
134 fObservables = observables_in_mc;
135
136 // Initialize the rest of the members
138
139 }
140
141
142 // CONSTRUCTOR
144 _minBinToPrint(-1), _maxBinToPrint(-1),
145 _label_print_width(20), _bin_print_width(12) {
146
147 // Save the model pointer
148 if( !model ) {
149 std::cout << "Error: The supplied pdf is nullptr" << std::endl;
150 throw hf_exc();
151 }
152
153 // Set the PDF member
154 fModel = model;
155 fObservables = observables;
156
157 // Get the observables
158 if( !observables ) {
159 std::cout << "Error: Supplied Observable set is nullptr" << std::endl;
160 throw hf_exc();
161 }
162 if( observables->empty() ) {
163 std::cout << "Error: Observable list: " << observables->GetName()
164 << " has no entries." << std::endl;
165 throw hf_exc();
166 }
167
168 // Initialize the rest of the members
170
171 }
172
173
174 void HistFactoryNavigation::PrintMultiDimHist(TH1* hist, int bin_print_width) {
175
176 // This is how ROOT makes us loop over histograms :(
177 int current_bin = 0;
178 int num_bins = hist->GetNbinsX()*hist->GetNbinsY()*hist->GetNbinsZ();
179 for(int i = 0; i < num_bins; ++i) {
180 // Avoid the overflow/underflow
181 current_bin++;
182 while( hist->IsBinUnderflow(current_bin) ||
183 hist->IsBinOverflow(current_bin) ) {
184 current_bin++;
185 }
186 // Check that we should print this bin
187 if( _minBinToPrint != -1 && i < _minBinToPrint) continue;
188 if( _maxBinToPrint != -1 && i > _maxBinToPrint) break;
189 std::cout << std::setw(bin_print_width) << hist->GetBinContent(current_bin);
190 }
191 std::cout << std::endl;
192
193 }
194
195
196
197 RooAbsPdf* HistFactoryNavigation::GetChannelPdf(const std::string& channel) {
198
199 std::map< std::string, RooAbsPdf* >::iterator itr;
200 itr = fChannelPdfMap.find(channel);
201
202 if( itr == fChannelPdfMap.end() ) {
203 std::cout << "Warning: Could not find channel: " << channel
204 << " in pdf: " << fModel->GetName() << std::endl;
205 return nullptr;
206 }
207
208 RooAbsPdf* pdf = itr->second;
209 if( pdf == nullptr ) {
210 std::cout << "Warning: Pdf associated with channel: " << channel
211 << " is nullptr" << std::endl;
212 return nullptr;
213 }
214
215 return pdf;
216
217 }
218
219 void HistFactoryNavigation::PrintState(const std::string& channel) {
220
221 //int label_print_width = 20;
222 //int bin_print_width = 12;
223 std::cout << std::endl << channel << ":" << std::endl;
224
225 // Get the map of Samples for this channel
226 std::map< std::string, RooAbsReal*> SampleFunctionMap = GetSampleFunctionMap(channel);
227
228 // Set the size of the print width if necessary
229 /*
230 for( std::map< std::string, RooAbsReal*>::iterator itr = SampleFunctionMap.begin();
231 itr != SampleFunctionMap.end(); ++itr) {
232 std::string sample_name = itr->first;
233 label_print_width = TMath::Max(label_print_width, (int)sample_name.size()+2);
234 }
235 */
236
237 // Loop over the SampleFunctionMap and print the individual histograms
238 // to get the total histogram for the channel
239 int num_bins = 0;
240 std::map< std::string, RooAbsReal*>::iterator itr = SampleFunctionMap.begin();
241 for( ; itr != SampleFunctionMap.end(); ++itr) {
242
243 std::string sample_name = itr->first;
244 std::string tmp_name = sample_name + channel + "_pretty_tmp";
245 std::unique_ptr<TH1> sample_hist{GetSampleHist(channel, sample_name, tmp_name)};
246 num_bins = sample_hist->GetNbinsX()*sample_hist->GetNbinsY()*sample_hist->GetNbinsZ();
247 std::cout << std::setw(_label_print_width) << sample_name;
248
249 // Print the content of the histogram
250 PrintMultiDimHist(sample_hist.get(), _bin_print_width);
251
252 }
253
254 // Make the line break as a set of "===============" ...
255 std::string line_break;
256 int high_bin = _maxBinToPrint==-1 ? num_bins : TMath::Min(_maxBinToPrint, (int)num_bins);
257 int low_bin = _minBinToPrint==-1 ? 1 : _minBinToPrint;
258 int break_length = (high_bin - low_bin + 1) * _bin_print_width;
259 break_length += _label_print_width;
260 for(int i = 0; i < break_length; ++i) {
261 line_break += "=";
262 }
263 std::cout << line_break << std::endl;
264
265 std::string tmp_name = channel + "_pretty_tmp";
266 std::unique_ptr<TH1> channel_hist{GetChannelHist(channel, tmp_name)};
267 std::cout << std::setw(_label_print_width) << "TOTAL:";
268
269 // Print the Histogram
270 PrintMultiDimHist(channel_hist.get(), _bin_print_width);
271
272 return;
273
274 }
275
276
278 // Loop over channels and print their states, one after another
279 for(unsigned int i = 0; i < fChannelNameVec.size(); ++i) {
281 }
282 }
283
284
285 void HistFactoryNavigation::SetPrintWidths(const std::string& channel) {
286
287 // Get the map of Samples for this channel
288 std::map< std::string, RooAbsReal*> SampleFunctionMap = GetSampleFunctionMap(channel);
289
290 // Get the max of the samples
291 for( std::map< std::string, RooAbsReal*>::iterator itr = SampleFunctionMap.begin();
292 itr != SampleFunctionMap.end(); ++itr) {
293 std::string sample_name = itr->first;
294 _label_print_width = TMath::Max(_label_print_width, (int)sample_name.size()+2);
295 }
296
297 _label_print_width = TMath::Max( _label_print_width, (int)channel.size() + 7);
298 }
299
300
302 const std::string& channel_to_print) {
303
304 // Print the contents of a 'HistFactory' RooDataset
305 // These are stored in a somewhat odd way that makes
306 // them difficult to inspect for humans.
307 // They have the following layout:
308 // =====================================================
309 // ChannelA ChannelB ChannelCat Weight
310 // -----------------------------------------------------
311 // bin_1_center 0 ChannelA bin_1_height
312 // bin_2_center 0 ChannelA bin_2_height
313 // 0 bin_1_center ChannelB bin_1_height
314 // 0 bin_2_center ChannelB bin_2_height
315 // ...etc...
316 // =====================================================
317
318 // int label_print_width = 20;
319 // int bin_print_width = 12;
320
321 // Get the Data Histogram for this channel
322 for( unsigned int i_chan=0; i_chan < fChannelNameVec.size(); ++i_chan) {
323
324 std::string channel_name = fChannelNameVec.at(i_chan);
325
326 // If we pass a channel string, we only print that one channel
327 if( !channel_to_print.empty() && channel_name != channel_to_print) continue;
328
329 std::unique_ptr<TH1> data_hist{GetDataHist(data, channel_name, channel_name+"_tmp")};
330 std::cout << std::setw(_label_print_width) << channel_name + " (data)";
331
332 // Print the Histogram
333 PrintMultiDimHist(data_hist.get(), _bin_print_width);
334 }
335 }
336
337
339 // Loop over all channels and print model
340 // (including all samples) and compare
341 // it to the supplied dataset
342
343 for( unsigned int i = 0; i < fChannelNameVec.size(); ++i) {
344 std::string channel = fChannelNameVec.at(i);
345 SetPrintWidths(channel);
346 PrintState(channel);
347 PrintDataSet(data, channel);
348 }
349
350 std::cout << std::endl;
351 }
352
353
354 void HistFactoryNavigation::PrintParameters(bool IncludeConstantParams) {
355
356 // Get the list of parameters
357 RooArgSet params;
359
360 std::cout << std::endl;
361
362 // Create the title row
363 std::cout << std::setw(30) << "Parameter";
364 std::cout << std::setw(15) << "Value"
365 << std::setw(15) << "Error Low"
366 << std::setw(15) << "Error High"
367 << std::endl;
368
369 // Loop over the parameters and print their values, etc
370 for (auto const *param : static_range_cast<RooRealVar *>(params)) {
371 if( !IncludeConstantParams && param->isConstant() ) continue;
372
373 std::cout << std::setw(30) << param->GetName();
374 std::cout << std::setw(15) << param->getVal();
375 if( !param->isConstant() ) {
376 std::cout << std::setw(15) << param->getErrorLo() << std::setw(15) << param->getErrorHi();
377 }
378 std::cout<< std::endl;
379 }
380 std::cout << std::endl;
381 }
382
383 void HistFactoryNavigation::PrintChannelParameters(const std::string& channel,
384 bool IncludeConstantParams) {
385
386 // Get the list of parameters
387 RooArgSet params;
389
390 // Get the pdf for this channel
391 RooAbsPdf* channel_pdf = GetChannelPdf(channel);
392
393 std::cout << std::endl;
394
395 // Create the title row
396 std::cout << std::setw(30) << "Parameter";
397 std::cout << std::setw(15) << "Value"
398 << std::setw(15) << "Error Low"
399 << std::setw(15) << "Error High"
400 << std::endl;
401
402 // Loop over the parameters and print their values, etc
403 for (auto const *param : static_range_cast<RooRealVar *>(params)) {
404 if( !IncludeConstantParams && param->isConstant() ) continue;
405 if( findChild(param->GetName(), channel_pdf)==nullptr ) continue;
406 std::cout << std::setw(30) << param->GetName();
407 std::cout << std::setw(15) << param->getVal();
408 if( !param->isConstant() ) {
409 std::cout << std::setw(15) << param->getErrorLo() << std::setw(15) << param->getErrorHi();
410 }
411 std::cout<< std::endl;
412 }
413 std::cout << std::endl;
414 }
415
416
417 void HistFactoryNavigation::PrintSampleParameters(const std::string& channel,
418 const std::string& sample,
419 bool IncludeConstantParams) {
420
421 // Get the list of parameters
422 RooArgSet params;
424
425 // Get the pdf for this channel
426 RooAbsReal* sample_func = SampleFunction(channel, sample);
427
428 std::cout << std::endl;
429
430 // Create the title row
431 std::cout << std::setw(30) << "Parameter";
432 std::cout << std::setw(15) << "Value"
433 << std::setw(15) << "Error Low"
434 << std::setw(15) << "Error High"
435 << std::endl;
436
437 // Loop over the parameters and print their values, etc
438 for (auto const *param : static_range_cast<RooRealVar *>(params)) {
439 if( !IncludeConstantParams && param->isConstant() ) continue;
440 if( findChild(param->GetName(), sample_func)==nullptr ) continue;
441 std::cout << std::setw(30) << param->GetName();
442 std::cout << std::setw(15) << param->getVal();
443 if( !param->isConstant() ) {
444 std::cout << std::setw(15) << param->getErrorLo() << std::setw(15) << param->getErrorHi();
445 }
446 std::cout<< std::endl;
447 }
448 std::cout << std::endl;
449 }
450
451
452 double HistFactoryNavigation::GetBinValue(int bin, const std::string& channel) {
453 // Get the total bin height for the ith bin (ROOT indexing convention)
454 // in channel 'channel'
455 // (Could be optimized, it uses an intermediate histogram for now...)
456
457 // Get the histogram, fetch the bin content, and return
458 std::unique_ptr<TH1> channel_hist_tmp{GetChannelHist(channel, channel+"_tmp")};
459 return channel_hist_tmp->GetBinContent(bin);
460 }
461
462
463 double HistFactoryNavigation::GetBinValue(int bin, const std::string& channel, const std::string& sample){
464 // Get the total bin height for the ith bin (ROOT indexing convention)
465 // in channel 'channel'
466 // (This will be slow if you plan on looping over it.
467 // Could be optimized, it uses an intermediate histogram for now...)
468
469 // Get the histogram, fetch the bin content, and return
470 std::unique_ptr<TH1> sample_hist_tmp{GetSampleHist(channel, sample, channel+"_tmp")};
471 return sample_hist_tmp->GetBinContent(bin);
472 }
473
474
475 std::map< std::string, RooAbsReal*> HistFactoryNavigation::GetSampleFunctionMap(const std::string& channel) {
476 // Get a map of strings to function pointers,
477 // which each function corresponds to a sample
478
479 std::map< std::string, std::map< std::string, RooAbsReal*> >::iterator channel_itr;
480 channel_itr = fChannelSampleFunctionMap.find(channel);
481 if( channel_itr==fChannelSampleFunctionMap.end() ){
482 std::cout << "Error: Channel: " << channel << " not found in Navigation" << std::endl;
483 throw hf_exc();
484 }
485
486 return channel_itr->second;
487 }
488
489
490 RooAbsReal* HistFactoryNavigation::SampleFunction(const std::string& channel, const std::string& sample){
491 // Return the function object pointer corresponding
492 // to a particular sample in a particular channel
493
494 std::map< std::string, std::map< std::string, RooAbsReal*> >::iterator channel_itr;
495 channel_itr = fChannelSampleFunctionMap.find(channel);
496 if( channel_itr==fChannelSampleFunctionMap.end() ){
497 std::cout << "Error: Channel: " << channel << " not found in Navigation" << std::endl;
498 throw hf_exc();
499 }
500
501 std::map< std::string, RooAbsReal*>& SampleMap = channel_itr->second;
502 std::map< std::string, RooAbsReal*>::iterator sample_itr;
503 sample_itr = SampleMap.find(sample);
504 if( sample_itr==SampleMap.end() ){
505 std::cout << "Error: Sample: " << sample << " not found in Navigation" << std::endl;
506 throw hf_exc();
507 }
508
509 return sample_itr->second;
510
511 }
512
513
515 // Get the observables for a particular channel
516
517 std::map< std::string, RooArgSet*>::iterator channel_itr;
518 channel_itr = fChannelObservMap.find(channel);
519 if( channel_itr==fChannelObservMap.end() ){
520 std::cout << "Error: Channel: " << channel << " not found in Navigation" << std::endl;
521 throw hf_exc();
522 }
523
524 return channel_itr->second;
525
526 }
527
528
529 TH1* HistFactoryNavigation::GetSampleHist(const std::string& channel, const std::string& sample,
530 const std::string& hist_name) {
531 // Get a histogram of the expected values for
532 // a particular sample in a particular channel
533 // Give a name, or a default one will be used
534
535 RooArgList observable_list( *GetObservableSet(channel) );
536
537 std::string name = hist_name;
538 if(hist_name.empty()) name = channel + "_" + sample + "_hist";
539
540 RooAbsReal* sample_function = SampleFunction(channel, sample);
541
542 return MakeHistFromRooFunction( sample_function, observable_list, name );
543
544 }
545
546
547 TH1* HistFactoryNavigation::GetChannelHist(const std::string& channel, const std::string& hist_name) {
548 // Get a histogram of the total expected value
549 // per bin for this channel
550 // Give a name, or a default one will be used
551
552 RooArgList observable_list( *GetObservableSet(channel) );
553
554 std::map< std::string, RooAbsReal*> SampleFunctionMap = GetSampleFunctionMap(channel);
555
556 // Okay, 'loop' once
557 TH1* total_hist = nullptr;
558 std::map< std::string, RooAbsReal*>::iterator itr = SampleFunctionMap.begin();
559 for( ; itr != SampleFunctionMap.end(); ++itr) {
560 std::string sample_name = itr->first;
561 std::string tmp_hist_name = sample_name + "_hist_tmp";
562 RooAbsReal* sample_function = itr->second;
563 std::unique_ptr<TH1> sample_hist{MakeHistFromRooFunction(sample_function, observable_list,
564 tmp_hist_name)};
565 total_hist = static_cast<TH1*>(sample_hist->Clone("TotalHist"));
566 break;
567 }
568 if (!total_hist)
569 return nullptr;
570
571 total_hist->Reset();
572
573 // Loop over the SampleFunctionMap and add up all the histograms
574 // to get the total histogram for the channel
575 itr = SampleFunctionMap.begin();
576 for( ; itr != SampleFunctionMap.end(); ++itr) {
577 std::string sample_name = itr->first;
578 std::string tmp_hist_name = sample_name + "_hist_tmp";
579 RooAbsReal* sample_function = itr->second;
580 std::unique_ptr<TH1> sample_hist{MakeHistFromRooFunction(sample_function, observable_list,
581 tmp_hist_name)};
582 total_hist->Add(sample_hist.get());
583 }
584
585 if(hist_name.empty()) total_hist->SetName(hist_name.c_str());
586 else total_hist->SetName( (channel + "_hist").c_str() );
587
588 return total_hist;
589
590 }
591
592
593 std::vector< std::string > HistFactoryNavigation::GetChannelSampleList(const std::string& channel) {
594
595 std::vector<std::string> sample_list;
596
597 std::map< std::string, RooAbsReal*> sample_map = fChannelSampleFunctionMap[channel];
598 std::map< std::string, RooAbsReal*>::iterator itr = sample_map.begin();;
599 for( ; itr != sample_map.end(); ++itr) {
600 sample_list.push_back( itr->first );
601 }
602
603 return sample_list;
604
605 }
606
607
609 const std::string& name) {
610
611 THStack* stack = new THStack(name.c_str(), "");
612
613 std::vector< std::string > samples = GetChannelSampleList(channel);
614
615 // Add the histograms
616 for( unsigned int i=0; i < samples.size(); ++i) {
617 std::string sample_name = samples.at(i);
618 TH1* hist = GetSampleHist(channel, sample_name, sample_name+"_tmp");
619 hist->SetLineColor(2+i);
620 hist->SetFillColor(2+i);
621 stack->Add(hist);
622 }
623
624 return stack;
625
626 }
627
628
630 const std::string& name) {
631
632 // TO DO:
633 // MAINTAIN THE ACTUAL RANGE, USING THE OBSERVABLES
634 // MAKE IT WORK FOR MULTI-DIMENSIONAL
635 //
636
637 std::unique_ptr<TList> dataset_list;
638
639 // If the dataset covers multiple categories,
640 // Split the dataset based on the categories
641 if(strcmp(fModel->ClassName(),"RooSimultaneous")==0){
642
643 // If so, get a list of the component pdf's:
644 RooSimultaneous* simPdf = static_cast<RooSimultaneous*>(fModel);
645 auto channelCat = static_cast<RooCategory const*>(&simPdf->indexCat());
646
647 dataset_list = std::unique_ptr<TList>{data->split(*channelCat)};
648
649 data = dynamic_cast<RooDataSet*>(dataset_list->FindObject(channel.c_str()) );
650
651 }
652
653 RooArgList vars( *GetObservableSet(channel) );
654
655 int dim = vars.size();
656
657 TH1* hist = nullptr;
658
659 if (!data) {
660 std::cout << "Error: To Create Histogram from RooDataSet" << std::endl;
661 if (dataset_list) {
662 dataset_list.reset();
663 }
664 throw hf_exc();
665 } else if( dim==1 ) {
666 RooRealVar* varX = static_cast<RooRealVar*>(vars.at(0));
667 hist = data->createHistogram( name.c_str(),*varX, RooFit::Binning(varX->getBinning()) );
668 }
669 else if( dim==2 ) {
670 RooRealVar* varX = static_cast<RooRealVar*>(vars.at(0));
671 RooRealVar* varY = static_cast<RooRealVar*>(vars.at(1));
672 hist = data->createHistogram( name.c_str(),*varX, RooFit::Binning(varX->getBinning()),
673 RooFit::YVar(*varY, RooFit::Binning(varY->getBinning())) );
674 }
675 else if( dim==3 ) {
676 RooRealVar* varX = static_cast<RooRealVar*>(vars.at(0));
677 RooRealVar* varY = static_cast<RooRealVar*>(vars.at(1));
678 RooRealVar* varZ = static_cast<RooRealVar*>(vars.at(2));
679 hist = data->createHistogram( name.c_str(),*varX, RooFit::Binning(varX->getBinning()),
680 RooFit::YVar(*varY, RooFit::Binning(varY->getBinning())),
681 RooFit::YVar(*varZ, RooFit::Binning(varZ->getBinning())) );
682 }
683 else {
684 std::cout << "Error: To Create Histogram from RooDataSet, Dimension must be 1, 2, or 3" << std::endl;
685 std::cout << "Observables: " << std::endl;
686 vars.Print("V");
687 if (dataset_list) {
688 dataset_list->Delete();
689 dataset_list.reset();
690 }
691 throw hf_exc();
692 }
693
694 if (dataset_list) {
695 dataset_list->Delete();
696 dataset_list.reset();
697 }
698
699 return hist;
700
701 }
702
703
704 void HistFactoryNavigation::DrawChannel(const std::string& channel, RooDataSet* data) {
705
706 // Get the stack
707 THStack* stack = GetChannelStack(channel, channel+"_stack_tmp");
708
709 stack->Draw();
710
711 if( data!=nullptr ) {
712 TH1* data_hist = GetDataHist(data, channel, channel+"_data_tmp");
713 data_hist->Draw("SAME");
714 }
715
716 }
717
718
719
721
722 // An internal method to recursively get all products,
723 // including if a RooProduct is a Product of RooProducts
724 // etc
725
726 RooArgSet allTerms;
727
728 // Get All Subnodes of this product
729 // Loop over the subnodes and add
730 for (auto *arg : node->components()) {
731 std::string ClassName = arg->ClassName();
732 if( ClassName == "RooProduct" ) {
733 RooProduct* prod = static_cast<RooProduct*>(arg);
734 allTerms.add( _GetAllProducts(prod) );
735 }
736 else {
737 allTerms.add(*arg);
738 }
739 }
740 return allTerms;
741 }
742
743
744
745
746 void HistFactoryNavigation::_GetNodes(RooAbsPdf* modelPdf, const RooArgSet* observables) {
747
748 // Get the pdf from the ModelConfig
749 //RooAbsPdf* modelPdf = mc->GetPdf();
750 //RooArgSet* observables = mc->GetObservables();
751
752 // Check if it is a simultaneous pdf or not
753 // (if it's an individual channel, it won't be, if it's
754 // combined, it's simultaneous)
755 // Fill the channel vectors based on the structure
756 // (Obviously, if it's not simultaneous, there will be
757 // only one entry in the vector for the single channel)
758 if(strcmp(modelPdf->ClassName(),"RooSimultaneous")==0){
759
760 // If so, get a list of the component pdf's:
761 auto simPdf = static_cast<RooSimultaneous*>(modelPdf);
762 auto channelCat = static_cast<RooCategory const*>(&simPdf->indexCat());
763
764 // Iterate over the categories and get the
765 // pdf and observables for each category
766 for (const auto& nameIdx : *channelCat) {
767 const std::string& ChannelName = nameIdx.first;
768 fChannelNameVec.push_back( ChannelName );
769 RooAbsPdf* pdftmp = simPdf->getPdf(ChannelName.c_str()) ;
770 RooArgSet* obstmp = std::unique_ptr<RooArgSet>{pdftmp->getObservables(*observables)}.release();
771 fChannelPdfMap[ChannelName] = pdftmp;
772 fChannelObservMap[ChannelName] = obstmp;
773 }
774
775 } else {
776 RooArgSet* obstmp = std::unique_ptr<RooArgSet>{modelPdf->getObservables(*observables)}.release();
777 // The channel name is model_CHANNEL
778 std::string ChannelName = modelPdf->GetName();
779 ChannelName = ChannelName.replace(0, 6, "");
780 fChannelNameVec.push_back(ChannelName);
781 fChannelPdfMap[ChannelName] = modelPdf;
782 fChannelObservMap[ChannelName] = obstmp;
783
784 }
785
786 // Okay, now we have maps of the pdfs
787 // and the observable list per channel
788 // We then loop over the channel pdfs:
789 // and find their RooRealSumPdfs
790 // std::map< std::string, RooRealSumPdf* > channelSumNodeMap;
791
792 for( unsigned int i = 0; i < fChannelNameVec.size(); ++i ) {
793
794 std::string ChannelName = fChannelNameVec.at(i);
795 RooAbsPdf* pdf = fChannelPdfMap[ChannelName];
796
797 // Loop over the pdf's components and find
798 // the (one) that is a RooRealSumPdf
799 // Based on the mode, we assume that node is
800 // the "unconstrained" pdf node for that channel
801 std::unique_ptr<RooArgSet> comps{pdf->getComponents()};
802 for (auto *arg : *comps) {
803 std::string ClassName = arg->ClassName();
804 if( ClassName == "RooRealSumPdf" ) {
805 fChannelSumNodeMap[ChannelName] = static_cast <RooRealSumPdf*>(arg);
806 break;
807 }
808 }
809 }
810
811 // Okay, now we have all necessary
812 // nodes filled for each channel.
813 for( unsigned int i = 0; i < fChannelNameVec.size(); ++i ) {
814
815 std::string ChannelName = fChannelNameVec.at(i);
816 RooRealSumPdf* sumPdf = dynamic_cast<RooRealSumPdf*>(fChannelSumNodeMap[ChannelName]);
817
818 // We now take the RooRealSumPdf and loop over
819 // its component functions. The RooRealSumPdf turns
820 // a list of functions (expected events or bin heights
821 // per sample) and turns it into a pdf.
822 // Therefore, we loop over it to find the expected
823 // height for the various samples
824
825 // First, create a map to store the function nodes
826 // for each sample in this channel
827 std::map< std::string, RooAbsReal*> sampleFunctionMap;
828
829 // Loop over the sample nodes in this
830 // channel's RooRealSumPdf
831 for (auto *func : static_range_cast<RooAbsReal *>(sumPdf->funcList())) {
832 // Do a bit of work to get the name of each sample
833 std::string SampleName = func->GetName();
834 if( SampleName.find("L_x_") != std::string::npos ) {
835 size_t index = SampleName.find("L_x_");
836 SampleName.replace( index, 4, "" );
837 }
838 if( SampleName.find(ChannelName) != std::string::npos ) {
839 size_t index = SampleName.find(ChannelName);
840 SampleName = SampleName.substr(0, index-1);
841 }
842
843 // And simply save this node into our map
844 sampleFunctionMap[SampleName] = func;
845 }
846
847 fChannelSampleFunctionMap[ChannelName] = sampleFunctionMap;
848
849 // Okay, now we have a list of histograms
850 // representing the samples for this channel.
851
852 }
853
854 }
855
856
857 RooAbsArg* HistFactoryNavigation::findChild(const std::string& name, RooAbsReal* parent) const {
858
859 RooAbsArg* term=nullptr;
860
861 // Check if it is a "component",
862 // ie a sub node:
863 std::unique_ptr<RooArgSet> comps{parent->getComponents()};
864 for (auto arg : *comps) {
865 std::string ArgName = arg->GetName();
866 if (ArgName == name) {
867 term = arg;
868 break;
869 }
870 }
871
872 if( term != nullptr ) return term;
873
874 // If that failed,
875 // Check if it's a Parameter
876 // (ie a RooRealVar)
877 RooArgSet args;
878 std::unique_ptr<RooArgSet> parameters{parent->getParameters(&args)};
879 for (auto *param : *parameters) {
880 std::string ParamName = param->GetName();
881 if( ParamName == name ) {
882 term = param;
883 break;
884 }
885 }
886
887 /* Not sure if we want to be silent
888 But since we're returning a pointer which can be nullptr,
889 I think it's the user's job to do checks on it.
890 A dereference will always cause a crash, so it won't
891 be silent for long...
892 if( term==nullptr ) {
893 std::cout << "Error: Failed to find node: " << name
894 << " as a child of: " << parent->GetName()
895 << std::endl;
896 }
897 */
898
899 return term;
900
901 }
902
903
905
906 std::string ConstraintTermName = parameter + "Constraint";
907
908 // First, as a sanity check, let's see if the parameter
909 // itself actually exists and if the model depends on it:
910 RooRealVar* param = dynamic_cast<RooRealVar*>(findChild(parameter, fModel));
911 if( param==nullptr ) {
912 std::cout << "Error: Couldn't Find parameter: " << parameter << " in model."
913 << std::endl;
914 return nullptr;
915 }
916
917 // The "gamma"'s use a different constraint term name
918 if( parameter.find("gamma_stat_") != std::string::npos ) {
919 ConstraintTermName = parameter + "_constraint";
920 }
921
922 // Now, get the constraint itself
923 RooAbsReal* term = dynamic_cast<RooAbsReal*>(findChild(ConstraintTermName, fModel));
924
925 if( term==nullptr ) {
926 std::cout << "Error: Couldn't Find constraint term for parameter: " << parameter
927 << " (Looked for '" << ConstraintTermName << "')" << std::endl;
928 return nullptr;
929 }
930
931 return term;
932
933 }
934
935
936 double HistFactoryNavigation::GetConstraintUncertainty(const std::string& parameter) {
937
938 RooAbsReal* constraintTerm = GetConstraintTerm(parameter);
939 if( constraintTerm==nullptr ) {
940 std::cout << "Error: Cannot get uncertainty because parameter: " << parameter
941 << " has no constraint term" << std::endl;
942 throw hf_exc();
943 }
944
945 // Get the type of constraint
946 std::string ConstraintType = constraintTerm->ClassName();
947
948 // Find its value
949 double sigma = 0.0;
950
951 if( ConstraintType.empty() ) {
952 std::cout << "Error: Constraint type is an empty string."
953 << " This simply should not be." << std::endl;
954 throw hf_exc();
955 }
956 else if( ConstraintType == "RooGaussian" ){
957
958 // Gaussian errors are the 'sigma' in the constraint term
959
960 // Get the name of the 'sigma' for the gaussian
961 // (I don't know of a way of doing RooGaussian::GetSigma() )
962 // For alpha's, the sigma points to a global RooConstVar
963 // with the name "1"
964 // For gamma_stat_*, the sigma is named *_sigma
965 std::string sigmaName;
966 if( parameter.find("alpha_")!=std::string::npos ) {
967 sigmaName = "1";;
968 }
969 else if( parameter.find("gamma_stat_")!=std::string::npos ) {
970 sigmaName = parameter + "_sigma";
971 }
972
973 // Get the sigma and its value
974 RooAbsReal* sigmaVar = dynamic_cast<RooAbsReal*>(constraintTerm->findServer(sigmaName.c_str()));
975 if( sigmaVar==nullptr ) {
976 std::cout << "Error: Failed to find the 'sigma' node: " << sigmaName
977 << " in the RooGaussian: " << constraintTerm->GetName() << std::endl;
978 throw hf_exc();
979 }
980 // If we find the uncertainty:
981 sigma = sigmaVar->getVal();
982 }
983 else if( ConstraintType == "RooPoisson" ){
984 // Poisson errors are given by inverting: tau = 1 / (sigma*sigma)
985 std::string tauName = "nom_" + parameter;
986 RooAbsReal* tauVar = dynamic_cast<RooAbsReal*>(constraintTerm->findServer(tauName.c_str()) );
987 if( tauVar==nullptr ) {
988 std::cout << "Error: Failed to find the nominal 'tau' node: " << tauName
989 << " for the RooPoisson: " << constraintTerm->GetName() << std::endl;
990 throw hf_exc();
991 }
992 double tau_val = tauVar->getVal();
993 sigma = 1.0 / TMath::Sqrt( tau_val );
994 }
995 else {
996 std::cout << "Error: Encountered unknown constraint type for Stat Uncertainties: "
997 << ConstraintType << std::endl;
998 throw hf_exc();
999 }
1000
1001 return sigma;
1002
1003 }
1004
1005 void HistFactoryNavigation::ReplaceNode(const std::string& ToReplace, RooAbsArg* ReplaceWith) {
1006
1007 // First, check that the node to replace is actually a node:
1008 RooAbsArg* nodeToReplace = findChild(ToReplace, fModel);
1009 if( nodeToReplace==nullptr ) {
1010 std::cout << "Error: Cannot replace node: " << ToReplace
1011 << " because this node wasn't found in: " << fModel->GetName()
1012 << std::endl;
1013 throw hf_exc();
1014 }
1015
1016
1017 const std::string attrib = "ORIGNAME:" + ToReplace;
1018 const bool oldAttrib = ReplaceWith->getAttribute(attrib.c_str());
1019 ReplaceWith->setAttribute(attrib.c_str());
1020 RooArgSet newServerSet{*ReplaceWith};
1021
1022 // Now that we have the node we want to replace, we have to
1023 // get its parent node
1024
1025 // Do this by looping over the clients and replacing their servers
1026 // (NOTE: This happens for ALL clients across the pdf)
1027 for (auto client : nodeToReplace->clients()) {
1028
1029 // Check if this client is a member of our pdf
1030 // (We probably don't want to mess with clients
1031 // if they aren't...)
1032 if( findChild(client->GetName(), fModel) == nullptr) continue;
1033
1034 // Now, do the replacement:
1035 client->redirectServers(newServerSet);
1036 std::cout << "Replaced: " << ToReplace << " with: " << ReplaceWith->GetName()
1037 << " in node: " << client->GetName() << std::endl;
1038
1039 }
1040
1041 // reset temporary attribute for server redirection
1042 ReplaceWith->setAttribute(attrib.c_str(), oldAttrib);
1043
1044 return;
1045
1046 }
1047
1048
1049 void HistFactoryNavigation::PrintSampleComponents(const std::string& channel,
1050 const std::string& sample) {
1051
1052 // Get the Sample Node
1053 RooAbsReal* sampleNode = SampleFunction(channel, sample);
1054
1055 // Get the observables for this channel
1056 RooArgList observable_list( *GetObservableSet(channel) );
1057
1058 // Make the total histogram for this sample
1059 std::string total_Name = sampleNode->GetName();
1060 std::unique_ptr<TH1> total_hist{MakeHistFromRooFunction( sampleNode, observable_list, total_Name + "_tmp")};
1061 unsigned int num_bins = total_hist->GetNbinsX()*total_hist->GetNbinsY()*total_hist->GetNbinsZ();
1062
1063 RooArgSet components;
1064
1065 // Let's see what it is...
1066 int label_print_width = 30;
1067 int bin_print_width = 12;
1068 if( strcmp(sampleNode->ClassName(),"RooProduct")==0){
1069 RooProduct* prod = dynamic_cast<RooProduct*>(sampleNode);
1070 components.add( _GetAllProducts(prod) );
1071 }
1072 else {
1073 components.add(*sampleNode);
1074 }
1075
1076 /////// NODE SIZE
1077 {
1078 for (auto *component : static_range_cast<RooAbsReal*>(components)) {
1079 std::string NodeName = component->GetName();
1080 label_print_width = TMath::Max(label_print_width, static_cast<int>(NodeName.size())+2);
1081 }
1082 }
1083
1084 // Now, loop over the components and print them out:
1085 std::cout << std::endl;
1086 std::cout << "Channel: " << channel << " Sample: " << sample << std::endl;
1087 std::cout << std::setw(label_print_width) << "Factor";
1088
1089 for(unsigned int i=0; i < num_bins; ++i) {
1090 if( _minBinToPrint != -1 && (int)i < _minBinToPrint) continue;
1091 if( _maxBinToPrint != -1 && (int)i > _maxBinToPrint) break;
1092 std::stringstream sstr;
1093 sstr << "Bin" << i;
1094 std::cout << std::setw(bin_print_width) << sstr.str();
1095 }
1096 std::cout << std::endl;
1097
1098
1099 for (auto *component : static_range_cast<RooAbsReal*>(components)) {
1100 std::string NodeName = component->GetName();
1101 // Make a histogram for this node
1102 // Do some horrible things to prevent some really
1103 // annoying messages from being printed
1106 std::unique_ptr<TH1> hist;
1107 hist.reset(MakeHistFromRooFunction( component, observable_list, NodeName+"_tmp"));
1109
1110 // Print the hist
1111 std::cout << std::setw(label_print_width) << NodeName;
1112
1113 // Print the Histogram
1114 PrintMultiDimHist(hist.get(), bin_print_width);
1115 }
1116 /////
1117 std::string line_break;
1118 int high_bin = _maxBinToPrint==-1 ? num_bins : TMath::Min(_maxBinToPrint, (int)num_bins);
1119 int low_bin = _minBinToPrint==-1 ? 1 : _minBinToPrint;
1120 int break_length = (high_bin - low_bin + 1) * bin_print_width;
1121 break_length += label_print_width;
1122 for(int i = 0; i < break_length; ++i) {
1123 line_break += "=";
1124 }
1125 std::cout << line_break << std::endl;
1126
1127 std::cout << std::setw(label_print_width) << "TOTAL:";
1128 PrintMultiDimHist(total_hist.get(), bin_print_width);
1129 }
1130
1131
1133 std::string name ) {
1134
1135 // Turn a RooAbsReal* into a TH1* based
1136 // on a template histogram.
1137 // The 'vars' arg list defines the x (and y and z variables)
1138 // Loop over the bins of the Template,
1139 // find the bin centers,
1140 // Scan the input Var over those bin centers,
1141 // and use the value of the function
1142 // to make the new histogram
1143
1144 // Make the new histogram
1145 // Cone and empty the template
1146 // TH1* hist = (TH1*) histTemplate.Clone( name.c_str() );
1147
1148 int dim = vars.size();
1149
1150 TH1* hist=nullptr;
1151
1152 if( dim==1 ) {
1153 RooRealVar* varX = static_cast<RooRealVar*>(vars.at(0));
1154 hist = func->createHistogram( name.c_str(),*varX, RooFit::Binning(varX->getBinning()), RooFit::Scaling(false) );
1155 }
1156 else if( dim==2 ) {
1157 RooRealVar* varX = static_cast<RooRealVar*>(vars.at(0));
1158 RooRealVar* varY = static_cast<RooRealVar*>(vars.at(1));
1159 hist = func->createHistogram( name.c_str(),*varX, RooFit::Binning(varX->getBinning()), RooFit::Scaling(false),
1160 RooFit::YVar(*varY, RooFit::Binning(varY->getBinning())) );
1161 }
1162 else if( dim==3 ) {
1163 RooRealVar* varX = static_cast<RooRealVar*>(vars.at(0));
1164 RooRealVar* varY = static_cast<RooRealVar*>(vars.at(1));
1165 RooRealVar* varZ = static_cast<RooRealVar*>(vars.at(2));
1166 hist = func->createHistogram( name.c_str(),*varX, RooFit::Binning(varX->getBinning()), RooFit::Scaling(false),
1167 RooFit::YVar(*varY, RooFit::Binning(varY->getBinning())),
1168 RooFit::YVar(*varZ, RooFit::Binning(varZ->getBinning())) );
1169 }
1170 else {
1171 std::cout << "Error: To Create Histogram from RooAbsReal function, Dimension must be 1, 2, or 3" << std::endl;
1172 throw hf_exc();
1173 }
1174
1175 return hist;
1176 }
1177
1178 // A simple wrapper to use a ModelConfig
1180 RooAbsPdf* modelPdf = mc->GetPdf();
1181 const RooArgSet* observables = mc->GetObservables();
1182 _GetNodes(modelPdf, observables);
1183 }
1184
1185
1186 void HistFactoryNavigation::SetConstant(const std::string& regExpr, bool constant) {
1187
1188 // Regex FTW
1189
1190 TString RegexTString(regExpr);
1191 TRegexp theRegExpr(RegexTString);
1192
1193 // Now, loop over all variables and
1194 // set the constant as
1195
1196 // Get the list of parameters
1197 RooArgSet params;
1199
1200 std::cout << std::endl;
1201
1202 // Create the title row
1203 std::cout << std::setw(30) << "Parameter";
1204 std::cout << std::setw(15) << "Value"
1205 << std::setw(15) << "Error Low"
1206 << std::setw(15) << "Error High"
1207 << std::endl;
1208
1209 // Loop over the parameters and print their values, etc
1210 for (auto *param : static_range_cast<RooRealVar *>(params)) {
1211
1212 std::string ParamName = param->GetName();
1213 TString ParamNameTString(ParamName);
1214
1215 // Use the Regex to skip all parameters that don't match
1216 //if( theRegExpr.Index(ParamNameTString, ParamName.size()) == -1 ) continue;
1217 Ssiz_t dummy;
1218 if( theRegExpr.Index(ParamNameTString, &dummy) == -1 ) continue;
1219
1220 param->setConstant( constant );
1221 std::cout << "Setting param: " << ParamName << " constant"
1222 << " (matches regex: " << regExpr << ")" << std::endl;
1223 }
1224 }
1225
1226 RooRealVar* HistFactoryNavigation::var(const std::string& varName) const {
1227
1228 RooAbsArg* arg = findChild(varName, fModel);
1229 if( !arg ) return nullptr;
1230
1231 RooRealVar* var_obj = dynamic_cast<RooRealVar*>(arg);
1232 return var_obj;
1233
1234 }
1235
1236 } // namespace HistFactory
1237} // namespace RooStats
1238
1239
1240
1241
#define ClassImp(name)
Definition Rtypes.h:377
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
char name[80]
Definition TGX11.cxx:110
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
bool redirectServers(const RooAbsCollection &newServerList, bool mustReplaceAll=false, bool nameChange=false, bool isRecursionStep=false)
Replace all direct servers of this object with the new servers in newServerList.
TIterator Use clients() and begin()
RooFit::OwningPtr< RooArgSet > getComponents() const
Create a RooArgSet with all components (branch nodes) of the expression tree headed by this object.
bool getAttribute(const Text_t *name) const
Check if a named attribute is set. By default, all attributes are unset.
void setAttribute(const Text_t *name, bool value=true)
Set (default) or clear a named boolean attribute of this object.
RooAbsArg * findServer(const char *name) const
Return server of this with name name. Returns nullptr if not found.
Definition RooAbsArg.h:210
const char * GetName() const override
Returns name of object.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
RooAbsArg * find(const char *name) const
Find object with given name in list.
void Print(Option_t *options=nullptr) const override
This method must be overridden when a class wants to print itself.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
TH1 * createHistogram(RooStringView varNameList, Int_t xbins=0, Int_t ybins=0, Int_t zbins=0) const
Create and fill a ROOT histogram TH1, TH2 or TH3 with the values of this function for the variables w...
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
Object to represent discrete states.
Definition RooCategory.h:28
Container class to hold unbinned data.
Definition RooDataSet.h:57
static RooMsgService & instance()
Return reference to singleton instance.
void setGlobalKillBelow(RooFit::MsgLevel level)
RooFit::MsgLevel globalKillBelow() const
Represents the product of a given set of RooAbsReal objects.
Definition RooProduct.h:29
RooArgList components()
Definition RooProduct.h:48
Implements a PDF constructed from a sum of functions:
const RooArgList & funcList() const
Variable that can be changed from the outside.
Definition RooRealVar.h:37
const RooAbsBinning & getBinning(const char *name=nullptr, bool verbose=true, bool createOnTheFly=false) const override
Return binning definition with name.
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
const RooAbsCategoryLValue & indexCat() const
void PrintParameters(bool IncludeConstantParams=false)
Print the current values and errors of pdf parameters.
std::vector< std::string > fChannelNameVec
The list of channels.
RooArgSet * GetObservableSet(const std::string &channel)
Get the set of observables for a given channel.
TH1 * GetDataHist(RooDataSet *data, const std::string &channel, const std::string &name="")
Get a histogram from the dataset for this channel.
void _GetNodes(ModelConfig *mc)
Fetch the node information for the pdf in question, and save it in the various collections in this cl...
std::map< std::string, RooAbsPdf * > fChannelSumNodeMap
Map of channel names to pdf without constraint.
void ReplaceNode(const std::string &ToReplace, RooAbsArg *ReplaceWith)
Find a node in the pdf and replace it with a new node These nodes can be functions,...
void PrintState()
Should pretty print all channels and the current values

void PrintSampleComponents(const std::string &channel, const std::string &sample)
Print the different components that make up a sample (NormFactors, Statistical Uncertainties,...
void SetPrintWidths(const std::string &channel)
Set the title and bin widths.
RooAbsPdf * fModel
The HistFactory Pdf Pointer.
TH1 * GetChannelHist(const std::string &channel, const std::string &name="")
Get the total channel histogram for this channel.
TH1 * GetSampleHist(const std::string &channel, const std::string &sample, const std::string &name="")
The (current) histogram for that sample This includes all parameters and interpolation.
TH1 * MakeHistFromRooFunction(RooAbsReal *func, RooArgList vars, std::string name="Hist")
Make a histogram from a function Edit so it can take a RooArgSet of parameters.
void PrintMultiDimHist(TH1 *hist, int bin_print_width)
Print a histogram's contents to the screen void PrettyPrintHistogram(TH1* hist);.
void PrintSampleParameters(const std::string &channel, const std::string &sample, bool IncludeConstantParams=false)
Print parameters that effect a particular sample.
std::map< std::string, RooAbsReal * > GetSampleFunctionMap(const std::string &channel)
Get a map of sample names to their functions for a particular channel.
RooAbsReal * SampleFunction(const std::string &channel, const std::string &sample)
Get the RooAbsReal function for a given sample in a given channel.
HistFactoryNavigation(ModelConfig *mc)
Initialize based on an already-created HistFactory Model.
double GetBinValue(int bin, const std::string &channel)
The value of the ith bin for the total in that channel.
RooAbsArg * findChild(const std::string &name, RooAbsReal *parent) const
Internal method implementation of finding a daughter node from a parent node (looping over all genera...
std::map< std::string, std::map< std::string, RooAbsReal * > > fChannelSampleFunctionMap
Map of Map of Channel, Sample names to Function Nodes Used by doing: fChannelSampleFunctionMap["MyCha...
std::map< std::string, RooAbsPdf * > fChannelPdfMap
Map of channel names to their full pdf's.
std::vector< std::string > GetChannelSampleList(const std::string &channel)
void PrintModelAndData(RooDataSet *data)
Print the model and the data, comparing channel by channel.
std::map< std::string, RooArgSet * > fChannelObservMap
Map of channel names to their set of ovservables.
RooAbsPdf * GetChannelPdf(const std::string &channel)
void PrintDataSet(RooDataSet *data, const std::string &channel="")
Print a "HistFactory style" RooDataSet in a readable way.
void PrintChannelParameters(const std::string &channel, bool IncludeConstantParams=false)
Print parameters that effect a particular channel.
void SetConstant(const std::string &regExpr=".*", bool constant=true)
RooArgSet _GetAllProducts(RooProduct *node)
Recursively get all products of products.
THStack * GetChannelStack(const std::string &channel, const std::string &name="")
Get a stack of all samples in a channel.
RooAbsReal * GetConstraintTerm(const std::string &parameter)
Get the constraint term for a given systematic (alpha or gamma)
void DrawChannel(const std::string &channel, RooDataSet *data=nullptr)
Draw a stack of the channel, and include data if the pointer is supplied.
RooRealVar * var(const std::string &varName) const
double GetConstraintUncertainty(const std::string &parameter)
Get the uncertainty based on the constraint term for a given systematic.
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition ModelConfig.h:35
const RooArgSet * GetObservables() const
get RooArgSet for observables (return nullptr if not existing)
RooAbsPdf * GetPdf() const
get model PDF (return nullptr if pdf has not been specified or does not exist)
Persistable container for RooFit projects.
TObject * obj(RooStringView name) const
Return any type of object (RooAbsArg, RooAbsData or generic object) with given name)
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
virtual Int_t GetNbinsY() const
Definition TH1.h:298
virtual Int_t GetNbinsZ() const
Definition TH1.h:299
virtual void Reset(Option_t *option="")
Reset this histogram: contents, errors, etc.
Definition TH1.cxx:7071
virtual Int_t GetNbinsX() const
Definition TH1.h:297
virtual Bool_t Add(TF1 *h1, Double_t c1=1, Option_t *option="")
Performs the operation: this = this + c1*f1 if errors are defined (see TH1::Sumw2),...
Definition TH1.cxx:826
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3066
Bool_t IsBinUnderflow(Int_t bin, Int_t axis=0) const
Return true if the bin is underflow.
Definition TH1.cxx:5182
Bool_t IsBinOverflow(Int_t bin, Int_t axis=0) const
Return true if the bin is overflow.
Definition TH1.cxx:5150
void SetName(const char *name) override
Change the name of this histogram.
Definition TH1.cxx:8928
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:5029
The Histogram stack class.
Definition THStack.h:40
virtual void Add(TH1 *h, Option_t *option="")
add a new histogram to the list Only 1-d and 2-d histograms currently supported.
Definition THStack.cxx:365
void Draw(Option_t *chopt="") override
Draw this multihist with its current attributes.
Definition THStack.cxx:450
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
Regular expression class.
Definition TRegexp.h:31
Ssiz_t Index(const TString &str, Ssiz_t *len, Ssiz_t start=0) const
Find the first occurrence of the regexp in string and return the position, or -1 if there is no match...
Definition TRegexp.cxx:213
Basic string class.
Definition TString.h:139
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg={})
RooCmdArg Scaling(bool flag)
RooCmdArg Binning(const RooAbsBinning &binning)
const Double_t sigma
MsgLevel
Verbosity level for RooMsgService::StreamConfig in RooMsgService.
Namespace for the RooStats classes.
Definition Asimov.h:19
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:250
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:662
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198