ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 /**
18 \file RooDecay.cxx
19 \class RooDecay
20 \ingroup Roofit
21 
22 Single or double sided decay function that can be analytically convolved
23 with any RooResolutionModel implementation
24 **/
25 
26 #include "RooFit.h"
27 
28 #include "Riostream.h"
29 #include "Riostream.h"
30 #include "RooDecay.h"
31 #include "RooRealVar.h"
32 #include "RooRandom.h"
33 
34 #include "TError.h"
35 
36 using namespace std;
37 
39 ;
40 
41 
42 
43 ////////////////////////////////////////////////////////////////////////////////
44 /// Constructor
45 
46 RooDecay::RooDecay(const char *name, const char *title,
47  RooRealVar& t, RooAbsReal& tau,
48  const RooResolutionModel& model, DecayType type) :
49  RooAbsAnaConvPdf(name,title,model,t),
50  _t("t","time",this,t),
51  _tau("tau","decay time",this,tau),
52  _type(type)
53 {
54  switch(type) {
55  case SingleSided:
56  _basisExp = declareBasis("exp(-@0/@1)",tau) ;
57  break ;
58  case Flipped:
59  _basisExp = declareBasis("exp(@0/@1)",tau) ;
60  break ;
61  case DoubleSided:
62  _basisExp = declareBasis("exp(-abs(@0)/@1)",tau) ;
63  break ;
64  }
65 }
66 
67 
68 
69 ////////////////////////////////////////////////////////////////////////////////
70 /// Copy constructor
71 
72 RooDecay::RooDecay(const RooDecay& other, const char* name) :
73  RooAbsAnaConvPdf(other,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 
83 ////////////////////////////////////////////////////////////////////////////////
84 /// Destructor
85 
87 {
88 }
89 
90 
91 
92 ////////////////////////////////////////////////////////////////////////////////
93 
94 Double_t RooDecay::coefficient(Int_t /*basisIndex*/) const
95 {
96  return 1 ;
97 }
98 
99 
100 
101 ////////////////////////////////////////////////////////////////////////////////
102 
103 Int_t RooDecay::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
104 {
105  if (matchArgs(directVars,generateVars,_t)) return 1 ;
106  return 0 ;
107 }
108 
109 
110 
111 ////////////////////////////////////////////////////////////////////////////////
112 
114 {
115  R__ASSERT(code==1) ;
116 
117  // Generate delta-t dependent
118  while(1) {
119  Double_t rand = RooRandom::uniform() ;
120  Double_t tval(0) ;
121 
122  switch(_type) {
123  case SingleSided:
124  tval = -_tau*log(rand);
125  break ;
126  case Flipped:
127  tval= +_tau*log(rand);
128  break ;
129  case DoubleSided:
130  tval = (rand<=0.5) ? -_tau*log(2*rand) : +_tau*log(2*(rand-0.5)) ;
131  break ;
132  }
133 
134  if (tval<_t.max() && tval>_t.min()) {
135  _t = tval ;
136  break ;
137  }
138  }
139 }
RooRealProxy _tau
Definition: RooDecay.h:43
virtual ~RooDecay()
Destructor.
Definition: RooDecay.cxx:86
RooRealProxy _t
Definition: RooDecay.h:42
RooDecay()
Definition: RooDecay.h:28
#define R__ASSERT(e)
Definition: TError.h:98
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
Definition: RooDecay.cxx:103
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
Int_t _basisExp
Definition: RooDecay.h:45
TThread * t[5]
Definition: threadsh1.C:13
TPaveLabel title(3, 27.1, 15, 28.7,"ROOT Environment and Tools")
static Double_t uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition: RooRandom.cxx:84
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
int type
Definition: TGX11.cxx:120
Int_t declareBasis(const char *expression, const RooArgList &params)
Declare a basis function for use in this physics model.
void generateEvent(Int_t code)
Interface for generation of anan event using the algorithm corresponding to the specified code...
Definition: RooDecay.cxx:113
DecayType _type
Definition: RooDecay.h:44
#define name(a, b)
Definition: linkTestLib0.cpp:5
ClassImp(RooDecay)
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
virtual Double_t coefficient(Int_t basisIndex) const
Definition: RooDecay.cxx:94
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
double log(double)
Single or double sided decay function that can be analytically convolved with any RooResolutionModel ...
Definition: RooDecay.h:22