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