RooResolutionModel.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitCore                                                       *
00004  * @(#)root/roofitcore:$Id: RooResolutionModel.cxx 30333 2009-09-21 15:39:17Z 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 //  RooResolutionModel is the base class of for PDFs that represent a
00020 //  resolution model that can be convoluted with physics a physics model of the form
00021 //
00022 //    Phys(x,a,b) = Sum_k coef_k(a) * basis_k(x,b)
00023 //  
00024 //  where basis_k are a limited number of functions in terms of the variable
00025 //  to be convoluted and coef_k are coefficients independent of the convolution
00026 //  variable.
00027 //  
00028 //  Classes derived from RooResolutionModel implement 
00029 //         _ _                        _                  _
00030 //   R_k(x,b,c) = Int(dx') basis_k(x',b) * resModel(x-x',c)
00031 // 
00032 //  which RooAbsAnaConvPdf uses to construct the pdf for [ Phys (x) R ] :
00033 //          _ _ _                 _          _ _
00034 //    PDF(x,a,b,c) = Sum_k coef_k(a) * R_k(x,b,c)
00035 //
00036 //  A minimal implementation of a RooResolutionModel consists of a
00037 //
00038 //    Int_t basisCode(const char* name)   
00039 //
00040 //  function indication which basis functions this resolution model supports, and
00041 //
00042 //    Double_t evaluate() 
00043 //
00044 //  Implementing the resolution model, optionally convoluted with one of the
00045 //  supported basis functions. RooResolutionModel objects can be used as regular
00046 //  PDFs (They inherit from RooAbsPdf), or as resolution model convoluted with
00047 //  a basis function. The implementation of evaluate() can identify the requested
00048 //  from of use from the basisCode() function. If zero, the regular PDF value
00049 //  should be calculated. If non-zero, the models value convoluted with the
00050 //  basis function identified by the code should be calculated.
00051 //
00052 //  Optionally, analytical integrals can be advertised and implemented, in the
00053 //  same way as done for regular PDFS (see RooAbsPdf for further details).
00054 //  Also in getAnalyticalIntegral()/analyticalIntegral() the implementation
00055 //  should use basisCode() to determine for which scenario the integral is
00056 //  requested.
00057 //
00058 //  The choice of basis returned by basisCode() is guaranteed not to change
00059 //  of the lifetime of a RooResolutionModel object.
00060 //
00061 
00062 #include "RooFit.h"
00063 
00064 #include "TClass.h"
00065 #include "TMath.h"
00066 #include "Riostream.h"
00067 #include "RooResolutionModel.h"
00068 #include "RooMsgService.h"
00069 #include "RooSentinel.h"
00070 
00071 ClassImp(RooResolutionModel) 
00072 ;
00073 
00074 RooFormulaVar* RooResolutionModel::_identity = 0;
00075 
00076 
00077 
00078 //_____________________________________________________________________________
00079 void RooResolutionModel::cleanup()
00080 {
00081   // Cleanup hook for RooSentinel atexit handler
00082   delete _identity ;
00083   _identity = 0 ;
00084 }
00085 
00086 
00087 //_____________________________________________________________________________
00088 RooResolutionModel::RooResolutionModel(const char *name, const char *title, RooRealVar& _x) : 
00089   RooAbsPdf(name,title), 
00090   x("x","Dependent or convolution variable",this,_x),
00091   _basisCode(0), _basis(0), 
00092   _ownBasis(kFALSE)
00093 {
00094   // Constructor with convolution variable 'x'
00095 
00096   if (!_identity) {
00097     _identity = identity() ; 
00098   }
00099 }
00100 
00101 
00102 
00103 //_____________________________________________________________________________
00104 RooResolutionModel::RooResolutionModel(const RooResolutionModel& other, const char* name) : 
00105   RooAbsPdf(other,name), 
00106   x("x",this,other.x),
00107   _basisCode(other._basisCode), _basis(0),
00108   _ownBasis(kFALSE)
00109 {
00110   // Copy constructor
00111 
00112   if (other._basis) {
00113     _basis = (RooFormulaVar*) other._basis->Clone() ;
00114     _ownBasis = kTRUE ;
00115     //_basis = other._basis ;
00116   }
00117 
00118   if (_basis) {
00119     TIterator* bsIter = _basis->serverIterator() ;
00120     RooAbsArg* basisServer ;
00121     while((basisServer = (RooAbsArg*)bsIter->Next())) {
00122       addServer(*basisServer,kTRUE,kFALSE) ;
00123     }
00124     delete bsIter ;
00125   }
00126 }
00127 
00128 
00129 
00130 //_____________________________________________________________________________
00131 RooResolutionModel::~RooResolutionModel()
00132 {
00133   // Destructor
00134 
00135   if (_ownBasis && _basis) {
00136     delete _basis ;
00137   }
00138 }
00139 
00140 
00141 
00142 //_____________________________________________________________________________
00143 RooFormulaVar* RooResolutionModel::identity() 
00144 { 
00145   // Return identity formula pointer
00146 
00147   if (!_identity) {
00148     _identity = new RooFormulaVar("identity","1",RooArgSet("")) ;  
00149     RooSentinel::activate() ;
00150   }
00151 
00152   return _identity ; 
00153 }
00154 
00155 
00156 
00157 //_____________________________________________________________________________
00158 RooResolutionModel* RooResolutionModel::convolution(RooFormulaVar* inBasis, RooAbsArg* owner) const
00159 {
00160   // Instantiate a clone of this resolution model representing a convolution with given
00161   // basis function. The owners object name is incorporated in the clones name
00162   // to avoid multiple convolution objects with the same name in complex PDF structures.
00163   // 
00164   // Note: The 'inBasis' formula expression must be a RooFormulaVar that encodes the formula
00165   // in the title of the object and this expression must be an exact match against the
00166   // implemented basis function strings (see derived class implementation of method basisCode()
00167   // for those strings
00168 
00169   // Check that primary variable of basis functions is our convolution variable  
00170   if (inBasis->getParameter(0) != x.absArg()) {
00171     coutE(InputArguments) << "RooResolutionModel::convolution(" << GetName() << "," << this  
00172                           << ") convolution parameter of basis function and PDF don't match" << endl 
00173                           << "basis->findServer(0) = " << inBasis->findServer(0) << endl 
00174                           << "x.absArg()           = " << x.absArg() << endl ;
00175     return 0 ;
00176   }
00177   
00178   if (basisCode(inBasis->GetTitle())==0) {
00179     coutE(InputArguments) << "RooResolutionModel::convolution(" << GetName() << "," << this  
00180                           << ") basis function '" << inBasis->GetTitle() << "' is not supported." << endl ;
00181     return 0 ;
00182   }
00183 
00184   TString newName(GetName()) ;
00185   newName.Append("_conv_") ;
00186   newName.Append(inBasis->GetName()) ;
00187   newName.Append("_[") ;
00188   newName.Append(owner->GetName()) ;
00189   newName.Append("]") ;
00190 
00191   RooResolutionModel* conv = (RooResolutionModel*) clone(newName) ;
00192   
00193   TString newTitle(conv->GetTitle()) ;
00194   newTitle.Append(" convoluted with basis function ") ;
00195   newTitle.Append(inBasis->GetName()) ;
00196   conv->SetTitle(newTitle.Data()) ;
00197 
00198   conv->changeBasis(inBasis) ;
00199 
00200   return conv ;
00201 }
00202 
00203 
00204 
00205 //_____________________________________________________________________________
00206 void RooResolutionModel::changeBasis(RooFormulaVar* inBasis) 
00207 {
00208   // Change the basis function we convolute with.
00209   // For one-time use by convolution() only.
00210 
00211   // Remove client-server link to old basis
00212   if (_basis) {
00213     TIterator* bsIter = _basis->serverIterator() ;
00214     RooAbsArg* basisServer ;
00215     while((basisServer = (RooAbsArg*)bsIter->Next())) {
00216       removeServer(*basisServer) ;
00217     }
00218     delete bsIter ;
00219 
00220     if (_ownBasis) {
00221       delete _basis ;
00222     }
00223   }
00224   _ownBasis = kFALSE ;
00225 
00226   // Change basis pointer and update client-server link
00227   _basis = inBasis ;
00228   if (_basis) {
00229     TIterator* bsIter = _basis->serverIterator() ;
00230     RooAbsArg* basisServer ;
00231     while((basisServer = (RooAbsArg*)bsIter->Next())) {
00232       addServer(*basisServer,kTRUE,kFALSE) ;
00233     }
00234     delete bsIter ;
00235   }
00236 
00237   _basisCode = inBasis?basisCode(inBasis->GetTitle()):0 ;
00238 }
00239 
00240 
00241 
00242 //_____________________________________________________________________________
00243 const RooRealVar& RooResolutionModel::basisConvVar() const 
00244 {
00245   // Return the convolution variable of the selection basis function.
00246   // This is, by definition, the first parameter of the basis function
00247 
00248   // Convolution variable is by definition first server of basis function
00249   TIterator* sIter = basis().serverIterator() ;
00250   RooRealVar* var = (RooRealVar*) sIter->Next() ;
00251   delete sIter ;
00252 
00253   return *var ;
00254 }
00255 
00256 
00257 
00258 //_____________________________________________________________________________
00259 RooRealVar& RooResolutionModel::convVar() const 
00260 {
00261   // Return the convolution variable of the resolution model
00262 
00263   return (RooRealVar&) x.arg() ;
00264 }
00265 
00266 
00267 
00268 //_____________________________________________________________________________
00269 Double_t RooResolutionModel::getVal(const RooArgSet* nset) const
00270 {
00271   // Modified version of RooAbsPdf::getVal(). If used as regular PDF, 
00272   // call RooAbsPdf::getVal(), otherwise return unnormalized value
00273   // regardless of specified normalization set
00274 
00275   if (!_basis) return RooAbsPdf::getVal(nset) ;
00276 
00277   // Return value of object. Calculated if dirty, otherwise cached value is returned.
00278   if (isValueDirty()) {
00279     _value = evaluate() ; 
00280 
00281     // WVE insert traceEval traceEval
00282     if (_verboseDirty) cxcoutD(Tracing) << "RooResolutionModel(" << GetName() << ") value = " << _value << endl ;
00283 
00284     clearValueDirty() ; 
00285     clearShapeDirty() ; 
00286   }
00287 
00288   return _value ;
00289 }
00290 
00291 
00292 
00293 //_____________________________________________________________________________
00294 Bool_t RooResolutionModel::redirectServersHook(const RooAbsCollection& newServerList, Bool_t mustReplaceAll, Bool_t nameChange, Bool_t /*isRecursive*/) 
00295 {
00296   // Forward redirectServers call to our basis function, which is not connected to either resolution
00297   // model or the physics model.
00298 
00299   if (!_basis) {
00300     _norm = 0 ;
00301     return kFALSE ; 
00302   } 
00303 
00304   RooFormulaVar* newBasis = (RooFormulaVar*) newServerList.find(_basis->GetName()) ;
00305   if (newBasis) {
00306 
00307     if (_ownBasis) {
00308       delete _basis ;
00309     }
00310 
00311     _basis = newBasis ;
00312     _ownBasis = kFALSE ;
00313   }
00314 
00315   _basis->redirectServers(newServerList,mustReplaceAll,nameChange) ;
00316     
00317   return (mustReplaceAll && !newBasis) ;
00318 }
00319 
00320 
00321 
00322 //_____________________________________________________________________________
00323 Bool_t RooResolutionModel::traceEvalHook(Double_t value) const 
00324 {
00325   // Floating point error checking and tracing for given float value
00326 
00327   // check for a math error or negative value
00328   return isnan(value) ;
00329 }
00330 
00331 
00332 
00333 //_____________________________________________________________________________
00334 void RooResolutionModel::normLeafServerList(RooArgSet& list) const 
00335 {
00336   // Return the list of servers used by our normalization integral
00337 
00338   _norm->leafNodeServerList(&list) ;
00339 }
00340 
00341 
00342 
00343 //_____________________________________________________________________________
00344 Double_t RooResolutionModel::getNorm(const RooArgSet* nset) const
00345 {
00346   // Return the integral of this PDF over all elements of 'nset'. 
00347 
00348   if (!nset) {
00349     return getVal() ;
00350   }
00351 
00352   syncNormalization(nset,kFALSE) ;
00353   if (_verboseEval>1) cxcoutD(Tracing) << IsA()->GetName() << "::getNorm(" << GetName() 
00354                                        << "): norm(" << _norm << ") = " << _norm->getVal() << endl ;
00355 
00356   Double_t ret = _norm->getVal() ;
00357   return ret ;
00358 }
00359 
00360 
00361 
00362 //_____________________________________________________________________________
00363 void RooResolutionModel::printMultiline(ostream& os, Int_t content, Bool_t verbose, TString indent) const
00364 {
00365   // Print info about this object to the specified stream. In addition to the info
00366   // from RooAbsArg::printStream() we add:
00367   //
00368   //     Shape : value, units, plot range
00369   //   Verbose : default binning and print label
00370 
00371   RooAbsPdf::printMultiline(os,content,verbose,indent) ;
00372 
00373   if(verbose) {
00374     os << indent << "--- RooResolutionModel ---" << endl;
00375     os << indent << "basis function = " ; 
00376     if (_basis) {
00377       _basis->printStream(os,kName|kAddress|kTitle,kSingleLine,indent) ;
00378     } else {
00379       os << "<none>" << endl ;
00380     }
00381   }
00382 }
00383 

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