/*****************************************************************************
 * Project: RooFit                                                           *
 * Package: RooFitCore                                                       *
 *    File: $Id: RooAbsReal.h,v 1.75 2007/07/13 21:50:24 wouter 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)             *
 *****************************************************************************/
#ifndef ROO_ABS_REAL
#define ROO_ABS_REAL

#include "RooAbsArg.h"
#include "RooCmdArg.h"
#include "RooCurve.h"
#include "RooArgSet.h"
#include "RooArgList.h"
#include "RooGlobalFunc.h"

class RooArgList ;
class RooDataSet ;
class RooPlot;
class RooRealVar;
class RooAbsFunc;
class RooAbsCategoryLValue ;
class RooCategory ;
class RooLinkedList ;
class RooNumIntConfig ;
class RooDataHist ;
class RooFunctor ;
class RooGenFunction ;
class RooMultiGenFunction ;
class RooFitResult ;
class RooMoment ;
class RooDerivative ;
class RooVectorDataStore ;

class TH1;
class TH1F;
class TH2F;
class TH3F;

#include <list>
#include <string>
#include <iostream>

class RooAbsReal : public RooAbsArg {
public:
  // Constructors, assignment etc
  RooAbsReal() ;
  RooAbsReal(const char *name, const char *title, const char *unit= "") ;
  RooAbsReal(const char *name, const char *title, Double_t minVal, Double_t maxVal, 
	     const char *unit= "") ;
  RooAbsReal(const RooAbsReal& other, const char* name=0);
  virtual ~RooAbsReal();

  // Return value and unit accessors
  inline Double_t getVal(const RooArgSet* set=0) const { 
/*     if (_fast && !_inhibitDirty && std::string("RooHistFunc")==IsA()->GetName()) std::cout << "RooAbsReal::getVal(" << GetName() << ") CLEAN value = " << _value << std::endl ;  */
#ifndef _WIN32
    return (_fast && !_inhibitDirty) ? _value : getValV(set) ; 
#else
    return (_fast && !inhibitDirty()) ? _value : getValV(set) ;     
#endif
  }
  inline  Double_t getVal(const RooArgSet& set) const { return _fast ? _value : getValV(&set) ; }

  virtual Double_t getValV(const RooArgSet* set=0) const ;

  Double_t getPropagatedError(const RooFitResult& fr) ;

  Bool_t operator==(Double_t value) const ;
  virtual Bool_t operator==(const RooAbsArg& other) ;
  virtual Bool_t isIdentical(const RooAbsArg& other, Bool_t assumeSameType=kFALSE)  ;


  inline const Text_t *getUnit() const { 
    // Return string with unit description
    return _unit.Data(); 
  }
  inline void setUnit(const char *unit) { 
    // Set unit description to given string
    _unit= unit; 
  }
  TString getTitle(Bool_t appendUnit= kFALSE) const;

  // Lightweight interface adaptors (caller takes ownership)
  RooAbsFunc *bindVars(const RooArgSet &vars, const RooArgSet* nset=0, Bool_t clipInvalid=kFALSE) const;

  // Create a fundamental-type object that can hold our value.
  RooAbsArg *createFundamental(const char* newname=0) const;

