CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Attributes | Private Member Functions | Friends
VerticalInterpHistPdf Class Reference

#include <VerticalInterpHistPdf.h>

Inheritance diagram for VerticalInterpHistPdf:

Public Member Functions

virtual TObject * clone (const char *newname) const
 
const RooArgList & coefList () const
 
Double_t evaluate () const
 
const RooArgList & funcList () const
 
Bool_t selfNormalized () const
 
 VerticalInterpHistPdf ()
 
 VerticalInterpHistPdf (const char *name, const char *title, const RooRealVar &x, const RooArgList &funcList, const RooArgList &coefList, Double_t smoothRegion=1., Int_t smoothAlgo=1)
 
 VerticalInterpHistPdf (const VerticalInterpHistPdf &other, const char *name=0)
 
const RooRealProxy & x () const
 
virtual ~VerticalInterpHistPdf ()
 

Protected Attributes

TH1 ** _cacheSingle
 not to be serialized More...
 
int * _cacheSingleGood
 not to be serialized More...
 
TH1 * _cacheTotal
 
TIterator * _coefIter
 Iterator over FUNC list. More...
 
RooListProxy _coefList
 
TIterator * _funcIter
 
RooListProxy _funcList
 
SimpleCacheSentry _sentry
 Iterator over coefficient list. More...
 
Int_t _smoothAlgo
 
Double_t _smoothRegion
 
RooRealProxy _x
 

Private Member Functions

double smoothStepFunc (double x) const
 
void syncComponent (int which) const
 
ClassDef(VerticalInterpHistPdf,
1) protected void 
syncTotal () const
 not to be serialized More...
 

Friends

class FastVerticalInterpHistPdf
 

Detailed Description

Definition at line 20 of file VerticalInterpHistPdf.h.

Constructor & Destructor Documentation

VerticalInterpHistPdf::VerticalInterpHistPdf ( )

Definition at line 25 of file VerticalInterpHistPdf.cc.

Referenced by clone().

25  :
26  _cacheTotal(0),
27  _cacheSingle(0),
29 {
30  // Default constructor
31  _funcIter = _funcList.createIterator() ;
32  _coefIter = _coefList.createIterator() ;
33 }
TIterator * _coefIter
Iterator over FUNC list.
int * _cacheSingleGood
not to be serialized
TH1 ** _cacheSingle
not to be serialized
VerticalInterpHistPdf::VerticalInterpHistPdf ( const char *  name,
const char *  title,
const RooRealVar &  x,
const RooArgList &  funcList,
const RooArgList &  coefList,
Double_t  smoothRegion = 1.,
Int_t  smoothAlgo = 1 
)

Definition at line 37 of file VerticalInterpHistPdf.cc.

References _coefIter, _coefList, _funcIter, and _funcList.

37  :
38  RooAbsPdf(name,title),
39  _x("x","Independent variable",this,const_cast<RooRealVar&>(x)),
40  _funcList("funcList","List of pdfs",this),
41  _coefList("coefList","List of coefficients",this), // we should get shapeDirty when coefficients change
42  _smoothRegion(smoothRegion),
43  _smoothAlgo(smoothAlgo),
44  _cacheTotal(0),
45  _cacheSingle(0),
47 {
48 
49  if (inFuncList.getSize()!=2*inCoefList.getSize()+1) {
50  coutE(InputArguments) << "VerticalInterpHistPdf::VerticalInterpHistPdf(" << GetName()
51  << ") number of pdfs and coefficients inconsistent, must have Nfunc=1+2*Ncoef" << endl ;
52  assert(0);
53  }
54 
55  TIterator* funcIter = inFuncList.createIterator() ;
56  RooAbsArg* func;
57  while((func = (RooAbsArg*)funcIter->Next())) {
58  RooAbsPdf *pdf = dynamic_cast<RooAbsPdf*>(func);
59  if (!pdf) {
60  coutE(InputArguments) << "ERROR: VerticalInterpHistPdf::VerticalInterpHistPdf(" << GetName() << ") function " << func->GetName() << " is not of type RooAbsPdf" << endl;
61  assert(0);
62  }
63  RooArgSet *params = pdf->getParameters(RooArgSet(x));
64  if (params->getSize() > 0) {
65  coutE(InputArguments) << "ERROR: VerticalInterpHistPdf::VerticalInterpHistPdf(" << GetName() << ") pdf " << func->GetName() << " has some parameters." << endl;
66  assert(0);
67  }
68  delete params;
69  _funcList.add(*func) ;
70  }
71  delete funcIter;
72 
73  TIterator* coefIter = inCoefList.createIterator() ;
74  RooAbsArg* coef;
75  while((coef = (RooAbsArg*)coefIter->Next())) {
76  if (!dynamic_cast<RooAbsReal*>(coef)) {
77  coutE(InputArguments) << "ERROR: VerticalInterpHistPdf::VerticalInterpHistPdf(" << GetName() << ") coefficient " << coef->GetName() << " is not of type RooAbsReal" << endl;
78  assert(0);
79  }
80  _coefList.add(*coef) ;
81  }
82  delete coefIter;
83 
84  _funcIter = _funcList.createIterator() ;
85  _coefIter = _coefList.createIterator() ;
86 
87 }
const RooRealProxy & x() const
TIterator * _coefIter
Iterator over FUNC list.
int * _cacheSingleGood
not to be serialized
TH1 ** _cacheSingle
not to be serialized
VerticalInterpHistPdf::VerticalInterpHistPdf ( const VerticalInterpHistPdf other,
const char *  name = 0 
)

