CMS 3D CMS Logo

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

#include <VerticalInterpHistPdf.h>

Inheritance diagram for FastVerticalInterpHistPdfBase:
FastVerticalInterpHistPdf FastVerticalInterpHistPdf2D

Classes

struct  Morph
 Must be public, for serialization. More...
 

Public Member Functions

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

Protected Member Functions

double smoothStepFunc (double x) const
 
void syncMorph (Morph &out, const FastTemplate &nominal, FastTemplate &lo, FastTemplate &hi) const
 not to be serialized More...
 
void syncTotal (FastTemplate &cache, const FastTemplate &cacheNominal, const FastTemplate &cacheNominalLog) const
 

Protected Attributes

TIterator * _coefIter
 Iterator over FUNC list. More...
 
RooListProxy _coefList
 
TIterator * _funcIter
 
RooListProxy _funcList
 
bool _init
 Iterator over coefficient list. More...
 
std::vector< RooAbsReal * > _morphParams
 not to be serialized More...
 
std::vector< Morph_morphs
 not to be serialized More...
 
SimpleCacheSentry _sentry
 not to be serialized More...
 
Int_t _smoothAlgo
 
Double_t _smoothRegion
 
RooRealProxy _x
 

Detailed Description

Definition at line 70 of file VerticalInterpHistPdf.h.

Constructor & Destructor Documentation

FastVerticalInterpHistPdfBase::FastVerticalInterpHistPdfBase ( )

Definition at line 230 of file VerticalInterpHistPdf.cc.

231 {
232  // Default constructor
233  _funcIter = _funcList.createIterator() ;
234  _coefIter = _coefList.createIterator() ;
235 }
TIterator * _coefIter
Iterator over FUNC list.
FastVerticalInterpHistPdfBase::FastVerticalInterpHistPdfBase ( const char *  name,
const char *  title,
const RooArgSet &  obs,
const RooArgList &  funcList,
const RooArgList &  coefList,
Double_t  smoothRegion = 1.,
Int_t  smoothAlgo = 1 
)

Definition at line 239 of file VerticalInterpHistPdf.cc.

References _coefIter, _coefList, _funcIter, _funcList, and TRACEME.

239  :
240  RooAbsPdf(name,title),
241  _funcList("funcList","List of pdfs",this),
242  _coefList("coefList","List of coefficients",this), // we should get shapeDirty when coefficients change
243  _smoothRegion(smoothRegion),
244  _smoothAlgo(smoothAlgo),
245  _init(false),
246  _morphs(), _morphParams()
247 {
248  TRACEME()
249 
250  if (inFuncList.getSize()!=2*inCoefList.getSize()+1) {
251  coutE(InputArguments) << "VerticalInterpHistPdf::VerticalInterpHistPdf(" << GetName()
252  << ") number of pdfs and coefficients inconsistent, must have Nfunc=1+2*Ncoef" << endl ;
253  assert(0);
254  }
255 
256  TIterator* funcIter = inFuncList.createIterator() ;
257  RooAbsArg* func;
258  while((func = (RooAbsArg*)funcIter->Next())) {
259  RooAbsPdf *pdf = dynamic_cast<RooAbsPdf*>(func);
260  if (!pdf) {
261  coutE(InputArguments) << "ERROR: VerticalInterpHistPdf::VerticalInterpHistPdf(" << GetName() << ") function " << func->GetName() << " is not of type RooAbsPdf" << endl;
262  assert(0);
263  }
264  RooArgSet *params = pdf->getParameters(obs);
265  if (params->getSize() > 0) {
266  coutE(InputArguments) << "ERROR: VerticalInterpHistPdf::VerticalInterpHistPdf(" << GetName() << ") pdf " << func->GetName() << " has some parameters." << endl;
267  assert(0);
268  }
269  delete params;
270  _funcList.add(*func) ;
271  }
272  delete funcIter;
273 
274  TIterator* coefIter = inCoefList.createIterator() ;
275  RooAbsArg* coef;
276  while((coef = (RooAbsArg*)coefIter->Next())) {
277  if (!dynamic_cast<RooAbsReal*>(coef)) {
278  coutE(InputArguments) << "ERROR: VerticalInterpHistPdf::VerticalInterpHistPdf(" << GetName() << ") coefficient " << coef->GetName() << " is not of type RooAbsReal" << endl;
279  assert(0);
280  }
281  _coefList.add(*coef) ;
282  }
283  delete coefIter;
284 
285  _funcIter = _funcList.createIterator() ;
286  _coefIter = _coefList.createIterator() ;
287 
288 }
#define TRACEME()
bool _init
Iterator over coefficient list.
TIterator * _coefIter
Iterator over FUNC list.
std::vector< RooAbsReal * > _morphParams
not to be serialized
std::vector< Morph > _morphs
not to be serialized
perl if(1 lt scalar(@::datatypes))
Definition: edlooper.cc:31
FastVerticalInterpHistPdfBase::FastVerticalInterpHistPdfBase ( const FastVerticalInterpHistPdfBase other,
const char *  name = 0 
)

