RooLinearVar.cxx

Go to the documentation of this file.
00001 /*****************************************************************************
00002  * Project: RooFit                                                           *
00003  * Package: RooFitCore                                                       *
00004  * @(#)root/roofitcore:$Id: RooLinearVar.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 // RooLinearVar is the most general form of a derived real-valued object that can
00020 // be used by RooRealIntegral to integrate over. The requirements for this are
00021 //
00022 //          - Can be modified directly (i.e. invertible formula)
00023 //          - Jacobian term in integral is constant (but not necessarily 1)
00024 //
00025 // This class implements the most general form that satisfy these requirement
00026 // 
00027 //    RLV = (slope)*x + (offset)
00028 //
00029 // X is required to be a RooRealVar to meet the invertibility criterium
00030 // (slope) and (offset) is are RooAbsReals, but may not overlap with x,
00031 // i.e. x may not be a server of (slope) and (offset)
00032 //
00033 // In the context of a dataset, (slope) may not contain any real-valued dependents
00034 // (satisfied constant Jacobian requirement). This check cannot be enforced at
00035 // construction time, but can be performed at run time through the isJacobianOK(depList)
00036 // member function.
00037 //
00038 //
00039 
00040 #include "RooFit.h"
00041 #include "Riostream.h"
00042 
00043 #include <math.h>
00044 #include <stdlib.h>
00045 #include <string.h>
00046 #include <ctype.h>
00047 #include "TClass.h"
00048 #include "TObjString.h"
00049 #include "TTree.h"
00050 #include "RooLinearVar.h"
00051 #include "RooStreamParser.h"
00052 #include "RooArgSet.h"
00053 #include "RooRealVar.h"
00054 #include "RooNumber.h"
00055 #include "RooBinning.h"
00056 #include "RooMsgService.h"
00057 
00058 
00059 
00060 ClassImp(RooLinearVar)
00061 
00062 
00063 //_____________________________________________________________________________
00064 RooLinearVar::RooLinearVar(const char *name, const char *title, RooAbsRealLValue& variable, 
00065                            const RooAbsReal& slope, const RooAbsReal& offset, const char *unit) :
00066   RooAbsRealLValue(name, title, unit), 
00067   _binning(variable.getBinning(),slope.getVal(),offset.getVal()),
00068   _var("var","variable",this,variable,kTRUE,kTRUE),
00069   _slope("slope","slope",this,(RooAbsReal&)slope),
00070   _offset("offset","offset",this,(RooAbsReal&)offset)
00071 {
00072   // Constructor with RooAbsRealLValue variable and RooAbsReal slope and offset
00073 
00074   // Slope and offset may not depend on variable
00075   if (slope.dependsOnValue(variable) || offset.dependsOnValue(variable)) {
00076     coutE(InputArguments) << "RooLinearVar::RooLinearVar(" << GetName() 
00077                           << "): ERROR, slope(" << slope.GetName() << ") and offset(" 
00078                           << offset.GetName() << ") may not depend on variable(" 
00079                           << variable.GetName() << ")" << endl ;
00080     assert(0) ;
00081   }
00082 
00083   // Initial plot range and number of bins from dependent variable
00084 //   setPlotRange(variable.getPlotMin()*_slope + _offset,
00085 //                variable.getPlotMax()*_slope + _offset) ;
00086 //   setPlotBins(variable.getPlotBins()) ;
00087                
00088 }  
00089 
00090 
00091 
00092 //_____________________________________________________________________________
00093 RooLinearVar::RooLinearVar(const RooLinearVar& other, const char* name) :
00094   RooAbsRealLValue(other,name), 
00095   _binning(other._binning),
00096   _var("var",this,other._var),
00097   _slope("slope",this,other._slope),
00098   _offset("offset",this,other._offset)
00099 {
00100   // Copy constructor
00101 }
00102 
00103 
00104 
00105 //_____________________________________________________________________________
00106 RooLinearVar::~RooLinearVar() 
00107 {
00108   // Destructor
00109 
00110   _altBinning.Delete() ;
00111 }
00112 
00113 
00114 
00115 //_____________________________________________________________________________
00116 Double_t RooLinearVar::evaluate() const
00117 {
00118   // Calculate current value of this object  
00119 
00120   return _offset + _var * _slope ;
00121 }
00122 
00123 
00124 
00125 //_____________________________________________________________________________
00126 void RooLinearVar::setVal(Double_t value) 
00127 {
00128   // Assign given value to linear transformation: sets input variable to (value-offset)/slope
00129   // If slope is zerom an error message is printed and no assignment is made
00130 
00131   //cout << "RooLinearVar::setVal(" << GetName() << "): new value = " << value << endl ;
00132 
00133   // Prevent DIV0 problems
00134   if (_slope == 0.) {
00135     coutE(Eval) << "RooLinearVar::setVal(" << GetName() << "): ERROR: slope is zero, cannot invert relation" << endl ;
00136     return ;
00137   }
00138 
00139   // Invert formula 'value = offset + slope*var'
00140   ((RooRealVar&)_var.arg()).setVal((value - _offset) / _slope) ;
00141 }
00142 
00143 
00144 
00145 //_____________________________________________________________________________
00146 Bool_t RooLinearVar::isJacobianOK(const RooArgSet& depList) const
00147 {
00148   // Returns true if Jacobian term associated with current
00149   // expression tree is indeed constant.
00150 
00151   if (!((RooAbsRealLValue&)_var.arg()).isJacobianOK(depList)) {
00152     return kFALSE ;
00153   }
00154 
00155   // Check if jacobian has no real-valued dependents
00156   RooAbsArg* arg ;
00157   TIterator* dIter = depList.createIterator() ;
00158   while ((arg=(RooAbsArg*)dIter->Next())) {
00159     if (arg->IsA()->InheritsFrom(RooAbsReal::Class())) {
00160       if (_slope.arg().dependsOnValue(*arg)) {
00161 //      cout << "RooLinearVar::isJacobianOK(" << GetName() << ") return kFALSE because slope depends on value of " << arg->GetName() << endl ;
00162         return kFALSE ;
00163       }
00164     }
00165   }
00166   delete dIter ;
00167 //   cout << "RooLinearVar::isJacobianOK(" << GetName() << ") return kTRUE" << endl ;
00168   return kTRUE ;
00169 }
00170 
00171 
00172 
00173 //_____________________________________________________________________________
00174 Double_t RooLinearVar::jacobian() const 
00175 {
00176   // Return value of Jacobian associated with the transformation
00177 
00178   return _slope*((RooAbsRealLValue&)_var.arg()).jacobian() ;
00179 }
00180 
00181 
00182 
00183 //_____________________________________________________________________________
00184 Bool_t RooLinearVar::readFromStream(istream& /*is*/, Bool_t /*compact*/, Bool_t /*verbose*/) 
00185 {
00186   // Read object contents from stream
00187   return kTRUE ;
00188 }
00189 
00190 
00191 
00192 //_____________________________________________________________________________
00193 void RooLinearVar::writeToStream(ostream& os, Bool_t compact) const
00194 {
00195   // Write object contents to stream
00196 
00197   if (compact) {
00198     os << getVal() ;
00199   } else {
00200     os << _slope.arg().GetName() << " * " << _var.arg().GetName() << " + " << _offset.arg().GetName() ;
00201   }
00202 }
00203 
00204 
00205 
00206 //_____________________________________________________________________________
00207  RooAbsBinning& RooLinearVar::getBinning(const char* name, Bool_t verbose, Bool_t createOnTheFly) 
00208 {
00209   // Retrieve binning of this linear transformation. A RooLinearVar does not have its own
00210   // binnings but uses linearly transformed binnings of teh input variable. If a given
00211   // binning exists on the input variable, it will also exists on this linear transformation
00212   // and a binning adaptor object is created on the fly.
00213 
00214   // Normalization binning
00215   if (name==0) {
00216     _binning.updateInput(((RooAbsRealLValue&)_var.arg()).getBinning(),_slope,_offset) ;
00217     return _binning ;
00218   } 
00219 
00220   // Alternative named range binnings, look for existing translator binning first
00221   RooLinTransBinning* altBinning = (RooLinTransBinning*) _altBinning.FindObject(name) ;
00222   if (altBinning) {
00223     altBinning->updateInput(((RooAbsRealLValue&)_var.arg()).getBinning(name,verbose),_slope,_offset) ;
00224     return *altBinning ;
00225   }
00226 
00227   // If binning is not found return default binning, if creation is not requested
00228   if (!createOnTheFly) {
00229     return _binning ;
00230   }
00231 
00232   // Create translator binning on the fly
00233   RooAbsBinning& sourceBinning = ((RooAbsRealLValue&)_var.arg()).getBinning(name,verbose) ;
00234   RooLinTransBinning* transBinning = new RooLinTransBinning(sourceBinning,_slope,_offset) ;
00235   _altBinning.Add(transBinning) ;
00236 
00237   return *transBinning ;
00238 }
00239 
00240 
00241 //_____________________________________________________________________________
00242 const RooAbsBinning& RooLinearVar::getBinning(const char* name, Bool_t verbose, Bool_t createOnTheFly) const
00243 {
00244   // Const version of getBinning()
00245 
00246   return const_cast<RooLinearVar*>(this)->getBinning(name,verbose,createOnTheFly) ;
00247 }
00248 
00249 
00250 //_____________________________________________________________________________
00251 Bool_t RooLinearVar::hasBinning(const char* name) const 
00252 {
00253   // Returns true if binning with given name exists.If a given binning
00254   // exists on the input variable, it will also exists on this linear
00255   // transformation.
00256 
00257   return ((RooAbsRealLValue&)_var.arg()).hasBinning(name) ;
00258 }

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