Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooNovosibirsk.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * File: $Id: RooNovosibirsk.h,v 1.7 2007/07/12 20:30:49 wouter Exp $
5 * Authors: *
6 * DB, Dieter Best, UC Irvine, best@slac.stanford.edu *
7 * HT, Hirohisa Tanaka SLAC tanaka@slac.stanford.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#ifndef ROO_NOVOSIBIRSK
17#define ROO_NOVOSIBIRSK
18
19#include "RooAbsPdf.h"
20#include "RooRealProxy.h"
21
22class RooRealVar;
23class RooAbsReal;
24
25class RooNovosibirsk : public RooAbsPdf {
26public:
27 // Your constructor needs a name and title and then a list of the
28 // dependent variables and parameters used by this PDF. Use an
29 // underscore in the variable names to distinguish them from your
30 // own local versions.
32 RooNovosibirsk(const char *name, const char *title,
33 RooAbsReal& _x, RooAbsReal& _peak,
34 RooAbsReal& _width, RooAbsReal& _tail);
35
36 RooNovosibirsk(const RooNovosibirsk& other,const char* name=0) ;
37
38 virtual TObject* clone(const char* newname) const { return new RooNovosibirsk(*this,newname); }
39
40 Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=0) const ;
41 Double_t analyticalIntegral(Int_t code, const char* rangeName=0) const ;
42
43 // An empty constructor is usually ok
44 inline virtual ~RooNovosibirsk() { }
45
46protected:
47 double evaluate() const;
48 void computeBatch(cudaStream_t*, double* output, size_t nEvents, RooFit::Detail::DataMap const&) const;
49 inline bool canComputeBatchWithCuda() const { return true; }
50
51private:
56
57 ClassDef(RooNovosibirsk,1) // Novosibirsk PDF
58};
59
60#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
RooNovosibirsk implements the Novosibirsk function.
RooRealProxy width
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
double evaluate() const
If tail=eta=0 the Belle distribution becomes gaussian.
virtual TObject * clone(const char *newname) const
virtual ~RooNovosibirsk()
RooRealProxy peak
RooRealProxy tail
void computeBatch(cudaStream_t *, double *output, size_t nEvents, RooFit::Detail::DataMap const &) const
Compute multiple values of Novosibirsk distribution.
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