Definition at line 291 of file VerticalInterpHistPdf.cc.

References _coefIter, _coefList, _funcIter, _funcList, _sentry, and SimpleCacheSentry::addVars().

291  :
292  RooAbsPdf(other,name),
293  _funcList("funcList",this,other._funcList),
294  _coefList("coefList",this,other._coefList),
296  _smoothAlgo(other._smoothAlgo),
297  _init(false),
298 #ifdef FastVerticalInterpHistPdf_CopyConstructor
300 #else
301  _morphs(), _morphParams()
302 #endif
303 {
304  // Copy constructor
305  _funcIter = _funcList.createIterator() ;
306  _coefIter = _coefList.createIterator() ;
307 #ifdef FastVerticalInterpHistPdf_CopyConstructor
309  _sentry.setValueDirty();
310 #endif
311 }
bool _init
Iterator over coefficient list.
TIterator * _coefIter
Iterator over FUNC list.
std::vector< RooAbsReal * > _morphParams
not to be serialized
std::vector< Morph > _morphs
not to be serialized
void addVars(const RooAbsCollection &vars)
SimpleCacheSentry _sentry
not to be serialized
FastVerticalInterpHistPdfBase::~FastVerticalInterpHistPdfBase ( )
virtual

Definition at line 315 of file VerticalInterpHistPdf.cc.

References _coefIter, and _funcIter.

316 {
317  // Destructor
318  delete _funcIter ;
319  delete _coefIter ;
320 }
TIterator * _coefIter
Iterator over FUNC list.

Member Function Documentation

virtual TObject* FastVerticalInterpHistPdfBase::clone ( const char *  newname) const
pure virtual
const RooArgList& FastVerticalInterpHistPdfBase::coefList ( ) const
inline

Definition at line 84 of file VerticalInterpHistPdf.h.

References _coefList.

84 { return _coefList ; }
Double_t FastVerticalInterpHistPdfBase::evaluate ( ) const
pure virtual
const RooArgList& FastVerticalInterpHistPdfBase::funcList ( ) const
inline

Definition at line 83 of file VerticalInterpHistPdf.h.

References _funcList.

83 { return _funcList ; }
Bool_t FastVerticalInterpHistPdfBase::selfNormalized ( ) const
inline

Definition at line 79 of file VerticalInterpHistPdf.h.

79 { return kTRUE; }
double FastVerticalInterpHistPdfBase::smoothStepFunc ( double  x) const
inlineprotected

Definition at line 114 of file VerticalInterpHistPdf.h.

References _smoothRegion.

Referenced by syncTotal().

114  {
115  if (fabs(x) >= _smoothRegion) return x > 0 ? +1 : -1;
116  double xnorm = x/_smoothRegion, xnorm2 = xnorm*xnorm;
117  return 0.125 * xnorm * (xnorm2 * (3.*xnorm2 - 10.) + 15);
118  }
Definition: DDAxes.h:10
void FastVerticalInterpHistPdfBase::syncMorph ( Morph out,
const FastTemplate nominal,
FastTemplate lo,
FastTemplate hi 
) const
protected

not to be serialized

Definition at line 380 of file VerticalInterpHistPdf.cc.

References _smoothAlgo, FastVerticalInterpHistPdfBase::Morph::diff, FastTemplate::LogRatio(), FastTemplate::Subtract(), FastVerticalInterpHistPdfBase::Morph::sum, and FastTemplate::SumDiff().

Referenced by FastVerticalInterpHistPdf::syncComponents(), and FastVerticalInterpHistPdf2D::syncComponents().

380  {
381  if (_smoothAlgo < 0) {
382  hi.LogRatio(nominal); lo.LogRatio(nominal);
383  //printf("Log-ratios for dimension %d: \n", dim); hi.Dump(); lo.Dump();
384  } else {
385  hi.Subtract(nominal); lo.Subtract(nominal);
386  //printf("Differences for dimension %d: \n", dim); hi.Dump(); lo.Dump();
387  }
388  FastTemplate::SumDiff(hi, lo, out.sum, out.diff);
389  //printf("Sum and diff for dimension %d: \n", dim); out.sum.Dump(); out.diff.Dump();
390 }
void Subtract(const FastTemplate &reference)
*this = *this - reference
tuple out
Definition: dbtoconf.py:99
static void SumDiff(const FastTemplate &h1, const FastTemplate &h2, FastTemplate &sum, FastTemplate &diff)
assigns sum and diff
void LogRatio(const FastTemplate &reference)
*this = log(*this)/(reference)
void FastVerticalInterpHistPdfBase::syncTotal ( FastTemplate cache,
const FastTemplate cacheNominal,
const FastTemplate cacheNominalLog 
) const
protected

