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