  // Analytical integration support
  virtual Int_t getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet, const char* rangeName=0) const ;
  virtual Double_t analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName=0) const ;
  virtual Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=0) const ;
  virtual Double_t analyticalIntegral(Int_t code, const char* rangeName=0) const ;
  virtual Bool_t forceAnalyticalInt(const RooAbsArg& /*dep*/) const { 
    // Interface to force RooRealIntegral to offer given observable for internal integration
    // even if this is deemed unsafe. This default implementation returns always flase
    return kFALSE ; 
  }
  virtual void forceNumInt(Bool_t flag=kTRUE) { 
    // If flag is true, all advertised analytical integrals will be ignored
    // and all integrals are calculated numerically
    _forceNumInt = flag ; 
  }
  Bool_t getForceNumInt() const { return _forceNumInt ; }

  // Chi^2 fits to histograms
  virtual RooFitResult* chi2FitTo(RooDataHist& data, const RooCmdArg& arg1=RooCmdArg::none(),  const RooCmdArg& arg2=RooCmdArg::none(),  
                              const RooCmdArg& arg3=RooCmdArg::none(),  const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),  
                              const RooCmdArg& arg6=RooCmdArg::none(),  const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
  virtual RooFitResult* chi2FitTo(RooDataHist& data, const RooLinkedList& cmdList) ;

  virtual RooAbsReal* createChi2(RooDataHist& data, const RooLinkedList& cmdList) ;
  virtual RooAbsReal* createChi2(RooDataHist& data, const RooCmdArg& arg1=RooCmdArg::none(),  const RooCmdArg& arg2=RooCmdArg::none(),  
				 const RooCmdArg& arg3=RooCmdArg::none(),  const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),  
				 const RooCmdArg& arg6=RooCmdArg::none(),  const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;

  // Chi^2 fits to X-Y datasets
  virtual RooFitResult* chi2FitTo(RooDataSet& xydata, const RooCmdArg& arg1=RooCmdArg::none(),  const RooCmdArg& arg2=RooCmdArg::none(),  
                              const RooCmdArg& arg3=RooCmdArg::none(),  const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),  
                              const RooCmdArg& arg6=RooCmdArg::none(),  const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
  virtual RooFitResult* chi2FitTo(RooDataSet& xydata, const RooLinkedList& cmdList) ;

  virtual RooAbsReal* createChi2(RooDataSet& data, const RooLinkedList& cmdList) ;
  virtual RooAbsReal* createChi2(RooDataSet& data, const RooCmdArg& arg1=RooCmdArg::none(),  const RooCmdArg& arg2=RooCmdArg::none(),  
				   const RooCmdArg& arg3=RooCmdArg::none(),  const RooCmdArg& arg4=RooCmdArg::none(), const RooCmdArg& arg5=RooCmdArg::none(),  
				   const RooCmdArg& arg6=RooCmdArg::none(),  const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;


  virtual RooAbsReal* createProfile(const RooArgSet& paramsOfInterest) ;


  RooAbsReal* createIntegral(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2=RooCmdArg::none(),
                             const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), 
			     const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(), 
			     const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) const ;

  RooAbsReal* createIntegral(const RooArgSet& iset, const char* rangeName) const { 
    // Create integral over observables in iset in range named rangeName
    return createIntegral(iset,0,0,rangeName) ; 
  }
  RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet& nset, const char* rangeName=0) const { 
    // Create integral over observables in iset in range named rangeName with integrand normalized over observables in nset
    return createIntegral(iset,&nset,0,rangeName) ; 
  }
  RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet& nset, const RooNumIntConfig& cfg, const char* rangeName=0) const { 
    // Create integral over observables in iset in range named rangeName with integrand normalized over observables in nset while
    // using specified configuration for any numeric integration
    return createIntegral(iset,&nset,&cfg,rangeName) ; 
  }
  RooAbsReal* createIntegral(const RooArgSet& iset, const RooNumIntConfig& cfg, const char* rangeName=0) const { 
    // Create integral over observables in iset in range named rangeName using specified configuration for any numeric integration
    return createIntegral(iset,0,&cfg,rangeName) ; 
  }
  virtual RooAbsReal* createIntegral(const RooArgSet& iset, const RooArgSet* nset=0, const RooNumIntConfig* cfg=0, const char* rangeName=0) const ;  


  void setParameterizeIntegral(const RooArgSet& paramVars) ;

  // Create running integrals
  RooAbsReal* createRunningIntegral(const RooArgSet& iset, const RooArgSet& nset=RooArgSet()) ;
  RooAbsReal* createRunningIntegral(const RooArgSet& iset, const RooCmdArg& arg1, const RooCmdArg& arg2=RooCmdArg::none(),
			const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), 
			const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(), 
			const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) ;
  RooAbsReal* createIntRI(const RooArgSet& iset, const RooArgSet& nset=RooArgSet()) ;
  RooAbsReal* createScanRI(const RooArgSet& iset, const RooArgSet& nset, Int_t numScanBins, Int_t intOrder) ;

  
  // Optimized accept/reject generator support
  virtual Int_t getMaxVal(const RooArgSet& vars) const ;
  virtual Double_t maxVal(Int_t code) const ;
  virtual Int_t minTrialSamples(const RooArgSet& /*arGenObs*/) const { return 0 ; }


  // Plotting options
  void setPlotLabel(const char *label);
  const char *getPlotLabel() const;

  virtual Double_t defaultErrorLevel() const { 
    // Return default level for MINUIT error analysis
    return 1.0 ; 
  }

  const RooNumIntConfig* getIntegratorConfig() const ;
  RooNumIntConfig* getIntegratorConfig() ;
  static RooNumIntConfig* defaultIntegratorConfig()  ;
  RooNumIntConfig* specialIntegratorConfig() const ;
  RooNumIntConfig* specialIntegratorConfig(Bool_t createOnTheFly) ;
  void setIntegratorConfig() ;
  void setIntegratorConfig(const RooNumIntConfig& config) ;

  virtual void fixAddCoefNormalization(const RooArgSet& addNormSet=RooArgSet(),Bool_t force=kTRUE) ;
  virtual void fixAddCoefRange(const char* rangeName=0,Bool_t force=kTRUE) ;

  virtual void preferredObservableScanOrder(const RooArgSet& obs, RooArgSet& orderedObs) const ;

  // User entry point for plotting
  virtual RooPlot* plotOn(RooPlot* frame, 
			  const RooCmdArg& arg1=RooCmdArg(), const RooCmdArg& arg2=RooCmdArg(),
			  const RooCmdArg& arg3=RooCmdArg(), const RooCmdArg& arg4=RooCmdArg(),
			  const RooCmdArg& arg5=RooCmdArg(), const RooCmdArg& arg6=RooCmdArg(),
			  const RooCmdArg& arg7=RooCmdArg(), const RooCmdArg& arg8=RooCmdArg(),
			  const RooCmdArg& arg9=RooCmdArg(), const RooCmdArg& arg10=RooCmdArg()
              ) const ;


  enum ScaleType { Raw, Relative, NumEvent, RelativeExpected } ;

  // Forwarder function for backward compatibility
  virtual RooPlot *plotSliceOn(RooPlot *frame, const RooArgSet& sliceSet, Option_t* drawOptions="L", 
			       Double_t scaleFactor=1.0, ScaleType stype=Relative, const RooAbsData* projData=0) const;

  // Fill an existing histogram
  TH1 *fillHistogram(TH1 *hist, const RooArgList &plotVars,
		     Double_t scaleFactor= 1, const RooArgSet *projectedVars= 0, Bool_t scaling=kTRUE,
		     const RooArgSet* condObs=0, Bool_t setError=kTRUE) const;

  // Create 1,2, and 3D histograms from and fill it
  TH1 *createHistogram(const char* varNameList, Int_t xbins=0, Int_t ybins=0, Int_t zbins=0) const ;
  TH1* createHistogram(const char *name, const RooAbsRealLValue& xvar, RooLinkedList& argList) const ;
  TH1 *createHistogram(const char *name, const RooAbsRealLValue& xvar,
                       const RooCmdArg& arg1=RooCmdArg::none(), const RooCmdArg& arg2=RooCmdArg::none(), 
                       const RooCmdArg& arg3=RooCmdArg::none(), const RooCmdArg& arg4=RooCmdArg::none(), 
                       const RooCmdArg& arg5=RooCmdArg::none(), const RooCmdArg& arg6=RooCmdArg::none(), 
                       const RooCmdArg& arg7=RooCmdArg::none(), const RooCmdArg& arg8=RooCmdArg::none()) const ;

  // Fill a RooDataHist
  RooDataHist* fillDataHist(RooDataHist *hist, const RooArgSet* nset, Double_t scaleFactor,
			    Bool_t correctForBinVolume=kFALSE, Bool_t showProgress=kFALSE) const ;

  // I/O streaming interface (machine readable)
  virtual Bool_t readFromStream(std::istream& is, Bool_t compact, Bool_t verbose=kFALSE) ;
  virtual void writeToStream(std::ostream& os, Bool_t compact) const ;

  // Printing interface (human readable)
  virtual void printValue(std::ostream& os) const ;
  virtual void printMultiline(std::ostream& os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const ;

  static void setCacheCheck(Bool_t flag) ;

  // Evaluation error logging 
  class EvalError {
  public:
    EvalError() { }
    EvalError(const EvalError& other) : _msg(other._msg), _srvval(other._srvval) { }
    void setMessage(const char* tmp) { std::string s(tmp); s.swap(_msg); }
    void setServerValues(const char* tmp) { std::string s(tmp); s.swap(_srvval); }
    std::string _msg;
    std::string _srvval;
  } ;

  enum ErrorLoggingMode { PrintErrors, CollectErrors, CountErrors, Ignore } ;
  static ErrorLoggingMode evalErrorLoggingMode() ;
  static void setEvalErrorLoggingMode(ErrorLoggingMode m) ;
  void logEvalError(const char* message, const char* serverValueString=0) const ;
  static void logEvalError(const RooAbsReal* originator, const char* origName, const char* message, const char* serverValueString=0) ;
  static void printEvalErrors(std::ostream&os=std::cout, Int_t maxPerNode=10000000) ;
  static Int_t numEvalErrors() ;
  static Int_t numEvalErrorItems() ;

   
  typedef std::map<const RooAbsArg*,std::pair<std::string,std::list<EvalError> > >::const_iterator EvalErrorIter ; 
  static EvalErrorIter evalErrorIter() ;

  static void clearEvalErrorLog() ;
  
  virtual Bool_t isBinnedDistribution(const RooArgSet& /*obs*/) const { return kFALSE ; }
  virtual std::list<Double_t>* binBoundaries(RooAbsRealLValue& /*obs*/, Double_t /*xlo*/, Double_t /*xhi*/) const { return 0 ; }
  virtual std::list<Double_t>* plotSamplingHint(RooAbsRealLValue& /*obs*/, Double_t /*xlo*/, Double_t /*xhi*/) const { 
    // Interface for returning an optional hint for initial sampling points when constructing a curve 
    // projected on observable.
    return 0 ; 
  }

  RooGenFunction* iGenFunction(RooRealVar& x, const RooArgSet& nset=RooArgSet()) ;
  RooMultiGenFunction* iGenFunction(const RooArgSet& observables, const RooArgSet& nset=RooArgSet()) ;

  RooFunctor* functor(const RooArgList& obs, const RooArgList& pars=RooArgList(), const RooArgSet& nset=RooArgSet()) const ;
  TF1* asTF(const RooArgList& obs, const RooArgList& pars=RooArgList(), const RooArgSet& nset=RooArgSet()) const ;

  RooDerivative* derivative(RooRealVar& obs, Int_t order=1, Double_t eps=0.001) ;
  RooDerivative* derivative(RooRealVar& obs, const RooArgSet& normSet, Int_t order, Double_t eps=0.001) ; 

  RooMoment* moment(RooRealVar& obs, Int_t order, Bool_t central, Bool_t takeRoot) ;
  RooMoment* moment(RooRealVar& obs, const RooArgSet& normObs, Int_t order, Bool_t central, Bool_t takeRoot, Bool_t intNormObs) ;

  RooMoment* mean(RooRealVar& obs) { return moment(obs,1,kFALSE,kFALSE) ; }
  RooMoment* mean(RooRealVar& obs, const RooArgSet& nset) { return moment(obs,nset,1,kFALSE,kFALSE,kTRUE) ; }
  RooMoment* sigma(RooRealVar& obs) { return moment(obs,2,kTRUE,kTRUE) ; }
  RooMoment* sigma(RooRealVar& obs, const RooArgSet& nset) { return moment(obs,nset,2,kTRUE,kTRUE,kTRUE) ; }

  Double_t findRoot(RooRealVar& x, Double_t xmin, Double_t xmax, Double_t yval) ;


  virtual Bool_t setData(RooAbsData& /*data*/, Bool_t /*cloneData*/=kTRUE) { return kTRUE ; }

  virtual void enableOffsetting(Bool_t) {} ;
  virtual Bool_t isOffsetting() const { return kFALSE ; }
  virtual Double_t offset() const { return 0 ; }
  
  static void setHideOffset(Bool_t flag);
  static Bool_t hideOffset() ;

protected:
  // Hook for objects with normalization-dependent parameters interperetation
  virtual void selectNormalization(const RooArgSet* depSet=0, Bool_t force=kFALSE) ;
  virtual void selectNormalizationRange(const char* rangeName=0, Bool_t force=kFALSE) ;

  // Helper functions for plotting
  Bool_t plotSanityChecks(RooPlot* frame) const ;
  void makeProjectionSet(const RooAbsArg* plotVar, const RooArgSet* allVars, 
			 RooArgSet& projectedVars, Bool_t silent) const ;

  TString integralNameSuffix(const RooArgSet& iset, const RooArgSet* nset=0, const char* rangeName=0, Bool_t omitEmpty=kFALSE) const ;


  Bool_t isSelectedComp() const ;

  
 public:
  const RooAbsReal* createPlotProjection(const RooArgSet& depVars, const RooArgSet& projVars) const ;
  const RooAbsReal* createPlotProjection(const RooArgSet& depVars, const RooArgSet& projVars, RooArgSet*& cloneSet) const ;
  const RooAbsReal *createPlotProjection(const RooArgSet &dependentVars, const RooArgSet *projectedVars,
				         RooArgSet *&cloneSet, const char* rangeName=0, const RooArgSet* condObs=0) const;
 protected:

  RooFitResult* chi2FitDriver(RooAbsReal& fcn, RooLinkedList& cmdList) ;

  RooPlot* plotOnWithErrorBand(RooPlot* frame,const RooFitResult& fr, Double_t Z, const RooArgSet* params, const RooLinkedList& argList, Bool_t method1) const ;

  // Support interface for subclasses to advertise their analytic integration
  // and generator capabilities in their analticalIntegral() and generateEvent()
  // implementations.
  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, 
		   const RooArgProxy& a) const ;
  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, 
		   const RooArgProxy& a, const RooArgProxy& b) const ;
  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, 
		   const RooArgProxy& a, const RooArgProxy& b, const RooArgProxy& c) const ;
  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, 
		   const RooArgProxy& a, const RooArgProxy& b, 		   
		   const RooArgProxy& c, const RooArgProxy& d) const ;

  Bool_t matchArgs(const RooArgSet& allDeps, RooArgSet& numDeps, 
		   const RooArgSet& set) const ;


  RooAbsReal* createIntObj(const RooArgSet& iset, const RooArgSet* nset, const RooNumIntConfig* cfg, const char* rangeName) const ;
  void findInnerMostIntegration(const RooArgSet& allObs, RooArgSet& innerObs, const char* rangeName) const ;


  // Internal consistency checking (needed by RooDataSet)
  virtual Bool_t isValid() const ;
  virtual Bool_t isValidReal(Double_t value, Bool_t printError=kFALSE) const ;

  // Function evaluation and error tracing
  Double_t traceEval(const RooArgSet* set) const ;
  virtual Bool_t traceEvalHook(Double_t /*value*/) const { 
    // Hook function to add functionality to evaluation tracing in derived classes
    return kFALSE ;
  }
  virtual Double_t evaluate() const = 0 ;

  // Hooks for RooDataSet interface
  friend class RooRealIntegral ;
  friend class RooVectorDataStore ;
  virtual void syncCache(const RooArgSet* set=0) { getVal(set) ; }
  virtual void copyCache(const RooAbsArg* source, Bool_t valueOnly=kFALSE, Bool_t setValDirty=kTRUE) ;
  virtual void attachToTree(TTree& t, Int_t bufSize=32000) ;
  virtual void attachToVStore(RooVectorDataStore& vstore) ;
  virtual void setTreeBranchStatus(TTree& t, Bool_t active) ;
  virtual void fillTreeBranch(TTree& t) ;

  friend class RooRealBinding ;
  Double_t _plotMin ;       // Minimum of plot range
  Double_t _plotMax ;       // Maximum of plot range
  Int_t    _plotBins ;      // Number of plot bins
  mutable Double_t _value ; // Cache for current value of object
  TString  _unit ;          // Unit for objects value
  TString  _label ;         // Plot label for objects value
  Bool_t   _forceNumInt ;   // Force numerical integration if flag set

  mutable Float_t _floatValue ; //! Transient cache for floating point values from tree branches 
  mutable Int_t   _intValue   ; //! Transient cache for integer values from tree branches 
  mutable Bool_t  _boolValue  ; //! Transient cache for bool values from tree branches 
  mutable UChar_t _byteValue  ; //! Transient cache for byte values from tree branches 
  mutable Char_t  _sbyteValue ; //! Transient cache for signed byte values from tree branches 
  mutable UInt_t  _uintValue  ; //! Transient cache for unsigned integer values from tree branches 

  friend class RooAbsPdf ;
  friend class RooAbsAnaConvPdf ;
  friend class RooRealProxy ;

  RooNumIntConfig* _specIntegratorConfig ; // Numeric integrator configuration specific for this object

  Bool_t   _treeVar ;       // !do not persist

  static Bool_t _cacheCheck ; // If true, always validate contents of clean which outcome of evaluate()

  friend class RooDataProjBinding ;
  friend class RooAbsOptGoodnessOfFit ;
  
  struct PlotOpt {
   PlotOpt() : drawOptions("L"), scaleFactor(1.0), stype(Relative), projData(0), binProjData(kFALSE), projSet(0), precision(1e-3), 
               shiftToZero(kFALSE),projDataSet(0),normRangeName(0),rangeLo(0),rangeHi(0),postRangeFracScale(kFALSE),wmode(RooCurve::Extended),
               projectionRangeName(0),curveInvisible(kFALSE), curveName(0),addToCurveName(0),addToWgtSelf(1.),addToWgtOther(1.),
               numCPU(1),interleave(RooFit::Interleave),curveNameSuffix(""), numee(10), eeval(0), doeeval(kFALSE), progress(kFALSE) {} ;
   Option_t* drawOptions ;
   Double_t scaleFactor ;	 
   ScaleType stype ;
   const RooAbsData* projData ;
   Bool_t binProjData ;
   const RooArgSet* projSet ;
   Double_t precision ;
   Bool_t shiftToZero ;
   const RooArgSet* projDataSet ;
   const char* normRangeName ;
   Double_t rangeLo ;
   Double_t rangeHi ;
   Bool_t postRangeFracScale ;
   RooCurve::WingMode wmode ;
   const char* projectionRangeName ;
   Bool_t curveInvisible ;
   const char* curveName ;
   const char* addToCurveName ;
   Double_t addToWgtSelf ;
   Double_t addToWgtOther ;
   Int_t    numCPU ;
   RooFit::MPSplit interleave ;
   const char* curveNameSuffix ; 
   Int_t    numee ;
   Double_t eeval ;
   Bool_t   doeeval ;
   Bool_t progress ;
  } ;

  // Plot implementation functions
  virtual RooPlot *plotOn(RooPlot* frame, PlotOpt o) const;

