a Code for the Combination of Indirect and Direct Constraints on High Energy Physics Models Logo
MVlnu.cpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2014 HEPfit Collaboration
3 *
4 *
5 * For the licensing terms see doc/COPYING.
6 */
7
8#include "StandardModel.h"
9#include "MVlnu.h"
10#include "std_make_vector.h"
11#include "gslpp_function_adapter.h"
12#include <boost/bind/bind.hpp>
13#include <limits>
14#include <gsl/gsl_sf_gegenbauer.h>
15#include <gsl/gsl_sf_expint.h>
16#include <limits>
17using namespace boost::placeholders;
18
19MVlnu::MVlnu(const StandardModel& SM_i, QCD::meson meson_i, QCD::meson vector_i, QCD::lepton lep_i)
20: mySM(SM_i), ig(ROOT::Math::Integration::kADAPTIVESINGULAR, ROOT::Math::Integration::kGAUSS51), wf()
21{
22 lep = lep_i;
23 meson = meson_i;
24 vectorM = vector_i;
25 CLNflag = false;
26 BGLflag = false;
27 DMflag = false;
28 btocNPpmflag = false;
29 NPanalysis = false;
30
31 checkcache_int_tau = false;
32 checkcache_int_mu = false;
33 checkcache_int_el = false;
34
35 double max_double = std::numeric_limits<double>::max();
36
37 hA1w1_cache = max_double;
38 rho2_cache = max_double;
39 R1w1_cache = max_double;
40 R2w1_cache = max_double;
41 N_A_cache = max_double;
42 N_1_cache =max_double;
43 N_2_cache = max_double;
44 j_A_cache= max_double;
45 j_0_cache = max_double;
46 j_1_cache = max_double;
47 j_2_cache = max_double;
48 k_A_cache = max_double;
49 k_0_cache = max_double;
50 k_1_cache = max_double;
51 k_2_cache = max_double;
52 l_A_cache = max_double;
53
54 af0_cache = max_double;
55 af1_cache = max_double;
56 af2_cache = max_double;
57 ag0_cache = max_double;
58 ag1_cache = max_double;
59 ag2_cache = max_double;
60 aF11_cache = max_double;
61 aF12_cache = max_double;
62 aF13_cache = max_double;
63 aF21_cache = max_double;
64 aF22_cache = max_double;
65 aF23_cache = max_double;
66
67 af_1_cache = max_double;
68 ag_1_cache = max_double;
69 aF1_1_cache = max_double;
70 aF2_1_cache = max_double;
71 af_2_cache = max_double;
72 ag_2_cache = max_double;
73 aF1_2_cache = max_double;
74 aF2_2_cache = max_double;
75 af_3_cache = max_double;
76 ag_3_cache = max_double;
77 aF1_3_cache = max_double;
78 aF2_3_cache = max_double;
79 af_4_cache = max_double;
80 ag_4_cache = max_double;
81 aF1_4_cache = max_double;
82 aF2_4_cache = max_double;
83 af_5_cache = max_double;
84 ag_5_cache = max_double;
85 aF1_5_cache = max_double;
86 aF2_5_cache = max_double;
87 af_6_cache = max_double;
88 ag_6_cache = max_double;
89 aF1_6_cache = max_double;
90 aF2_6_cache = max_double;
91 af_7_cache = max_double;
92 ag_7_cache = max_double;
93 aF1_7_cache = max_double;
94 aF2_7_cache = max_double;
95 af_8_cache = max_double;
96 ag_8_cache = max_double;
97 aF1_8_cache = max_double;
98 aF2_8_cache = max_double;
99 af_9_cache = max_double;
100 ag_9_cache = max_double;
101 aF1_9_cache = max_double;
102 aF2_9_cache = max_double;
103 af_10_cache = max_double;
104 ag_10_cache = max_double;
105 aF1_10_cache = max_double;
106 aF2_10_cache = max_double;
107 bf_1_cache = max_double;
108 bg_1_cache = max_double;
109 bF1_1_cache = max_double;
110 bF2_1_cache = max_double;
111 bf_2_cache = max_double;
112 bg_2_cache = max_double;
113 bF1_2_cache = max_double;
114 bF2_2_cache = max_double;
115 bf_3_cache = max_double;
116 bg_3_cache = max_double;
117 bF1_3_cache = max_double;
118 bF2_3_cache = max_double;
119 bf_4_cache = max_double;
120 bg_4_cache = max_double;
121 bF1_4_cache = max_double;
122 bF2_4_cache = max_double;
123 bf_5_cache = max_double;
124 bg_5_cache = max_double;
125 bF1_5_cache = max_double;
126 bF2_5_cache = max_double;
127 bf_6_cache = max_double;
128 bg_6_cache = max_double;
129 bF1_6_cache = max_double;
130 bF2_6_cache = max_double;
131 bf_7_cache = max_double;
132 bg_7_cache = max_double;
133 bF1_7_cache = max_double;
134 bF2_7_cache = max_double;
135 bf_8_cache = max_double;
136 bg_8_cache = max_double;
137 bF1_8_cache = max_double;
138 bF2_8_cache = max_double;
139 bf_9_cache = max_double;
140 bg_9_cache = max_double;
141 bF1_9_cache = max_double;
142 bF2_9_cache = max_double;
143 bf_10_cache = max_double;
144 bg_10_cache = max_double;
145 bF1_10_cache = max_double;
146 bF2_10_cache = max_double;
147
148 CS_cache = max_double;
149 CSp_cache = max_double;
150 CP_cache = max_double;
151 CPp_cache = max_double;
152 CV_cache = max_double;
153 CVp_cache = max_double;
154 CA_cache = max_double;
155 CAp_cache = max_double;
156 CT_cache = max_double;
157 CTp_cache = max_double;
158
159 ig.SetRelTolerance(1.E-6); // set relative tolerance
160 ig.SetAbsTolerance(1.E-6); // set absoulte tolerance
161
162}
163
165}
166
167std::vector<std::string> MVlnu::initializeMVlnuParameters()
168{
172 btocNPpmflag = (mySM.getModelName().compare("RealWeakEFTCCPM") == 0);
173 NPanalysis = (mySM.getModelName().compare("RealWeakEFTCCPM") == 0 || mySM.getModelName().compare("RealWeakEFTCC") == 0);
174
175 if (CLNflag + BGLflag + DMflag != true) throw std::runtime_error("MVlnu: Set only one among CLNflag, BGLflag, DMflag to true");
176 else mvlnuParameters = make_vector<std::string>();
177 if (CLNflag) {
178 mvlnuParameters.clear();
179 if (vectorM == StandardModel::D_star_P) mvlnuParameters = make_vector<std::string>()
180 << "hA1w1" << "rho2" << "R1w1" << "R2w1"
181 << "N_A" << "N_1" << "N_2" << "j_A" << "j_0" << "j_1" << "j_2"
182 << "k_A" << "k_0" << "k_1" << "k_2" << "l_A";
183 }
184 else if (BGLflag) {
185 mvlnuParameters.clear();
186 if (vectorM == StandardModel::D_star_P) mvlnuParameters = make_vector<std::string>()
187 << "af0" << "af1" << "af2" << "ag0" << "ag1" << "ag2"
188 << "aF11" << "aF12" << "aF13" << "aF21" << "aF22" << "aF23"
189 << "mBcstV1" << "mBcstV2" << "mBcstV3" << "mBcstV4"
190 << "mBcstA1" << "mBcstA2" << "mBcstA3" << "mBcstA4"
191 << "mBcstP1" << "mBcstP2" << "mBcstP3"
192 << "chiTV" << "chiTA" << "chiTP" << "nI";
193 }
194 else if (DMflag){
195 mvlnuParameters.clear();
196 if (vectorM == StandardModel::D_star_P) mvlnuParameters = make_vector<std::string>()
197 << "af_1" << "ag_1" << "aF1_1" << "aF2_1"
198 << "af_2" << "ag_2" << "aF1_2" << "aF2_2"
199 << "af_3" << "ag_3" << "aF1_3" << "aF2_3"
200 << "af_4" << "ag_4" << "aF1_4" << "aF2_4"
201 << "af_5" << "ag_5" << "aF1_5" << "aF2_5"
202 << "af_6" << "ag_6" << "aF1_6" << "aF2_6"
203 << "af_7" << "ag_7" << "aF1_7" << "aF2_7"
204 << "af_8" << "ag_8" << "aF1_8" << "aF2_8"
205 << "af_9" << "ag_9" << "aF1_9" << "aF2_9"
206 << "af_10" << "ag_10" << "aF1_10" << "aF2_10"
207 << "bf_1" << "bg_1" << "bF1_1" << "bF2_1"
208 << "bf_2" << "bg_2" << "bF1_2" << "bF2_2"
209 << "bf_3" << "bg_3" << "bF1_3" << "bF2_3"
210 << "bf_4" << "bg_4" << "bF1_4" << "bF2_4"
211 << "bf_5" << "bg_5" << "bF1_5" << "bF2_5"
212 << "bf_6" << "bg_6" << "bF1_6" << "bF2_6"
213 << "bf_7" << "bg_7" << "bF1_7" << "bF2_7"
214 << "bf_8" << "bg_8" << "bF1_8" << "bF2_8"
215 << "bf_9" << "bg_9" << "bF1_9" << "bF2_9"
216 << "bf_10" << "bg_10" << "bF1_10" << "bF2_10";
217 }
218 else {
219 std::stringstream out;
220 out << vectorM;
221 throw std::runtime_error("MVlnu: vector " + out.str() + " not implemented");
222 }
223
226 return mvlnuParameters;
227}
228
229void MVlnu::updateParameters()
230{
231 if (!mySM.getFlavour().getUpdateFlag(meson, vectorM, lep)) return;
232
234 Mnu = 0.; // neutrinos assumed to be massless
238 w0 = (MM*MM+MV*MV)/(2.*MM*MV);
239 RV = 2.*sqrt(MM*MV)/(MM+MV);
240 mu_b = MM; // mySM.getMub();
241 Mb = mySM.getQuarks(QCD::BOTTOM).getMass(); // add the PS b mass
242 Mc = mySM.getQuarks(QCD::CHARM).getMass(); // add the PS b mass
243 ale_mub = mySM.Ale(mu_b,FULLNLO);
244 /* Amplitude propto 4*GF*Vij/sqrt(2) & kinematics requires 1/(2^9 pi^3 MB^3) */
245 amplsq_factor = 1./(64.*M_PI*M_PI*M_PI*MM*MM*MM);
246 q2min = Mlep*Mlep;
247 q2max = (MM-MV)*(MM-MV);
248
249 MV_o_MM = MV / MM;
250 sqrtMV_o_MM = sqrt(MV_o_MM);
251
252 /* SM + NP Wilson coefficients */
253 gslpp::complex norm = 4./sqrt(2.);
254 gslpp::vector<gslpp::complex> ** allcoeff_bclnu = mySM.getFlavour().ComputeCoeffdiujlknu(2,1,0,mu_b);
255 CV = (*(allcoeff_bclnu[LO]))(0)/norm*(1.+ale_mub/M_PI*log(mySM.getMz()/mu_b))/2.;
256 CA = -CV;
257 CVp = (*(allcoeff_bclnu[LO]))(1)/norm/2.;
258 CAp = -CVp;
259 CS = (*(allcoeff_bclnu[LO]))(2)/norm/2.;
260 CSp = (*(allcoeff_bclnu[LO]))(3)/norm/2.;
261 CP = -CS;
262 CPp = -CSp;
263 C7 = 0.;
264 C7p = 0.;
265 CT = (*(allcoeff_bclnu[LO]))(4)/norm/2.;
266 CTp = 0.;
267
268 /* SM + NP Wilson coefficients */
269 if (NPanalysis) {
270 if (lep == StandardModel::TAU) {
271 if (btocNPpmflag) {
272 CV += (mySM.getCCC3() - mySM.getCCC4()) / 2. / M_SQRT2;
273 CVp = (mySM.getCCC3() + mySM.getCCC4()) / 2. / M_SQRT2;
274 CA -= (mySM.getCCC3() - mySM.getCCC4()) / 2. / M_SQRT2;
275 CAp = -(mySM.getCCC3() + mySM.getCCC4()) / 2. / M_SQRT2;
276 CS = (mySM.getCCC1() - mySM.getCCC2()) / 2. / M_SQRT2;
277 CSp = (mySM.getCCC1() + mySM.getCCC2()) / 2. / M_SQRT2;
278 CP = -(mySM.getCCC1() - mySM.getCCC2()) / 2. / M_SQRT2;
279 CPp = -(mySM.getCCC1() + mySM.getCCC2()) / 2. / M_SQRT2;
280 CTp = mySM.getOptionalParameter("CT_NP");
281 } else {
282 CV += mySM.getCCC3() / 2.;
283 CVp = mySM.getCCC4() / 2.;
284 CA -= mySM.getCCC3() / 2.;
285 CAp = -mySM.getCCC4() / 2.;
286 CS = mySM.getCCC1() / 2.;
287 CSp = mySM.getCCC2() / 2.;
288 CP = -mySM.getCCC1() / 2.;
289 CPp = -mySM.getCCC2() / 2.;
290 CTp = mySM.getCCC5();
291 }
292 }
293 }
294
295 switch (vectorM) {
297 if (CLNflag) {
298 hA1w1 = mySM.getOptionalParameter("hA1w1");
299 rho2 = mySM.getOptionalParameter("rho2");
300 R1w1 = mySM.getOptionalParameter("R1w1");
301 R2w1 = mySM.getOptionalParameter("R2w1");
302 N_A = mySM.getOptionalParameter("N_A");
303 N_1 = mySM.getOptionalParameter("N_1");
304 N_2 = mySM.getOptionalParameter("N_2");
305 j_A = mySM.getOptionalParameter("j_A");
306 j_0 = mySM.getOptionalParameter("j_0");
307 j_1 = mySM.getOptionalParameter("j_1");
308 j_2 = mySM.getOptionalParameter("j_2");
309 k_A = mySM.getOptionalParameter("k_A");
310 k_0 = mySM.getOptionalParameter("k_0");
311 k_1 = mySM.getOptionalParameter("k_1");
312 k_2 = mySM.getOptionalParameter("k_2");
313 l_A = mySM.getOptionalParameter("l_A");
314 }
315 else if (BGLflag) {
316 af0 = mySM.getOptionalParameter("af0");
317 af1 = mySM.getOptionalParameter("af1");
318 af2 = mySM.getOptionalParameter("af2");
319 ag0 = mySM.getOptionalParameter("ag0");
320 ag1 = mySM.getOptionalParameter("ag1");
321 ag2 = mySM.getOptionalParameter("ag2");
322 aF11 = mySM.getOptionalParameter("aF11");
323 aF12 = mySM.getOptionalParameter("aF12");
324 aF13 = mySM.getOptionalParameter("aF13");
325 aF21 = mySM.getOptionalParameter("aF21");
326 aF22 = mySM.getOptionalParameter("aF22");
327 aF23 = mySM.getOptionalParameter("aF23");
328 mBcstV1 = mySM.getOptionalParameter("mBcstV1");
329 mBcstV2 = mySM.getOptionalParameter("mBcstV2");
330 mBcstV3 = mySM.getOptionalParameter("mBcstV3");
331 mBcstV4 = mySM.getOptionalParameter("mBcstV4");
332 mBcstA1 = mySM.getOptionalParameter("mBcstA1");
333 mBcstA2 = mySM.getOptionalParameter("mBcstA2");
334 mBcstA3 = mySM.getOptionalParameter("mBcstA3");
335 mBcstA4 = mySM.getOptionalParameter("mBcstA4");
336 mBcstP1 = mySM.getOptionalParameter("mBcstP1");
337 mBcstP2 = mySM.getOptionalParameter("mBcstP2");
338 mBcstP3 = mySM.getOptionalParameter("mBcstP3");
339 chiTV = mySM.getOptionalParameter("chiTV");
340 chiTA = mySM.getOptionalParameter("chiTA");
341 chiTP = mySM.getOptionalParameter("chiTP");
342 nI = mySM.getOptionalParameter("nI");
343
344 zV1 = sqrt((MM+MV)*(MM+MV)-mBcstV1*mBcstV1)-sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV));
345 zV1 /= (sqrt((MM+MV)*(MM+MV)-mBcstV1*mBcstV1)+sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV)));
346 zV2 = sqrt((MM+MV)*(MM+MV)-mBcstV2*mBcstV2)-sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV));
347 zV2 /= (sqrt((MM+MV)*(MM+MV)-mBcstV2*mBcstV2)+sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV)));
348 zV3 = sqrt((MM+MV)*(MM+MV)-mBcstV3*mBcstV3)-sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV));
349 zV3 /= (sqrt((MM+MV)*(MM+MV)-mBcstV3*mBcstV3)+sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV)));
350 zV4 = sqrt((MM+MV)*(MM+MV)-mBcstV4*mBcstV4)-sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV));
351 zV4 /= (sqrt((MM+MV)*(MM+MV)-mBcstV4*mBcstV4)+sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV)));
352
353 zA1 = sqrt((MM+MV)*(MM+MV)-mBcstA1*mBcstA1)-sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV));
354 zA1 /= (sqrt((MM+MV)*(MM+MV)-mBcstA1*mBcstA1)+sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV)));
355 zA2 = sqrt((MM+MV)*(MM+MV)-mBcstA2*mBcstA2)-sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV));
356 zA2 /= (sqrt((MM+MV)*(MM+MV)-mBcstA2*mBcstA2)+sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV)));
357 zA3 = sqrt((MM+MV)*(MM+MV)-mBcstA3*mBcstA3)-sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV));
358 zA3 /= (sqrt((MM+MV)*(MM+MV)-mBcstA3*mBcstA3)+sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV)));
359 zA4 = sqrt((MM+MV)*(MM+MV)-mBcstA4*mBcstA4)-sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV));
360 zA4 /= (sqrt((MM+MV)*(MM+MV)-mBcstA4*mBcstA4)+sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV)));
361
362 zP1 = sqrt((MM+MV)*(MM+MV)-mBcstP1*mBcstP1)-sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV));
363 zP1 /= (sqrt((MM+MV)*(MM+MV)-mBcstP1*mBcstP1)+sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV)));
364 zP2 = sqrt((MM+MV)*(MM+MV)-mBcstP2*mBcstP2)-sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV));
365 zP2 /= (sqrt((MM+MV)*(MM+MV)-mBcstP2*mBcstP2)+sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV)));
366 zP3 = sqrt((MM+MV)*(MM+MV)-mBcstP3*mBcstP3)-sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV));
367 zP3 /= (sqrt((MM+MV)*(MM+MV)-mBcstP3*mBcstP3)+sqrt((MM+MV)*(MM+MV)-(MM-MV)*(MM-MV)));
368 }
369 else if (DMflag) {
370 af_1 = mySM.getOptionalParameter("af_1");
371 ag_1 = mySM.getOptionalParameter("ag_1");
372 aF1_1 = mySM.getOptionalParameter("aF1_1");
373 aF2_1 = mySM.getOptionalParameter("aF2_1");
374 af_2 = mySM.getOptionalParameter("af_2");
375 ag_2 = mySM.getOptionalParameter("ag_2");
376 aF1_2 = mySM.getOptionalParameter("aF1_2");
377 aF2_2 = mySM.getOptionalParameter("aF2_2");
378 af_3 = mySM.getOptionalParameter("af_3");
379 ag_3 = mySM.getOptionalParameter("ag_3");
380 aF1_3 = mySM.getOptionalParameter("aF1_3");
381 aF2_3 = mySM.getOptionalParameter("aF2_3");
382 af_4 = mySM.getOptionalParameter("af_4");
383 ag_4 = mySM.getOptionalParameter("ag_4");
384 aF1_4 = mySM.getOptionalParameter("aF1_4");
385 aF2_4 = mySM.getOptionalParameter("aF2_4");
386 af_5 = mySM.getOptionalParameter("af_5");
387 ag_5 = mySM.getOptionalParameter("ag_5");
388 aF1_5 = mySM.getOptionalParameter("aF1_5");
389 aF2_5 = mySM.getOptionalParameter("aF2_5");
390 af_6 = mySM.getOptionalParameter("af_6");
391 ag_6 = mySM.getOptionalParameter("ag_6");
392 aF1_6 = mySM.getOptionalParameter("aF1_6");
393 aF2_6 = mySM.getOptionalParameter("aF2_6");
394 af_7 = mySM.getOptionalParameter("af_7");
395 ag_7 = mySM.getOptionalParameter("ag_7");
396 aF1_7 = mySM.getOptionalParameter("aF1_7");
397 aF2_7 = mySM.getOptionalParameter("aF2_7");
398 af_8 = mySM.getOptionalParameter("af_8");
399 ag_8 = mySM.getOptionalParameter("ag_8");
400 aF1_8 = mySM.getOptionalParameter("aF1_8");
401 aF2_8 = mySM.getOptionalParameter("aF2_8");
402 af_9 = mySM.getOptionalParameter("af_9");
403 ag_9 = mySM.getOptionalParameter("ag_9");
404 aF1_9 = mySM.getOptionalParameter("aF1_9");
405 aF2_9 = mySM.getOptionalParameter("aF2_9");
406 af_10 = mySM.getOptionalParameter("af_10");
407 ag_10 = mySM.getOptionalParameter("ag_10");
408 aF1_10 = mySM.getOptionalParameter("aF1_10");
409 aF2_10 = mySM.getOptionalParameter("aF2_10");
410 bf_1 = mySM.getOptionalParameter("bf_1");
411 bg_1 = mySM.getOptionalParameter("bg_1");
412 bF1_1 = mySM.getOptionalParameter("bF1_1");
413 bF2_1 = mySM.getOptionalParameter("bF2_1");
414 bf_2 = mySM.getOptionalParameter("bf_2");
415 bg_2 = mySM.getOptionalParameter("bg_2");
416 bF1_2 = mySM.getOptionalParameter("bF1_2");
417 bF2_2 = mySM.getOptionalParameter("bF2_2");
418 bf_3 = mySM.getOptionalParameter("bf_3");
419 bg_3 = mySM.getOptionalParameter("bg_3");
420 bF1_3 = mySM.getOptionalParameter("bF1_3");
421 bF2_3 = mySM.getOptionalParameter("bF2_3");
422 bf_4 = mySM.getOptionalParameter("bf_4");
423 bg_4 = mySM.getOptionalParameter("bg_4");
424 bF1_4 = mySM.getOptionalParameter("bF1_4");
425 bF2_4 = mySM.getOptionalParameter("bF2_4");
426 bf_5 = mySM.getOptionalParameter("bf_5");
427 bg_5 = mySM.getOptionalParameter("bg_5");
428 bF1_5 = mySM.getOptionalParameter("bF1_5");
429 bF2_5 = mySM.getOptionalParameter("bF2_5");
430 bf_6 = mySM.getOptionalParameter("bf_6");
431 bg_6 = mySM.getOptionalParameter("bg_6");
432 bF1_6 = mySM.getOptionalParameter("bF1_6");
433 bF2_6 = mySM.getOptionalParameter("bF2_6");
434 bf_7 = mySM.getOptionalParameter("bf_7");
435 bg_7 = mySM.getOptionalParameter("bg_7");
436 bF1_7 = mySM.getOptionalParameter("bF1_7");
437 bF2_7 = mySM.getOptionalParameter("bF2_7");
438 bf_8 = mySM.getOptionalParameter("bf_8");
439 bg_8 = mySM.getOptionalParameter("bg_8");
440 bF1_8 = mySM.getOptionalParameter("bF1_8");
441 bF2_8 = mySM.getOptionalParameter("bF2_8");
442 bf_9 = mySM.getOptionalParameter("bf_9");
443 bg_9 = mySM.getOptionalParameter("bg_9");
444 bF1_9 = mySM.getOptionalParameter("bF1_9");
445 bF2_9 = mySM.getOptionalParameter("bF2_9");
446 bf_10 = mySM.getOptionalParameter("bf_10");
447 bg_10 = mySM.getOptionalParameter("bg_10");
448 bF1_10 = mySM.getOptionalParameter("bF1_10");
449 bF2_10 = mySM.getOptionalParameter("bF2_10");
450 }
451 else{};
452 break;
453 default:
454 std::stringstream out;
455 out << vectorM;
456 throw std::runtime_error("MVlnu: vector " + out.str() + " not implemented");
457 }
458
459 if ((hA1w1 != hA1w1_cache) || (rho2 != rho2_cache) || (R1w1 != R1w1_cache) || (R2w1 != R2w1_cache)
460 || (N_A != N_A_cache) || (N_1 != N_1_cache) || (N_2 != N_2_cache) || (l_A != l_A_cache)
461 || (j_A != j_A_cache) || (j_0 != j_0_cache) || (j_1 != j_1_cache) || (j_2 != j_2_cache)
462 || (k_A != k_A_cache) || (k_0 != k_0_cache) || (k_1 != k_1_cache) || (k_2 != k_2_cache)
463 || (af0 != af0_cache) || (af1 != af1_cache) || (af2 != af2_cache)
464 || (ag0 != ag0_cache) || (ag1 != af1_cache) || (ag2 != af2_cache)
465 || (aF11 != aF11_cache) || (aF12 != aF12_cache) || (aF13 != aF13_cache)
466 || (aF21 != aF21_cache) || (aF22 != aF22_cache) || (aF23 != aF23_cache)
467 || (af_1 != af_1_cache) || (ag_1 != ag_1_cache) || (aF1_1 != aF1_1_cache) || (aF2_1 != aF2_1_cache)
468 || (af_2 != af_2_cache) || (ag_2 != ag_2_cache) || (aF1_2 != aF1_2_cache) || (aF2_2 != aF2_2_cache)
469 || (af_3 != af_3_cache) || (ag_3 != ag_3_cache) || (aF1_3 != aF1_3_cache) || (aF2_3 != aF2_3_cache)
470 || (af_4 != af_4_cache) || (ag_4 != ag_4_cache) || (aF1_4 != aF1_4_cache) || (aF2_4 != aF2_4_cache)
471 || (af_5 != af_5_cache) || (ag_5 != ag_5_cache) || (aF1_5 != aF1_5_cache) || (aF2_5 != aF2_5_cache)
472 || (af_6 != af_6_cache) || (ag_6 != ag_6_cache) || (aF1_6 != aF1_6_cache) || (aF2_6 != aF2_6_cache)
473 || (af_7 != af_7_cache) || (ag_7 != ag_7_cache) || (aF1_7 != aF1_7_cache) || (aF2_7 != aF2_7_cache)
474 || (af_8 != af_8_cache) || (ag_8 != ag_8_cache) || (aF1_8 != aF1_8_cache) || (aF2_8 != aF2_8_cache)
475 || (af_9 != af_9_cache) || (ag_9 != ag_9_cache) || (aF1_9 != aF1_9_cache) || (aF2_9 != aF2_9_cache)
476 || (af_10 != af_10_cache) || (ag_10 != ag_10_cache) || (aF1_10 != aF1_10_cache) || (aF2_10 != aF2_10_cache)
477 || (bf_1 != bf_1_cache) || (bg_1 != bg_1_cache) || (bF1_1 != bF1_1_cache) || (bF2_1 != bF2_1_cache)
478 || (bf_2 != bf_2_cache) || (bg_2 != bg_2_cache) || (bF1_2 != bF1_2_cache) || (bF2_2 != bF2_2_cache)
479 || (bf_3 != bf_3_cache) || (bg_3 != bg_3_cache) || (bF1_3 != bF1_3_cache) || (bF2_3 != bF2_3_cache)
480 || (bf_4 != bf_4_cache) || (bg_4 != bg_4_cache) || (bF1_4 != bF1_4_cache) || (bF2_4 != bF2_4_cache)
481 || (bf_5 != bf_5_cache) || (bg_5 != bg_5_cache) || (bF1_5 != bF1_5_cache) || (bF2_5 != bF2_5_cache)
482 || (bf_6 != bf_6_cache) || (bg_6 != bg_6_cache) || (bF1_6 != bF1_6_cache) || (bF2_6 != bF2_6_cache)
483 || (bf_7 != bf_7_cache) || (bg_7 != bg_7_cache) || (bF1_7 != bF1_7_cache) || (bF2_7 != bF2_7_cache)
484 || (bf_8 != bf_8_cache) || (bg_8 != bg_8_cache) || (bF1_8 != bF1_8_cache) || (bF2_8 != bF2_8_cache)
485 || (bf_9 != bf_9_cache) || (bg_9 != bg_9_cache) || (bF1_9 != bF1_9_cache) || (bF2_9 != bF2_9_cache)
486 || (bf_10 != bf_10_cache) || (bg_10 != bg_10_cache) || (bF1_10 != bF1_10_cache) || (bF2_10 != bF2_10_cache)
487 || (CS != CS_cache) || (CSp != CSp_cache)
488 || (CP != CP_cache) || (CPp != CPp_cache)
489 || (CV != CV_cache) || (CVp != CVp_cache)
490 || (CA != CA_cache) || (CAp != CAp_cache)
491 || (CT != CT_cache) || (CTp != CTp_cache)) {
492 checkcache_int_tau = false;
493 checkcache_int_mu = false;
494 checkcache_int_el = false;
495 }
496
497 if (!checkcache_int_tau || !checkcache_int_mu || !checkcache_int_el) {
498 if (lep == StandardModel::TAU) {
499 cached_intJ1s_tau = integrateJ(1, q2min, q2max);
500 cached_intJ1c_tau = integrateJ(2, q2min, q2max);
501 cached_intJ2s_tau = integrateJ(3, q2min, q2max);
502 cached_intJ2c_tau = integrateJ(4, q2min, q2max);
503 cached_intJ3_tau = integrateJ(5, q2min, q2max);
504 cached_intJ6s_tau = integrateJ(8, q2min, q2max);
505 cached_intJ6c_tau = integrateJ(9, q2min, q2max);
506 cached_intJ9_tau = integrateJ(12, q2min, q2max);
507 cached_intJ4_tau = 0.;
508 cached_intJ5_tau = 0.;
509 cached_intJ7_tau = 0.;
510 cached_intJ8_tau = 0.;
511 // not needed at present
512 /*
513 cached_intJ4_tau = integrateJ(6,q2min,q2max);
514 cached_intJ5_tau = integrateJ(7,q2min,q2max);
515 cached_intJ7_tau = integrateJ(10,q2min,q2max);
516 cached_intJ8_tau = integrateJ(11,q2min,q2max);
517 */
518 checkcache_int_tau = true;
519 }
520 if (lep == StandardModel::MU) {
521 cached_intJ1s_mu = integrateJ(1, q2min, q2max);
522 cached_intJ1c_mu = integrateJ(2, q2min, q2max);
523 cached_intJ2s_mu = integrateJ(3, q2min, q2max);
524 cached_intJ2c_mu = integrateJ(4, q2min, q2max);
525 cached_intJ3_mu = integrateJ(5, q2min, q2max);
526 cached_intJ6s_mu = integrateJ(8, q2min, q2max);
527 cached_intJ6c_mu = integrateJ(9, q2min, q2max);
528 cached_intJ9_mu = integrateJ(12, q2min, q2max);
529 cached_intJ4_mu = 0.;
530 cached_intJ5_mu = 0.;
531 cached_intJ7_mu = 0.;
532 cached_intJ8_mu = 0.;
533 // not needed at present
534 /*
535 cached_intJ4_mu = integrateJ(6,q2min,q2max);
536 cached_intJ5_mu = integrateJ(7,q2min,q2max);
537 cached_intJ7_mu = integrateJ(10,q2min,q2max);
538 cached_intJ8_mu = integrateJ(11,q2min,q2max);
539 */
540 checkcache_int_mu = true;
541 }
543 cached_intJ1s_el = integrateJ(1, q2min, q2max);
544 cached_intJ1c_el = integrateJ(2, q2min, q2max);
545 cached_intJ2s_el = integrateJ(3, q2min, q2max);
546 cached_intJ2c_el = integrateJ(4, q2min, q2max);
547 cached_intJ3_el = integrateJ(5, q2min, q2max);
548 cached_intJ6s_el = integrateJ(8, q2min, q2max);
549 cached_intJ6c_el = integrateJ(9, q2min, q2max);
550 cached_intJ9_el = integrateJ(12, q2min, q2max);
551 cached_intJ4_el = 0.;
552 cached_intJ5_el = 0.;
553 cached_intJ7_el = 0.;
554 cached_intJ8_el = 0.;
555 // not needed at present
556 /*
557 cached_intJ4_el = integrateJ(6,q2min,q2max);
558 cached_intJ5_el = integrateJ(7,q2min,q2max);
559 cached_intJ7_el = integrateJ(10,q2min,q2max);
560 cached_intJ8_el = integrateJ(11,q2min,q2max);
561 */
562 checkcache_int_el = true;
563 }
564 }
565 if (CLNflag) {
566 hA1w1_cache = hA1w1;
567 rho2_cache = rho2;
568 R1w1_cache = R1w1;
569 R2w1_cache = R2w1;
570 N_A_cache = N_A;
571 N_1_cache =N_1;
572 N_2_cache = N_2;
573 j_A_cache= j_A;
574 j_0_cache = j_0;
575 j_1_cache = j_1;
576 j_2_cache = j_2;
577 k_A_cache = k_A;
578 k_0_cache = k_0;
579 k_1_cache = k_1;
580 k_2_cache = k_2;
581 l_A_cache = l_A;
582 }
583 else if (BGLflag){
584 af0_cache = af0;
585 af1_cache = af1;
586 af2_cache = af2;
587 ag0_cache = ag0;
588 ag1_cache = ag1;
589 ag2_cache = ag2;
590 aF11_cache = aF11;
591 aF12_cache = aF12;
592 aF13_cache = aF13;
593 aF21_cache = aF21;
594 aF22_cache = aF22;
595 aF23_cache = aF23;
596 }
597 else if (DMflag){
598 af_1_cache = af_1;
599 ag_1_cache = ag_1;
600 aF1_1_cache = aF1_1;
601 aF2_1_cache = aF2_1;
602 af_2_cache = af_2;
603 ag_2_cache = ag_2;
604 aF1_2_cache = aF1_2;
605 aF2_2_cache = aF2_2;
606 af_3_cache = af_3;
607 ag_3_cache = ag_3;
608 aF1_3_cache = aF1_3;
609 aF2_3_cache = aF2_3;
610 af_4_cache = af_4;
611 ag_4_cache = ag_4;
612 aF1_4_cache = aF1_4;
613 aF2_4_cache = aF2_4;
614 af_5_cache = af_5;
615 ag_5_cache = ag_5;
616 aF1_5_cache = aF1_5;
617 aF2_5_cache = aF2_5;
618 af_6_cache = af_6;
619 ag_6_cache = ag_6;
620 aF1_6_cache = aF1_6;
621 aF2_6_cache = aF2_6;
622 af_7_cache = af_7;
623 ag_7_cache = ag_7;
624 aF1_7_cache = aF1_7;
625 aF2_7_cache = aF2_7;
626 af_8_cache = af_8;
627 ag_8_cache = ag_8;
628 aF1_8_cache = aF1_8;
629 aF2_8_cache = aF2_8;
630 af_9_cache = af_9;
631 ag_9_cache = ag_9;
632 aF1_9_cache = aF1_9;
633 aF2_9_cache = aF2_9;
634 af_10_cache = af_10;
635 ag_10_cache = ag_10;
636 aF1_10_cache = aF1_10;
637 aF2_10_cache = aF2_10;
638 bf_1_cache = bf_1;
639 bg_1_cache = bg_1;
640 bF1_1_cache = bF1_1;
641 bF2_1_cache = bF2_1;
642 bf_2_cache = bf_2;
643 bg_2_cache = bg_2;
644 bF1_2_cache = bF1_2;
645 bF2_2_cache = bF2_2;
646 bf_3_cache = bf_3;
647 bg_3_cache = bg_3;
648 bF1_3_cache = bF1_3;
649 bF2_3_cache = bF2_3;
650 bf_4_cache = bf_4;
651 bg_4_cache = bg_4;
652 bF1_4_cache = bF1_4;
653 bF2_4_cache = bF2_4;
654 bf_5_cache = bf_5;
655 bg_5_cache = bg_5;
656 bF1_5_cache = bF1_5;
657 bF2_5_cache = bF2_5;
658 bf_6_cache = bf_6;
659 bg_6_cache = bg_6;
660 bF1_6_cache = bF1_6;
661 bF2_6_cache = bF2_6;
662 bf_7_cache = bf_7;
663 bg_7_cache = bg_7;
664 bF1_7_cache = bF1_7;
665 bF2_7_cache = bF2_7;
666 bf_8_cache = bf_8;
667 bg_8_cache = bg_8;
668 bF1_8_cache = bF1_8;
669 bF2_8_cache = bF2_8;
670 bf_9_cache = bf_9;
671 bg_9_cache = bg_9;
672 bF1_9_cache = bF1_9;
673 bF2_9_cache = bF2_9;
674 bf_10_cache = bf_10;
675 bg_10_cache = bg_10;
676 bF1_10_cache = bF1_10;
677 bF2_10_cache = bF2_10;
678 }
679 else{};
680
681 CS_cache = CS;
682 CSp_cache = CSp;
683 CP_cache = CP;
684 CPp_cache = CPp;
685 CV_cache = CV;
686 CVp_cache = CVp;
687 CA_cache = CA;
688 CAp_cache = CAp;
689 CT_cache = CT;
690 CTp_cache = CTp;
691
693
694 return;
695
696}
697
698/*******************************************************************************
699 * Kinematic functions *
700 * ****************************************************************************/
701
702double MVlnu::lambda_half(double a, double b, double c)
703{
704 return sqrt(a*a+b*b+c*c-2.*(a*b+a*c+b*c));
705}
706/*******************************************************************************
707 * Form factors *
708 * ****************************************************************************/
709
710double MVlnu::phi_f(double z)
711{
712 double prefac = 4. * (MV_o_MM) / MM / MM * sqrt(nI / (3. * M_PI * chiTA));
713 double num = (1. + z) * sqrt((1. - z)*(1. - z)*(1. - z));
714 double den = (1. + MV_o_MM)*(1. - z) + 2. * sqrtMV_o_MM * (1. + z);
715 double den4 = den * den * den*den;
716 return prefac * num / den4;
717}
718
719double MVlnu::f_BGL(double q2)
720{
721 double w = w0 - q2 / (2. * MM * MV);
722 double z = (sqrt(w + 1.) - M_SQRT2) / (sqrt(w + 1.) + M_SQRT2);
723 double Pfacf = (z - zA1) / (1. - z * zA1)*(z - zA2) / (1. - z * zA2)*(z - zA3) / (1. - z * zA3)*(z - zA4) / (1. - z * zA4);
724 double phif = phi_f(z);
725 return (af0 + af1 * z + af2 * z * z) / phif / Pfacf;
726}
727
728double MVlnu::phi_g(double z)
729{
730 double prefac = sqrt(nI / (3. * M_PI * chiTV));
731 double num = 16. * (MV_o_MM)*(MV_o_MM)*(1. + z)*(1. + z) / sqrt(1. - z);
732 double den = (1. + MV_o_MM)*(1. - z) + 2. * sqrtMV_o_MM * (1. + z);
733 double den4 = den * den * den*den;
734 return prefac * num / den4;
735}
736
737double MVlnu::g_BGL(double q2)
738{
739 double w = w0 - q2 / (2. * MM * MV);
740 double z = (sqrt(w + 1.) - M_SQRT2) / (sqrt(w + 1.) + M_SQRT2);
741 double Pfacg = (z - zV1) / (1. - z * zV1)*(z - zV2) / (1. - z * zV2)*(z - zV3) / (1. - z * zV3)*(z - zV4) / (1. - z * zV4);
742 double phig = phi_g(z);
743 return (ag0 + ag1 * z + ag2 * z * z) / phig / Pfacg;
744}
745
746double MVlnu::phi_F1(double z)
747{
748 double prefac = 4. * (MV_o_MM) / MM / MM / MM * sqrt(nI / (6. * M_PI * chiTA));
749 double num = (1. + z) * sqrt((1. - z)*(1. - z)*(1. - z)*(1. - z)*(1. - z));
750 double den = (1. + MV_o_MM)*(1. - z) + 2. * sqrtMV_o_MM * (1. + z);
751 double den5 = den * den * den * den*den;
752 return prefac * num / den5;
753}
754
755double MVlnu::F1_BGL(double q2)
756{
757 double w = w0 - q2 / (2. * MM * MV);
758 double z = (sqrt(w + 1.) - M_SQRT2) / (sqrt(w + 1.) + M_SQRT2);
759 double PfacF1 = (z - zA1) / (1. - z * zA1)*(z - zA2) / (1. - z * zA2)*(z - zA3) / (1. - z * zA3)*(z - zA4) / (1. - z * zA4);
760 double phiF1 = phi_F1(z);
761 double aF10 = (MM - MV)*(phi_F1(0.) / phi_f(0.)) * af0; // F1(z=0) = (MM-MV)*f(z=0)
762 return (aF10 + aF11 * z + aF12 * z * z + aF13 * z * z * z) / phiF1 / PfacF1;
763}
764
765double MVlnu::phi_F2(double z)
766{
767 double prefac = 8. * sqrt(2.)*(MV_o_MM)*(MV_o_MM) * sqrt(nI / (M_PI * chiTP));
768 double num = (1. + z)*(1. + z) / sqrt(1. - z);
769 double den = (1. + MV_o_MM)*(1. - z) + 2. * sqrtMV_o_MM * (1. + z);
770 double den4 = den * den * den*den;
771 return prefac * num / den4;
772}
773
774double MVlnu::F2_BGL(double q2)
775{
776 double w = w0 - q2 / (2. * MM * MV);
777 double z = (sqrt(w + 1.) - M_SQRT2) / (sqrt(w + 1.) + M_SQRT2);
778 double z0 = (sqrt(w0 + 1.) - M_SQRT2) / (sqrt(w0 + 1.) + M_SQRT2);
779 double PfacF2 = (z - zP1) / (1. - z * zP1)*(z - zP2) / (1. - z * zP2)*(z - zP3) / (1. - z * zP3);
780 double PfacF2z0 = (z0 - zP1) / (1. - z0 * zP1)*(z0 - zP2) / (1. - z0 * zP2)*(z0 - zP3) / (1. - z0 * zP3);
781 double phiF2 = phi_F2(z);
782 double phiF2z0 = phi_F2(z0);
783 double aF20 = PfacF2z0 * phiF2z0 * 2. * F1_BGL(0.) / (MM * MM - MV * MV) - aF21 * z0 - aF22 * z0*z0; // F2(q2=0) = 2.*F1(q2=0)/(MM*MM-MV*MV)
784 return (aF20 + aF21 * z + aF22 * z * z + aF23 * z * z * z) / phiF2 / PfacF2;
785}
786
787double MVlnu::hA1(double q2)
788{
789 double w = w0 - q2 / (2. * MM * MV);
790 double z = (sqrt(w + 1.) - M_SQRT2) / (sqrt(w + 1.) + M_SQRT2);
791 if (CLNflag) return hA1w1 * N_A * (1. - j_A * 8. * rho2 * z + k_A * (53. * rho2 - 15.) * z * z - l_A * (231. * rho2 - 91.) * z * z * z);
792 else if (BGLflag) return(f_BGL(q2) / sqrt(MM * MV) / (1. + w));
793 else if (DMflag) {
794 double w = w0 - q2 / (2. * MM * MV);
795 double f_fac = 0.;
796 if (w<1.05) {f_fac = af_1 + bf_1 * w;}
797 else if (w<1.10) {f_fac = af_2 + bf_2 * w;}
798 else if (w<1.15) {f_fac = af_3 + bf_3 * w;}
799 else if (w<1.20) {f_fac = af_4 + bf_4 * w;}
800 else if (w<1.25) {f_fac = af_5 + bf_5 * w;}
801 else if (w<1.30) {f_fac = af_6 + bf_6 * w;}
802 else if (w<1.35) {f_fac = af_7 + bf_7 * w;}
803 else if (w<1.40) {f_fac = af_8 + bf_8 * w;}
804 else if (w<1.45) {f_fac = af_9 + bf_9 * w;}
805 else {f_fac = af_10 + bf_10 * w;}
806 return f_fac / sqrt(MM * MV) / (1. + w);
807 }
808 else return 0.;
809}
810
811double MVlnu::R1(double q2)
812{
813 double w = w0 - q2 / (2. * MM * MV);
814 return N_1 * R1w1 - j_1 * 0.12 * (w - 1.) + k_1 * 0.05 * (w - 1.)*(w - 1.);
815}
816
817double MVlnu::R2(double q2)
818{
819 double w = w0 - q2 / (2. * MM * MV);
820 return N_2 * R2w1 + j_2 * 0.11 * (w - 1.) - k_2 * 0.06 * (w - 1.)*(w - 1.);
821}
822
823double MVlnu::R0(double q2)
824{
825 double w = w0 - q2 / (2. * MM * MV);
826 /* form factor relation among A0, A1 and A2 at q2=0 */
827 double R2q2at0 = R2(0.);
828 double R0q2at0 = (MM + MV - (MM - MV) * R2q2at0) / (2. * MV);
829 // caveat: HQET rel at the kinematic endpoint, q2 = 0 ...
830 double R0w1 = R0q2at0 + j_0 * 0.11 * (w0 - 1.) - k_0 * 0.01 * (w0 - 1.)*(w0 - 1.);
831 // one may consider "lattice" R0w1 = 1.14 +- O(10%) + consistency rel at q2 = 0 ...
832 return R0w1 - j_0 * 0.11 * (w - 1.) + k_0 * 0.01 * (w - 1.)*(w - 1.);
833}
834
835double MVlnu::V(double q2)
836{
837 if (CLNflag) return R1(q2) / RV * hA1(q2);
838 else if (BGLflag) return (MM + MV) * g_BGL(q2) / 2.;
839 else if (DMflag) {
840 double w = w0 - q2 / (2. * MM * MV);
841 double g_fac = 0.;
842 if (w<1.05) {g_fac = ag_1 + bg_1 * w;}
843 else if (w<1.10) {g_fac = ag_2 + bg_2 * w;}
844 else if (w<1.15) {g_fac = ag_3 + bg_3 * w;}
845 else if (w<1.20) {g_fac = ag_4 + bg_4 * w;}
846 else if (w<1.25) {g_fac = ag_5 + bg_5 * w;}
847 else if (w<1.30) {g_fac = ag_6 + bg_6 * w;}
848 else if (w<1.35) {g_fac = ag_7 + bg_7 * w;}
849 else if (w<1.40) {g_fac = ag_8 + bg_8 * w;}
850 else if (w<1.45) {g_fac = ag_9 + bg_9 * w;}
851 else {g_fac = ag_10 + bg_10 * w;}
852 return (MM + MV) * g_fac / 2.;
853 }
854 else return 0.;
855}
856
857double MVlnu::A0(double q2)
858{
859 double w = w0 - q2 / (2. * MM * MV);
860 // A0 = RV * P1 = R0 * A1
861 if (CLNflag) return R0(q2) * (w + 1.) * RV / 2. * hA1(q2);
862 // F2 = P1 / (RV / 2.)
863 else if (BGLflag) return F2_BGL(q2) / 2.;
864 else if (DMflag) {
865 double w = w0 - q2 / (2. * MM * MV);
866 double F2_fac = 0.;
867 if (w<1.05) {F2_fac = aF2_1 + bF2_1 * w;}
868 else if (w<1.10) {F2_fac = aF2_2 + bF2_2 * w;}
869 else if (w<1.15) {F2_fac = aF2_3 + bF2_3 * w;}
870 else if (w<1.20) {F2_fac = aF2_4 + bF2_4 * w;}
871 else if (w<1.25) {F2_fac = aF2_5 + bF2_5 * w;}
872 else if (w<1.30) {F2_fac = aF2_6 + bF2_6 * w;}
873 else if (w<1.35) {F2_fac = aF2_7 + bF2_7 * w;}
874 else if (w<1.40) {F2_fac = aF2_8 + bF2_8 * w;}
875 else if (w<1.45) {F2_fac = aF2_9 + bF2_9 * w;}
876 else {F2_fac = aF2_10 + bF2_10 * w;}
877 return F2_fac / 2.;
878 }
879 else return 0.;
880}
881
882double MVlnu::A1(double q2)
883{
884 /* form factor in 1:1 with hA1 */
885 double w = w0 - q2 / (2. * MM * MV);
886 return (w + 1.) * RV / 2. * hA1(q2);
887}
888
889double MVlnu::A2(double q2)
890{
891 double w = w0 - q2 / (2. * MM * MV);
892 if (CLNflag) return R2(q2) / RV * hA1(q2);
893 else if (BGLflag) return (MM + MV) / 2. / (w * w - 1.) / MM / MV * ((w - MV_o_MM) * f_BGL(q2) - F1_BGL(q2) / MM);
894 else if (DMflag) {
895 double w = w0 - q2 / (2. * MM * MV);
896 double f_fac = 0.;
897 double F1_fac = 0.;
898 if (w<1.05) {f_fac = af_1 + bf_1 * w; F1_fac = aF1_1 + bF1_1 * w;}
899 else if (w<1.10) {f_fac = af_2 + bf_2 * w; F1_fac = aF1_2 + bF1_2 * w;}
900 else if (w<1.15) {f_fac = af_3 + bf_3 * w; F1_fac = aF1_3 + bF1_3 * w;}
901 else if (w<1.20) {f_fac = af_4 + bf_4 * w; F1_fac = aF1_4 + bF1_4 * w;}
902 else if (w<1.25) {f_fac = af_5 + bf_5 * w; F1_fac = aF1_5 + bF1_5 * w;}
903 else if (w<1.30) {f_fac = af_6 + bf_6 * w; F1_fac = aF1_6 + bF1_6 * w;}
904 else if (w<1.35) {f_fac = af_7 + bf_7 * w; F1_fac = aF1_7 + bF1_7 * w;}
905 else if (w<1.40) {f_fac = af_8 + bf_8 * w; F1_fac = aF1_8 + bF1_8 * w;}
906 else if (w<1.45) {f_fac = af_9 + bf_9 * w; F1_fac = aF1_9 + bF1_9 * w;}
907 else {f_fac = af_10 + bf_10 * w; F1_fac = aF1_10 + bF1_10 * w;}
908 return (MM + MV) / 2. / (w * w - 1.) / MM / MV * ((w - MV_o_MM) * f_fac - F1_fac / MM);
909 }
910 else return 0.;
911}
912
913double MVlnu::A12(double q2)
914{
915 return (A1(q2)*(MM + MV)*(MM + MV)*(MM * MM - MV * MV - q2) -
916 A2(q2)*(MM * MM * MM * MM + (MV * MV - q2)*(MV * MV - q2) -
917 2. * MM * MM * (MV * MV + q2))) / (16. * MM * MV * MV * (MM + MV));
918}
919
920double MVlnu::T1(double q2)
921{
922 double delta_T1 = 0.;
923 return (Mb + Mc) / (MM + MV) * V(q2)*(1. + delta_T1);
924}
925
926double MVlnu::T2(double q2)
927{
928 double delta_T2 = 0.;
929 return (Mb - Mc) / (MM - MV) * A1(q2)*(1. + delta_T2);
930}
931
932double MVlnu::T23(double q2)
933{
934 double delta_T23 = 0.;
935 return ((Mb - Mc)*((MM - MV)*(MM - MV) - q2)*((MM + MV)*(MM + MV) - q2) * A0(q2) +
936 8 * MM * MV * (MV * MV - MM * MM) * A12(q2)) / (4. * MM * (MV - MM) * MV * q2)*(1. + delta_T23);
937}
938/********************************************************************************
939 * Helicity amplitudes (normalization such that all H \propto (mass scale)^-1) *
940 * *****************************************************************************/
941
942gslpp::complex MVlnu::HV0(double q2)
943{
944 return 4. * gslpp::complex::i() * MM * MV / (sqrt(q2)*(MM + MV))*((CV - CVp)*(MM + MV) * A12(q2) + Mb * (C7 - C7p) * T23(q2));
945}
946
947gslpp::complex MVlnu::HVp(double q2)
948{
949 return gslpp::complex::i()*((((CV + CVp) * lambda_half(MM*MM, MV*MV, q2) * V(q2)-(MM + MV)*(MM + MV)*(CV - CVp) * A1(q2))) / (2. * (MM + MV))
950 + (Mb / q2)*((C7 + C7p) * lambda_half(MM*MM, MV*MV, q2) * T1(q2)-(C7 - C7p)*(MM * MM - MV * MV) * T2(q2)));
951}
952
953gslpp::complex MVlnu::HVm(double q2)
954{
955 return gslpp::complex::i()*(((-(CV + CVp) * lambda_half(MM*MM, MV*MV, q2) * V(q2)-(MM + MV)*(MM + MV)*(CV - CVp) * A1(q2))) / (2. * (MM + MV))
956 + (Mb / q2)*(-(C7 + C7p) * lambda_half(MM*MM, MV*MV, q2) * T1(q2)-(C7 - C7p)*(MM * MM - MV * MV) * T2(q2)));
957}
958
959gslpp::complex MVlnu::HAp(double q2)
960{
961 return gslpp::complex::i()*((CA + CAp) * lambda_half(MM*MM, MV*MV, q2) * V(q2)-(MM + MV)*(MM + MV)*(CA - CAp) * A1(q2)) / (2. * (MM + MV));
962}
963
964gslpp::complex MVlnu::HAm(double q2)
965{
966 return gslpp::complex::i()*(-(CA + CAp) * lambda_half(MM*MM, MV*MV, q2) * V(q2)-(MM + MV)*(MM + MV)*(CA - CAp) * A1(q2)) / (2. * (MM + MV));
967}
968
969gslpp::complex MVlnu::HA0(double q2)
970{
971 return 4. * gslpp::complex::i() * MV * MM / (sqrt(q2))*(CA - CAp) * A12(q2);
972}
973
974gslpp::complex MVlnu::HP(double q2)
975{
976 return gslpp::complex::i() * lambda_half(MM*MM, MV*MV, q2) / 2. * ((CP - CPp) / (Mb + Mc)+(Mlep + Mnu) / q2 * (CA - CAp)) * A0(q2);
977}
978
979gslpp::complex MVlnu::HS(double q2)
980{
981 return gslpp::complex::i() * lambda_half(MM*MM, MV*MV, q2) / 2. * ((CS - CSp) / (Mb + Mc)+(Mlep - Mnu) / q2 * (CV - CVp)) * A0(q2);
982}
983
984gslpp::complex MVlnu::HT0(double q2)
985{
986 return 2. * M_SQRT2 * (MM * MV) / (MM + MV)*(CT + CTp) * T23(q2);
987}
988
989gslpp::complex MVlnu::HT0t(double q2)
990{
991 return 2. * (MM * MV) / (MM + MV)*(CT - CTp) * T23(q2);
992}
993
994gslpp::complex MVlnu::HTp(double q2)
995{
996 return ((CT - CTp) * lambda_half(MM*MM, MV*MV, q2) * T1(q2)-(CT + CTp)*(MM * MM - MV * MV) * T2(q2)) / (sqrt(2. * q2));
997}
998
999gslpp::complex MVlnu::HTpt(double q2)
1000{
1001 return ((CT + CTp) * lambda_half(MM*MM, MV*MV, q2) * T1(q2)-(CT - CTp)*(MM * MM - MV * MV) * T2(q2)) / (2. * sqrt(q2));
1002}
1003
1004gslpp::complex MVlnu::HTm(double q2)
1005{
1006 return (-(CT - CTp) * lambda_half(MM*MM, MV*MV, q2) * T1(q2)-(CT + CTp)*(MM * MM - MV * MV) * T2(q2)) / (sqrt(2. * q2));
1007}
1008
1009gslpp::complex MVlnu::HTmt(double q2)
1010{
1011 return (-(CT + CTp) * lambda_half(MM*MM, MV*MV, q2) * T1(q2)-(CT - CTp)*(MM * MM - MV * MV) * T2(q2)) / (2. * sqrt(q2));
1012}
1013/*******************************************************************************
1014 * Generalized angular coefficients (see 1506.03970) *
1015 * ****************************************************************************/
1016
1017gslpp::complex MVlnu::G000(double q2)
1018{
1019 double lambda_MM = lambda_half(MM*MM, MV*MV, q2);
1020 double lambda_lep = lambda_half(Mlep*Mlep, Mnu*Mnu, q2);
1021 double lambda_lep2 = lambda_lep*lambda_lep;
1022 double Elep = sqrt(Mlep * Mlep + lambda_lep2 / (4. * q2));
1023 double Enu = sqrt(Mnu * Mnu + lambda_lep2 / (4. * q2));
1024 double Gprefactor = lambda_MM * lambda_lep / q2;
1025
1026 return Gprefactor * (4. / 9. * (3. * Elep * Enu + lambda_lep2 / (4. * q2))*(HVp(q2).abs2() + HVm(q2).abs2() + HV0(q2).abs2() + HAp(q2).abs2() + HAm(q2).abs2() + HA0(q2).abs2()) +
1027 4. * Mlep * Mnu / 3. * (HVp(q2).abs2() + HVm(q2).abs2() + HV0(q2).abs2() - HAp(q2).abs2() - HAm(q2).abs2() - HA0(q2).abs2()) +
1028 4. / 3. * ((Elep * Enu - Mlep * Mnu + lambda_lep2 / (4. * q2)) * HS(q2).abs2()+(Elep * Enu + Mlep * Mnu + lambda_lep2 / (4. * q2)) * HP(q2).abs2()) +
1029 16. / 9. * (3. * (Elep * Enu + Mlep * Mnu) - lambda_lep2 / (4. * q2))*(HTpt(q2).abs2() + HTmt(q2).abs2() + HT0t(q2).abs2()) +
1030 8. / 9. * (3. * (Elep * Enu - Mlep * Mnu) - lambda_lep2 / (4. * q2))*(HTp(q2).abs2() + HTm(q2).abs2() + HT0(q2).abs2()) +
1031 16. / 3. * (Mlep * Enu + Mnu * Elep)*(HVp(q2) * HTpt(q2).conjugate() + HVm(q2) * HTmt(q2).conjugate() + HV0(q2) * HT0t(q2).conjugate()).imag() +
1032 8. * M_SQRT2 / 3. * (Mlep * Enu - Mnu * Elep)*(HAp(q2) * HTp(q2).conjugate() + HAm(q2) * HTm(q2).conjugate() + HA0(q2) * HT0(q2).conjugate()).imag());
1033}
1034
1035gslpp::complex MVlnu::G010(double q2)
1036{
1037 double lambda_MM = lambda_half(MM*MM, MV*MV, q2);
1038 double lambda_lep = lambda_half(Mlep*Mlep, Mnu*Mnu, q2);
1039 double Gprefactor = lambda_MM * lambda_lep / q2;
1040
1041 return Gprefactor * 4. / 3. * lambda_lep * ((HVp(q2) * HAp(q2).conjugate() - HVm(q2) * HAm(q2).conjugate()).real()+
1042 (2. * M_SQRT2) / q2 * (Mlep * Mlep - Mnu * Mnu)*(HTp(q2) * HTpt(q2).conjugate() - HTm(q2) * HTmt(q2).conjugate()).real() +
1043 2. * (Mlep + Mnu) / sqrt(q2)*(HAp(q2) * HTpt(q2).conjugate() - HAm(q2) * HTmt(q2).conjugate()).imag() +
1044 M_SQRT2 * (Mlep - Mnu) / sqrt(q2)*(HVp(q2) * HTp(q2).conjugate() - HVm(q2) * HTm(q2).conjugate()).imag()-
1045 (Mlep - Mnu) / sqrt(q2)*(HA0(q2) * HP(q2).conjugate()).real()-(Mlep + Mnu) / sqrt(q2)*(HV0(q2) * HS(q2).conjugate()).real()+
1046 (M_SQRT2 * HT0(q2) * HP(q2).conjugate() + 2. * HT0t(q2) * HS(q2).conjugate()).imag());
1047
1048}
1049
1050gslpp::complex MVlnu::G020(double q2)
1051{
1052 double lambda_MM = lambda_half(MM*MM, MV*MV, q2);
1053 double lambda_lep = lambda_half(Mlep*Mlep, Mnu*Mnu, q2);
1054 double lambda_lep2 = lambda_lep*lambda_lep;
1055 double Gprefactor = lambda_MM * lambda_lep / q2;
1056
1057 return -Gprefactor * 2. / 9. * lambda_lep2 / q2 * ((-HVp(q2).abs2() - HVm(q2).abs2() + 2. * HV0(q2).abs2() - HAp(q2).abs2() - HAm(q2).abs2() + 2. * HA0(q2).abs2()) -
1058 2. * (2. * HT0(q2).abs2() - HTp(q2).abs2() - HTm(q2).abs2()) - 4. * (2. * HT0t(q2).abs2() - HTpt(q2).abs2() - HTmt(q2).abs2()));
1059}
1060
1061gslpp::complex MVlnu::G200(double q2)
1062{
1063 double lambda_MM = lambda_half(MM*MM, MV*MV, q2);
1064 double lambda_lep = lambda_half(Mlep*Mlep, Mnu*Mnu, q2);
1065 double lambda_lep2 = lambda_lep*lambda_lep;
1066 double Elep = sqrt(Mlep * Mlep + lambda_lep2 / (4. * q2));
1067 double Enu = sqrt(Mnu * Mnu + lambda_lep2 / (4. * q2));
1068 double Gprefactor = lambda_MM * lambda_lep / q2;
1069
1070 return Gprefactor * (-4. / 9. * (3. * Elep * Enu + lambda_lep2 / (4. * q2))*(HVp(q2).abs2() + HVm(q2).abs2() - 2. * HV0(q2).abs2() + HAp(q2).abs2() + HAm(q2).abs2() - 2. * HA0(q2).abs2()) -
1071 4. / 3. * Mlep * Mnu * (HVp(q2).abs2() + HVm(q2).abs2() - 2. * HV0(q2).abs2() - HAp(q2).abs2() - HAm(q2).abs2() + 2. * HA0(q2).abs2()) +
1072 8. / 3. * (Elep * Enu - Mlep * Mnu + lambda_lep2 / (4. * q2)) * HS(q2).abs2() + 8. / 3. * (Elep * Enu + Mlep * Mnu + lambda_lep2 / (4. * q2)) * HP(q2).abs2() -
1073 16. / 9. * (3. * (Elep * Enu + Mlep * Mnu) - lambda_lep2 / (4. * q2))*(HTpt(q2).abs2() + HTmt(q2).abs2() - 2. * HT0t(q2).abs2()) - 8. / 9. * (3. * (Elep * Enu - Mlep * Mnu) - lambda_lep2 / (4. * q2))*
1074 (HTp(q2).abs2() + HTm(q2).abs2() - 2 * HT0(q2).abs2()) - 16. / 3. * (Mlep * Enu + Mnu * Elep)*(HVp(q2) * HTpt(q2).conjugate() + HVm(q2) * HTmt(q2).conjugate() - 2. * HV0(q2) * HT0t(q2).conjugate()).imag() -
1075 8. * M_SQRT2 / 3. * (Mlep * Enu - Mnu * Elep)*(HAp(q2) * HTp(q2).conjugate() + HAm(q2) * HTm(q2).conjugate() - 2. * HA0(q2) * HT0(q2).conjugate()).imag());
1076}
1077
1078gslpp::complex MVlnu::G210(double q2)
1079{
1080 double lambda_MM = lambda_half(MM*MM, MV*MV, q2);
1081 double lambda_lep = lambda_half(Mlep*Mlep, Mnu*Mnu, q2);
1082 double Gprefactor = lambda_MM * lambda_lep / q2;
1083
1084 return -Gprefactor * 4. * lambda_lep / 3. * ((HVp(q2) * HAp(q2).conjugate() - HVm(q2) * HAm(q2).conjugate()).real() +
1085 2. * M_SQRT2 * (Mlep * Mlep - Mnu * Mnu) / q2 * (HTp(q2) * HTpt(q2).conjugate() - HTm(q2) * HTmt(q2).conjugate()).real() +
1086 2. * (Mlep + Mnu) / sqrt(q2)*(HAp(q2) * HTpt(q2).conjugate() - HAm(q2) * HTmt(q2).conjugate()).imag() +
1087 M_SQRT2 * (Mlep - Mnu) / sqrt(q2)*(HVp(q2) * HTp(q2).conjugate() - HVm(q2) * HTm(q2).conjugate()).imag() +
1088 2. * (Mlep - Mnu) / sqrt(q2)*(HA0(q2) * HP(q2).conjugate()).real() + 2. * (Mlep + Mnu) / sqrt(q2)*(HV0(q2) * HS(q2).conjugate()).real() -
1089 2. * (M_SQRT2 * HT0(q2) * HP(q2).conjugate() + 2. * HT0t(q2) * HS(q2).conjugate()).imag());
1090}
1091
1092gslpp::complex MVlnu::G220(double q2)
1093{
1094 double lambda_MM = lambda_half(MM*MM, MV*MV, q2);
1095 double lambda_lep = lambda_half(Mlep*Mlep, Mnu*Mnu, q2);
1096 double lambda_lep2 = lambda_lep*lambda_lep;
1097 double Gprefactor = lambda_MM * lambda_lep / q2;
1098
1099 return -Gprefactor * 2. / 9. * lambda_lep2 / q2 * (HVp(q2).abs2() + HVm(q2).abs2() + 4. * HV0(q2).abs2() + HAp(q2).abs2() +
1100 HAm(q2).abs2() + 4. * HA0(q2).abs2() - 2. * (HTp(q2).abs2() + HTm(q2).abs2() + 4. * HT0(q2).abs2()) -
1101 4. * (HTpt(q2).abs2() + HTmt(q2).abs2() + 4. * HT0t(q2).abs2()));
1102}
1103
1104gslpp::complex MVlnu::G211(double q2)
1105{
1106 double lambda_MM = lambda_half(MM*MM, MV*MV, q2);
1107 double lambda_lep = lambda_half(Mlep*Mlep, Mnu*Mnu, q2);
1108 double Gprefactor = lambda_MM * lambda_lep / q2;
1109
1110 return Gprefactor * 4. / sqrt(3.) * lambda_lep * (HVp(q2) * HA0(q2).conjugate() + HAp(q2) * HV0(q2).conjugate() -
1111 HV0(q2) * HAm(q2).conjugate() - HA0(q2) * HVm(q2).conjugate()+(Mlep + Mnu) / sqrt(q2)*(HVp(q2) * HS(q2).conjugate() + HS(q2) * HVm(q2).conjugate()) -
1112 gslpp::complex::i() * M_SQRT2 * (HP(q2) * HTm(q2).conjugate() - HTp(q2) * HP(q2).conjugate() + M_SQRT2 * (HS(q2) * HTmt(q2).conjugate() - HTpt(q2) * HS(q2).conjugate()))+
1113 (Mlep - Mnu) / sqrt(q2)*(HAp(q2) * HP(q2).conjugate() + HP(q2) * HAm(q2).conjugate()) -
1114 gslpp::complex::i()*2. * (Mlep + Mnu) / sqrt(q2)*(HAp(q2) * HT0t(q2).conjugate() + HT0t(q2) * HAm(q2).conjugate() - HTpt(q2) * HA0(q2).conjugate() - HA0(q2) * HTmt(q2).conjugate()) -
1115 gslpp::complex::i() * M_SQRT2 * (Mlep - Mnu) / sqrt(q2)*(HVp(q2) * HT0(q2).conjugate() + HT0(q2) * HVm(q2).conjugate() - HTp(q2) * HV0(q2).conjugate() - HV0(q2) * HTm(q2).conjugate()) +
1116 2. * M_SQRT2 * (Mlep * Mlep - Mnu * Mnu) / q2 * (HTp(q2) * HT0t(q2).conjugate() + HTpt(q2) * HT0(q2).conjugate() - HT0(q2) * HTmt(q2).conjugate() - HT0t(q2) * HTm(q2).conjugate()));
1117}
1118
1119gslpp::complex MVlnu::G221(double q2)
1120{
1121 double lambda_MM = lambda_half(MM*MM, MV*MV, q2);
1122 double lambda_lep = lambda_half(Mlep*Mlep, Mnu*Mnu, q2);
1123 double lambda_lep2 = lambda_lep*lambda_lep;
1124 double Gprefactor = lambda_MM * lambda_lep / q2;
1125
1126 return Gprefactor * 4. / 3. * lambda_lep2 / q2 * (HVp(q2) * HV0(q2).conjugate() + HV0(q2) * HVm(q2).conjugate() +
1127 HAp(q2) * HA0(q2).conjugate() + HA0(q2) * HAm(q2).conjugate() - 2. * (HTp(q2) * HT0(q2).conjugate() +
1128 HT0(q2) * HTm(q2).conjugate() + 2. * (HTpt(q2) * HT0t(q2).conjugate() + HT0t(q2) * HTmt(q2).conjugate())));
1129}
1130
1131gslpp::complex MVlnu::G222(double q2)
1132{
1133 double lambda_MM = lambda_half(MM*MM, MV*MV, q2);
1134 double lambda_lep = lambda_half(Mlep*Mlep, Mnu*Mnu, q2);
1135 double lambda_lep2 = lambda_lep*lambda_lep;
1136 double Gprefactor = lambda_MM * lambda_lep / q2;
1137
1138 return -Gprefactor * 8. / 3. * lambda_lep2 / q2 * (HVp(q2) * HVm(q2).conjugate() + HAp(q2) * HAm(q2).conjugate() -
1139 2. * (HTp(q2) * HTm(q2).conjugate() + 2. * HTpt(q2) * HTmt(q2).conjugate()));
1140}
1141/***************************************************************************
1142 * 12 independent J angular coefficients (see again 1506.03970) *
1143 * ************************************************************************/
1144
1145double MVlnu::J1s(double q2)
1146{
1147 if ((q2 < Mlep * Mlep) or (q2 > (MM - MV)*(MM - MV))) return 0.;
1148 return amplsq_factor * (8. * G000(q2) + 2. * G020(q2) - 4. * G200(q2) - G220(q2)).real() / 3.;
1149}
1150
1151double MVlnu::J1c(double q2)
1152{
1153 if ((q2 < Mlep * Mlep) or (q2 > (MM - MV)*(MM - MV))) return 0.;
1154 return amplsq_factor * (8. * G000(q2) + 2. * G020(q2) + 8. * G200(q2) + 2. * G220(q2)).real() / 3.;
1155}
1156
1157double MVlnu::J2s(double q2)
1158{
1159 if ((q2 < Mlep * Mlep) or (q2 > (MM - MV)*(MM - MV))) return 0.;
1160 return amplsq_factor * (2. * G020(q2) - G220(q2)).real();
1161}
1162
1163double MVlnu::J2c(double q2)
1164{
1165 if ((q2 < Mlep * Mlep) or (q2 > (MM - MV)*(MM - MV))) return 0.;
1166 return amplsq_factor * (2. * (G020(q2) + G220(q2))).real();
1167}
1168
1169double MVlnu::J3(double q2)
1170{
1171 if ((q2 < Mlep * Mlep) or (q2 > (MM - MV)*(MM - MV))) return 0.;
1172 return amplsq_factor * (G222(q2).real());
1173}
1174
1175double MVlnu::J4(double q2)
1176{
1177 if ((q2 < Mlep * Mlep) or (q2 > (MM - MV)*(MM - MV))) return 0.;
1178 return -amplsq_factor * (G221(q2).real());
1179}
1180
1181double MVlnu::J5(double q2)
1182{
1183 if ((q2 < Mlep * Mlep) or (q2 > (MM - MV)*(MM - MV))) return 0.;
1184 return amplsq_factor * (2. * G211(q2).real() / sqrt(3.));
1185}
1186
1187double MVlnu::J6s(double q2)
1188{
1189 if ((q2 < Mlep * Mlep) or (q2 > (MM - MV)*(MM - MV))) return 0.;
1190 return -amplsq_factor * (4. * (2. * G010(q2) - G210(q2)).real() / 3.);
1191}
1192
1193double MVlnu::J6c(double q2)
1194{
1195 if (q2 < Mlep * Mlep) return 0.;
1196 if (q2 > (MM - MV)*(MM - MV)) return 0.;
1197 return -amplsq_factor * (8. * (G010(q2) + G210(q2)).real() / 3.);
1198}
1199
1200double MVlnu::J7(double q2)
1201{
1202 if ((q2 < Mlep * Mlep) or (q2 > (MM - MV)*(MM - MV))) return 0.;
1203 return -amplsq_factor * (2. * sqrt(3.)*(G211(q2).imag()) / 3.);
1204}
1205
1206double MVlnu::J8(double q2)
1207{
1208 if ((q2 < Mlep * Mlep) or (q2 > (MM - MV)*(MM - MV))) return 0.;
1209 return amplsq_factor * (G221(q2).imag());
1210}
1211
1212double MVlnu::J9(double q2)
1213{
1214 if ((q2 < Mlep * Mlep) or (q2 > (MM - MV)*(MM - MV))) return 0.;
1215 return -amplsq_factor * (G222(q2).imag());
1216}
1217/***************************************************************************
1218 * Integration of angular coefficients Js *
1219 * ************************************************************************/
1220
1221double MVlnu::integrateJ(int i, double q2_min, double q2_max)
1222{
1223 switch (i) {
1224 case 1:
1225 if (lep == StandardModel::TAU) if (checkcache_int_tau && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ1s_tau;
1226 if (lep == StandardModel::MU) if (checkcache_int_mu && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ1s_mu;
1227 if (lep == StandardModel::ELECTRON) if (checkcache_int_el && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ1s_el;
1228 wf=ROOT::Math::Functor1D(&(*this),&MVlnu::J1s);
1229 ig.SetFunction(wf);
1230 J_res = ig.Integral(q2_min, q2_max);
1231 return J_res;
1232 break;
1233 case 2:
1234 if (lep == StandardModel::TAU) if (checkcache_int_tau && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ1c_tau;
1235 if (lep == StandardModel::MU) if (checkcache_int_mu && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ1c_mu;
1236 if (lep == StandardModel::ELECTRON) if (checkcache_int_el && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ1c_el;
1237 wf=ROOT::Math::Functor1D(&(*this),&MVlnu::J1c);
1238 ig.SetFunction(wf);
1239 J_res = ig.Integral(q2_min, q2_max);
1240 return J_res;
1241 break;
1242 case 3:
1243 if (lep == StandardModel::TAU) if (checkcache_int_tau && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ2s_tau;
1244 if (lep == StandardModel::MU) if (checkcache_int_mu && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ2s_mu;
1245 if (lep == StandardModel::ELECTRON) if (checkcache_int_el && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ2s_el;
1246 wf=ROOT::Math::Functor1D(&(*this),&MVlnu::J2s);
1247 ig.SetFunction(wf);
1248 J_res = ig.Integral(q2_min, q2_max);
1249 return J_res;
1250 break;
1251 case 4:
1252 if (lep == StandardModel::TAU) if (checkcache_int_tau && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ2c_tau;
1253 if (lep == StandardModel::MU) if (checkcache_int_mu && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ2c_mu;
1254 if (lep == StandardModel::ELECTRON) if (checkcache_int_el && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ2c_el;
1255 wf=ROOT::Math::Functor1D(&(*this),&MVlnu::J2c);
1256 ig.SetFunction(wf);
1257 J_res = ig.Integral(q2_min, q2_max);
1258 return J_res;
1259 break;
1260 case 5:
1261 if (lep == StandardModel::TAU) if (checkcache_int_tau && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ3_tau;
1262 if (lep == StandardModel::MU) if (checkcache_int_mu && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ3_mu;
1263 if (lep == StandardModel::ELECTRON) if (checkcache_int_el && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ3_el;
1264 wf=ROOT::Math::Functor1D(&(*this),&MVlnu::J3);
1265 ig.SetFunction(wf);
1266 J_res = ig.Integral(q2_min, q2_max);
1267 return J_res;
1268 break;
1269 case 6:
1270 if (lep == StandardModel::TAU) if (checkcache_int_tau && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ4_tau;
1271 if (lep == StandardModel::MU) if (checkcache_int_mu && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ4_mu;
1272 if (lep == StandardModel::ELECTRON) if (checkcache_int_el && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ4_el;
1273 wf=ROOT::Math::Functor1D(&(*this),&MVlnu::J4);
1274 ig.SetFunction(wf);
1275 J_res = ig.Integral(q2_min, q2_max);
1276 return J_res;
1277 break;
1278 case 7:
1279 if (lep == StandardModel::TAU) if (checkcache_int_tau && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ5_tau;
1280 if (lep == StandardModel::MU) if (checkcache_int_mu && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ5_mu;
1281 if (lep == StandardModel::ELECTRON) if (checkcache_int_el && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ5_el;
1282 wf=ROOT::Math::Functor1D(&(*this),&MVlnu::J5);
1283 ig.SetFunction(wf);
1284 J_res = ig.Integral(q2_min, q2_max);
1285 return J_res;
1286 break;
1287 case 8:
1288 if (lep == StandardModel::TAU) if (checkcache_int_tau && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ6s_tau;
1289 if (lep == StandardModel::MU) if (checkcache_int_mu && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ6s_mu;
1290 if (lep == StandardModel::ELECTRON) if (checkcache_int_el && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ6s_el;
1291 wf=ROOT::Math::Functor1D(&(*this),&MVlnu::J6s);
1292 ig.SetFunction(wf);
1293 J_res = ig.Integral(q2_min, q2_max);
1294 return J_res;
1295 break;
1296 case 9:
1297 if (lep == StandardModel::TAU) if (checkcache_int_tau && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ6c_mu;
1298 if (lep == StandardModel::MU) if (checkcache_int_mu && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ6c_mu;
1299 if (lep == StandardModel::ELECTRON) if (checkcache_int_el && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ6c_el;
1300 wf=ROOT::Math::Functor1D(&(*this),&MVlnu::J6c);
1301 ig.SetFunction(wf);
1302 J_res = ig.Integral(q2_min, q2_max);
1303 return J_res;
1304 break;
1305 case 10:
1306 if (lep == StandardModel::TAU) if (checkcache_int_tau && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ7_tau;
1307 if (lep == StandardModel::MU) if (checkcache_int_mu && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ7_mu;
1308 if (lep == StandardModel::ELECTRON) if (checkcache_int_el && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ7_el;
1309 wf=ROOT::Math::Functor1D(&(*this),&MVlnu::J7);
1310 ig.SetFunction(wf);
1311 J_res = ig.Integral(q2_min, q2_max);
1312 return J_res;
1313 break;
1314 case 11:
1315 if (lep == StandardModel::TAU) if (checkcache_int_tau && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ8_tau;
1316 if (lep == StandardModel::MU) if (checkcache_int_mu && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ8_mu;
1317 if (lep == StandardModel::ELECTRON) if (checkcache_int_el && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ8_el;
1318 wf=ROOT::Math::Functor1D(&(*this),&MVlnu::J8);
1319 ig.SetFunction(wf);
1320 J_res = ig.Integral(q2_min, q2_max);
1321 return J_res;
1322 break;
1323 case 12:
1324 if (lep == StandardModel::TAU) if (checkcache_int_tau && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ9_tau;
1325 if (lep == StandardModel::MU) if (checkcache_int_mu && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ9_mu;
1326 if (lep == StandardModel::ELECTRON) if (checkcache_int_el && (q2_min == q2min) && (q2_max == q2max)) return cached_intJ9_el;
1327 wf=ROOT::Math::Functor1D(&(*this),&MVlnu::J9);
1328 ig.SetFunction(wf);
1329 J_res = ig.Integral(q2_min, q2_max);
1330 return J_res;
1331 break;
1332 default:
1333 std::stringstream out;
1334 out << i;
1335 throw std::runtime_error("MVlnu::integrateJ: index " + out.str() + " not implemented");
1336 }
1337}
1338
1339double MVlnu::dGammadw(double q2)
1340{
1341 updateParameters();
1342
1343 return 3. / 4. * (2. * J1s(q2) + J1c(q2)) - 1. / 4. * (2. * J2s(q2) + J2c(q2));
1344}
1345
1346double MVlnu::getDeltaGammaDeltaw(double w_min, double w_max)
1347{
1348 updateParameters();
1349
1350 double q2_min = (2. * MM * MV)*(w0 - w_max); // min is Mlep*Mlep;
1351 double q2_max = (2. * MM * MV)*(w0 - w_min); // max is (MM-MV)*(MM-MV);
1352
1353 double intJ1s = integrateJ(1, q2_min, q2_max);
1354 double intJ1c = integrateJ(2, q2_min, q2_max);
1355 double intJ2s = integrateJ(3, q2_min, q2_max);
1356 double intJ2c = integrateJ(4, q2_min, q2_max);
1357
1358 return 3. / 4. * (2. * intJ1s + intJ1c) - 1. / 4. * (2. * intJ2s + intJ2c);
1359}
1360
1361double MVlnu::dGammadcldq2(double q2, double cl)
1362{
1363 updateParameters();
1364
1365 return 3. / 8. * ((J1s(q2) + 2. * J1c(q2)) + cl * (J6s(q2) + 2. * J6c(q2))+(2. * cl * cl - 1.)*(J2s(q2) + 2. * J2c(q2)));
1366}
1367
1368double MVlnu::dGammadcl(double cl)
1369{
1370 updateParameters();
1371
1372 double intJ1s = integrateJ(1, q2min, q2max);
1373 double intJ1c = integrateJ(2, q2min, q2max);
1374 double intJ2s = integrateJ(3, q2min, q2max);
1375 double intJ2c = integrateJ(4, q2min, q2max);
1376 double intJ6s = integrateJ(8, q2min, q2max);
1377 double intJ6c = integrateJ(9, q2min, q2max);
1378
1379 return 3. / 8. * ((intJ1s + 2. * intJ1c) + cl * (intJ6s + 2. * intJ6c)+(2. * cl * cl - 1.)*(intJ2s + 2. * intJ2c));
1380}
1381
1382double MVlnu::getDeltaGammaDeltacl(double cl_min, double cl_max)
1383{
1384 updateParameters();
1385
1386 double intJ1s = integrateJ(1, q2min, q2max);
1387 double intJ1c = integrateJ(2, q2min, q2max);
1388 double intJ2s = integrateJ(3, q2min, q2max);
1389 double intJ2c = integrateJ(4, q2min, q2max);
1390 double intJ6s = integrateJ(8, q2min, q2max);
1391 double intJ6c = integrateJ(9, q2min, q2max);
1392
1393 return 3. / 8. * ((cl_max - cl_min)*(intJ1c + 2. * intJ1s)+
1394 (cl_max * cl_max - cl_min * cl_min) / 2. * (intJ6c + 2. * intJ6s)+
1395 (2. / 3. * (cl_max * cl_max * cl_max - cl_min * cl_min * cl_min)-(cl_max - cl_min))*(intJ2c + 2. * intJ2s));
1396}
1397
1398double MVlnu::dGammadcVdq2(double q2, double cl)
1399{
1400 updateParameters();
1401
1402 return 3. / 8. * ((J1s(q2) + 2. * J1c(q2)) + cl * (J6s(q2) + 2. * J6c(q2))+(2. * cl * cl - 1.)*(J2s(q2) + 2. * J2c(q2)));
1403}
1404
1405double MVlnu::dGammadcV(double cV)
1406{
1407 updateParameters();
1408
1409 double intJ1s = integrateJ(1, q2min, q2max);
1410 double intJ1c = integrateJ(2, q2min, q2max);
1411 double intJ2s = integrateJ(3, q2min, q2max);
1412 double intJ2c = integrateJ(4, q2min, q2max);
1413
1414 return 3. / 8. * (cV * cV * (3. * intJ1c - intJ2c)+(1. - cV * cV)*(3. * intJ1s - intJ2s));
1415}
1416
1417double MVlnu::getDeltaGammaDeltacV(double cV_min, double cV_max)
1418{
1419 updateParameters();
1420
1421 double intJ1s = integrateJ(1, q2min, q2max);
1422 double intJ1c = integrateJ(2, q2min, q2max);
1423 double intJ2s = integrateJ(3, q2min, q2max);
1424 double intJ2c = integrateJ(4, q2min, q2max);
1425
1426 return 3. / 8. * ((cV_max * cV_max * cV_max - cV_min * cV_min * cV_min) / 3. * (3. * intJ1c - intJ2c)+
1427 ((cV_max - cV_min)-(cV_max * cV_max * cV_max - cV_min * cV_min * cV_min) / 3.)*(3. * intJ1s - intJ2s));
1428}
1429
1430double MVlnu::dGammadchidq2(double q2, double chi)
1431{
1432 updateParameters();
1433
1434 return (3. * J1c(q2) + 6. * J1s(q2) - J2c(q2) - 2. * J2s(q2)) / 8. / M_PI +
1435 cos(2. * chi) / 2. / M_PI * J3(q2) + sin(2. * chi) / 2. / M_PI * J9(q2);
1436}
1437
1438double MVlnu::dGammadchi(double chi)
1439{
1440 updateParameters();
1441
1442 double intJ1s = integrateJ(1, q2min, q2max);
1443 double intJ1c = integrateJ(2, q2min, q2max);
1444 double intJ2s = integrateJ(3, q2min, q2max);
1445 double intJ2c = integrateJ(4, q2min, q2max);
1446 double intJ3 = integrateJ(5, q2min, q2max);
1447 double intJ9 = integrateJ(12, q2min, q2max);
1448
1449 return ((3. * intJ1c + 6. * intJ1s - intJ2c - 2. * intJ2s) / 4. +
1450 cos(2. * chi) * intJ3 + sin(2. * chi) * intJ9) / 2. / M_PI;
1451}
1452
1453double MVlnu::getDeltaGammaDeltachi(double chi_min, double chi_max)
1454{
1455 updateParameters();
1456
1457 double intJ1s = integrateJ(1, q2min, q2max);
1458 double intJ1c = integrateJ(2, q2min, q2max);
1459 double intJ2s = integrateJ(3, q2min, q2max);
1460 double intJ2c = integrateJ(4, q2min, q2max);
1461 double intJ3 = integrateJ(5, q2min, q2max);
1462 double intJ9 = integrateJ(12, q2min, q2max);
1463
1464 return ((chi_max - chi_min)*(3. * intJ1c + 6. * intJ1s - intJ2c - 2. * intJ2s) / 4. +
1465 (sin(2. * chi_max) - sin(2. * chi_min)) / 2. * intJ3 -
1466 (cos(2. * chi_max) - cos(2. * chi_min)) / 2. * intJ9) / (2. * M_PI);
1467}
1468
1470{
1471 updateParameters();
1472
1473 double intJ1s = integrateJ(1, q2min, q2max);
1474 double intJ1c = integrateJ(2, q2min, q2max);
1475 double intJ2s = integrateJ(3, q2min, q2max);
1476 double intJ2c = integrateJ(4, q2min, q2max);
1477
1478 double DeltaJL = (3. * intJ1c - intJ2c) / 4.;
1479 double DeltaJ = 3. / 4. * (2. * intJ1s + intJ1c) - 1. / 4. * (2. * intJ2s + intJ2c);
1480 return DeltaJL/DeltaJ;
1481
1482}
1483
1485{
1486 updateParameters();
1487
1488 return ag0 * ag0 + ag1 * ag1 + ag2*ag2;
1489
1490}
1491
1493{
1494 updateParameters();
1495
1496 double aF10 = (MM - MV)*(phi_F1(0.) / phi_f(0.)) * af0;
1497 return af0 * af0 + af1 * af1 + af2 * af2 + aF10 * aF10 + aF11 * aF11 + aF12*aF12 + aF13*aF13;
1498}
1499
1501{
1502 updateParameters();
1503
1504 double z0 = (sqrt(w0 + 1.) - M_SQRT2) / (sqrt(w0 + 1.) + M_SQRT2);
1505 double PfacF2z0 = (z0 - zP1) / (1. - z0 * zP1)*(z0 - zP2) / (1. - z0 * zP2)*(z0 - zP3) / (1. - z0 * zP3);
1506 double phiF2z0 = phi_F2(z0);
1507 double aF20 = PfacF2z0 * phiF2z0 * 2. * F1_BGL(0.) / (MM * MM - MV * MV) - aF21 * z0 - aF22 * z0*z0 - aF23 * z0*z0*z0;
1508 return aF20 * aF20 + aF21 * aF21 + aF22*aF22 + aF23*aF23;
1509}
1510
1512{
1513 updateParameters();
1514
1515 return hA1(q2max);
1516}
1517
1518double MVlnu::get_hA1(double w)
1519{
1520 updateParameters();
1521 double q2 = (2. * MM * MV)*(w0 - w);
1522
1523 return hA1(q2);
1524}
1525
1526double MVlnu::get_hA2(double w)
1527{
1528 updateParameters();
1529 double q2 = (2. * MM * MV)*(w0 - w);
1530
1531 return (hA1(q2) * (1. + w) - RV * (A0(q2) * (1. + MV_o_MM) - A2(q2) * (MV_o_MM - w))) / (1. + MV_o_MM*MV_o_MM - 2.* MV_o_MM * w);
1532 // return (hA1(q2) * (1. + w) - A0(q2) * RV * (1. + MV_o_MM) - A2(q2) * RV * (w - MV_o_MM)) / (1. + MV_o_MM*MV_o_MM - 2.*MV_o_MM*w) ;
1533}
1534
1535double MVlnu::get_hA3(double w)
1536{
1537 updateParameters();
1538 double q2 = (2. * MM * MV)*(w0 - w);
1539
1540 return A2(q2) * RV - MV_o_MM * get_hA2(w) ;
1541}
1542
1543double MVlnu::get_hV(double w)
1544{
1545 updateParameters();
1546 double q2 = (2. * MM * MV)*(w0 - w);
1547
1548 return V(q2) * RV;
1549}
1550
1551double MVlnu::get_R1(double w)
1552{
1553 updateParameters();
1554 double q2 = (2. * MM * MV)*(w0 - w);
1555
1556 if (CLNflag) return R1(q2);
1557 else if (BGLflag) return V(q2) * RV / hA1(q2);
1558 return 0.;
1559}
1560
1561double MVlnu::get_R2(double w)
1562{
1563 updateParameters();
1564 double q2 = (2. * MM * MV)*(w0 - w);
1565
1566 if (CLNflag) return R2(q2);
1567 else if (BGLflag) return A2(q2) * RV / hA1(q2);
1568 return 0.;
1569}
1570
1571double MVlnu::get_R0(double w)
1572{
1573 updateParameters();
1574 double q2 = (2. * MM * MV)*(w0 - w);
1575
1576 if (CLNflag) return R0(q2);
1577 else if (BGLflag) return A0(q2) / A1(q2);
1578 return 0.;
1579}
1580
1581
1582/***************************************************************************
1583 * SM computation ... lep hel asymmetry, see 1203.2654, 1707.09509 *
1584 * ************************************************************************/
1585
1586double MVlnu::Hplus(double q2){
1587 double abs_p = lambda_half(MM*MM,MV*MV,q2);
1588 return (MM+MV)*A1(q2)-2.*MM/(MM+MV)*abs_p*V(q2);
1589}
1590
1591double MVlnu::Hminus(double q2){
1592 double abs_p = lambda_half(MM*MM,MV*MV,q2);
1593 return (MM+MV)*A1(q2)+2.*MM/(MM+MV)*abs_p*V(q2);
1594}
1595
1596double MVlnu::H0(double q2){
1597 double abs_p = lambda_half(MM*MM,MV*MV,q2);
1598 return ((MM*MM-MV*MV-q2)*(MM+MV)*A1(q2)-4.*MM*MM*abs_p*abs_p/(MM+MV)*A2(q2))/(2.*MV*sqrt(q2));
1599}
1600
1601double MVlnu::H0t(double q2){
1602 double abs_p = lambda_half(MM*MM,MV*MV,q2);
1603 return 2.*MM*abs_p*A0(q2)/sqrt(q2);
1604}
1605
1606double MVlnu::dGpdq2(double q2){
1607
1608 updateParameters();
1609
1610 double abs_p = lambda_half(MM*MM,MV*MV,q2);
1611 double lep_factor = (1.-Mlep*Mlep/q2)*(1.-Mlep*Mlep/q2)*Mlep*Mlep/(2.*q2);
1612 return 2./3.*amplsq_factor*abs_p*q2*lep_factor*(Hplus(q2)*Hplus(q2)+Hminus(q2)*Hminus(q2)+H0(q2)*H0(q2)
1613 +3.*H0t(q2)*H0t(q2));
1614}
1615
1616double MVlnu::dGmdq2(double q2){
1617
1618 updateParameters();
1619
1620 double abs_p = lambda_half(MM*MM,MV*MV,q2);
1621 double lep_factor = (1.-Mlep*Mlep/q2)*(1.-Mlep*Mlep/q2);
1622 return 2./3.*amplsq_factor*abs_p*q2*lep_factor*(Hplus(q2)*Hplus(q2)+Hminus(q2)*Hminus(q2)+H0(q2)*H0(q2));
1623}
1624
1625double MVlnu::integrateGpm(int i, double q2_min, double q2_max)
1626{
1627 switch (i) {
1628 case 1:
1629 wf=ROOT::Math::Functor1D(&(*this),&MVlnu::dGpdq2);
1630 ig.SetFunction(wf);
1631 J_res = ig.Integral(q2_min, q2_max);
1632 return J_res;
1633 break;
1634 case 2:
1635 wf=ROOT::Math::Functor1D(&(*this),&MVlnu::dGmdq2);
1636 ig.SetFunction(wf);
1637 J_res = ig.Integral(q2_min, q2_max);
1638 return J_res;
1639 break;
1640 default:
1641 std::stringstream out;
1642 out << i;
1643 throw std::runtime_error("MVlnu::integrateGpm: index " + out.str() + " not implemented");
1644 }
1645}
1646
1648{
1649 updateParameters();
1650
1651 double DeltaGammaPlus = integrateGpm(1,q2min,q2max);
1652 double DeltaGammaMinus = integrateGpm(2,q2min,q2max);
1653
1654 return (DeltaGammaPlus-DeltaGammaMinus)/(DeltaGammaPlus+DeltaGammaMinus);
1655
1656}
@ LO
Definition: OrderScheme.h:34
@ FULLNLO
Definition: OrderScheme.h:38
bool getFlagBGL() const
Definition: Flavour.h:344
bool getFlagCLN() const
Definition: Flavour.h:340
void setUpdateFlag(QCD::meson meson_i, QCD::meson meson_j, QCD::lepton lep_i, bool updated_i) const
sets the update flag for the initial and final state dependent object for .
Definition: Flavour.cpp:309
bool getUpdateFlag(QCD::meson meson_i, QCD::meson meson_j, QCD::lepton lep_i) const
gets the update flag for the initial and final state dependent object for .
Definition: Flavour.cpp:334
bool getFlagDM() const
Definition: Flavour.h:348
gslpp::vector< gslpp::complex > ** ComputeCoeffdiujlknu(int i, int j, int k, double mu) const
Computes the Wilson coefficient for the Hamiltonian transitions in the JMS basis ordered as CnueduVL...
Definition: Flavour.cpp:152
double get_hA1w1()
return A1 form factor at
Definition: MVlnu.cpp:1511
double get_R0(double w)
return at
Definition: MVlnu.cpp:1571
double MV_o_MM
Definition: MVlnu.h:169
double get_unitarity_V_BGL()
Vector unitarity constraint for BGL parameters.
Definition: MVlnu.cpp:1484
double sqrtMV_o_MM
Definition: MVlnu.h:170
double w0
Definition: MVlnu.h:177
double getDeltaGammaDeltacl(double cl_min, double cl_max)
The integral of from to .
Definition: MVlnu.cpp:1382
double get_unitarity_P_BGL()
Pseudoscalar unitarity constraint for BGL parameters.
Definition: MVlnu.cpp:1500
double getPlep()
Binned lepton helicity asymmetry .
Definition: MVlnu.cpp:1647
double mu_b
Definition: MVlnu.h:179
bool btocNPpmflag
Definition: MVlnu.h:166
double get_hA1(double w)
return at
Definition: MVlnu.cpp:1518
std::vector< std::string > initializeMVlnuParameters()
Definition: MVlnu.cpp:167
double get_hV(double w)
return at
Definition: MVlnu.cpp:1543
QCD::meson meson
Definition: MVlnu.h:160
bool BGLflag
Definition: MVlnu.h:164
double get_hA2(double w)
return at
Definition: MVlnu.cpp:1526
const StandardModel & mySM
Definition: MVlnu.h:158
double get_R2(double w)
return at
Definition: MVlnu.cpp:1561
bool NPanalysis
Definition: MVlnu.h:167
double RV
Definition: MVlnu.h:178
double MM
Definition: MVlnu.h:175
double getDeltaGammaDeltaw(double w_min, double w_max)
The integral of from to .
Definition: MVlnu.cpp:1346
double Mb
Definition: MVlnu.h:180
double get_R1(double w)
return at
Definition: MVlnu.cpp:1551
double Mc
Definition: MVlnu.h:181
double getFL()
Binned D* polarization fraction .
Definition: MVlnu.cpp:1469
double Mnu
Definition: MVlnu.h:174
double width
Definition: MVlnu.h:182
QCD::meson vectorM
Definition: MVlnu.h:161
bool DMflag
Definition: MVlnu.h:165
virtual ~MVlnu()
Destructor.
Definition: MVlnu.cpp:164
bool CLNflag
Definition: MVlnu.h:163
QCD::lepton lep
Definition: MVlnu.h:159
double get_hA3(double w)
return at
Definition: MVlnu.cpp:1535
double MV
Definition: MVlnu.h:176
MVlnu(const StandardModel &SM_i, QCD::meson meson_i, QCD::meson vector_i, QCD::lepton lep_i)
Constructor.
Definition: MVlnu.cpp:19
double get_unitarity_A_BGL()
Axial unitarity constraint for BGL parameters.
Definition: MVlnu.cpp:1492
double getDeltaGammaDeltachi(double chi_min, double chi_max)
The integral of from to .
Definition: MVlnu.cpp:1453
double Mlep
Definition: MVlnu.h:173
double getDeltaGammaDeltacV(double cV_min, double cV_max)
The integral of from to .
Definition: MVlnu.cpp:1417
std::vector< std::string > mvlnuParameters
Definition: MVlnu.h:162
double computeWidth() const
A method to compute the width of the meson from its lifetime.
Definition: Meson.cpp:521
std::string getModelName() const
A method to fetch the name of the model.
Definition: Model.h:59
const double & getMass() const
A get method to access the particle mass.
Definition: Particle.h:61
meson
An enum type for mesons.
Definition: QCD.h:336
@ D_star_P
Definition: QCD.h:352
const double getOptionalParameter(std::string name) const
A method to get parameters that are specific to only one set of observables.
Definition: QCD.h:450
const Meson & getMesons(const QCD::meson m) const
A get method to access a meson as an object of the type Meson.
Definition: QCD.h:526
@ BOTTOM
Definition: QCD.h:329
@ CHARM
Definition: QCD.h:326
lepton
An enum type for leptons.
Definition: QCD.h:310
@ MU
Definition: QCD.h:314
@ ELECTRON
Definition: QCD.h:312
@ TAU
Definition: QCD.h:316
const Particle & getQuarks(const QCD::quark q) const
A get method to access a quark as an object of the type Particle.
Definition: QCD.h:536
void initializeMeson(QCD::meson meson_i) const
A method to initialize a meson.
Definition: QCD.cpp:280
A model class for the Standard Model.
virtual const double getCCC1() const
A virtual implementation for the RealWeakEFTCC class.
const Particle & getLeptons(const QCD::lepton p) const
A get method to retrieve the member object of a lepton.
virtual const double getCCC3() const
A virtual implementation for the RealWeakEFTCC class.
const double getMz() const
A get method to access the mass of the boson .
const double Ale(double mu, orders order, bool Nf_thr=true) const
The running electromagnetic coupling in the scheme.
virtual const double getCCC4() const
A virtual implementation for the RealWeakEFTCC class.
const Flavour & getFlavour() const
virtual const double getCCC2() const
A virtual implementation for the RealWeakEFTCC class.
virtual const double getCCC5() const
A virtual implementation for the RealWeakEFTCC class.