Definition at line 93 of file VerticalInterpHistPdf.cc.

References _coefIter, _coefList, _funcIter, and _funcList.

93  :
94  RooAbsPdf(other,name),
95  _x("x",this,other._x),
96  _funcList("funcList",this,other._funcList),
97  _coefList("coefList",this,other._coefList),
99  _smoothAlgo(other._smoothAlgo),
100  _cacheTotal(0),
101  _cacheSingle(0),
102  _cacheSingleGood(0)
103 {
104  // Copy constructor
105 
106  _funcIter = _funcList.createIterator() ;
107  _coefIter = _coefList.createIterator() ;
108 }
TIterator * _coefIter
Iterator over FUNC list.
int * _cacheSingleGood
not to be serialized
TH1 ** _cacheSingle
not to be serialized
VerticalInterpHistPdf::~VerticalInterpHistPdf ( )
virtual

Definition at line 113 of file VerticalInterpHistPdf.cc.

References _cacheSingle, _cacheSingleGood, _cacheTotal, _coefIter, _funcIter, _funcList, and i.

114 {
115  // Destructor
116  delete _funcIter ;
117  delete _coefIter ;
118  if (_cacheTotal) {
119  delete _cacheTotal;
120  for (int i = 0; i < _funcList.getSize(); ++i) delete _cacheSingle[i];
121  delete [] _cacheSingle;
122  delete [] _cacheSingleGood;
123  }
124 }
int i
Definition: DBlmapReader.cc:9
TIterator * _coefIter
Iterator over FUNC list.
int * _cacheSingleGood
not to be serialized
TH1 ** _cacheSingle
not to be serialized

Member Function Documentation

virtual TObject* VerticalInterpHistPdf::clone ( const char *  newname) const
inlinevirtual

Definition at line 26 of file VerticalInterpHistPdf.h.

References VerticalInterpHistPdf().

26 { return new VerticalInterpHistPdf(*this,newname) ; }
const RooArgList& VerticalInterpHistPdf::coefList ( ) const
inline

Definition at line 34 of file VerticalInterpHistPdf.h.

References _coefList.

34 { return _coefList ; }
Double_t VerticalInterpHistPdf::evaluate ( ) const

Definition at line 127 of file VerticalInterpHistPdf.cc.

References _cacheTotal, _coefList, _sentry, _x, gather_cfg::cout, SimpleCacheSentry::good(), i, and syncTotal().