public:
  // PlotOn with command list
  virtual RooPlot* plotOn(RooPlot* frame, RooLinkedList& cmdList) const ;

 protected:
  virtual RooPlot *plotAsymOn(RooPlot *frame, const RooAbsCategoryLValue& asymCat, PlotOpt o) const;


private:

  static ErrorLoggingMode _evalErrorMode ;
  static std::map<const RooAbsArg*,std::pair<std::string,std::list<EvalError> > > _evalErrorList ;
  static Int_t _evalErrorCount ;

  Bool_t matchArgsByName(const RooArgSet &allArgs, RooArgSet &matchedArgs, const TList &nameList) const;

protected:


  friend class RooRealSumPdf ;
  friend class RooAddPdf ;
  friend class RooAddModel ;
  void selectComp(Bool_t flag) { 
    // If flag is true, only selected component will be included in evaluates of RooAddPdf components
    _selectComp = flag ; 
  }
  static void globalSelectComp(Bool_t flag) ;
  Bool_t _selectComp ;               //! Component selection flag for RooAbsPdf::plotCompOn
  static Bool_t _globalSelectComp ;  // Global activation switch for component selection

  mutable RooArgSet* _lastNSet ; //!
  static Bool_t _hideOffset ; // Offset hiding flag

  ClassDef(RooAbsReal,2) // Abstract real-valued variable
};

