/***************************************************************************** * Project: RooFit * * Package: RooFitCore * * File: $Id: RooSimPdfBuilder.cc,v 1.32 2005/06/20 15:45:14 wverkerke Exp $ * Authors: * * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu * * DK, David Kirkby, UC Irvine, dkirkby@uci.edu * * * * Copyright (c) 2000-2005, Regents of the University of California * * and Stanford University. All rights reserved. * * * * Redistribution and use in source and binary forms, * * with or without modification, are permitted according to the terms * * listed in LICENSE (http://roofit.sourceforge.net/license.txt) * *****************************************************************************/ // -- CLASS DESCRIPTION [PDF] -- // //// //
// RooSimPdfBuilder is a powerful tool to build RooSimultaneous // PDFs that are defined in terms component PDFs that are identical in // structure, but have different parameters. //
// //
// The following example demonstrates the essence of RooSimPdfBuilder:
// Given a dataset D with a RooRealVar X and a RooCategory C that has
// state C1 and C2.
//
//
//
// Coding this example directly with RooFit classes gives // (we assume dataset D and variables C and X have been declared previously) //
// // RooRealVar m("m","mean of gaussian",-10,10) ; // RooRealVar s_C1("s_C1","sigma of gaussian C1",0,20) ; // RooRealVar s_C2("s_C2","sigma of gaussian C2",0,20) ; // RooGaussian gauss_C1("gauss_C1","gaussian C1",X,m,s_C1) ; // RooGaussian gauss_C2("gauss_C2","gaussian C2",X,m,s_C2) ; // // RooRealVar k_C1("k_C1","ArgusBG kappa parameter C1",-50,0) ; // RooRealVar k_C2("k_C2","ArgusBG kappa parameter C2",-50,0) ; // RooRealVar xm("xm","ArgusBG cutoff point",5.29) ; // RooArgusBG argus_C1("argus_C1","argus background C1",X,k_C1,xm) ; // RooArgusBG argus_C2("argus_C2","argus background C2",X,k_C2,xm) ; // // RooRealVar gfrac("gfrac","fraction of gaussian",0.,1.) ; // RooAddPdf pdf_C1("pdf_C1","gauss+argus_C1",RooArgList(gauss_C1,argus_C1),gfrac) ; // RooAddPdf pdf_C2("pdf_C2","gauss+argus_C2",RooArgList(gauss_C2,argus_C2),gfrac) ; // // RooSimultaneous simPdf("simPdf","simPdf",C) ; // simPdf.addPdf(pdf_C1,"C1") ; // simPdf.addPdf(pdf_C2,"C2") ; // //// Coding this example with RooSimPdfBuilder gives //
// // RooRealVar m("m","mean of gaussian",-10,10) ; // RooRealVar s("s","sigma of gaussian",0,20) ; // RooGaussian gauss("gauss","gaussian",X,m,s) ; // // RooRealVar k("k","ArgusBG kappa parameter",-50,0) ; // RooRealVar xm("xm","ArgusBG cutoff point",5.29) ; // RooArgusBG argus("argus","argus background",X,k,xm) ; // // RooRealVar gfrac("gfrac","fraction of gaussian",0.,1.) ; // RooAddPdf pdf("pdf","gauss+argus",RooArgList(gauss,argus),gfrac) ; // // RooSimPdfBuilder builder(pdf) ; // RooArgSet* config = builder.createProtoBuildConfig() ; // (*config)["physModels"] = "pdf" ; // Name of the PDF we are going to work with // (*config)["splitCats"] = "C" ; // Category used to differentiate sub-datasets // (*config)["pdf"] = "C : k,s" ; // Prescription to taylor PDF parameters k and s // // for each data subset designated by C states // RooSimultaneous* simPdf = builder.buildPdf(*config,&D) ; // //// The above snippet of code demonstrates the concept of RooSimPdfBuilder: // the user defines a single 'prototype' PDF that defines the structure of all // PDF components of the RooSimultaneous PDF to be built. RooSimPdfBuilder // then takes this prototype and replicates it as a component // PDF for each state of the C index category. //
//// In the above example RooSimPdfBuilder // will first replicate k and s into // k_C1,k_C2 and s_C1,s_C2, as prescribed in the // configuration. Then it will recursively replicate all PDF nodes that depend on // the 'split' parameter nodes: gauss into gauss_C1,C2, argus // into argus_C1,C2 and finally pdf into pdf_C1,pdf_C2. // When PDFs for all states of C have been replicated // they are assembled into a RooSimultaneous PDF, which is returned by the buildPdf() // method. //
//// Although in this very simple example the use of RooSimPdfBuilder doesn't // reduce the amount of code much, it is already easier to read and maintain // because there is no duplicate code. As the complexity of the RooSimultaneous // to be built increases, the advantages of RooSimPdfBuilder will become more and // more apparent. //
// // //// Each builder configuration needs at minumum two lines, physModels and splitCats, which identify // the ingredients of the build. In this section we only explain the building rules for // builds from a single prototype PDF. In that case the physModels line always reads //
// // physModels = {pdfName} // //// The second line, splitCats, indicates which categories are going to be used to // differentiate the various subsets of the 'master' input data set. You can enter // a single category here, or multiple if necessary: //
// // splitCats = {catName} [{catName} ...] // //// All listed splitcats must be RooCategories that appear in the dataset provided to // RooSimPdfBuilder::buildPdf() //
//// The parameter splitting prescriptions, the essence of each build configuration // can be supplied in a third line carrying the name of the pdf listed in physModels //
// // pdfName = {splitCat} : {parameter} [,{parameter},....] // //// Each pdf can have only one line with splitting rules, but multiple rules can be // supplied in each line, e.g. //
// // pdfName = {splitCat} : {parameter} [,{parameter},....] // {splitCat} : {parameter} [,{parameter},....] // //// Conversely, each parameter can only have one splitting prescription, but it may be split // by multiple categories, e.g. //
// // pdfName = {splitCat1},{splitCat2} : {parameter} // //// instructs RooSimPdfBuilder to build a RooSuperCategory // of {splitCat1} and {splitCat2} // and split {parameter} with that RooSuperCategory //
//// Here is an example of a builder configuration that uses several of the options discussed // above: //
// // physModels = pdf // splitCats = tagCat runBlock // pdf = tagCat : signalRes,bkgRes // runBlock : fudgeFactor // tagCat,runBlock : kludgeParam // // //// The prototype builder configuration returned by // RooSimPdfBuilder::createProtoBuildConfig() is a pointer to a RooArgSet filled with // initially blank RooStringVars named physModels,splitCats and one additional for each // PDF supplied to the RooSimPdfBuilders constructor (with the same name) //
//// In macro code, the easiest way to assign new values to these RooStringVars // is to use RooArgSets array operator and the RooStringVars assignment operator, e.g. //
// // (*config)["physModels"] = "Blah" ; // //// To enter multiple splitting rules simply separate consecutive rules by whitespace // (not newlines), e.g. //
// // (*config)["physModels"] = "Blah " // << note trailing space here // "Blah 2" ; // //// In this example, the C++ compiler will concatenate the two string literals (without inserting // any whitespace), so the extra space after 'Blah' is important here. //
//// Alternatively, you can read the configuration from an ASCII file, as you can // for any RooArgSet using RooArgSet::readFromFile(). In that case the ASCII file // can follow the syntax of the examples above and the '\\' line continuation // sequence can be used to fold a long splitting rule over multiple lines. //
// // RooArgSet* config = builder.createProtoBuildConfig() ; // config->readFromFile("config.txt") ; // // --- config.txt ---------------- // physModels = pdf // splitCats = tagCat // pdf = tagCat : bogusPar // ------------------------------- // // // //// It is also possible to build a RooSimultaneous PDF from multiple PDF prototypes. // This is appropriate for cases where the input prototype PDF would otherwise be // a RooSimultaneous PDF by itself. In such cases we don't feed a single // RooSimultaneous PDF into RooSimPdfBuilder, instead we feed it its ingredients and // add a prescription to the builder configuration that corresponds to the // PDF-category state mapping of the prototype RooSimultaneous. //
//// The constructor of the RooSimPdfBuilder will look as follows: //
// // RooSimPdfBuilder builder(RooArgSet(pdfA,pdfB,...)) ; // //// The physModels line is now expanded to carry the pdf->state mapping information // that the prototype RooSimultaneous would have. I.e. //
// // physModels = mode : pdfA=modeA pdfB=modeB // //// is equivalent to a prototype RooSimultaneous constructed as //
// // RooSimultanous simPdf("simPdf","simPdf",mode); // simPdf.addPdf(pdfA,"modeA") ; // simPdf.addPdf(pdfB,"modeB") ; // //// The rest of the builder configuration works the same, except that // each prototype PDF now has its own set of splitting rules, e.g. //
// // physModels = mode : pdfA=modeA pdfB=modeB // splitCats = tagCat // pdfA = tagCat : bogusPar // pdfB = tagCat : fudgeFactor // //
// Please note that
//
//
// there are still only 2 sets of splitting rules: one for pdfA and one // for pdfB. However, you can differentiate between modeA and modeC in // the above example. The technique is to use mode as splitting category, e.g. //
// // physModels = mode : pdfA=modeA pdfB=modeB pdfA=modeC pdfB=modeD // splitCats = tagCat // pdfA = tagCat : bogusPar // mode : funnyPar // pdfB = mode : kludgeFactor // //// will result in an individual set of funnyPar parameters for modeA and modeC // labeled funnyPar_modeA and funnyPar_modeB and an individual set of // kludgeFactor parameters for pdfB, kludgeFactor_modeB and kludgeFactor_modeD. // Please note that for splits in the master index category (mode) only the // applicable states are built (A,C for pdfA, B,D for pdfB) //
// // //// You can request to limit the list of states of each splitCat that // will be considered in the build. This limitation is requested in the // each build as follows: //
// // splitCats = tagCat(Lep,Kao) RunBlock(Run1) // //// In this example the splitting of tagCat is limited to states Lep,Kao // and the splitting of runBlock is limited to Run1. The splits apply // globally to each build, i.e. every parameter split requested in this // build will be limited according to these specifications. //
//// NB: Partial builds have no pdf associated with the unbuilt states of the // limited splits. Running such a pdf on a dataset that contains data with // unbuilt states will result in this data being ignored completely. //
// // //// It is possible to make non-trivial parameter splits with RooSimPdfBuilder. // Trivial splits are considered simple splits in one (fundamental) category // in the dataset or a split in a RooSuperCategory 'product' of multiple // fundamental categories in the dataset. Non-trivial splits can be performed // using an intermediate 'category function' (RooMappedCategory, // RooGenericCategory,RooThresholdCategory etc), i.e. any RooAbsCategory // derived objects that calculates its output as function of one or more // input RooRealVars and/or RooCategories. //
//// Such 'function categories' objects must be constructed by the user prior // to building the PDF. In the RooSimPdfBuilder::buildPdf() function these // objects can be passed in an optional RooArgSet called 'auxiliary categories': //
// // const RooSimultaneous* buildPdf(const RooArgSet& buildConfig, const RooAbsData* dataSet, // const RooArgSet& auxSplitCats, Bool_t verbose=kFALSE) { // ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^ // //// Objects passed in this argset can subsequently be used in the build configuration, e.g. //
// // RooMappedCategory tagMap("tagMap","Mapped tagging category",tagCat,"CutBased") ; // tagMap.map("Lep","CutBased") ; // tagMap.map("Kao","CutBased") ; // tagMap.map("NT*","NeuralNet") ; // ... // builder.buildPdf(config,D,tagMap) ; // ^^^^^^ //// In the above example signalRes will be split in signalRes_Kao,signalRes_Lep, // signalRes_NT1,signalRes_NT2, while fudgeFactor will be split in fudgeFactor_CutBased // and fudgeFactor_NeuralNet. //
//// Category functions passed in the auxSplitCats RooArgSet can be used regularly // in the splitting configuration. They should not be listed in splitCats, // but must be able to be expressed completely in terms of the splitCats that // are listed. //
// // //// Sometimes you want to build multiple PDFs for independent consecutive fits // that share some of their parameters. For example, we have two prototype PDFs // pdfA(x;p,q) and pdfB(x;p,r) that have a common parameter p. // We want to build a RooSimultaneous for both pdfA and B, // which involves a split of parameter p and we would like to build the // simultaneous pdfs simA and simB such that still share their (now split) parameters // p_XXX. This is accomplished by letting a single instance of RooSimPdfBuilder handle // the builds of both pdfA and pdfB, as illustrated in this example: //
// // RooSimPdfBuilder builder(RooArgSet(pdfA,pdfB)) ; // // RooArgSet* configA = builder.createProtoBuildConfig() ; // (*configA)["physModels"] = "pdfA" ; // (*configA)["splitCats"] = "C" ; // (*configA)["pdf"] = "C : p" ; // RooSimultaneous* simA = builder.buildPdf(*configA,&D) ; // // RooArgSet* configB = builder.createProtoBuildConfig() ; // (*configA)["physModels"] = "pdfB" ; // (*configA)["splitCats"] = "C" ; // (*configA)["pdf"] = "C : p" ; // RooSimultaneous* simB = builder.buildPdf(*configB,&D) ; // // //// The RooSimPdfBuilder instance owns all the objects it creates, including the top-level // RooSimultaneous returned by buildPdf(). Therefore the builder instance should // exist as long as the constructed PDFs needs to exist. //
// // // #ifndef _REENTRANT #define _REENTRANT #endif #include "RooFit.h" #include <string.h> #include <string.h> // Matthew D. Langston <langston@SLAC.Stanford.EDU> // // Microsoft doesn't supply strings.h. However, strings.h is only // required for strtok_r (a reentrant version of strtok), and // Microsoft's version of strtok is already thread safe. #ifndef _WIN32 #include <strings.h> #endif #include "RooSimPdfBuilder.h" #include "RooRealVar.h" #include "RooFormulaVar.h" #include "RooAbsCategory.h" #include "RooCategory.h" #include "RooStringVar.h" #include "RooMappedCategory.h" #include "RooRealIntegral.h" #include "RooDataSet.h" #include "RooArgSet.h" #include "RooPlot.h" #include "RooAddPdf.h" #include "RooLinearVar.h" #include "RooTruthModel.h" #include "RooAddModel.h" #include "RooProdPdf.h" #include "RooCustomizer.h" #include "RooThresholdCategory.h" #include "RooMultiCategory.h" #include "RooSuperCategory.h" #include "RooSimultaneous.h" #include "RooTrace.h" #include "RooFitResult.h" #include "RooDataHist.h" #include "RooGenericPdf.h" ClassImp(RooSimPdfBuilder) ; RooSimPdfBuilder::RooSimPdfBuilder(const RooArgSet& protoPdfSet) : _protoPdfSet(protoPdfSet) { _compSplitCatSet.setHashTableSize(1000) ; _splitNodeList.setHashTableSize(10000) ; } RooArgSet* RooSimPdfBuilder::createProtoBuildConfig() { // Make RooArgSet of configuration objects RooArgSet* buildConfig = new RooArgSet ; buildConfig->addOwned(* new RooStringVar("physModels","List and mapping of physics models to include in build","",4096)) ; buildConfig->addOwned(* new RooStringVar("splitCats","List of categories used for splitting","",1024)) ; TIterator* iter = _protoPdfSet.createIterator() ; RooAbsPdf* proto ; while ((proto=(RooAbsPdf*)iter->Next())) { buildConfig->addOwned(* new RooStringVar(proto->GetName(),proto->GetName(),"",4096)) ; } delete iter ; return buildConfig ; } void RooSimPdfBuilder::addSpecializations(const RooArgSet& specSet) { _splitNodeList.addOwned(specSet) ; } const RooSimultaneous* RooSimPdfBuilder::buildPdf(const RooArgSet& buildConfig, const RooArgSet& dependents, const RooArgSet* auxSplitCats, Bool_t verbose) { // Initialize needed components const char* spaceChars = " \t" ; // Retrieve physics index category char buf[2048] ; strcpy(buf,((RooStringVar*)buildConfig.find("physModels"))->getVal()) ; RooAbsCategoryLValue* physCat(0) ; if (strstr(buf," : ")) { const char* physCatName = strtok(buf,spaceChars) ; physCat = dynamic_cast<RooAbsCategoryLValue*>(dependents.find(physCatName)) ; if (!physCat) { cout << "RooSimPdfBuilder::buildPdf: ERROR physics index category " << physCatName << " not found in dataset variables" << endl ; return 0 ; } cout << "RooSimPdfBuilder::buildPdf: category indexing physics model: " << physCatName << endl ; } // Create list of physics models to be built char *physName ; RooArgSet physModelSet ; if (physCat) { // Absorb colon token strtok(0,spaceChars) ; physName = strtok(0,spaceChars) ; } else { physName = strtok(buf,spaceChars) ; } if (!physName) { cout << "RooSimPdfBuilder::buildPdf: ERROR: No models specified, nothing to do!" << endl ; return 0 ; } Bool_t first(kTRUE) ; RooArgSet stateMap ; while(physName) { char *stateName(0) ; // physName may be <state>=<pdfName> or just <pdfName> is state and pdf have identical names if (strchr(physName,'=')) { // Must have a physics category for mapping to make sense if (!physCat) { cout << "RooSimPdfBuilder::buildPdf: WARNING: without physCat specification " << "<physCatState>=<pdfProtoName> association is meaningless" << endl ; } stateName = physName ; physName = strchr(stateName,'=') ; *(physName++) = 0 ; } else { stateName = physName ; } RooAbsPdf* physModel = (RooAbsPdf*) _protoPdfSet.find(physName) ; if (!physModel) { cout << "RooSimPdfBuilder::buildPdf: ERROR requested physics model " << physName << " is not defined" << endl ; return 0 ; } // Check if state mapping has already been defined if (stateMap.find(stateName)) { cout << "RooSimPdfBuilder::buildPdf: WARNING: multiple PDFs specified for state " << stateName << ", only first will be used" << endl ; continue ; } // Add pdf to list of models to be processed physModelSet.add(*physModel,kTRUE) ; // silence duplicate insertion warnings // Store state->pdf mapping stateMap.addOwned(* new RooStringVar(stateName,stateName,physName)) ; // Continue with next mapping physName = strtok(0,spaceChars) ; if (first) { first = kFALSE ; } else if (physCat==0) { cout << "RooSimPdfBuilder::buildPdf: WARNING: without physCat specification, only the first model will be used" << endl ; break ; } } cout << "RooSimPdfBuilder::buildPdf: list of physics models " ; physModelSet.Print("1") ; // Create list of dataset categories to be used in splitting TList splitStateList ; RooArgSet splitCatSet ; strcpy(buf,((RooStringVar*)buildConfig.find("splitCats"))->getVal()) ; char *catName = strtok(buf,spaceChars) ; char *stateList(0) ; while(catName) { // Chop off optional list of selected states char* tokenPtr(0) ; if (strchr(catName,'(')) { #ifndef _WIN32 catName = strtok_r(catName,"(",&tokenPtr) ; stateList = strtok_r(0,")",&tokenPtr) ; #else catName = strtok(catName,"(") ; stateList = strtok(0,")") ; #endif } else { stateList = 0 ; } RooCategory* splitCat = dynamic_cast<RooCategory*>(dependents.find(catName)) ; if (!splitCat) { cout << "RooSimPdfBuilder::buildPdf: ERROR requested split category " << catName << " is not a RooCategory in the dataset" << endl ; return 0 ; } splitCatSet.add(*splitCat) ; // Process optional state list if (stateList) { cout << "RooSimPdfBuilder::buildPdf: splitting of category " << catName << " restricted to states (" << stateList << ")" << endl ; // Create list named after this splitCat holding its selected states TList* slist = new TList ; slist->SetName(catName) ; splitStateList.Add(slist) ; #ifndef _WIN32 char* stateLabel = strtok_r(stateList,",",&tokenPtr) ; #else char* stateLabel = strtok(stateList,",") ; #endif while(stateLabel) { // Lookup state label and require it exists const RooCatType* type = splitCat->lookupType(stateLabel) ; if (!type) { cout << "RooSimPdfBuilder::buildPdf: ERROR splitCat " << splitCat->GetName() << " doesn't have a state named " << stateLabel << endl ; splitStateList.Delete() ; return 0 ; } slist->Add((TObject*)type) ; #ifndef _WIN32 stateLabel = strtok_r(0,",",&tokenPtr) ; #else stateLabel = strtok(0,",") ; #endif } } catName = strtok(0,spaceChars) ; } if (physCat) splitCatSet.add(*physCat) ; RooSuperCategory masterSplitCat("masterSplitCat","Master splitting category",splitCatSet) ; cout << "RooSimPdfBuilder::buildPdf: list of splitting categories " ; splitCatSet.Print("1") ; // Clone auxiliary split cats and attach to splitCatSet RooArgSet auxSplitSet ; RooArgSet* auxSplitCloneSet(0) ; if (auxSplitCats) { // Deep clone auxililary split cats auxSplitCloneSet = (RooArgSet*) auxSplitCats->snapshot(kTRUE) ; if (!auxSplitCloneSet) { cout << "RooSimPdfBuilder::buildPdf(" << GetName() << ") Couldn't deep-clone set auxiliary splitcats, abort." << endl ; return 0 ; } TIterator* iter = auxSplitCats->createIterator() ; RooAbsArg* arg ; while((arg=(RooAbsArg*)iter->Next())) { // Find counterpart in cloned set RooAbsArg* aux = auxSplitCats->find(arg->GetName()) ; // Check that there is no fundamental splitCat in the dataset with the bane of the auxiliary split if (splitCatSet.find(aux->GetName())) { cout << "RooSimPdfBuilder::buildPdf: WARNING: dataset contains a fundamental splitting category " << endl << " with the same name as an auxiliary split function (" << aux->GetName() << "). " << endl << " Auxiliary split function will be ignored" << endl ; continue ; } // Check that all servers of this aux cat are contained in splitCatSet RooArgSet* parSet = aux->getParameters(splitCatSet) ; if (parSet->getSize()>0) { cout << "RooSimPdfBuilder::buildPdf: WARNING: ignoring auxiliary category " << aux->GetName() << " because it has servers that are not listed in splitCatSet: " ; parSet->Print("1") ; delete parSet ; continue ; } // Redirect servers to splitCatSet aux->recursiveRedirectServers(splitCatSet) ; // Add top level nodes to auxSplitSet auxSplitSet.add(*aux) ; } delete iter ; cout << "RooSimPdfBuilder::buildPdf: list of auxiliary splitting categories " ; auxSplitSet.Print("1") ; } TList* customizerList = new TList ; // Loop over requested physics models and build components TIterator* physIter = physModelSet.createIterator() ; RooAbsPdf* physModel ; while((physModel=(RooAbsPdf*)physIter->Next())) { cout << "RooSimPdfBuilder::buildPdf: processing physics model " << physModel->GetName() << endl ; RooCustomizer* physCustomizer = new RooCustomizer(*physModel,masterSplitCat,_splitNodeList) ; customizerList->Add(physCustomizer) ; // Parse the splitting rules for this physics model RooStringVar* ruleStr = (RooStringVar*) buildConfig.find(physModel->GetName()) ; if (ruleStr) { strcpy(buf,ruleStr->getVal()) ; char *tokenPtr(0) ; #ifndef _WIN32 char* token = strtok_r(buf,spaceChars,&tokenPtr) ; #else char* token = strtok(buf,spaceChars) ; #endif enum Mode { SplitCat, Colon, ParamList } ; Mode mode(SplitCat) ; char* splitCatName ; RooAbsCategory* splitCat(0) ; while(token) { switch (mode) { case SplitCat: { splitCatName = token ; if (strchr(splitCatName,',')) { // Composite splitting category // Check if already instantiated splitCat = (RooAbsCategory*) _compSplitCatSet.find(splitCatName) ; TString origCompCatName(splitCatName) ; if (!splitCat) { // Build now #ifndef _WIN32 char *tokptr = 0; char *catName = strtok_r(token,",",&tokptr) ; #else char *catName = strtok(token,",") ; #endif RooArgSet compCatSet ; while(catName) { RooAbsArg* cat = splitCatSet.find(catName) ; // If not, check if it is an auxiliary splitcat if (!cat) { cat = (RooAbsCategory*) auxSplitSet.find(catName) ; } if (!cat) { cout << "RooSimPdfBuilder::buildPdf: ERROR " << catName << " not found in the primary or auxilary splitcat list" << endl ; customizerList->Delete() ; delete customizerList ; splitStateList.Delete() ; return 0 ; } compCatSet.add(*cat) ; #ifndef _WIN32 catName = strtok_r(0,",",&tokptr) ; #else catName = strtok(0,",") ; #endif } splitCat = new RooMultiCategory(origCompCatName,origCompCatName,compCatSet) ; _compSplitCatSet.addOwned(*splitCat) ; //cout << "composite splitcat: " << splitCat->GetName() ; } } else { // Simple splitting category // First see if it is a simple splitting category splitCat = (RooAbsCategory*) splitCatSet.find(splitCatName) ; // If not, check if it is an auxiliary splitcat if (!splitCat) { splitCat = (RooAbsCategory*) auxSplitSet.find(splitCatName) ; } if (!splitCat) { cout << "RooSimPdfBuilder::buildPdf: ERROR splitting category " << splitCatName << " not found in the primary or auxiliary splitcat list" << endl ; customizerList->Delete() ; delete customizerList ; splitStateList.Delete() ; return 0 ; } } mode = Colon ; break ; } case Colon: { if (strcmp(token,":")) { cout << "RooSimPdfBuilder::buildPdf: ERROR in parsing, expected ':' after " << splitCat << ", found " << token << endl ; customizerList->Delete() ; delete customizerList ; splitStateList.Delete() ; return 0 ; } mode = ParamList ; break ; } case ParamList: { // Verify the validity of the parameter list and build the corresponding argset RooArgSet splitParamList ; RooArgSet* paramList = physModel->getParameters(dependents) ; // wve -- add nodes to parameter list RooArgSet* compList = physModel->getComponents() ; paramList->add(*compList) ; delete compList ; Bool_t lastCharIsComma = (token[strlen(token)-1]==',') ; #ifndef _WIN32 char *tokptr = 0 ; char *paramName = strtok_r(token,",",&tokptr) ; #else char *tokptr(0) ; char *paramName = strtok(token,",") ; #endif // Check for fractional split option 'param_name[remainder_state]' char *remainderState = 0 ; #ifndef _WIN32 char *tokptr2 = 0 ; if (paramName && strtok_r(paramName,"[",&tokptr2)) { remainderState = strtok_r(0,"]",&tokptr2) ; } #else if (paramName && strtok(paramName,"[")) { remainderState = strtok(0,"]") ; } #endif while(paramName) { // If fractional split is specified, check that remainder state is a valid state of this split cat if (remainderState) { if (!splitCat->lookupType(remainderState)) { cout << "RooSimPdfBuilder::buildPdf: ERROR fraction split of parameter " << paramName << " has invalid remainder state name: " << remainderState << endl ; delete paramList ; customizerList->Delete() ; delete customizerList ; splitStateList.Delete() ; return 0 ; } } RooAbsArg* param = paramList->find(paramName) ; if (!param) { cout << "RooSimPdfBuilder::buildPdf: ERROR " << paramName << " is not a parameter of physics model " << physModel->GetName() << endl ; delete paramList ; customizerList->Delete() ; delete customizerList ; splitStateList.Delete() ; return 0 ; } splitParamList.add(*param) ; // Build split leaf of fraction splits here if (remainderState) { // Check if we are splitting a real-valued parameter if (!dynamic_cast<RooAbsReal*>(param)) { cout << "RooSimPdfBuilder::buildPdf: ERROR fraction split requested of non-real valued parameter " << param->GetName() << endl ; delete paramList ; customizerList->Delete() ; delete customizerList ; splitStateList.Delete() ; return 0 ; } // Check if we are doing a restricted build TList* remStateSplitList = static_cast<TList*>(splitStateList.FindObject(splitCat->GetName())) ; // If so, check if remainder state is actually being built. if (remStateSplitList && !remStateSplitList->FindObject(remainderState)) { cout << "RooSimPdfBuilder::buildPdf: ERROR " << paramName << " remainder state " << remainderState << " in parameter split " << param->GetName() << " is not actually being built" << endl ; delete paramList ; customizerList->Delete() ; delete customizerList ; splitStateList.Delete() ; return 0 ; } TIterator* iter = splitCat->typeIterator() ; RooCatType* type ; RooArgList fracLeafList ; TString formExpr("1") ; Int_t i(0) ; while((type=(RooCatType*)iter->Next())) { // Skip remainder state if (!TString(type->GetName()).CompareTo(remainderState)) continue ; // If restricted build is requested, skip states of splitcat that are not built if (remStateSplitList && !remStateSplitList->FindObject(type->GetName())) { continue ; } // Construct name of split leaf TString splitLeafName(param->GetName()) ; splitLeafName.Append("_") ; splitLeafName.Append(type->GetName()) ; // Check if split leaf already exists RooAbsArg* splitLeaf = _splitNodeList.find(splitLeafName) ; if (!splitLeaf) { // If not create it now splitLeaf = (RooAbsArg*) param->clone(splitLeafName) ; _splitNodeList.addOwned(*splitLeaf) ; } fracLeafList.add(*splitLeaf) ; formExpr.Append(Form("-@%d",i++)) ; } delete iter ; // Construct RooFormulaVar expresssing remainder of fraction TString remLeafName(param->GetName()) ; remLeafName.Append("_") ; remLeafName.Append(remainderState) ; // Check if no specialization was already specified for remainder state if (!_splitNodeList.find(remLeafName)) { RooAbsArg* remLeaf = new RooFormulaVar(remLeafName,formExpr,fracLeafList) ; _splitNodeList.addOwned(*remLeaf) ; cout << "RooSimPdfBuilder::buildPdf: creating remainder fraction formula for " << remainderState << " specialization of split parameter " << param->GetName() << " " << formExpr << endl ; } } // Parse next parameter name #ifndef _WIN32 paramName = strtok_r(0,",",&tokptr) ; if (paramName && strtok_r(paramName,"[",&tokptr2)) { remainderState = strtok_r(0,"]",&tokptr2) ; } #else paramName = strtok(0,",") ; if (paramName && strtok(paramName,"[")) { remainderState = strtok(0,"]") ; } #endif } // Add the rule to the appropriate customizer ; physCustomizer->splitArgs(splitParamList,*splitCat) ; delete paramList ; if (!lastCharIsComma) mode = SplitCat ; break ; } } #ifndef _WIN32 token = strtok_r(0,spaceChars,&tokenPtr) ; #else token = strtok(0,spaceChars) ; #endif } if (mode!=SplitCat) { cout << "RooSimPdfBuilder::buildPdf: ERROR in parsing, expected " << (mode==Colon?":":"parameter list") << " after " << token << endl ; } //RooArgSet* paramSet = physModel->getParameters(dependents) ; } else { cout << "RooSimPdfBuilder::buildPdf: no splitting rules for " << physModel->GetName() << endl ; } } cout << "RooSimPdfBuilder::buildPdf: configured customizers for all physics models" << endl ; customizerList->Print() ; // Create fit category from physCat and splitCatList ; RooArgSet fitCatList ; if (physCat) fitCatList.add(*physCat) ; fitCatList.add(splitCatSet) ; TIterator* fclIter = fitCatList.createIterator() ; RooSuperCategory *fitCat = new RooSuperCategory("fitCat","fitCat",fitCatList) ; // Create master PDF RooSimultaneous* simPdf = new RooSimultaneous("simPdf","simPdf",*fitCat) ; // Add component PDFs to master PDF TIterator* fcIter = fitCat->typeIterator() ; RooCatType* fcState ; while((fcState=(RooCatType*)fcIter->Next())) { // Select fitCat state fitCat->setLabel(fcState->GetName()) ; // Check if this fitCat state is selected fclIter->Reset() ; RooAbsCategory* splitCat ; Bool_t select(kTRUE) ; while((splitCat=(RooAbsCategory*)fclIter->Next())) { // Find selected state list TList* slist = (TList*) splitStateList.FindObject(splitCat->GetName()) ; if (!slist) continue ; RooCatType* type = (RooCatType*) slist->FindObject(splitCat->getLabel()) ; if (!type) { select = kFALSE ; } } if (!select) continue ; // Select appropriate PDF for this physCat state RooCustomizer* physCustomizer ; if (physCat) { RooStringVar* physNameVar = (RooStringVar*) stateMap.find(physCat->getLabel()) ; if (!physNameVar) continue ; physCustomizer = (RooCustomizer*) customizerList->FindObject(physNameVar->getVal()); } else { physCustomizer = (RooCustomizer*) customizerList->First() ; } cout << "RooSimPdfBuilder::buildPdf: Customizing physics model " << physCustomizer->GetName() << " for mode " << fcState->GetName() << endl ; // Customizer PDF for current state and add to master simPdf RooAbsPdf* fcPdf = (RooAbsPdf*) physCustomizer->build(masterSplitCat.getLabel(),verbose) ; simPdf->addPdf(*fcPdf,fcState->GetName()) ; } delete fcIter ; // Move customizers (owning the cloned branch node components) to the attic _retiredCustomizerList.AddAll(customizerList) ; delete customizerList ; delete fclIter ; splitStateList.Delete() ; if (auxSplitCloneSet) delete auxSplitCloneSet ; delete physIter ; return simPdf ; } RooSimPdfBuilder::~RooSimPdfBuilder() { _retiredCustomizerList.Delete() ; }