Logo ROOT   6.10/09
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 #include <string.h>
422 
423 #ifndef _WIN32
424 #include <strings.h>
425 #else
426 
427 
428 static char *strtok_r(char *s1, const char *s2, char **lasts)
429 {
430  char *ret;
431 
432  if (s1 == NULL)
433  s1 = *lasts;
434  while(*s1 && strchr(s2, *s1))
435  ++s1;
436  if(*s1 == '\0')
437  return NULL;
438  ret = s1;
439  while(*s1 && !strchr(s2, *s1))
440  ++s1;
441  if(*s1)
442  *s1++ = '\0';
443  *lasts = s1;
444  return ret;
445 }
446 
447 #endif
448 
449 #include "Riostream.h"
450 #include "RooSimPdfBuilder.h"
451 
452 #include "RooRealVar.h"
453 #include "RooFormulaVar.h"
454 #include "RooAbsCategory.h"
455 #include "RooCategory.h"
456 #include "RooStringVar.h"
457 #include "RooMappedCategory.h"
458 #include "RooRealIntegral.h"
459 #include "RooDataSet.h"
460 #include "RooArgSet.h"
461 #include "RooPlot.h"
462 #include "RooAddPdf.h"
463 #include "RooLinearVar.h"
464 #include "RooTruthModel.h"
465 #include "RooAddModel.h"
466 #include "RooProdPdf.h"
467 #include "RooCustomizer.h"
468 #include "RooThresholdCategory.h"
469 #include "RooMultiCategory.h"
470 #include "RooSuperCategory.h"
471 #include "RooSimultaneous.h"
472 #include "RooTrace.h"
473 #include "RooFitResult.h"
474 #include "RooDataHist.h"
475 #include "RooGenericPdf.h"
476 #include "RooMsgService.h"
477 
478 using namespace std ;
479 
481 ;
482 
483 
484 
485 ////////////////////////////////////////////////////////////////////////////////
486 
488  _protoPdfSet(protoPdfSet)
489 {
493 }
494 
495 
496 
497 
498 ////////////////////////////////////////////////////////////////////////////////
499 /// Make RooArgSet of configuration objects
500 
502 {
503  RooArgSet* buildConfig = new RooArgSet ;
504  buildConfig->addOwned(* new RooStringVar("physModels","List and mapping of physics models to include in build","",4096)) ;
505  buildConfig->addOwned(* new RooStringVar("splitCats","List of categories used for splitting","",1024)) ;
506 
508  RooAbsPdf* proto ;
509  while ((proto=(RooAbsPdf*)iter->Next())) {
510  buildConfig->addOwned(* new RooStringVar(proto->GetName(),proto->GetName(),"",4096)) ;
511  }
512  delete iter ;
513 
514  return buildConfig ;
515 }
516 
517 
518 
519 ////////////////////////////////////////////////////////////////////////////////
520 
522 {
523  _splitNodeList.add(specSet) ;
524 }
525 
526 
527 
528 ////////////////////////////////////////////////////////////////////////////////
529 /// Initialize needed components
530 
531 RooSimultaneous* RooSimPdfBuilder::buildPdf(const RooArgSet& buildConfig, const RooArgSet& dependents,
532  const RooArgSet* auxSplitCats, Bool_t verbose)
533 {
534  const char* spaceChars = " \t" ;
535 
536  // Retrieve physics index category
537  Int_t buflen = strlen(((RooStringVar*)buildConfig.find("physModels"))->getVal())+1 ;
538  char *buf = new char[buflen] ;
539 
540  strlcpy(buf,((RooStringVar*)buildConfig.find("physModels"))->getVal(),buflen) ;
541  RooAbsCategoryLValue* physCat(0) ;
542  if (strstr(buf," : ")) {
543  const char* physCatName = strtok(buf,spaceChars) ;
544  physCat = dynamic_cast<RooAbsCategoryLValue*>(dependents.find(physCatName)) ;
545  if (!physCat) {
546  coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR physics index category " << physCatName
547  << " not found in dataset variables" << endl ;
548  delete[] buf ;
549  return 0 ;
550  }
551  coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: category indexing physics model: " << physCatName << endl ;
552  }
553 
554  // Create list of physics models to be built
555  char *physName ;
556  RooArgSet physModelSet ;
557  if (physCat) {
558  // Absorb colon token
559  strtok(0,spaceChars) ;
560  physName = strtok(0,spaceChars) ;
561  } else {
562  physName = strtok(buf,spaceChars) ;
563  }
564 
565  if (!physName) {
566  coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR: No models specified, nothing to do!" << endl ;
567  delete[] buf ;
568  return 0 ;
569  }
570 
571  Bool_t first(kTRUE) ;
572  RooArgSet stateMap ;
573  while(physName) {
574 
575  char *stateName(0) ;
576 
577  // physName may be <state>=<pdfName> or just <pdfName> is state and pdf have identical names
578  if (strchr(physName,'=')) {
579  // Must have a physics category for mapping to make sense
580  if (!physCat) {
581  coutW(ObjectHandling) << "RooSimPdfBuilder::buildPdf: WARNING: without physCat specification "
582  << "<physCatState>=<pdfProtoName> association is meaningless" << endl ;
583  }
584  stateName = physName ;
585  physName = strchr(stateName,'=') ;
586  if (physName) {
587  *(physName++) = 0 ;
588  }
589  } else {
590  stateName = physName ;
591  }
592 
593  RooAbsPdf* physModel = (RooAbsPdf*) (physName ? _protoPdfSet.find(physName) : 0 );
594  if (!physModel) {
595  coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR requested physics model "
596  << (physName?physName:"(null)") << " is not defined" << endl ;
597  delete[] buf ;
598  return 0 ;
599  }
600 
601  // Check if state mapping has already been defined
602  if (stateMap.find(stateName)) {
603  coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: multiple PDFs specified for state "
604  << stateName << ", only first will be used" << endl ;
605  continue ;
606  }
607 
608  // Add pdf to list of models to be processed
609  physModelSet.add(*physModel,kTRUE) ; // silence duplicate insertion warnings
610 
611  // Store state->pdf mapping
612  stateMap.addOwned(* new RooStringVar(stateName,stateName,physName)) ;
613 
614  // Continue with next mapping
615  physName = strtok(0,spaceChars) ;
616  if (first) {
617  first = kFALSE ;
618  } else if (physCat==0) {
619  coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: without physCat specification, only the first model will be used" << endl ;
620  break ;
621  }
622  }
623  coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: list of physics models " << physModelSet << endl ;
624 
625 
626 
627  // Create list of dataset categories to be used in splitting
628  TList splitStateList ;
629  RooArgSet splitCatSet ;
630 
631  delete[] buf ;
632  buflen = strlen(((RooStringVar*)buildConfig.find("splitCats"))->getVal())+1 ;
633  buf = new char[buflen] ;
634  strlcpy(buf,((RooStringVar*)buildConfig.find("splitCats"))->getVal(),buflen) ;
635 
636  char *catName = strtok(buf,spaceChars) ;
637  char *stateList(0) ;
638  while(catName) {
639 
640  // Chop off optional list of selected states
641  char* tokenPtr(0) ;
642  if (strchr(catName,'(')) {
643 
644  catName = strtok_r(catName,"(",&tokenPtr) ;
645  stateList = strtok_r(0,")",&tokenPtr) ;
646 
647  } else {
648  stateList = 0 ;
649  }
650 
651  RooCategory* splitCat = catName ? dynamic_cast<RooCategory*>(dependents.find(catName)) : 0 ;
652  if (!splitCat) {
653  coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR requested split category " << (catName?catName:"(null)")
654  << " is not a RooCategory in the dataset" << endl ;
655  delete[] buf ;
656  return 0 ;
657  }
658  splitCatSet.add(*splitCat) ;
659 
660  // Process optional state list
661  if (stateList) {
662  coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: splitting of category " << catName
663  << " restricted to states (" << stateList << ")" << endl ;
664 
665  // Create list named after this splitCat holding its selected states
666  TList* slist = new TList ;
667  slist->SetName(catName) ;
668  splitStateList.Add(slist) ;
669 
670  char* stateLabel = strtok_r(stateList,",",&tokenPtr) ;
671 
672  while(stateLabel) {
673  // Lookup state label and require it exists
674  const RooCatType* type = splitCat->lookupType(stateLabel) ;
675  if (!type) {
676  coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR splitCat " << splitCat->GetName()
677  << " doesn't have a state named " << stateLabel << endl ;
678  splitStateList.Delete() ;
679  delete[] buf ;
680  return 0 ;
681  }
682  slist->Add((TObject*)type) ;
683 
684  stateLabel = strtok_r(0,",",&tokenPtr) ;
685  }
686  }
687 
688  catName = strtok(0,spaceChars) ;
689  }
690  if (physCat) splitCatSet.add(*physCat) ;
691  RooSuperCategory masterSplitCat("masterSplitCat","Master splitting category",splitCatSet) ;
692 
693  coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: list of splitting categories " << splitCatSet << endl ;
694 
695  // Clone auxiliary split cats and attach to splitCatSet
696  RooArgSet auxSplitSet ;
697  RooArgSet* auxSplitCloneSet(0) ;
698  if (auxSplitCats) {
699  // Deep clone auxililary split cats
700  auxSplitCloneSet = (RooArgSet*) auxSplitCats->snapshot(kTRUE) ;
701  if (!auxSplitCloneSet) {
702  coutE(InputArguments) << "RooSimPdfBuilder::buildPdf(" << GetName() << ") Couldn't deep-clone set auxiliary splitcats, abort." << endl ;
703  delete[] buf ;
704  return 0 ;
705  }
706 
707  TIterator* iter = auxSplitCats->createIterator() ;
708  RooAbsArg* arg ;
709  while((arg=(RooAbsArg*)iter->Next())) {
710  // Find counterpart in cloned set
711  RooAbsArg* aux = auxSplitCats->find(arg->GetName()) ;
712 
713  // Check that there is no fundamental splitCat in the dataset with the bane of the auxiliary split
714  if (splitCatSet.find(aux->GetName())) {
715  coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: dataset contains a fundamental splitting category " << endl
716  << " with the same name as an auxiliary split function (" << aux->GetName() << "). " << endl
717  << " Auxiliary split function will be ignored" << endl ;
718  continue ;
719  }
720 
721  // Check that all servers of this aux cat are contained in splitCatSet
722  RooArgSet* parSet = aux->getParameters(splitCatSet) ;
723  if (parSet->getSize()>0) {
724  coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: ignoring auxiliary category " << aux->GetName()
725  << " because it has servers that are not listed in splitCatSet: " << *parSet << endl ;
726  delete parSet ;
727  continue ;
728  }
729 
730  // Redirect servers to splitCatSet
731  aux->recursiveRedirectServers(splitCatSet) ;
732 
733  // Add top level nodes to auxSplitSet
734  auxSplitSet.add(*aux) ;
735  }
736  delete iter ;
737 
738  coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: list of auxiliary splitting categories " << auxSplitSet << endl ;
739  }
740 
741 
742  TList* customizerList = new TList ;
743 
744  // Loop over requested physics models and build components
745  TIterator* physIter = physModelSet.createIterator() ;
746  RooAbsPdf* physModel ;
747  while((physModel=(RooAbsPdf*)physIter->Next())) {
748  coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: processing physics model " << physModel->GetName() << endl ;
749 
750  RooCustomizer* physCustomizer = new RooCustomizer(*physModel,masterSplitCat,_splitNodeList) ;
751  customizerList->Add(physCustomizer) ;
752 
753  // Parse the splitting rules for this physics model
754  RooStringVar* ruleStr = (RooStringVar*) buildConfig.find(physModel->GetName()) ;
755  if (ruleStr) {
756 
757  delete[] buf ;
758  buflen = strlen(ruleStr->getVal())+1 ;
759  buf = new char[buflen] ;
760 
761  strlcpy(buf,ruleStr->getVal(),buflen) ;
762  char *tokenPtr(0) ;
763 
764  char* token = strtok_r(buf,spaceChars,&tokenPtr) ;
765 
766  enum Mode { SplitCat, Colon, ParamList } ;
767  Mode mode(SplitCat) ;
768 
769  char* splitCatName ;
770  RooAbsCategory* splitCat(0) ;
771 
772  while(token) {
773 
774  switch (mode) {
775  case SplitCat:
776  {
777  splitCatName = token ;
778 
779  if (strchr(splitCatName,',')) {
780  // Composite splitting category
781 
782  // Check if already instantiated
783  splitCat = (RooAbsCategory*) _compSplitCatSet.find(splitCatName) ;
784  TString origCompCatName(splitCatName) ;
785  if (!splitCat) {
786  // Build now
787 
788  char *tokptr = 0;
789  char *catName2 = strtok_r(token,",",&tokptr) ;
790 
791  RooArgSet compCatSet ;
792  while(catName2) {
793  RooAbsArg* cat = splitCatSet.find(catName2) ;
794 
795  // If not, check if it is an auxiliary splitcat
796  if (!cat) {
797  cat = (RooAbsCategory*) auxSplitSet.find(catName2) ;
798  }
799 
800  if (!cat) {
801  coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR " << catName2
802  << " not found in the primary or auxilary splitcat list" << endl ;
803  customizerList->Delete() ;
804  delete customizerList ;
805 
806  splitStateList.Delete() ;
807  delete[] buf ;
808  return 0 ;
809  }
810  compCatSet.add(*cat) ;
811 
812  catName2 = strtok_r(0,",",&tokptr) ;
813  }
814 
815 
816  // check that any auxSplitCats in compCatSet do not depend on any other
817  // fundamental or auxiliary splitting categories in the composite set.
818  TIterator* iter = compCatSet.createIterator() ;
819  RooAbsArg* arg ;
820  while((arg=(RooAbsArg*)iter->Next())) {
821  RooArgSet tmp(compCatSet) ;
822  tmp.remove(*arg) ;
823  if (arg->dependsOnValue(tmp)) {
824  coutE(InputArguments) << "RooSimPdfBuilder::buildPDF: ERROR: Ill defined split: auxiliary splitting category " << arg->GetName()
825  << " used in composite split " << compCatSet << " depends on one or more of the other splitting categories in the composite split" << endl ;
826 
827  // Cleanup and axit
828  customizerList->Delete() ;
829  delete customizerList ;
830  splitStateList.Delete() ;
831  delete[] buf ;
832  return 0 ;
833  }
834  }
835  delete iter ;
836 
837  splitCat = new RooMultiCategory(origCompCatName,origCompCatName,compCatSet) ;
838  _compSplitCatSet.addOwned(*splitCat) ;
839  //cout << "composite splitcat: " << splitCat->GetName() ;
840  }
841  } else {
842  // Simple splitting category
843 
844  // First see if it is a simple splitting category
845  splitCat = (RooAbsCategory*) splitCatSet.find(splitCatName) ;
846 
847  // If not, check if it is an auxiliary splitcat
848  if (!splitCat) {
849  splitCat = (RooAbsCategory*) auxSplitSet.find(splitCatName) ;
850  }
851 
852  if (!splitCat) {
853  coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR splitting category "
854  << splitCatName << " not found in the primary or auxiliary splitcat list" << endl ;
855  customizerList->Delete() ;
856  delete customizerList ;
857  splitStateList.Delete() ;
858  delete[] buf ;
859  return 0 ;
860  }
861  }
862 
863  mode = Colon ;
864  break ;
865  }
866  case Colon:
867  {
868  if (strcmp(token,":")) {
869  coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR in parsing, expected ':' after "
870  << splitCat << ", found " << token << endl ;
871  customizerList->Delete() ;
872  delete customizerList ;
873  splitStateList.Delete() ;
874  delete[] buf ;
875  return 0 ;
876  }
877  mode = ParamList ;
878  break ;
879  }
880  case ParamList:
881  {
882  // Verify the validity of the parameter list and build the corresponding argset
883  RooArgSet splitParamList ;
884  RooArgSet* paramList = physModel->getParameters(dependents) ;
885 
886  // wve -- add nodes to parameter list
887  RooArgSet* compList = physModel->getComponents() ;
888  paramList->add(*compList) ;
889  delete compList ;
890 
891  Bool_t lastCharIsComma = (token[strlen(token)-1]==',') ;
892 
893  char *tokptr = 0 ;
894  char *paramName = strtok_r(token,",",&tokptr) ;
895 
896  // Check for fractional split option 'param_name[remainder_state]'
897  char *remainderState = 0 ;
898  char *tokptr2 = 0 ;
899  if (paramName && strtok_r(paramName,"[",&tokptr2)) {
900  remainderState = strtok_r(0,"]",&tokptr2) ;
901  }
902 
903  while(paramName) {
904 
905  // If fractional split is specified, check that remainder state is a valid state of this split cat
906  if (remainderState) {
907  if (!splitCat->lookupType(remainderState)) {
908  coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR fraction split of parameter "
909  << paramName << " has invalid remainder state name: " << remainderState << endl ;
910  delete paramList ;
911  customizerList->Delete() ;
912  delete customizerList ;
913  splitStateList.Delete() ;
914  delete[] buf ;
915  return 0 ;
916  }
917  }
918 
919  RooAbsArg* param = paramList->find(paramName) ;
920  if (!param) {
921  coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR " << paramName
922  << " is not a parameter of physics model " << physModel->GetName() << endl ;
923  delete paramList ;
924  customizerList->Delete() ;
925  delete customizerList ;
926  splitStateList.Delete() ;
927  delete[] buf ;
928  return 0 ;
929  }
930  splitParamList.add(*param) ;
931 
932  // Build split leaf of fraction splits here
933  if (remainderState) {
934 
935  // Check if we are splitting a real-valued parameter
936  if (!dynamic_cast<RooAbsReal*>(param)) {
937  coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR fraction split requested of non-real valued parameter "
938  << param->GetName() << endl ;
939  delete paramList ;
940  customizerList->Delete() ;
941  delete customizerList ;
942  splitStateList.Delete() ;
943  delete[] buf ;
944  return 0 ;
945  }
946 
947  // Check if we are doing a restricted build
948  TList* remStateSplitList = static_cast<TList*>(splitStateList.FindObject(splitCat->GetName())) ;
949 
950  // If so, check if remainder state is actually being built.
951  if (remStateSplitList && !remStateSplitList->FindObject(remainderState)) {
952  coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR " << paramName
953  << " remainder state " << remainderState << " in parameter split "
954  << param->GetName() << " is not actually being built" << endl ;
955  delete paramList ;
956  customizerList->Delete() ;
957  delete customizerList ;
958  splitStateList.Delete() ;
959  delete[] buf ;
960  return 0 ;
961  }
962 
963  TIterator* iter = splitCat->typeIterator() ;
964  RooCatType* type ;
965  RooArgList fracLeafList ;
966  TString formExpr("1") ;
967  Int_t i(0) ;
968 
969  while((type=(RooCatType*)iter->Next())) {
970 
971  // Skip remainder state
972  if (!TString(type->GetName()).CompareTo(remainderState)) continue ;
973 
974  // If restricted build is requested, skip states of splitcat that are not built
975  if (remStateSplitList && !remStateSplitList->FindObject(type->GetName())) {
976  continue ;
977  }
978 
979  // Construct name of split leaf
980  TString splitLeafName(param->GetName()) ;
981  splitLeafName.Append("_") ;
982  splitLeafName.Append(type->GetName()) ;
983 
984  // Check if split leaf already exists
985  RooAbsArg* splitLeaf = _splitNodeList.find(splitLeafName) ;
986  if (!splitLeaf) {
987  // If not create it now
988  splitLeaf = (RooAbsArg*) param->clone(splitLeafName) ;
989  _splitNodeList.add(*splitLeaf) ;
990  _splitNodeListOwned.addOwned(*splitLeaf) ;
991  }
992  fracLeafList.add(*splitLeaf) ;
993  formExpr.Append(Form("-@%d",i++)) ;
994  }
995  delete iter ;
996 
997  // Construct RooFormulaVar expresssing remainder of fraction
998  TString remLeafName(param->GetName()) ;
999  remLeafName.Append("_") ;
1000  remLeafName.Append(remainderState) ;
1001 
1002  // Check if no specialization was already specified for remainder state
1003  if (!_splitNodeList.find(remLeafName)) {
1004  RooAbsArg* remLeaf = new RooFormulaVar(remLeafName,formExpr,fracLeafList) ;
1005  _splitNodeList.add(*remLeaf) ;
1006  _splitNodeListOwned.addOwned(*remLeaf) ;
1007  coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: creating remainder fraction formula for " << remainderState
1008  << " specialization of split parameter " << param->GetName() << " " << formExpr << endl ;
1009  }
1010  }
1011 
1012  // Parse next parameter name
1013  paramName = strtok_r(0,",",&tokptr) ;
1014  if (paramName && strtok_r(paramName,"[",&tokptr2)) {
1015  remainderState = strtok_r(0,"]",&tokptr2) ;
1016  }
1017  }
1018 
1019  // Add the rule to the appropriate customizer ;
1020  physCustomizer->splitArgs(splitParamList,*splitCat) ;
1021 
1022  delete paramList ;
1023 
1024  if (!lastCharIsComma) mode = SplitCat ;
1025  break ;
1026  }
1027  }
1028 
1029  token = strtok_r(0,spaceChars,&tokenPtr) ;
1030 
1031  }
1032  if (mode!=SplitCat) {
1033  coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR in parsing, expected "
1034  << (mode==Colon?":":"parameter list") << " after " << (token?token:"(null)") << endl ;
1035  }
1036 
1037  //RooArgSet* paramSet = physModel->getParameters(dependents) ;
1038  } else {
1039  coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: no splitting rules for " << physModel->GetName() << endl ;
1040  }
1041  }
1042 
1043  coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: configured customizers for all physics models" << endl ;
1044  if (oodologI((TObject*)0,ObjectHandling)) {
1045  customizerList->Print() ;
1046  }
1047 
1048  // Create fit category from physCat and splitCatList ;
1049  RooArgSet fitCatList ;
1050  if (physCat) fitCatList.add(*physCat) ;
1051  fitCatList.add(splitCatSet) ;
1052  TIterator* fclIter = fitCatList.createIterator() ;
1053  RooSuperCategory *fitCat = new RooSuperCategory("fitCat","fitCat",fitCatList) ;
1054 
1055  // Create master PDF
1056  RooSimultaneous* simPdf = new RooSimultaneous("simPdf","simPdf",*fitCat) ;
1057 
1058  // Add component PDFs to master PDF
1059  TIterator* fcIter = fitCat->typeIterator() ;
1060 
1061  RooCatType* fcState ;
1062  while((fcState=(RooCatType*)fcIter->Next())) {
1063  // Select fitCat state
1064  fitCat->setLabel(fcState->GetName()) ;
1065 
1066  // Check if this fitCat state is selected
1067  fclIter->Reset() ;
1068  RooAbsCategory* splitCat ;
1069  Bool_t select(kTRUE) ;
1070  while((splitCat=(RooAbsCategory*)fclIter->Next())) {
1071  // Find selected state list
1072  TList* slist = (TList*) splitStateList.FindObject(splitCat->GetName()) ;
1073  if (!slist) continue ;
1074  RooCatType* type = (RooCatType*) slist->FindObject(splitCat->getLabel()) ;
1075  if (!type) {
1076  select = kFALSE ;
1077  }
1078  }
1079  if (!select) continue ;
1080 
1081 
1082  // Select appropriate PDF for this physCat state
1083  RooCustomizer* physCustomizer ;
1084  if (physCat) {
1085  RooStringVar* physNameVar = (RooStringVar*) stateMap.find(physCat->getLabel()) ;
1086  if (!physNameVar) continue ;
1087  physCustomizer = (RooCustomizer*) customizerList->FindObject(physNameVar->getVal());
1088  } else {
1089  physCustomizer = (RooCustomizer*) customizerList->First() ;
1090  }
1091 
1092  coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: Customizing physics model " << physCustomizer->GetName()
1093  << " for mode " << fcState->GetName() << endl ;
1094 
1095  // Customizer PDF for current state and add to master simPdf
1096  RooAbsPdf* fcPdf = (RooAbsPdf*) physCustomizer->build(masterSplitCat.getLabel(),verbose) ;
1097  simPdf->addPdf(*fcPdf,fcState->GetName()) ;
1098  }
1099  delete fcIter ;
1100 
1101  // Move customizers (owning the cloned branch node components) to the attic
1102  _retiredCustomizerList.AddAll(customizerList) ;
1103  delete customizerList ;
1104 
1105  delete fclIter ;
1106  splitStateList.Delete() ;
1107 
1108  if (auxSplitCloneSet) delete auxSplitCloneSet ;
1109  delete physIter ;
1110 
1111  delete[] buf ;
1112  _simPdfList.push_back(simPdf) ;
1113  _fitCatList.push_back(fitCat) ;
1114  return simPdf ;
1115 }
1116 
1117 
1118 
1119 
1120 
1121 ////////////////////////////////////////////////////////////////////////////////
1122 
1124 {
1126 
1127  std::list<RooSimultaneous*>::iterator iter = _simPdfList.begin() ;
1128  while(iter != _simPdfList.end()) {
1129  delete *iter ;
1130  ++iter ;
1131  }
1132 
1133  std::list<RooSuperCategory*>::iterator iter2 = _fitCatList.begin() ;
1134  while(iter2 != _fitCatList.end()) {
1135  delete *iter2 ;
1136  ++iter2 ;
1137  }
1138 
1139 }
1140 
1141 
virtual TObject * clone(const char *newname=0) const =0
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
TIterator * createIterator(Bool_t dir=kIterForward) const
#define coutE(a)
Definition: RooMsgService.h:34
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
virtual const char * getVal() const
Return value of object. Calculated if dirty, otherwise cached value is returned.
Definition: RooStringVar.h:34
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:409
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:86
virtual void Reset()=0
std::list< RooSimultaneous * > _simPdfList
#define coutI(a)
Definition: RooMsgService.h:31
virtual void AddAll(const TCollection *col)
Add all objects from collection col to this collection.
Definition: TCollection.cxx:58
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
STL namespace.
#define coutW(a)
Definition: RooMsgService.h:33
#define NULL
Definition: RtypesCore.h:88
Iterator abstract base class.
Definition: TIterator.h:30
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:501
Bool_t addPdf(const RooAbsPdf &pdf, const char *catLabel)
Associate given PDF with index category state label &#39;catLabel&#39;.
RooCatType is an auxilary class for RooAbsCategory and defines a a single category state...
Definition: RooCatType.h:22
RooAbsCategoryLValue is the common abstract base class for objects that represent a discrete value th...
RooArgSet _splitNodeListOwned
A doubly linked list.
Definition: TList.h:43
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 TObject * First() const
Return the first object in the list. Returns 0 when list is empty.
Definition: TList.cxx:561
RooSimultaneous * buildPdf(const RooArgSet &buildConfig, const RooArgSet &dependents, const RooArgSet *auxSplitCats=0, Bool_t verbose=kFALSE)
Initialize needed components.
TIterator * typeIterator() const
Return iterator over all defined states.
virtual const Text_t * GetName() const
Returns name of object.
Definition: RooCatType.h:44
RooArgSet * getComponents() const
Definition: RooAbsArg.cxx:688
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.
bool verbose
char * Form(const char *fmt,...)
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:90
virtual const char * getLabel() const
Return label string of current state.
void SetName(const char *name)
Definition: TCollection.h:111
RooArgSet _splitNodeList
RooCategory represents a fundamental (non-derived) discrete value object.
Definition: RooCategory.h:24
std::list< RooSuperCategory * > _fitCatList
const Bool_t kFALSE
Definition: RtypesCore.h:92
RooMultiCategory consolidates several RooAbsCategory objects into a single category.
#define ClassImp(name)
Definition: Rtypes.h:336
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooArgSet _compSplitCatSet
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&#39;t match any of...
Definition: RooAbsArg.cxx:560
#define oodologI(o, a)
Definition: RooMsgService.h:71
int type
Definition: TGX11.cxx:120
void addSpecializations(const RooArgSet &specSet)
Bool_t dependsOnValue(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0) const
Definition: RooAbsArg.h:88
Mother of all ROOT objects.
Definition: TObject.h:37
void splitArgs(const RooArgSet &argSet, const RooAbsCategory &splitCat)
Split all arguments in &#39;set&#39; into individualized clones for each defined state of &#39;splitCat&#39;...
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
RooArgSet _protoPdfSet
RooAbsArg * build(const char *masterCatState, Bool_t verbose=kFALSE)
Build a clone of the prototype executing all registered &#39;replace&#39; rules and &#39;split&#39; rules for the mas...
virtual void Add(TObject *obj)
Definition: TList.h:77
virtual TObject * Next()=0
RooAbsCategory is the common abstract base class for objects that represent a discrete value with a f...
const char * proto
Definition: civetweb.c:11652
RooSimPdfBuilder(const RooArgSet &pdfProtoList)
RooSuperCategory consolidates several RooAbsCategoryLValue objects into a single category.
Definition: first.py:1
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
void setHashTableSize(Int_t i)
virtual void Print(Option_t *option="") const
Default print for collections, calls Print(option, 1).
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:364
Bool_t recursiveRedirectServers(const RooAbsCollection &newServerList, Bool_t mustReplaceAll=kFALSE, Bool_t nameChange=kFALSE, Bool_t recurseInNewSet=kTRUE)
Definition: RooAbsArg.cxx:1088
RooStringVar implements a string values RooAbsArg.
Definition: RooStringVar.h:23
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...
RooSimultaneous facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset...
const Bool_t kTRUE
Definition: RtypesCore.h:91
RooArgSet * createProtoBuildConfig()
Make RooArgSet of configuration objects.