Logo ROOT   6.16/01
Reference Guide
RooSimPdfBuilder.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofitcore:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17//////////////////////////////////////////////////////////////////////////////
18//
19// Begin_Html
20//
21// <b>This tool has now been superceded by RooSimWSTool</b>
22//
23// <p>
24// <tt>RooSimPdfBuilder</tt> is a powerful tool to build <tt>RooSimultaneous</tt>
25// PDFs that are defined in terms component PDFs that are identical in
26// structure, but have different parameters.
27// </p>
28//
29// <h2>Example</h2>
30//
31// <p>
32// The following example demonstrates the essence of <tt>RooSimPdfBuilder</tt>:
33// Given a dataset D with a <tt>RooRealVar X</tt> and a <tt>RooCategory C</tt> that has
34// state C1 and C2.
35// <ul>
36// <li> We want to fit the distribution of <tt>X</tt> with a Gaussian+ArgusBG PDF,
37// <li> We want to fit the data subsets <tt>D(C==C1)</tt> and <tt>D(C==C2)</tt> separately and simultaneously.
38// <li> The PDFs to fit data subsets D_C1 and D_C2 are identical except for
39// <ul>
40// <li> the kappa parameter of the ArgusBG PDF and
41// <li> the sigma parameter of the gaussian PDF
42// </ul>
43// where each PDF will have its own copy of the parameter
44// </ul>
45// </p>
46// <p>
47// Coding this example directly with RooFit classes gives
48// (we assume dataset D and variables C and X have been declared previously)
49// </p>
50// <pre>
51// RooRealVar m("m","mean of gaussian",-10,10) ;
52// RooRealVar s_C1("s_C1","sigma of gaussian C1",0,20) ;
53// RooRealVar s_C2("s_C2","sigma of gaussian C2",0,20) ;
54// RooGaussian gauss_C1("gauss_C1","gaussian C1",X,m,s_C1) ;
55// RooGaussian gauss_C2("gauss_C2","gaussian C2",X,m,s_C2) ;
56//
57// RooRealVar k_C1("k_C1","ArgusBG kappa parameter C1",-50,0) ;
58// RooRealVar k_C2("k_C2","ArgusBG kappa parameter C2",-50,0) ;
59// RooRealVar xm("xm","ArgusBG cutoff point",5.29) ;
60// RooArgusBG argus_C1("argus_C1","argus background C1",X,k_C1,xm) ;
61// RooArgusBG argus_C2("argus_C2","argus background C2",X,k_C2,xm) ;
62//
63// RooRealVar gfrac("gfrac","fraction of gaussian",0.,1.) ;
64// RooAddPdf pdf_C1("pdf_C1","gauss+argus_C1",RooArgList(gauss_C1,argus_C1),gfrac) ;
65// RooAddPdf pdf_C2("pdf_C2","gauss+argus_C2",RooArgList(gauss_C2,argus_C2),gfrac) ;
66//
67// RooSimultaneous simPdf("simPdf","simPdf",C) ;
68// simPdf.addPdf(pdf_C1,"C1") ;
69// simPdf.addPdf(pdf_C2,"C2") ;
70// </pre>
71// <p>
72// Coding this example with RooSimPdfBuilder gives
73// </p>
74// <pre>
75// RooRealVar m("m","mean of gaussian",-10,10) ;
76// RooRealVar s("s","sigma of gaussian",0,20) ;
77// RooGaussian gauss("gauss","gaussian",X,m,s) ;
78//
79// RooRealVar k("k","ArgusBG kappa parameter",-50,0) ;
80// RooRealVar xm("xm","ArgusBG cutoff point",5.29) ;
81// RooArgusBG argus("argus","argus background",X,k,xm) ;
82//
83// RooRealVar gfrac("gfrac","fraction of gaussian",0.,1.) ;
84// RooAddPdf pdf("pdf","gauss+argus",RooArgList(gauss,argus),gfrac) ;
85//
86// RooSimPdfBuilder builder(pdf) ;
87// RooArgSet* config = builder.createProtoBuildConfig() ;
88// (*config)["physModels"] = "pdf" ; // Name of the PDF we are going to work with
89// (*config)["splitCats"] = "C" ; // Category used to differentiate sub-datasets
90// (*config)["pdf"] = "C : k,s" ; // Prescription to taylor PDF parameters k and s
91// // for each data subset designated by C states
92// RooSimultaneous* simPdf = builder.buildPdf(*config,&D) ;
93// </pre>
94// <p>
95// The above snippet of code demonstrates the concept of <tt>RooSimPdfBuilder</tt>:
96// the user defines a single <i>'prototype' PDF</i> that defines the structure of all
97// PDF components of the <tt>RooSimultaneous</tt> PDF to be built. <tt>RooSimPdfBuilder</tt>
98// then takes this prototype and replicates it as a component
99// PDF for each state of the C index category.
100// </p>
101// <p>
102// In the above example </tt>RooSimPdfBuilder</tt>
103// will first replicate <tt>k</tt> and <tt>s</tt> into
104// <tt>k_C1,k_C2</tt> and <tt>s_C1,s_C2</tt>, as prescribed in the
105// configuration. Then it will recursively replicate all PDF nodes that depend on
106// the 'split' parameter nodes: <tt>gauss</tt> into <tt>gauss_C1,C2</tt>, <tt>argus</tt>
107// into <tt>argus_C1,C2</tt> and finally <tt>pdf</tt> into <tt>pdf_C1,pdf_C2</tt>.
108// When PDFs for all states of C have been replicated
109// they are assembled into a <tt>RooSimultaneous</tt> PDF, which is returned by the <tt>buildPdf()</tt>
110// method.
111// </p>
112// <p>
113// Although in this very simple example the use of <tt>RooSimPdfBuilder</tt> doesn't
114// reduce the amount of code much, it is already easier to read and maintain
115// because there is no duplicate code. As the complexity of the <tt>RooSimultaneous</tt>
116// to be built increases, the advantages of <tt>RooSimPdfBuilder</tt> will become more and
117// more apparent.
118// </p>
119//
120//
121// <h2>Builder configuration rules for a single prototype PDF</h2>
122// <p>
123// Each builder configuration needs at minumum two lines, <tt>physModels</tt> and <tt>splitCats</tt>, which identify
124// the ingredients of the build. In this section we only explain the building rules for
125// builds from a single prototype PDF. In that case the <tt>physModels</tt> line always reads
126// </p>
127// <pre>
128// physModels = {pdfName}
129// </pre>
130// <p>
131// The second line, <tt>splitCats</tt>, indicates which categories are going to be used to
132// differentiate the various subsets of the 'master' input data set. You can enter
133// a single category here, or multiple if necessary:
134// </p>
135// <pre>
136// splitCats = {catName} [{catName} ...]
137// </pre>
138// <p>
139// All listed splitcats must be <tt>RooCategories</tt> that appear in the dataset provided to
140// <tt>RooSimPdfBuilder::buildPdf()</tt>
141// </p>
142// <p>
143// The parameter splitting prescriptions, the essence of each build configuration
144// can be supplied in a third line carrying the name of the pdf listed in <tt>physModels</tt>
145// </p>
146// <pre>
147// pdfName = {splitCat} : {parameter} [,{parameter},....]
148// </pre>
149// <p>
150// Each pdf can have only one line with splitting rules, but multiple rules can be
151// supplied in each line, e.g.
152// </p>
153// <pre>
154// pdfName = {splitCat} : {parameter} [,{parameter},....]
155// {splitCat} : {parameter} [,{parameter},....]
156// </pre>
157// <p>
158// Conversely, each parameter can only have one splitting prescription, but it may be split
159// by multiple categories, e.g.
160// </p>
161// <pre>
162// pdfName = {splitCat1},{splitCat2} : {parameter}
163// </pre>
164// <p>
165// instructs <tt>RooSimPdfBuilder</tt> to build a <tt>RooSuperCategory</tt>
166// of <tt>{splitCat1}</tt> and <tt>{splitCat2}</tt>
167// and split <tt>{parameter}</tt> with that <tt>RooSuperCategory</tt>
168// </p>
169// <p>
170// Here is an example of a builder configuration that uses several of the options discussed
171// above:
172// </p>
173// <pre>
174// physModels = pdf
175// splitCats = tagCat runBlock
176// pdf = tagCat : signalRes,bkgRes
177// runBlock : fudgeFactor
178// tagCat,runBlock : kludgeParam
179// </pre>
180//
181// <h2>How to enter configuration data</h2>
182//
183// <p>
184// The prototype builder configuration returned by
185// <tt>RooSimPdfBuilder::createProtoBuildConfig()</tt> is a pointer to a <tt>RooArgSet</tt> filled with
186// initially blank <tt>RooStringVars</tt> named <tt>physModels,splitCats</tt> and one additional for each
187// PDF supplied to the <tt>RooSimPdfBuilders</tt> constructor (with the same name)
188// </p>
189// <p>
190// In macro code, the easiest way to assign new values to these <tt>RooStringVars</tt>
191// is to use <tt>RooArgSet</tt>s array operator and the <tt>RooStringVar</tt>s assignment operator, e.g.
192// </p>
193// <pre>
194// (*config)["physModels"] = "Blah" ;
195// </pre>
196// <p>
197// To enter multiple splitting rules simply separate consecutive rules by whitespace
198// (not newlines), e.g.
199// </p>
200// <pre>
201// (*config)["physModels"] = "Blah " // << note trailing space here
202// "Blah 2" ;
203// </pre>
204// <p>
205// In this example, the C++ compiler will concatenate the two string literals (without inserting
206// any whitespace), so the extra space after 'Blah' is important here.
207// </p>
208// <p>
209// Alternatively, you can read the configuration from an ASCII file, as you can
210// for any <tt>RooArgSet</tt> using <tt>RooArgSet::readFromFile()</tt>. In that case the ASCII file
211// can follow the syntax of the examples above and the '<tt>\</tt>' line continuation
212// sequence can be used to fold a long splitting rule over multiple lines.
213// </p>
214// <pre>
215// RooArgSet* config = builder.createProtoBuildConfig() ;
216// config->readFromFile("config.txt") ;
217//
218// --- config.txt ----------------
219// physModels = pdf
220// splitCats = tagCat
221// pdf = tagCat : bogusPar
222// -------------------------------
223// </pre>
224//
225//
226// <h2>Working with multiple prototype PDFs</h2>
227// <p>
228// It is also possible to build a <tt>RooSimultaneous</tt> PDF from multiple PDF prototypes.
229// This is appropriate for cases where the input prototype PDF would otherwise be
230// a <tt>RooSimultaneous</tt> PDF by itself. In such cases we don't feed a single
231// <tt>RooSimultaneous</tt> PDF into <tt>RooSimPdfBuilder</tt>, instead we feed it its ingredients and
232// add a prescription to the builder configuration that corresponds to the
233// PDF-category state mapping of the prototype <tt>RooSimultaneous</tt>.
234// </p>
235// <p>
236// The constructor of the <tt>RooSimPdfBuilder</tt> will look as follows:
237// </p>
238// <pre>
239// RooSimPdfBuilder builder(RooArgSet(pdfA,pdfB,...)) ;
240// </pre>
241// <p>
242// The <tt>physModels</tt> line is now expanded to carry the pdf->state mapping information
243// that the prototype <tt>RooSimultaneous</tt> would have. I.e.
244// </p>
245// <pre>
246// physModels = mode : pdfA=modeA pdfB=modeB
247// </pre>
248// <p>
249// is equivalent to a prototype <tt>RooSimultaneous</tt> constructed as
250// </p>
251// <pre>
252// RooSimultanous simPdf("simPdf","simPdf",mode);
253// simPdf.addPdf(pdfA,"modeA") ;
254// simPdf.addPdf(pdfB,"modeB") ;
255// </pre>
256// <p>
257// The rest of the builder configuration works the same, except that
258// each prototype PDF now has its own set of splitting rules, e.g.
259// </p>
260// <pre>
261// physModels = mode : pdfA=modeA pdfB=modeB
262// splitCats = tagCat
263// pdfA = tagCat : bogusPar
264// pdfB = tagCat : fudgeFactor
265// </pre>
266// <p>
267// Please note that
268// <ul>
269// <li> The master index category ('mode' above) doesn't have to be listed in
270// <tt>splitCats</tt>, this is implicit.
271//
272// <li> The number of splitting prescriptions goes by the
273// number of prototype PDFs and not by the number of states of the
274// master index category (mode in the above and below example).
275// </ul>
276//
277// In the following case:
278//</p>
279// <pre>
280// physModels = mode : pdfA=modeA pdfB=modeB pdfA=modeC pdfB=modeD
281// </pre>
282// <p>
283// there are still only 2 sets of splitting rules: one for <tt>pdfA</tt> and one
284// for <tt>pdfB</tt>. However, you <i>can</i> differentiate between <tt>modeA</tt> and <tt>modeC</tt> in
285// the above example. The technique is to use <tt>mode</tt> as splitting category, e.g.
286// </p>
287// <pre>
288// physModels = mode : pdfA=modeA pdfB=modeB pdfA=modeC pdfB=modeD
289// splitCats = tagCat
290// pdfA = tagCat : bogusPar
291// mode : funnyPar
292// pdfB = mode : kludgeFactor
293// </pre>
294// <p>
295// will result in an individual set of <tt>funnyPar</tt> parameters for <tt>modeA</tt> and <tt>modeC</tt>
296// labeled <tt>funnyPar_modeA</tt> and <tt>funnyPar_modeB</tt> and an individual set of
297// kludgeFactor parameters for <tt>pdfB</tt>, <tt>kludgeFactor_modeB</tt> and <tt>kludgeFactor_modeD</tt>.
298// Please note that for splits in the master index category (mode) only the
299// applicable states are built (A,C for <tt>pdfA</tt>, B,D for <tt>pdfB</tt>)
300// </p>
301//
302//
303// <h2>Advanced options</h2>
304//
305// <h4>Partial splits</h4>
306// <p>
307// You can request to limit the list of states of each splitCat that
308// will be considered in the build. This limitation is requested in the
309// each build as follows:
310// </p>
311// <pre>
312// splitCats = tagCat(Lep,Kao) RunBlock(Run1)
313// </pre>
314// <p>
315// In this example the splitting of <tt>tagCat</tt> is limited to states <tt>Lep,Kao</tt>
316// and the splitting of <tt>runBlock</tt> is limited to <tt>Run1</tt>. The splits apply
317// globally to each build, i.e. every parameter split requested in this
318// build will be limited according to these specifications.
319// </p>
320// <p>
321// NB: Partial builds have no pdf associated with the unbuilt states of the
322// limited splits. Running such a pdf on a dataset that contains data with
323// unbuilt states will result in this data being ignored completely.
324// </p>
325//
326//
327// <h4>Non-trivial splits</h4>
328// <p>
329// It is possible to make non-trivial parameter splits with <tt>RooSimPdfBuilder</tt>.
330// Trivial splits are considered simple splits in one (fundamental) category
331// in the dataset or a split in a <tt>RooSuperCategory</tt> 'product' of multiple
332// fundamental categories in the dataset. Non-trivial splits can be performed
333// using an intermediate 'category function' (<tt>RooMappedCategory,
334// RooGenericCategory,RooThresholdCategory</tt> etc), i.e. any <tt>RooAbsCategory</tt>
335// derived objects that calculates its output as function of one or more
336// input <tt>RooRealVars</tt> and/or <tt>RooCategories</tt>.
337// </p>
338// <p>
339// Such 'function categories' objects must be constructed by the user prior
340// to building the PDF. In the <tt>RooSimPdfBuilder::buildPdf()</tt> function these
341// objects can be passed in an optional <tt>RooArgSet</tt> called 'auxiliary categories':
342// </p>
343// <pre>
344// const <tt>RooSimultaneous</tt>* buildPdf(const RooArgSet& buildConfig, const RooAbsData* dataSet,
345// const RooArgSet& auxSplitCats, Bool_t verbose=kFALSE) {
346// ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
347// </pre>
348// <p>
349// Objects passed in this argset can subsequently be used in the build configuration, e.g.
350// </p>
351// <pre>
352// RooMappedCategory tagMap("tagMap","Mapped tagging category",tagCat,"CutBased") ;
353// tagMap.map("Lep","CutBased") ;
354// tagMap.map("Kao","CutBased") ;
355// tagMap.map("NT*","NeuralNet") ;
356// ...
357// builder.buildPdf(config,D,tagMap) ;
358// ^^^^^^
359//<Contents of config>
360// physModels = pdf
361// splitCats = tagCat runBlock
362// pdf = tagCat : signalRes
363// tagMap : fudgeFactor
364// ^^^^^^
365// </pre>
366// <p>
367// In the above example <tt>signalRes</tt> will be split in <tt>signalRes_Kao,signalRes_Lep,
368// signalRes_NT1,signalRes_NT2</tt>, while <tt>fudgeFactor</tt> will be split in <tt>fudgeFactor_CutBased</tt>
369// and <tt>fudgeFactor_NeuralNet</tt>.
370// </p>
371// <p>
372// Category functions passed in the auxSplitCats <tt>RooArgSet</tt> can be used regularly
373// in the splitting configuration. They should not be listed in <tt>splitCats</tt>,
374// but must be able to be expressed <i>completely</i> in terms of the <tt>splitCats</tt> that
375// are listed.
376// </p>
377//
378//
379// <h4>Multiple connected builds</h4>
380// <p>
381// Sometimes you want to build multiple PDFs for independent consecutive fits
382// that share some of their parameters. For example, we have two prototype PDFs
383// <tt>pdfA(x;p,q)</tt> and <tt>pdfB(x;p,r)</tt> that have a common parameter <tt>p</tt>.
384// We want to build a <tt>RooSimultaneous</tt> for both <tt>pdfA</tt> and <tt>B</tt>,
385// which involves a split of parameter <tt>p</tt> and we would like to build the
386// simultaneous pdfs </tt>simA</tt> and <tt>simB</tt> such that still share their (now split) parameters
387// <tt>p_XXX</tt>. This is accomplished by letting a single instance of <tt>RooSimPdfBuilder</tt> handle
388// the builds of both <tt>pdfA</tt> and <tt>pdfB</tt>, as illustrated in this example:
389// </p>
390// <pre>
391// RooSimPdfBuilder builder(RooArgSet(pdfA,pdfB)) ;
392//
393// RooArgSet* configA = builder.createProtoBuildConfig() ;
394// (*configA)["physModels"] = "pdfA" ;
395// (*configA)["splitCats"] = "C" ;
396// (*configA)["pdf"] = "C : p" ;
397// RooSimultaneous* simA = builder.buildPdf(*configA,&D) ;
398//
399// RooArgSet* configB = builder.createProtoBuildConfig() ;
400// (*configA)["physModels"] = "pdfB" ;
401// (*configA)["splitCats"] = "C" ;
402// (*configA)["pdf"] = "C : p" ;
403// RooSimultaneous* simB = builder.buildPdf(*configB,&D) ;
404// </pre>
405//
406// <h2>Ownership of constructed PDFs</h2>
407// <p>
408// The <tt>RooSimPdfBuilder</tt> instance owns all the objects it creates, including the top-level
409// <tt>RooSimultaneous</tt> returned by <tt>buildPdf()</tt>. Therefore the builder instance should
410// exist as long as the constructed PDFs needs to exist.
411// </p>
412//
413// End_Html
414//
415
416
417
418#include "RooFit.h"
419
420#include <string.h>
421
422#ifndef _WIN32
423#include <strings.h>
424#endif
425
426#include "Riostream.h"
427#include "RooSimPdfBuilder.h"
428
429#include "RooRealVar.h"
430#include "RooFormulaVar.h"
431#include "RooAbsCategory.h"
432#include "RooCategory.h"
433#include "RooStringVar.h"
434#include "RooMappedCategory.h"
435#include "RooRealIntegral.h"
436#include "RooDataSet.h"
437#include "RooArgSet.h"
438#include "RooPlot.h"
439#include "RooAddPdf.h"
440#include "RooLinearVar.h"
441#include "RooTruthModel.h"
442#include "RooAddModel.h"
443#include "RooProdPdf.h"
444#include "RooCustomizer.h"
445#include "RooThresholdCategory.h"
446#include "RooMultiCategory.h"
447#include "RooSuperCategory.h"
448#include "RooSimultaneous.h"
449#include "RooTrace.h"
450#include "RooFitResult.h"
451#include "RooDataHist.h"
452#include "RooGenericPdf.h"
453#include "RooMsgService.h"
454
455using namespace std ;
456
458;
459
460
461
462////////////////////////////////////////////////////////////////////////////////
463
465 _protoPdfSet(protoPdfSet)
466{
470}
471
472
473
474
475////////////////////////////////////////////////////////////////////////////////
476/// Make RooArgSet of configuration objects
477
479{
480 RooArgSet* buildConfig = new RooArgSet ;
481 buildConfig->addOwned(* new RooStringVar("physModels","List and mapping of physics models to include in build","",4096)) ;
482 buildConfig->addOwned(* new RooStringVar("splitCats","List of categories used for splitting","",1024)) ;
483
486 while ((proto=(RooAbsPdf*)iter->Next())) {
487 buildConfig->addOwned(* new RooStringVar(proto->GetName(),proto->GetName(),"",4096)) ;
488 }
489 delete iter ;
490
491 return buildConfig ;
492}
493
494
495
496////////////////////////////////////////////////////////////////////////////////
497
499{
500 _splitNodeList.add(specSet) ;
501}
502
503
504
505////////////////////////////////////////////////////////////////////////////////
506/// Initialize needed components
507
508RooSimultaneous* RooSimPdfBuilder::buildPdf(const RooArgSet& buildConfig, const RooArgSet& dependents,
509 const RooArgSet* auxSplitCats, Bool_t verbose)
510{
511 const char* spaceChars = " \t" ;
512
513 // Retrieve physics index category
514 Int_t buflen = strlen(((RooStringVar*)buildConfig.find("physModels"))->getVal())+1 ;
515 char *buf = new char[buflen] ;
516
517 strlcpy(buf,((RooStringVar*)buildConfig.find("physModels"))->getVal(),buflen) ;
518 RooAbsCategoryLValue* physCat(0) ;
519 if (strstr(buf," : ")) {
520 const char* physCatName = strtok(buf,spaceChars) ;
521 physCat = dynamic_cast<RooAbsCategoryLValue*>(dependents.find(physCatName)) ;
522 if (!physCat) {
523 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR physics index category " << physCatName
524 << " not found in dataset variables" << endl ;
525 delete[] buf ;
526 return 0 ;
527 }
528 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: category indexing physics model: " << physCatName << endl ;
529 }
530
531 // Create list of physics models to be built
532 char *physName ;
533 RooArgSet physModelSet ;
534 if (physCat) {
535 // Absorb colon token
536 strtok(0,spaceChars) ;
537 physName = strtok(0,spaceChars) ;
538 } else {
539 physName = strtok(buf,spaceChars) ;
540 }
541
542 if (!physName) {
543 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR: No models specified, nothing to do!" << endl ;
544 delete[] buf ;
545 return 0 ;
546 }
547
549 RooArgSet stateMap ;
550 while(physName) {
551
552 char *stateName(0) ;
553
554 // physName may be <state>=<pdfName> or just <pdfName> is state and pdf have identical names
555 if (strchr(physName,'=')) {
556 // Must have a physics category for mapping to make sense
557 if (!physCat) {
558 coutW(ObjectHandling) << "RooSimPdfBuilder::buildPdf: WARNING: without physCat specification "
559 << "<physCatState>=<pdfProtoName> association is meaningless" << endl ;
560 }
561 stateName = physName ;
562 physName = strchr(stateName,'=') ;
563 if (physName) {
564 *(physName++) = 0 ;
565 }
566 } else {
567 stateName = physName ;
568 }
569
570 RooAbsPdf* physModel = (RooAbsPdf*) (physName ? _protoPdfSet.find(physName) : 0 );
571 if (!physModel) {
572 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR requested physics model "
573 << (physName?physName:"(null)") << " is not defined" << endl ;
574 delete[] buf ;
575 return 0 ;
576 }
577
578 // Check if state mapping has already been defined
579 if (stateMap.find(stateName)) {
580 coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: multiple PDFs specified for state "
581 << stateName << ", only first will be used" << endl ;
582 continue ;
583 }
584
585 // Add pdf to list of models to be processed
586 physModelSet.add(*physModel,kTRUE) ; // silence duplicate insertion warnings
587
588 // Store state->pdf mapping
589 stateMap.addOwned(* new RooStringVar(stateName,stateName,physName)) ;
590
591 // Continue with next mapping
592 physName = strtok(0,spaceChars) ;
593 if (first) {
594 first = kFALSE ;
595 } else if (physCat==0) {
596 coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: without physCat specification, only the first model will be used" << endl ;
597 break ;
598 }
599 }
600 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: list of physics models " << physModelSet << endl ;
601
602
603
604 // Create list of dataset categories to be used in splitting
605 TList splitStateList ;
606 RooArgSet splitCatSet ;
607
608 delete[] buf ;
609 buflen = strlen(((RooStringVar*)buildConfig.find("splitCats"))->getVal())+1 ;
610 buf = new char[buflen] ;
611 strlcpy(buf,((RooStringVar*)buildConfig.find("splitCats"))->getVal(),buflen) ;
612
613 char *catName = strtok(buf,spaceChars) ;
614 char *stateList(0) ;
615 while(catName) {
616
617 // Chop off optional list of selected states
618 char* tokenPtr(0) ;
619 if (strchr(catName,'(')) {
620
621 catName = R__STRTOK_R(catName,"(",&tokenPtr) ;
622 stateList = R__STRTOK_R(0,")",&tokenPtr) ;
623
624 } else {
625 stateList = 0 ;
626 }
627
628 RooCategory* splitCat = catName ? dynamic_cast<RooCategory*>(dependents.find(catName)) : 0 ;
629 if (!splitCat) {
630 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR requested split category " << (catName?catName:"(null)")
631 << " is not a RooCategory in the dataset" << endl ;
632 delete[] buf ;
633 return 0 ;
634 }
635 splitCatSet.add(*splitCat) ;
636
637 // Process optional state list
638 if (stateList) {
639 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: splitting of category " << catName
640 << " restricted to states (" << stateList << ")" << endl ;
641
642 // Create list named after this splitCat holding its selected states
643 TList* slist = new TList ;
644 slist->SetName(catName) ;
645 splitStateList.Add(slist) ;
646
647 char* stateLabel = R__STRTOK_R(stateList,",",&tokenPtr) ;
648
649 while(stateLabel) {
650 // Lookup state label and require it exists
651 const RooCatType* type = splitCat->lookupType(stateLabel) ;
652 if (!type) {
653 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR splitCat " << splitCat->GetName()
654 << " doesn't have a state named " << stateLabel << endl ;
655 splitStateList.Delete() ;
656 delete[] buf ;
657 return 0 ;
658 }
659 slist->Add((TObject*)type) ;
660
661 stateLabel = R__STRTOK_R(0,",",&tokenPtr) ;
662 }
663 }
664
665 catName = strtok(0,spaceChars) ;
666 }
667 if (physCat) splitCatSet.add(*physCat) ;
668 RooSuperCategory masterSplitCat("masterSplitCat","Master splitting category",splitCatSet) ;
669
670 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: list of splitting categories " << splitCatSet << endl ;
671
672 // Clone auxiliary split cats and attach to splitCatSet
673 RooArgSet auxSplitSet ;
674 RooArgSet* auxSplitCloneSet(0) ;
675 if (auxSplitCats) {
676 // Deep clone auxililary split cats
677 auxSplitCloneSet = (RooArgSet*) auxSplitCats->snapshot(kTRUE) ;
678 if (!auxSplitCloneSet) {
679 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf(" << GetName() << ") Couldn't deep-clone set auxiliary splitcats, abort." << endl ;
680 delete[] buf ;
681 return 0 ;
682 }
683
684 TIterator* iter = auxSplitCats->createIterator() ;
685 RooAbsArg* arg ;
686 while((arg=(RooAbsArg*)iter->Next())) {
687 // Find counterpart in cloned set
688 RooAbsArg* aux = auxSplitCats->find(arg->GetName()) ;
689
690 // Check that there is no fundamental splitCat in the dataset with the bane of the auxiliary split
691 if (splitCatSet.find(aux->GetName())) {
692 coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: dataset contains a fundamental splitting category " << endl
693 << " with the same name as an auxiliary split function (" << aux->GetName() << "). " << endl
694 << " Auxiliary split function will be ignored" << endl ;
695 continue ;
696 }
697
698 // Check that all servers of this aux cat are contained in splitCatSet
699 RooArgSet* parSet = aux->getParameters(splitCatSet) ;
700 if (parSet->getSize()>0) {
701 coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: ignoring auxiliary category " << aux->GetName()
702 << " because it has servers that are not listed in splitCatSet: " << *parSet << endl ;
703 delete parSet ;
704 continue ;
705 }
706
707 // Redirect servers to splitCatSet
708 aux->recursiveRedirectServers(splitCatSet) ;
709
710 // Add top level nodes to auxSplitSet
711 auxSplitSet.add(*aux) ;
712 }
713 delete iter ;
714
715 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: list of auxiliary splitting categories " << auxSplitSet << endl ;
716 }
717
718
719 TList* customizerList = new TList ;
720
721 // Loop over requested physics models and build components
722 TIterator* physIter = physModelSet.createIterator() ;
723 RooAbsPdf* physModel ;
724 while((physModel=(RooAbsPdf*)physIter->Next())) {
725 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: processing physics model " << physModel->GetName() << endl ;
726
727 RooCustomizer* physCustomizer = new RooCustomizer(*physModel,masterSplitCat,_splitNodeList) ;
728 customizerList->Add(physCustomizer) ;
729
730 // Parse the splitting rules for this physics model
731 RooStringVar* ruleStr = (RooStringVar*) buildConfig.find(physModel->GetName()) ;
732 if (ruleStr) {
733
734 delete[] buf ;
735 buflen = strlen(ruleStr->getVal())+1 ;
736 buf = new char[buflen] ;
737
738 strlcpy(buf,ruleStr->getVal(),buflen) ;
739 char *tokenPtr(0) ;
740
741 char* token = R__STRTOK_R(buf,spaceChars,&tokenPtr) ;
742
743 enum Mode { SplitCat, Colon, ParamList } ;
744 Mode mode(SplitCat) ;
745
746 char* splitCatName ;
747 RooAbsCategory* splitCat(0) ;
748
749 while(token) {
750
751 switch (mode) {
752 case SplitCat:
753 {
754 splitCatName = token ;
755
756 if (strchr(splitCatName,',')) {
757 // Composite splitting category
758
759 // Check if already instantiated
760 splitCat = (RooAbsCategory*) _compSplitCatSet.find(splitCatName) ;
761 TString origCompCatName(splitCatName) ;
762 if (!splitCat) {
763 // Build now
764
765 char *tokptr = 0;
766 char *catName2 = R__STRTOK_R(token,",",&tokptr) ;
767
768 RooArgSet compCatSet ;
769 while(catName2) {
770 RooAbsArg* cat = splitCatSet.find(catName2) ;
771
772 // If not, check if it is an auxiliary splitcat
773 if (!cat) {
774 cat = (RooAbsCategory*) auxSplitSet.find(catName2) ;
775 }
776
777 if (!cat) {
778 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR " << catName2
779 << " not found in the primary or auxilary splitcat list" << endl ;
780 customizerList->Delete() ;
781 delete customizerList ;
782
783 splitStateList.Delete() ;
784 delete[] buf ;
785 return 0 ;
786 }
787 compCatSet.add(*cat) ;
788
789 catName2 = R__STRTOK_R(0,",",&tokptr) ;
790 }
791
792
793 // check that any auxSplitCats in compCatSet do not depend on any other
794 // fundamental or auxiliary splitting categories in the composite set.
795 TIterator* iter = compCatSet.createIterator() ;
796 RooAbsArg* arg ;
797 while((arg=(RooAbsArg*)iter->Next())) {
798 RooArgSet tmp(compCatSet) ;
799 tmp.remove(*arg) ;
800 if (arg->dependsOnValue(tmp)) {
801 coutE(InputArguments) << "RooSimPdfBuilder::buildPDF: ERROR: Ill defined split: auxiliary splitting category " << arg->GetName()
802 << " used in composite split " << compCatSet << " depends on one or more of the other splitting categories in the composite split" << endl ;
803
804 // Cleanup and axit
805 customizerList->Delete() ;
806 delete customizerList ;
807 splitStateList.Delete() ;
808 delete[] buf ;
809 return 0 ;
810 }
811 }
812 delete iter ;
813
814 splitCat = new RooMultiCategory(origCompCatName,origCompCatName,compCatSet) ;
815 _compSplitCatSet.addOwned(*splitCat) ;
816 //cout << "composite splitcat: " << splitCat->GetName() ;
817 }
818 } else {
819 // Simple splitting category
820
821 // First see if it is a simple splitting category
822 splitCat = (RooAbsCategory*) splitCatSet.find(splitCatName) ;
823
824 // If not, check if it is an auxiliary splitcat
825 if (!splitCat) {
826 splitCat = (RooAbsCategory*) auxSplitSet.find(splitCatName) ;
827 }
828
829 if (!splitCat) {
830 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR splitting category "
831 << splitCatName << " not found in the primary or auxiliary splitcat list" << endl ;
832 customizerList->Delete() ;
833 delete customizerList ;
834 splitStateList.Delete() ;
835 delete[] buf ;
836 return 0 ;
837 }
838 }
839
840 mode = Colon ;
841 break ;
842 }
843 case Colon:
844 {
845 if (strcmp(token,":")) {
846 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR in parsing, expected ':' after "
847 << splitCat << ", found " << token << endl ;
848 customizerList->Delete() ;
849 delete customizerList ;
850 splitStateList.Delete() ;
851 delete[] buf ;
852 return 0 ;
853 }
854 mode = ParamList ;
855 break ;
856 }
857 case ParamList:
858 {
859 // Verify the validity of the parameter list and build the corresponding argset
860 RooArgSet splitParamList ;
861 RooArgSet* paramList = physModel->getParameters(dependents) ;
862
863 // wve -- add nodes to parameter list
864 RooArgSet* compList = physModel->getComponents() ;
865 paramList->add(*compList) ;
866 delete compList ;
867
868 Bool_t lastCharIsComma = (token[strlen(token)-1]==',') ;
869
870 char *tokptr = 0 ;
871 char *paramName = R__STRTOK_R(token,",",&tokptr) ;
872
873 // Check for fractional split option 'param_name[remainder_state]'
874 char *remainderState = 0 ;
875 char *tokptr2 = 0 ;
876 if (paramName && R__STRTOK_R(paramName,"[",&tokptr2)) {
877 remainderState = R__STRTOK_R(0,"]",&tokptr2) ;
878 }
879
880 while(paramName) {
881
882 // If fractional split is specified, check that remainder state is a valid state of this split cat
883 if (remainderState) {
884 if (!splitCat->lookupType(remainderState)) {
885 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR fraction split of parameter "
886 << paramName << " has invalid remainder state name: " << remainderState << endl ;
887 delete paramList ;
888 customizerList->Delete() ;
889 delete customizerList ;
890 splitStateList.Delete() ;
891 delete[] buf ;
892 return 0 ;
893 }
894 }
895
896 RooAbsArg* param = paramList->find(paramName) ;
897 if (!param) {
898 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR " << paramName
899 << " is not a parameter of physics model " << physModel->GetName() << endl ;
900 delete paramList ;
901 customizerList->Delete() ;
902 delete customizerList ;
903 splitStateList.Delete() ;
904 delete[] buf ;
905 return 0 ;
906 }
907 splitParamList.add(*param) ;
908
909 // Build split leaf of fraction splits here
910 if (remainderState) {
911
912 // Check if we are splitting a real-valued parameter
913 if (!dynamic_cast<RooAbsReal*>(param)) {
914 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR fraction split requested of non-real valued parameter "
915 << param->GetName() << endl ;
916 delete paramList ;
917 customizerList->Delete() ;
918 delete customizerList ;
919 splitStateList.Delete() ;
920 delete[] buf ;
921 return 0 ;
922 }
923
924 // Check if we are doing a restricted build
925 TList* remStateSplitList = static_cast<TList*>(splitStateList.FindObject(splitCat->GetName())) ;
926
927 // If so, check if remainder state is actually being built.
928 if (remStateSplitList && !remStateSplitList->FindObject(remainderState)) {
929 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR " << paramName
930 << " remainder state " << remainderState << " in parameter split "
931 << param->GetName() << " is not actually being built" << endl ;
932 delete paramList ;
933 customizerList->Delete() ;
934 delete customizerList ;
935 splitStateList.Delete() ;
936 delete[] buf ;
937 return 0 ;
938 }
939
940 TIterator* iter = splitCat->typeIterator() ;
942 RooArgList fracLeafList ;
943 TString formExpr("1") ;
944 Int_t i(0) ;
945
946 while((type=(RooCatType*)iter->Next())) {
947
948 // Skip remainder state
949 if (!TString(type->GetName()).CompareTo(remainderState)) continue ;
950
951 // If restricted build is requested, skip states of splitcat that are not built
952 if (remStateSplitList && !remStateSplitList->FindObject(type->GetName())) {
953 continue ;
954 }
955
956 // Construct name of split leaf
957 TString splitLeafName(param->GetName()) ;
958 splitLeafName.Append("_") ;
959 splitLeafName.Append(type->GetName()) ;
960
961 // Check if split leaf already exists
962 RooAbsArg* splitLeaf = _splitNodeList.find(splitLeafName) ;
963 if (!splitLeaf) {
964 // If not create it now
965 splitLeaf = (RooAbsArg*) param->clone(splitLeafName) ;
966 _splitNodeList.add(*splitLeaf) ;
967 _splitNodeListOwned.addOwned(*splitLeaf) ;
968 }
969 fracLeafList.add(*splitLeaf) ;
970 formExpr.Append(Form("-@%d",i++)) ;
971 }
972 delete iter ;
973
974 // Construct RooFormulaVar expresssing remainder of fraction
975 TString remLeafName(param->GetName()) ;
976 remLeafName.Append("_") ;
977 remLeafName.Append(remainderState) ;
978
979 // Check if no specialization was already specified for remainder state
980 if (!_splitNodeList.find(remLeafName)) {
981 RooAbsArg* remLeaf = new RooFormulaVar(remLeafName,formExpr,fracLeafList) ;
982 _splitNodeList.add(*remLeaf) ;
983 _splitNodeListOwned.addOwned(*remLeaf) ;
984 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: creating remainder fraction formula for " << remainderState
985 << " specialization of split parameter " << param->GetName() << " " << formExpr << endl ;
986 }
987 }
988
989 // Parse next parameter name
990 paramName = R__STRTOK_R(0,",",&tokptr) ;
991 if (paramName && R__STRTOK_R(paramName,"[",&tokptr2)) {
992 remainderState = R__STRTOK_R(0,"]",&tokptr2) ;
993 }
994 }
995
996 // Add the rule to the appropriate customizer ;
997 physCustomizer->splitArgs(splitParamList,*splitCat) ;
998
999 delete paramList ;
1000
1001 if (!lastCharIsComma) mode = SplitCat ;
1002 break ;
1003 }
1004 }
1005
1006 token = R__STRTOK_R(0,spaceChars,&tokenPtr) ;
1007
1008 }
1009 if (mode!=SplitCat) {
1010 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR in parsing, expected "
1011 << (mode==Colon?":":"parameter list") << " after " << (token?token:"(null)") << endl ;
1012 }
1013
1014 //RooArgSet* paramSet = physModel->getParameters(dependents) ;
1015 } else {
1016 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: no splitting rules for " << physModel->GetName() << endl ;
1017 }
1018 }
1019
1020 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: configured customizers for all physics models" << endl ;
1021 if (oodologI((TObject*)0,ObjectHandling)) {
1022 customizerList->Print() ;
1023 }
1024
1025 // Create fit category from physCat and splitCatList ;
1026 RooArgSet fitCatList ;
1027 if (physCat) fitCatList.add(*physCat) ;
1028 fitCatList.add(splitCatSet) ;
1029 TIterator* fclIter = fitCatList.createIterator() ;
1030 RooSuperCategory *fitCat = new RooSuperCategory("fitCat","fitCat",fitCatList) ;
1031
1032 // Create master PDF
1033 RooSimultaneous* simPdf = new RooSimultaneous("simPdf","simPdf",*fitCat) ;
1034
1035 // Add component PDFs to master PDF
1036 TIterator* fcIter = fitCat->typeIterator() ;
1037
1038 RooCatType* fcState ;
1039 while((fcState=(RooCatType*)fcIter->Next())) {
1040 // Select fitCat state
1041 fitCat->setLabel(fcState->GetName()) ;
1042
1043 // Check if this fitCat state is selected
1044 fclIter->Reset() ;
1045 RooAbsCategory* splitCat ;
1046 Bool_t select(kTRUE) ;
1047 while((splitCat=(RooAbsCategory*)fclIter->Next())) {
1048 // Find selected state list
1049 TList* slist = (TList*) splitStateList.FindObject(splitCat->GetName()) ;
1050 if (!slist) continue ;
1051 RooCatType* type = (RooCatType*) slist->FindObject(splitCat->getLabel()) ;
1052 if (!type) {
1053 select = kFALSE ;
1054 }
1055 }
1056 if (!select) continue ;
1057
1058
1059 // Select appropriate PDF for this physCat state
1060 RooCustomizer* physCustomizer ;
1061 if (physCat) {
1062 RooStringVar* physNameVar = (RooStringVar*) stateMap.find(physCat->getLabel()) ;
1063 if (!physNameVar) continue ;
1064 physCustomizer = (RooCustomizer*) customizerList->FindObject(physNameVar->getVal());
1065 } else {
1066 physCustomizer = (RooCustomizer*) customizerList->First() ;
1067 }
1068
1069 coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: Customizing physics model " << physCustomizer->GetName()
1070 << " for mode " << fcState->GetName() << endl ;
1071
1072 // Customizer PDF for current state and add to master simPdf
1073 RooAbsPdf* fcPdf = (RooAbsPdf*) physCustomizer->build(masterSplitCat.getLabel(),verbose) ;
1074 simPdf->addPdf(*fcPdf,fcState->GetName()) ;
1075 }
1076 delete fcIter ;
1077
1078 // Move customizers (owning the cloned branch node components) to the attic
1079 _retiredCustomizerList.AddAll(customizerList) ;
1080 delete customizerList ;
1081
1082 delete fclIter ;
1083 splitStateList.Delete() ;
1084
1085 if (auxSplitCloneSet) delete auxSplitCloneSet ;
1086 delete physIter ;
1087
1088 delete[] buf ;
1089 _simPdfList.push_back(simPdf) ;
1090 _fitCatList.push_back(fitCat) ;
1091 return simPdf ;
1092}
1093
1094
1095
1096
1097
1098////////////////////////////////////////////////////////////////////////////////
1099
1101{
1103
1104 std::list<RooSimultaneous*>::iterator iter = _simPdfList.begin() ;
1105 while(iter != _simPdfList.end()) {
1106 delete *iter ;
1107 ++iter ;
1108 }
1109
1110 std::list<RooSuperCategory*>::iterator iter2 = _fitCatList.begin() ;
1111 while(iter2 != _fitCatList.end()) {
1112 delete *iter2 ;
1113 ++iter2 ;
1114 }
1115
1116}
1117
1118
#define coutI(a)
Definition: RooMsgService.h:31
#define coutW(a)
Definition: RooMsgService.h:33
#define coutE(a)
Definition: RooMsgService.h:34
#define oodologI(o, a)
Definition: RooMsgService.h:71
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:363
char * R__STRTOK_R(char *str, const char *delim, char **saveptr)
Definition: Rtypes.h:484
int type
Definition: TGX11.cxx:120
char * Form(const char *fmt,...)
const char * proto
Definition: civetweb.c:16604
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
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:533
Bool_t dependsOnValue(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0) const
Definition: RooAbsArg.h:88
RooArgSet * getComponents() const
Definition: RooAbsArg.cxx:680
Bool_t recursiveRedirectServers(const RooAbsCollection &newServerList, Bool_t mustReplaceAll=kFALSE, Bool_t nameChange=kFALSE, Bool_t recurseInNewSet=kTRUE)
Definition: RooAbsArg.cxx:1080
virtual TObject * clone(const char *newname=0) const =0
RooAbsCategoryLValue is the common abstract base class for objects that represent a discrete value th...
RooAbsCategory is the common abstract base class for objects that represent a discrete value with a f...
TIterator * typeIterator() const
Return iterator over all defined states.
virtual const char * getLabel() const
Return label string of current state.
const RooCatType * lookupType(Int_t index, Bool_t printError=kFALSE) const
Find our type corresponding to the specified index, or return 0 for no match.
Int_t getSize() const
RooAbsCollection * snapshot(Bool_t deepCopy=kTRUE) const
Take a snap shot of current collection contents: An owning collection is returned containing clones o...
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
void setHashTableSize(Int_t i)
TIterator * createIterator(Bool_t dir=kIterForward) const
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
virtual Bool_t addOwned(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling addOwned() for each element in the source...
Definition: RooArgSet.h:92
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
virtual const Text_t * GetName() const
Returns name of object.
Definition: RooCatType.h:44
RooCategory represents a fundamental (non-derived) discrete value object.
Definition: RooCategory.h:24
void splitArgs(const RooArgSet &argSet, const RooAbsCategory &splitCat)
Split all arguments in 'set' into individualized clones for each defined state of 'splitCat'.
RooAbsArg * build(const char *masterCatState, Bool_t verbose=kFALSE)
Build a clone of the prototype executing all registered 'replace' rules and 'split' rules for the mas...
RooMultiCategory consolidates several RooAbsCategory objects into a single category.
RooArgSet _compSplitCatSet
RooArgSet _protoPdfSet
std::list< RooSuperCategory * > _fitCatList
RooArgSet _splitNodeList
std::list< RooSimultaneous * > _simPdfList
RooSimultaneous * buildPdf(const RooArgSet &buildConfig, const RooArgSet &dependents, const RooArgSet *auxSplitCats=0, Bool_t verbose=kFALSE)
Initialize needed components.
RooArgSet * createProtoBuildConfig()
Make RooArgSet of configuration objects.
void addSpecializations(const RooArgSet &specSet)
RooArgSet _splitNodeListOwned
RooSimPdfBuilder(const RooArgSet &pdfProtoList)
RooSimultaneous facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
Bool_t addPdf(const RooAbsPdf &pdf, const char *catLabel)
Associate given PDF with index category state label 'catLabel'.
RooStringVar implements a string values RooAbsArg.
Definition: RooStringVar.h:23
virtual const char * getVal() const
Return value of object. Calculated if dirty, otherwise cached value is returned.
Definition: RooStringVar.h:34
RooSuperCategory consolidates several RooAbsCategoryLValue objects into a single category.
virtual Bool_t setLabel(const char *label, Bool_t printError=kTRUE)
Set the value of the super category by specifying the state name by setting the state names of the co...
virtual void Print(Option_t *option="") const
Default print for collections, calls Print(option, 1).
void SetName(const char *name)
Definition: TCollection.h:204
virtual void AddAll(const TCollection *col)
Add all objects from collection col to this collection.
Iterator abstract base class.
Definition: TIterator.h:30
virtual void Reset()=0
virtual TObject * Next()=0
A doubly linked list.
Definition: TList.h:44
virtual void Add(TObject *obj)
Definition: TList.h:87
virtual TObject * FindObject(const char *name) const
Delete a TObjLink object.
Definition: TList.cxx:574
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:467
virtual TObject * First() const
Return the first object in the list. Returns 0 when list is empty.
Definition: TList.cxx:655
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Mother of all ROOT objects.
Definition: TObject.h:37
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:357
@ InputArguments
Definition: RooGlobalFunc.h:58
@ ObjectHandling
Definition: RooGlobalFunc.h:58
Definition: first.py:1
STL namespace.