Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooBreitWigner.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * File: $Id: RooBreitWigner.h,v 1.9 2007/07/12 20:30:49 wouter Exp $
5 * Authors: *
6 * AS, Abi Soffer, Colorado State University, abi@slac.stanford.edu *
7 * TS, Thomas Schietinger, SLAC, schieti@slac.stanford.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * Colorado State University *
11 * and Stanford University. All rights reserved. *
12 * *
13 * Redistribution and use in source and binary forms, *
14 * with or without modification, are permitted according to the terms *
15 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
16 *****************************************************************************/
17#ifndef ROO_BREITWIGNER
18#define ROO_BREITWIGNER
19
20#include "RooAbsPdf.h"
21#include "RooRealProxy.h"
22
23class RooRealVar;
24
25class RooBreitWigner : public RooAbsPdf {
26public:
28 RooBreitWigner(const char *name, const char *title,
29 RooAbsReal& _x, RooAbsReal& _mean, RooAbsReal& _width);
30 RooBreitWigner(const RooBreitWigner& other, const char* name=0) ;
31 virtual TObject* clone(const char* newname) const { return new RooBreitWigner(*this,newname); }
32 inline virtual ~RooBreitWigner() { }
33
34 Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=0) const ;
35 Double_t analyticalIntegral(Int_t code, const char* rangeName=0) const ;
36
37protected:
38
42
43 double evaluate() const;
44 void computeBatch(cudaStream_t*, double* output, size_t nEvents, RooFit::Detail::DataMap const&) const;
45 inline bool canComputeBatchWithCuda() const { return true; }
46
47// void initGenerator();
48// Int_t generateDependents();
49
50private:
51
52 ClassDef(RooBreitWigner,1) // Breit Wigner PDF
53};
54
55#endif
double Double_t
Definition RtypesCore.h:59
#define ClassDef(name, id)
Definition Rtypes.h:325
char name[80]
Definition TGX11.cxx:110
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:35
Class RooBreitWigner is a RooAbsPdf implementation that models a non-relativistic Breit-Wigner shape.
double evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
void computeBatch(cudaStream_t *, double *output, size_t nEvents, RooFit::Detail::DataMap const &) const
Compute multiple values of BreitWigner distribution.
RooRealProxy width
virtual TObject * clone(const char *newname) const
virtual ~RooBreitWigner()
RooRealProxy mean
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
bool canComputeBatchWithCuda() const
RooRealProxy x
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
Mother of all ROOT objects.
Definition TObject.h:41
static void output(int code)
Definition gifencode.c:226