JHUGen MELA  JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
RooSpinOne_7D Class Reference

#include <RooSpinOne_7D.h>

Inheritance diagram for RooSpinOne_7D:
Collaboration diagram for RooSpinOne_7D:

Public Member Functions

 RooSpinOne_7D ()
 
 RooSpinOne_7D (const char *name, const char *title, RooAbsReal &_mzz, RooAbsReal &_m1, RooAbsReal &_m2, RooAbsReal &_h1, RooAbsReal &_h2, RooAbsReal &_hs, RooAbsReal &_Phi, RooAbsReal &_Phi1, RooAbsReal &_g1Val, RooAbsReal &_g2Val, RooAbsReal &_R1Val, RooAbsReal &_R2Val, RooAbsReal &_aParam, RooAbsReal &_mZ, RooAbsReal &_gamZ)
 
 RooSpinOne_7D (const RooSpinOne_7D &other, const char *name=0)
 
virtual TObject * clone (const char *newname) const
 
virtual ~RooSpinOne_7D ()
 
Int_t getAnalyticalIntegral (RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
 
Double_t analyticalIntegral (Int_t code, const char *rangeName=0) const
 
void setZZ4fOrdering (Bool_t flag=true)
 

Protected Member Functions

Double_t evaluate () const
 

Protected Attributes

RooRealProxy mzz
 
RooRealProxy m1
 
RooRealProxy m2
 
RooRealProxy h1
 
RooRealProxy h2
 
RooRealProxy hs
 
RooRealProxy Phi
 
RooRealProxy Phi1
 
RooRealProxy g1Val
 
RooRealProxy g2Val
 
RooRealProxy R1Val
 
RooRealProxy R2Val
 
RooRealProxy aParam
 
RooRealProxy mZ
 
RooRealProxy gamZ
 
Bool_t ZZ4fOrdering
 

Detailed Description

Definition at line 18 of file RooSpinOne_7D.h.

Constructor & Destructor Documentation

◆ RooSpinOne_7D() [1/3]

RooSpinOne_7D::RooSpinOne_7D ( )
inline

Definition at line 20 of file RooSpinOne_7D.h.

20 {};

◆ RooSpinOne_7D() [2/3]

RooSpinOne_7D::RooSpinOne_7D ( const char *  name,
const char *  title,
RooAbsReal &  _mzz,
RooAbsReal &  _m1,
RooAbsReal &  _m2,
RooAbsReal &  _h1,
RooAbsReal &  _h2,
RooAbsReal &  _hs,
RooAbsReal &  _Phi,
RooAbsReal &  _Phi1,
RooAbsReal &  _g1Val,
RooAbsReal &  _g2Val,
RooAbsReal &  _R1Val,
RooAbsReal &  _R2Val,
RooAbsReal &  _aParam,
RooAbsReal &  _mZ,
RooAbsReal &  _gamZ 
)

Definition at line 20 of file RooSpinOne_7D.cc.

36  :
37  RooAbsPdf(name, title),
38  mzz("mzz", "mzz", this, _mzz),
39  m1("m1", "m1", this, _m1),
40  m2("m2", "m2", this, _m2),
41  h1("h1", "h1", this, _h1),
42  h2("h2", "h2", this, _h2),
43  hs("hs", "hs", this, _hs),
44  Phi("Phi", "Phi", this, _Phi),
45  Phi1("Phi1", "Phi1", this, _Phi1),
46  g1Val("g1Val", "g1Val", this, _g1Val),
47  g2Val("g2Val", "g2Val", this, _g2Val),
48  R1Val("R1Val", "R1Val", this, _R1Val),
49  R2Val("R2Val", "R2Val", this, _R2Val),
50  aParam("aParam", "aParam", this, _aParam),
51  mZ("mZ", "mZ", this, _mZ),
52  gamZ("gamZ", "gamZ", this, _gamZ),
53  ZZ4fOrdering(true)
54 {}

◆ RooSpinOne_7D() [3/3]

RooSpinOne_7D::RooSpinOne_7D ( const RooSpinOne_7D other,
const char *  name = 0 
)

Definition at line 57 of file RooSpinOne_7D.cc.

57  :
58 RooAbsPdf(other, name),
59 mzz("mzz", this, other.mzz),
60 m1("m1", this, other.m1),
61 m2("m2", this, other.m2),
62 h1("h1", this, other.h1),
63 h2("h2", this, other.h2),
64 hs("hs", this, other.hs),
65 Phi("Phi", this, other.Phi),
66 Phi1("Phi1", this, other.Phi1),
67 g1Val("g1Val", this, other.g1Val),
68 g2Val("g2Val", this, other.g2Val),
69 R1Val("R1Val", this, other.R1Val),
70 R2Val("R2Val", this, other.R2Val),
71 aParam("aParam", this, other.aParam),
72 mZ("mZ", this, other.mZ),
73 gamZ("gamZ", this, other.gamZ),
75 {}

◆ ~RooSpinOne_7D()

virtual RooSpinOne_7D::~RooSpinOne_7D ( )
inlinevirtual

Definition at line 39 of file RooSpinOne_7D.h.

39 {}

Member Function Documentation

◆ analyticalIntegral()

Double_t RooSpinOne_7D::analyticalIntegral ( Int_t  code,
const char *  rangeName = 0 
) const

Definition at line 224 of file RooSpinOne_7D.cc.

224  {
225  bool isZZ = true;
226  if (mZ < 90.) isZZ = false;
227  if (isZZ) {
228  if ((m1+m2) > mzz || m2>m1) return 1e-9;
229  }
230  else {
231  if ((m1+m2) > mzz) return 1e-9;
232  }
233  double nanval = sqrt((1 - TMath::Power(m1 - m2, 2)/TMath::Power(mzz, 2))*(1 - TMath::Power(m1 + m2, 2)/TMath::Power(mzz, 2)));
234  if (nanval != nanval) return 1e-9;
235 
236  //-----------------------------------------------------------------------
237  // propagator
238  //-----------------------------------------------------------------------
239 
240  Double_t betaValSquared = (1.-(pow(m1-m2, 2)/pow(mzz, 2)))*(1.-(pow(m1+m2, 2)/pow(mzz, 2)));
241  Double_t betaVal = sqrt(betaValSquared);
242 
243  Double_t term1Coeff = (pow(m1, 3))/((pow(m1, 2)-pow(mZ, 2))*(pow(m1, 2)-pow(mZ, 2))+pow(mZ, 2)*pow(gamZ, 2));
244  Double_t term2Coeff = (pow(m2, 3))/((pow(m2, 2)-pow(mZ, 2))*(pow(m2, 2)-pow(mZ, 2))+pow(mZ, 2)*pow(gamZ, 2));
245 
246  //-----------------------------------------------------------------------
247  // Helicity Amplitudes
248  //-----------------------------------------------------------------------
249  // calculating the angular parameters from the coupling constants
250  // See http://www.pha.jhu.edu/~gritsan/FORM/result_spin1.txt
251 
252  Double_t x = (mzz*mzz-m1*m1-m2*m2)/(2.0*m1*m2);
253 
254  Double_t hs_ = -hs;
255 
256  Double_t f00Real = TMath::Power(mzz, -1)*g1Val*(sqrt(x*x-1))*(m1*m1-m2*m2);
257  Double_t f00Imag = 0.;
258 
259  Double_t fppReal = 0.;
260  Double_t fppImag = TMath::Power(mzz, -1)*g2Val*(m1*m1-m2*m2);
261 
262  Double_t fmmReal = 0.;
263  Double_t fmmImag = -1. * fppImag;
264 
265 
266  Double_t fp0Real = m1*g1Val * (sqrt(x*x-1));
267  Double_t fp0Imag = Power(mzz, -2.)*Power(m2, 3)*g2Val * (-1./2.)
268  + Power(mzz, -2.)*m1*m1*m2*g2Val * (1 + 2*(x*x-1))
269  + Power(mzz, -2.)*Power(m1, 4)*Power(m2, -1)*g2Val * (-1./2.)
270  + m2*g2Val * (-1./2.)
271  + m1*m1*Power(m2, -1)*g2Val * (1./2.);
272 
273 
274  Double_t f0pReal = m2*g1Val * (-sqrt(x*x-1));
275  Double_t f0pImag = Power(mzz, -2.)*Power(m1, -1)*Power(m2, 4)*g2Val * (1./2.)
276  + Power(mzz, -2.)*m1*Power(m2, 2)*g2Val * (-1 - 2*(x*x-1))
277  + Power(mzz, -2.)*Power(m1, 3)*g2Val * (1./2.)
278  + Power(m1, -1)*Power(m2, 2)*g2Val * (-1./2.)
279  + m1*g2Val * (1./2.);
280 
281  Double_t f0mReal = m2*g1Val * (-sqrt(x*x-1));
282  Double_t f0mImag = Power(mzz, -2.)*Power(m1, -1)*Power(m2, 4)*g2Val * (-1./2.)
283  + Power(mzz, -2.)*m1*Power(m2, 2)*g2Val * (1 + 2*(x*x-1))
284  + Power(mzz, -2.)*Power(m1, 3)*g2Val * (-1./2.)
285  + Power(m1, -1)*Power(m2, 2)*g2Val * (1./2.)
286  + m1*g2Val * (-1./2.);
287 
288  Double_t fm0Real = m1*g1Val * (sqrt(x*x-1));
289  Double_t fm0Imag = Power(mzz, -2.)*Power(m2, 3)*g2Val * (1./2.)
290  + Power(mzz, -2.)*m1*m1*m2*g2Val * (-1 - 2*(x*x-1))
291  + Power(mzz, -2.)*Power(m1, 4)*Power(m2, -1)*g2Val * (1./2.)
292  + m2*g2Val * (1./2.)
293  + m1*m1*Power(m2, -1)*g2Val * (-1./2.);
294 
295  Double_t f00 = f00Imag*f00Imag + f00Real*f00Real;
296  Double_t fpp = fppImag*fppImag + fppReal*fppReal;
297  Double_t fmm = fmmImag*fmmImag + fmmReal*fmmReal;
298  Double_t fp0 = fp0Imag*fp0Imag + fp0Real*fp0Real;
299  Double_t f0p = f0pImag*f0pImag + f0pReal*f0pReal;
300  Double_t fm0 = fm0Imag*fm0Imag + fm0Real*fm0Real;
301  Double_t f0m = f0mImag*f0mImag + f0mReal*f0mReal;
302 
303  Double_t phi00=atan2(f00Imag, f00Real);
304  Double_t phipp=atan2(fppImag, fppReal)-phi00;
305  Double_t phimm=atan2(fmmImag, fmmReal)-phi00;
306  Double_t phip0=atan2(fp0Imag, fp0Real)-phi00;
307  Double_t phi0p=atan2(f0pImag, f0pReal)-phi00;
308  Double_t phim0=atan2(fm0Imag, fm0Real)-phi00;
309  Double_t phi0m=atan2(f0mImag, f0mReal)-phi00;
310 
311  Double_t integral=0;
312 
313  switch (code)
314  {
315  // Integrate out all angles
316  case 6:
317  {
318  integral = 0.;
319  integral+=
320  (32.*f00*TMath::Power(Pi(), 2))/27.;
321  integral+=
322  (32.*fpp*TMath::Power(Pi(), 2))/27.;
323  integral+=
324  (32.*fmm*TMath::Power(Pi(), 2))/27.;
325  integral+=
326  (32.*fp0*TMath::Power(Pi(), 2))/27.;
327  integral+=
328  (32.*f0m*TMath::Power(Pi(), 2))/27.;
329  integral+=
330  (32.*f0p*TMath::Power(Pi(), 2))/27.;
331  integral+=
332  (32.*fm0*TMath::Power(Pi(), 2))/27.;
333  return term1Coeff*term2Coeff*betaVal*integral;
334 
335  }
336 
337  // projections onto Phi1, integrate all other angles
338  case 5:
339  {
340  integral = 0.;
341  integral +=
342  (16*f00*Pi())/27.;
343  integral +=
344  (16*fpp*Pi())/27.;
345  integral +=
346  (16*fmm*Pi())/27.;
347  integral +=
348  (16*fp0*Pi())/27.;
349  integral +=
350  (16*f0m*Pi())/27.;
351  integral +=
352  (16*f0p*Pi())/27.;
353  integral +=
354  (16*fm0*Pi())/27.;
355  integral +=
356  (8*Sqrt(fm0)*Sqrt(fp0)*Pi()*Cos(2.*Phi1 + phim0 - phip0))/27.;
357  return term1Coeff*term2Coeff*betaVal*integral;
358 
359  }
360  // projection to Phi, integrate all other angles
361  case 4:
362  {
363  integral = 0.;
364  integral +=
365  (16*f00*Pi())/27.;
366  integral +=
367  (16*fmm*Pi())/27.;
368  integral +=
369  (16*fpp*Pi())/27.;
370  integral +=
371  (16*fp0*Pi())/27.;
372  integral +=
373  (16*f0m*Pi())/27.;
374  integral +=
375  (16*f0p*Pi())/27.;
376  integral +=
377  (16*fm0*Pi())/27.;
378  integral +=
379  (Sqrt(f00)*Sqrt(fmm)*TMath::Power(TMath::Pi(), 3)*R1Val*R2Val*Cos(Phi - phimm))/12.;
380  integral +=
381  (Sqrt(f00)*Sqrt(fpp)*TMath::Power(TMath::Pi(), 3)*R1Val*R2Val*Cos(Phi + phipp))/12.;
382  integral +=
383  (8*Sqrt(fmm)*Sqrt(fpp)*TMath::Pi()*Cos(2*Phi - phimm + phipp))/27.;
384  integral +=
385  (Sqrt(f0m)*Sqrt(fp0)*TMath::Power(Pi(), 3)*R1Val*R2Val*
386  Cos(Phi - phi0m + phip0))/12.;
387  integral +=
388  (Sqrt(f0p)*Sqrt(fm0)*TMath::Power(Pi(), 3)*R1Val*R2Val*
389  Cos(Phi + phi0p - phim0))/12.;
390  return term1Coeff*term2Coeff*betaVal*integral;
391  }
392 
393  // projections to h2, integrate over all others
394  case 3:
395  {
396  integral = 0.;
397  integral +=
398  (-8*f00*(-1 + TMath::Power(h2, 2))*TMath::Power(TMath::Pi(), 2))/9.;
399  integral +=
400  (4*fmm*TMath::Power(TMath::Pi(), 2)*(1 + TMath::Power(h2, 2) - 2*h2*R2Val))/9.;
401  integral +=
402  (4*fpp*TMath::Power(TMath::Pi(), 2)*(1 + TMath::Power(h2, 2) + 2*h2*R2Val))/9.;
403  integral +=
404  (-8*fp0*(-1 + TMath::Power(h2, 2))*TMath::Power(TMath::Pi(), 2))/9.;
405  integral +=
406  (4*f0m*TMath::Power(TMath::Pi(), 2)*(1 + TMath::Power(h2, 2) - 2*h2*R2Val))/9.;
407  integral +=
408  (4*f0p*TMath::Power(TMath::Pi(), 2)*(1 + TMath::Power(h2, 2) + 2*h2*R2Val))/9.;
409  integral +=
410  (-8*fm0*(-1 + TMath::Power(h2, 2))*TMath::Power(TMath::Pi(), 2))/9.;
411  return term1Coeff*term2Coeff*betaVal*integral;
412  }
413  // projections to h1, integrate all others
414  case 2:
415  {
416  integral=0.;
417  integral +=
418  (-8*f00*(-1 + TMath::Power(h1, 2))*TMath::Power(TMath::Pi(), 2))/9.;
419  integral +=
420  (4*fmm*TMath::Power(TMath::Pi(), 2)*(1 + TMath::Power(h1, 2) - 2*h1*R1Val))/9.;
421  integral +=
422  (4*fpp*TMath::Power(TMath::Pi(), 2)*(1 + TMath::Power(h1, 2) + 2*h1*R1Val))/9.;
423  integral +=
424  (4*fp0*TMath::Power(TMath::Pi(), 2)*(1 + TMath::Power(h1, 2) + 2*h1*R1Val))/9.;
425  integral +=
426  (-8*f0m*(-1 + TMath::Power(h1, 2))*TMath::Power(TMath::Pi(), 2))/9.;
427  integral +=
428  (-8*f0p*(-1 + TMath::Power(h1, 2))*TMath::Power(TMath::Pi(), 2))/9.;
429  integral +=
430  (4*fm0*TMath::Power(TMath::Pi(), 2)*(1 + TMath::Power(h1, 2) - 2*h1*R1Val))/9.;
431  return term1Coeff*term2Coeff*betaVal*integral;
432  }
433 
434  // projections to hs, integrate all others
435  case 1:
436  {
437  integral = 0.;
438  integral +=
439  (-8*f00*(-1 + TMath::Power(hs_, 2))*TMath::Power(TMath::Pi(), 2))/9.;
440  integral +=
441  (-8*fmm*(-1 + TMath::Power(hs_, 2))*TMath::Power(TMath::Pi(), 2))/9.;
442  integral +=
443  (-8*fpp*(-1 + TMath::Power(hs_, 2))*TMath::Power(TMath::Pi(), 2))/9.;
444  integral +=
445  (4.*fp0*TMath::Power(Pi(), 2)*(1. + TMath::Power(hs_, 2)))/9.;
446  integral +=
447  (4.*f0m*TMath::Power(Pi(), 2)*(1. + TMath::Power(hs_, 2)))/9.;
448  integral +=
449  (4.*f0p*TMath::Power(Pi(), 2)*(1. + TMath::Power(hs_, 2)))/9.;
450  integral +=
451  (4.*fm0*TMath::Power(Pi(), 2)*(1. + TMath::Power(hs_, 2)))/9.;
452  return term1Coeff*term2Coeff*betaVal*integral;
453 
454  }
455  }
456  assert(0);
457  return 0;
458 }

◆ clone()

virtual TObject* RooSpinOne_7D::clone ( const char *  newname) const
inlinevirtual

Definition at line 38 of file RooSpinOne_7D.h.

38 { return new RooSpinOne_7D(*this, newname); }

◆ evaluate()

Double_t RooSpinOne_7D::evaluate ( ) const
protected

Definition at line 78 of file RooSpinOne_7D.cc.

79 {
80 
81 
82  bool isZZ = true;
83  if ( mZ < 90.) isZZ = false;
84  if ( isZZ ) {
85  if( (m1+m2) > mzz || m2>m1 ) return 1e-9;
86  } else {
87  if( (m1+m2) > mzz ) return 1e-9;
88  }
89  double nanval = sqrt((1 - TMath::Power(m1 - m2,2)/TMath::Power(mzz,2))*(1 - TMath::Power(m1 + m2,2)/TMath::Power(mzz,2)));
90  if (nanval != nanval) return 1e-9;
91 
92  //-----------------------------------------------------------------------
93  // propagator
94  //-----------------------------------------------------------------------
95 
96  Double_t betaValSquared = (1.-(pow(m1-m2,2)/pow(mzz,2)))*(1.-(pow(m1+m2,2)/pow(mzz,2)));
97  Double_t betaVal = sqrt(betaValSquared);
98 
99  Double_t term1Coeff = (pow(m1,3))/( (pow(m1,2)-pow(mZ,2))*(pow(m1,2)-pow(mZ,2))+pow(mZ,2)*pow(gamZ,2) );
100  Double_t term2Coeff = (pow(m2,3))/( (pow(m2,2)-pow(mZ,2))*(pow(m2,2)-pow(mZ,2))+pow(mZ,2)*pow(gamZ,2) );
101 
102  //-----------------------------------------------------------------------
103  // Helicity Amplitudes
104  //-----------------------------------------------------------------------
105 
106  // calculating the angular parameters from the coupling constants
107  // See http://www.pha.jhu.edu/~gritsan/FORM/result_spin1.txt
108  Double_t x = (mzz*mzz-m1*m1-m2*m2)/(2.0*m1*m2);
109 
110  Double_t hs_ = -hs;
111 
112  Double_t f00Real = TMath::Power(mzz,-1)*g1Val*( sqrt(x*x-1) )*(m1*m1-m2*m2);
113  Double_t f00Imag = 0.;
114 
115  Double_t fppReal = 0.;
116  Double_t fppImag = TMath::Power(mzz,-1)*g2Val*(m1*m1-m2*m2);
117 
118  Double_t fmmReal = 0.;
119  Double_t fmmImag = -1. * fppImag;
120 
121 
122  Double_t fp0Real = m1*g1Val * ( sqrt(x*x-1) );
123  Double_t fp0Imag = Power(mzz,-2.)*Power(m2,3)*g2Val * ( - 1./2. )
124  + Power(mzz,-2.)*m1*m1*m2*g2Val * ( 1 + 2*(x*x-1) )
125  + Power(mzz,-2.)*Power(m1,4)*Power(m2,-1)*g2Val * ( - 1./2. )
126  + m2*g2Val * ( - 1./2. )
127  + m1*m1*Power(m2,-1)*g2Val * ( 1./2. );
128 
129 
130  Double_t f0pReal = m2*g1Val * ( - sqrt(x*x-1) );
131  Double_t f0pImag = Power(mzz,-2.)*Power(m1,-1)*Power(m2,4)*g2Val * ( 1./2. )
132  + Power(mzz,-2.)*m1*Power(m2,2)*g2Val * ( - 1 - 2*(x*x-1) )
133  + Power(mzz,-2.)*Power(m1,3)*g2Val * ( 1./2. )
134  + Power(m1,-1)*Power(m2,2)*g2Val * ( - 1./2. )
135  + m1*g2Val * ( 1./2. );
136 
137  Double_t f0mReal = m2*g1Val * ( - sqrt(x*x-1) );
138  Double_t f0mImag = Power(mzz,-2.)*Power(m1,-1)*Power(m2,4)*g2Val * ( - 1./2. )
139  + Power(mzz,-2.)*m1*Power(m2,2)*g2Val * ( 1 + 2*(x*x-1) )
140  + Power(mzz,-2.)*Power(m1,3)*g2Val * ( - 1./2. )
141  + Power(m1,-1)*Power(m2,2)*g2Val * ( 1./2. )
142  + m1*g2Val * ( - 1./2. );
143 
144  Double_t fm0Real = m1*g1Val * ( sqrt(x*x-1) );
145  Double_t fm0Imag = Power(mzz,-2.)*Power(m2,3)*g2Val * ( 1./2. )
146  + Power(mzz,-2.)*m1*m1*m2*g2Val * ( - 1 - 2*(x*x-1) )
147  + Power(mzz,-2.)*Power(m1,4)*Power(m2,-1)*g2Val * ( 1./2. )
148  + m2*g2Val * ( 1./2. )
149  + m1*m1*Power(m2,-1)*g2Val * ( - 1./2. );
150 
151  Double_t f00 = f00Imag*f00Imag + f00Real*f00Real;
152  Double_t fpp = fppImag*fppImag + fppReal*fppReal;
153  Double_t fmm = fmmImag*fmmImag + fmmReal*fmmReal;
154  Double_t fp0 = fp0Imag*fp0Imag + fp0Real*fp0Real;
155  Double_t f0p = f0pImag*f0pImag + f0pReal*f0pReal;
156  Double_t fm0 = fm0Imag*fm0Imag + fm0Real*fm0Real;
157  Double_t f0m = f0mImag*f0mImag + f0mReal*f0mReal;
158 
159  Double_t phi00=atan2(f00Imag,f00Real);
160  Double_t phipp=atan2(fppImag,fppReal)-phi00;
161  Double_t phimm=atan2(fmmImag,fmmReal)-phi00;
162  Double_t phip0=atan2(fp0Imag,fp0Real)-phi00;
163  Double_t phi0p=atan2(f0pImag,f0pReal)-phi00;
164  Double_t phim0=atan2(fm0Imag,fm0Real)-phi00;
165  Double_t phi0m=atan2(f0mImag,f0mReal)-phi00;
166 
167  Double_t value=0;
168 
169  value += -(f00*(-1 + TMath::Power(h1,2))*(-1 + TMath::Power(h2,2))*(-1 + TMath::Power(hs_,2)))/8.;
170  value += -(fmm*(-1 + TMath::Power(hs_,2))*(1 + TMath::Power(h1,2) - 2*h1*R1Val)*(1 + TMath::Power(h2,2) - 2*h2*R2Val))/32.;
171  value += -(fpp*(-1 + TMath::Power(hs_,2))*(1 + TMath::Power(h1,2) + 2*h1*R1Val)*(1 + TMath::Power(h2,2) + 2*h2*R2Val))/32.;
172 
173  value += -(fp0*(-1 + TMath::Power(h2,2))*(1 + TMath::Power(hs_,2))*(1 + TMath::Power(h1,2) + 2*h1*R1Val))/32.;
174  value += -(f0m*(-1 + TMath::Power(h1,2))*(1 + TMath::Power(hs_,2))*(1 + TMath::Power(h2,2) - 2*h2*R2Val))/32.;
175  value += -(f0p*(-1 + TMath::Power(h1,2))*(1 + TMath::Power(hs_,2))*(1 + TMath::Power(h2,2) + 2*h2*R2Val))/32.;
176  value += -(fm0*(-1 + TMath::Power(h2,2))*(1 + TMath::Power(hs_,2))*(1 + TMath::Power(h1,2) - 2*h1*R1Val))/32.; // fix
177 
178  value += -(Sqrt(f00)*Sqrt(fmm)*Sqrt(1 - TMath::Power(h1,2))*Sqrt(1 - TMath::Power(h2,2))*(-1 + TMath::Power(hs_,2))*(h1 - R1Val)*(h2 - R2Val)*Cos(Phi - phimm))/8.;
179  value += -(Sqrt(f00)*Sqrt(fpp)*Sqrt(1 - TMath::Power(h1,2))*Sqrt(1 - TMath::Power(h2,2))*(-1 + TMath::Power(hs_,2))*(h1 + R1Val)*(h2 + R2Val)*Cos(Phi + phipp))/8.;
180  value += (Sqrt(f00)*Sqrt(fp0)*Sqrt(1 - TMath::Power(h1,2))*(-1 + TMath::Power(h2,2))*hs_*Sqrt(1 - TMath::Power(hs_,2))*(h1 + R1Val)*Cos(Phi1 - phip0))/8.;
181  value += -(Sqrt(f00)*Sqrt(f0m)*(-1 + TMath::Power(h1,2))*Sqrt(1 - TMath::Power(h2,2))*hs_*Sqrt(1 - TMath::Power(hs_,2))*(h2 - R2Val)*Cos(Phi - phi0m + Phi1))/8.;
182  value += -(Sqrt(f00)*Sqrt(f0p)*(-1 + TMath::Power(h1,2))*Sqrt(1 - TMath::Power(h2,2))*hs_*Sqrt(1 - TMath::Power(hs_,2))*(h2 + R2Val)*Cos(Phi + phi0p + Phi1))/8.;
183  value += (Sqrt(f00)*Sqrt(fm0)*Sqrt(1 - TMath::Power(h1,2))*(-1 + TMath::Power(h2,2))*hs_*Sqrt(1 - TMath::Power(hs_,2))*(h1 - R1Val)*Cos(Phi1 + phim0))/8.; // fix
184 
185  value += -(Sqrt(fmm)*Sqrt(fpp)*(-1 + TMath::Power(h1,2))*(-1 + TMath::Power(h2,2))*(-1 + TMath::Power(hs_,2))*Cos(2*Phi - phimm + phipp))/16.;
186  value += -(Sqrt(fmm)*Sqrt(fp0)*(-1 + TMath::Power(h1,2))*Sqrt(1 - TMath::Power(h2,2))*hs_*Sqrt(1 - TMath::Power(hs_,2))*(h2 - R2Val)*Cos(Phi - Phi1 - phimm + phip0))/16.;
187  value += (Sqrt(f0m)*Sqrt(fmm)*Sqrt(1 - TMath::Power(h1,2))*hs_*Sqrt(1 - TMath::Power(hs_,2))*(h1 - R1Val)*(1 + TMath::Power(h2,2) - 2*h2*R2Val)*Cos(phi0m - Phi1 - phimm))/16.;
188  value += (Sqrt(f0p)*Sqrt(fmm)*Sqrt(1 - TMath::Power(h1,2))*(-1 + TMath::Power(h2,2))*hs_*Sqrt(1 - TMath::Power(hs_,2))*(h1 - R1Val)*Cos(2*Phi + phi0p + Phi1 - phimm))/16.; // fix
189  value += -(Sqrt(fm0)*Sqrt(fmm)*Sqrt(1 - TMath::Power(h2,2))*hs_*Sqrt(1 - TMath::Power(hs_,2))*(1 + TMath::Power(h1,2) - 2*h1*R1Val)*(h2 - R2Val)*Cos(Phi + Phi1 + phim0 - phimm))/16.;
190  value += -(Sqrt(fp0)*Sqrt(fpp)*Sqrt(1 - TMath::Power(h2,2))*hs_*Sqrt(1 - TMath::Power(hs_,2))*(1 + TMath::Power(h1,2) + 2*h1*R1Val)*(h2 + R2Val)*Cos(Phi + Phi1 - phip0 + phipp))/16.;
191  value += (Sqrt(f0m)*Sqrt(fpp)*Sqrt(1 - TMath::Power(h1,2))*(-1 + TMath::Power(h2,2))*hs_*Sqrt(1 - TMath::Power(hs_,2))*(h1 + R1Val)*Cos(2*Phi - phi0m + Phi1 + phipp))/16.;
192  value += (Sqrt(f0p)*Sqrt(fpp)*Sqrt(1 - TMath::Power(h1,2))*hs_*Sqrt(1 - TMath::Power(hs_,2))*(h1 + R1Val)*(1 + TMath::Power(h2,2) + 2*h2*R2Val)*Cos(phi0p + Phi1 - phipp))/16.;
193  value += -(Sqrt(fm0)*Sqrt(fpp)*(-1 + TMath::Power(h1,2))*Sqrt(1 - TMath::Power(h2,2))*hs_*Sqrt(1 - TMath::Power(hs_,2))*(h2 + R2Val)*Cos(Phi - Phi1 - phim0 + phipp))/16.;
194 
195  value += -(Sqrt(f0m)*Sqrt(fp0)*Sqrt(1 - TMath::Power(h1,2))*Sqrt(1 - TMath::Power(h2,2))*
196  (1 + TMath::Power(hs_,2))*(h1 + R1Val)*(h2 - R2Val)*Cos(Phi - phi0m + phip0))/16.;
197  value += -(Sqrt(f0p)*Sqrt(fp0)*Sqrt(1 - TMath::Power(h1,2))*Sqrt(1 - TMath::Power(h2,2))*
198  (-1 + TMath::Power(hs_,2))*(h1 + R1Val)*(h2 + R2Val)*
199  Cos(Phi + phi0p + 2*Phi1 - phip0))/16.;
200  value += -(Sqrt(fm0)*Sqrt(fp0)*(-1 + TMath::Power(h1,2))*(-1 + TMath::Power(h2,2))*
201  (-1 + TMath::Power(hs_,2))*Cos(2*Phi1 + phim0 - phip0))/16.;
202  value += -(Sqrt(f0m)*Sqrt(f0p)*(-1 + TMath::Power(h1,2))*(-1 + TMath::Power(h2,2))*
203  (-1 + TMath::Power(hs_,2))*Cos(2*Phi - phi0m + phi0p + 2*Phi1))/16.;
204  value += -(Sqrt(f0m)*Sqrt(fm0)*Sqrt(1 - TMath::Power(h1,2))*Sqrt(1 - TMath::Power(h2,2))*
205  (-1 + TMath::Power(hs_,2))*(h1 - R1Val)*(h2 - R2Val)*
206  Cos(Phi - phi0m + 2*Phi1 + phim0))/16.;
207  value += -(Sqrt(f0p)*Sqrt(fm0)*Sqrt(1 - TMath::Power(h1,2))*Sqrt(1 - TMath::Power(h2,2))*
208  (1 + TMath::Power(hs_,2))*(h1 - R1Val)*(h2 + R2Val)*Cos(Phi + phi0p - phim0))/16.;
209 
210  return value*term1Coeff*term2Coeff*betaVal;
211 
212 }

◆ getAnalyticalIntegral()

Int_t RooSpinOne_7D::getAnalyticalIntegral ( RooArgSet &  allVars,
RooArgSet &  analVars,
const char *  rangeName = 0 
) const

Definition at line 214 of file RooSpinOne_7D.cc.

214  {
215  if (matchArgs(allVars, analVars, RooArgSet(*hs.absArg(), *h1.absArg(), *h2.absArg(), *Phi.absArg(), *Phi1.absArg()))) return 6;
216  if (matchArgs(allVars, analVars, hs, h1, h2, Phi)) return 5;
217  if (matchArgs(allVars, analVars, hs, h1, h2, Phi1)) return 4;
218  if (matchArgs(allVars, analVars, hs, h1, Phi, Phi1)) return 3;
219  if (matchArgs(allVars, analVars, hs, h2, Phi, Phi1)) return 2;
220  if (matchArgs(allVars, analVars, h1, h2, Phi, Phi1)) return 1;
221  return 0;
222 }

◆ setZZ4fOrdering()

void RooSpinOne_7D::setZZ4fOrdering ( Bool_t  flag = true)

Definition at line 460 of file RooSpinOne_7D.cc.

460 { ZZ4fOrdering=flag; }

Member Data Documentation

◆ aParam

RooRealProxy RooSpinOne_7D::aParam
protected

Definition at line 59 of file RooSpinOne_7D.h.

◆ g1Val

RooRealProxy RooSpinOne_7D::g1Val
protected

Definition at line 55 of file RooSpinOne_7D.h.

◆ g2Val

RooRealProxy RooSpinOne_7D::g2Val
protected

Definition at line 56 of file RooSpinOne_7D.h.

◆ gamZ

RooRealProxy RooSpinOne_7D::gamZ
protected

Definition at line 61 of file RooSpinOne_7D.h.

◆ h1

RooRealProxy RooSpinOne_7D::h1
protected

Definition at line 50 of file RooSpinOne_7D.h.

◆ h2

RooRealProxy RooSpinOne_7D::h2
protected

Definition at line 51 of file RooSpinOne_7D.h.

◆ hs

RooRealProxy RooSpinOne_7D::hs
protected

Definition at line 52 of file RooSpinOne_7D.h.

◆ m1

RooRealProxy RooSpinOne_7D::m1
protected

Definition at line 48 of file RooSpinOne_7D.h.

◆ m2

RooRealProxy RooSpinOne_7D::m2
protected

Definition at line 49 of file RooSpinOne_7D.h.

◆ mZ

RooRealProxy RooSpinOne_7D::mZ
protected

Definition at line 60 of file RooSpinOne_7D.h.

◆ mzz

RooRealProxy RooSpinOne_7D::mzz
protected

Definition at line 47 of file RooSpinOne_7D.h.

◆ Phi

RooRealProxy RooSpinOne_7D::Phi
protected

Definition at line 53 of file RooSpinOne_7D.h.

◆ Phi1

RooRealProxy RooSpinOne_7D::Phi1
protected

Definition at line 54 of file RooSpinOne_7D.h.

◆ R1Val

RooRealProxy RooSpinOne_7D::R1Val
protected

Definition at line 57 of file RooSpinOne_7D.h.

◆ R2Val

RooRealProxy RooSpinOne_7D::R2Val
protected

Definition at line 58 of file RooSpinOne_7D.h.

◆ ZZ4fOrdering

Bool_t RooSpinOne_7D::ZZ4fOrdering
protected

Definition at line 63 of file RooSpinOne_7D.h.


The documentation for this class was generated from the following files:
value
pymela::gHIGGS_KAPPA value("gHIGGS_KAPPA_TILDE", pymela::gHIGGS_KAPPA_TILDE) .value("SIZE_HQQ"
RooSpinOne_7D::aParam
RooRealProxy aParam
Definition: RooSpinOne_7D.h:59
RooSpinOne_7D::m2
RooRealProxy m2
Definition: RooSpinOne_7D.h:49
RooSpinOne_7D::ZZ4fOrdering
Bool_t ZZ4fOrdering
Definition: RooSpinOne_7D.h:63
RooSpinOne_7D::h2
RooRealProxy h2
Definition: RooSpinOne_7D.h:51
RooSpinOne_7D::g2Val
RooRealProxy g2Val
Definition: RooSpinOne_7D.h:56
RooSpinOne_7D::gamZ
RooRealProxy gamZ
Definition: RooSpinOne_7D.h:61
RooSpinOne_7D::R2Val
RooRealProxy R2Val
Definition: RooSpinOne_7D.h:58
RooSpinOne_7D::R1Val
RooRealProxy R1Val
Definition: RooSpinOne_7D.h:57
RooSpinOne_7D::h1
RooRealProxy h1
Definition: RooSpinOne_7D.h:50
RooSpinOne_7D::mzz
RooRealProxy mzz
Definition: RooSpinOne_7D.h:47
RooSpinOne_7D::m1
RooRealProxy m1
Definition: RooSpinOne_7D.h:48
RooSpinOne_7D::hs
RooRealProxy hs
Definition: RooSpinOne_7D.h:52
RooSpinOne_7D::Phi
RooRealProxy Phi
Definition: RooSpinOne_7D.h:53
RooSpinOne_7D::g1Val
RooRealProxy g1Val
Definition: RooSpinOne_7D.h:55
RooSpinOne_7D::mZ
RooRealProxy mZ
Definition: RooSpinOne_7D.h:60
RooSpinOne_7D::Phi1
RooRealProxy Phi1
Definition: RooSpinOne_7D.h:54
RooSpinOne_7D::RooSpinOne_7D
RooSpinOne_7D()
Definition: RooSpinOne_7D.h:20