Logo ROOT  
Reference Guide
Loading...
Searching...
No Matches
RooDecay.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17/** \class RooDecay
18 \ingroup Roofit
19
20Single or double sided decay function that can be analytically convolved
21with any RooResolutionModel implementation. It declares the basis functions
22for the analytical convolution with a RooResolutionModel. See RooAbsAnaConvPdf.
23\f[
24 \mathrm{basis} = \begin{cases}
25 \exp\left(-\frac{t}{\tau}\right) & \mathrm{SingleSided} \\
26 \exp\left( \frac{t}{\tau}\right) & \mathrm{Flipped} \\
27 \exp\left(-\frac{|t|}{\tau}\right) & \mathrm{DoubleSided}
28 \end{cases}
29\f]
30**/
31
32#include "RooDecay.h"
33
34#include "RooRealVar.h"
35#include "RooRandom.h"
36
37#include "TError.h"
38
39
40////////////////////////////////////////////////////////////////////////////////
41/// Create a new RooDecay.
42/// \param[in] name Name of this object.
43/// \param[in] title Title (for *e.g.* plotting)
44/// \param[in] t Convolution variable (*e.g.* time).
45/// \param[in] tau Decay constant.
46/// \param[in] model Resolution model for the convolution.
47/// \param[in] type One of the decays types `SingleSided, Flipped, DoubleSided`
48RooDecay::RooDecay(const char *name, const char *title,
49 RooRealVar& t, RooAbsReal& tau,
50 const RooResolutionModel& model, DecayType type) :
51 RooAbsAnaConvPdf(name,title,model,t),
52 _t("t","time",this,t),
53 _tau("tau","decay time",this,tau),
54 _type(type)
55{
56 switch(type) {
57 case SingleSided:
58 _basisExp = declareBasis("exp(-@0/@1)",tau) ;
59 break ;
60 case Flipped:
61 _basisExp = declareBasis("exp(@0/@1)",tau) ;
62 break ;
63 case DoubleSided:
64 _basisExp = declareBasis("exp(-abs(@0)/@1)",tau) ;
65 break ;
66 }
67}
68
69////////////////////////////////////////////////////////////////////////////////
70/// Copy constructor
71
72RooDecay::RooDecay(const RooDecay& other, const char* name) :
74 _t("t",this,other._t),
75 _tau("tau",this,other._tau),
76 _type(other._type),
77 _basisExp(other._basisExp)
78{
79}
80
81////////////////////////////////////////////////////////////////////////////////
82
83double RooDecay::coefficient(Int_t /*basisIndex*/) const
84{
85 return 1 ;
86}
87
88////////////////////////////////////////////////////////////////////////////////
89
90Int_t RooDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool /*staticInitOK*/) const
91{
92 if (matchArgs(directVars,generateVars,_t)) return 1 ;
93 return 0 ;
94}
95
96////////////////////////////////////////////////////////////////////////////////
97
99{
100 R__ASSERT(code==1) ;
101
102 // Generate delta-t dependent
103 while(true) {
104 double rand = RooRandom::uniform() ;
105 double tval(0) ;
106
107 switch(_type) {
108 case SingleSided:
109 tval = -_tau*log(rand);
110 break ;
111 case Flipped:
112 tval= +_tau*log(rand);
113 break ;
114 case DoubleSided:
115 tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
116 break ;
117 }
118
119 if (tval<_t.max() && tval>_t.min()) {
120 _t = tval ;
121 break ;
122 }
123 }
124}
int Int_t
Signed integer 4 bytes (int).
Definition RtypesCore.h:59
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
char name[80]
Definition TGX11.cxx:148
RooAbsAnaConvPdf()
Default constructor, required for persistence.
Int_t declareBasis(const char *expression, const RooArgList &params)
Declare a basis function for use in this physics model.
friend class RooAbsReal
Definition RooAbsPdf.h:342
bool matchArgs(const RooArgSet &allDeps, RooArgSet &analDeps, const RooArgProxy &a, const Proxies &... proxies) const
Definition RooAbsReal.h:425
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
double coefficient(Int_t basisIndex) const override
Definition RooDecay.cxx:83
RooDecay()
Definition RooDecay.h:28
@ DoubleSided
Definition RooDecay.h:25
@ SingleSided
Definition RooDecay.h:25
@ Flipped
Definition RooDecay.h:25
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const override
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
Definition RooDecay.cxx:90
Int_t _basisExp
Definition RooDecay.h:53
RooRealProxy _t
Definition RooDecay.h:50
RooRealProxy _tau
Definition RooDecay.h:51
void generateEvent(Int_t code) override
Interface for generation of an event using the algorithm corresponding to the specified code.
Definition RooDecay.cxx:98
DecayType _type
Definition RooDecay.h:52
static double uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1).
Definition RooRandom.cxx:77
Variable that can be changed from the outside.
Definition RooRealVar.h:37
RooResolutionModel is the base class for PDFs that represent a resolution model that can be convolute...