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