JHUGen MELA  v2.4.1
Matrix element calculations as used in JHUGen. MELA is an important tool that was used for the Higgs boson discovery and for precise measurements of its structure and interactions. Please see the website https://spin.pha.jhu.edu/ and papers cited there for more details, and kindly cite those papers when using this code.
RooSpinZero_7DComplex_withAccep_HVV.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  accepParameters const& _accepParams,
15  RooSpin::VdecayType _Vdecay1, RooSpin::VdecayType _Vdecay2,
16  TVar::VerbosityLevel verbosity_
17 ) : RooSpinZero(
18  name, title,
19  _measurables,
20  _parameters,
21  _couplings,
22  _Vdecay1, _Vdecay2,
23  verbosity_
24 ),
25 aPhi("aPhi", "aPhi", this, (RooAbsReal&)*(_accepParams.aPhi)),
26 bPhi("bPhi", "bPhi", this, (RooAbsReal&)*(_accepParams.bPhi)),
27 cPhi("cPhi", "cPhi", this, (RooAbsReal&)*(_accepParams.cPhi)),
28 dPhi("dPhi", "dPhi", this, (RooAbsReal&)*(_accepParams.dPhi)),
29 ePhi("ePhi", "ePhi", this, (RooAbsReal&)*(_accepParams.ePhi)),
30 aPhi1("aPhi1", "aPhi1", this, (RooAbsReal&)*(_accepParams.aPhi1)),
31 bPhi1("bPhi1", "bPhi1", this, (RooAbsReal&)*(_accepParams.bPhi1)),
32 cPhi1("cPhi1", "cPhi1", this, (RooAbsReal&)*(_accepParams.cPhi1)),
33 dPhi1("dPhi1", "dPhi1", this, (RooAbsReal&)*(_accepParams.dPhi1)),
34 ePhi1("ePhi1", "ePhi1", this, (RooAbsReal&)*(_accepParams.ePhi1)),
35 aH1("aH1", "aH1", this, (RooAbsReal&)*(_accepParams.aH1)),
36 bH1("bH1", "bH1", this, (RooAbsReal&)*(_accepParams.bH1)),
37 cH1("cH1", "cH1", this, (RooAbsReal&)*(_accepParams.cH1)),
38 dH1("dH1", "dH1", this, (RooAbsReal&)*(_accepParams.dH1)),
39 eH1("eH1", "eH1", this, (RooAbsReal&)*(_accepParams.eH1)),
40 aH2("aH2", "aH2", this, (RooAbsReal&)*(_accepParams.aH2)),
41 bH2("bH2", "bH2", this, (RooAbsReal&)*(_accepParams.bH2)),
42 cH2("cH2", "cH2", this, (RooAbsReal&)*(_accepParams.cH2)),
43 dH2("dH2", "dH2", this, (RooAbsReal&)*(_accepParams.dH2)),
44 eH2("eH2", "eH2", this, (RooAbsReal&)*(_accepParams.eH2)),
45 aHs("aHs", "aHs", this, (RooAbsReal&)*(_accepParams.aHs)),
46 bHs("bHs", "bHs", this, (RooAbsReal&)*(_accepParams.bHs)),
47 cHs("cHs", "cHs", this, (RooAbsReal&)*(_accepParams.cHs)),
48 dHs("dHs", "dHs", this, (RooAbsReal&)*(_accepParams.dHs)),
49 eHs("eHs", "eHs", this, (RooAbsReal&)*(_accepParams.eHs)),
50 aM1("aM1", "aM1", this, (RooAbsReal&)*(_accepParams.aM1)),
51 bM1("bM1", "bM1", this, (RooAbsReal&)*(_accepParams.bM1)),
52 cM1("cM1", "cM1", this, (RooAbsReal&)*(_accepParams.cM1)),
53 dM1("dM1", "dM1", this, (RooAbsReal&)*(_accepParams.dM1)),
54 aM2("aM2", "aM2", this, (RooAbsReal&)*(_accepParams.aM2)),
55 bM2("bM2", "bM2", this, (RooAbsReal&)*(_accepParams.bM2)),
56 cM2("cM2", "cM2", this, (RooAbsReal&)*(_accepParams.cM2)),
57 dM2("dM2", "dM2", this, (RooAbsReal&)*(_accepParams.dM2)),
58 ZZ4fOrdering(true)
59 {}
60 
61 
63  const RooSpinZero_7DComplex_withAccep_HVV& other, const char* name
64 ) : RooSpinZero(other, name),
65 aPhi("aPhi", this, other.aPhi),
66 bPhi("bPhi", this, other.bPhi),
67 cPhi("cPhi", this, other.cPhi),
68 dPhi("dPhi", this, other.dPhi),
69 ePhi("ePhi", this, other.ePhi),
70 aPhi1("aPhi1", this, other.aPhi1),
71 bPhi1("bPhi1", this, other.bPhi1),
72 cPhi1("cPhi1", this, other.cPhi1),
73 dPhi1("dPhi1", this, other.dPhi1),
74 ePhi1("ePhi1", this, other.ePhi1),
75 aH1("aH1", this, other.aH1),
76 bH1("bH1", this, other.bH1),
77 cH1("cH1", this, other.cH1),
78 dH1("dH1", this, other.dH1),
79 eH1("eH1", this, other.eH1),
80 aH2("aH2", this, other.aH2),
81 bH2("bH2", this, other.bH2),
82 cH2("cH2", this, other.cH2),
83 dH2("dH2", this, other.dH2),
84 eH2("eH2", this, other.eH2),
85 aHs("aHs", this, other.aHs),
86 bHs("bHs", this, other.bHs),
87 cHs("cHs", this, other.cHs),
88 dHs("dHs", this, other.dHs),
89 eHs("eHs", this, other.eHs),
90 aM1("aM1", this, other.aM1),
91 bM1("bM1", this, other.bM1),
92 cM1("cM1", this, other.cM1),
93 dM1("dM1", this, other.dM1),
94 aM2("aM2", this, other.aM2),
95 bM2("bM2", this, other.bM2),
96 cM2("cM2", this, other.cM2),
97 dM2("dM2", this, other.dM2),
98 ZZ4fOrdering(other.ZZ4fOrdering)
99 {}
100 
102  Double_t& A00term, Double_t& Appterm, Double_t& Ammterm,
103  Double_t& A00ppterm, Double_t& A00mmterm, Double_t& Appmmterm,
104  const Int_t code,
105  int VGammaVpmode1, int VGammaVpmode2
106 ) const{
107  const Double_t Pi = TMath::Pi();
108 
109  Double_t R1Val, R2Val;
110  calculateVffR1R2(R1Val, R2Val, VGammaVpmode1==1, VGammaVpmode2==1);
111  if (VGammaVpmode1==2 || VGammaVpmode2==2){
112  Double_t RVp1Val=0, RVp2Val=0;
113  calculateVprimeffR1R2(RVp1Val, RVp2Val);
114  if (VGammaVpmode1==2) R1Val=RVp1Val;
115  if (VGammaVpmode2==2) R1Val=RVp2Val;
116  }
117 
118  Double_t A00Re, A00Im, AppRe, AppIm, AmmRe, AmmIm;
119  calculateAmplitudes(A00Re, A00Im, AppRe, AppIm, AmmRe, AmmIm, VGammaVpmode1, VGammaVpmode2);
120 
121  Double_t A00 = A00Im*A00Im + A00Re*A00Re;
122  Double_t App = AppIm*AppIm + AppRe*AppRe;
123  Double_t Amm = AmmIm*AmmIm + AmmRe*AmmRe;
124 
125  Double_t phi00=atan2(A00Im, A00Re);
126  Double_t phipp=atan2(AppIm, AppRe)-phi00;
127  Double_t phimm=atan2(AmmIm, AmmRe)-phi00;
128 
129  Double_t A00_prefactor = 1.;
130  Double_t Amm_pp_prefactor = 1.;
131  Double_t A00mm_prefactor = 2.;
132  Double_t A00pp_prefactor = 2.;
133  Double_t Ammpp_prefactor = 2.;
134 
135  A00term = A00*A00_prefactor;
136  Appterm = App*Amm_pp_prefactor;
137  Ammterm = Amm*Amm_pp_prefactor;
138  A00ppterm = sqrt(A00*App)*A00pp_prefactor;
139  A00mmterm = sqrt(A00*Amm)*A00mm_prefactor;
140  Appmmterm = sqrt(Amm*App)*Ammpp_prefactor;
141 
142  if ((code % prime_h1)==0){
143  Double_t A00_h1int = 4.*aH1/3. + 4.*bH1/15. + 4.*cH1/35. + 4.*dH1/63. + 4.*eH1/99.;
144  Double_t Ammpp_h1int = 8.*aH1/3. + 16.*bH1/15. + 24.*cH1/35. + 32.*dH1/63. + 40.*eH1/99.;
145  Double_t A0m_h1int = (128.*aH1 + 32.*bH1 + 16.*cH1 + 10.*dH1 + 7.*eH1)*Pi/256.*R1Val;
146  Double_t A0p_h1int = (128.*aH1 + 32.*bH1 + 16.*cH1 + 10.*dH1 + 7.*eH1)*Pi/256.*R1Val;
147  Double_t Amp_h1int = 4.*aH1/3. + 4.*bH1/15. + 4.*cH1/35. + 4.*dH1/63. + 4.*eH1/99.;
148 
149  if (A00term!=0) A00term *= A00_h1int;
150  if (Appterm!=0) Appterm *= Ammpp_h1int;
151  if (Ammterm!=0) Ammterm *= Ammpp_h1int;
152  if (A00ppterm!=0) A00ppterm *= A0p_h1int;
153  if (A00mmterm!=0) A00mmterm *= A0m_h1int;
154  if (Appmmterm!=0) Appmmterm *= Amp_h1int;
155  }
156  else{
157  Double_t common_fac = (aH1 + bH1*pow(h1, 2) + cH1*pow(h1, 4) +dH1*pow(h1, 6) + eH1*pow(h1, 8));
158 
159  if (A00term!=0) A00term *= (1. - pow(h1, 2))*common_fac;
160  if (Appterm!=0) Appterm *= (1. + pow(h1, 2) - 2.*h1*R1Val)*common_fac;
161  if (Ammterm!=0) Ammterm *= (1. + pow(h1, 2) + 2.*h1*R1Val)*common_fac;
162  if (A00ppterm!=0) A00ppterm *= sqrt(1. - pow(h1, 2))*(R1Val - h1)*common_fac;
163  if (A00mmterm!=0) A00mmterm *= sqrt(1. - pow(h1, 2))*(R1Val + h1)*common_fac;
164  if (Appmmterm!=0) Appmmterm *= (1. - pow(h1, 2))*common_fac;
165  }
166 
167  if ((code % prime_h2)==0){
168  Double_t A00_h2int = 4.*aH2/3. + 4.*bH2/15. + 4.*cH2/35. + 4.*dH2/63. + 4.*eH2/99.;
169  Double_t Ammpp_h2int = 8.*aH2/3. + 16.*bH2/15. + 24.*cH2/35. + 32.*dH2/63. + 40.*eH2/99.;
170  Double_t A0m_h2int = (128.*aH2 + 32.*bH2 + 16.*cH2 + 10.*dH2 + 7.*eH2)*Pi/256.*R2Val;
171  Double_t A0p_h2int = (128.*aH2 + 32.*bH2 + 16.*cH2 + 10.*dH2 + 7.*eH2)*Pi/256.*R2Val;
172  Double_t Amp_h2int = 4.*aH2/3. + 4.*bH2/15. + 4.*cH2/35. + 4.*dH2/63. + 4.*eH2/99.;
173 
174  if (A00term!=0) A00term *= A00_h2int;
175  if (Appterm!=0) Appterm *= Ammpp_h2int;
176  if (Ammterm!=0) Ammterm *= Ammpp_h2int;
177  if (A00ppterm!=0) A00ppterm *= A0p_h2int;
178  if (A00mmterm!=0) A00mmterm *= A0m_h2int;
179  if (Appmmterm!=0) Appmmterm *= Amp_h2int;
180  }
181  else{
182  Double_t common_fac = (aH2 + bH2*pow(h2, 2) + cH2*pow(h2, 4) +dH2*pow(h2, 6) + eH2*pow(h2, 8));
183 
184  if (A00term!=0) A00term *= (1. - pow(h2, 2))*common_fac;
185  if (Appterm!=0) Appterm *= (1. + pow(h2, 2) - 2.*h2*R2Val)*common_fac;
186  if (Ammterm!=0) Ammterm *= (1. + pow(h2, 2) + 2.*h2*R2Val)*common_fac;
187  if (A00ppterm!=0) A00ppterm *= sqrt(1. - pow(h2, 2))*(R2Val - h2)*common_fac;
188  if (A00mmterm!=0) A00mmterm *= sqrt(1. - pow(h2, 2))*(R2Val + h2)*common_fac;
189  if (Appmmterm!=0) Appmmterm *= (1. - pow(h2, 2))*common_fac;
190  }
191 
192  if ((code % prime_hs)==0){
193  Double_t A_hsint = 2.*aHs + 2.*bHs/3. + 2.*cHs/5. + 2.*dHs/7. + 2.*eHs/9.; // All of them
194 
195  if (A00term!=0) A00term *= A_hsint;
196  if (Appterm!=0) Appterm *= A_hsint;
197  if (Ammterm!=0) Ammterm *= A_hsint;
198  if (A00ppterm!=0) A00ppterm *= A_hsint;
199  if (A00mmterm!=0) A00mmterm *= A_hsint;
200  if (Appmmterm!=0) Appmmterm *= A_hsint;
201  }
202  else{
203  Double_t common_fac = (aHs + bHs*pow(hs, 2) + cHs*pow(hs, 4) + dHs*pow(hs, 6) +eHs*pow(hs, 8));
204 
205  if (A00term!=0) A00term *= common_fac;
206  if (Appterm!=0) Appterm *= common_fac;
207  if (Ammterm!=0) Ammterm *= common_fac;
208  if (A00ppterm!=0) A00ppterm *= common_fac;
209  if (A00mmterm!=0) A00mmterm *= common_fac;
210  if (Appmmterm!=0) Appmmterm *= common_fac;
211  }
212 
213  if ((code % prime_Phi)==0){
214  Double_t A00mmpp_phiint = 2.*aPhi*Pi;
215  Double_t A0p_phiint = bPhi*Pi*cos(phipp);
216  Double_t A0m_phiint = bPhi*Pi*cos(phimm);
217  Double_t Amp_phiint = cPhi*Pi*cos(phimm - phipp);
218 
219  if (A00term!=0) A00term *= A00mmpp_phiint;
220  if (Appterm!=0) Appterm *= A00mmpp_phiint;
221  if (Ammterm!=0) Ammterm *= A00mmpp_phiint;
222  if (A00ppterm!=0) A00ppterm *= A0p_phiint;
223  if (A00mmterm!=0) A00mmterm *= A0m_phiint;
224  if (Appmmterm!=0) Appmmterm *= Amp_phiint;
225  }
226  else{
227  Double_t common_fac = (aPhi + bPhi*cos(Phi) + cPhi*cos(2*Phi) + dPhi*cos(3*Phi) +ePhi*cos(4*Phi));
228 
229  if (A00term!=0) A00term *= common_fac;
230  if (Appterm!=0) Appterm *= common_fac;
231  if (Ammterm!=0) Ammterm *= common_fac;
232  if (A00ppterm!=0) A00ppterm *= cos(Phi + phipp)*common_fac;
233  if (A00mmterm!=0) A00mmterm *= cos(Phi - phimm)*common_fac;
234  if (Appmmterm!=0) Appmmterm *= cos(2*Phi - phimm + phipp)*common_fac;
235  }
236 
237  if ((code % prime_Phi1)==0){
238  Double_t A_phi1int = 2.*aPhi1*Pi;
239 
240  if (A00term!=0) A00term *= A_phi1int;
241  if (Appterm!=0) Appterm *= A_phi1int;
242  if (Ammterm!=0) Ammterm *= A_phi1int;
243  if (A00ppterm!=0) A00ppterm *= A_phi1int;
244  if (A00mmterm!=0) A00mmterm *= A_phi1int;
245  if (Appmmterm!=0) Appmmterm *= A_phi1int;
246  }
247  else{
248  Double_t common_fac = (aPhi1 + bPhi1*cos(Phi1) + cPhi1*cos(2*Phi1) + dPhi1*cos(3*Phi1) +ePhi1*cos(4*Phi1));
249 
250  if (A00term!=0) A00term *= common_fac;
251  if (Appterm!=0) Appterm *= common_fac;
252  if (Ammterm!=0) Ammterm *= common_fac;
253  if (A00ppterm!=0) A00ppterm *= common_fac;
254  if (A00mmterm!=0) A00mmterm *= common_fac;
255  if (Appmmterm!=0) Appmmterm *= common_fac;
256  }
257 
258  if (verbosity>=TVar::DEBUG){
259  MELAout << "RooSpinZero_7DComplex_withAccep_HVV::evaluatePolarizationTerms( " << code << " , " << VGammaVpmode1 << " , " << VGammaVpmode2 << " ):" << endl;
260  MELAout << "\t- |A00|**2 = " << A00term << endl;
261  MELAout << "\t- |A++|**2 = " << Appterm << endl;
262  MELAout << "\t- |A--|**2 = " << Ammterm << endl;
263  MELAout << "\t- |A00||A++| = " << A00ppterm << endl;
264  MELAout << "\t- |A00||A--| = " << A00mmterm << endl;
265  MELAout << "\t- |A++||A--| = " << Appmmterm << endl;
266  }
267 }
268 
270  Double_t mV;
271  getMVGamV(&mV);
272  bool isZZ = (mV >= 90.);
273  Double_t epsilon=1e-15;
274  Double_t m1_=m1; if (Vdecay1==RooSpin::kVdecayType_GammaOnshell) m1_=0;
275  Double_t m2_=m2; if (Vdecay2==RooSpin::kVdecayType_GammaOnshell) m2_=0;
276  if (
277  (m1_+m2_)>m12 ||
278  (isZZ && Vdecay1==Vdecay2 && ZZ4fOrdering && fabs(m2_-mV)<fabs(m1_-mV) && Vdecay2!=RooSpin::kVdecayType_GammaOnshell) ||
281  ) return epsilon;
282 
283  Int_t code = intCodeStart;
285  code *= prime_Phi;
289  }
290 
291  Double_t betaValSq = (1.-(pow(m1_-m2_, 2)/pow(m12, 2)))*(1.-(pow(m1_+m2_, 2)/pow(m12, 2)));
292  if (betaValSq<=0.) return epsilon;
293  Double_t betaVal = sqrt(betaValSq);
294 
295  Double_t term1Coeff = 1;
296  Double_t term2Coeff = 1;
297  if (Vdecay1!=RooSpin::kVdecayType_GammaOnshell) term1Coeff = 2.*m1_*GeVunit; // dm**2 = 2m dm
298  if (Vdecay2!=RooSpin::kVdecayType_GammaOnshell) term2Coeff = 2.*m2_*GeVunit;
299 
300  Double_t value = 0;
301  for (int VGammaVpmode1=0; VGammaVpmode1<=2; VGammaVpmode1++){
302  for (int VGammaVpmode2=0; VGammaVpmode2<=2; VGammaVpmode2++){
303  if (!(
304  (VGammaVpmode1==1 || Vdecay1!=RooSpin::kVdecayType_GammaOnshell)
305  &&
306  (VGammaVpmode2==1 || Vdecay2!=RooSpin::kVdecayType_GammaOnshell)
307  )
308  ||
309  (VGammaVpmode1==1 && VGammaVpmode2==2) || (VGammaVpmode1==2 && VGammaVpmode2==1)
310  ||
311  !computeNeededAmplitude(VGammaVpmode1, VGammaVpmode2)
312  ) continue;
313  Double_t val_A00=0, val_App=0, val_Amm=0, val_A0p=0, val_A0m=0, val_Amp=0;
314  evaluatePolarizationTerms(val_A00, val_App, val_Amm, val_A0p, val_A0m, val_Amp, code, VGammaVpmode1, VGammaVpmode2);
315  value += val_A00 + val_App + val_Amm + val_A0p + val_A0m + val_Amp;
316  }
317  }
318  value *= betaVal*term1Coeff*term2Coeff
319  *(1+aM1*m1_+bM1*m1_*m1_+cM1*m1_*m1_*m1_+dM1*m1_*m1_*m1_*m1_)
320  *(1+aM2*m2_+bM2*m2_*m2_+cM2*m2_*m2_*m2_+dM2*m2_*m2_*m2_*m2_);
321 
322  if (!(value==value)){
323  MELAout << "Evaluate NaN=" << value << " at "
324  << "h1=" << h1 << '\t'
325  << "h2=" << h2 << '\t'
326  << "hs=" << hs << '\t'
327  << "Phi1=" << Phi1 << '\t'
328  << "Phi=" << Phi << '\t'
329  << "m1=" << m1 << '\t'
330  << "m2=" << m2 << '\t'
331  << "m12=" << m12 << '\t'
332  << endl;
333  MELAout << "Possible sources:\n"
334  << "betaVal=" << betaVal << '\t'
335  << "term1Coeff=" << term1Coeff << '\t'
336  << "term2Coeff=" << term2Coeff
337  << endl;
338  }
339  return value;
340 }
341 
342 Int_t RooSpinZero_7DComplex_withAccep_HVV::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const{
343  Int_t code = intCodeStart;
344  if (checkFundamentalType(h1)){ if (matchArgs(allVars, analVars, h1) || Vdecay1==RooSpin::kVdecayType_GammaOnshell) code *= prime_h1; }
345  if (checkFundamentalType(h2)){ if (matchArgs(allVars, analVars, h2) || Vdecay2==RooSpin::kVdecayType_GammaOnshell) code *= prime_h2; }
346  if (checkFundamentalType(hs)){ if (matchArgs(allVars, analVars, hs)) code *= prime_hs; }
349  if (code==1) code=0;
350  return code;
351 }
352 Double_t RooSpinZero_7DComplex_withAccep_HVV::analyticalIntegral(Int_t code, const char* /*rangeName*/) const{
353  Double_t mV;
354  getMVGamV(&mV);
355  bool isZZ = (mV >= 90.);
356  Double_t epsilon=1e-10;
357  Double_t m1_=m1; if (Vdecay1==RooSpin::kVdecayType_GammaOnshell) m1_=0;
358  Double_t m2_=m2; if (Vdecay2==RooSpin::kVdecayType_GammaOnshell) m2_=0;
359  if (
360  (m1_+m2_)>m12 ||
361  (isZZ && Vdecay1==Vdecay2 && ZZ4fOrdering && fabs(m2_-mV)<fabs(m1_-mV) && Vdecay2!=RooSpin::kVdecayType_GammaOnshell) ||
364  ) return epsilon;
365 
366  Double_t betaValSq = (1.-(pow(m1_-m2_, 2)/pow(m12, 2)))*(1.-(pow(m1_+m2_, 2)/pow(m12, 2)));
367  if (betaValSq<=0.) return epsilon;
368  Double_t betaVal = sqrt(betaValSq);
369 
370  Double_t term1Coeff = 1;
371  Double_t term2Coeff = 1;
372  if (Vdecay1!=RooSpin::kVdecayType_GammaOnshell) term1Coeff = 2.*m1_*GeVunit; // dm**2 = 2m dm
373  if (Vdecay2!=RooSpin::kVdecayType_GammaOnshell) term2Coeff = 2.*m2_*GeVunit;
374 
375  Double_t value = 0;
376  for (int VGammaVpmode1=0; VGammaVpmode1<=2; VGammaVpmode1++){
377  for (int VGammaVpmode2=0; VGammaVpmode2<=2; VGammaVpmode2++){
378  if (!(
379  (VGammaVpmode1==0 && VGammaVpmode2==0 && Vdecay1!=RooSpin::kVdecayType_GammaOnshell && Vdecay2!=RooSpin::kVdecayType_GammaOnshell)
380  ||
381  (VGammaVpmode1==0 && VGammaVpmode2==1 && Vdecay1!=RooSpin::kVdecayType_GammaOnshell)
382  ||
383  (VGammaVpmode1==1 && VGammaVpmode2==0 && Vdecay2!=RooSpin::kVdecayType_GammaOnshell)
384  ||
385  (VGammaVpmode1==1 && VGammaVpmode2==1)
386  ||
387  (VGammaVpmode1==0 && VGammaVpmode2==2 && Vdecay1!=RooSpin::kVdecayType_GammaOnshell && Vdecay2!=RooSpin::kVdecayType_GammaOnshell)
388  ||
389  (VGammaVpmode1==2 && VGammaVpmode2==0 && Vdecay1!=RooSpin::kVdecayType_GammaOnshell && Vdecay2!=RooSpin::kVdecayType_GammaOnshell)
390  ||
391  (VGammaVpmode1==2 && VGammaVpmode2==2 && Vdecay1!=RooSpin::kVdecayType_GammaOnshell && Vdecay2!=RooSpin::kVdecayType_GammaOnshell)
392  )) continue;
393  Double_t val_A00=0, val_App=0, val_Amm=0, val_A0p=0, val_A0m=0, val_Amp=0;
394  evaluatePolarizationTerms(val_A00, val_App, val_Amm, val_A0p, val_A0m, val_Amp, code, VGammaVpmode1, VGammaVpmode2);
395  value += val_A00 + val_App + val_Amm + val_A0p + val_A0m + val_Amp;
396  }
397  }
398  value *= betaVal*term1Coeff*term2Coeff
399  *(1+aM1*m1_+bM1*m1_*m1_+cM1*m1_*m1_*m1_+dM1*m1_*m1_*m1_*m1_)
400  *(1+aM2*m2_+bM2*m2_*m2_+cM2*m2_*m2_*m2_+dM2*m2_*m2_*m2_*m2_);
401 
402  if (!(value==value)){
403  MELAout << "Integral NaN=" << value << " at "
404  << "h1=" << h1 << '\t'
405  << "h2=" << h2 << '\t'
406  << "hs=" << hs << '\t'
407  << "Phi1=" << Phi1 << '\t'
408  << "Phi=" << Phi << '\t'
409  << "m1=" << m1 << '\t'
410  << "m2=" << m2 << '\t'
411  << "m12=" << m12 << '\t'
412  << endl;
413  MELAout << "Possible sources:\n"
414  << "betaVal=" << betaVal << '\t'
415  << "term1Coeff=" << term1Coeff << '\t'
416  << "term2Coeff=" << term2Coeff
417  << endl;
418  }
419  return value;
420 }
421 
RooSpinZero_7DComplex_withAccep_HVV::dPhi1
RooRealProxy dPhi1
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:78
RooSpinZero_7DComplex_withAccep_HVV::cHs
RooRealProxy cHs
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:92
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
RooSpinZero_7DComplex_withAccep_HVV::dPhi
RooRealProxy dPhi
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:73
RooSpinZero_7DComplex_withAccep_HVV::aM2
RooRealProxy aM2
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:101
RooSpinZero_7DComplex_withAccep_HVV::ZZ4fOrdering
Bool_t ZZ4fOrdering
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:106
RooSpinZero_7DComplex_withAccep_HVV::aPhi1
RooRealProxy aPhi1
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:75
RooSpinZero_7DComplex_withAccep_HVV::evaluatePolarizationTerms
void evaluatePolarizationTerms(Double_t &A00term, Double_t &Appterm, Double_t &Ammterm, Double_t &A00ppterm, Double_t &A00mmterm, Double_t &Appmmterm, const Int_t code, int VGammaVpmode1=0, int VGammaVpmode2=0) const
Definition: RooSpinZero_7DComplex_withAccep_HVV.cc:101
RooSpinZero_7DComplex_withAccep_HVV::ePhi
RooRealProxy ePhi
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:74
RooSpinZero_7DComplex_withAccep_HVV::aH1
RooRealProxy aH1
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:80
RooSpin::m12
RooRealProxy m12
Definition: RooSpin.h:113
RooSpin::prime_h1
@ prime_h1
Definition: RooSpin.h:39
RooSpinZero_7DComplex_withAccep_HVV::getAnalyticalIntegral
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Definition: RooSpinZero_7DComplex_withAccep_HVV.cc:342
RooSpin::prime_Phi
@ prime_Phi
Definition: RooSpin.h:42
RooSpinZero_7DComplex_withAccep_HVV::dM1
RooRealProxy dM1
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:99
RooSpinZero_7DComplex_withAccep_HVV::aH2
RooRealProxy aH2
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:85
RooSpinZero_7DComplex_withAccep_HVV::evaluate
Double_t evaluate() const
Definition: RooSpinZero_7DComplex_withAccep_HVV.cc:269
RooSpin::prime_Phi1
@ prime_Phi1
Definition: RooSpin.h:43
RooSpinZero_7DComplex_withAccep_HVV::aM1
RooRealProxy aM1
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:96
RooSpinZero_7DComplex_withAccep_HVV::bM1
RooRealProxy bM1
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:97
RooSpin::m1
RooRealProxy m1
Definition: RooSpin.h:111
RooSpinZero::computeNeededAmplitude
virtual Bool_t computeNeededAmplitude(int VGammaVpmode1, int VGammaVpmode2) const final
Definition: RooSpinZero.cc:583
RooSpinZero_7DComplex_withAccep_HVV::cPhi
RooRealProxy cPhi
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:72
RooSpin::calculateVprimeffR1R2
virtual void calculateVprimeffR1R2(Double_t &R1Val, Double_t &R2Val) const
Definition: RooSpin.cc:343
MELAStreamHelpers::MELAout
MELAOutputStreamer MELAout
RooSpin::checkFundamentalType
virtual Bool_t checkFundamentalType(const RooRealProxy &proxy) const
Definition: RooSpin.cc:392
RooSpinZero_7DComplex_withAccep_HVV::bPhi
RooRealProxy bPhi
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:71
RooSpinZero_7DComplex_withAccep_HVV::ePhi1
RooRealProxy ePhi1
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:79
RooSpinZero_7DComplex_withAccep_HVV::cH2
RooRealProxy cH2
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:87
RooSpinZero_7DComplex_withAccep_HVV::setZZ4fOrdering
void setZZ4fOrdering(Bool_t flag=true)
Definition: RooSpinZero_7DComplex_withAccep_HVV.cc:422
RooSpinZero_7DComplex_withAccep_HVV::eH1
RooRealProxy eH1
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:84
RooSpinZero_7DComplex_withAccep_HVV::analyticalIntegral
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Definition: RooSpinZero_7DComplex_withAccep_HVV.cc:352
RooSpin::verbosity
TVar::VerbosityLevel verbosity
Definition: RooSpin.h:139
RooSpinZero_7DComplex_withAccep_HVV::bH1
RooRealProxy bH1
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:81
RooSpinZero_7DComplex_withAccep_HVV::RooSpinZero_7DComplex_withAccep_HVV
RooSpinZero_7DComplex_withAccep_HVV()
Definition: RooSpinZero_7DComplex_withAccep_HVV.cc:8
RooSpin::Vdecay2
RooSpin::VdecayType Vdecay2
Definition: RooSpin.h:136
RooSpinZero
Definition: RooSpinZero.h:7
RooSpinZero_7DComplex_withAccep_HVV::dM2
RooRealProxy dM2
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:104
RooSpinZero_7DComplex_withAccep_HVV::accepParameters
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:11
RooSpin::intCodeStart
Int_t intCodeStart
Definition: RooSpin.h:138
RooSpin::kVdecayType_GammaOnshell
@ kVdecayType_GammaOnshell
Definition: RooSpin.h:30
RooSpin::prime_h2
@ prime_h2
Definition: RooSpin.h:40
RooSpinZero_7DComplex_withAccep_HVV::cM2
RooRealProxy cM2
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:103
RooSpinZero_7DComplex_withAccep_HVV::bHs
RooRealProxy bHs
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:91
RooSpinZero_7DComplex_withAccep_HVV.h
RooSpinZero_7DComplex_withAccep_HVV::dHs
RooRealProxy dHs
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:93
RooSpinZero_7DComplex_withAccep_HVV::eHs
RooRealProxy eHs
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:94
RooSpinZero_7DComplex_withAccep_HVV::bPhi1
RooRealProxy bPhi1
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:76
RooSpinZero_7DComplex_withAccep_HVV::cM1
RooRealProxy cM1
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:98
RooSpinZero_7DComplex_withAccep_HVV::cPhi1
RooRealProxy cPhi1
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:77
RooSpinZero_7DComplex_withAccep_HVV::bH2
RooRealProxy bH2
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:86
RooSpinZero::modelCouplings
Definition: RooSpinZero.h:10
RooSpin::m2
RooRealProxy m2
Definition: RooSpin.h:112
RooSpin::h2
RooRealProxy h2
Definition: RooSpin.h:109
RooSpinZero::calculateAmplitudes
virtual void calculateAmplitudes(Double_t &A00Re, Double_t &A00Im, Double_t &AppRe, Double_t &AppIm, Double_t &AmmRe, Double_t &AmmIm, int VGammaVpmode1=0, int VGammaVpmode2=0) const
Definition: RooSpinZero.cc:477
RooSpinZero_7DComplex_withAccep_HVV::dH1
RooRealProxy dH1
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:83
RooSpin::hs
RooRealProxy hs
Definition: RooSpin.h:114
MELAStreamHelpers
Definition: MELAStreamHelpers.hh:7
RooSpinZero_7DComplex_withAccep_HVV
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:7
RooSpinZero_7DComplex_withAccep_HVV::bM2
RooRealProxy bM2
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:102
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
RooSpinZero_7DComplex_withAccep_HVV::eH2
RooRealProxy eH2
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:89
TVar::DEBUG
@ DEBUG
Definition: TVar.hh:51
RooSpin::Vdecay1
RooSpin::VdecayType Vdecay1
Definition: RooSpin.h:135
RooSpin::modelMeasurables
Definition: RooSpin.h:50
RooSpinZero_7DComplex_withAccep_HVV::dH2
RooRealProxy dH2
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:88
RooSpin::h1
RooRealProxy h1
Definition: RooSpin.h:108
RooSpin::modelParameters
Definition: RooSpin.h:61
RooSpinZero_7DComplex_withAccep_HVV::aHs
RooRealProxy aHs
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:90
RooSpinZero_7DComplex_withAccep_HVV::aPhi
RooRealProxy aPhi
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:70
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
RooSpinZero_7DComplex_withAccep_HVV::cH1
RooRealProxy cH1
Definition: RooSpinZero_7DComplex_withAccep_HVV.h:82
RooSpin::prime_hs
@ prime_hs
Definition: RooSpin.h:41
RooSpin::Phi1
RooRealProxy Phi1
Definition: RooSpin.h:115