#endif
 RooAbsReal.h:1
 RooAbsReal.h:2
 RooAbsReal.h:3
 RooAbsReal.h:4
 RooAbsReal.h:5
 RooAbsReal.h:6
 RooAbsReal.h:7
 RooAbsReal.h:8
 RooAbsReal.h:9
 RooAbsReal.h:10
 RooAbsReal.h:11
 RooAbsReal.h:12
 RooAbsReal.h:13
 RooAbsReal.h:14
 RooAbsReal.h:15
 RooAbsReal.h:16
 RooAbsReal.h:17
 RooAbsReal.h:18
 RooAbsReal.h:19
 RooAbsReal.h:20
 RooAbsReal.h:21
 RooAbsReal.h:22
 RooAbsReal.h:23
 RooAbsReal.h:24
 RooAbsReal.h:25
 RooAbsReal.h:26
 RooAbsReal.h:27
 RooAbsReal.h:28
 RooAbsReal.h:29
 RooAbsReal.h:30
 RooAbsReal.h:31
 RooAbsReal.h:32
 RooAbsReal.h:33
 RooAbsReal.h:34
 RooAbsReal.h:35
 RooAbsReal.h:36
 RooAbsReal.h:37
 RooAbsReal.h:38
 RooAbsReal.h:39
 RooAbsReal.h:40
 RooAbsReal.h:41
 RooAbsReal.h:42
 RooAbsReal.h:43
 RooAbsReal.h:44
 RooAbsReal.h:45
 RooAbsReal.h:46
 RooAbsReal.h:47
 RooAbsReal.h:48
 RooAbsReal.h:49
 RooAbsReal.h:50
 RooAbsReal.h:51
 RooAbsReal.h:52
 RooAbsReal.h:53
 RooAbsReal.h:54
 RooAbsReal.h:55
 RooAbsReal.h:56
 RooAbsReal.h:57
 RooAbsReal.h:58
 RooAbsReal.h:59
 RooAbsReal.h:60
 RooAbsReal.h:61
 RooAbsReal.h:62
 RooAbsReal.h:63
 RooAbsReal.h:64
 RooAbsReal.h:65
 RooAbsReal.h:66
 RooAbsReal.h:67
 RooAbsReal.h:68
 RooAbsReal.h:69
 RooAbsReal.h:70
 RooAbsReal.h:71
 RooAbsReal.h:72
 RooAbsReal.h:73
 RooAbsReal.h:74
 RooAbsReal.h:75
 RooAbsReal.h:76
 RooAbsReal.h:77
 RooAbsReal.h:78
 RooAbsReal.h:79
 RooAbsReal.h:80
 RooAbsReal.h:81
 RooAbsReal.h:82
 RooAbsReal.h:83
 RooAbsReal.h:84
 RooAbsReal.h:85
 RooAbsReal.h:86
 RooAbsReal.h:87
 RooAbsReal.h:88
 RooAbsReal.h:89
 RooAbsReal.h:90
 RooAbsReal.h:91
 RooAbsReal.h:92
 RooAbsReal.h:93
 RooAbsReal.h:94
 RooAbsReal.h:95
 RooAbsReal.h:96
 RooAbsReal.h:97
 RooAbsReal.h:98
 RooAbsReal.h:99
 RooAbsReal.h:100
 RooAbsReal.h:101
 RooAbsReal.h:102
 RooAbsReal.h:103
 RooAbsReal.h:104
 RooAbsReal.h:105
 RooAbsReal.h:106
 RooAbsReal.h:107
 RooAbsReal.h:108
 RooAbsReal.h:109
 RooAbsReal.h:110
 RooAbsReal.h:111
 RooAbsReal.h:112
 RooAbsReal.h:113
 RooAbsReal.h:114
 RooAbsReal.h:115
 RooAbsReal.h:116
 RooAbsReal.h:117
 RooAbsReal.h:118
 RooAbsReal.h:119
 RooAbsReal.h:120
 RooAbsReal.h:121
 RooAbsReal.h:122
 RooAbsReal.h:123
 RooAbsReal.h:124
 RooAbsReal.h:125
 RooAbsReal.h:126
 RooAbsReal.h:127
 RooAbsReal.h:128
 RooAbsReal.h:129
 RooAbsReal.h:130
 RooAbsReal.h:131
 RooAbsReal.h:132
 RooAbsReal.h:133
 RooAbsReal.h:134
 RooAbsReal.h:135
 RooAbsReal.h:136
 RooAbsReal.h:137
 RooAbsReal.h:138
 RooAbsReal.h:139
 RooAbsReal.h:140
 RooAbsReal.h:141
 RooAbsReal.h:142
 RooAbsReal.h:143
 RooAbsReal.h:144
 RooAbsReal.h:145
 RooAbsReal.h:146
 RooAbsReal.h:147
 RooAbsReal.h:148
 RooAbsReal.h:149
 RooAbsReal.h:150
 RooAbsReal.h:151
 RooAbsReal.h:152
 RooAbsReal.h:153
 RooAbsReal.h:154
 RooAbsReal.h:155
 RooAbsReal.h:156
 RooAbsReal.h:157
 RooAbsReal.h:158
 RooAbsReal.h:159
 RooAbsReal.h:160
 RooAbsReal.h:161
 RooAbsReal.h:162
 RooAbsReal.h:163
 RooAbsReal.h:164
 RooAbsReal.h:165
 RooAbsReal.h:166
 RooAbsReal.h:167
 RooAbsReal.h:168
 RooAbsReal.h:169
 RooAbsReal.h:170
 RooAbsReal.h:171
 RooAbsReal.h:172
 RooAbsReal.h:173
 RooAbsReal.h:174
 RooAbsReal.h:175
 RooAbsReal.h:176
 RooAbsReal.h:177
 RooAbsReal.h:178
 RooAbsReal.h:179
 RooAbsReal.h:180
 RooAbsReal.h:181
 RooAbsReal.h:182
 RooAbsReal.h:183
 RooAbsReal.h:184
 RooAbsReal.h:185
 RooAbsReal.h:186
 RooAbsReal.h:187
 RooAbsReal.h:188
 RooAbsReal.h:189
 RooAbsReal.h:190
 RooAbsReal.h:191
 RooAbsReal.h:192
 RooAbsReal.h:193
 RooAbsReal.h:194
 RooAbsReal.h:195
 RooAbsReal.h:196
 RooAbsReal.h:197
 RooAbsReal.h:198
 RooAbsReal.h:199
 RooAbsReal.h:200
 RooAbsReal.h:201
 RooAbsReal.h:202
 RooAbsReal.h:203
 RooAbsReal.h:204
 RooAbsReal.h:205
 RooAbsReal.h:206
 RooAbsReal.h:207
 RooAbsReal.h:208
 RooAbsReal.h:209
 RooAbsReal.h:210
 RooAbsReal.h:211
 RooAbsReal.h:212
 RooAbsReal.h:213
 RooAbsReal.h:214
 RooAbsReal.h:215
 RooAbsReal.h:216
 RooAbsReal.h:217
 RooAbsReal.h:218
 RooAbsReal.h:219
 RooAbsReal.h:220
 RooAbsReal.h:221
 RooAbsReal.h:222
 RooAbsReal.h:223
 RooAbsReal.h:224
 RooAbsReal.h:225
 RooAbsReal.h:226
 RooAbsReal.h:227
 RooAbsReal.h:228
 RooAbsReal.h:229
 RooAbsReal.h:230
 RooAbsReal.h:231
 RooAbsReal.h:232
 RooAbsReal.h:233
 RooAbsReal.h:234
 RooAbsReal.h:235
 RooAbsReal.h:236
 RooAbsReal.h:237
 RooAbsReal.h:238
 RooAbsReal.h:239
 RooAbsReal.h:240
 RooAbsReal.h:241
 RooAbsReal.h:242
 RooAbsReal.h:243
 RooAbsReal.h:244
 RooAbsReal.h:245
 RooAbsReal.h:246
 RooAbsReal.h:247
 RooAbsReal.h:248
 RooAbsReal.h:249
 RooAbsReal.h:250
 RooAbsReal.h:251
 RooAbsReal.h:252
 RooAbsReal.h:253
 RooAbsReal.h:254
 RooAbsReal.h:255
 RooAbsReal.h:256
 RooAbsReal.h:257
 RooAbsReal.h:258
 RooAbsReal.h:259
 RooAbsReal.h:260
 RooAbsReal.h:261
 RooAbsReal.h:262
 RooAbsReal.h:263
 RooAbsReal.h:264
 RooAbsReal.h:265
 RooAbsReal.h:266
 RooAbsReal.h:267
 RooAbsReal.h:268
 RooAbsReal.h:269
 RooAbsReal.h:270
 RooAbsReal.h:271
 RooAbsReal.h:272
 RooAbsReal.h:273
 RooAbsReal.h:274
 RooAbsReal.h:275
 RooAbsReal.h:276
 RooAbsReal.h:277
 RooAbsReal.h:278
 RooAbsReal.h:279
 RooAbsReal.h:280
 RooAbsReal.h:281
 RooAbsReal.h:282
 RooAbsReal.h:283
 RooAbsReal.h:284
 RooAbsReal.h:285
 RooAbsReal.h:286
 RooAbsReal.h:287
 RooAbsReal.h:288
 RooAbsReal.h:289
 RooAbsReal.h:290
 RooAbsReal.h:291
 RooAbsReal.h:292
 RooAbsReal.h:293
 RooAbsReal.h:294
 RooAbsReal.h:295
 RooAbsReal.h:296
 RooAbsReal.h:297
 RooAbsReal.h:298
 RooAbsReal.h:299
 RooAbsReal.h:300
 RooAbsReal.h:301
 RooAbsReal.h:302
 RooAbsReal.h:303
 RooAbsReal.h:304
 RooAbsReal.h:305
 RooAbsReal.h:306
 RooAbsReal.h:307
 RooAbsReal.h:308
 RooAbsReal.h:309
 RooAbsReal.h:310
 RooAbsReal.h:311
 RooAbsReal.h:312
 RooAbsReal.h:313
 RooAbsReal.h:314
 RooAbsReal.h:315
 RooAbsReal.h:316
 RooAbsReal.h:317
 RooAbsReal.h:318
 RooAbsReal.h:319
 RooAbsReal.h:320
 RooAbsReal.h:321
 RooAbsReal.h:322
 RooAbsReal.h:323
 RooAbsReal.h:324
 RooAbsReal.h:325
 RooAbsReal.h:326
 RooAbsReal.h:327
 RooAbsReal.h:328
 RooAbsReal.h:329
 RooAbsReal.h:330
 RooAbsReal.h:331
 RooAbsReal.h:332
 RooAbsReal.h:333
 RooAbsReal.h:334
 RooAbsReal.h:335
 RooAbsReal.h:336
 RooAbsReal.h:337
 RooAbsReal.h:338
 RooAbsReal.h:339
 RooAbsReal.h:340
 RooAbsReal.h:341
 RooAbsReal.h:342
 RooAbsReal.h:343
 RooAbsReal.h:344
 RooAbsReal.h:345
 RooAbsReal.h:346
 RooAbsReal.h:347
 RooAbsReal.h:348
 RooAbsReal.h:349
 RooAbsReal.h:350
 RooAbsReal.h:351
 RooAbsReal.h:352
 RooAbsReal.h:353
 RooAbsReal.h:354
 RooAbsReal.h:355
 RooAbsReal.h:356
 RooAbsReal.h:357
 RooAbsReal.h:358
 RooAbsReal.h:359
 RooAbsReal.h:360
 RooAbsReal.h:361
 RooAbsReal.h:362
 RooAbsReal.h:363
 RooAbsReal.h:364
 RooAbsReal.h:365
 RooAbsReal.h:366
 RooAbsReal.h:367
 RooAbsReal.h:368
 RooAbsReal.h:369
 RooAbsReal.h:370
 RooAbsReal.h:371
 RooAbsReal.h:372
 RooAbsReal.h:373
 RooAbsReal.h:374
 RooAbsReal.h:375
 RooAbsReal.h:376
 RooAbsReal.h:377
 RooAbsReal.h:378
 RooAbsReal.h:379
 RooAbsReal.h:380
 RooAbsReal.h:381
 RooAbsReal.h:382
 RooAbsReal.h:383
 RooAbsReal.h:384
 RooAbsReal.h:385
 RooAbsReal.h:386
 RooAbsReal.h:387
 RooAbsReal.h:388
 RooAbsReal.h:389
 RooAbsReal.h:390
 RooAbsReal.h:391
 RooAbsReal.h:392
 RooAbsReal.h:393
 RooAbsReal.h:394
 RooAbsReal.h:395
 RooAbsReal.h:396
 RooAbsReal.h:397
 RooAbsReal.h:398
 RooAbsReal.h:399
 RooAbsReal.h:400
 RooAbsReal.h:401
 RooAbsReal.h:402
 RooAbsReal.h:403
 RooAbsReal.h:404
 RooAbsReal.h:405
 RooAbsReal.h:406
 RooAbsReal.h:407
 RooAbsReal.h:408
 RooAbsReal.h:409
 RooAbsReal.h:410
 RooAbsReal.h:411
 RooAbsReal.h:412
 RooAbsReal.h:413
 RooAbsReal.h:414
 RooAbsReal.h:415
 RooAbsReal.h:416
 RooAbsReal.h:417
 RooAbsReal.h:418
 RooAbsReal.h:419
 RooAbsReal.h:420
 RooAbsReal.h:421
 RooAbsReal.h:422
 RooAbsReal.h:423
 RooAbsReal.h:424
 RooAbsReal.h:425
 RooAbsReal.h:426
 RooAbsReal.h:427
 RooAbsReal.h:428
 RooAbsReal.h:429
 RooAbsReal.h:430
 RooAbsReal.h:431
 RooAbsReal.h:432
 RooAbsReal.h:433
 RooAbsReal.h:434
 RooAbsReal.h:435
 RooAbsReal.h:436
 RooAbsReal.h:437
 RooAbsReal.h:438
 RooAbsReal.h:439
 RooAbsReal.h:440
 RooAbsReal.h:441
 RooAbsReal.h:442
 RooAbsReal.h:443
 RooAbsReal.h:444
 RooAbsReal.h:445
 RooAbsReal.h:446
 RooAbsReal.h:447
 RooAbsReal.h:448
 RooAbsReal.h:449
 RooAbsReal.h:450
 RooAbsReal.h:451
 RooAbsReal.h:452
 RooAbsReal.h:453
 RooAbsReal.h:454
 RooAbsReal.h:455
 RooAbsReal.h:456
 RooAbsReal.h:457
 RooAbsReal.h:458
 RooAbsReal.h:459
 RooAbsReal.h:460
 RooAbsReal.h:461
 RooAbsReal.h:462
 RooAbsReal.h:463
 RooAbsReal.h:464
 RooAbsReal.h:465
 RooAbsReal.h:466
 RooAbsReal.h:467
 RooAbsReal.h:468
 RooAbsReal.h:469
 RooAbsReal.h:470
 RooAbsReal.h:471
 RooAbsReal.h:472
 RooAbsReal.h:473
 RooAbsReal.h:474
 RooAbsReal.h:475
 RooAbsReal.h:476
 RooAbsReal.h:477
 RooAbsReal.h:478
 RooAbsReal.h:479
 RooAbsReal.h:480
 RooAbsReal.h:481
 RooAbsReal.h:482
 RooAbsReal.h:483
 RooAbsReal.h:484
 RooAbsReal.h:485
 RooAbsReal.h:486