JHUGen MELA  JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
RooSpinTwo_7DComplex_ppHVV.cc
Go to the documentation of this file.
2 
3 
4 using namespace std;
5 using namespace MELAStreamHelpers;
6 
7 
10  const char *name, const char *title,
11  modelMeasurables const& _measurables,
12  modelParameters const& _parameters,
13  modelCouplings const& _couplings,
14  RooSpin::VdecayType _Vdecay1, RooSpin::VdecayType _Vdecay2,
15  TVar::VerbosityLevel verbosity_
16 ) : RooSpinTwo(
17  name, title,
18  _measurables,
19  _parameters,
20  _couplings,
21  _Vdecay1, _Vdecay2,
22  verbosity_
23 ),
24 ZZ4fOrdering(true)
25 {}
26 
27 
29  const RooSpinTwo_7DComplex_ppHVV& other, const char* name
30 ) : RooSpinTwo(other, name),
31 ZZ4fOrdering(other.ZZ4fOrdering)
32 {}
33 
34 
35 Double_t RooSpinTwo_7DComplex_ppHVV::evaluateH1Factor(Int_t i1, Int_t j1, Int_t helicity, Int_t code) const{
36  const Double_t Pi = TMath::Pi();
37  Double_t dHel = (Double_t)helicity;
38  Double_t result = 0;
39  if ((code % prime_h1)==0){
40  if ((i1==0 && j1==0) || (i1==-1 && j1==1) || (i1==1 && j1==-1)) result = 4./3.; // 15 amps
41  else if (i1==1 && j1==1) result = 8./3.; // 6 amps
42  else if (i1==-1 && j1==-1) result = 8./3.; // 6 amps
43  else if ((i1==0 && j1==1) || (i1==1 && j1==0)) result = Pi*dHel/2.; // 9 amps
44  else if ((i1==0 && j1==-1) || (i1==-1 && j1==0)) result = Pi*dHel/2.; // 9 amps
45  }
46  else{
47  if ((i1==0 && j1==0) || (i1==-1 && j1==1) || (i1==1 && j1==-1)) result = 1.-pow(h1, 2); // 15 amps
48  else if (i1==1 && j1==1) result = 1.+pow(h1, 2)-2.*h1*dHel; // 6 amps
49  else if (i1==-1 && j1==-1) result = 1.+pow(h1, 2)+2.*h1*dHel; // 6 amps
50  else if ((i1==0 && j1==1) || (i1==1 && j1==0)) result = sqrt(fabs(1.-pow(h1, 2)))*(dHel-h1); // 9 amps
51  else if ((i1==0 && j1==-1) || (i1==-1 && j1==0)) result = sqrt(fabs(1.-pow(h1, 2)))*(dHel+h1); // 9 amps
52  }
53  return result;
54 }
55 Double_t RooSpinTwo_7DComplex_ppHVV::evaluateH2Factor(Int_t i2, Int_t j2, Int_t helicity, Int_t code) const{
56  const Double_t Pi = TMath::Pi();
57  Double_t dHel = (Double_t)helicity;
58  Double_t result = 0;
59  if ((code % prime_h2)==0){
60  if ((i2==0 && j2==0) || (i2==-1 && j2==1) || (i2==1 && j2==-1)) result = 4./3.; // 15 amps
61  else if (i2==1 && j2==1) result = 8./3.; // 6 amps
62  else if (i2==-1 && j2==-1) result = 8./3.; // 6 amps
63  else if ((i2==0 && j2==1) || (i2==1 && j2==0)) result = Pi*dHel/2.; // 9 amps
64  else if ((i2==0 && j2==-1) || (i2==-1 && j2==0)) result = Pi*dHel/2.; // 9 amps
65  }
66  else{
67  if ((i2==0 && j2==0) || (i2==-1 && j2==1) || (i2==1 && j2==-1)) result = 1.-pow(h2, 2); // 15 amps
68  else if (i2==1 && j2==1) result = 1.+pow(h2, 2)-2.*h2*dHel; // 6 amps
69  else if (i2==-1 && j2==-1) result = 1.+pow(h2, 2)+2.*h2*dHel; // 6 amps
70  else if ((i2==0 && j2==1) || (i2==1 && j2==0)) result = sqrt(fabs(1.-pow(h2, 2)))*(dHel-h2); // 9 amps
71  else if ((i2==0 && j2==-1) || (i2==-1 && j2==0)) result = sqrt(fabs(1.-pow(h2, 2)))*(dHel+h2); // 9 amps
72  }
73  return result;
74 }
75 Double_t RooSpinTwo_7DComplex_ppHVV::evaluateHSFactor(Int_t di, Int_t dj, Int_t code) const{
76  Double_t f_spinz0 = 1. - f_spinz1 - f_spinz2;
77  if (f_spinz0<0.) f_spinz0=0;
78  Double_t hsneg = -hs; if (fabs(hsneg)>1.) hsneg *= 1./fabs(hsneg);
79 
80  Double_t AF200 = 0;
81  Double_t AF201 = 0;
82  Double_t AF202 = 0;
83  Double_t AF211 = 0;
84  Double_t AF2m11 = 0;
85  Double_t AF212 = 0;
86  Double_t AF2m12 = 0;
87  Double_t AF222 = 0;
88  Double_t AF2m22 = 0;
89  if ((code % prime_hs)==0){
90  AF200 = 2.*((2.*f_spinz0 + 3.*f_spinz2) - 6.*(2.*f_spinz0 - 2.*f_spinz1 + f_spinz2)/3. +3.*(6.*f_spinz0 - 4.*f_spinz1 + f_spinz2)/5.)/8.; // 6
91  //AF201 = 0; // 12
92  AF202 = ((2.*f_spinz0 - f_spinz2)*(4./3.) - (6.*f_spinz0 - 4.*f_spinz1 + f_spinz2)*(4./15.))*(-sqrt(3./2.)/8.); // 6
93 
94  AF211 = 2.*((f_spinz1 + f_spinz2) + 3.*(2.*f_spinz0 - f_spinz1)/3. - (6.*f_spinz0 - 4.*f_spinz1 + f_spinz2)/5.)/4.; // 6
95  AF2m11 = ((f_spinz1 - f_spinz2)*(4./3.) + (6.*f_spinz0 - 4.*f_spinz1 + f_spinz2)*(4./15.))/(-4.); // 4
96 
97  //AF212 = 0; // 4
98  //AF2m12 = 0; // 4
99 
100  AF222 = 2.*((6.*f_spinz0 + 4.*f_spinz1 + f_spinz2) - 6.*(2.*f_spinz0 - f_spinz2)/3. + (6.*f_spinz0 - 4.*f_spinz1 + f_spinz2)/5.)/16.; // 2
101  AF2m22 = 16./15.*(6*f_spinz0 - 4*f_spinz1 + f_spinz2)/16.; // 1
102  }
103  else{
104  AF200 = ((2.*f_spinz0 + 3.*f_spinz2) - 6.*(2.*f_spinz0 - 2.*f_spinz1 + f_spinz2)*pow(hsneg, 2) +3.*(6.*f_spinz0 - 4.*f_spinz1 + f_spinz2)*pow(hsneg, 4))/8.; // 6 // F(2)00
105  AF201 = (hsneg*sqrt(1 - pow(hsneg, 2))*((2.*f_spinz0 - 2.*f_spinz1 + f_spinz2) - (6.*f_spinz0 - 4.*f_spinz1 + f_spinz2)*pow(hsneg, 2)))*(-sqrt(6.)/8.); // 12 // F(2)01
106  AF202 = ((1. - pow(hsneg, 2))*((2.*f_spinz0 - f_spinz2) - (6.*f_spinz0 - 4.*f_spinz1 + f_spinz2)*pow(hsneg, 2)))*(-sqrt(3./2.)/8.); // 6 // F(2)02
107 
108  AF211 = ((f_spinz1 + f_spinz2) + 3.*(2.*f_spinz0 - f_spinz1)*pow(hsneg, 2) - (6.*f_spinz0 - 4.*f_spinz1 + f_spinz2)*pow(hsneg, 4))/4.; // 6 // F(2)11
109  AF2m11 = ((1. - pow(hsneg, 2))*((f_spinz1 - f_spinz2) + (6.*f_spinz0 - 4.*f_spinz1 + f_spinz2)*pow(hsneg, 2)))/(-4.); // 4 // F(2)-11
110 
111  AF212 = (hsneg*sqrt(1 - pow(hsneg, 2))*((6.*f_spinz0 - 3.*f_spinz2) - (6.*f_spinz0 - 4.*f_spinz1 + f_spinz2)*pow(hsneg, 2)))/(-8.); // 4 // F(2)12
112  AF2m12 = ((6*f_spinz0 - 4*f_spinz1 + f_spinz2)*hsneg*pow(1. - pow(hsneg, 2), 1.5))/(-8.); // 4 // F(2)-12
113 
114  AF222 = ((6.*f_spinz0 + 4.*f_spinz1 + f_spinz2) - 6.*(2.*f_spinz0 - f_spinz2)*pow(hsneg, 2) + (6.*f_spinz0 - 4.*f_spinz1 + f_spinz2)*pow(hsneg, 4))/16.; // 2 // F(2)22
115  AF2m22 = ((6*f_spinz0 - 4*f_spinz1 + f_spinz2)*pow(1. - pow(hsneg, 2), 2))/16.; // 1 // F(2)-22
116  }
117 
118  Double_t result = 0;
119  if (di==0 && dj==0) result = AF200;
120  else if ((di==1 && dj==1) || (di==-1 && dj==-1)) result = AF211;
121  else if ((di==2 && dj==2) || (di==-2 && dj==-2)) result = AF222;
122 
123  else if ((di==1 && dj==0) || (di==0 && dj==1)) result = AF201;
124  else if ((di==-1 && dj==0) || (di==0 && dj==-1)) result = -AF201;
125  else if ((di==1 && dj==2) || (di==2 && dj==1)) result = AF212;
126  else if ((di==-1 && dj==-2) || (di==-2 && dj==-1)) result = -AF212;
127 
128  else if ((di==2 && dj==0) || (di==0 && dj==2) || (di==-2 && dj==0) || (di==0 && dj==-2)) result = AF202;
129  else if ((di==1 && dj==-1) || (di==-1 && dj==1)) result = AF2m11;
130 
131  else if ((di==-1 && dj==2) || (di==2 && dj==-1)) result = AF2m12;
132  else if ((di==-2 && dj==1) || (di==1 && dj==-2)) result = -AF2m12;
133 
134  else if ((di==-2 && dj==2) || (di==2 && dj==-2)) result = AF2m22;
135  return result;
136 }
137 Double_t RooSpinTwo_7DComplex_ppHVV::evaluatePhi1PhiFactor(Int_t i1, Int_t i2, Int_t j1, Int_t j2, Int_t code, Double_t extraPhase1, Double_t extraPhase2) const{
138  const Double_t Pi = TMath::Pi();
139 
140  Double_t result = 0;
141  Double_t phase = 0;
142  Double_t phasePhi = 0;
143  Double_t phasePhi1 = 0;
144 
145  phasePhi1 += -i1;
146  phasePhi1 += i2;
147  phasePhi += i2;
148 
149  phasePhi1 += j1;
150  phasePhi1 += -j2;
151  phasePhi += -j2;
152 
153  if ((code % prime_Phi)==0 && (code % prime_Phi1)==0){
154  if (i1==j1 && i2==j2) result = 4.*pow(Pi, 2);
155  // Everything else is 0!
156  }
157  else if ((code % prime_Phi)==0){
158  if (i1==j1 && i2==j2) result = 2.*Pi;
159  else if (phasePhi==0.){
160  phase = Phi1*phasePhi1+extraPhase1-extraPhase2;
161  result = cos(phase)*2.*Pi;
162  }
163  }
164  else if ((code % prime_Phi1)==0){
165  if (i1==j1 && i2==j2) result = 2.*Pi;
166  else if (phasePhi1==0.){
167  phase = Phi*phasePhi+extraPhase1-extraPhase2;
168  result = cos(phase)*2.*Pi;
169  }
170  }
171  else{
172  phase = Phi1*phasePhi1+Phi*phasePhi+extraPhase1-extraPhase2;
173  result = cos(phase);
174  }
175  return result;
176 }
177 
178 void RooSpinTwo_7DComplex_ppHVV::evaluatePolarizationTerms(std::vector<Double_t>& Axxyyterm, const Int_t code, bool isGammaV1, bool isGammaV2) const{
179  Double_t R1Val, R2Val;
180  calculateVffR1R2(R1Val, R2Val, isGammaV1, isGammaV2);
181 
182  Double_t
183  A00Re, A00Im,
184  AppRe, AppIm,
185  A0pRe, A0pIm, Ap0Re, Ap0Im,
186  AmmRe, AmmIm, A0mRe, A0mIm, Am0Re, Am0Im,
187  ApmRe, ApmIm, AmpRe, AmpIm;
188 
190  A00Re, A00Im,
191  AppRe, AppIm,
192  A0pRe, A0pIm, Ap0Re, Ap0Im,
193  AmmRe, AmmIm, A0mRe, A0mIm, Am0Re, Am0Im,
194  ApmRe, ApmIm, AmpRe, AmpIm,
195  isGammaV1, isGammaV2
196  );
197 
198  std::vector<Double_t> ARexy, AImxy;
199  std::vector<Int_t> index_x,index_y;
200  ARexy.push_back(AmmRe);
201  ARexy.push_back(Am0Re);
202  ARexy.push_back(AmpRe);
203  ARexy.push_back(A0mRe);
204  ARexy.push_back(A00Re);
205  ARexy.push_back(A0pRe);
206  ARexy.push_back(ApmRe);
207  ARexy.push_back(Ap0Re);
208  ARexy.push_back(AppRe);
209  AImxy.push_back(AmmIm);
210  AImxy.push_back(Am0Im);
211  AImxy.push_back(AmpIm);
212  AImxy.push_back(A0mIm);
213  AImxy.push_back(A00Im);
214  AImxy.push_back(A0pIm);
215  AImxy.push_back(ApmIm);
216  AImxy.push_back(Ap0Im);
217  AImxy.push_back(AppIm);
218  index_x.push_back(-1);
219  index_x.push_back(-1);
220  index_x.push_back(-1);
221  index_x.push_back(0);
222  index_x.push_back(0);
223  index_x.push_back(0);
224  index_x.push_back(1);
225  index_x.push_back(1);
226  index_x.push_back(1);
227  index_y.push_back(-1);
228  index_y.push_back(0);
229  index_y.push_back(1);
230  index_y.push_back(-1);
231  index_y.push_back(0);
232  index_y.push_back(1);
233  index_y.push_back(-1);
234  index_y.push_back(0);
235  index_y.push_back(1);
236 
237  for (unsigned int ii=0; ii<index_x.size(); ii++){
238  Int_t i1 = index_x.at(ii);
239  Int_t i2 = index_y.at(ii);
240  Int_t i12 = i1-i2;
241  Double_t ARexyi = ARexy.at(ii);
242  Double_t AImxyi = AImxy.at(ii);
243 
244  Double_t Axyi = sqrt(ARexyi*ARexyi + AImxyi*AImxyi);
245  Double_t phixyi = atan2(AImxyi, ARexyi);
246 
247  for (unsigned int jj=ii; jj<index_x.size(); jj++){
248  Int_t j1 = index_x.at(jj);
249  Int_t j2 = index_y.at(jj);
250  Int_t j12 = j1-j2;
251  Double_t ARexyj = ARexy.at(jj);
252  Double_t AImxyj = AImxy.at(jj);
253 
254  Double_t Axyj = sqrt(ARexyj*ARexyj + AImxyj*AImxyj);
255  Double_t phixyj = atan2(AImxyj, ARexyj);
256 
257  Double_t globalFactor = 1;
258  if (ii!=jj) globalFactor = 2;
259  Double_t phifactor = evaluatePhi1PhiFactor(i1, i2, j1, j2, code, phixyi, phixyj);
260  if (phifactor==0) continue;
261  Double_t h1h2factor[4]={ 0 };
262  for (int ih1=0; ih1<2; ih1++){ // (LL, LR), (RL, RR)
263  Int_t hh1 = 1-2*ih1; // L/R
264  Double_t h1factor = evaluateH1Factor(i1, j1, hh1, code);
265  for (int ih2=0; ih2<2; ih2++){
266  Int_t hh2 = 1-2*ih2; // L/R
267  Double_t h2factor = evaluateH2Factor(i2, j2, hh2, code);
268 
269  Double_t fraction = (1.+((Double_t)hh1)*R1Val+((Double_t)hh2)*R2Val+((Double_t)hh1)*R1Val*((Double_t)hh2)*R2Val)/4.;
270  h1h2factor[2*ih1+ih2] = fraction*h1factor*h2factor;
271  }
272  }
273  Double_t hsfactor=evaluateHSFactor(i12, j12, code);
274 
275  globalFactor *= Axyi*Axyj*phifactor*hsfactor;
276  Double_t result = 0;
277  for (int ih1=0; ih1<2; ih1++){ // (LL, LR), (RL, RR)
278  for (int ih2=0; ih2<2; ih2++) result += globalFactor*h1h2factor[2*ih1+ih2];
279  }
280  if (result!=0.) Axxyyterm.push_back(result);
281  }
282  }
283  return;
284 }
285 
287  Double_t mV;
288  getMVGamV(&mV);
289  bool isZZ = (mV >= 90.);
290  Double_t epsilon=1e-15;
291  Double_t m1_=m1; if (Vdecay1==RooSpin::kVdecayType_GammaOnshell) m1_=0;
292  Double_t m2_=m2; if (Vdecay2==RooSpin::kVdecayType_GammaOnshell) m2_=0;
293  if (
294  (m1_+m2_)>m12 ||
295  (isZZ && Vdecay1==Vdecay2 && ZZ4fOrdering && fabs(m2_-mV)<fabs(m1_-mV) && Vdecay2!=RooSpin::kVdecayType_GammaOnshell) ||
298  ) return epsilon;
299 
300  Int_t code = intCodeStart;
302  code *= prime_Phi;
306  }
307 
308  Double_t betaValSq = (1.-(pow(m1_-m2_, 2)/pow(m12, 2)))*(1.-(pow(m1_+m2_, 2)/pow(m12, 2)));
309  if (betaValSq<=0.) return epsilon;
310  Double_t betaVal = sqrt(betaValSq);
311 
312  Double_t term1Coeff = 1;
313  Double_t term2Coeff = 1;
314  if (Vdecay1!=RooSpin::kVdecayType_GammaOnshell) term1Coeff = 2.*m1_*GeVunit;
315  if (Vdecay2!=RooSpin::kVdecayType_GammaOnshell) term2Coeff = 2.*m2_*GeVunit;
316 
317  std::vector<Double_t> Axxyyterm;
318  evaluatePolarizationTerms(Axxyyterm, code);
319  Double_t value = 0;
320  for (unsigned int s=0; s<Axxyyterm.size(); s++) value += Axxyyterm.at(s);
321  value *= term1Coeff*term2Coeff*betaVal;
322 
323  if (!(value==value)) MELAout << "Evaluate NaN=" << value << endl;
324  if (value<=0.){
325  MELAout << "Evaluated value<=0: " << value << endl;
326  value=epsilon;
327  }
328  return value;
329 }
330 
331 Int_t RooSpinTwo_7DComplex_ppHVV::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const{
332  Int_t code = intCodeStart;
333  if (checkFundamentalType(h1)){ if (matchArgs(allVars, analVars, h1) || Vdecay1==RooSpin::kVdecayType_GammaOnshell) code *= prime_h1; }
334  if (checkFundamentalType(h2)){ if (matchArgs(allVars, analVars, h2) || Vdecay2==RooSpin::kVdecayType_GammaOnshell) code *= prime_h2; }
335  if (checkFundamentalType(hs)){ if (matchArgs(allVars, analVars, hs)) code *= prime_hs; }
338  if (code==1) code=0;
339  return code;
340 }
341 Double_t RooSpinTwo_7DComplex_ppHVV::analyticalIntegral(Int_t code, const char* /*rangeName*/) const{
342  Double_t mV;
343  getMVGamV(&mV);
344  bool isZZ = (mV >= 90.);
345  Double_t epsilon=1e-10;
346  Double_t m1_=m1; if (Vdecay1==RooSpin::kVdecayType_GammaOnshell) m1_=0;
347  Double_t m2_=m2; if (Vdecay2==RooSpin::kVdecayType_GammaOnshell) m2_=0;
348  if (
349  (m1_+m2_)>m12 ||
350  (isZZ && Vdecay1==Vdecay2 && ZZ4fOrdering && fabs(m2_-mV)<fabs(m1_-mV) && Vdecay2!=RooSpin::kVdecayType_GammaOnshell) ||
353  ) return epsilon;
354 
355  Double_t betaValSq = (1.-(pow(m1_-m2_, 2)/pow(m12, 2)))*(1.-(pow(m1_+m2_, 2)/pow(m12, 2)));
356  if (betaValSq<=0.) return epsilon;
357  Double_t betaVal = sqrt(betaValSq);
358 
359  Double_t term1Coeff = 1;
360  Double_t term2Coeff = 1;
361  if (Vdecay1!=RooSpin::kVdecayType_GammaOnshell) term1Coeff = 2.*m1_*GeVunit;
362  if (Vdecay2!=RooSpin::kVdecayType_GammaOnshell) term2Coeff = 2.*m2_*GeVunit;
363 
364  std::vector<Double_t> Axxyyterm;
365  evaluatePolarizationTerms(Axxyyterm, code);
366  Double_t value = 0;
367  for (unsigned int s=0; s<Axxyyterm.size(); s++) value += Axxyyterm.at(s);
368  value *= term1Coeff*term2Coeff*betaVal;
369 
370  if (!(value==value)){
371  MELAout << "Integral NaN=" << value << " at "
372  << "h1=" << h1 << '\t'
373  << "h2=" << h2 << '\t'
374  << "hs=" << hs << '\t'
375  << "Phi1=" << Phi1 << '\t'
376  << "Phi=" << Phi << '\t'
377  << "m1=" << m1_ << '\t'
378  << "m2=" << m2_ << '\t'
379  << "m12=" << m12 << '\t'
380  << endl;
381  MELAout << "Possible sources:\n"
382  << "betaVal=" << betaVal << '\t'
383  << "term1Coeff=" << term1Coeff << '\t'
384  << "term2Coeff=" << term2Coeff << '\t'
385  << endl;
386  }
387  if (value<=0.){
388  MELAout << "Evaluated integral<=0: " << value << endl;
389  value=epsilon;
390  }
391  return value;
392 }
393 
RooSpinTwo_7DComplex_ppHVV::evaluateHSFactor
Double_t evaluateHSFactor(Int_t di, Int_t dj, Int_t code) const
Definition: RooSpinTwo_7DComplex_ppHVV.cc:75
value
pymela::gHIGGS_KAPPA value("gHIGGS_KAPPA_TILDE", pymela::gHIGGS_KAPPA_TILDE) .value("SIZE_HQQ"
RooSpin::Phi
RooRealProxy Phi
Definition: RooSpin.h:110
TVar::VerbosityLevel
VerbosityLevel
Definition: TVar.hh:47
RooSpinTwo_7DComplex_ppHVV::evaluatePhi1PhiFactor
Double_t evaluatePhi1PhiFactor(Int_t i1, Int_t i2, Int_t j1, Int_t j2, Int_t code, Double_t extraPhase1, Double_t extraPhase2) const
Definition: RooSpinTwo_7DComplex_ppHVV.cc:137
RooSpin::prime_hs
@ prime_hs
Definition: RooSpin.h:41
RooSpinTwo::f_spinz1
RooRealProxy f_spinz1
Definition: RooSpinTwo.h:63
ih1
int ih1
Definition: TMCFM.hh:51
RooSpinTwo::f_spinz2
RooRealProxy f_spinz2
Definition: RooSpinTwo.h:64
RooSpin::m12
RooRealProxy m12
Definition: RooSpin.h:113
RooSpinTwo::calculateAmplitudes
virtual void calculateAmplitudes(Double_t &A00Re, Double_t &A00Im, Double_t &AppRe, Double_t &AppIm, Double_t &A0pRe, Double_t &A0pIm, Double_t &Ap0Re, Double_t &Ap0Im, Double_t &AmmRe, Double_t &AmmIm, Double_t &A0mRe, Double_t &A0mIm, Double_t &Am0Re, Double_t &Am0Im, Double_t &ApmRe, Double_t &ApmIm, Double_t &AmpRe, Double_t &AmpIm, bool isGammaV1=false, bool isGammaV2=false) const
Definition: RooSpinTwo.cc:139
RooSpin::m1
RooRealProxy m1
Definition: RooSpin.h:111
RooSpin::prime_Phi
@ prime_Phi
Definition: RooSpin.h:42
MELAStreamHelpers::MELAout
MELAOutputStreamer MELAout
RooSpin::checkFundamentalType
virtual Bool_t checkFundamentalType(const RooRealProxy &proxy) const
Definition: RooSpin.cc:392
RooSpinTwo_7DComplex_ppHVV::analyticalIntegral
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Definition: RooSpinTwo_7DComplex_ppHVV.cc:341
RooSpin::prime_h2
@ prime_h2
Definition: RooSpin.h:40
RooSpinTwo_7DComplex_ppHVV::evaluate
Double_t evaluate() const
Definition: RooSpinTwo_7DComplex_ppHVV.cc:286
RooSpinTwo_7DComplex_ppHVV::setZZ4fOrdering
void setZZ4fOrdering(Bool_t flag=true)
Definition: RooSpinTwo_7DComplex_ppHVV.cc:394
RooSpin::Vdecay2
RooSpin::VdecayType Vdecay2
Definition: RooSpin.h:136
RooSpinTwo_7DComplex_ppHVV::evaluateH2Factor
Double_t evaluateH2Factor(Int_t i2, Int_t j2, Int_t helicity, Int_t code) const
Definition: RooSpinTwo_7DComplex_ppHVV.cc:55
RooSpinTwo_7DComplex_ppHVV::RooSpinTwo_7DComplex_ppHVV
RooSpinTwo_7DComplex_ppHVV()
Definition: RooSpinTwo_7DComplex_ppHVV.cc:8
RooSpin::intCodeStart
Int_t intCodeStart
Definition: RooSpin.h:138
RooSpinTwo_7DComplex_ppHVV
Definition: RooSpinTwo_7DComplex_ppHVV.h:7
RooSpin::kVdecayType_GammaOnshell
@ kVdecayType_GammaOnshell
Definition: RooSpin.h:30
RooSpinTwo::modelCouplings
Definition: RooSpinTwo.h:10
RooSpin::m2
RooRealProxy m2
Definition: RooSpin.h:112
RooSpin::h2
RooRealProxy h2
Definition: RooSpin.h:109
RooSpinTwo_7DComplex_ppHVV::evaluateH1Factor
Double_t evaluateH1Factor(Int_t i1, Int_t j1, Int_t helicity, Int_t code) const
Definition: RooSpinTwo_7DComplex_ppHVV.cc:35
RooSpin::hs
RooRealProxy hs
Definition: RooSpin.h:114
MELAStreamHelpers
Definition: MELAStreamHelpers.hh:7
RooSpinTwo_7DComplex_ppHVV.h
ih2
int ih2
Definition: TMCFM.hh:51
RooSpin::prime_h1
@ prime_h1
Definition: RooSpin.h:39
RooSpinTwo_7DComplex_ppHVV::ZZ4fOrdering
Bool_t ZZ4fOrdering
Definition: RooSpinTwo_7DComplex_ppHVV.h:38
RooSpin::GeVunit
static constexpr Double_t GeVunit
Definition: RooSpin.h:141
RooSpin::getMVGamV
virtual void getMVGamV(Double_t *mV=0, Double_t *gamV=0) const
Definition: RooSpin.cc:305
RooSpinTwo_7DComplex_ppHVV::getAnalyticalIntegral
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Definition: RooSpinTwo_7DComplex_ppHVV.cc:331
RooSpin::prime_Phi1
@ prime_Phi1
Definition: RooSpin.h:43
RooSpin::Vdecay1
RooSpin::VdecayType Vdecay1
Definition: RooSpin.h:135
RooSpin::modelMeasurables
Definition: RooSpin.h:50
RooSpin::h1
RooRealProxy h1
Definition: RooSpin.h:108
RooSpinTwo_7DComplex_ppHVV::evaluatePolarizationTerms
void evaluatePolarizationTerms(std::vector< Double_t > &Axxyyterm, const Int_t code, bool isGammaV1=false, bool isGammaV2=false) const
Definition: RooSpinTwo_7DComplex_ppHVV.cc:178
RooSpin::modelParameters
Definition: RooSpin.h:61
RooSpin::calculateVffR1R2
virtual void calculateVffR1R2(Double_t &R1Val, Double_t &R2Val, bool isGammaV1=false, bool isGammaV2=false) const
Definition: RooSpin.cc:297
RooSpin::VdecayType
VdecayType
Definition: RooSpin.h:28
RooSpinTwo
Definition: RooSpinTwo.h:7
RooSpin::Phi1
RooRealProxy Phi1
Definition: RooSpin.h:115