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
40
41////////////////////////////////////////////////////////////////////////////////
42/// Create a new RooDecay.
43/// \param[in] name Name of this object.
44/// \param[in] title Title (for *e.g.* plotting)
45/// \param[in] t Convolution variable (*e.g.* time).
46/// \param[in] tau Decay constant.
47/// \param[in] model Resolution model for the convolution.
48/// \param[in] type One of the decays types `SingleSided, Flipped, DoubleSided`
49RooDecay::RooDecay(const char *name, const char *title,
50 RooRealVar& t, RooAbsReal& tau,
51 const RooResolutionModel& model, DecayType type) :
52 RooAbsAnaConvPdf(name,title,model,t),
53 _t("t","time",this,t),
54 _tau("tau","decay time",this,tau),
55 _type(type)
56{
57 switch(type) {
58 case SingleSided:
59 _basisExp = declareBasis("exp(-@0/@1)",tau) ;
60 break ;
61 case Flipped:
62 _basisExp = declareBasis("exp(@0/@1)",tau) ;
63 break ;
64 case DoubleSided:
65 _basisExp = declareBasis("exp(-abs(@0)/@1)",tau) ;
66 break ;
67 }
68}
69
70////////////////////////////////////////////////////////////////////////////////
71/// Copy constructor
72
73RooDecay::RooDecay(const RooDecay& other, const char* name) :
75 _t("t",this,other._t),
76 _tau("tau",this,other._tau),
77 _type(other._type),
78 _basisExp(other._basisExp)
79{
80}
81
82////////////////////////////////////////////////////////////////////////////////
83
84double RooDecay::coefficient(Int_t /*basisIndex*/) const
85{
86 return 1 ;
87}
88
89////////////////////////////////////////////////////////////////////////////////
90
91Int_t RooDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool /*staticInitOK*/) const
92{
93 if (matchArgs(directVars,generateVars,_t)) return 1 ;
94 return 0 ;
95}
96
97////////////////////////////////////////////////////////////////////////////////
98
100{
101 R__ASSERT(code==1) ;
102
103 // Generate delta-t dependent
104 while(true) {
105 double rand = RooRandom::uniform() ;
106 double tval(0) ;
107
108 switch(_type) {
109 case SingleSided:
110 tval = -_tau*log(rand);
111 break ;
112 case Flipped:
113 tval= +_tau*log(rand);
114 break ;
115 case DoubleSided:
116 tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
117 break ;
118 }
119
120 if (tval<_t.max() && tval>_t.min()) {
121 _t = tval ;
122 break ;
123 }
124 }
125}
#define ClassImp(name)
Definition Rtypes.h:382
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
Definition TGX11.cxx:110
Base class for PDFs that represent a physics model that can be analytically convolved with a resoluti...
Int_t declareBasis(const char *expression, const RooArgList &params)
Declare a basis function for use in this physics model.
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
bool matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Single or double sided decay function that can be analytically convolved with any RooResolutionModel ...
Definition RooDecay.h:22
double coefficient(Int_t basisIndex) const override
Definition RooDecay.cxx:84
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:91
Int_t _basisExp
Definition RooDecay.h:44
RooRealProxy _t
Definition RooDecay.h:41
RooRealProxy _tau
Definition RooDecay.h:42
void generateEvent(Int_t code) override
Interface for generation of an event using the algorithm corresponding to the specified code.
Definition RooDecay.cxx:99
DecayType _type
Definition RooDecay.h:43
static double uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition RooRandom.cxx:78
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...