RooSimPdfBuilder.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitCore                                                       *
00004  * @(#)root/roofitcore:$Id: RooSimPdfBuilder.cxx 36274 2010-10-11 08:05:03Z wouter $
00005  * Authors:                                                                  *
00006  *   WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu       *
00007  *   DK, David Kirkby,    UC Irvine,         dkirkby@uci.edu                 *
00008  *                                                                           *
00009  * Copyright (c) 2000-2005, Regents of the University of California          *
00010  *                          and Stanford University. All rights reserved.    *
00011  *                                                                           *
00012  * Redistribution and use in source and binary forms,                        *
00013  * with or without modification, are permitted according to the terms        *
00014  * listed in LICENSE (http://roofit.sourceforge.net/license.txt)             *
00015  *****************************************************************************/
00016 
00017 //////////////////////////////////////////////////////////////////////////////
00018 //
00019 // Begin_Html
00020 //
00021 // <b>This tool has now been superceded by RooSimWSTool</b>
00022 // 
00023 //  <p>
00024 //    <tt>RooSimPdfBuilder</tt> is a powerful tool to build <tt>RooSimultaneous</tt>
00025 //    PDFs that are defined in terms component PDFs that are identical in
00026 //    structure, but have different parameters. 
00027 //  </p>
00028 //
00029 //  <h2>Example</h2>
00030 //
00031 //  <p>
00032 //    The following example demonstrates the essence of <tt>RooSimPdfBuilder</tt>:
00033 //    Given a dataset D with a <tt>RooRealVar X</tt> and a <tt>RooCategory C</tt> that has
00034 //    state C1 and C2. 
00035 //    <ul>
00036 //    <li> We want to fit the distribution of <tt>X</tt> with a Gaussian+ArgusBG PDF, 
00037 //    <li> We want to fit the data subsets <tt>D(C==C1)</tt> and <tt>D(C==C2)</tt> separately and simultaneously. 
00038 //    <li> The PDFs to fit data subsets D_C1 and D_C2 are identical except for 
00039 //         <ul>
00040 //         <li> the kappa parameter of the ArgusBG PDF and 
00041 //         <li> the sigma parameter of the gaussian PDF
00042 //         </ul>
00043 //         where each PDF will have its own copy of the parameter
00044 //    </ul>
00045 //  </p>
00046 //  <p>
00047 //    Coding this example directly with RooFit classes gives
00048 //    (we assume dataset D and variables C and X have been declared previously)
00049 //  </p>
00050 //  <pre>
00051 // RooRealVar m("m","mean of gaussian",-10,10) ;
00052 // RooRealVar s_C1("s_C1","sigma of gaussian C1",0,20) ;
00053 // RooRealVar s_C2("s_C2","sigma of gaussian C2",0,20) ;
00054 // RooGaussian gauss_C1("gauss_C1","gaussian C1",X,m,s_C1) ;
00055 // RooGaussian gauss_C2("gauss_C2","gaussian C2",X,m,s_C2) ;
00056 //
00057 // RooRealVar k_C1("k_C1","ArgusBG kappa parameter C1",-50,0) ;
00058 // RooRealVar k_C2("k_C2","ArgusBG kappa parameter C2",-50,0) ;
00059 // RooRealVar xm("xm","ArgusBG cutoff point",5.29) ;
00060 // RooArgusBG argus_C1("argus_C1","argus background C1",X,k_C1,xm) ;
00061 // RooArgusBG argus_C2("argus_C2","argus background C2",X,k_C2,xm) ;
00062 //
00063 // RooRealVar gfrac("gfrac","fraction of gaussian",0.,1.) ;
00064 // RooAddPdf pdf_C1("pdf_C1","gauss+argus_C1",RooArgList(gauss_C1,argus_C1),gfrac) ;
00065 // RooAddPdf pdf_C2("pdf_C2","gauss+argus_C2",RooArgList(gauss_C2,argus_C2),gfrac) ;
00066 //
00067 // RooSimultaneous simPdf("simPdf","simPdf",C) ;   
00068 // simPdf.addPdf(pdf_C1,"C1") ;
00069 // simPdf.addPdf(pdf_C2,"C2") ;
00070 //  </pre>
00071 //  <p>
00072 //    Coding this example with RooSimPdfBuilder gives
00073 //  </p>
00074 //  <pre>
00075 // RooRealVar m("m","mean of gaussian",-10,10) ;
00076 // RooRealVar s("s","sigma of gaussian",0,20) ;
00077 // RooGaussian gauss("gauss","gaussian",X,m,s) ;
00078 //
00079 // RooRealVar k("k","ArgusBG kappa parameter",-50,0) ;
00080 // RooRealVar xm("xm","ArgusBG cutoff point",5.29) ;
00081 // RooArgusBG argus("argus","argus background",X,k,xm) ;
00082 //
00083 // RooRealVar gfrac("gfrac","fraction of gaussian",0.,1.) ;
00084 // RooAddPdf pdf("pdf","gauss+argus",RooArgList(gauss,argus),gfrac) ;
00085 //
00086 // RooSimPdfBuilder builder(pdf) ;
00087 // RooArgSet* config = builder.createProtoBuildConfig() ;
00088 // (*config)["physModels"] = "pdf" ;      // Name of the PDF we are going to work with
00089 // (*config)["splitCats"]  = "C" ;        // Category used to differentiate sub-datasets
00090 // (*config)["pdf"]        = "C : k,s" ;  // Prescription to taylor PDF parameters k and s 
00091 //                                        // for each data subset designated by C states
00092 // RooSimultaneous* simPdf = builder.buildPdf(*config,&D) ;
00093 //  </pre>
00094 //  <p>
00095 //    The above snippet of code demonstrates the concept of <tt>RooSimPdfBuilder</tt>:
00096 //    the user defines a single <i>'prototype' PDF</i> that defines the structure of all
00097 //    PDF components of the <tt>RooSimultaneous</tt> PDF to be built. <tt>RooSimPdfBuilder</tt> 
00098 //    then takes this prototype and replicates it as a component 
00099 //    PDF for each state of the C index category.
00100 //  </p>
00101 //  <p>
00102 //    In the above example </tt>RooSimPdfBuilder</tt>
00103 //    will first replicate <tt>k</tt> and <tt>s</tt> into 
00104 //    <tt>k_C1,k_C2</tt> and <tt>s_C1,s_C2</tt>, as prescribed in the
00105 //    configuration. Then it will recursively replicate all PDF nodes that depend on
00106 //    the 'split' parameter nodes: <tt>gauss</tt> into <tt>gauss_C1,C2</tt>, <tt>argus</tt> 
00107 //    into <tt>argus_C1,C2</tt> and finally <tt>pdf</tt> into <tt>pdf_C1,pdf_C2</tt>. 
00108 //     When PDFs for all states of C have been replicated
00109 //    they are assembled into a <tt>RooSimultaneous</tt> PDF, which is returned by the <tt>buildPdf()</tt>
00110 //    method.
00111 //  </p>
00112 //  <p>
00113 //    Although in this very simple example the use of <tt>RooSimPdfBuilder</tt> doesn't
00114 //    reduce the amount of code much, it is already easier to read and maintain
00115 //    because there is no duplicate code. As the complexity of the <tt>RooSimultaneous</tt>
00116 //    to be built increases, the advantages of <tt>RooSimPdfBuilder</tt> will become more and
00117 //    more apparent.
00118 //  </p>
00119 //
00120 //
00121 //  <h2>Builder configuration rules for a single prototype PDF</h2>
00122 //  <p>
00123 //    Each builder configuration needs at minumum two lines, <tt>physModels</tt> and <tt>splitCats</tt>, which identify
00124 //    the ingredients of the build. In this section we only explain the building rules for
00125 //    builds from a single prototype PDF. In that case the <tt>physModels</tt> line always reads
00126 //  </p>
00127 //  <pre>
00128 //  physModels = {pdfName}
00129 //  </pre>
00130 //  <p>
00131 //    The second line, <tt>splitCats</tt>, indicates which categories are going to be used to 
00132 //    differentiate the various subsets of the 'master' input data set. You can enter
00133 //    a single category here, or multiple if necessary:
00134 //  </p>
00135 //  <pre>
00136 // splitCats = {catName} [{catName} ...]
00137 //  </pre>
00138 //  <p>
00139 //    All listed splitcats must be <tt>RooCategories</tt> that appear in the dataset provided to
00140 //    <tt>RooSimPdfBuilder::buildPdf()</tt>
00141 //  </p>
00142 //  <p>
00143 //    The parameter splitting prescriptions, the essence of each build configuration
00144 //    can be supplied in a third line carrying the name of the pdf listed in <tt>physModels</tt>
00145 //  </p>
00146 //  <pre>
00147 // pdfName = {splitCat} : {parameter} [,{parameter},....]
00148 //  </pre>
00149 //  <p>     
00150 //    Each pdf can have only one line with splitting rules, but multiple rules can be
00151 //    supplied in each line, e.g.
00152 //  </p>
00153 //  <pre>
00154 // pdfName = {splitCat} : {parameter} [,{parameter},....] 
00155 //           {splitCat} : {parameter} [,{parameter},....]
00156 //  </pre>
00157 //  <p>
00158 //    Conversely, each parameter can only have one splitting prescription, but it may be split
00159 //    by multiple categories, e.g.
00160 //  </p>
00161 //  <pre>
00162 // pdfName = {splitCat1},{splitCat2} : {parameter}
00163 //  </pre>
00164 //  <p>
00165 //    instructs <tt>RooSimPdfBuilder</tt> to build a <tt>RooSuperCategory</tt> 
00166 //    of <tt>{splitCat1}</tt> and <tt>{splitCat2}</tt>
00167 //    and split <tt>{parameter}</tt> with that <tt>RooSuperCategory</tt>
00168 //  </p>
00169 //  <p>
00170 //    Here is an example of a builder configuration that uses several of the options discussed
00171 //    above:
00172 //  </p>
00173 //  <pre>
00174 //   physModels = pdf
00175 //   splitCats  = tagCat runBlock
00176 //   pdf        = tagCat          : signalRes,bkgRes 
00177 //                runBlock        : fudgeFactor      
00178 //                tagCat,runBlock : kludgeParam
00179 //  </pre>
00180 //
00181 //  <h2>How to enter configuration data</h2>
00182 //
00183 //  <p>
00184 //    The prototype builder configuration returned by 
00185 //    <tt>RooSimPdfBuilder::createProtoBuildConfig()</tt> is a pointer to a <tt>RooArgSet</tt> filled with
00186 //    initially blank <tt>RooStringVars</tt> named <tt>physModels,splitCats</tt> and one additional for each
00187 //    PDF supplied to the <tt>RooSimPdfBuilders</tt> constructor (with the same name)
00188 //  </p>
00189 //  <p>
00190 //    In macro code, the easiest way to assign new values to these <tt>RooStringVars</tt>
00191 //    is to use <tt>RooArgSet</tt>s array operator and the <tt>RooStringVar</tt>s assignment operator, e.g.
00192 //  </p>
00193 //  <pre>
00194 // (*config)["physModels"] = "Blah" ;
00195 //  </pre>
00196 //  <p>
00197 //    To enter multiple splitting rules simply separate consecutive rules by whitespace
00198 //    (not newlines), e.g.           
00199 //  </p>
00200 //  <pre>
00201 // (*config)["physModels"] = "Blah " // << note trailing space here
00202 //                          "Blah 2" ;
00203 //  </pre>
00204 //  <p>
00205 //    In this example, the C++ compiler will concatenate the two string literals (without inserting
00206 //    any whitespace), so the extra space after 'Blah' is important here.
00207 //  </p>
00208 //  <p>      
00209 //    Alternatively, you can read the configuration from an ASCII file, as you can
00210 //    for any <tt>RooArgSet</tt> using <tt>RooArgSet::readFromFile()</tt>. In that case the ASCII file
00211 //    can follow the syntax of the examples above and the '<tt>\</tt>' line continuation 
00212 //    sequence can be used to fold a long splitting rule over multiple lines.
00213 //  </p>
00214 //  <pre>
00215 // RooArgSet* config = builder.createProtoBuildConfig() ;
00216 // config->readFromFile("config.txt") ;
00217 //
00218 // --- config.txt ----------------    
00219 // physModels = pdf
00220 // splitCats  = tagCat
00221 // pdf        = tagCat : bogusPar
00222 // -------------------------------
00223 //  </pre>
00224 //
00225 //
00226 //  <h2>Working with multiple prototype PDFs</h2>
00227 //  <p>
00228 //    It is also possible to build a <tt>RooSimultaneous</tt> PDF from multiple PDF prototypes.
00229 //    This is appropriate for cases where the input prototype PDF would otherwise be 
00230 //    a <tt>RooSimultaneous</tt> PDF by itself. In such cases we don't feed a single
00231 //    <tt>RooSimultaneous</tt> PDF into <tt>RooSimPdfBuilder</tt>, instead we feed it its ingredients and
00232 //    add a prescription to the builder configuration that corresponds to the 
00233 //    PDF-category state mapping of the prototype <tt>RooSimultaneous</tt>.
00234 //  </p>
00235 //  <p>
00236 //    The constructor of the <tt>RooSimPdfBuilder</tt> will look as follows:
00237 //  </p>
00238 //  <pre>
00239 //  RooSimPdfBuilder builder(RooArgSet(pdfA,pdfB,...)) ;
00240 //  </pre>
00241 //  <p>
00242 //    The <tt>physModels</tt> line is now expanded to carry the pdf->state mapping information
00243 //    that the prototype <tt>RooSimultaneous</tt> would have. I.e.
00244 //  </p>
00245 //  <pre>
00246 // physModels = mode : pdfA=modeA  pdfB=modeB
00247 //  </pre>
00248 //  <p>
00249 //    is equivalent to a prototype <tt>RooSimultaneous</tt> constructed as
00250 //  </p>
00251 //  <pre>
00252 // RooSimultanous simPdf("simPdf","simPdf",mode);
00253 // simPdf.addPdf(pdfA,"modeA") ;
00254 // simPdf.addPdf(pdfB,"modeB") ;
00255 //  </pre>
00256 //  <p>
00257 //    The rest of the builder configuration works the same, except that
00258 //    each prototype PDF now has its own set of splitting rules, e.g.
00259 //  </p>
00260 //  <pre>
00261 // physModels = mode : pdfA=modeA  pdfB=modeB
00262 // splitCats  = tagCat
00263 // pdfA       = tagCat : bogusPar
00264 // pdfB       = tagCat : fudgeFactor   
00265 //  </pre>
00266 //  <p>
00267 //    Please note that 
00268 //    <ul>
00269 //    <li> The master index category ('mode' above) doesn't have to be listed in 
00270 //         <tt>splitCats</tt>, this is implicit.
00271 //
00272 //    <li> The number of splitting prescriptions goes by the
00273 //          number of prototype PDFs and not by the number of states of the
00274 //          master index category (mode in the above and below example). 
00275 //  </ul>
00276 //
00277 //  In the following case:
00278 //</p>
00279 //  <pre>
00280 //    physModels = mode : pdfA=modeA  pdfB=modeB  pdfA=modeC  pdfB=modeD
00281 //  </pre>
00282 //  <p>
00283 //    there are still only 2 sets of splitting rules: one for <tt>pdfA</tt> and one
00284 //    for <tt>pdfB</tt>. However, you <i>can</i> differentiate between <tt>modeA</tt> and <tt>modeC</tt> in 
00285 //    the above example. The technique is to use <tt>mode</tt> as splitting category, e.g.
00286 //  </p>
00287 //  <pre>
00288 //    physModels = mode : pdfA=modeA  pdfB=modeB  pdfA=modeC  pdfB=modeD
00289 //    splitCats = tagCat
00290 //    pdfA      = tagCat : bogusPar 
00291 //                mode   : funnyPar
00292 //    pdfB      = mode   : kludgeFactor
00293 //  </pre>
00294 //  <p>
00295 //    will result in an individual set of <tt>funnyPar</tt> parameters for <tt>modeA</tt> and <tt>modeC</tt>
00296 //    labeled <tt>funnyPar_modeA</tt> and <tt>funnyPar_modeB</tt> and an individual set of
00297 //    kludgeFactor parameters for <tt>pdfB</tt>, <tt>kludgeFactor_modeB</tt> and <tt>kludgeFactor_modeD</tt>. 
00298 //    Please note that for splits in the master index category (mode) only the 
00299 //    applicable states are built (A,C for <tt>pdfA</tt>, B,D for <tt>pdfB</tt>)
00300 //  </p>
00301 //
00302 //
00303 //  <h2>Advanced options</h2>
00304 //
00305 //  <h4>Partial splits</h4>
00306 //  <p>
00307 //    You can request to limit the list of states of each splitCat that
00308 //    will be considered in the build. This limitation is requested in the 
00309 //    each build as follows:
00310 //  </p>
00311 //  <pre>
00312 // splitCats = tagCat(Lep,Kao) RunBlock(Run1)
00313 //  </pre>
00314 //  <p>
00315 //    In this example the splitting of <tt>tagCat</tt> is limited to states <tt>Lep,Kao</tt>
00316 //    and the splitting of <tt>runBlock</tt> is limited to <tt>Run1</tt>. The splits apply
00317 //    globally to each build, i.e. every parameter split requested in this
00318 //    build will be limited according to these specifications. 
00319 //  </p>
00320 //  <p>      
00321 //    NB: Partial builds have no pdf associated with the unbuilt states of the 
00322 //    limited splits. Running such a pdf on a dataset that contains data with 
00323 //    unbuilt states will result in this data being ignored completely.
00324 //  </p>
00325 //
00326 //
00327 //  <h4>Non-trivial splits</h4>
00328 //  <p>
00329 //    It is possible to make non-trivial parameter splits with <tt>RooSimPdfBuilder</tt>.
00330 //    Trivial splits are considered simple splits in one (fundamental) category
00331 //    in the dataset or a split in a <tt>RooSuperCategory</tt> 'product' of multiple
00332 //    fundamental categories in the dataset. Non-trivial splits can be performed
00333 //    using an intermediate 'category function' (<tt>RooMappedCategory,
00334 //    RooGenericCategory,RooThresholdCategory</tt> etc), i.e. any <tt>RooAbsCategory</tt>
00335 //    derived objects that calculates its output as function of one or more
00336 //    input <tt>RooRealVars</tt> and/or <tt>RooCategories</tt>.
00337 //  </p>
00338 //  <p>
00339 //    Such 'function categories' objects must be constructed by the user prior
00340 //    to building the PDF. In the <tt>RooSimPdfBuilder::buildPdf()</tt> function these
00341 //    objects can be passed in an optional <tt>RooArgSet</tt> called 'auxiliary categories':
00342 //  </p>
00343 //  <pre>
00344 //   const <tt>RooSimultaneous</tt>* buildPdf(const RooArgSet& buildConfig, const RooAbsData* dataSet, 
00345 //                                   const RooArgSet& auxSplitCats, Bool_t verbose=kFALSE) {
00346 //                                   ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
00347 //  </pre>
00348 //  <p>
00349 //    Objects passed in this argset can subsequently be used in the build configuration, e.g.
00350 //  </p>
00351 //  <pre>
00352 // RooMappedCategory tagMap("tagMap","Mapped tagging category",tagCat,"CutBased") ;
00353 // tagMap.map("Lep","CutBased") ;
00354 // tagMap.map("Kao","CutBased") ;
00355 // tagMap.map("NT*","NeuralNet") ;                                                                                          
00356 // ...
00357 // builder.buildPdf(config,D,tagMap) ;
00358 //                          ^^^^^^
00359 //<Contents of config>
00360 //   physModels = pdf
00361 //   splitCats  = tagCat runBlock
00362 //   pdf        = tagCat          : signalRes 
00363 //                tagMap          : fudgeFactor      
00364 //                ^^^^^^
00365 //  </pre>
00366 //  <p>
00367 //    In the above example <tt>signalRes</tt> will be split in <tt>signalRes_Kao,signalRes_Lep,
00368 //    signalRes_NT1,signalRes_NT2</tt>, while <tt>fudgeFactor</tt> will be split in <tt>fudgeFactor_CutBased</tt>
00369 //    and <tt>fudgeFactor_NeuralNet</tt>.
00370 //  </p>
00371 //  <p>
00372 //    Category functions passed in the auxSplitCats <tt>RooArgSet</tt> can be used regularly
00373 //    in the splitting configuration. They should not be listed in <tt>splitCats</tt>,
00374 //    but must be able to be expressed <i>completely</i> in terms of the <tt>splitCats</tt> that 
00375 //    are listed.
00376 //  </p>
00377 //
00378 //  
00379 //  <h4>Multiple connected builds</h4>
00380 //  <p>
00381 //    Sometimes you want to build multiple PDFs for independent consecutive fits 
00382 //    that share some of their parameters. For example, we have two prototype PDFs 
00383 //    <tt>pdfA(x;p,q)</tt> and <tt>pdfB(x;p,r)</tt> that have a common parameter <tt>p</tt>. 
00384 //    We want to build a <tt>RooSimultaneous</tt> for both <tt>pdfA</tt> and <tt>B</tt>, 
00385 //    which involves a split of parameter <tt>p</tt> and we would like to build the
00386 //    simultaneous pdfs </tt>simA</tt> and <tt>simB</tt> such that still share their (now split) parameters
00387 //    <tt>p_XXX</tt>. This is accomplished by letting a single instance of <tt>RooSimPdfBuilder</tt> handle
00388 //    the builds of both <tt>pdfA</tt> and <tt>pdfB</tt>, as illustrated in this example:
00389 //  </p>
00390 //  <pre>
00391 // RooSimPdfBuilder builder(RooArgSet(pdfA,pdfB)) ;
00392 //
00393 // RooArgSet* configA = builder.createProtoBuildConfig() ;
00394 // (*configA)["physModels"] = "pdfA" ;     
00395 // (*configA)["splitCats"]  = "C" ;        
00396 // (*configA)["pdf"]        = "C : p" ;  
00397 // RooSimultaneous* simA = builder.buildPdf(*configA,&D) ;
00398 //
00399 // RooArgSet* configB = builder.createProtoBuildConfig() ;
00400 // (*configA)["physModels"] = "pdfB" ;     
00401 // (*configA)["splitCats"]  = "C" ;        
00402 // (*configA)["pdf"]        = "C : p" ;  
00403 // RooSimultaneous* simB = builder.buildPdf(*configB,&D) ;
00404 //  </pre>
00405 //
00406 //  <h2>Ownership of constructed PDFs</h2>
00407 //  <p>
00408 //    The <tt>RooSimPdfBuilder</tt> instance owns all the objects it creates, including the top-level
00409 //    <tt>RooSimultaneous</tt> returned by <tt>buildPdf()</tt>. Therefore the builder instance should 
00410 //    exist as long as the constructed PDFs needs to exist.
00411 //  </p>
00412 //
00413 // End_Html
00414 //
00415 
00416 
00417 
00418 #include "RooFit.h"
00419 
00420 #include <string.h>
00421 #include <string.h>
00422 
00423 #ifndef _WIN32
00424 #include <strings.h>
00425 #else
00426 
00427 
00428 static char *strtok_r(char *s1, const char *s2, char **lasts)
00429 {
00430   char *ret;
00431   
00432   if (s1 == NULL)
00433     s1 = *lasts;
00434   while(*s1 && strchr(s2, *s1))
00435     ++s1;
00436   if(*s1 == '\0')
00437     return NULL;
00438   ret = s1;
00439   while(*s1 && !strchr(s2, *s1))
00440     ++s1;
00441   if(*s1)
00442     *s1++ = '\0';
00443   *lasts = s1;
00444   return ret;
00445 }
00446 
00447 #endif
00448 
00449 #include "Riostream.h"
00450 #include "RooSimPdfBuilder.h"
00451 
00452 #include "RooRealVar.h"
00453 #include "RooFormulaVar.h"
00454 #include "RooAbsCategory.h"
00455 #include "RooCategory.h"
00456 #include "RooStringVar.h"
00457 #include "RooMappedCategory.h"
00458 #include "RooRealIntegral.h"
00459 #include "RooDataSet.h"
00460 #include "RooArgSet.h"
00461 #include "RooPlot.h"
00462 #include "RooAddPdf.h"
00463 #include "RooLinearVar.h"
00464 #include "RooTruthModel.h"
00465 #include "RooAddModel.h"
00466 #include "RooProdPdf.h"
00467 #include "RooCustomizer.h"
00468 #include "RooThresholdCategory.h"
00469 #include "RooMultiCategory.h"
00470 #include "RooSuperCategory.h"
00471 #include "RooSimultaneous.h"
00472 #include "RooTrace.h"
00473 #include "RooFitResult.h"
00474 #include "RooDataHist.h"
00475 #include "RooGenericPdf.h"
00476 #include "RooMsgService.h"
00477 
00478 using namespace std ;
00479 
00480 ClassImp(RooSimPdfBuilder)
00481 ;
00482 
00483 
00484 
00485 //_____________________________________________________________________________
00486 RooSimPdfBuilder::RooSimPdfBuilder(const RooArgSet& protoPdfSet) :
00487   _protoPdfSet(protoPdfSet)
00488 {
00489   _compSplitCatSet.setHashTableSize(1000) ;
00490   _splitNodeList.setHashTableSize(10000) ;
00491   _splitNodeListOwned.setHashTableSize(10000) ;
00492 }
00493 
00494 
00495 
00496 
00497 //_____________________________________________________________________________
00498 RooArgSet* RooSimPdfBuilder::createProtoBuildConfig()
00499 {
00500   // Make RooArgSet of configuration objects
00501   RooArgSet* buildConfig = new RooArgSet ;
00502   buildConfig->addOwned(* new RooStringVar("physModels","List and mapping of physics models to include in build","",4096)) ;
00503   buildConfig->addOwned(* new RooStringVar("splitCats","List of categories used for splitting","",1024)) ;
00504 
00505   TIterator* iter = _protoPdfSet.createIterator() ;
00506   RooAbsPdf* proto ;
00507   while ((proto=(RooAbsPdf*)iter->Next())) {
00508     buildConfig->addOwned(* new RooStringVar(proto->GetName(),proto->GetName(),"",4096)) ;
00509   }
00510   delete iter ;
00511 
00512   return buildConfig ;
00513 }
00514 
00515 
00516 
00517 //_____________________________________________________________________________
00518 void RooSimPdfBuilder::addSpecializations(const RooArgSet& specSet) 
00519 {
00520   _splitNodeList.add(specSet) ;
00521 }
00522 
00523 
00524 
00525 //_____________________________________________________________________________
00526 RooSimultaneous* RooSimPdfBuilder::buildPdf(const RooArgSet& buildConfig, const RooArgSet& dependents,
00527                                             const RooArgSet* auxSplitCats, Bool_t verbose)
00528 {
00529   // Initialize needed components
00530   const char* spaceChars = " \t" ;
00531 
00532   // Retrieve physics index category
00533   Int_t buflen = strlen(((RooStringVar*)buildConfig.find("physModels"))->getVal())+1 ;
00534   char *buf = new char[buflen] ;
00535 
00536   strlcpy(buf,((RooStringVar*)buildConfig.find("physModels"))->getVal(),buflen) ;
00537   RooAbsCategoryLValue* physCat(0) ;
00538   if (strstr(buf," : ")) {
00539     const char* physCatName = strtok(buf,spaceChars) ;
00540     physCat = dynamic_cast<RooAbsCategoryLValue*>(dependents.find(physCatName)) ;
00541     if (!physCat) {
00542       coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR physics index category " << physCatName 
00543                             << " not found in dataset variables" << endl ;
00544       delete[] buf ;
00545       return 0 ;      
00546     }
00547     coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: category indexing physics model: " << physCatName << endl ;
00548   }
00549 
00550   // Create list of physics models to be built
00551   char *physName ;
00552   RooArgSet physModelSet ;
00553   if (physCat) {
00554     // Absorb colon token
00555     strtok(0,spaceChars) ;
00556     physName = strtok(0,spaceChars) ;
00557   } else {
00558     physName = strtok(buf,spaceChars) ;
00559   }
00560 
00561   if (!physName) {
00562     coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR: No models specified, nothing to do!" << endl ;
00563     delete[] buf ;
00564     return 0 ;
00565   }
00566 
00567   Bool_t first(kTRUE) ;
00568   RooArgSet stateMap ;
00569   while(physName) {
00570 
00571     char *stateName(0) ;
00572 
00573     // physName may be <state>=<pdfName> or just <pdfName> is state and pdf have identical names
00574     if (strchr(physName,'=')) {
00575       // Must have a physics category for mapping to make sense
00576       if (!physCat) {
00577         coutW(ObjectHandling) << "RooSimPdfBuilder::buildPdf: WARNING: without physCat specification "
00578                          << "<physCatState>=<pdfProtoName> association is meaningless" << endl ;
00579       }
00580       stateName = physName ;
00581       physName = strchr(stateName,'=') ;
00582       if (physName) {
00583         *(physName++) = 0 ;      
00584       }
00585     } else {
00586       stateName = physName ;
00587     }
00588 
00589     RooAbsPdf* physModel = (RooAbsPdf*) (physName ? _protoPdfSet.find(physName) : 0 );
00590     if (!physModel) {
00591       coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR requested physics model " 
00592                             << (physName?physName:"(null)") << " is not defined" << endl ;
00593       delete[] buf ;
00594       return 0 ;
00595     }    
00596 
00597     // Check if state mapping has already been defined
00598     if (stateMap.find(stateName)) {
00599       coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: multiple PDFs specified for state " 
00600                             << stateName << ", only first will be used" << endl ;
00601       continue ;
00602     }
00603 
00604     // Add pdf to list of models to be processed
00605     physModelSet.add(*physModel,kTRUE) ; // silence duplicate insertion warnings
00606 
00607     // Store state->pdf mapping    
00608     stateMap.addOwned(* new RooStringVar(stateName,stateName,physName)) ;
00609 
00610     // Continue with next mapping
00611     physName = strtok(0,spaceChars) ;
00612     if (first) {
00613       first = kFALSE ;
00614     } else if (physCat==0) {
00615       coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: without physCat specification, only the first model will be used" << endl ;
00616       break ;
00617     }
00618   }
00619   coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: list of physics models " << physModelSet << endl ;
00620 
00621 
00622 
00623   // Create list of dataset categories to be used in splitting
00624   TList splitStateList ;
00625   RooArgSet splitCatSet ;
00626 
00627   delete[] buf ; 
00628   buflen = strlen(((RooStringVar*)buildConfig.find("splitCats"))->getVal())+1 ;
00629   buf = new char[buflen] ;
00630   strlcpy(buf,((RooStringVar*)buildConfig.find("splitCats"))->getVal(),buflen) ;
00631 
00632   char *catName = strtok(buf,spaceChars) ;
00633   char *stateList(0) ;
00634   while(catName) {
00635 
00636     // Chop off optional list of selected states
00637     char* tokenPtr(0) ;
00638     if (strchr(catName,'(')) {
00639 
00640       catName = strtok_r(catName,"(",&tokenPtr) ;
00641       stateList = strtok_r(0,")",&tokenPtr) ;
00642 
00643     } else {
00644       stateList = 0 ;
00645     }
00646 
00647     RooCategory* splitCat = catName ? dynamic_cast<RooCategory*>(dependents.find(catName)) : 0 ;
00648     if (!splitCat) {
00649       coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR requested split category " << (catName?catName:"(null)") 
00650                             << " is not a RooCategory in the dataset" << endl ;
00651       delete[] buf ;
00652       return 0 ;
00653     }
00654     splitCatSet.add(*splitCat) ;
00655 
00656     // Process optional state list
00657     if (stateList) {
00658       coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: splitting of category " << catName 
00659                        << " restricted to states (" << stateList << ")" << endl ;
00660 
00661       // Create list named after this splitCat holding its selected states
00662       TList* slist = new TList ;
00663       slist->SetName(catName) ;
00664       splitStateList.Add(slist) ;
00665 
00666       char* stateLabel = strtok_r(stateList,",",&tokenPtr) ;
00667 
00668       while(stateLabel) {
00669         // Lookup state label and require it exists
00670         const RooCatType* type = splitCat->lookupType(stateLabel) ;
00671         if (!type) {
00672           coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR splitCat " << splitCat->GetName() 
00673                                 << " doesn't have a state named " << stateLabel << endl ;
00674           splitStateList.Delete() ;
00675           delete[] buf ;
00676           return 0 ;
00677         }
00678         slist->Add((TObject*)type) ;
00679 
00680         stateLabel = strtok_r(0,",",&tokenPtr) ;
00681       }
00682     }
00683     
00684     catName = strtok(0,spaceChars) ;
00685   }
00686   if (physCat) splitCatSet.add(*physCat) ;
00687   RooSuperCategory masterSplitCat("masterSplitCat","Master splitting category",splitCatSet) ;
00688   
00689   coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: list of splitting categories " << splitCatSet << endl ;
00690 
00691   // Clone auxiliary split cats and attach to splitCatSet
00692   RooArgSet auxSplitSet ;
00693   RooArgSet* auxSplitCloneSet(0) ;
00694   if (auxSplitCats) {
00695     // Deep clone auxililary split cats
00696     auxSplitCloneSet = (RooArgSet*) auxSplitCats->snapshot(kTRUE) ;
00697     if (!auxSplitCloneSet) {
00698       coutE(InputArguments) << "RooSimPdfBuilder::buildPdf(" << GetName() << ") Couldn't deep-clone set auxiliary splitcats, abort." << endl ;
00699       delete[] buf ;
00700       return 0 ;
00701     }
00702 
00703     TIterator* iter = auxSplitCats->createIterator() ;
00704     RooAbsArg* arg ;
00705     while((arg=(RooAbsArg*)iter->Next())) {
00706       // Find counterpart in cloned set
00707       RooAbsArg* aux = auxSplitCats->find(arg->GetName()) ;
00708 
00709       // Check that there is no fundamental splitCat in the dataset with the bane of the auxiliary split
00710       if (splitCatSet.find(aux->GetName())) {
00711         coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: dataset contains a fundamental splitting category " << endl 
00712                               << " with the same name as an auxiliary split function (" << aux->GetName() << "). " << endl 
00713                               << " Auxiliary split function will be ignored" << endl ;
00714         continue ;
00715       }
00716 
00717       // Check that all servers of this aux cat are contained in splitCatSet
00718       RooArgSet* parSet = aux->getParameters(splitCatSet) ;
00719       if (parSet->getSize()>0) {
00720         coutW(InputArguments) << "RooSimPdfBuilder::buildPdf: WARNING: ignoring auxiliary category " << aux->GetName() 
00721                               << " because it has servers that are not listed in splitCatSet: " << *parSet << endl ;
00722         delete parSet ;
00723         continue ;
00724       }
00725 
00726       // Redirect servers to splitCatSet
00727       aux->recursiveRedirectServers(splitCatSet) ;
00728 
00729       // Add top level nodes to auxSplitSet
00730       auxSplitSet.add(*aux) ;
00731     }
00732     delete iter ;
00733 
00734     coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: list of auxiliary splitting categories " << auxSplitSet << endl ;
00735   }
00736 
00737 
00738   TList* customizerList = new TList ;
00739 
00740   // Loop over requested physics models and build components
00741   TIterator* physIter = physModelSet.createIterator() ;
00742   RooAbsPdf* physModel ;
00743   while((physModel=(RooAbsPdf*)physIter->Next())) {
00744     coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: processing physics model " << physModel->GetName() << endl ;
00745 
00746     RooCustomizer* physCustomizer = new RooCustomizer(*physModel,masterSplitCat,_splitNodeList) ;
00747     customizerList->Add(physCustomizer) ;
00748 
00749     // Parse the splitting rules for this physics model
00750     RooStringVar* ruleStr = (RooStringVar*) buildConfig.find(physModel->GetName()) ;
00751     if (ruleStr) {
00752 
00753       delete[] buf ; 
00754       buflen = strlen(ruleStr->getVal())+1 ;
00755       buf = new char[buflen] ;
00756 
00757       strlcpy(buf,ruleStr->getVal(),buflen) ;
00758       char *tokenPtr(0) ;
00759 
00760       char* token = strtok_r(buf,spaceChars,&tokenPtr) ;
00761       
00762       enum Mode { SplitCat, Colon, ParamList } ;
00763       Mode mode(SplitCat) ;
00764 
00765       char* splitCatName ;
00766       RooAbsCategory* splitCat(0) ;
00767 
00768       while(token) {
00769 
00770         switch (mode) {
00771         case SplitCat:
00772           {
00773             splitCatName = token ;
00774                     
00775             if (strchr(splitCatName,',')) {
00776               // Composite splitting category
00777               
00778               // Check if already instantiated
00779               splitCat = (RooAbsCategory*) _compSplitCatSet.find(splitCatName) ;              
00780               TString origCompCatName(splitCatName) ;
00781               if (!splitCat) {
00782                 // Build now
00783 
00784                 char *tokptr = 0;
00785                 char *catName2 = strtok_r(token,",",&tokptr) ;
00786 
00787                 RooArgSet compCatSet ;
00788                 while(catName2) {
00789                   RooAbsArg* cat = splitCatSet.find(catName2) ;
00790                   
00791                   // If not, check if it is an auxiliary splitcat
00792                   if (!cat) {
00793                     cat = (RooAbsCategory*) auxSplitSet.find(catName2) ;
00794                   }
00795 
00796                   if (!cat) {
00797                     coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR " << catName2
00798                                           << " not found in the primary or auxilary splitcat list" << endl ;
00799                     customizerList->Delete() ;
00800                     delete customizerList ;
00801 
00802                     splitStateList.Delete() ;
00803                     delete[] buf ;
00804                     return 0 ;
00805                   }
00806                   compCatSet.add(*cat) ;
00807 
00808                   catName2 = strtok_r(0,",",&tokptr) ;
00809                 }               
00810 
00811 
00812                 // check that any auxSplitCats in compCatSet do not depend on any other
00813                 // fundamental or auxiliary splitting categories in the composite set.
00814                 TIterator* iter = compCatSet.createIterator() ;
00815                 RooAbsArg* arg ;
00816                 while((arg=(RooAbsArg*)iter->Next())) {
00817                   RooArgSet tmp(compCatSet) ;
00818                   tmp.remove(*arg) ;
00819                   if (arg->dependsOnValue(tmp)) {
00820                     coutE(InputArguments) << "RooSimPdfBuilder::buildPDF: ERROR: Ill defined split: auxiliary splitting category " << arg->GetName() 
00821                                           << " used in composite split " << compCatSet << " depends on one or more of the other splitting categories in the composite split" << endl ;
00822                     
00823                     // Cleanup and axit
00824                     customizerList->Delete() ;
00825                     delete customizerList ;
00826                     splitStateList.Delete() ;
00827                     delete[] buf ;
00828                     return 0 ;
00829                   }
00830                 }
00831                 delete iter ;
00832 
00833                 splitCat = new RooMultiCategory(origCompCatName,origCompCatName,compCatSet) ;
00834                 _compSplitCatSet.addOwned(*splitCat) ;
00835                 //cout << "composite splitcat: " << splitCat->GetName() ;
00836               }
00837             } else {
00838               // Simple splitting category
00839               
00840               // First see if it is a simple splitting category
00841               splitCat = (RooAbsCategory*) splitCatSet.find(splitCatName) ;
00842 
00843               // If not, check if it is an auxiliary splitcat
00844               if (!splitCat) {
00845                 splitCat = (RooAbsCategory*) auxSplitSet.find(splitCatName) ;
00846               }
00847 
00848               if (!splitCat) {
00849                 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR splitting category " 
00850                                       << splitCatName << " not found in the primary or auxiliary splitcat list" << endl ;
00851                 customizerList->Delete() ;
00852                 delete customizerList ;
00853                 splitStateList.Delete() ;
00854                 delete[] buf ;
00855                 return 0 ;
00856               }
00857             }
00858             
00859             mode = Colon ;
00860             break ;
00861           }
00862         case Colon:
00863           {
00864             if (strcmp(token,":")) {
00865               coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR in parsing, expected ':' after " 
00866                                     << splitCat << ", found " << token << endl ;
00867               customizerList->Delete() ;
00868               delete customizerList ;
00869               splitStateList.Delete() ;
00870               delete[] buf ;
00871               return 0 ;            
00872             }
00873             mode = ParamList ;
00874             break ;
00875           }
00876         case ParamList:
00877           {
00878             // Verify the validity of the parameter list and build the corresponding argset
00879             RooArgSet splitParamList ;
00880             RooArgSet* paramList = physModel->getParameters(dependents) ;
00881 
00882             // wve -- add nodes to parameter list
00883             RooArgSet* compList = physModel->getComponents() ;
00884             paramList->add(*compList) ;
00885             delete compList ;
00886 
00887             Bool_t lastCharIsComma = (token[strlen(token)-1]==',') ;
00888 
00889             char *tokptr = 0 ;
00890             char *paramName = strtok_r(token,",",&tokptr) ;
00891 
00892             // Check for fractional split option 'param_name[remainder_state]'
00893             char *remainderState = 0 ;
00894             char *tokptr2 = 0 ;
00895             if (paramName && strtok_r(paramName,"[",&tokptr2)) {
00896               remainderState = strtok_r(0,"]",&tokptr2) ;
00897             }
00898 
00899             while(paramName) {
00900 
00901               // If fractional split is specified, check that remainder state is a valid state of this split cat
00902               if (remainderState) {
00903                 if (!splitCat->lookupType(remainderState)) {
00904                   coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR fraction split of parameter " 
00905                                         << paramName << " has invalid remainder state name: " << remainderState << endl ;
00906                   delete paramList ;
00907                   customizerList->Delete() ;
00908                   delete customizerList ;
00909                   splitStateList.Delete() ;
00910                   delete[] buf ;
00911                   return 0 ;
00912                 }
00913               }       
00914 
00915               RooAbsArg* param = paramList->find(paramName) ;
00916               if (!param) {
00917                 coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR " << paramName 
00918                                       << " is not a parameter of physics model " << physModel->GetName() << endl ;
00919                 delete paramList ;
00920                 customizerList->Delete() ;
00921                 delete customizerList ;
00922                 splitStateList.Delete() ;
00923                 delete[] buf ;
00924                 return 0 ;
00925               }
00926               splitParamList.add(*param) ;
00927               
00928               // Build split leaf of fraction splits here
00929               if (remainderState) {
00930 
00931                 // Check if we are splitting a real-valued parameter
00932                 if (!dynamic_cast<RooAbsReal*>(param)) {
00933                   coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR fraction split requested of non-real valued parameter " 
00934                                         << param->GetName() << endl ;
00935                   delete paramList ;
00936                   customizerList->Delete() ;
00937                   delete customizerList ;
00938                   splitStateList.Delete() ;
00939                   delete[] buf ;
00940                   return 0 ;
00941                 }
00942 
00943                 // Check if we are doing a restricted build
00944                 TList* remStateSplitList = static_cast<TList*>(splitStateList.FindObject(splitCat->GetName())) ;
00945 
00946                 // If so, check if remainder state is actually being built.
00947                 if (remStateSplitList && !remStateSplitList->FindObject(remainderState)) {
00948                   coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR " << paramName 
00949                                         << " remainder state " << remainderState << " in parameter split " 
00950                                         << param->GetName() << " is not actually being built" << endl ;
00951                   delete paramList ;
00952                   customizerList->Delete() ;
00953                   delete customizerList ;
00954                   splitStateList.Delete() ;
00955                   delete[] buf ;
00956                   return 0 ;              
00957                 }
00958 
00959                 TIterator* iter = splitCat->typeIterator() ;
00960                 RooCatType* type ;
00961                 RooArgList fracLeafList ;
00962                 TString formExpr("1") ;
00963                 Int_t i(0) ;
00964 
00965                 while((type=(RooCatType*)iter->Next())) {
00966 
00967                   // Skip remainder state
00968                   if (!TString(type->GetName()).CompareTo(remainderState)) continue ;
00969 
00970                   // If restricted build is requested, skip states of splitcat that are not built
00971                   if (remStateSplitList && !remStateSplitList->FindObject(type->GetName())) {
00972                     continue ;
00973                   }
00974                   
00975                   // Construct name of split leaf
00976                   TString splitLeafName(param->GetName()) ;
00977                   splitLeafName.Append("_") ;
00978                   splitLeafName.Append(type->GetName()) ;
00979                   
00980                   // Check if split leaf already exists
00981                   RooAbsArg* splitLeaf = _splitNodeList.find(splitLeafName) ;
00982                   if (!splitLeaf) {
00983                     // If not create it now
00984                     splitLeaf = (RooAbsArg*) param->clone(splitLeafName) ;
00985                     _splitNodeList.add(*splitLeaf) ;
00986                     _splitNodeListOwned.addOwned(*splitLeaf) ;
00987                   }
00988                   fracLeafList.add(*splitLeaf) ;
00989                   formExpr.Append(Form("-@%d",i++)) ;
00990                 }
00991                 delete iter ;           
00992 
00993                 // Construct RooFormulaVar expresssing remainder of fraction            
00994                 TString remLeafName(param->GetName()) ;
00995                 remLeafName.Append("_") ;
00996                 remLeafName.Append(remainderState) ;
00997 
00998                 // Check if no specialization was already specified for remainder state
00999                 if (!_splitNodeList.find(remLeafName)) {
01000                   RooAbsArg* remLeaf = new RooFormulaVar(remLeafName,formExpr,fracLeafList) ;
01001                   _splitNodeList.add(*remLeaf) ;
01002                   _splitNodeListOwned.addOwned(*remLeaf) ;
01003                   coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: creating remainder fraction formula for " << remainderState 
01004                                    << " specialization of split parameter " << param->GetName() << " " << formExpr << endl ;
01005                 }
01006               }
01007 
01008               // Parse next parameter name
01009               paramName = strtok_r(0,",",&tokptr) ;
01010               if (paramName && strtok_r(paramName,"[",&tokptr2)) {
01011                 remainderState = strtok_r(0,"]",&tokptr2) ;
01012               }
01013             }
01014 
01015             // Add the rule to the appropriate customizer ;
01016             physCustomizer->splitArgs(splitParamList,*splitCat) ;
01017 
01018             delete paramList ;
01019 
01020             if (!lastCharIsComma) mode = SplitCat ;
01021             break ;
01022           }
01023         }
01024 
01025         token = strtok_r(0,spaceChars,&tokenPtr) ;
01026 
01027       }
01028       if (mode!=SplitCat) {
01029         coutE(InputArguments) << "RooSimPdfBuilder::buildPdf: ERROR in parsing, expected " 
01030                               << (mode==Colon?":":"parameter list") << " after " << (token?token:"(null)") << endl ;
01031       }
01032 
01033       //RooArgSet* paramSet = physModel->getParameters(dependents) ;
01034     } else {
01035       coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: no splitting rules for " << physModel->GetName() << endl ;
01036     }    
01037   }
01038   
01039   coutI(ObjectHandling)  << "RooSimPdfBuilder::buildPdf: configured customizers for all physics models" << endl ;
01040   if (oodologI((TObject*)0,ObjectHandling)) {
01041     customizerList->Print() ;
01042   }
01043 
01044   // Create fit category from physCat and splitCatList ;
01045   RooArgSet fitCatList ;
01046   if (physCat) fitCatList.add(*physCat) ;
01047   fitCatList.add(splitCatSet) ;
01048   TIterator* fclIter = fitCatList.createIterator() ;
01049   RooSuperCategory *fitCat = new RooSuperCategory("fitCat","fitCat",fitCatList) ;
01050 
01051   // Create master PDF 
01052   RooSimultaneous* simPdf = new RooSimultaneous("simPdf","simPdf",*fitCat) ;
01053 
01054   // Add component PDFs to master PDF
01055   TIterator* fcIter = fitCat->typeIterator() ;
01056 
01057   RooCatType* fcState ;  
01058   while((fcState=(RooCatType*)fcIter->Next())) {
01059     // Select fitCat state
01060     fitCat->setLabel(fcState->GetName()) ;
01061 
01062     // Check if this fitCat state is selected
01063     fclIter->Reset() ;
01064     RooAbsCategory* splitCat ;
01065     Bool_t select(kTRUE) ;
01066     while((splitCat=(RooAbsCategory*)fclIter->Next())) {
01067       // Find selected state list 
01068       TList* slist = (TList*) splitStateList.FindObject(splitCat->GetName()) ;
01069       if (!slist) continue ;
01070       RooCatType* type = (RooCatType*) slist->FindObject(splitCat->getLabel()) ;
01071       if (!type) {
01072         select = kFALSE ;
01073       }
01074     }
01075     if (!select) continue ;
01076 
01077     
01078     // Select appropriate PDF for this physCat state
01079     RooCustomizer* physCustomizer ;
01080     if (physCat) {      
01081       RooStringVar* physNameVar = (RooStringVar*) stateMap.find(physCat->getLabel()) ;
01082       if (!physNameVar) continue ;
01083       physCustomizer = (RooCustomizer*) customizerList->FindObject(physNameVar->getVal());  
01084     } else {
01085       physCustomizer = (RooCustomizer*) customizerList->First() ;
01086     }
01087 
01088     coutI(ObjectHandling) << "RooSimPdfBuilder::buildPdf: Customizing physics model " << physCustomizer->GetName() 
01089                      << " for mode " << fcState->GetName() << endl ;    
01090 
01091     // Customizer PDF for current state and add to master simPdf
01092     RooAbsPdf* fcPdf = (RooAbsPdf*) physCustomizer->build(masterSplitCat.getLabel(),verbose) ;
01093     simPdf->addPdf(*fcPdf,fcState->GetName()) ;
01094   }
01095   delete fcIter ;
01096 
01097   // Move customizers (owning the cloned branch node components) to the attic
01098   _retiredCustomizerList.AddAll(customizerList) ;
01099   delete customizerList ;
01100 
01101   delete fclIter ;
01102   splitStateList.Delete() ;
01103 
01104   if (auxSplitCloneSet) delete auxSplitCloneSet ;
01105   delete physIter ;
01106 
01107   delete[] buf ;
01108   _simPdfList.push_back(simPdf) ;
01109   _fitCatList.push_back(fitCat) ;
01110   return simPdf ;
01111 }
01112 
01113 
01114 
01115 
01116 
01117 //_____________________________________________________________________________
01118 RooSimPdfBuilder::~RooSimPdfBuilder() 
01119 {
01120   _retiredCustomizerList.Delete() ;
01121 
01122   std::list<RooSimultaneous*>::iterator iter = _simPdfList.begin() ;
01123   while(iter != _simPdfList.end()) {
01124     delete *iter ;
01125     ++iter ;
01126   }
01127 
01128   std::list<RooSuperCategory*>::iterator iter2 = _fitCatList.begin() ;
01129   while(iter2 != _fitCatList.end()) {
01130     delete *iter2 ;
01131     ++iter2 ;
01132   }
01133   
01134 }
01135  
01136 

Generated on Tue Jul 5 15:07:31 2011 for ROOT_528-00b_version by  doxygen 1.5.1