|
JHUGen MELA
JHUGen v7.5.6, MELA v2.4.2
Matrix element calculations as used in JHUGen.
|
This is the main Mela object!
More...
#include <Mela.h>
|
| Mela (double LHCsqrts_=13., double mh_=125., TVar::VerbosityLevel verbosity_=TVar::ERROR) |
| the MELA constructor More...
|
|
| Mela (const Mela &other) |
| Copy constructor for MELA. More...
|
|
| ~Mela () |
| MELA destructor. More...
|
|
void | build (double mh_) |
| This is the actual building of the tool that occurs in each instance of the Mela::Mela constructor. More...
|
|
void | setProcess (TVar::Process myModel, TVar::MatrixElement myME, TVar::Production myProduction) |
| Sets the process, matrix element, and production that MELA is to use for this event. Calls ZZMatrixElement::set_Process, which calls TEvtProb::SetProcess. More...
|
|
void | setVerbosity (TVar::VerbosityLevel verbosity_=TVar::ERROR) |
| Sets the verbosity for MELA outside of the initial constructor. More...
|
|
TVar::VerbosityLevel | getVerbosity () |
| Gets the current verbosity level for MELA. More...
|
|
void | setMelaLeptonInterference (TVar::LeptonInterference myLepInterf=TVar::DefaultLeptonInterf) |
| Sets the MELA Lepton Interference. More...
|
|
void | setRemoveLeptonMasses (bool MasslessLeptonSwitch=true) |
| either permits or forbids massive leptons. More...
|
|
void | setRemoveJetMasses (bool MasslessLeptonSwitch=true) |
| either permits or forbids massive jets. More...
|
|
void | setMelaPrimaryHiggsMass (double myHiggsMass) |
| Sets the mass of the "primary" higgs. More...
|
|
void | setMelaHiggsMass (double myHiggsMass, int index=0) |
| Sets the mass of your chosen Higgs. More...
|
|
void | setMelaHiggsWidth (double myHiggsWidth=-1, int index=0) |
| Sets the width of your chosen Higgs. More...
|
|
void | setMelaHiggsMassWidth (double myHiggsMass, double myHiggsWidth, int index) |
| a combination of setMelaHiggsMass and setMelaHiggsWidth. More...
|
|
void | setRenFacScaleMode (TVar::EventScaleScheme renormalizationSch, TVar::EventScaleScheme factorizationSch, double ren_sf, double fac_sf) |
| Sets the renormalization and the factorization schemes. More...
|
|
void | setCandidateDecayMode (TVar::CandidateDecayMode mode) |
| Sets the decay mode for your event. More...
|
|
void | setCurrentCandidateFromIndex (unsigned int icand) |
| Switches the candidate that you are working on to another candidate based off of an index. More...
|
|
void | setCurrentCandidate (MELACandidate *cand) |
| Switches the candidate that you are working on to another candidate object specified. More...
|
|
void | setInputEvent (SimpleParticleCollection_t *pDaughters, SimpleParticleCollection_t *pAssociated=0, SimpleParticleCollection_t *pMothers=0, bool isGen=false) |
| Sets the input event for MELA. MELA cannot run without this. More...
|
|
void | resetInputEvent () |
| Resets the event in preparation for the next iteration of the event loop. More...
|
|
void | setTempCandidate (SimpleParticleCollection_t *pDaughters, SimpleParticleCollection_t *pAssociated=0, SimpleParticleCollection_t *pMothers=0, bool isGen=false) |
| Sets a temporary MELA candidate, by setting melaCand in Xcal2 to a temporary candidate without pushing this candidate to the candList of Xcal2. More...
|
|
void | appendTopCandidate (SimpleParticleCollection_t *TopDaughters) |
| Adds a top quark as a MELA candidate. More...
|
|
void | resetMass (double inmass, int ipart) |
| Resets the mass for a particle that is an electroweak parameter according to its id. More...
|
|
void | resetWidth (double inwidth, int ipart) |
| Resets the width for a particle that is an electroweak parameter according to its id. More...
|
|
void | resetYukawaMass (double inmass, int ipart) |
| Resets the Yukawa coupling (or "Yukawa mass") for a specific particle. Does not change its "intrinsic" mass. More...
|
|
void | resetQuarkMasses () |
| Resets the masses of each quark to their original values. More...
|
|
void | resetMCFM_EWKParameters (double ext_Gf, double ext_aemmz, double ext_mW, double ext_mZ, double ext_xW, int ext_ewscheme=3) |
| Resets the electroweak parameters back to their defaults. More...
|
|
double | getPrimaryMass (int ipart) |
| A function to get the current primary EW/QCD parameters from MELA. More...
|
|
double | getPrimaryWidth (int ipart) |
| A function to get the current primary EW/QCD parameters from MELA. More...
|
|
double | getHiggsWidthAtPoleMass (double mass) |
| Returns the width of the Higgs at a given pole mass as a calculation. More...
|
|
void | SetMadgraphCKMElements (double ckmlambda=0.2265, double ckma=0.79, double ckmrho=0.141, double ckmeta=0.357, bool force_refresh=false) |
|
std::complex< double > | GetMadgraphCKMElement (int iquark, int jquark) |
|
MelaIO * | getIORecord () |
| Returns the MELAIO object, and by consequence, the entire parton-by-parton matrix element record. More...
|
|
MELACandidate * | getCurrentCandidate () |
| Gets the current MELA top-level (input) candList object. More...
|
|
int | getCurrentCandidateIndex () |
| Returns the index of the current MELA candidate - returns -1 if there is no candidate to be found. More...
|
|
int | getNCandidates () |
| Returns the size of the candidate list TEvtProb::candList. More...
|
|
std::vector< MELATopCandidate_t * > * | getTopCandidateCollection () |
| Same as getNCandidates but specifically for Top Quark Candidates. More...
|
|
void | getConstant (float &prob) |
| This function returns a multiplicative constant to the matrix element calculation. More...
|
|
void | getPAux (float &prob) |
|
RooSpin::modelMeasurables | getMeasurablesRRV () |
| Returns a RooSpin::modelMeasureables object containing many observable quantities. More...
|
|
const TVar::event_scales_type & | getRenFacScaleMode () const |
| Returns the current RenFac Scale Mode for MELA. More...
|
|
void | computeDecayAngles (float &qH, float &m1, float &m2, float &costheta1, float &costheta2, float &Phi, float &costhetastar, float &Phi1) |
| computes the decay angles for gg -> H -> ZZ as defined in Figure 1 of arXiv:1001.3396 More...
|
|
void | computeVBFAngles (float &Q2V1, float &Q2V2, float &costheta1, float &costheta2, float &Phi, float &costhetastar, float &Phi1) |
| computes the production angles for Vector Boson Fusion (VBF) -> H -> ZZ as defined in Figure 1 of arXiv:1208.4018 More...
|
|
void | computeVBFAngles_ComplexBoost (float &Q2V1, float &Q2V2, float &costheta1_real, float &costheta1_imag, float &costheta2_real, float &costheta2_imag, float &Phi, float &costhetastar, float &Phi1) |
|
void | computeVHAngles (float &mVstar, float &mV, float &costheta1, float &costheta2, float &Phi, float &costhetastar, float &Phi1) |
| computes the production for Vector Boson Fusion (VBF) -> H -> ZZ as defined in Figure 1 of arXiv:1208.4018 More...
|
|
void | computeTTHAngles (int topDecay, float &mT1, float &mW1, float &mT2, float &mW2, float &costheta1, float &costheta2, float &Phi, float &costhetastar, float &Phi1, float &hbb, float &hWW, float &Phibb, float &Phi1bb, float &hWplusf, float &PhiWplusf, float &hWminusf, float &PhiWminusf) |
| computes the angles for a ttH system More...
|
|
void | computeP_selfDspin0 (double selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2], float &prob, bool useConstant=true) |
| This is a function that calls Mela::computeP with preset quark, gluon, and vector boson couplings for spin 0. More...
|
|
void | computeP_selfDspin1 (double selfDZqqcoupl_input[SIZE_ZQQ][2], double selfDZvvcoupl_input[SIZE_ZVV][2], float &prob, bool useConstant=true) |
| This is a function that calls Mela::computeP with preset quark, gluon, and vector boson couplings for spin 1. More...
|
|
void | computeP_selfDspin1 (double selfDZvvcoupl_input[SIZE_ZVV][2], float &prob, bool useConstant=true) |
| This is a function that calls Mela::computeP with preset quark, gluon, and vector boson couplings for spin 1, but sets default values for Zqq couplings. More...
|
|
void | computeP_selfDspin2 (double selfDGggcoupl_input[SIZE_GGG][2], double selfDGqqcoupl_input[SIZE_GQQ][2], double selfDGvvcoupl_input[SIZE_GVV][2], float &prob, bool useConstant=true) |
| This is a function that calls Mela::computeP with preset quark, gluon, and vector boson couplings for spin 2. More...
|
|
void | computeP_selfDspin2 (double selfDGggcoupl_input[SIZE_GGG][2], double selfDGvvcoupl_input[SIZE_GVV][2], float &prob, bool useConstant=true) |
| This is a function that calls Mela::computeP with preset quark, gluon, and vector boson couplings for spin 1, but sets default values for Zqq couplings. More...
|
|
void | computeP (float &prob, bool useConstant=true) |
| Computes the probability for the probabilities on the decay side of things using the constituent daughter 4 vectors along with any jets that could exist. More...
|
|
void | computeD_CP (TVar::MatrixElement myME, TVar::Process myType, float &prob) |
| computes the value of D_CP More...
|
|
void | computeProdDecP (double selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2], double selfDHwwcoupl_input[nSupportedHiggses][SIZE_HVV][2], double selfDaTQGCcoupl_input[SIZE_ATQGC][2], double selfDAZffcoupl_input[SIZE_AZff][2], float &prob, bool useConstant=true) |
| computes the combined production and decay probability while taking in coupling arrays More...
|
|
void | computeProdDecP (float &prob, bool useConstant=true) |
| computes the combined production and decay probability, with couplings set beforehand More...
|
|
void | computeProdP (double selfDHggcoupl_input[SIZE_HGG][2], double selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2], double selfDHwwcoupl_input[nSupportedHiggses][SIZE_HVV][2], float &prob, bool useConstant=true) |
| computes Production side probabilities while taking in coupling arrays More...
|
|
void | computeProdP (float &prob, bool useConstant=true) |
| computes Production side probabilities using couplings set beforehand More...
|
|
void | computeProdP_VH (double selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2], float &prob, bool includeHiggsDecay=false, bool useConstant=true) |
|
void | computeProdP_VH (float &prob, bool includeHiggsDecay=false, bool useConstant=true) |
|
void | computeProdP_ttH (float &prob, int topProcess=2, int topDecay=0, bool useConstant=true) |
|
void | compute4FermionWeight (float &w) |
|
void | getXPropagator (TVar::ResonancePropagatorScheme scheme, float &prop) |
|
void | computePM4l (TVar::SuperMelaSyst syst, float &prob) |
|
void | computeDijetConvBW (float &prob, bool useTrueBW=false) |
|
void | computeD_gg (TVar::MatrixElement myME, TVar::Process myType, float &prob) |
|
std::vector< TLorentzVector > | calculate4Momentum (double Mx, double M1, double M2, double theta, double theta1, double theta2, double Phi1, double Phi) |
|
|
void | printLogo () const |
|
void | setSpinZeroCouplings () |
|
void | setSpinOneCouplings () |
|
void | setSpinTwoCouplings () |
|
void | setATQGCCouplings () |
|
void | setAZffCouplings () |
|
bool | configureAnalyticalPDFs () |
|
void | reset_SelfDCouplings () |
|
void | reset_PAux () |
|
void | reset_CandRef () |
|
void | constructDggr (float bkg_VAMCFM_noscale, float ggzz_VAMCFM_noscale, float ggHZZ_prob_pure_noscale, float ggHZZ_prob_int_noscale, float widthScale, float &myDggr) |
|
void | getPConstantHandles () |
|
void | deletePConstantHandles () |
|
MelaPConstant * | getPConstantHandle (TVar::MatrixElement me_, TVar::Production prod_, TVar::Process proc_, TString relpath, TString spname, const bool useSqrts=false) |
|
void | deletePConstantHandle (MelaPConstant *&handle) |
|
void | computeConstant (float &prob) |
|
void | setConstant () |
|
float | getConstant_JHUGenUndecayed () |
|
float | getConstant_4l () |
|
float | getConstant_2l2q () |
|
float | getConstant_4q () |
|
float | getConstant_FourFermionDecay (const int &decid) |
|
This is the main Mela object!
Definition at line 54 of file Mela.h.
◆ Mela() [1/2]
the MELA constructor
- Parameters
-
[in] | LHCsqrts_ | The luminosity of your collider in TeV. Default is LHC at 13 TeV. |
[in] | mh_ | The mass of the Higgs in GeV. Default is 125 GeV |
[in] | verbosity_ | The verbosity of MELA that you desire, as defined in TVar::ERROR |
Definition at line 53 of file Mela.cc.
◆ Mela() [2/2]
Mela::Mela |
( |
const Mela & |
other | ) |
|
Copy constructor for MELA.
- Parameters
-
[in] | other | another MELA instance |
Definition at line 70 of file Mela.cc.
◆ ~Mela()
MELA destructor.
Definition at line 81 of file Mela.cc.
◆ appendTopCandidate()
◆ build()
void Mela::build |
( |
double |
mh_ | ) |
|
This is the actual building of the tool that occurs in each instance of the Mela::Mela constructor.
- Parameters
-
mh_ | This is the mass of the Higgs in GeV |
Definition at line 139 of file Mela.cc.
145 const double maxSqrts = 8.;
153 const string mcfmWarning = MELAPKGPATH +
"data/ffwarn.dat"; symlink(mcfmWarning.c_str(),
"ffwarn.dat");
154 const string mcfm_brsm_o = MELAPKGPATH +
"data/br.sm1"; symlink(mcfm_brsm_o.c_str(),
"br.sm1");
155 const string mcfm_brsm_t = MELAPKGPATH +
"data/br.sm2"; symlink(mcfm_brsm_t.c_str(),
"br.sm2");
156 const string mcfmInput1 = MELAPKGPATH +
"data/input.DAT"; symlink(mcfmInput1.c_str(),
"input.DAT");
157 const string mcfmInput2 = MELAPKGPATH +
"data/process.DAT"; symlink(mcfmInput2.c_str(),
"process.DAT");
159 mkdir(
"Pdfdata", S_IRWXU | S_IRWXG | S_IROTH | S_IXOTH);
160 const string mcfmInput3 = MELAPKGPATH +
"data/Pdfdata/cteq6l1.tbl"; symlink(mcfmInput3.c_str(),
"Pdfdata/cteq6l1.tbl");
161 const string mcfmInput4 = MELAPKGPATH +
"data/Pdfdata/cteq6l.tbl"; symlink(mcfmInput4.c_str(),
"Pdfdata/cteq6l.tbl");
165 mzz_rrv =
new RooRealVar(
"mzz",
"m_{ZZ}", mh_, 0., 1000.);
166 z1mass_rrv =
new RooRealVar(
"z1mass",
"m_{Z1}", 0., 160.);
167 z2mass_rrv =
new RooRealVar(
"z2mass",
"m_{Z2}", 0., 200.);
168 costhetastar_rrv =
new RooRealVar(
"costhetastar",
"cos#theta^{*}", -1., 1.);
169 costheta1_rrv =
new RooRealVar(
"costheta1",
"cos#theta_{1}", -1., 1.);
170 costheta2_rrv =
new RooRealVar(
"costheta2",
"cos#theta_{2}", -1., 1.);
171 phi_rrv=
new RooRealVar(
"phi",
"#Phi", -TMath::Pi(), TMath::Pi());
172 phi1_rrv=
new RooRealVar(
"phi1",
"#Phi_{1}", -TMath::Pi(), TMath::Pi());
173 Y_rrv =
new RooRealVar(
"Yzz",
"#Y_{ZZ}", 0, -4, 4);
174 upFrac_rrv =
new RooRealVar(
"upFrac",
"fraction up-quarks", .5, 0., 1.);
191 qqZZmodel =
new RooqqZZ_JHU_ZgammaZZ_fast(
"qqZZmodel",
"qqZZmodel", *
z1mass_rrv, *
z2mass_rrv, *
costheta1_rrv, *
costheta2_rrv, *
phi_rrv, *
costhetastar_rrv, *
phi1_rrv, *
mzz_rrv, *
upFrac_rrv);
194 const string path_HiggsWidthFile = MELAPKGPATH +
"data/HiggsTotalWidth_YR3.txt";
196 const string path_nnpdf = MELAPKGPATH +
"data/Pdfdata/NNPDF30_lo_as_0130.LHgrid";
197 char path_nnpdf_c[] =
"Pdfdata/NNPDF30_lo_as_0130.LHgrid";
200 symlink(path_nnpdf.c_str(), path_nnpdf_c);
217 RooMsgService::instance().getStream(1).removeTopic(NumIntegration);
218 RooMsgService::instance().setStreamStatus(1, kFALSE);
219 RooMsgService::instance().setStreamStatus(0, kFALSE);
223 if (superMELA_LHCsqrts > maxSqrts) superMELA_LHCsqrts = maxSqrts;
226 sprintf(cardpath,
"data/CombinationInputs/SM_inputs_%dTeV/inputs_4mu.txt", superMELA_LHCsqrts);
227 string cardfile = MELAPKGPATH + cardpath;
235 if (superDijetSqrts<13.) superDijetSqrts=13.;
◆ calculate4Momentum()
std::vector< TLorentzVector > Mela::calculate4Momentum |
( |
double |
Mx, |
|
|
double |
M1, |
|
|
double |
M2, |
|
|
double |
theta, |
|
|
double |
theta1, |
|
|
double |
theta2, |
|
|
double |
Phi1, |
|
|
double |
Phi |
|
) |
| |
◆ cleanLinkedFiles()
void Mela::cleanLinkedFiles |
( |
| ) |
|
|
static |
Definition at line 126 of file Mela.cc.
127 remove(
"ffwarn.dat");
131 remove(
"process.DAT");
132 remove(
"Pdfdata/cteq6l1.tbl");
133 remove(
"Pdfdata/cteq6l.tbl");
134 remove(
"Pdfdata/NNPDF30_lo_as_0130.LHgrid");
◆ compute4FermionWeight()
void Mela::compute4FermionWeight |
( |
float & |
w | ) |
|
Definition at line 1846 of file Mela.cc.
1851 bool hasFailed=
false;
1852 int id_original[2][2];
1853 for (
int iv=0; iv<2; iv++){
1857 for (
int ivd=0; ivd<2; ivd++) id_original[iv][ivd]=Vi->
getDaughter(ivd)->
id;
1884 float dXsec_HZZ_JHU, dXsec_HZZ_JHU_interf;
1889 computeP(dXsec_HZZ_JHU_interf,
false);
1892 w = dXsec_HZZ_JHU_interf / dXsec_HZZ_JHU;
◆ computeConstant()
void Mela::computeConstant |
( |
float & |
prob | ) |
|
|
protected |
◆ computeD_CP()
computes the value of D_CP
- See also
- calls Mela::computeP_selfDspin0
- Parameters
-
[in] | myME | the matrix element you want to use, set as a TVar::MatrixElement |
[in] | myType | the process you want to use, set as a TVar::Process |
[out] | prob | This is the value of the output probability - edited by reference |
Definition at line 1378 of file Mela.cc.
1390 coupl_mix[0][0][0] =1.;
1391 coupl_mix[0][3][0] =2.521;
1392 coupl_1[0][0][0] =1.;
1393 coupl_2[0][3][0] =2.521;
1396 coupl_mix[0][0][0] =1.;
1397 coupl_mix[0][3][1] =2.521;
1398 coupl_1[0][0][0] =1.;
1399 coupl_2[0][3][1] =2.521;
1402 coupl_mix[0][0][0] =1.;
1403 coupl_mix[0][1][0] = 1.638;
1404 coupl_1[0][0][0] =1.;
1405 coupl_2[0][1][0] = 1.638;
1408 coupl_mix[0][0][0] =1.;
1409 coupl_mix[0][1][1] = 1.638;
1410 coupl_1[0][0][0] =1.;
1411 coupl_2[0][1][1] = 1.638;
1414 coupl_mix[0][0][0] =1.;
1415 coupl_mix[0][11][0] = 12046.01;
1416 coupl_1[0][0][0] =1.;
1417 coupl_2[0][11][0] = 12046.01;
1420 coupl_mix[0][0][0] =1.;
1421 coupl_mix[0][4][0] = 0.0688;
1422 coupl_1[0][0][0] =1.;
1423 coupl_2[0][4][0] = 0.0688;
1426 coupl_mix[0][0][0] =1.;
1427 coupl_mix[0][7][0] = -0.0898;
1428 coupl_1[0][0][0] =1.;
1429 coupl_2[0][7][0] = -0.0898;
1432 coupl_mix[0][0][0] =1.;
1433 coupl_mix[0][6][0] = 0.0855;
1434 coupl_1[0][0][0] =1.;
1435 coupl_2[0][6][0] = 0.0855;
1438 coupl_mix[0][0][0] =1.;
1439 coupl_mix[0][9][0] = -0.0907;
1440 coupl_1[0][0][0] =1.;
1441 coupl_2[0][9][0] = -0.0907;
1444 coupl_mix[0][0][0] =1.;
1445 coupl_mix[0][30][0] = -7591.914;
1446 coupl_1[0][0][0] =1.;
1447 coupl_2[0][30][0] = -7591.914;
1450 coupl_mix[0][0][0] =1.;
1451 coupl_mix[0][30][1] = -7591.914;
1452 coupl_1[0][0][0] =1.;
1453 coupl_2[0][30][1] = -7591.914;
1456 MELAout <<
"Error: Not supported!"<<endl;
1464 prob = pMix- p1- p2;
◆ computeD_gg()
Definition at line 1973 of file Mela.cc.
1980 MELAout <<
"Only support MCFM and D_gg10"<<endl;
1986 float bkg_VAMCFM, ggzz_VAMCFM_noscale, ggHZZ_prob_pure_noscale, ggHZZ_prob_int_noscale, bkgHZZ_prob_noscale;
1992 bkgHZZ_prob_noscale /= ggScale;
1993 ggHZZ_prob_pure_noscale /= ggScale;
1994 ggzz_VAMCFM_noscale /= ggScale;
1996 ggHZZ_prob_int_noscale = bkgHZZ_prob_noscale - ggHZZ_prob_pure_noscale - ggzz_VAMCFM_noscale;
2000 constructDggr(bkg_VAMCFM, ggzz_VAMCFM_noscale, ggHZZ_prob_pure_noscale, ggHZZ_prob_int_noscale, 10., prob);
◆ computeDecayAngles()
void Mela::computeDecayAngles |
( |
float & |
qH, |
|
|
float & |
m1, |
|
|
float & |
m2, |
|
|
float & |
costheta1, |
|
|
float & |
costheta2, |
|
|
float & |
Phi, |
|
|
float & |
costhetastar, |
|
|
float & |
Phi1 |
|
) |
| |
computes the decay angles for gg -> H -> ZZ as defined in Figure 1 of arXiv:1001.3396
- See also
- Calls TUtil::computeAngles to fully calculate angles
-
Use TUtil::computeAngles if you would like to calculate the angles through inputting your own 4-vectors and quantities rather than through the MELA event loop
- Attention
- This function requires there to be an input event already defined
- Warning
- This function will edit all values inplace. Every value is technically an output.
-
The selection for m1 and m2 are based upon which of the leptons is closer to the z mass. $m_1 \leq m_2$.
- Parameters
-
[out] | qH | The mass of the "Higgs" as reconstructed by the constituent 4 leptons |
[out] | m1 | The mass of the first decay particle as reconstructed by 2 of the leptons |
[out] | m2 | The mass of the first decay particle as reconstructed by 2 of the leptons |
[out] | costheta1 | $\cos(\theta_1)$ |
[out] | costheta2 | $\cos(\theta_2)$ |
[out] | Phi | $\phi$ |
[out] | costhetastar | $\cos(\theta*)$ |
[out] | Phi1 | $\phi_1$ |
Definition at line 598 of file Mela.cc.
611 qH=0; m1=0; m2=0; costheta1=0; costheta2=0; Phi=0; costhetastar=0; Phi1=0;
615 TLorentzVector
const nullVector(0, 0, 0, 0);
618 simple_event_record mela_event;
619 mela_event.AssociationCode=partIncCode;
623 if (daughters.size()<2 || daughters.size()>4 || mela_event.intermediateVid.size()!=2){
624 if (
myVerbosity_>=
TVar::ERROR)
MELAerr <<
"Mela::computeDecayAngles: Number of daughters " << daughters.size() <<
" or number of intermediate Vs " << mela_event.intermediateVid <<
" not supported!" << endl;
629 if (daughters.size()%2==1){
for (
unsigned int ipar=daughters.size(); ipar<4; ipar++) daughters.push_back(
SimpleParticle_t(-9000, nullVector)); }
630 else if (daughters.size()==2){
634 qH = (daughters.at(0).second+daughters.at(1).second+daughters.at(2).second+daughters.at(3).second).M();
635 m1 = (daughters.at(0).second+daughters.at(1).second).M();
636 m2 = (daughters.at(2).second+daughters.at(3).second).M();
639 costhetastar, costheta1, costheta2, Phi, Phi1,
640 daughters.at(0).second, daughters.at(0).first,
641 daughters.at(1).second, daughters.at(1).first,
642 daughters.at(2).second, daughters.at(2).first,
643 daughters.at(3).second, daughters.at(3).first
647 if (!std::isfinite(costhetastar)) costhetastar=0;
648 if (!std::isfinite(costheta1)) costheta1=0;
649 if (!std::isfinite(costheta2)) costheta2=0;
650 if (!std::isfinite(Phi)) Phi=0;
651 if (!std::isfinite(Phi1)) Phi1=0;
◆ computeDijetConvBW()
void Mela::computeDijetConvBW |
( |
float & |
prob, |
|
|
bool |
useTrueBW = false |
|
) |
| |
◆ computeP()
void Mela::computeP |
( |
float & |
prob, |
|
|
bool |
useConstant = true |
|
) |
| |
Computes the probability for the probabilities on the decay side of things using the constituent daughter 4 vectors along with any jets that could exist.
- See also
- Mela::computeDecayAngles
- Attention
- All JHUGen matrix elements, productions, and processes can be used here
-
The following production modes are supported using the MCFM Matrix Element (as named in TVar::Production) TVar::ZZGG TVar::ZZQQB TVar::ZZQQB_STU TVar::ZZQQB_S TVar::ZZQQB_TU TVar::ZZINDEPENDENT TVar::JJQCD && process == TVar::bkgZJets
- Parameters
-
[out] | prob | This is the value of the output probability - edited by reference |
[in] | useConstant | This turns on the calculation of a corrective constant to different probabilities through Mela::getConstant. If you would like the "pure" MELA calculation to be run, set useConstant to false. By default true. |
Definition at line 1223 of file Mela.cc.
1232 TLorentzVector
const nullVector(0, 0, 0, 0);
1233 float mZZ=0, mZ1=0, mZ2=0, costheta1=0, costheta2=0, Phi=0, costhetastar=0, Phi1=0;
1238 costheta1, costheta2, Phi,
1249 Y_rrv->setConstant(
true);
1252 if (computeAnaMELA){
1255 prob = integral->getVal();
1258 else prob =
pdf->getVal();
1262 Y_rrv->setConstant(
false);
1275 costheta1, costheta2, Phi,
1280 MELAout <<
"Mela::computeP: Condition (myME_ == TVar::MCFM && myProduction_ == TVar::ZZINDEPENDENT && myModel_ == TVar::bkgZZ/WW/ZGamma/ZJJ/GammaGamma)." << endl;
1281 vector<TLorentzVector> pDauVec =
calculate4Momentum(mZZ, mZ1, mZ1, acos(costhetastar), acos(costheta1), acos(costheta2), Phi1, Phi);
1283 <<
"\tOriginal mZZ=" << mZZ <<
" "
1284 <<
"m1=" << mZ1 <<
" "
1285 <<
"m2=" << mZ2 <<
" "
1286 <<
"h1=" << costheta1 <<
" "
1287 <<
"h2=" << costheta2 <<
" "
1288 <<
"Phi=" << Phi <<
" "
1289 <<
"hs=" << costhetastar <<
" "
1290 <<
"Phi1=" << Phi1 << endl;
1291 MELAout <<
"\tfor daughters:" << endl;
1292 for (
int iv=0; iv<2; iv++){
1296 <<
"x=" << pDauVec.at(2*iv+idau).X() <<
" "
1297 <<
"y=" << pDauVec.at(2*iv+idau).Y() <<
" "
1298 <<
"z=" << pDauVec.at(2*iv+idau).Z() <<
" "
1299 <<
"t=" << pDauVec.at(2*iv+idau).T() << endl;
1306 int gridsize_hs = 5;
1309 double hs_step = (hs_max - hs_min) /
double(gridsize_hs);
1311 int gridsize_phi1 = 5;
1312 double phi1_min = 0;
1313 double phi1_max = TMath::Pi();
1314 double phi1_step = (phi1_max - phi1_min) /
double(gridsize_phi1);
1316 for (
int i_hs = 0; i_hs < gridsize_hs + 1; i_hs++){
1317 double hs_val = hs_min + i_hs * hs_step;
1318 for (
int i_phi1 = 0; i_phi1 < gridsize_phi1 + 1; i_phi1++){
1319 double phi1_val = phi1_min + i_phi1 * phi1_step;
1324 vector<TLorentzVector> pDauVec =
calculate4Momentum(mZZ, mZ1, mZ2, acos(hs_val), acos(costheta1), acos(costheta2), phi1_val, Phi);
1325 for (
int iv=0; iv<2; iv++){
1328 daughters.push_back(tmpPair);
1332 MELAout <<
"Mela::computeP: hs, Phi1 are now " << hs_val <<
" " << phi1_val << endl;
1333 unsigned int idau=1;
1335 MELAout <<
"Dau " << idau <<
" "
1336 <<
"id=" << tmpPart.first <<
" "
1337 <<
"x=" << tmpPart.second.X() <<
" "
1338 <<
"y=" << tmpPart.second.Y() <<
" "
1339 <<
"z=" << tmpPart.second.Z() <<
" "
1340 <<
"t=" << tmpPart.second.T() << endl;
1344 vector<MELAParticle*> partList_tmp;
1345 vector<MELACandidate*> candList_tmp;
1361 for (
MELAParticle* tmpPart:partList_tmp)
delete tmpPart;
1366 prob = prob / float((gridsize_hs + 1) * (gridsize_phi1 +1));
◆ computeP_selfDspin0()
void Mela::computeP_selfDspin0 |
( |
double |
selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2], |
|
|
float & |
prob, |
|
|
bool |
useConstant = true |
|
) |
| |
This is a function that calls Mela::computeP with preset quark, gluon, and vector boson couplings for spin 0.
- See also
- Mela::computeP, Mela::selfDHqqcoupl, Mela::selfDHggcoupl, Mela::selfDHttcoupl, Mela::selfDHbbcoupl, Mela::selfDHzzcoupl, Mela::selfDHwwcoupl
- Warning
- I (Mohit Srivastav) would personally avoid using this function, as it obfuscates the couplings you are really using .
-
It is always a good idea to set couplings yourself so that you know what you're really doing.
- Parameters
-
[in] | selfDHvvcoupl_input | an array corresponding to the couplings of Mela::selfDHzzcoupl (between the Higgs and the Z Boson), Mela::selfDHwwcoupl using the indices in TCouplingsBase.hh |
[out] | prob | This is the value of the output probability - edited by reference |
[in] | useConstant | This turns on the calculation of a corrective constant to different probabilities through Mela::getConstant. If you would like the "pure" MELA calculation to be run, set useConstant to false. By default true. |
Definition at line 1117 of file Mela.cc.
1139 for (
int im=0; im<2; im++){
1141 selfDHzzcoupl[jh][ic][im] = selfDHvvcoupl_input[jh][ic][im];
1142 selfDHwwcoupl[jh][ic][im] = selfDHvvcoupl_input[jh][ic][im];
◆ computeP_selfDspin1() [1/2]
void Mela::computeP_selfDspin1 |
( |
double |
selfDZqqcoupl_input[SIZE_ZQQ][2], |
|
|
double |
selfDZvvcoupl_input[SIZE_ZVV][2], |
|
|
float & |
prob, |
|
|
bool |
useConstant = true |
|
) |
| |
This is a function that calls Mela::computeP with preset quark, gluon, and vector boson couplings for spin 1.
- See also
- Mela::computeP, Mela::selfDHqqcoupl, Mela::selfDHggcoupl, Mela::selfDHttcoupl, Mela::selfDHbbcoupl, Mela::selfDHzzcoupl, Mela::selfDHwwcoupl
- Warning
- I (Mohit Srivastav) would personally avoid using this function, as it obfuscates the couplings you are really using .
-
It is always a good idea to set couplings yourself so that you know what you're really doing.
- Parameters
-
[in] | selfDZqqcoupl_input | an array corresponding to the couplings of Mela::selfDZqqcoupl (between the Higgs and quarks) using the indices in TCouplingsBase.hh |
[in] | selfDZvvcoupl_input | an array corresponding to the couplings of Mela::selfDZvvcoupl (between the Higgs and the Z/W Bosons) using the indices in TCouplingsBase.hh |
[out] | prob | This is the value of the output probability - edited by reference |
[in] | useConstant | This turns on the calculation of a corrective constant to different probabilities through Mela::getConstant. If you would like the "pure" MELA calculation to be run, set useConstant to false. By default true. |
Definition at line 1152 of file Mela.cc.
1158 for (
int im=0; im<2; im++){
◆ computeP_selfDspin1() [2/2]
void Mela::computeP_selfDspin1 |
( |
double |
selfDZvvcoupl_input[SIZE_ZVV][2], |
|
|
float & |
prob, |
|
|
bool |
useConstant = true |
|
) |
| |
This is a function that calls Mela::computeP with preset quark, gluon, and vector boson couplings for spin 1, but sets default values for Zqq couplings.
- See also
- Mela::computeP, Mela::selfDHqqcoupl, Mela::selfDHggcoupl, Mela::selfDHttcoupl, Mela::selfDHbbcoupl, Mela::selfDHzzcoupl, Mela::selfDHwwcoupl
- Warning
- I (Mohit Srivastav) would personally avoid using this function, as it obfuscates the couplings you are really using .
-
It is always a good idea to set couplings yourself so that you know what you're really doing.
- Parameters
-
[in] | selfDZvvcoupl_input | an array corresponding to the couplings of Mela::selfDZvvcoupl (between the Higgs and the Z/W Bosons) using the indices in TCouplingsBase.hh |
[out] | prob | This is the value of the output probability - edited by reference |
[in] | useConstant | This turns on the calculation of a corrective constant to different probabilities through Mela::getConstant. If you would like the "pure" MELA calculation to be run, set useConstant to false. By default true. |
Definition at line 1168 of file Mela.cc.
1177 for (
int im=0; im<2; im++){
◆ computeP_selfDspin2() [1/2]
void Mela::computeP_selfDspin2 |
( |
double |
selfDGggcoupl_input[SIZE_GGG][2], |
|
|
double |
selfDGqqcoupl_input[SIZE_GQQ][2], |
|
|
double |
selfDGvvcoupl_input[SIZE_GVV][2], |
|
|
float & |
prob, |
|
|
bool |
useConstant = true |
|
) |
| |
This is a function that calls Mela::computeP with preset quark, gluon, and vector boson couplings for spin 2.
- See also
- Mela::computeP, Mela::selfDGggcoupl, Mela::selfDGqqcoupl, Mela::selfDGvvcoupl
- Warning
- I (Mohit Srivastav) would personally avoid using this function, as it obfuscates the couplings you are really using .
-
It is always a good idea to set couplings yourself so that you know what you're really doing.
- Parameters
-
[in] | selfDGggcoupl_input | an array corresponding to the couplings of Mela::selfDGggcoupl (between the "graviton" and gluons) using the indices in TCouplingsBase.hh |
[in] | selfDGqqcoupl_input | an array corresponding to the couplings of Mela::selfDGqqcoupl (between the "graviton" and quarks) using the indices in TCouplingsBase.hh |
[in] | selfDGvvcoupl_input | an array corresponding to the couplings of Mela::selfDGvvcoupl (between the "graviton" and the W/Z bosons) using the indices in TCouplingsBase.hh |
[out] | prob | This is the value of the output probability - edited by reference |
[in] | useConstant | This turns on the calculation of a corrective constant to different probabilities through Mela::getConstant. If you would like the "pure" MELA calculation to be run, set useConstant to false. By default true. |
Definition at line 1186 of file Mela.cc.
1193 for (
int im=0; im<2; im++){
◆ computeP_selfDspin2() [2/2]
void Mela::computeP_selfDspin2 |
( |
double |
selfDGggcoupl_input[SIZE_GGG][2], |
|
|
double |
selfDGvvcoupl_input[SIZE_GVV][2], |
|
|
float & |
prob, |
|
|
bool |
useConstant = true |
|
) |
| |
This is a function that calls Mela::computeP with preset quark, gluon, and vector boson couplings for spin 1, but sets default values for Zqq couplings.
- See also
- Mela::computeP, Mela::selfDGggcoupl, Mela::selfDGqqcoupl, Mela::selfDGvvcoupl
- Warning
- I (Mohit Srivastav) would personally avoid using this function, as it obfuscates the couplings you are really using .
-
It is always a good idea to set couplings yourself so that you know what you're really doing.
- Parameters
-
[in] | selfDGggcoupl_input | an array corresponding to the couplings of Mela::selfDGggcoupl (between the "graviton" and gluons) using the indices in TCouplingsBase.hh |
[in] | selfDGvvcoupl_input | an array corresponding to the couplings of Mela::selfDGvvcoupl (between the "graviton" and the W/Z bosons) using the indices in TCouplingsBase.hh |
[out] | prob | This is the value of the output probability - edited by reference |
[in] | useConstant | This turns on the calculation of a corrective constant to different probabilities through Mela::getConstant. If you would like the "pure" MELA calculation to be run, set useConstant to false. By default true. |
Definition at line 1204 of file Mela.cc.
1213 for (
int im=0; im<2; im++){
◆ computePM4l()
Definition at line 1903 of file Mela.cc.
1909 bool hasFailed=
false;
1910 int id_original[2][2];
1911 for (
int iv=0; iv<2; iv++){
1915 for (
int ivd=0; ivd<2; ivd++) id_original[iv][ivd]=Vi->
getDaughter(ivd)->
id;
1919 if (abs(id_original[0][0])==11 && abs(id_original[1][0])==11 && abs(id_original[0][1])==11 && abs(id_original[1][1])==11)
super->
SetDecayChannel(
"4e");
1920 else if (abs(id_original[0][0])==13 && abs(id_original[1][0])==13 && abs(id_original[0][1])==13 && abs(id_original[1][1])==13)
super->
SetDecayChannel(
"4mu");
1922 (abs(id_original[0][0])==11 && abs(id_original[0][1])==11 && abs(id_original[1][0])==13 && abs(id_original[1][1])==13)
1924 (abs(id_original[0][0])==13 && abs(id_original[0][1])==13 && abs(id_original[1][0])==11 && abs(id_original[1][1])==11)
1926 else{
if (
myVerbosity_>=
TVar::ERROR)
MELAerr <<
"Mela::computePM4l: SuperMELA is currently not implemented for decay states other than 4e. 4mu, 2e2mu." << endl; hasFailed=
true; }
1947 std::pair<double, double> m4lP =
super->
M4lProb(mZZtmp);
◆ computeProdDecP() [1/2]
void Mela::computeProdDecP |
( |
double |
selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2], |
|
|
double |
selfDHwwcoupl_input[nSupportedHiggses][SIZE_HVV][2], |
|
|
double |
selfDaTQGCcoupl_input[SIZE_ATQGC][2], |
|
|
double |
selfDAZffcoupl_input[SIZE_AZff][2], |
|
|
float & |
prob, |
|
|
bool |
useConstant = true |
|
) |
| |
computes the combined production and decay probability while taking in coupling arrays
- Attention
- The JHUGen Matrix element is not supported
- See also
- Calls Mela::computeProdDecP(prob, useConstant)
-
Mela::selfDHvvcoupl, Mela::selfDHwwcoupl, Mela::selfDaTQGCcoupl, Mela::selfDAZffcoupl
- Parameters
-
[in] | selfDHvvcoupl_input | input for the couplings set in Mela::selfDHvvcoupl |
[in] | selfDHwwcoupl_input | input for the couplings set in Mela::selfDHwwcoupl |
[in] | selfDaTQGCcoupl_input | input for the couplings set in Mela::selfDaTQGCcoupl |
[in] | selfDAZffcoupl_input | input for the couplings set in Mela::selfDAZffcoupl |
[out] | prob | This is the value of the output probability - edited by reference |
[in] | useConstant | This turns on the calculation of a corrective constant to different probabilities through Mela::getConstant. If you would like the "pure" MELA calculation to be run, set useConstant to false. By default true. |
Definition at line 1469 of file Mela.cc.
1478 for (
int im=0; im<2; im++){
1480 selfDHzzcoupl[jh][ic][im] = selfDHvvcoupl_input[jh][ic][im];
1481 selfDHwwcoupl[jh][ic][im] = selfDHwwcoupl_input[jh][ic][im];
1485 for (
int im=0; im<2; im++){
1488 for (
int im=0; im<2; im++){
◆ computeProdDecP() [2/2]
void Mela::computeProdDecP |
( |
float & |
prob, |
|
|
bool |
useConstant = true |
|
) |
| |
◆ computeProdP() [1/2]
void Mela::computeProdP |
( |
double |
selfDHggcoupl_input[SIZE_HGG][2], |
|
|
double |
selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2], |
|
|
double |
selfDHwwcoupl_input[nSupportedHiggses][SIZE_HVV][2], |
|
|
float & |
prob, |
|
|
bool |
useConstant = true |
|
) |
| |
computes Production side probabilities while taking in coupling arrays
- Attention
- The MCFM Matrix element is not supported
- See also
- Calls Mela::computeProdP(prob, useConstant)
-
Mela:selfDHggcoupl, Mela::selfDHvvcoupl, Mela::selfDHwwcoupl
- Parameters
-
[in] | selfDHggcoupl_input | input for the couplings set in Mela::setlfDHggcoupl |
[in] | selfDHvvcoupl_input | input for the couplings set in Mela::selfDHvvcoupl |
[in] | selfDHwwcoupl_input | input for the couplings set in Mela::selfDHwwcoupl |
[out] | prob | This is the value of the output probability - edited by reference |
[in] | useConstant | This turns on the calculation of a corrective constant to different probabilities through Mela::getConstant. If you would like the "pure" MELA calculation to be run, set useConstant to false. By default true. |
Definition at line 1540 of file Mela.cc.
1547 for (
int im=0; im<2; im++){
for (
int ic=0; ic<
SIZE_HGG; ic++)
selfDHggcoupl[0][ic][im] = selfDHggcoupl_input[ic][im]; }
1549 for (
int im=0; im<2; im++){
1551 selfDHzzcoupl[jh][ic][im] = selfDHvvcoupl_input[jh][ic][im];
1552 selfDHwwcoupl[jh][ic][im] = selfDHwwcoupl_input[jh][ic][im];
◆ computeProdP() [2/2]
void Mela::computeProdP |
( |
float & |
prob, |
|
|
bool |
useConstant = true |
|
) |
| |
computes Production side probabilities using couplings set beforehand
- Attention
- The MCFM Matrix element is not supported
- Parameters
-
[out] | prob | This is the value of the output probability - edited by reference |
[in] | useConstant | This turns on the calculation of a corrective constant to different probabilities through Mela::getConstant. If you would like the "pure" MELA calculation to be run, set useConstant to false. By default true. |
Definition at line 1561 of file Mela.cc.
1573 TLorentzVector nullFourVector(0, 0, 0, 0);
1574 bool isJet2Fake =
false;
1577 unsigned int firstJetIndex=0;
1578 TLorentzVector jet1, higgs;
1579 TLorentzVector jet1massless(0, 0, 0, 0);
1580 TLorentzVector jet2massless(0, 0, 0, 0);
1587 if (tmpPart->passSelection){
1603 const double threshold = 1000.*
LHCsqrts/2.;
1604 TLorentzVector pTotal = higgs+jet1massless+jet2massless;
1605 double sysZ = pTotal.Z();
1606 if (fabs(sysZ)>threshold){
1607 double maxpz2 = threshold - higgs.Z() - jet1massless.Z();
1608 if (fabs(maxpz2)>0.){
1609 double ratio = jet2massless.Z()/maxpz2;
1610 double absp=sqrt(pow(jet2massless.Pt(), 2)+pow(jet2massless.Z()*ratio, 2));
1611 if (
myVerbosity_>=
TVar::INFO)
MELAout <<
"Mela::computeProdP, isJet2Fake=true case: Rescaling pz of fake jet by " << ratio <<
" and energy = " << absp <<
"." << endl;
1612 jet2massless.SetXYZT(jet2massless.X(), jet2massless.Y(), jet2massless.Z()*ratio, absp);
1615 if (
myVerbosity_>=
TVar::INFO)
MELAout <<
"Mela::computeProdP, isJet2Fake=true case: Unable to rescaling pz of fake jet since max(|pz|)<0. Setting to 0 with appropriate energy = pT = " << jet2massless.Pt() <<
"." << endl;
1616 jet2massless.SetXYZT(jet2massless.X(), jet2massless.Y(), 0., jet2massless.Pt());
1625 firstJet->
p4.SetXYZT(jet1massless.X(), jet1massless.Y(), jet1massless.Z(), jet1massless.T());
1634 std::vector<double> etaArray;
1635 std::vector<double> pArray;
1636 double eta_max = 10;
1637 if (jet2massless.Pt()>0.) eta_max = max(eta_max, 1.2*fabs(jet2massless.Eta()));
1638 double eta_min = -eta_max;
1640 for (
int iter=0; iter<nGrid; iter++){
1643 double jet2temp_eta = ((double)iter)*(eta_max-eta_min) / (((
double)nGrid) - 1.) + eta_min;
1644 etaArray.push_back(jet2temp_eta);
1645 double jet2temp_sinh_eta = TMath::SinH(jet2temp_eta);
1646 double jet2temp_pz = jet2massless.Pt()*jet2temp_sinh_eta;
1647 fakeJet.p4.SetZ(jet2temp_pz);
1648 fakeJet.p4.SetX(jet2massless.X()); fakeJet.p4.SetY(jet2massless.Y()); fakeJet.p4.SetT(fakeJet.p4.P());
1651 const double threshold = 1000.*
LHCsqrts/2.;
1652 TLorentzVector pTotal = higgs+jet1massless+fakeJet.p4;
1653 double sys = (pTotal.T()+fabs(pTotal.Z()))/2.;
1654 if (fabs(sys)<threshold){
1658 pArray.push_back((
double)prob_temp);
1661 const double grid_precision = 0.15;
1663 for (
int iG=0; iG<nGrid-1; iG++){
1664 if (pArray[iG]==pArray[iG+1])
continue;
1665 if (etaArray[iG]==etaArray[iG+1])
continue;
1669 TGraph interpolator(nGrid, etaArray.data(), pArray.data());
1670 double derivative_first = (pArray[1]-pArray[0])/(etaArray[1]-etaArray[0]);
1671 double derivative_last = (pArray[nGrid-1]-pArray[nGrid-2])/(etaArray[nGrid-1]-etaArray[nGrid-2]);
1672 TSpline3 spline(
"spline", &interpolator,
"b1e1", derivative_first, derivative_last);
1673 double x_middle = (etaArray[iG]+etaArray[iG+1])*0.5;
1674 double y_middle = (pArray[iG]+pArray[iG+1])*0.5;
1675 double y_sp = spline.Eval(x_middle);
1676 if (y_sp<0) y_sp = 0;
1678 std::vector<double>::iterator gridIt;
1680 if (fabs(y_sp-y_middle)<grid_precision*fabs(y_middle) || fabs(etaArray[iG+1]-etaArray[iG])<1e-3){
1681 gridIt = pArray.begin()+iG+1;
1682 pArray.insert(gridIt, y_sp);
1683 gridIt = etaArray.begin()+iG+1;
1684 etaArray.insert(gridIt, x_middle);
1690 double jet2temp_eta = x_middle;
1691 gridIt = etaArray.begin()+iG+1;
1692 etaArray.insert(gridIt, x_middle);
1693 double jet2temp_sinh_eta = TMath::SinH(jet2temp_eta);
1694 double jet2temp_pz = jet2massless.Pt()*jet2temp_sinh_eta;
1695 fakeJet.p4.SetZ(jet2temp_pz);
1696 fakeJet.p4.SetX(jet2massless.X()); fakeJet.p4.SetY(jet2massless.Y()); fakeJet.p4.SetT(fakeJet.p4.P());
1699 const double threshold = 1000.*
LHCsqrts/2.;
1700 TLorentzVector pTotal = higgs+jet1massless+fakeJet.p4;
1701 double sys = (pTotal.T()+fabs(pTotal.Z()))/2.;
1702 if (fabs(sys)<threshold){
1706 gridIt = pArray.begin()+iG+1;
1707 pArray.insert(gridIt, (
double)prob_temp);
1716 int iGFirst=0, iGLast=nGrid-1;
1717 for (
int iG=1; iG<nGrid; iG++){
1718 if (pArray[iG]>0 && pArray[iG-1]==0){
1723 for (
int iG=nGrid-2; iG>=0; iG--){
1724 if (pArray[iG]>0 && pArray[iG+1]==0){
1729 double dEtaGrid = etaArray[iGLast] - etaArray[iGFirst];
1730 for (
int iG=iGFirst; iG<iGLast-1; iG++){
1731 double dEta = etaArray[iG+1] - etaArray[iG];
1732 double sumProb = pArray[iG]+pArray[iG+1];
1734 dEta = dEta/dEtaGrid;
1735 double addProb = sumProb*dEta;
1739 firstJet->
p4.SetXYZT(jet1.X(), jet1.Y(), jet1.Z(), jet1.T());
1744 if (
melaCand!=candOriginal)
MELAerr <<
"Mela::computeProdP: melaCand!=candOriginal at the end of the fake jet scenario!" << endl;
◆ computeProdP_ttH()
void Mela::computeProdP_ttH |
( |
float & |
prob, |
|
|
int |
topProcess = 2 , |
|
|
int |
topDecay = 0 , |
|
|
bool |
useConstant = true |
|
) |
| |
◆ computeProdP_VH() [1/2]
void Mela::computeProdP_VH |
( |
double |
selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2], |
|
|
float & |
prob, |
|
|
bool |
includeHiggsDecay = false , |
|
|
bool |
useConstant = true |
|
) |
| |
Definition at line 1769 of file Mela.cc.
1778 for (
int im=0; im<2; im++){
1781 selfDHzzcoupl[jh][ic][im] = selfDHvvcoupl_input[jh][ic][im];
1782 selfDHwwcoupl[jh][ic][im] = selfDHvvcoupl_input[jh][ic][im];
◆ computeProdP_VH() [2/2]
void Mela::computeProdP_VH |
( |
float & |
prob, |
|
|
bool |
includeHiggsDecay = false , |
|
|
bool |
useConstant = true |
|
) |
| |
◆ computeTTHAngles()
void Mela::computeTTHAngles |
( |
int |
topDecay, |
|
|
float & |
mT1, |
|
|
float & |
mW1, |
|
|
float & |
mT2, |
|
|
float & |
mW2, |
|
|
float & |
costheta1, |
|
|
float & |
costheta2, |
|
|
float & |
Phi, |
|
|
float & |
costhetastar, |
|
|
float & |
Phi1, |
|
|
float & |
hbb, |
|
|
float & |
hWW, |
|
|
float & |
Phibb, |
|
|
float & |
Phi1bb, |
|
|
float & |
hWplusf, |
|
|
float & |
PhiWplusf, |
|
|
float & |
hWminusf, |
|
|
float & |
PhiWminusf |
|
) |
| |
computes the angles for a ttH system
- See also
- Calls TUtil::computeTTHAngles to fully calculate angles
- Attention
- This function is ill-documented. Use at your own discretion.
- Parameters
-
[in] | topDecay | If topDecay > 0, use unstable tops, otherwise use stable top quarks in the matrix element |
[out] | mT1 | |
[out] | mW1 | |
[out] | mT2 | |
[out] | mW2 | |
[out] | costheta1 | |
[out] | costheta2 | |
[out] | Phi | |
[out] | costhetastar | |
[out] | Phi1 | |
[out] | hbb | |
[out] | hWW | |
[out] | Phibb | |
[out] | Phi1bb | |
[out] | hWplusf | |
[out] | PhiWplusf | |
[out] | hWminusf | |
[out] | PhiWminusf | |
Definition at line 935 of file Mela.cc.
968 =costheta1=costheta2=Phi=costhetastar=Phi1
969 =hbb=hWW=Phibb=Phi1bb
971 =hWminusf=PhiWminusf=0;
975 TLorentzVector
const nullVector(0, 0, 0, 0);
978 int nRequested_Tops=1;
979 int nRequested_Antitops=1;
983 simple_event_record mela_event;
984 mela_event.AssociationCode=partIncCode;
985 mela_event.nRequested_Tops=nRequested_Tops;
986 mela_event.nRequested_Antitops=nRequested_Antitops;
997 for (
unsigned int itd=0; itd<mela_event.pTopDaughters.at(0).size(); itd++) topDaughters.push_back(mela_event.pTopDaughters.at(0).at(itd));
998 for (
unsigned int itd=0; itd<mela_event.pAntitopDaughters.at(0).size(); itd++) antitopDaughters.push_back(mela_event.pAntitopDaughters.at(0).at(itd));
1001 for (
unsigned int itop=0; itop<mela_event.pStableTops.size(); itop++) topDaughters.push_back(mela_event.pStableTops.at(itop));
1002 for (
unsigned int itop=0; itop<mela_event.pStableAntitops.size(); itop++) antitopDaughters.push_back(mela_event.pStableAntitops.at(itop));
1005 if (topDecay==0 && (mela_event.pStableTops.size()<1 || mela_event.pStableAntitops.size()<1)){
1007 <<
"TUtil::TTHiggsMatEl: Number of stable tops (" << mela_event.pStableTops.size() <<
")"
1008 <<
" and number of stable antitops (" << mela_event.pStableAntitops.size() <<
")"
1009 <<
" in ttH process are not 1!" << endl;
1012 else if (topDecay>0 && (mela_event.pTopDaughters.size()<1 || mela_event.pAntitopDaughters.size()<1)){
1014 <<
"Mela::computeTTHAngles: Number of set of top daughters (" << mela_event.pTopDaughters.size() <<
")"
1015 <<
" and number of set of antitop daughters (" << mela_event.pAntitopDaughters.size() <<
")"
1016 <<
" in ttH process are not 1!" << endl;
1019 else if (topDecay>0 && (mela_event.pTopDaughters.at(0).size()!=3 || mela_event.pAntitopDaughters.at(0).size()!=3)){
1021 <<
"Mela::computeTTHAngles: Number of top daughters (" << mela_event.pTopDaughters.at(0).size() <<
")"
1022 <<
" and number of antitop daughters (" << mela_event.pAntitopDaughters.at(0).size() <<
")"
1023 <<
" in ttH process are not 3!" << endl;
1026 if (topDaughters.size()<3){
for (
size_t ip=topDaughters.size(); ip<3; ip++) topDaughters.emplace_back(-9000, nullVector); }
1027 if (antitopDaughters.size()<3){
for (
size_t ip=antitopDaughters.size(); ip<3; ip++) antitopDaughters.emplace_back(-9000, nullVector); }
1030 if (daughters.size()>4){
1033 for (
auto it=daughters.cbegin()+1; it!=daughters.cend(); it++){ firstPart.second = firstPart.second + it->second; }
1034 daughters.erase(daughters.begin()+4, daughters.end());
1036 if (daughters.size()%2==1){
for (
unsigned int ipar=daughters.size(); ipar<4; ipar++) daughters.push_back(
SimpleParticle_t(-9000, nullVector)); }
1037 else if (daughters.size()==2){
1039 daughters.insert(daughters.begin()+1,
SimpleParticle_t(-9000, nullVector));
1046 for (
size_t ip=0; ip<topDaughters.size(); ip++){
1047 if (ip>0){ pT += topDaughters.at(ip).second; pW += topDaughters.at(ip).second; }
1048 else pT += topDaughters.at(ip).second;
1056 for (
size_t ip=0; ip<antitopDaughters.size(); ip++){
1057 if (ip>0){ pT += antitopDaughters.at(ip).second; pW += antitopDaughters.at(ip).second; }
1058 else pT += antitopDaughters.at(ip).second;
1065 costhetastar, costheta1, costheta2, Phi, Phi1,
1066 hbb, hWW, Phibb, Phi1bb,
1068 hWminusf, PhiWminusf,
1070 daughters.at(0).second, daughters.at(0).first,
1071 daughters.at(1).second, daughters.at(1).first,
1072 daughters.at(2).second, daughters.at(2).first,
1073 daughters.at(3).second, daughters.at(3).first,
1075 topDaughters.at(0).second, topDaughters.at(0).first,
1076 topDaughters.at(1).second, topDaughters.at(1).first,
1077 topDaughters.at(2).second, topDaughters.at(2).first,
1079 antitopDaughters.at(0).second, antitopDaughters.at(0).first,
1080 antitopDaughters.at(1).second, antitopDaughters.at(1).first,
1081 antitopDaughters.at(2).second, antitopDaughters.at(2).first,
1083 &(mothers.at(0).second), mothers.at(0).first,
1084 &(mothers.at(1).second), mothers.at(1).first
1088 if (!std::isfinite(costhetastar)) costhetastar=0;
1089 if (!std::isfinite(costheta1)) costheta1=0;
1090 if (!std::isfinite(costheta2)) costheta2=0;
1091 if (!std::isfinite(Phi)) Phi=0;
1092 if (!std::isfinite(Phi1)) Phi1=0;
1094 if (!std::isfinite(hbb)) hbb=0;
1095 if (!std::isfinite(hWW)) hWW=0;
1096 if (!std::isfinite(Phibb)) Phibb=0;
1097 if (!std::isfinite(Phi1bb)) Phi1bb=0;
1099 if (!std::isfinite(hWplusf)) hWplusf=0;
1100 if (!std::isfinite(PhiWplusf)) PhiWplusf=0;
1102 if (!std::isfinite(hWminusf)) hWminusf=0;
1103 if (!std::isfinite(PhiWminusf)) PhiWminusf=0;
1106 <<
"Mela::computeTTHAngles: (h1, h2, Phi, hs, Phi1) = "
1107 << costheta1 <<
", " << costheta2 <<
", " << Phi <<
", "
1108 << costhetastar <<
", " << Phi1 << endl;
◆ computeVBFAngles()
void Mela::computeVBFAngles |
( |
float & |
Q2V1, |
|
|
float & |
Q2V2, |
|
|
float & |
costheta1, |
|
|
float & |
costheta2, |
|
|
float & |
Phi, |
|
|
float & |
costhetastar, |
|
|
float & |
Phi1 |
|
) |
| |
computes the production angles for Vector Boson Fusion (VBF) -> H -> ZZ as defined in Figure 1 of arXiv:1208.4018
- See also
- Calls TUtil::computeVBFAngles to fully calculate angles
-
Use TUtil::computeVBFAngles if you would like to calculate the angles through inputting your own 4-vectors and quantities rather than through the MELA event loop
- Attention
- This function requires there to be an input event already defined
- Warning
- TO COMPUTE VBF ANGLES CORRECTLY DO NOT INPUT A GEN EVENT (an event with mother particles). This will cause the calculator to make the angles relative to the initial-state particles, making the angle a useless discriminator.
-
This function will edit all values inplace. Every value is technically an output.
-
The selection for m1 and m2 are based upon which of the leptons is closer to the z mass. $m_1 \leq m_2$.
- Parameters
-
[out] | Q2V1 | The mass of the first decay particle as reconstructed by 2 of the leptons |
[out] | Q2V2 | The mass of the first decay particle as reconstructed by 2 of the leptons |
[out] | costheta1 | $\cos(\theta_1)$ |
[out] | costheta2 | $\cos(\theta_2)$ |
[out] | Phi | $\phi$ |
[out] | costhetastar | $\cos(\theta*)$ |
[out] | Phi1 | $\phi_1$ |
Definition at line 660 of file Mela.cc.
672 Q2V1=0; Q2V2=0; costheta1=0; costheta2=0; Phi=0; costhetastar=0; Phi1=0;
676 TLorentzVector
const nullVector(0, 0, 0, 0);
678 int nRequested_AssociatedJets=2;
680 simple_event_record mela_event;
681 mela_event.AssociationCode=partIncCode;
682 mela_event.nRequested_AssociatedJets=nRequested_AssociatedJets;
688 if ((
int)aparts.size()!=nRequested_AssociatedJets){
if (
myVerbosity_>=
TVar::ERROR)
MELAerr <<
"Mela::computeVBFAngles: Number of associated particles is not 2!" << endl;
return; }
691 if (daughters.size()>4){
694 for (
auto it=daughters.cbegin()+1; it!=daughters.cend(); it++){ firstPart.second = firstPart.second + it->second; }
695 daughters.erase(daughters.begin()+4, daughters.end());
697 if (daughters.size()%2==1){
for (
unsigned int ipar=daughters.size(); ipar<4; ipar++) daughters.push_back(
SimpleParticle_t(-9000, nullVector)); }
698 else if (daughters.size()==2){
704 MELAout <<
"Mela::computeVBFAngles: Giving the following particles to TUtil::computeVBFAngles:" << endl;
705 for (
unsigned int i=0; i<std::min(daughters.size(), (SimpleParticleCollection_t::size_type) 4); i++)
MELAout << daughters.at(i) << endl;
706 for (
unsigned int i=0; i<std::min(aparts.size(), (SimpleParticleCollection_t::size_type) 2); i++)
MELAout << aparts.at(i) << endl;
707 for (
unsigned int i=0; i<std::min(mothers.size(), (SimpleParticleCollection_t::size_type) 2); i++)
MELAout << mothers.at(i) << endl;
711 costhetastar, costheta1, costheta2, Phi, Phi1, Q2V1, Q2V2,
712 daughters.at(0).second, daughters.at(0).first,
713 daughters.at(1).second, daughters.at(1).first,
714 daughters.at(2).second, daughters.at(2).first,
715 daughters.at(3).second, daughters.at(3).first,
716 aparts.at(0).second, aparts.at(0).first,
717 aparts.at(1).second, aparts.at(1).first,
718 &(mothers.at(0).second), mothers.at(0).first,
719 &(mothers.at(1).second), mothers.at(1).first
723 if (!std::isfinite(costhetastar)) costhetastar=0;
724 if (!std::isfinite(costheta1)) costheta1=0;
725 if (!std::isfinite(costheta2)) costheta2=0;
726 if (!std::isfinite(Phi)) Phi=0;
727 if (!std::isfinite(Phi1)) Phi1=0;
728 if (!std::isfinite(Q2V1)) Q2V1=0;
729 if (!std::isfinite(Q2V2)) Q2V2=0;
732 <<
"Mela::computeVBFAngles: (Q2_1, Q2_2, h1, h2, Phi, hs, Phi1) = "
733 << Q2V1 <<
", " << Q2V2 <<
", "
734 << costheta1 <<
", " << costheta2 <<
", " << Phi <<
", "
735 << costhetastar <<
", " << Phi1 << endl;
◆ computeVBFAngles_ComplexBoost()
void Mela::computeVBFAngles_ComplexBoost |
( |
float & |
Q2V1, |
|
|
float & |
Q2V2, |
|
|
float & |
costheta1_real, |
|
|
float & |
costheta1_imag, |
|
|
float & |
costheta2_real, |
|
|
float & |
costheta2_imag, |
|
|
float & |
Phi, |
|
|
float & |
costhetastar, |
|
|
float & |
Phi1 |
|
) |
| |
Definition at line 742 of file Mela.cc.
754 Q2V1=0; Q2V2=0; costheta1_real=0; costheta2_real=0; costheta1_imag=0; costheta2_imag=0; Phi=0; costhetastar=0; Phi1=0;
758 TLorentzVector
const nullVector(0, 0, 0, 0);
760 int nRequested_AssociatedJets=2;
762 simple_event_record mela_event;
763 mela_event.AssociationCode=partIncCode;
764 mela_event.nRequested_AssociatedJets=nRequested_AssociatedJets;
770 if ((
int) aparts.size()!=nRequested_AssociatedJets){
if (
myVerbosity_>=
TVar::ERROR)
MELAerr <<
"Mela::computeVBFAngles_ComplexBoost: Number of associated particles is not 2!" << endl;
return; }
773 if (daughters.size()>4){
776 for (
auto it=daughters.cbegin()+1; it!=daughters.cend(); it++){ firstPart.second = firstPart.second + it->second; }
777 daughters.erase(daughters.begin()+4, daughters.end());
779 if (daughters.size()%2==1){
for (
unsigned int ipar=daughters.size(); ipar<4; ipar++) daughters.push_back(
SimpleParticle_t(-9000, nullVector)); }
780 else if (daughters.size()==2){
786 costhetastar, costheta1_real, costheta1_imag, costheta2_real, costheta2_imag, Phi, Phi1, Q2V1, Q2V2,
787 daughters.at(0).second, daughters.at(0).first,
788 daughters.at(1).second, daughters.at(1).first,
789 daughters.at(2).second, daughters.at(2).first,
790 daughters.at(3).second, daughters.at(3).first,
791 aparts.at(0).second, aparts.at(0).first,
792 aparts.at(1).second, aparts.at(1).first,
793 &(mothers.at(0).second), mothers.at(0).first,
794 &(mothers.at(1).second), mothers.at(1).first
798 if (!std::isfinite(costhetastar)) costhetastar=0;
799 if (!std::isfinite(costheta1_real)) costheta1_real=0;
800 if (!std::isfinite(costheta2_real)) costheta2_real=0;
801 if (!std::isfinite(costheta1_imag)) costheta1_imag=0;
802 if (!std::isfinite(costheta2_imag)) costheta2_imag=0;
803 if (!std::isfinite(Phi)) Phi=0;
804 if (!std::isfinite(Phi1)) Phi1=0;
805 if (!std::isfinite(Q2V1)) Q2V1=0;
806 if (!std::isfinite(Q2V2)) Q2V2=0;
809 <<
", " << costheta1_real <<
" + " << costheta1_imag <<
"i, "
810 << costheta2_real <<
" + " << costheta2_imag <<
"i, " << Phi <<
", "
811 << costhetastar <<
", " << Phi1 << endl;
813 <<
"Mela::computeVBFAngles_ComplexBoost: (Q2_1, Q2_2, h1, h2, Phi, hs, Phi1) = "
814 << Q2V1 <<
", " << Q2V2 <<
", "
815 << costheta1_real <<
" + " << costheta1_imag <<
"i, "
816 << costheta2_real <<
" + " << costheta2_imag <<
"i, "
817 << Phi <<
", " << costhetastar <<
", " << Phi1 << endl;
819 else if (
myVerbosity_>=
TVar::DEBUG)
MELAerr <<
"Mela::computeVBFAngles_ComplexBoost: No possible melaCand in TEvtProb to compute angles." << endl;
◆ computeVHAngles()
void Mela::computeVHAngles |
( |
float & |
mVstar, |
|
|
float & |
mV, |
|
|
float & |
costheta1, |
|
|
float & |
costheta2, |
|
|
float & |
Phi, |
|
|
float & |
costhetastar, |
|
|
float & |
Phi1 |
|
) |
| |
computes the production for Vector Boson Fusion (VBF) -> H -> ZZ as defined in Figure 1 of arXiv:1208.4018
- See also
- Calls TUtil::computeVHAngles to fully calculate angles
-
Use TUtil::computeVHAngles if you would like to calculate the angles through inputting your own 4-vectors and quantities rather than through the MELA event loop
- Attention
- This function requires there to be an input event already defined
- Warning
- This function will edit all values inplace. Every value is technically an output.
-
The selection for m1 and m2 are based upon which of the leptons is closer to the z mass. $m_1 \leq m_2$.
- Parameters
-
[out] | mVstar | The mass of the virtual Z boson in ZH production |
[out] | mV | The mass of the Z in decay |
[out] | costheta1 | $\cos(\theta_1)$ |
[out] | costheta2 | $\cos(\theta_2)$ |
[out] | Phi | $\phi$ |
[out] | costhetastar | $\cos(\theta*)$ |
[out] | Phi1 | $\phi_1$ |
Definition at line 826 of file Mela.cc.
838 mVstar = 0; mV = 0; costheta1=0; costheta2=0; Phi=0; costhetastar=0; Phi1=0;
842 TLorentzVector
const nullVector(0, 0, 0, 0);
849 int nRequested_AssociatedJets=0;
850 int nRequested_AssociatedLeptons=0;
851 int nRequested_AssociatedPhotons=0;
852 int AssociationVCompatibility=0;
856 nRequested_AssociatedJets=2;
860 nRequested_AssociatedLeptons=2;
864 nRequested_AssociatedPhotons=1;
869 simple_event_record mela_event;
870 mela_event.AssociationCode=partIncCode;
871 mela_event.AssociationVCompatibility=AssociationVCompatibility;
872 mela_event.nRequested_AssociatedJets=nRequested_AssociatedJets;
873 mela_event.nRequested_AssociatedLeptons=nRequested_AssociatedLeptons;
874 mela_event.nRequested_AssociatedPhotons=nRequested_AssociatedPhotons;
882 MELAerr <<
"Mela::computeVHAngles: Number of associated particles (" << aparts.size() <<
") is less than ";
884 else MELAerr << nRequested_AssociatedPhotons;
891 if (daughters.size()>4){
894 for (
auto it=daughters.cbegin()+1; it!=daughters.cend(); it++){ firstPart.second = firstPart.second + it->second; }
895 daughters.erase(daughters.begin()+4, daughters.end());
897 if (daughters.size()%2==1){
for (
unsigned int ipar=daughters.size(); ipar<4; ipar++) daughters.push_back(
SimpleParticle_t(-9000, nullVector)); }
898 else if (daughters.size()==2){
904 costhetastar, costheta1, costheta2, Phi, Phi1, mVstar, mV,
905 daughters.at(0).second, daughters.at(0).first,
906 daughters.at(1).second, daughters.at(1).first,
907 daughters.at(2).second, daughters.at(2).first,
908 daughters.at(3).second, daughters.at(3).first,
909 aparts.at(0).second, aparts.at(0).first,
910 aparts.at(1).second, aparts.at(1).first,
911 &(mothers.at(0).second), mothers.at(0).first,
912 &(mothers.at(1).second), mothers.at(1).first
916 if (!std::isfinite(costhetastar)) costhetastar=0;
917 if (!std::isfinite(costheta1)) costheta1=0;
918 if (!std::isfinite(costheta2)) costheta2=0;
919 if (!std::isfinite(Phi)) Phi=0;
920 if (!std::isfinite(Phi1)) Phi1=0;
921 if (!std::isfinite(mVstar)) mVstar=0;
922 if (!std::isfinite(mV)) mV=0;
925 <<
"Mela::computeVHAngles: (mVstar, mV, h1, h2, Phi, hs, Phi1) = "
926 << mVstar <<
", " << mV <<
", " << costheta1 <<
", " << costheta2 <<
", " << Phi <<
", "
927 << costhetastar <<
", " << Phi1 << endl;
◆ configureAnalyticalPDFs()
bool Mela::configureAnalyticalPDFs |
( |
| ) |
|
|
protected |
Definition at line 2008 of file Mela.cc.
2051 for (
int im=0; im<2; im++){
2175 Double_t fppReal = 1./sqrt(6.) * (
c1/4.*2. + 2.*
c2);
2176 Double_t fppImag = 1./sqrt(6.) *
c5;
2177 Double_t fmmReal = 1./sqrt(6.) * (
c1/4.*2. + 2.*
c2);
2178 Double_t fmmImag = 1./sqrt(6.)*
c5;
2179 Double_t fmpReal = 1./4.*
c1*2.;
2180 Double_t fmpImag = 0;
2181 Double_t fpp = fppImag*fppImag + fppReal*fppReal;
2182 Double_t fmm = fmmImag*fmmImag + fmmReal*fmmReal;
2183 Double_t fmp = fmpImag*fmpImag + fmpReal*fmpReal;
2196 MELAout <<
"Mela::configureAnalyticalPDFs -> ERROR TVar::Process not applicable!!! ME: " <<
myME_ <<
", model: " <<
myModel_ << endl;
2200 if (
pdf==0) noPass=
true;
◆ constructDggr()
void Mela::constructDggr |
( |
float |
bkg_VAMCFM_noscale, |
|
|
float |
ggzz_VAMCFM_noscale, |
|
|
float |
ggHZZ_prob_pure_noscale, |
|
|
float |
ggHZZ_prob_int_noscale, |
|
|
float |
widthScale, |
|
|
float & |
myDggr |
|
) |
| |
|
protected |
Definition at line 1959 of file Mela.cc.
1967 float total_sig_ME = (widthScale * ggHZZ_prob_pure_noscale + sqrt(widthScale) * ggHZZ_prob_int_noscale + ggzz_VAMCFM_noscale);
1968 float total_bkg_ME = bkg_VAMCFM_noscale;
1969 float kd_denominator = (total_sig_ME+total_bkg_ME);
1970 if (kd_denominator>0.) myDggr = total_sig_ME/(total_sig_ME+total_bkg_ME);
◆ deletePConstantHandle()
Definition at line 2824 of file Mela.cc.
2826 delete handle; handle=0;
◆ deletePConstantHandles()
void Mela::deletePConstantHandles |
( |
| ) |
|
|
protected |
◆ getConstant()
void Mela::getConstant |
( |
float & |
prob | ) |
|
This function returns a multiplicative constant to the matrix element calculation.
- Attention
- if useConstant=True in the following functions the constant will be applied to the matrix element:
- computeP_selfDspin0
- computeP_selfDspin1
- computeP_selfDspin1
- computeP_selfDspin2
- computeP_selfDspin2
- computeP
- computeProdDecP
- computeProdDecP
- computeProdP
- computeProdP
- computeProdP_VH
- computeProdP_VH
- computeProdP_ttH
- computePM4l
- Parameters
-
[out] | prob | This is the number that will be multiplied inplace by the constant. |
Definition at line 2207 of file Mela.cc.
◆ getConstant_2l2q()
float Mela::getConstant_2l2q |
( |
| ) |
|
|
protected |
◆ getConstant_4l()
float Mela::getConstant_4l |
( |
| ) |
|
|
protected |
◆ getConstant_4q()
float Mela::getConstant_4q |
( |
| ) |
|
|
protected |
Definition at line 2345 of file Mela.cc.
2348 const int decid = 121;
◆ getConstant_FourFermionDecay()
float Mela::getConstant_FourFermionDecay |
( |
const int & |
decid | ) |
|
|
protected |
Definition at line 2352 of file Mela.cc.
2355 const bool is4mu = (decid==28561);
2356 const bool is4e = (decid==14641 || decid==50625);
2357 const bool is2mu2e = (decid==20449 || decid==27225 || decid==38025 || decid==169 || decid==121);
2359 const unsigned int nPossibleHandles=6;
2362 float constant_tmp=0;
2472 if (hasVBF) pchandle[0]=hvbf;
2473 if (hasZH) pchandle[1]=hzh;
2474 if (hasWH) pchandle[2]=hwh;
2477 if (hasVBF) pchandle[3]=hvbs;
2478 if (hasZH) pchandle[4]=hzzz;
2479 if (hasWH) pchandle[5]=hwzz;
2500 bool hasNullHandle=
true;
2501 for (
unsigned int ihandle=0; ihandle<nPossibleHandles; ihandle++){
if (pchandle[ihandle]!=0){ constant_tmp += pchandle[ihandle]->
Eval(
getIORecord(),
myVerbosity_); hasNullHandle=
false; } }
2502 if (hasNullHandle)
return constant;
2504 constant = constant_tmp;
◆ getConstant_JHUGenUndecayed()
float Mela::getConstant_JHUGenUndecayed |
( |
| ) |
|
|
protected |
Definition at line 2290 of file Mela.cc.
2295 unsigned int iarray=0;
2320 else return constant;
◆ getCurrentCandidate()
◆ getCurrentCandidateIndex()
int Mela::getCurrentCandidateIndex |
( |
| ) |
|
◆ getHiggsWidthAtPoleMass()
double Mela::getHiggsWidthAtPoleMass |
( |
double |
mass | ) |
|
◆ getIORecord()
MelaIO * Mela::getIORecord |
( |
| ) |
|
◆ GetMadgraphCKMElement()
std::complex< double > Mela::GetMadgraphCKMElement |
( |
int |
iquark, |
|
|
int |
jquark |
|
) |
| |
◆ getMeasurablesRRV()
Returns a RooSpin::modelMeasureables object containing many observable quantities.
- Attention
- The parameters in question are:
- $\cos(\theta_1)$
- $\cos(\theta_2)$
- $\phi$
- mass of the first Z boson ($M_{Z1}$)
- mass of the second Z boson ($M_{Z2}$)
- mass of the full particle ($M_{ZZ}$)
- $\cos(\theta*)$
- $\phi_1$
- Y
Definition at line 568 of file Mela.cc.
◆ getNCandidates()
int Mela::getNCandidates |
( |
| ) |
|
◆ getPAux()
void Mela::getPAux |
( |
float & |
prob | ) |
|
◆ getPConstantHandle()
Definition at line 2704 of file Mela.cc.
2715 string cfile_fullpath;
2720 const string path = MELAPKGPATH +
"data/";
2724 const unsigned int npossiblesqrts=3;
2725 const double possible_sqrts[npossiblesqrts]={ 7, 8, 13 };
2726 vector<double> trysqrts;
2727 for (
unsigned isq=0; isq<npossiblesqrts; isq++){
2728 double val = possible_sqrts[isq];
2730 bool inserted=
false;
2731 for (std::vector<double>::iterator it = trysqrts.begin(); it<trysqrts.end(); it++){
2734 trysqrts.insert(it, val);
2738 if (!inserted) trysqrts.push_back(val);
2740 for (
auto& dsqrts:trysqrts){
2741 TString strsqrts = Form(
"%s_%.0f%s", relpath.Data(), dsqrts,
"TeV");
2742 cfile_fullpath = path;
2743 cfile_fullpath.append(strsqrts.Data());
2744 cfile_fullpath.append(
".root");
2745 pchandle =
new MelaPConstant(me_, prod_, proc_, cfile_fullpath.c_str(), spname.Data());
2747 if (
myVerbosity_>=
TVar::DEBUG)
MELAout <<
"Mela::getPConstantHandle: Full path and spline name: " << cfile_fullpath <<
", " << spname <<
" is valid." << endl;
2751 if (
myVerbosity_>=
TVar::DEBUG)
MELAout <<
"Mela::getPConstantHandle: Full path and spline name: " << cfile_fullpath <<
", " << spname <<
" is invalid." << endl;
2757 cfile_fullpath = path;
2758 cfile_fullpath.append(relpath.Data());
2759 cfile_fullpath.append(
".root");
2760 pchandle =
new MelaPConstant(me_, prod_, proc_, cfile_fullpath.c_str(), spname.Data());
2764 if (
myVerbosity_>=
TVar::DEBUG && pchandle==0)
MELAerr <<
"Mela::getPConstantHandle: Handle of " << spname <<
" from " << cfile_fullpath <<
" is invalid!" << endl;
◆ getPConstantHandles()
void Mela::getPConstantHandles |
( |
| ) |
|
|
protected |
Definition at line 2509 of file Mela.cc.
2568 TString filename, spname;
2571 filename =
"pAvgSmooth_JHUGen_JJQCD_HSMHiggs";
2572 spname =
"P_ConserveDifermionMass";
2574 spname =
"P_MomentumToEnergy";
2578 filename =
"pAvgSmooth_JHUGen_JQCD_HSMHiggs";
2579 spname =
"P_ConserveDifermionMass";
2581 spname =
"P_MomentumToEnergy";
2585 filename =
"pAvgSmooth_JHUGen_JJVBF_HSMHiggs";
2586 spname =
"P_ConserveDifermionMass";
2588 spname =
"P_MomentumToEnergy";
2592 filename =
"pAvgSmooth_JHUGen_Had_ZH_HSMHiggs";
2593 spname =
"P_ConserveDifermionMass";
2595 spname =
"P_MomentumToEnergy";
2599 filename =
"pAvgSmooth_JHUGen_Had_WH_HSMHiggs";
2600 spname =
"P_ConserveDifermionMass";
2602 spname =
"P_MomentumToEnergy";
2607 filename =
"pAvgSmooth_MCFM_JJQCD_bkgZJets_13TeV_2l2q";
2608 spname =
"P_ConserveDifermionMass";
2612 filename =
"pAvgSmooth_JHUGen_ZZGG_HSMHiggs";
2613 spname =
"P_ConserveDifermionMass_4mu";
2615 spname =
"P_ConserveDifermionMass_4e";
2617 spname =
"P_ConserveDifermionMass_2mu2e";
2620 filename =
"pAvgSmooth_MCFM_ZZGG_HSMHiggs";
2621 spname =
"P_ConserveDifermionMass_4mu";
2623 spname =
"P_ConserveDifermionMass_4e";
2625 spname =
"P_ConserveDifermionMass_2mu2e";
2628 filename =
"pAvgSmooth_MCFM_JJVBF_S_HSMHiggs";
2629 spname =
"P_ConserveDifermionMass_4mu";
2631 spname =
"P_ConserveDifermionMass_4e";
2633 spname =
"P_ConserveDifermionMass_2mu2e";
2636 filename =
"pAvgSmooth_MCFM_Had_ZH_S_HSMHiggs";
2637 spname =
"P_ConserveDifermionMass_4mu";
2639 spname =
"P_ConserveDifermionMass_4e";
2641 spname =
"P_ConserveDifermionMass_2mu2e";
2644 filename =
"pAvgSmooth_MCFM_Had_WH_S_HSMHiggs";
2645 spname =
"P_ConserveDifermionMass_4mu";
2647 spname =
"P_ConserveDifermionMass_4e";
2649 spname =
"P_ConserveDifermionMass_2mu2e";
2653 filename =
"pAvgSmooth_MCFM_ZZGG_bkgZZ";
2654 spname =
"P_ConserveDifermionMass_4mu";
2656 spname =
"P_ConserveDifermionMass_4e";
2658 spname =
"P_ConserveDifermionMass_2mu2e";
2661 filename =
"pAvgSmooth_MCFM_ZZQQB_bkgZZ";
2662 spname =
"P_ConserveDifermionMass_4mu";
2664 spname =
"P_ConserveDifermionMass_4e";
2666 spname =
"P_ConserveDifermionMass_2mu2e";
2669 filename =
"pAvgSmooth_MCFM_JJVBF_bkgZZ";
2670 spname =
"P_ConserveDifermionMass_4mu";
2672 spname =
"P_ConserveDifermionMass_4e";
2674 spname =
"P_ConserveDifermionMass_2mu2e";
2677 filename =
"pAvgSmooth_MCFM_Had_ZH_bkgZZ";
2678 spname =
"P_ConserveDifermionMass_4mu";
2680 spname =
"P_ConserveDifermionMass_4e";
2682 spname =
"P_ConserveDifermionMass_2mu2e";
2685 filename =
"pAvgSmooth_MCFM_Had_WH_bkgZZ";
2686 spname =
"P_ConserveDifermionMass_4mu";
2688 spname =
"P_ConserveDifermionMass_4e";
2690 spname =
"P_ConserveDifermionMass_2mu2e";
2693 filename =
"pAvgSmooth_MCFM_JJQCD_bkgZZ";
2694 spname =
"P_ConserveDifermionMass_4mu";
2696 spname =
"P_ConserveDifermionMass_4e";
2698 spname =
"P_ConserveDifermionMass_2mu2e";
◆ getPrimaryMass()
double Mela::getPrimaryMass |
( |
int |
ipart | ) |
|
A function to get the current primary EW/QCD parameters from MELA.
- See also
- Wrapper for ZZMatrixElement::get_PrimaryMass, which is a Wrapper for TEvtProb::GetPrimaryMass, which calls TUtil::GetMass.
- Warning
- (Higgs mass/width used in the Matrix Element could be different). Refer to setMelaHiggsMassWidth or setMELAHiggsMass for explicit references to those.
- Attention
- You can input negative ids to edit antiparticles (i.e. 24 -> W+, -24 -> W-)
-
In MCFM both the particle and its antiparticle's mass are changed at the same time. In JHUGen this is not the case - they are separate entries.
-
The particles that you can change when using the MCFM matrix element are as as follow:
absolute value of ipart | Particle Name |
8 | t Prime quark (4th generation) |
7 | b Prime quark (4th generation) |
6 | Top Quark |
5 | Bottom Quark |
4 | Charm Quark |
3 | Strange Quark |
2 | Up Quark |
1 | Down Quark |
11 | Electron |
13 | Muon |
15 | Tau |
23 | Z Boson |
24 | W Boson |
25 | Higgs Boson |
-
The particles that you can change when using the JHUGen matrix element are as follow:
ipart | Particle Name | Sign Sensitive |
0 OR 21 | Gluon | No |
1 | Down Quark | Yes |
2 | Up Quark | Yes |
3 | Strange Quark | Yes |
4 | Charm Quark | Yes |
5 | Bottom Quark | Yes |
6 | Top Quark | Yes |
11 | Electron | Yes |
22 | Photon | No |
23 | Z Boson | No |
24 | W Boson | Yes |
13 | Muon | Yes |
15 | Tau | Yes |
12 | Electron Neutrino | Yes |
14 | Muon Neutrino | Yes |
16 | Tau Neutrino | Yes |
25 | Higgs Boson | No |
32 | Z Prime | No |
33 | Z Prime 2 | No |
34 | W Prime | Yes |
- Parameters
-
[in] | ipart | the particle you would like to get the mass of. |
Definition at line 544 of file Mela.cc.
◆ getPrimaryWidth()
double Mela::getPrimaryWidth |
( |
int |
ipart | ) |
|
A function to get the current primary EW/QCD parameters from MELA.
- See also
- Wrapper for ZZMatrixElement::get_PrimaryMass, which is a Wrapper for TEvtProb::GetPrimaryMass, which calls TUtil::GetMass.
- Warning
- (Higgs mass/width used in the Matrix Element could be different). Refer to setMelaHiggsMassWidth or setMELAHiggsMass for explicit references to those.
- Attention
- You can input negative ids to edit antiparticles (i.e. 24 -> W+, -24 -> W-)
-
In MCFM both the particle and its antiparticle's mass are changed at the same time. In JHUGen this is not the case - they are separate entries.
-
The particles that you can change when using the MCFM matrix element are as as follow:
absolute value of ipart | Particle Name |
15 | Tau |
23 | Z Boson |
24 | W Boson |
25 | Higgs Boson |
-
The particles that you can change when using the JHUGen matrix element are as follow:
ipart | Particle Name | Sign Sensitive |
0 OR 21 | Gluon | No |
1 | Down Quark | Yes |
2 | Up Quark | Yes |
3 | Strange Quark | Yes |
4 | Charm Quark | Yes |
5 | Bottom Quark | Yes |
6 | Top Quark | Yes |
11 | Electron | Yes |
22 | Photon | No |
23 | Z Boson | No |
24 | W Boson | Yes |
13 | Muon | Yes |
15 | Tau | Yes |
12 | Electron Neutrino | Yes |
14 | Muon Neutrino | Yes |
16 | Tau Neutrino | Yes |
25 | Higgs Boson | No |
32 | Z Prime | No |
33 | Z Prime 2 | No |
34 | W Prime | Yes |
- Parameters
-
[in] | ipart | the particle you would like to get the width of. |
Definition at line 545 of file Mela.cc.
◆ getRenFacScaleMode()
◆ getTopCandidateCollection()
◆ getVerbosity()
◆ getXPropagator()
◆ printLogo()
void Mela::printLogo |
( |
| ) |
const |
|
protected |
Definition at line 246 of file Mela.cc.
247 vector<string> logolines;
248 logolines.push_back(
"MELA (Matrix Element Likelihood Approach)");
249 logolines.push_back(
"");
250 logolines.push_back(
"Data analysis and Monte Carlo weights package");
251 logolines.push_back(
"for analyses of resonances produced at pp, ppbar, and e+e- colliders, featuring:");
252 logolines.push_back(
"");
253 logolines.push_back(
"* JHUGenMELA *");
254 logolines.push_back(
"Signal calculations based on analytical pdf.s, and JHU Generator (JHUGen) matrix elements");
255 logolines.push_back(
"(See JHUGen credits below)");
256 logolines.push_back(
"");
257 logolines.push_back(
"* MCFM *");
258 logolines.push_back(
"Signal, background, and interference calculations, modified based on JHUGen matrix elements");
259 logolines.push_back(
"(See MCFM credits below)");
260 logolines.push_back(
"");
261 logolines.push_back(
"* MADGRAPH *");
262 logolines.push_back(
"SMEFTsim calculations for specified processes");
263 logolines.push_back(
"(See MADGRAPH credits below)");
264 logolines.push_back(
"For more details: http://spin.pha.jhu.edu");
265 logolines.push_back(
"");
266 size_t maxlinesize = 0;
267 for (
auto const&
l:logolines) maxlinesize = std::max(maxlinesize,
l.length());
269 unsigned int iline=0;
270 for (
auto const&
l:logolines){
280 logolines.push_back(
"MadGraph5_aMC@NLO");
281 logolines.push_back(
"");
282 logolines.push_back(
"Going Beyond");
283 logolines.push_back(
"");
284 logolines.push_back(
"http://madgraph.hep.uiuc.edu");
285 logolines.push_back(
"http://madgraph.phys.ucl.ac.be");
286 logolines.push_back(
"http://amcatnlo.cern.ch");
287 logolines.push_back(
"");
288 logolines.push_back(
"The MadGraph5_aMC@NLO team");
289 logolines.push_back(
"");
290 logolines.push_back(
"Utilizing the SMEFTsim Package");
291 logolines.push_back(
"I. Brivio, Y. Jiang, M. Trott, arXiv: 1709.06492");
292 logolines.push_back(
"I. Brivio, arXiv: 2012.11343");
293 logolines.push_back(
"");
295 for (
auto const&
l:logolines) maxlinesize = std::max(maxlinesize,
l.length());
298 for (
auto const&
l:logolines){
◆ reset_CandRef()
void Mela::reset_CandRef |
( |
| ) |
|
|
protected |
◆ reset_PAux()
void Mela::reset_PAux |
( |
| ) |
|
|
protected |
◆ reset_SelfDCouplings()
void Mela::reset_SelfDCouplings |
( |
| ) |
|
|
protected |
Definition at line 436 of file Mela.cc.
443 for (
int im=0; im<2; im++){
468 for (
int im=0; im<2; im++){
482 for (
int im=0; im<2; im++){
488 for (
int im=0; im<2; im++){
499 for (
int im=0; im<2; im++){
511 for (
int im=0; im<2; im++){
516 for (
int im=0; im<2; im++){
◆ resetInputEvent()
void Mela::resetInputEvent |
( |
| ) |
|
Resets the event in preparation for the next iteration of the event loop.
- See also
- Wrapper for ZZMatrixElement::reset_InputEvent, which is a wrapper for TEvtProb::ResetInputEvent
- Attention
- Without resetting the input event at the end of each for loop, behavior could be unexpected!
-
It is important to call this at the end of every event loop iteration to clean up TEvtProb!
Definition at line 362 of file Mela.cc.
◆ resetMass()
void Mela::resetMass |
( |
double |
inmass, |
|
|
int |
ipart |
|
) |
| |
Resets the mass for a particle that is an electroweak parameter according to its id.
- See also
- Wrapper for ZZMatrixElement::reset_Mass, which is a wrapper for TEvtProb::ResetMass, which is a wrapper for TUtil::SetMass
-
This finally interfaces with either the function SetMass in mod_parameters.F90 for JHUGen, accessed through .TModParameters.hh
-
in mod_parameters.F90 the conversion from ipart to the values in JHUGen's internal code happens with the function convertLHEreverse(Part)
-
or the masses in the JHUGen-MCFM library at /src/Inc/masses.F or at TMCFM::spinzerohiggs_anomcoupl (for tPrime and bPrime)
- Warning
- It is not recommended to use this function to edit the mass of the Higgs in MCFM. Please use either setMelaHiggsMass or setMelaHiggsMassWidth for that.
-
This function changes BOTH the Yukawa mass as well as the nominal mass in Madgraph! If you would like to decouple the Yukawa mass from the nominal mass please use Mela::resetYukawaMass.
- Attention
- You can input negative ids to edit antiparticles (i.e. 24 -> W+, -24 -> W-)
-
In MCFM/Madgraph both the particle and its antiparticle's mass are changed at the same time. In JHUGen this is not the case - they are separate entries.
-
The particles that you can change when using the MCFM and Madgraph matrix elements are as as follow:
absolute value of ipart | Particle Name |
8 | t Prime quark (4th generation) |
7 | b Prime quark (4th generation) |
6 | Top Quark |
5 | Bottom Quark |
4 | Charm Quark |
3 | Strange Quark |
2 | Up Quark |
1 | Down Quark |
11 | Electron |
13 | Muon |
15 | Tau |
23 | Z Boson |
24 | W Boson |
25 | Higgs Boson |
-
The particles that you can change when using the JHUGen matrix element are as follow:
ipart | Particle Name | Sign Sensitive |
0 OR 21 | Gluon | No |
1 | Down Quark | Yes |
2 | Up Quark | Yes |
3 | Strange Quark | Yes |
4 | Charm Quark | Yes |
5 | Bottom Quark | Yes |
6 | Top Quark | Yes |
11 | Electron | Yes |
22 | Photon | No |
23 | Z Boson | No |
24 | W Boson | Yes |
13 | Muon | Yes |
15 | Tau | Yes |
12 | Electron Neutrino | Yes |
14 | Muon Neutrino | Yes |
16 | Tau Neutrino | Yes |
25 | Higgs Boson | No |
32 | Z Prime | No |
33 | Z Prime 2 | No |
34 | W Prime | Yes |
- Parameters
-
[in] | inmass | the mass that you want in GeV |
[in] | ipart | the particle whose mass you wish to change |
Definition at line 522 of file Mela.cc.
◆ resetMCFM_EWKParameters()
void Mela::resetMCFM_EWKParameters |
( |
double |
ext_Gf, |
|
|
double |
ext_aemmz, |
|
|
double |
ext_mW, |
|
|
double |
ext_mZ, |
|
|
double |
ext_xW, |
|
|
int |
ext_ewscheme = 3 |
|
) |
| |
◆ resetQuarkMasses()
void Mela::resetQuarkMasses |
( |
| ) |
|
Resets the masses of each quark to their original values.
- See also
- Wrapper for ZZMatrixElement::reset_QuarkMasses, which is a wrapper for TEvtProb::ResetQuarkMasses, which calls TEvtProb::ResetMass
- Warning
- You should run this command each time you would like to reset the mass of a quark to mitigate unexpected behavior
- Attention
- The quarks that are reset are as follow, with their default masses in GeV listed as well:
Quark | Mass |
Down | 0.001 |
Up | 0.005 |
Strange | 0.1 |
Charm | 1.275 |
Bottom | 4.75 |
Top | 173.2 |
B Prime | 1e5 |
T Prime | 1e5 |
Definition at line 539 of file Mela.cc.
◆ resetWidth()
void Mela::resetWidth |
( |
double |
inwidth, |
|
|
int |
ipart |
|
) |
| |
Resets the width for a particle that is an electroweak parameter according to its id.
- See also
- Wrapper for ZZMatrixElement::reset_Width, which is a wrapper for TEvtProb::ResetWidth, which is a wrapper for TUtil::SetDecayWidth
-
This finally interfaces with either the function SetDecayWidth in mod_parameters.F90 for JHUGen, accessed through .TModParameters.hh
-
in mod_parameters.F90 the conversion from ipart to the values in JHUGen's internal code happens with the function convertLHEreverse(Part)
-
or the masses in the JHUGen-MCFM library at /src/Inc/masses.F or at TMCFM::spinzerohiggs_anomcoupl (for tPrime and bPrime)
- Warning
- It is not recommended to use this function to edit the mass of the Higgs in MCFM. Please use either setMelaHiggsWidth or setMelaHiggsMassWidth for that.
- Attention
- You can input negative ids to edit antiparticles (i.e. 24 -> W+, -24 -> W-)
-
In MCFM both the particle and its antiparticle's mass are changed at the same time. In JHUGen this is not the case - they are separate entries.
-
The particles that you can change when using the MCFM matrix element are as as follow:
absolute value of ipart | Particle Name |
15 | Tau |
23 | Z Boson |
24 | W Boson |
25 | Higgs Boson |
-
The particles that you can change when using the JHUGen matrix element are as follow:
ipart | Particle Name | Sign Sensitive |
0 OR 21 | Gluon | No |
1 | Down Quark | Yes |
2 | Up Quark | Yes |
3 | Strange Quark | Yes |
4 | Charm Quark | Yes |
5 | Bottom Quark | Yes |
6 | Top Quark | Yes |
11 | Electron | Yes |
22 | Photon | No |
23 | Z Boson | No |
24 | W Boson | Yes |
13 | Muon | Yes |
15 | Tau | Yes |
12 | Electron Neutrino | Yes |
14 | Muon Neutrino | Yes |
16 | Tau Neutrino | Yes |
25 | Higgs Boson | No |
32 | Z Prime | No |
33 | Z Prime 2 | No |
34 | W Prime | Yes |
- Parameters
-
[in] | inwidth | the width that you want in GeV |
[in] | ipart | the particle whose mass you wish to change |
Definition at line 523 of file Mela.cc.
◆ resetYukawaMass()
void Mela::resetYukawaMass |
( |
double |
inmass, |
|
|
int |
ipart |
|
) |
| |
Resets the Yukawa coupling (or "Yukawa mass") for a specific particle. Does not change its "intrinsic" mass.
- See also
- Resets the values
mdl_ym{b,c,s,up,do,e,m,tau}
in the params_r_
common block for Madgraph
- Attention
- The particles whose Yukawa mass you can change are as follow:
ipart | Particle Name |
1 | Down Quark |
2 | Up Quark |
3 | Strange Quark |
4 | Charm Quark |
5 | Bottom Quark |
6 | Top Quark |
11 | Electron |
13 | Muon |
15 | Tau |
- Parameters
-
[in] | inmass | the Yukawa mass you want in GeV |
[in] | ipart | the particle whose Yukawa mass you wish to change |
Definition at line 524 of file Mela.cc.
525 const int ipartabs = abs(ipart);
536 MELAerr <<
"Particle with id " << ipart <<
" does not have supported Yukawa Couplings!" << endl;
◆ setATQGCCouplings()
void Mela::setATQGCCouplings |
( |
| ) |
|
|
protected |
◆ setAZffCouplings()
void Mela::setAZffCouplings |
( |
| ) |
|
|
protected |
◆ setCandidateDecayMode()
◆ setConstant()
void Mela::setConstant |
( |
| ) |
|
|
protected |
Definition at line 2214 of file Mela.cc.
2285 if (
std::isnan(constant) || std::isinf(constant) || constant<=0.) constant=0;
2286 else constant=1./constant;
◆ setCurrentCandidate()
◆ setCurrentCandidateFromIndex()
void Mela::setCurrentCandidateFromIndex |
( |
unsigned int |
icand | ) |
|
◆ setInputEvent()
Sets the input event for MELA. MELA cannot run without this.
- See also
- Wrapper for ZZMatrixElement::set_InputEvent, which is a wrapper for TEvtProb::SetInputEvent, which calls TUtil::ConvertVectorFormat
- Attention
- An input event must be set for each event in an event loop, otherwise MELA will throw a segmentation error
- Parameters
-
[in] | pDaughters | A SimpleParticleCollection_t of particle daughters (decay products) |
[in] | pAssociated | A SimpleParticleCollection_t of associated particles (i.e. jets), by default 0 (no jets) |
[in] | pMothers | A SimpleParticleCollection_t of particle mothers (i.e. gluons), by default 0 (reco data contains no mother information) |
[in] | isGen | A boolean signifying whether the event in question is a Gen event or a reco event, by default false (reco) |
Definition at line 349 of file Mela.cc.
◆ SetMadgraphCKMElements()
void Mela::SetMadgraphCKMElements |
( |
double |
ckmlambda = 0.2265 , |
|
|
double |
ckma = 0.79 , |
|
|
double |
ckmrho = 0.141 , |
|
|
double |
ckmeta = 0.357 , |
|
|
bool |
force_refresh = false |
|
) |
| |
◆ setMelaHiggsMass()
void Mela::setMelaHiggsMass |
( |
double |
myHiggsMass, |
|
|
int |
index = 0 |
|
) |
| |
Sets the mass of your chosen Higgs.
- See also
- Wrapper for ZZMatrixElement::set_mHiggs.
- Parameters
-
[in] | myHiggsMass | This is the mass of the Higgs that you would like, in GeV |
[in] | index | This is either 0 or 1, depending on if you want to change the first or second resonance. By default it is 0. |
Definition at line 342 of file Mela.cc.
◆ setMelaHiggsMassWidth()
void Mela::setMelaHiggsMassWidth |
( |
double |
myHiggsMass, |
|
|
double |
myHiggsWidth, |
|
|
int |
index |
|
) |
| |
a combination of setMelaHiggsMass and setMelaHiggsWidth.
- See also
- Wrapper for ZZMatrixElement::set_mHiggs_wHiggs.
- Parameters
-
[in] | myHiggsMass | This is the mass of the Higgs that you would like, in GeV |
[in] | myHiggsWidth | This is the mass of the Higgs that you would like, in GeV |
[in] | index | This is either 0 or 1, depending on if you want to change the first or second resonance. By default it is 0. |
Definition at line 344 of file Mela.cc.
◆ setMelaHiggsWidth()
void Mela::setMelaHiggsWidth |
( |
double |
myHiggsWidth = -1 , |
|
|
int |
index = 0 |
|
) |
| |
Sets the width of your chosen Higgs.
- See also
- Wrapper for ZZMatrixElement::set_wHiggs
- Parameters
-
[in] | myHiggsWidth | This is the mass of the Higgs that you would like, in GeV |
[in] | index | This is either 0 or 1, depending on if you want to change the first or second resonance. By default it is 0. |
Definition at line 343 of file Mela.cc.
◆ setMelaLeptonInterference()
Sets the MELA Lepton Interference.
- Parameters
-
[in] | myLepInterf | sets the myLepInterf_ variable to one found in TVar::LeptonInterface through ZZMatrixElement::set_LeptonInterface. TVar::DefaultLeptonInterf by default. |
Definition at line 345 of file Mela.cc.
◆ setMelaPrimaryHiggsMass()
void Mela::setMelaPrimaryHiggsMass |
( |
double |
myHiggsMass | ) |
|
Sets the mass of the "primary" higgs.
- See also
- Wrapper for the function ZZMatrixElement::set_PrimaryHiggsMass, which is a wrapper for TEvtProb::SetPrimaryHiggsMass
- Parameters
-
[in] | myHiggsMass | This is the mass of the Higgs that you would like, in GeV |
- Attention
- The primary Higgs is the first resonance - nominally MELA can have 2 resonances when working through JHUGen-MCFM
Definition at line 341 of file Mela.cc.
◆ setProcess()
Sets the process, matrix element, and production that MELA is to use for this event. Calls ZZMatrixElement::set_Process, which calls TEvtProb::SetProcess.
- Attention
- Remember to set the process for each event, otherwise the MELA event loop will throw a segmentation error.
- Parameters
-
Definition at line 309 of file Mela.cc.
326 MELAout <<
"Production mode " <<
myProduction_ <<
" is not currently supported by MADMELA!" << endl;
◆ setRemoveJetMasses()
void Mela::setRemoveJetMasses |
( |
bool |
MasslessLeptonSwitch = true | ) |
|
either permits or forbids massive jets.
- See also
- Wrapper for the function TUtil::applyJetMassCorrection
- Parameters
-
[in] | MasslessLeptonSwitch | Whether you would like your jets to be massive, by default true. |
Definition at line 557 of file Mela.cc.
◆ setRemoveLeptonMasses()
void Mela::setRemoveLeptonMasses |
( |
bool |
MasslessLeptonSwitch = true | ) |
|
either permits or forbids massive leptons.
- See also
- Wrapper for the function TUtil::applyLeptonMassCorrection
- Parameters
-
[in] | MasslessLeptonSwitch | Whether you would like your leptons to be massive, by default true. |
Definition at line 556 of file Mela.cc.
◆ setRenFacScaleMode()
◆ setSpinOneCouplings()
void Mela::setSpinOneCouplings |
( |
| ) |
|
|
protected |
◆ setSpinTwoCouplings()
void Mela::setSpinTwoCouplings |
( |
| ) |
|
|
protected |
◆ setSpinZeroCouplings()
void Mela::setSpinZeroCouplings |
( |
| ) |
|
|
protected |
◆ setTempCandidate()
Sets a temporary MELA candidate, by setting melaCand in Xcal2 to a temporary candidate without pushing this candidate to the candList of Xcal2.
- See also
- Wrapper for ZZMatrixElement::set_TempCandidate
- Parameters
-
[in] | pDaughters | A SimpleParticleCollection_t of particle daughters (decay products) |
[in] | pAssociated | A SimpleParticleCollection_t of associated particles (i.e. jets), by default 0 (no jets) |
[in] | pMothers | A SimpleParticleCollection_t of particle mothers (i.e. gluons), by default 0 (reco data contains no mother information) |
[in] | isGen | A boolean signifying whether the event in question is a Gen event or a reco event, by default false (reco) |
Definition at line 366 of file Mela.cc.
◆ setVerbosity()
Sets the verbosity for MELA outside of the initial constructor.
- Parameters
-
[in] | verbosity_ | The verbosity of MELA that you desire, as defined in TVar::ERROR. |
Definition at line 330 of file Mela.cc.
◆ auxiliaryProb
float Mela::auxiliaryProb |
|
protected |
◆ costheta1_rrv
RooRealVar* Mela::costheta1_rrv |
◆ costheta2_rrv
RooRealVar* Mela::costheta2_rrv |
◆ costhetastar_rrv
RooRealVar* Mela::costhetastar_rrv |
◆ ggSpin0Model
◆ LHCsqrts
◆ melaCand
◆ melaRandomNumber
TRandom3 Mela::melaRandomNumber |
◆ myLepInterf_
◆ myME_
◆ myModel_
◆ myProduction_
◆ myVerbosity_
◆ mzz_rrv
RooRealVar* Mela::mzz_rrv |
◆ pAvgSmooth_JHUGen_Had_WH_HSMHiggs
◆ pAvgSmooth_JHUGen_Had_ZH_HSMHiggs
◆ pAvgSmooth_JHUGen_JJQCD_HSMHiggs
◆ pAvgSmooth_JHUGen_JJVBF_HSMHiggs
◆ pAvgSmooth_JHUGen_JQCD_HSMHiggs
◆ pAvgSmooth_JHUGen_ZZGG_HSMHiggs_2mu2e
◆ pAvgSmooth_JHUGen_ZZGG_HSMHiggs_4e
◆ pAvgSmooth_JHUGen_ZZGG_HSMHiggs_4mu
◆ pAvgSmooth_MCFM_Had_WH_bkgZZ_2mu2e
◆ pAvgSmooth_MCFM_Had_WH_bkgZZ_4e
◆ pAvgSmooth_MCFM_Had_WH_bkgZZ_4mu
◆ pAvgSmooth_MCFM_Had_WH_S_HSMHiggs_2mu2e
◆ pAvgSmooth_MCFM_Had_WH_S_HSMHiggs_4e
◆ pAvgSmooth_MCFM_Had_WH_S_HSMHiggs_4mu
◆ pAvgSmooth_MCFM_Had_ZH_bkgZZ_2mu2e
◆ pAvgSmooth_MCFM_Had_ZH_bkgZZ_4e
◆ pAvgSmooth_MCFM_Had_ZH_bkgZZ_4mu
◆ pAvgSmooth_MCFM_Had_ZH_S_HSMHiggs_2mu2e
◆ pAvgSmooth_MCFM_Had_ZH_S_HSMHiggs_4e
◆ pAvgSmooth_MCFM_Had_ZH_S_HSMHiggs_4mu
◆ pAvgSmooth_MCFM_JJQCD_bkgZJets_2l2q
◆ pAvgSmooth_MCFM_JJQCD_bkgZZ_2mu2e
◆ pAvgSmooth_MCFM_JJQCD_bkgZZ_4e
◆ pAvgSmooth_MCFM_JJQCD_bkgZZ_4mu
◆ pAvgSmooth_MCFM_JJVBF_bkgZZ_2mu2e
◆ pAvgSmooth_MCFM_JJVBF_bkgZZ_4e
◆ pAvgSmooth_MCFM_JJVBF_bkgZZ_4mu
◆ pAvgSmooth_MCFM_JJVBF_S_HSMHiggs_2mu2e
◆ pAvgSmooth_MCFM_JJVBF_S_HSMHiggs_4e
◆ pAvgSmooth_MCFM_JJVBF_S_HSMHiggs_4mu
◆ pAvgSmooth_MCFM_ZZGG_bkgZZ_2mu2e
◆ pAvgSmooth_MCFM_ZZGG_bkgZZ_4e
◆ pAvgSmooth_MCFM_ZZGG_bkgZZ_4mu
◆ pAvgSmooth_MCFM_ZZGG_HSMHiggs_2mu2e
◆ pAvgSmooth_MCFM_ZZGG_HSMHiggs_4e
◆ pAvgSmooth_MCFM_ZZGG_HSMHiggs_4mu
◆ pAvgSmooth_MCFM_ZZQQB_bkgZZ_2mu2e
◆ pAvgSmooth_MCFM_ZZQQB_bkgZZ_4e
◆ pAvgSmooth_MCFM_ZZQQB_bkgZZ_4mu
◆ pdf
◆ phi1_rrv
RooRealVar* Mela::phi1_rrv |
◆ phi_rrv
RooRealVar* Mela::phi_rrv |
◆ qqZZmodel
◆ spin1Model
◆ spin2Model
◆ super
◆ superDijet
◆ upFrac_rrv
RooRealVar* Mela::upFrac_rrv |
◆ Y_rrv
◆ z1mass_rrv
RooRealVar* Mela::z1mass_rrv |
◆ z2mass_rrv
RooRealVar* Mela::z2mass_rrv |
◆ ZZME
The documentation for this class was generated from the following files:
RooAbsReal * gzgs4List[1][2]
void reset_Width(double inmass, int ipart)
void getConstant(float &prob)
This function returns a multiplicative constant to the matrix element calculation.
std::pair< int, TLorentzVector > SimpleParticle_t
void scaleMomentumToEnergy(const TLorentzVector &massiveJet, TLorentzVector &masslessJet, double mass=0)
RooSpinTwo::modelCouplings couplings
MelaPConstant * pAvgSmooth_MCFM_ZZGG_bkgZZ_4e
bool isALepton(const int id)
MelaPConstant * pAvgSmooth_MCFM_ZZGG_HSMHiggs_4e
MelaPConstant * pAvgSmooth_MCFM_ZZQQB_bkgZZ_4mu
void setVerbosity(TVar::VerbosityLevel verbosity)
MELAParticle * getDaughter(int index) const
virtual void resetHypotheses()
MelaPConstant * pAvgSmooth_JHUGen_ZZGG_HSMHiggs_2mu2e
double selfDGqqcoupl[SIZE_GQQ][2]
double selfDGvpvpcoupl[SIZE_GVV][2]
MelaPConstant * pAvgSmooth_MCFM_ZZQQB_bkgZZ_4e
void initialize_madMELA()
Initializes all of the values in Madgraph for proper usage. These include values like alpha,...
virtual void setTensorPolarization(int ig, double initval)
MelaPConstant * pAvgSmooth_MCFM_Had_ZH_bkgZZ_4e
TVar::VerbosityLevel myVerbosity_
void build(double mh_)
This is the actual building of the tool that occurs in each instance of the Mela::Mela constructor.
MelaPConstant * pAvgSmooth_MCFM_JJVBF_S_HSMHiggs_2mu2e
void addAssociatedJet(MELAParticle *myParticle)
void get_XPropagator(TVar::ResonancePropagatorScheme scheme, float &prop)
void computeDecayAngles(float &qH, float &m1, float &m2, float &costheta1, float &costheta2, float &Phi, float &costhetastar, float &Phi1)
computes the decay angles for gg -> H -> ZZ as defined in Figure 1 of arXiv:1001.3396
MelaPConstant * pAvgSmooth_MCFM_JJQCD_bkgZZ_4e
void computeProdP_VH(double selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2], float &prob, bool includeHiggsDecay=false, bool useConstant=true)
void set_mHiggs_wHiggs(double mh_, double gah_, int index)
MelaPConstant * pAvgSmooth_MCFM_Had_ZH_bkgZZ_4mu
virtual void addHypothesis(int ig, double initval, double iphase=0)
int getNDaughters() const
std::map< std::pair< TVar::Process, TVar::Production >, MG_process_double > * updateMap
std::vector< MELATopCandidate_t * > * get_TopCandidateCollection()
int selfDHwwCLambda_qsq[nSupportedHiggses][SIZE_HVV_CQSQ]
void reset_Mass(double inmass, int ipart)
@ nFermionMassRemovalSchemes
MelaIO * getIORecord()
Returns the MELAIO object, and by consequence, the entire parton-by-parton matrix element record.
MELAParticle * getSortedV(int index) const
TString ProductionName(TVar::Production temp)
MelaPConstant * pAvgSmooth_JHUGen_ZZGG_HSMHiggs_4mu
double selfDHttcoupl[nSupportedHiggses][SIZE_HQQ][2]
void setRemoveLeptonMasses(bool MasslessLeptonSwitch=true)
either permits or forbids massive leptons.
void SetVerbosity(bool verb=true)
MelaPConstant * pAvgSmooth_MCFM_JJQCD_bkgZZ_4mu
void setCurrentCandidate(MELACandidate *cand)
Switches the candidate that you are working on to another candidate object specified.
MelaPConstant * pAvgSmooth_MCFM_ZZGG_bkgZZ_2mu2e
RooAbsReal * g2List[8][2]
RooAbsReal * gzgs3List[1][2]
void setMelaHiggsWidth(double myHiggsWidth=-1, int index=0)
Sets the width of your chosen Higgs.
void computeProdXS_JH(float &mevalue)
void computeTTHAngles(float &hs, float &hincoming, float &hTT, float &PhiTT, float &Phi1, float &hbb, float &hWW, float &Phibb, float &Phi1bb, float &hWplusf, float &PhiWplusf, float &hWminusf, float &PhiWminusf, TLorentzVector p4M11, int Z1_lept1Id, TLorentzVector p4M12, int Z1_lept2Id, TLorentzVector p4M21, int Z2_lept1Id, TLorentzVector p4M22, int Z2_lept2Id, TLorentzVector b, int bId, TLorentzVector Wplusf, int WplusfId, TLorentzVector Wplusfb, int WplusfbId, TLorentzVector bbar, int bbarId, TLorentzVector Wminusf, int WminusfId, TLorentzVector Wminusfb, int WminusfbId, TLorentzVector *injet1=0, int injet1Id=0, TLorentzVector *injet2=0, int injet2Id=0)
MelaPConstant * getPConstantHandle(TVar::MatrixElement me_, TVar::Production prod_, TVar::Process proc_, TString relpath, TString spname, const bool useSqrts=false)
void computeProdP(double selfDHggcoupl_input[SIZE_HGG][2], double selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2], double selfDHwwcoupl_input[nSupportedHiggses][SIZE_HVV][2], float &prob, bool useConstant=true)
computes Production side probabilities while taking in coupling arrays
double selfDHzpzpcoupl[SIZE_HVV][2]
MelaPConstant * pAvgSmooth_MCFM_Had_WH_S_HSMHiggs_4e
void applyLeptonMassCorrection(bool flag=true)
double getMEConst() const
void SetVerbosity(const TVar::VerbosityLevel verbosity_)
@ kUseAssociated_StableTops
double selfDWpffcoupl[SIZE_Vpff][2]
MelaPConstant * pAvgSmooth_MCFM_JJVBF_S_HSMHiggs_4mu
double selfDZpffcoupl[SIZE_Vpff][2]
double selfDaTQGCcoupl[SIZE_ATQGC][2]
TVar::Production myProduction_
MelaPConstant * pAvgSmooth_MCFM_Had_ZH_S_HSMHiggs_4mu
void computeConstant(float &prob)
void set_Verbosity(TVar::VerbosityLevel verbosity_)
void set_SpinZeroContact(double selfDHzzpcoupl[SIZE_HVV][2], double selfDHzpzpcoupl[SIZE_HVV][2], double selfDHwwpcoupl[SIZE_HVV][2], double selfDHwpwpcoupl[SIZE_HVV][2])
MelaPConstant * pAvgSmooth_MCFM_Had_WH_S_HSMHiggs_2mu2e
void set_SpinOneCouplings(double selfDZqqcoupl[SIZE_ZQQ][2], double selfDZvvcoupl[SIZE_ZVV][2])
void PrintCandidateSummary(MELACandidate *cand)
bool isAZBoson(const int id)
void setProcess(TVar::Process myModel, TVar::MatrixElement myME, TVar::Production myProduction)
Sets the process, matrix element, and production that MELA is to use for this event....
double get_HiggsWidthAtPoleMass(double mass)
void deletePConstantHandles()
void setSpinZeroCouplings()
MelaPConstant * pAvgSmooth_JHUGen_JJVBF_HSMHiggs[TVar::nFermionMassRemovalSchemes-1]
logical function isnan(x)
std::vector< TLorentzVector > Calculate4Momentum(double Mx, double M1, double M2, double theta, double theta1, double theta2, double Phi1, double Phi)
RooqqZZ_JHU_ZgammaZZ_fast * qqZZmodel
double selfDHqqcoupl[nSupportedHiggses][SIZE_HQQ][2]
MelaPConstant * pAvgSmooth_MCFM_JJVBF_bkgZZ_2mu2e
ScalarPdfFactory_HVV * ggSpin0Model
RooRealVar * costheta1_rrv
void set_SpinZeroCouplings(double selfDHggcoupl[nSupportedHiggses][SIZE_HGG][2], double selfDHg4g4coupl[nSupportedHiggses][SIZE_HGG][2], double selfDHqqcoupl[nSupportedHiggses][SIZE_HQQ][2], double selfDHbbcoupl[nSupportedHiggses][SIZE_HQQ][2], double selfDHttcoupl[nSupportedHiggses][SIZE_HQQ][2], double selfDHb4b4coupl[nSupportedHiggses][SIZE_HQQ][2], double selfDHt4t4coupl[nSupportedHiggses][SIZE_HQQ][2], double selfDHzzcoupl[nSupportedHiggses][SIZE_HVV][2], double selfDHwwcoupl[nSupportedHiggses][SIZE_HVV][2], double selfDHzzLambda_qsq[nSupportedHiggses][SIZE_HVV_LAMBDAQSQ][SIZE_HVV_CQSQ], double selfDHwwLambda_qsq[nSupportedHiggses][SIZE_HVV_LAMBDAQSQ][SIZE_HVV_CQSQ], int selfDHzzCLambda_qsq[nSupportedHiggses][SIZE_HVV_CQSQ], int selfDHwwCLambda_qsq[nSupportedHiggses][SIZE_HVV_CQSQ], double selfDSMEFTSimcoupl[SIZE_SMEFT], bool diffHWW=false)
MelaPConstant * pAvgSmooth_MCFM_Had_WH_bkgZZ_2mu2e
void append_TopCandidate(SimpleParticleCollection_t *TopDaughters)
double selfDGvvcoupl[SIZE_GVV][2]
MELACandidate * shallowCopy()
bool isAPhoton(const int id)
void set_TempCandidate(SimpleParticleCollection_t *pDaughters, SimpleParticleCollection_t *pAssociated=0, SimpleParticleCollection_t *pMothers=0, bool isGen=false)
void set_CandidateDecayMode(TVar::CandidateDecayMode mode)
double get_PrimaryMass(int ipart)
void set_AZffCouplings(double selfDAZffcoupl[SIZE_AZff][2])
void set_SpinTwoCouplings(double selfDGqqcoupl[SIZE_GQQ][2], double selfDGggcoupl[SIZE_GGG][2], double selfDGvvcoupl[SIZE_GVV][2])
void computeP_selfDspin0(double selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2], float &prob, bool useConstant=true)
This is a function that calls Mela::computeP with preset quark, gluon, and vector boson couplings for...
int get_CurrentCandidateIndex()
void set_Process(TVar::Process process_, TVar::MatrixElement me_, TVar::Production production_)
MELACandidate * get_CurrentCandidate()
MELAOutputStreamer MELAout
void setDefaultMadgraphValues()
Sets the Madgraph values to their defaults. This function is called after every call to a compute fun...
MelaPConstant * pAvgSmooth_JHUGen_JJQCD_HSMHiggs[TVar::nFermionMassRemovalSchemes-1]
MelaPConstant * pAvgSmooth_MCFM_ZZGG_HSMHiggs_4mu
void set_mHiggs(double mh_, int index)
bool differentiate_HWW_HZZ
MelaPConstant * pAvgSmooth_MCFM_JJVBF_bkgZZ_4e
void set_VprimeContactCouplings(double selfDZpffcoupl[SIZE_Vpff][2], double selfDWpffcoupl[SIZE_Vpff][2], double M_Zprime, double Ga_Zprime, double M_Wprime, double Ga_Wprime)
void SetPathToCards(std::string dirToCards)
virtual void addHypothesis(int ig, int ilam, double iphase=0, double altparam_fracval=0)
MELACandidate * ConvertVectorFormat(SimpleParticleCollection_t *pDaughters, SimpleParticleCollection_t *pAssociated, SimpleParticleCollection_t *pMothers, bool isGen, std::vector< MELAParticle * > *particleList, std::vector< MELACandidate * > *candList)
MelaPConstant * pAvgSmooth_MCFM_JJVBF_bkgZZ_4mu
MelaPConstant * pAvgSmooth_MCFM_Had_WH_bkgZZ_4e
std::pair< double, double > M4lProb(double m4l)
complex(8), parameter, public c5
double Eval(const MelaIO *RcdME, TVar::VerbosityLevel verbosity) const
RooAbsReal * Lambda_z4qsq[SIZE_HVV_CQSQ]
double selfDHg4g4coupl[nSupportedHiggses][SIZE_HGG][2]
MelaPConstant * pAvgSmooth_JHUGen_JQCD_HSMHiggs[TVar::nFermionMassRemovalSchemes-1]
void computeP(float &prob, bool useConstant=true)
Computes the probability for the probabilities on the decay side of things using the constituent daug...
int configure(TVar::Process model_)
const TVar::event_scales_type & get_RenFacScaleMode() const
RooAbsReal * Lambda_z1qsq[SIZE_HVV_CQSQ]
std::string GetMELAPath()
float getConstant_JHUGenUndecayed()
void set_InputEvent(SimpleParticleCollection_t *pDaughters, SimpleParticleCollection_t *pAssociated=0, SimpleParticleCollection_t *pMothers=0, bool isGen=false)
void setCandidateDecayMode(TVar::CandidateDecayMode mode)
Sets the decay mode for your event.
void set_CurrentCandidateFromIndex(unsigned int icand)
void setMelaHiggsMass(double myHiggsMass, int index=0)
Sets the mass of your chosen Higgs.
TensorPdfFactory_ppHVV * spin2Model
void makeParamsConst(bool yesNo=true)
double selfDZqqcoupl[SIZE_ZQQ][2]
RooAbsReal * ggsgs4List[1][2]
double selfDHb4b4coupl[nSupportedHiggses][SIZE_HQQ][2]
double selfDHwpwpcoupl[SIZE_HVV][2]
void computeVHAngles(float &costhetastar, float &costheta1, float &costheta2, float &Phi, float &Phi1, float &m1, float &m2, TLorentzVector p4M11, int Z1_lept1Id, TLorentzVector p4M12, int Z1_lept2Id, TLorentzVector p4M21, int Z2_lept1Id, TLorentzVector p4M22, int Z2_lept2Id, TLorentzVector jet1, int jet1Id, TLorentzVector jet2, int jet2Id, TLorentzVector *injet1=0, int injet1Id=0, TLorentzVector *injet2=0, int injet2Id=0)
double GetSigShapeSystematic(std::string parName)
void computeProdXS_VH(float &mevalue, bool includeHiggsDecay=false)
RooAbsReal * gzgs2List[1][2]
TVar::LeptonInterference myLepInterf_
double selfDHzzLambda_qsq[nSupportedHiggses][SIZE_HVV_LAMBDAQSQ][SIZE_HVV_CQSQ]
void computeProdXS_JJH(float &mevalue)
void computeProdDecP(double selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2], double selfDHwwcoupl_input[nSupportedHiggses][SIZE_HVV][2], double selfDaTQGCcoupl_input[SIZE_ATQGC][2], double selfDAZffcoupl_input[SIZE_AZff][2], float &prob, bool useConstant=true)
computes the combined production and decay probability while taking in coupling arrays
void reset_SelfDCouplings()
virtual void resetHypotheses()
void computeAngles(float &costhetastar, float &costheta1, float &costheta2, float &Phi, float &Phi1, TLorentzVector Z1_lept1, int Z1_lept1Id, TLorentzVector Z1_lept2, int Z1_lept2Id, TLorentzVector Z2_lept1, int Z2_lept1Id, TLorentzVector Z2_lept2, int Z2_lept2Id)
double selfDHwwcoupl[nSupportedHiggses][SIZE_HVV][2]
void computeProdXS_VVHVV(float &mevalue)
TVar::MatrixElement myME_
double GetSigShapeParameter(std::string parName)
double selfDGvvpcoupl[SIZE_GVV][2]
void constructDggr(float bkg_VAMCFM_noscale, float ggzz_VAMCFM_noscale, float ggHZZ_prob_pure_noscale, float ggHZZ_prob_int_noscale, float widthScale, float &myDggr)
void set_RenFacScaleMode(TVar::EventScaleScheme renormalizationSch, TVar::EventScaleScheme factorizationSch, double ren_sf, double fac_sf)
double selfDAZffcoupl[SIZE_AZff][2]
MelaPConstant * pAvgSmooth_MCFM_JJQCD_bkgZJets_2l2q
int selfDHzzCLambda_qsq[nSupportedHiggses][SIZE_HVV_CQSQ]
double get_PrimaryHiggsMass()
void setSpinOneCouplings()
MELAParticle * getAssociatedJet(int index) const
void computeXS(float &mevalue)
double selfDZvvcoupl[SIZE_ZVV][2]
MelaPConstant * pAvgSmooth_MCFM_Had_WH_S_HSMHiggs_4mu
MelaPConstant * pAvgSmooth_MCFM_JJQCD_bkgZZ_2mu2e
MelaPConstant * pAvgSmooth_MCFM_ZZGG_bkgZZ_4mu
virtual void makeParamsConst(bool yesNo)
MelaPConstant * pAvgSmooth_MCFM_Had_ZH_S_HSMHiggs_2mu2e
SuperDijetMela * superDijet
MELAOutputStreamer MELAerr
float getConstant_FourFermionDecay(const int &decid)
TRandom3 melaRandomNumber
void computeProdP_ttH(float &prob, int topProcess=2, int topDecay=0, bool useConstant=true)
std::complex< double > GetMadgraphCKMElement(int iquark, int jquark)
void writeCentered(const T &val, char fillch=' ', std::streamsize gapsize=0)
double selfDHt4t4coupl[nSupportedHiggses][SIZE_HQQ][2]
double selfDHwwLambda_qsq[nSupportedHiggses][SIZE_HVV_LAMBDAQSQ][SIZE_HVV_CQSQ]
MelaPConstant * pAvgSmooth_MCFM_Had_ZH_bkgZZ_2mu2e
void setMelaPrimaryHiggsMass(double myHiggsMass)
Sets the mass of the "primary" higgs.
void computeFakeJet(TLorentzVector const &realJet, TLorentzVector const &others, TLorentzVector &fakeJet)
struct madMela::@2 params_r_
void setMEConst(const double &val)
double selfDSMEFTSimcoupl[SIZE_SMEFT]
MelaPConstant * pAvgSmooth_JHUGen_Had_ZH_HSMHiggs[TVar::nFermionMassRemovalSchemes-1]
complex(8), parameter, public c1
RooRealVar * costhetastar_rrv
std::vector< MELAParticle * > & getAssociatedJets()
RooRealVar * costheta2_rrv
double selfDHwwpcoupl[SIZE_HVV][2]
RooAbsReal * Lambda_z3qsq[SIZE_HVV_CQSQ]
RooAbsReal * gzgs1List[1][2]
bool isAGluon(const int id)
bool configureAnalyticalPDFs()
RooAbsReal * g3List[8][2]
double selfDHzzcoupl[nSupportedHiggses][SIZE_HVV][2]
std::vector< TLorentzVector > calculate4Momentum(double Mx, double M1, double M2, double theta, double theta1, double theta2, double Phi1, double Phi)
void set_CurrentCandidate(MELACandidate *cand)
void getPConstantHandles()
RooAbsReal * g1List[8][2]
void setSpinTwoCouplings()
TVar::FermionMassRemoval JetMassScheme
RooSpinZero::modelCouplings couplings
@ kUseAssociated_UnstableTops
RooAbsReal * ggsgs3List[1][2]
VectorPdfFactory * spin1Model
MelaPConstant * pAvgSmooth_MCFM_ZZGG_HSMHiggs_2mu2e
double selfDHzzpcoupl[SIZE_HVV][2]
void set_LeptonInterference(TVar::LeptonInterference myLepInterf)
void set_aTQGCCouplings(double selfDaTQGCcoupl[SIZE_ATQGC][2])
std::vector< SimpleParticle_t > SimpleParticleCollection_t
MelaPConstant * pAvgSmooth_MCFM_ZZQQB_bkgZZ_2mu2e
void SetMadgraphCKMElements(double ckmlambda=0.2265, double ckma=0.79, double ckmrho=0.141, double ckmeta=0)
double selfDHggcoupl[nSupportedHiggses][SIZE_HGG][2]
RooAbsReal * Lambda_z2qsq[SIZE_HVV_CQSQ]
MelaPConstant * pAvgSmooth_JHUGen_Had_WH_HSMHiggs[TVar::nFermionMassRemovalSchemes-1]
void deletePConstantHandle(MelaPConstant *&handle)
void reset_MCFM_EWKParameters(double ext_Gf, double ext_aemmz, double ext_mW, double ext_mZ, double ext_xW, int ext_ewscheme=3)
void setMelaLeptonInterference(TVar::LeptonInterference myLepInterf=TVar::DefaultLeptonInterf)
Sets the MELA Lepton Interference.
RooAbsReal * g4List[8][2]
void computeProdXS_ttH(float &mevalue, int topProcess, int topDecay=0)
MelaPConstant * pAvgSmooth_MCFM_JJVBF_S_HSMHiggs_4e
void set_SpinTwoContact(double selfDGvvpcoupl[SIZE_GVV][2], double selfDGvpvpcoupl[SIZE_GVV][2])
double get_PrimaryWidth(int ipart)
RooAbsReal * cLambda_qsq[SIZE_HVV_CQSQ]
void SetDecayChannel(std::string myChan)
double selfDHbbcoupl[nSupportedHiggses][SIZE_HQQ][2]
void computeVBFAngles(float &costhetastar, float &costheta1, float &costheta2, float &Phi, float &Phi1, float &Q2V1, float &Q2V2, TLorentzVector p4M11, int Z1_lept1Id, TLorentzVector p4M12, int Z1_lept2Id, TLorentzVector p4M21, int Z2_lept1Id, TLorentzVector p4M22, int Z2_lept2Id, TLorentzVector jet1, int jet1Id, TLorentzVector jet2, int jet2Id, TLorentzVector *injet1=0, int injet1Id=0, TLorentzVector *injet2=0, int injet2Id=0)
MelaPConstant * pAvgSmooth_MCFM_Had_WH_bkgZZ_4mu
void GetBoostedParticleVectors(MELACandidate *melaCand, TVar::simple_event_record &mela_event, TVar::VerbosityLevel verbosity=TVar::DEBUG)
RooAbsReal * bList[SIZE_GVV][2]
double selfDGggcoupl[SIZE_GGG][2]
MelaPConstant * pAvgSmooth_JHUGen_ZZGG_HSMHiggs_4e
RooAbsReal * ggsgs2List[1][2]
MELACandidate * getCurrentCandidate()
Gets the current MELA top-level (input) candList object.
bool isAJet(const int id)
void computeVBFAngles_ComplexBoost(float &costhetastar, float &costheta1_real, float &costheta1_imag, float &costheta2_real, float &costheta2_imag, float &Phi, float &Phi1, float &Q2V1, float &Q2V2, TLorentzVector p4M11, int Z1_lept1Id, TLorentzVector p4M12, int Z1_lept2Id, TLorentzVector p4M21, int Z2_lept1Id, TLorentzVector p4M22, int Z2_lept2Id, TLorentzVector jet1, int jet1Id, TLorentzVector jet2, int jet2Id, TLorentzVector *injet1=0, int injet1Id=0, TLorentzVector *injet2=0, int injet2Id=0)
complex(8), parameter, public c2
void set_PrimaryHiggsMass(double mh)
void set_wHiggs(double gah_, int index)
float GetConvBW(TVar::Production prod, MELACandidate *cand, bool useTrueBW)
MelaPConstant * pAvgSmooth_MCFM_Had_ZH_S_HSMHiggs_4e
void applyJetMassCorrection(bool flag=true)