128 {
129  if (_cacheTotal == 0) setupCaches();
130 #if 0
131  printf("Evaluate called for x %f, cache status %d: ", _x.arg().getVal(), _sentry.good());
132  int ndim = _coefList.getSize();
133  for (int i = 0; i < ndim; ++i) {
134  RooAbsReal *rar = dynamic_cast<RooAbsReal *>(_coefList.at(i));
135  printf("%s = %+6.4f ", rar->GetName(), rar->getVal());
136  }
137  std::cout << std::endl;
138 #endif
139 
140  if (!_sentry.good()) syncTotal();
141  int nbin = _cacheTotal->GetNbinsX(), ibin = _cacheTotal->FindBin(_x);
142  if (ibin < 1) ibin = 1;
143  else if (ibin > nbin) ibin = nbin;
144  return _cacheTotal->GetBinContent(ibin);
145 }
int i
Definition: DBlmapReader.cc:9
ClassDef(VerticalInterpHistPdf, 1) protected void syncTotal() const
not to be serialized
bool good() const
SimpleCacheSentry _sentry
Iterator over coefficient list.
tuple cout
Definition: gather_cfg.py:121
const RooArgList& VerticalInterpHistPdf::funcList ( ) const
inline

Definition at line 33 of file VerticalInterpHistPdf.h.

References _funcList.

33 { return _funcList ; }
Bool_t VerticalInterpHistPdf::selfNormalized ( ) const
inline

Definition at line 29 of file VerticalInterpHistPdf.h.

29 { return kTRUE; }
double VerticalInterpHistPdf::smoothStepFunc ( double  x) const
inlineprivate

Definition at line 63 of file VerticalInterpHistPdf.h.

References _smoothRegion.

Referenced by syncTotal().

63  {
64  if (fabs(x) >= _smoothRegion) return x > 0 ? +1 : -1;
65  double xnorm = x/_smoothRegion, xnorm2 = xnorm*xnorm;
66  return 0.125 * xnorm * (xnorm2 * (3.*xnorm2 - 10.) + 15);
67  }
const RooRealProxy & x() const
void VerticalInterpHistPdf::syncComponent ( int  which) const
private

Definition at line 148 of file VerticalInterpHistPdf.cc.

References _cacheSingle, _cacheSingleGood, _funcList, _smoothAlgo, _x, b, i, create_public_lumi_plots::log, and detailsBasic3DVector::y.

Referenced by syncTotal().

148  {
149  RooAbsPdf *pdfi = dynamic_cast<RooAbsPdf *>(_funcList.at(i));
150  if (_cacheSingle[i] != 0) delete _cacheSingle[i];
151  _cacheSingle[i] = pdfi->createHistogram("",dynamic_cast<const RooRealVar &>(_x.arg()));
152  _cacheSingle[i]->SetDirectory(0);
153  if (_cacheSingle[i]->Integral("width")) { _cacheSingle[i]->Scale(1.0/_cacheSingle[i]->Integral("width")); }
154  if (i > 0) {
155  for (int b = 1, nb = _cacheSingle[i]->GetNbinsX(); b <= nb; ++b) {
156  double y = _cacheSingle[i]->GetBinContent(b);
157  double y0 = _cacheSingle[0]->GetBinContent(b);
158  if (_smoothAlgo < 0) {
159  if (y > 0 && y0 > 0) {
160  double logk = log(y/y0);
161  // odd numbers correspond to up variations, even numbers to down variations,
162  // and down variations need -log(kappa) instead of log(kappa)
163  _cacheSingle[i]->SetBinContent(b, logk);
164  } else {
165  _cacheSingle[i]->SetBinContent(b, 0);
166  }
167  } else {
168  _cacheSingle[i]->SetBinContent(b, y - y0);
169  }
170  }
171  }
172  _cacheSingleGood[i] = true;
173 }
int i
Definition: DBlmapReader.cc:9
int * _cacheSingleGood
not to be serialized
double b
Definition: hdecay.h:120
TH1 ** _cacheSingle
not to be serialized
void VerticalInterpHistPdf::syncTotal ( ) const
private

not to be serialized

Definition at line 175 of file VerticalInterpHistPdf.cc.

References _cacheSingle, _cacheSingleGood, _cacheTotal, _coefIter, _coefList, _sentry, _smoothAlgo, alpha, b, alignCSCRings::e, create_public_lumi_plots::exp, i, SimpleCacheSentry::reset(), smoothStepFunc(), syncComponent(), and x().

Referenced by evaluate().