Definition at line 427 of file VerticalInterpHistPdf.cc.

References _coefList, _init, _morphParams, _morphs, _sentry, _smoothAlgo, a, b, FastTemplate::CopyValues(), FastTemplate::CropUnderflows(), diffTreeTool::diff, FastTemplate::Exp(), i, FastTemplate::Meld(), SimpleCacheSentry::reset(), smoothStepFunc(), TRACEME, and x.

Referenced by FastVerticalInterpHistPdf::syncTotal(), and FastVerticalInterpHistPdf2D::syncTotal().

427  {
428  TRACEME()
429  /* === how the algorithm works, in theory ===
430  * let dhi = h_hi - h_nominal
431  * dlo = h_lo - h_nominal
432  * and x be the morphing parameter
433  * we define alpha = x * 0.5 * ((dhi-dlo) + (dhi+dlo)*smoothStepFunc(x));
434  * which satisfies:
435  * alpha(0) = 0
436  * alpha(+1) = dhi
437  * alpha(-1) = dlo
438  * alpha(x >= +1) = |x|*dhi
439  * alpha(x <= -1) = |x|*dlo
440  * alpha is continuous and has continuous first and second derivative, as smoothStepFunc has them
441  * === and in practice ===
442  * we already have computed the histogram for diff=(dhi-dlo) and sum=(dhi+dlo)
443  * so we just do template += (0.5 * x) * (diff + smoothStepFunc(x) * sum)
444  * ========================================== */
445 
446  // start from nominal
447  cache.CopyValues(_smoothAlgo < 0 ? cacheNominalLog : cacheNominal);
448  //printf("Cache initialized to nominal template: \n"); cacheNominal.Dump();
449 
450  // apply all morphs one by one
451  for (int i = 0, ndim = _coefList.getSize(); i < ndim; ++i) {
452  double x = _morphParams[i]->getVal();
453  double a = 0.5*x, b = smoothStepFunc(x);
454  cache.Meld(_morphs[i].diff, _morphs[i].sum, a, b);
455  //printf("Merged transformation for dimension %d, x = %+5.3f, step = %.3f: \n", i, x, b); cache.Dump();
456  }
457 
458  // if necessary go back to linear scale
459  if (_smoothAlgo < 0) {
460  cache.Exp();
461  //printf("Done exponential tranformation\n"); cache.Dump();
462  } else {
463  cache.CropUnderflows();
464  }
465 
466  // mark as done
467  _sentry.reset();
468  _init = true;
469 }
int i
Definition: DBlmapReader.cc:9
#define TRACEME()
void CropUnderflows(T minimum=1e-9)
protect from underflows (*this = max(*this, minimum));
bool _init
Iterator over coefficient list.
std::vector< RooAbsReal * > _morphParams
not to be serialized
std::vector< Morph > _morphs
not to be serialized
void Exp()
*this = exp(*this)
double b
Definition: hdecay.h:120
double smoothStepFunc(double x) const
SimpleCacheSentry _sentry
not to be serialized
double a
Definition: hdecay.h:121
Definition: DDAxes.h:10
void Meld(const FastTemplate &diff, const FastTemplate &sum, T x, T y)
Does this += x * (diff + (sum)*y)

Member Data Documentation

TIterator* FastVerticalInterpHistPdfBase::_coefIter
protected
RooListProxy FastVerticalInterpHistPdfBase::_coefList
protected
TIterator* FastVerticalInterpHistPdfBase::_funcIter
protected
RooListProxy FastVerticalInterpHistPdfBase::_funcList
protected
bool FastVerticalInterpHistPdfBase::_init
mutableprotected
std::vector<RooAbsReal *> FastVerticalInterpHistPdfBase::_morphParams
mutableprotected
std::vector<Morph> FastVerticalInterpHistPdfBase::_morphs
mutableprotected
SimpleCacheSentry FastVerticalInterpHistPdfBase::_sentry
mutableprotected
Int_t FastVerticalInterpHistPdfBase::_smoothAlgo
protected
Double_t FastVerticalInterpHistPdfBase::_smoothRegion
protected

Definition at line 92 of file VerticalInterpHistPdf.h.

Referenced by smoothStepFunc().

RooRealProxy FastVerticalInterpHistPdfBase::_x
protected

Definition at line 89 of file VerticalInterpHistPdf.h.