#include "RooFit.h"
#include "RooMsgService.h"
#include "TIterator.h"
#include "TIterator.h"
#include "TList.h"
#include "RooAddModel.h"
#include "RooDataSet.h"
#include "RooRealProxy.h"
#include "RooPlot.h"
#include "RooRealVar.h"
#include "RooAddGenContext.h"
#include "RooRealConstant.h"
#include "RooNameReg.h"
#include "RooMsgService.h"
#include "Riostream.h"
ClassImp(RooAddModel)
;
RooAddModel::RooAddModel() :
  _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,kFALSE,kFALSE),
  _refCoefRangeName(0),
  _codeReg(10),
  _snormList(0)
{
  _pdfIter   = _pdfList.createIterator() ;
  _coefIter  = _coefList.createIterator() ;
  _coefCache = new Double_t[10] ;
  _coefErrCount = _errorCount ;
}
RooAddModel::RooAddModel(const char *name, const char *title, const RooArgList& inPdfList, const RooArgList& inCoefList, Bool_t ownPdfList) :
  RooResolutionModel(name,title,((RooResolutionModel*)inPdfList.at(0))->convVar()),
  _refCoefNorm("!refCoefNorm","Reference coefficient normalization set",this,kFALSE,kFALSE),
  _refCoefRangeName(0),
  _projectCoefs(kFALSE),
  _projCacheMgr(this,10),
  _intCacheMgr(this,10),
  _codeReg(10),
  _pdfList("pdfs","List of PDFs",this),
  _coefList("coefficients","List of coefficients",this),
  _haveLastCoef(kFALSE),
  _allExtendable(kFALSE)
{ 
  
  
  
  
  
  
  if (inPdfList.getSize()>inCoefList.getSize()+1) {
    coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName() 
			  << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1" << endl ;
    assert(0) ;
  }
  _pdfIter  = _pdfList.createIterator() ;
  _coefIter = _coefList.createIterator() ;
 
  
  TIterator* pdfIter = inPdfList.createIterator() ;
  TIterator* coefIter = inCoefList.createIterator() ;
  RooAbsPdf* pdf ;
  RooAbsReal* coef ;
  while((coef = (RooAbsPdf*)coefIter->Next())) {
    pdf = (RooAbsPdf*) pdfIter->Next() ;
    if (!pdf) {
      coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName() 
			    << ") number of pdfs and coefficients inconsistent, must have Npdf=Ncoef or Npdf=Ncoef+1" << endl ;
      assert(0) ;
    }
    if (!dynamic_cast<RooAbsReal*>(coef)) {
      coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName() << ") coefficient " << coef->GetName() << " is not of type RooAbsReal, ignored" << endl ;
      continue ;
    }
    if (!dynamic_cast<RooAbsReal*>(pdf)) {
      coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName() << ") pdf " << pdf->GetName() << " is not of type RooAbsPdf, ignored" << endl ;
      continue ;
    }
    _pdfList.add(*pdf) ;
    _coefList.add(*coef) ;    
  }
  pdf = (RooAbsPdf*) pdfIter->Next() ;
  if (pdf) {
    if (!dynamic_cast<RooAbsReal*>(pdf)) {
      coutE(InputArguments) << "RooAddModel::RooAddModel(" << GetName() << ") last pdf " << coef->GetName() << " is not of type RooAbsPdf, fatal error" << endl ;
      assert(0) ;
    }
    _pdfList.add(*pdf) ;  
  } else {
    _haveLastCoef=kTRUE ;
  }
  delete pdfIter ;
  delete coefIter  ;
  _coefCache = new Double_t[_pdfList.getSize()] ;
  _coefErrCount = _errorCount ;
  if (ownPdfList) {
    _ownedComps.addOwned(_pdfList) ;
  }
}
RooAddModel::RooAddModel(const RooAddModel& other, const char* name) :
  RooResolutionModel(other,name),
  _refCoefNorm("!refCoefNorm",this,other._refCoefNorm),
  _refCoefRangeName((TNamed*)other._refCoefRangeName),
  _projectCoefs(other._projectCoefs),
  _projCacheMgr(other._projCacheMgr,this),
  _intCacheMgr(other._intCacheMgr,this),
  _codeReg(other._codeReg),
  _pdfList("pdfs",this,other._pdfList),
  _coefList("coefficients",this,other._coefList),
  _haveLastCoef(other._haveLastCoef),
  _allExtendable(other._allExtendable)
{
  
  _pdfIter  = _pdfList.createIterator() ;
  _coefIter = _coefList.createIterator() ;
  _coefCache = new Double_t[_pdfList.getSize()] ;
  _coefErrCount = _errorCount ;
}
RooAddModel::~RooAddModel()
{
  
  delete _pdfIter ;
  delete _coefIter ;
  if (_coefCache) delete[] _coefCache ;
}
void RooAddModel::fixCoefNormalization(const RooArgSet& refCoefNorm) 
{
  if (refCoefNorm.getSize()==0) {
    _projectCoefs = kFALSE ;
    return ;
  }
  _projectCoefs = kTRUE ;  
  _refCoefNorm.removeAll() ;
  _refCoefNorm.add(refCoefNorm) ;
  _projCacheMgr.reset() ;
}
void RooAddModel::fixCoefRange(const char* rangeName)
{
  _refCoefRangeName = (TNamed*)RooNameReg::ptr(rangeName) ;
  if (_refCoefRangeName) _projectCoefs = kTRUE ;
}
RooResolutionModel* RooAddModel::convolution(RooFormulaVar* inBasis, RooAbsArg* owner) const
{
  
  
  
  
  
  
  if (inBasis->findServer(0) != x.absArg()) {
    coutE(InputArguments) << "RooAddModel::convolution(" << GetName() 
			  << ") convolution parameter of basis function and PDF don't match" << endl ;
    ccoutE(InputArguments) << "basis->findServer(0) = " << inBasis->findServer(0) << " " << inBasis->findServer(0)->GetName() << endl ;
    ccoutE(InputArguments) << "x.absArg()           = " << x.absArg() << " " << x.absArg()->GetName() << endl ;
    inBasis->Print("v") ;
    return 0 ;
  }
  TString newName(GetName()) ;
  newName.Append("_conv_") ;
  newName.Append(inBasis->GetName()) ;
  newName.Append("_[") ;
  newName.Append(owner->GetName()) ;
  newName.Append("]") ;
  TString newTitle(GetTitle()) ;
  newTitle.Append(" convoluted with basis function ") ;
  newTitle.Append(inBasis->GetName()) ;
  _pdfIter->Reset() ;
  RooResolutionModel* model ;
  RooArgList modelList ;
  while((model = (RooResolutionModel*)_pdfIter->Next())) {       
    
    RooResolutionModel* conv = model->convolution(inBasis,owner) ;    
    modelList.add(*conv) ;
  }
  _coefIter->Reset() ;
  RooAbsReal* coef ;
  RooArgList theCoefList ;  
  while((coef = (RooAbsReal*)_coefIter->Next())) {
    theCoefList.add(*coef) ;
  }
    
  RooAddModel* convSum = new RooAddModel(newName,newTitle,modelList,theCoefList,kTRUE) ;
  convSum->changeBasis(inBasis) ;
  return convSum ;
}
Int_t RooAddModel::basisCode(const char* name) const 
{
  
  
  
  
  TIterator* mIter = _pdfList.createIterator() ;
  RooResolutionModel* model ;
  Bool_t first(kTRUE), code(0) ;
    while((model = (RooResolutionModel*)mIter->Next())) {
      Int_t subCode = model->basisCode(name) ;
      if (first) {
	code = subCode ;
      } else if (subCode==0) {
	code = 0 ;
      }
  }
  delete mIter ;
  return code ;
}
RooAddModel::CacheElem* RooAddModel::getProjCache(const RooArgSet* nset, const RooArgSet* iset, const char* rangeName) const
{
  
  CacheElem* cache = (CacheElem*) _projCacheMgr.getObj(nset,iset,0,RooNameReg::ptr(rangeName)) ;
  if (cache) {
    return cache ;
  }
  
  cache = new CacheElem ;
  
  
  RooArgSet *fullDepList = getObservables(nset) ;
  if (iset) {
    fullDepList->remove(*iset,kTRUE,kTRUE) ;
  }    
  
  _pdfIter->Reset() ;
  _coefIter->Reset() ;
  RooAbsPdf* pdf ;
  RooAbsReal* coef ;
  while((pdf=(RooAbsPdf*)_pdfIter->Next())) {    
    coef=(RooAbsPdf*)_coefIter->Next() ;
    
    RooArgSet supNSet(*fullDepList) ;
    
    RooArgSet* pdfDeps = pdf->getObservables(nset) ;
    if (pdfDeps) {
      supNSet.remove(*pdfDeps,kTRUE,kTRUE) ;
      delete pdfDeps ; 
    }
    
    RooArgSet* coefDeps = coef ? coef->getObservables(nset) : 0 ;
    if (coefDeps) {
      supNSet.remove(*coefDeps,kTRUE,kTRUE) ;
      delete coefDeps ;
    }
    
    RooAbsReal* snorm ;
    TString name(GetName()) ;
    name.Append("_") ;
    name.Append(pdf->GetName()) ;
    name.Append("_SupNorm") ;
    if (supNSet.getSize()>0) {
      snorm = new RooRealIntegral(name,"Supplemental normalization integral",RooRealConstant::value(1.0),supNSet) ;
    } else {
      snorm = new RooRealVar(name,"Unit Supplemental normalization integral",1.0) ;
    }
    cache->_suppNormList.addOwned(*snorm) ;
  }
  delete fullDepList ;
    
  if (_verboseEval>1) {
    cxcoutD(Caching) << "RooAddModel::syncSuppNormList(" << GetName() << ") synching supplemental normalization list for norm" << (nset?*nset:RooArgSet()) << endl ;
    if (dologD(Caching)) {
      cache->_suppNormList.Print("v") ;
    }
  }
  
  
  if (!_projectCoefs || _basis!=0 ) {
    _projCacheMgr.setObj(nset,iset,cache,RooNameReg::ptr(rangeName)) ;
    return cache ;
  }
  
  RooArgSet* nset2 = nset ? getObservables(nset) : new RooArgSet() ;
  
  if (!nset2->equals(_refCoefNorm) || _refCoefRangeName !=0 || rangeName !=0 ) {
    
    coutI(Caching)  << "RooAddModel::syncCoefProjList(" << GetName() << ") creating coefficient projection integrals" << endl ;
    ccoutI(Caching) << "  from current normalization: "  ; nset2->Print("1") ;
    ccoutI(Caching) << "          with current range: " << (rangeName?rangeName:"<none>") << endl ;
    ccoutI(Caching) << "  to reference normalization: "  ; _refCoefNorm.Print("1") ; 
    ccoutI(Caching) << "        with reference range: " << (_refCoefRangeName?RooNameReg::str(_refCoefRangeName):"<none>") << endl ;
    
    
    _pdfIter->Reset() ;
    RooAbsPdf* thePdf ;
    while((thePdf=(RooAbsPdf*)_pdfIter->Next())) {
      
      RooAbsReal* pdfProj ;
      if (!nset2->equals(_refCoefNorm)) {
	pdfProj = thePdf->createIntegral(*nset2,_refCoefNorm) ;
	pdfProj->setOperMode(operMode()) ;
      } else {
	TString name(GetName()) ;
	name.Append("_") ;
	name.Append(thePdf->GetName()) ;
	name.Append("_ProjectNorm") ;
	pdfProj = new RooRealVar(name,"Unit Projection normalization integral",1.0) ;
      }
      cache->_projList.addOwned(*pdfProj) ;
      
      RooArgSet supNormSet(_refCoefNorm) ;
      RooArgSet* deps = thePdf->getParameters(RooArgSet()) ;
      supNormSet.remove(*deps,kTRUE,kTRUE) ;
      delete deps ;
      RooAbsReal* snorm ;
      TString name(GetName()) ;
      name.Append("_") ;
      name.Append(thePdf->GetName()) ;
      name.Append("_ProjSupNorm") ;
      if (supNormSet.getSize()>0) {
	snorm = new RooRealIntegral(name,"Projection Supplemental normalization integral",
				    RooRealConstant::value(1.0),supNormSet) ;
      } else {
	snorm = new RooRealVar(name,"Unit Projection Supplemental normalization integral",1.0) ;
      }
      cache->_suppProjList.addOwned(*snorm) ;
      
      RooAbsReal* rangeProj1 ;
      if (_refCoefRangeName && _refCoefNorm.getSize()>0) {
	rangeProj1 = thePdf->createIntegral(_refCoefNorm,_refCoefNorm,RooNameReg::str(_refCoefRangeName)) ;
	rangeProj1->setOperMode(operMode()) ;
      } else {
	TString theName(GetName()) ;
	theName.Append("_") ;
	theName.Append(thePdf->GetName()) ;
	theName.Append("_RangeNorm1") ;
	rangeProj1 = new RooRealVar(theName,"Unit range normalization integral",1.0) ;
      }
      cache->_refRangeProjList.addOwned(*rangeProj1) ;
      
      
      RooAbsReal* rangeProj2 ;
      if (rangeName && _refCoefNorm.getSize()>0) {
	rangeProj2 = thePdf->createIntegral(_refCoefNorm,_refCoefNorm,rangeName) ;
	rangeProj2->setOperMode(operMode()) ;
      } else {
	TString theName(GetName()) ;
	theName.Append("_") ;
	theName.Append(thePdf->GetName()) ;
	theName.Append("_RangeNorm2") ;
	rangeProj2 = new RooRealVar(theName,"Unit range normalization integral",1.0) ;
      }
      cache->_rangeProjList.addOwned(*rangeProj2) ;
    }               
  }
  delete nset2 ;
  _projCacheMgr.setObj(nset,iset,cache,RooNameReg::ptr(rangeName)) ;
  return cache ;
}
void RooAddModel::updateCoefficients(CacheElem& cache, const RooArgSet* nset) const 
{
  
  
  Int_t i ;
  
  if (_allExtendable) {
    
    
    Double_t coefSum(0) ;
    for (i=0 ; i<_pdfList.getSize() ; i++) {
      _coefCache[i] = ((RooAbsPdf*)_pdfList.at(i))->expectedEvents(_refCoefNorm.getSize()>0?&_refCoefNorm:nset) ;
      coefSum += _coefCache[i] ;
    }
    if (coefSum==0.) {
      coutW(Eval) << "RooAddModel::updateCoefCache(" << GetName() << ") WARNING: total number of expected events is 0" << endl ;
    } else {
      for (i=0 ; i<_pdfList.getSize() ; i++) {
	_coefCache[i] /= coefSum ;
      }			            
    }
    
  } else {
    if (_haveLastCoef) {
      
      
      Double_t coefSum(0) ;
      for (i=0 ; i<_coefList.getSize() ; i++) {
	_coefCache[i] = ((RooAbsPdf*)_coefList.at(i))->getVal(nset) ;
	coefSum += _coefCache[i] ;
      }		
      for (i=0 ; i<_coefList.getSize() ; i++) {
	_coefCache[i] /= coefSum ;
      }			
    } else {
      
      
      Double_t lastCoef(1) ;
      for (i=0 ; i<_coefList.getSize() ; i++) {
	_coefCache[i] = ((RooAbsPdf*)_coefList.at(i))->getVal(nset) ;
 	cxcoutD(Caching) << "SYNC: orig coef[" << i << "] = " << _coefCache[i] << endl ;
	lastCoef -= _coefCache[i] ;
      }			
      _coefCache[_coefList.getSize()] = lastCoef ;
      cxcoutD(Caching) << "SYNC: orig coef[" << _coefList.getSize() << "] = " << _coefCache[_coefList.getSize()] << endl ;
      
      
      
      if ((lastCoef<-1e-05 || (lastCoef-1)>1e-5) && _coefErrCount-->0) {
	coutW(Eval) << "RooAddModel::updateCoefCache(" << GetName() 
		    << " WARNING: sum of PDF coefficients not in range [0-1], value=" 
		    << 1-lastCoef << endl ;
	if (_coefErrCount==0) {
	  coutW(Eval) << " (no more will be printed)" << endl  ;
	}
      } 
    }
  }
  
  
  if ((!_projectCoefs) || cache._projList.getSize()==0) {
    
    return ;
  }
  
  Double_t coefSum(0) ;
  for (i=0 ; i<_pdfList.getSize() ; i++) {
    RooAbsPdf::globalSelectComp(kTRUE) ;    
    RooAbsReal* pp = ((RooAbsReal*)cache._projList.at(i)) ; 
    RooAbsReal* sn = ((RooAbsReal*)cache._suppProjList.at(i)) ; 
    RooAbsReal* r1 = ((RooAbsReal*)cache._refRangeProjList.at(i)) ;
    RooAbsReal* r2 = ((RooAbsReal*)cache._rangeProjList.at(i)) ;
    
    if (dologD(Eval)) {
      cxcoutD(Eval) << "pp = " << pp->GetName() << endl 
		    << "sn = " << sn->GetName() << endl 
		    << "r1 = " << r1->GetName() << endl 
		    << "r2 = " << r2->GetName() << endl ;
      r1->printStream(ccoutD(Eval),kName|kArgs|kValue,kSingleLine) ;
      r1->printCompactTree(ccoutD(Eval)) ;
    }
    Double_t proj = pp->getVal()/sn->getVal()*(r2->getVal()/r1->getVal()) ;  
    
    RooAbsPdf::globalSelectComp(kFALSE) ;
    _coefCache[i] *= proj ;
    coefSum += _coefCache[i] ;
  }
  for (i=0 ; i<_pdfList.getSize() ; i++) {
    _coefCache[i] /= coefSum ;
  }
   
  
}
Double_t RooAddModel::evaluate() const 
{
  
  const RooArgSet* nset = _normSet ; 
  CacheElem* cache = getProjCache(nset) ;
  updateCoefficients(*cache,nset) ;
  
  
  _pdfIter->Reset() ;
  _coefIter->Reset() ;
  RooAbsPdf* pdf ;
  Double_t snormVal ;
  Double_t value(0) ;
  Int_t i(0) ;
  while((pdf = (RooAbsPdf*)_pdfIter->Next())) {
    if (_coefCache[i]!=0.) {
      snormVal = nset ? ((RooAbsReal*)cache->_suppNormList.at(i))->getVal() : 1.0 ;
      Double_t pdfVal = pdf->getVal(nset) ;
      
      if (pdf->isSelectedComp()) {
	value += pdfVal*_coefCache[i]/snormVal ;
 	cxcoutD(Eval) << "RooAddModel::evaluate(" << GetName() << ")  value += [" 
 			<< pdf->GetName() << "] " << pdfVal << " * " << _coefCache[i] << " / " << snormVal << endl ;
      }
    }
    i++ ;
  }
  return value ;
}
void RooAddModel::resetErrorCounters(Int_t resetValue)
{
  
  
  RooAbsPdf::resetErrorCounters(resetValue) ;
  _coefErrCount = resetValue ;
}
Bool_t RooAddModel::checkObservables(const RooArgSet* nset) const 
{
  
  
  
  Bool_t ret(kFALSE) ;
  _pdfIter->Reset() ;
  _coefIter->Reset() ;
  RooAbsReal* coef ;
  RooAbsReal* pdf ;
  while((coef=(RooAbsReal*)_coefIter->Next())) {
    pdf = (RooAbsReal*)_pdfIter->Next() ;
    if (pdf->observableOverlaps(nset,*coef)) {
      coutE(InputArguments) << "RooAddModel::checkObservables(" << GetName() << "): ERROR: coefficient " << coef->GetName() 
			    << " and PDF " << pdf->GetName() << " have one or more dependents in common" << endl ;
      ret = kTRUE ;
    }
  }
  
  return ret ;
}
Int_t RooAddModel::getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, 
					 const RooArgSet* normSet, const char* rangeName) const 
{
  if (_forceNumInt) return 0 ;
  
  analVars.add(allVars) ;
  
  Int_t code ;
  RooArgList *cilist ;
  getCompIntList(normSet,&allVars,cilist,code,rangeName) ;
  
  return code+1 ;
  
}
void RooAddModel::getCompIntList(const RooArgSet* nset, const RooArgSet* iset, pRooArgList& compIntList, Int_t& code, const char* isetRangeName) const 
{
  
  Int_t sterileIdx(-1) ;
  IntCacheElem* cache = (IntCacheElem*) _intCacheMgr.getObj(nset,iset,&sterileIdx,RooNameReg::ptr(isetRangeName)) ;
  if (cache) {
    code = _intCacheMgr.lastIndex() ;
    compIntList = &cache->_intList ;
    
    return ;
  }
  
  cache = new IntCacheElem ;
  
  _pdfIter->Reset() ;
  RooResolutionModel* model ;
  while ((model=(RooResolutionModel*)_pdfIter->Next())) {
    RooAbsReal* intPdf = model->createIntegral(*iset,nset,0,isetRangeName) ;
    cache->_intList.addOwned(*intPdf) ;
  }
  
  code = _intCacheMgr.setObj(nset,iset,(RooAbsCacheElement*)cache,RooNameReg::ptr(isetRangeName)) ;
  
  compIntList = &cache->_intList ;
}
Double_t RooAddModel::analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const 
{
  
  
  if (code==0) {
    return getVal(normSet) ;
  }
  
  IntCacheElem* cache = (IntCacheElem*) _intCacheMgr.getObjByIndex(code-1) ;
  
  RooArgList* compIntList ;
  
  if (cache==0) {
    RooArgSet* vars = getParameters(RooArgSet()) ;
    RooArgSet* nset = _intCacheMgr.nameSet1ByIndex(code-1)->select(*vars) ;
    RooArgSet* iset = _intCacheMgr.nameSet2ByIndex(code-1)->select(*vars) ;
    Int_t code2(-1) ;
    getCompIntList(nset,iset,compIntList,code2,rangeName) ;
    delete vars ;
    delete nset ;
    delete iset ;
  } else {
    compIntList = &cache->_intList ;
  }
  
  const RooArgSet* nset = _normSet ; 
  CacheElem* pcache = getProjCache(nset) ;
  updateCoefficients(*pcache,nset) ;
  
  
  TIterator* compIntIter = compIntList->createIterator() ;
  _coefIter->Reset() ;
  RooAbsReal* pdfInt ;
  Double_t snormVal ;
  Double_t value(0) ;
  Int_t i(0) ;
  while((pdfInt = (RooAbsReal*)compIntIter->Next())) {
    if (_coefCache[i]!=0.) {
      snormVal = nset ? ((RooAbsReal*)pcache->_suppNormList.at(i))->getVal() : 1.0 ;
      Double_t intVal = pdfInt->getVal(nset) ;
      value += intVal*_coefCache[i]/snormVal ;
      cxcoutD(Eval) << "RooAddModel::evaluate(" << GetName() << ")  value += [" 
		      << pdfInt->GetName() << "] " << intVal << " * " << _coefCache[i] << " / " << snormVal << endl ;
    }
    i++ ;
  }
  delete compIntIter ;
  
  return value ;
  
}
Double_t RooAddModel::expectedEvents(const RooArgSet* nset) const 
{  
  
  
  Double_t expectedTotal(0.0);
  RooAbsPdf* pdf ;
    
  if (_allExtendable) {
    
    
    _pdfIter->Reset() ;
    while((pdf = (RooAbsPdf*)_pdfIter->Next())) {      
      expectedTotal += pdf->expectedEvents(nset) ;
    }   
    
  } else {
    
    
    _coefIter->Reset() ;
    RooAbsReal* coef ;
    while((coef=(RooAbsReal*)_coefIter->Next())) {
      expectedTotal += coef->getVal() ;
    }   
  }
  return expectedTotal;
}
void RooAddModel::selectNormalization(const RooArgSet* depSet, Bool_t force) 
{
  
  if (!force && _refCoefNorm.getSize()!=0) {
    return ;
  }
  if (!depSet) {
    fixCoefNormalization(RooArgSet()) ;
    return ;
  }
  RooArgSet* myDepSet = getObservables(depSet) ;
  fixCoefNormalization(*myDepSet) ;
  delete myDepSet ;
}
void RooAddModel::selectNormalizationRange(const char* rangeName, Bool_t force) 
{
  
  if (!force && _refCoefRangeName) {
    return ;
  }
  fixCoefRange(rangeName) ;
}
RooArgList RooAddModel::CacheElem::containedArgs(Action) 
{
  RooArgList allNodes;
  allNodes.add(_projList) ;
  allNodes.add(_suppProjList) ;
  allNodes.add(_refRangeProjList) ;
  allNodes.add(_rangeProjList) ;
  return allNodes ;
}
RooArgList RooAddModel::IntCacheElem::containedArgs(Action) 
{
  RooArgList allNodes(_intList) ;
  return allNodes ;
}
Last change: Tue May 13 17:03:46 2008
Last generated: 2008-05-13 17:03
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.