Logo ROOT  
Reference Guide
RooBreitWigner.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
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
18/** \class RooBreitWigner
19 \ingroup Roofit
20
21Class RooBreitWigner is a RooAbsPdf implementation
22that models a non-relativistic Breit-Wigner shape
23**/
24
25#include "RooFit.h"
26
27#include "Riostream.h"
28#include "Riostream.h"
29#include <math.h>
30
31#include "RooBreitWigner.h"
32#include "RooAbsReal.h"
33#include "RooRealVar.h"
34#include "BatchHelpers.h"
35// #include "RooFitTools/RooRandom.h"
36
37using namespace std;
38
40
41////////////////////////////////////////////////////////////////////////////////
42
43RooBreitWigner::RooBreitWigner(const char *name, const char *title,
44 RooAbsReal& _x, RooAbsReal& _mean,
45 RooAbsReal& _width) :
46 RooAbsPdf(name,title),
47 x("x","Dependent",this,_x),
48 mean("mean","Mean",this,_mean),
49 width("width","Width",this,_width)
50{
51}
52
53////////////////////////////////////////////////////////////////////////////////
54
56 RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
57 width("width",this,other.width)
58{
59}
60
61////////////////////////////////////////////////////////////////////////////////
62
64{
65 Double_t arg= x - mean;
66 return 1. / (arg*arg + 0.25*width*width);
67}
68
69////////////////////////////////////////////////////////////////////////////////
70
71namespace {
72//Author: Emmanouil Michalainas, CERN 21 August 2019
73
74template<class Tx, class Tmean, class Twidth>
75void compute( size_t batchSize,
76 double * __restrict output,
77 Tx X, Tmean M, Twidth W)
78{
79 for (size_t i=0; i<batchSize; i++) {
80 const double arg = X[i]-M[i];
81 output[i] = 1 / (arg*arg + 0.25*W[i]*W[i]);
82 }
83}
84};
85
86RooSpan<double> RooBreitWigner::evaluateBatch(std::size_t begin, std::size_t batchSize) const {
87 using namespace BatchHelpers;
88 auto xData = x.getValBatch(begin, batchSize);
89 auto meanData = mean.getValBatch(begin, batchSize);
90 auto widthData = width.getValBatch(begin, batchSize);
91 const bool batchX = !xData.empty();
92 const bool batchMean = !meanData.empty();
93 const bool batchWidth = !widthData.empty();
94
95 if (!batchX && !batchMean && !batchWidth) {
96 return {};
97 }
98 batchSize = findSize({ xData, meanData, widthData });
99 auto output = _batchData.makeWritableBatchUnInit(begin, batchSize);
100
101 if (batchX && !batchMean && !batchWidth ) {
102 compute(batchSize, output.data(), xData, BracketAdapter<double>(mean), BracketAdapter<double>(width));
103 }
104 else if (!batchX && batchMean && !batchWidth ) {
105 compute(batchSize, output.data(), BracketAdapter<double>(x), meanData, BracketAdapter<double>(width));
106 }
107 else if (batchX && batchMean && !batchWidth ) {
108 compute(batchSize, output.data(), xData, meanData, BracketAdapter<double>(width));
109 }
110 else if (!batchX && !batchMean && batchWidth ) {
111 compute(batchSize, output.data(), BracketAdapter<double>(x), BracketAdapter<double>(mean), widthData);
112 }
113 else if (batchX && !batchMean && batchWidth ) {
114 compute(batchSize, output.data(), xData, BracketAdapter<double>(mean), widthData);
115 }
116 else if (!batchX && batchMean && batchWidth ) {
117 compute(batchSize, output.data(), BracketAdapter<double>(x), meanData, widthData);
118 }
119 else if (batchX && batchMean && batchWidth ) {
120 compute(batchSize, output.data(), xData, meanData, widthData);
121 }
122 return output;
123}
124
125////////////////////////////////////////////////////////////////////////////////
126
127Int_t RooBreitWigner::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
128{
129 if (matchArgs(allVars,analVars,x)) return 1 ;
130 return 0 ;
131}
132
133////////////////////////////////////////////////////////////////////////////////
134
135Double_t RooBreitWigner::analyticalIntegral(Int_t code, const char* rangeName) const
136{
137 switch(code) {
138 case 1:
139 {
140 Double_t c = 2./width;
141 return c*(atan(c*(x.max(rangeName)-mean)) - atan(c*(x.min(rangeName)-mean)));
142 }
143 }
144
145 assert(0) ;
146 return 0 ;
147}
#define c(i)
Definition: RSha256.hxx:101
int Int_t
Definition: RtypesCore.h:41
double Double_t
Definition: RtypesCore.h:55
#define ClassImp(name)
Definition: Rtypes.h:365
include TDocParser_001 C image html pict1_TDocParser_001 png width
Definition: TDocParser.cxx:121
char name[80]
Definition: TGX11.cxx:109
double atan(double)
RooSpan< double > makeWritableBatchUnInit(std::size_t begin, std::size_t batchSize)
Make a batch and return a span pointing to the pdf-local memory.
Definition: BatchData.cxx:58
Little adapter that gives a bracket operator to types that don't have one.
Definition: BatchHelpers.h:58
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:59
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
BatchHelpers::BatchData _batchData
Definition: RooAbsReal.h:447
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
Class RooBreitWigner is a RooAbsPdf implementation that models a non-relativistic Breit-Wigner shape.
Double_t 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.
RooSpan< double > evaluateBatch(std::size_t begin, std::size_t batchSize) const
Evaluate function for a batch of input data points.
RooRealProxy width
RooRealProxy mean
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooRealProxy x
A simple container to hold a batch of data values.
Definition: RooSpan.h:32
Double_t min(const char *rname=0) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
RooSpan< const double > getValBatch(std::size_t begin, std::size_t batchSize) const
Double_t max(const char *rname=0) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
Double_t x[n]
Definition: legend1.C:17
size_t findSize(std::vector< RooSpan< const double > > parameters)
This function returns the minimum size of the non-zero-sized batches.
static void output(int code)
Definition: gifencode.c:226