175  {
176  int ndim = _coefList.getSize();
177  for (int i = 0; i < 2*ndim+1; ++i) {
178  if (!_cacheSingleGood[i]) syncComponent(i);
179  }
180  for (int b = 1, nb = _cacheTotal->GetNbinsX(); b <= nb; ++b) {
181  double val = _cacheSingle[0]->GetBinContent(b);
182  _coefIter->Reset();
183  for (int i = 0; i < ndim; ++i) {
184  double dhi = _cacheSingle[2*i+1]->GetBinContent(b);
185  double dlo = _cacheSingle[2*i+2]->GetBinContent(b);
186  double x = (dynamic_cast<RooAbsReal *>(_coefIter->Next()))->getVal();
187  double alpha = x * 0.5 * ((dhi-dlo) + (dhi+dlo)*smoothStepFunc(x));
188  // alpha(0) = 0
189  // alpha(+1) = dhi
190  // alpha(-1) = dlo
191  // alpha(x >= +1) = |x|*dhi
192  // alpha(x <= -1) = |x|*dlo
193  // alpha is continuous and has continuous first and second derivative, as smoothStepFunc has them
194  if (_smoothAlgo < 0) {
195  val *= exp(alpha);
196  } else {
197  val += alpha;
198  }
199  }
200  if (val <= 0) val = 1e-9;
201  _cacheTotal->SetBinContent(b, val);
202  }
203  double norm = _cacheTotal->Integral("width");
204  if (norm > 0) _cacheTotal->Scale(1.0/norm);
205  _sentry.reset();
206 }
int i
Definition: DBlmapReader.cc:9
float alpha
Definition: AMPTWrapper.h:95
const RooRealProxy & x() const
TIterator * _coefIter
Iterator over FUNC list.
int * _cacheSingleGood
not to be serialized
double smoothStepFunc(double x) const
void syncComponent(int which) const
SimpleCacheSentry _sentry
Iterator over coefficient list.
double b
Definition: hdecay.h:120
TH1 ** _cacheSingle
not to be serialized
const RooRealProxy& VerticalInterpHistPdf::x ( ) const
inline

Friends And Related Function Documentation

friend class FastVerticalInterpHistPdf
friend

Definition at line 37 of file VerticalInterpHistPdf.h.

Member Data Documentation

TH1** VerticalInterpHistPdf::_cacheSingle
mutableprotected

not to be serialized

Definition at line 52 of file VerticalInterpHistPdf.h.

Referenced by syncComponent(), syncTotal(), and ~VerticalInterpHistPdf().

int* VerticalInterpHistPdf::_cacheSingleGood
mutableprotected

not to be serialized

Definition at line 53 of file VerticalInterpHistPdf.h.

Referenced by syncComponent(), syncTotal(), and ~VerticalInterpHistPdf().

TH1* VerticalInterpHistPdf::_cacheTotal
mutableprotected

Definition at line 49 of file VerticalInterpHistPdf.h.

Referenced by evaluate(), syncTotal(), and ~VerticalInterpHistPdf().

TIterator* VerticalInterpHistPdf::_coefIter
protected

Iterator over FUNC list.

Definition at line 45 of file VerticalInterpHistPdf.h.

Referenced by syncTotal(), VerticalInterpHistPdf(), and ~VerticalInterpHistPdf().

RooListProxy VerticalInterpHistPdf::_coefList
protected

Definition at line 41 of file VerticalInterpHistPdf.h.

Referenced by coefList(), evaluate(), syncTotal(), and VerticalInterpHistPdf().

TIterator* VerticalInterpHistPdf::_funcIter
protected

Definition at line 44 of file VerticalInterpHistPdf.h.

Referenced by VerticalInterpHistPdf(), and ~VerticalInterpHistPdf().

RooListProxy VerticalInterpHistPdf::_funcList
protected
SimpleCacheSentry VerticalInterpHistPdf::_sentry
mutableprotected

Iterator over coefficient list.

Definition at line 48 of file VerticalInterpHistPdf.h.

Referenced by evaluate(), and syncTotal().

Int_t VerticalInterpHistPdf::_smoothAlgo
protected

Definition at line 43 of file VerticalInterpHistPdf.h.

Referenced by syncComponent(), and syncTotal().

Double_t VerticalInterpHistPdf::_smoothRegion
protected

Definition at line 42 of file VerticalInterpHistPdf.h.

Referenced by smoothStepFunc().

RooRealProxy VerticalInterpHistPdf::_x
protected

Definition at line 39 of file VerticalInterpHistPdf.h.

Referenced by evaluate(), syncComponent(), and x().