JHUGen MELA
v2.4.1
Matrix element calculations as used in JHUGen. MELA is an important tool that was used for the Higgs boson discovery and for precise measurements of its structure and interactions. Please see the website https://spin.pha.jhu.edu/ and papers cited there for more details, and kindly cite those papers when using this code.
|
Go to the documentation of this file.
38 #include "TLorentzVector.h"
62 void scaleMomentumToEnergy(
const TLorentzVector& massiveJet, TLorentzVector& masslessJet,
double mass=0);
65 TLorentzVector
const& jet1,
int const& jet1Id,
66 TLorentzVector
const& jet2,
int const& jet2Id,
67 double m1=0,
double m2=0
72 void computeFakeJet(TLorentzVector
const& realJet, TLorentzVector
const& others, TLorentzVector& fakeJet);
75 std::pair<TLorentzVector, TLorentzVector>
ComplexBoost(TVector3
const& beta, TLorentzVector
const& p4);
88 TLorentzVector Z1_lept1,
int Z1_lept1Id,
89 TLorentzVector Z1_lept2,
int Z1_lept2Id,
90 TLorentzVector Z2_lept1,
int Z2_lept1Id,
91 TLorentzVector Z2_lept2,
int Z2_lept2Id
100 TLorentzVector Z1_lept1,
int Z1_lept1Id,
101 TLorentzVector Z1_lept2,
int Z1_lept2Id,
102 TLorentzVector Z2_lept1,
int Z2_lept1Id,
103 TLorentzVector Z2_lept2,
int Z2_lept2Id
114 TLorentzVector p4M11,
int Z1_lept1Id,
115 TLorentzVector p4M12,
int Z1_lept2Id,
116 TLorentzVector p4M21,
int Z2_lept1Id,
117 TLorentzVector p4M22,
int Z2_lept2Id,
118 TLorentzVector jet1,
int jet1Id,
119 TLorentzVector jet2,
int jet2Id,
120 TLorentzVector* injet1=0,
int injet1Id=0,
121 TLorentzVector* injet2=0,
int injet2Id=0
125 float& costheta1_real,
float& costheta1_imag,
126 float& costheta2_real,
float& costheta2_imag,
131 TLorentzVector p4M11,
int Z1_lept1Id,
132 TLorentzVector p4M12,
int Z1_lept2Id,
133 TLorentzVector p4M21,
int Z2_lept1Id,
134 TLorentzVector p4M22,
int Z2_lept2Id,
135 TLorentzVector jet1,
int jet1Id,
136 TLorentzVector jet2,
int jet2Id,
137 TLorentzVector* injet1=0,
int injet1Id=0,
138 TLorentzVector* injet2=0,
int injet2Id=0
148 TLorentzVector p4M11,
int Z1_lept1Id,
149 TLorentzVector p4M12,
int Z1_lept2Id,
150 TLorentzVector p4M21,
int Z2_lept1Id,
151 TLorentzVector p4M22,
int Z2_lept2Id,
152 TLorentzVector jet1,
int jet1Id,
153 TLorentzVector jet2,
int jet2Id,
154 TLorentzVector* injet1=0,
int injet1Id=0,
155 TLorentzVector* injet2=0,
int injet2Id=0
183 TLorentzVector p4M11,
int Z1_lept1Id,
184 TLorentzVector p4M12,
int Z1_lept2Id,
186 TLorentzVector p4M21,
int Z2_lept1Id,
187 TLorentzVector p4M22,
int Z2_lept2Id,
189 TLorentzVector b,
int bId,
190 TLorentzVector Wplusf,
int WplusfId,
191 TLorentzVector Wplusfb,
int WplusfbId,
192 TLorentzVector bbar,
int bbarId,
193 TLorentzVector Wminusf,
int WminusfId,
194 TLorentzVector Wminusfb,
int WminusfbId,
195 TLorentzVector* injet1=0,
int injet1Id=0,
196 TLorentzVector* injet2=0,
int injet2Id=0
200 void SetEwkCouplingParameters(
double ext_Gf,
double ext_aemmz,
double ext_mW,
double ext_mZ,
double ext_xW,
int ext_ewscheme);
201 void SetMass(
double inmass,
int ipart);
209 void SetCKMElements(
double* invckm_ud,
double* invckm_us,
double* invckm_cd,
double* invckm_cs,
double* invckm_ts,
double* invckm_tb,
double* invckm_ub=0,
double* invckm_cb=0,
double* invckm_td=0);
212 void SetAlphaS(
double& Q_ren,
double& Q_fac,
double multiplier_ren,
double multiplier_fac,
int mynloop,
int mynflav, std::string mypartons);
225 std::vector<int>* partOrder, std::vector<int>* apartOrder
230 void InitJHUGenMELA(
const char* pathtoPDFSet,
int PDFMember,
double collider_sqrts);
278 bool includeHiggsDecay,
285 int topDecay,
int topProcess,
320 std::vector<MELAParticle*>* particleList,
321 std::vector<MELACandidate*>* candList
329 std::vector<MELAParticle*>* particleList,
330 std::vector<MELAThreeBodyDecayCandidate*>* tbdCandList
void scaleMomentumToEnergy(const TLorentzVector &massiveJet, TLorentzVector &masslessJet, double mass=0)
bool MCFM_chooser(const TVar::Process &process, const TVar::Production &production, const TVar::LeptonInterference &leptonInterf, const TVar::VerbosityLevel &verbosity, const TVar::simple_event_record &mela_event)
double ResonancePropagator(double const &sqrts, TVar::ResonancePropagatorScheme scheme)
void SetJHUGenSpinZeroGGCouplings(double Hggcoupl[SIZE_HGG][2])
int WipeMEArray(const TVar::Process &process, const TVar::Production &production, const int id[mxpart], double msq[nmsq][nmsq], const TVar::VerbosityLevel &verbosity)
bool MCFM_smalls(double s[][mxpart], int npart)
double SumMEPDF(const TLorentzVector &p0, const TLorentzVector &p1, double msq[nmsq][nmsq], MelaIO *RcdME, const double &EBEAM, const TVar::VerbosityLevel &verbosity)
double BBHiggsMatEl(const TVar::Process &process, const TVar::Production &production, const TVar::MatrixElement &matrixElement, TVar::event_scales_type *event_scales, MelaIO *RcdME, const double &EBEAM, int botProcess, TVar::VerbosityLevel verbosity)
double InterpretScaleScheme(const TVar::Production &production, const TVar::MatrixElement &matrixElement, const TVar::EventScaleScheme &scheme, TLorentzVector p[mxpart])
void SetJHUGenHiggsMassWidth(double MReso, double GaReso)
double HJJMatEl(const TVar::Process &process, const TVar::Production &production, const TVar::MatrixElement &matrixElement, TVar::event_scales_type *event_scales, MelaIO *RcdME, const double &EBEAM, TVar::VerbosityLevel verbosity)
void SetCKMElements(double *invckm_ud, double *invckm_us, double *invckm_cd, double *invckm_cs, double *invckm_ts, double *invckm_tb, double *invckm_ub=0, double *invckm_cb=0, double *invckm_td=0)
void SetMCFMAZffCouplings(bool useBSM, AZffCouplings const *Zcouplings)
void SetAlphaS(double &Q_ren, double &Q_fac, double multiplier_ren, double multiplier_fac, int mynloop, int mynflav, std::string mypartons)
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)
void applyLeptonMassCorrection(bool flag=true)
void SetMass(double inmass, int ipart)
void SetJHUGenVprimeContactCouplings(double Zpffcoupl[SIZE_Vpff][2], double Wpffcoupl[SIZE_Vpff][2])
void PrintCandidateSummary(MELACandidate *cand)
void SetMCFMSpinZeroCouplings(bool useBSM, SpinZeroCouplings const *Hcouplings, bool forceZZ)
void SetJHUGenSpinOneCouplings(double Zqqcoupl[SIZE_ZQQ][2], double Zvvcoupl[SIZE_ZVV][2])
bool forbidMassiveLeptons
Remove fermion mass if the flag is set to true.
void GetAlphaS(double *alphas_, double *alphasmz_)
void SetDecayWidth(double inwidth, int ipart)
bool MCFM_masscuts(double s[][mxpart], const TVar::Process &process)
void setLeptonMassScheme(TVar::FermionMassRemoval scheme=TVar::ConserveDifermionMass)
void SetJHUGenSpinTwoCouplings(double Gacoupl[SIZE_GGG][2], double Gvvcoupl[SIZE_GVV][2], double Gvvpcoupl[SIZE_GVV][2], double Gvpvpcoupl[SIZE_GVV][2], double qLeftRightcoupl[SIZE_GQQ][2])
MELACandidate * ConvertVectorFormat(SimpleParticleCollection_t *pDaughters, SimpleParticleCollection_t *pAssociated, SimpleParticleCollection_t *pMothers, bool isGen, std::vector< MELAParticle * > *particleList, std::vector< MELACandidate * > *candList)
void SetJHUGenSpinZeroVVCouplings(double Hvvcoupl[SIZE_HVV][2], double Hvvpcoupl[SIZE_HVV][2], double Hvpvpcoupl[SIZE_HVV][2], int Hvvcoupl_cqsq[SIZE_HVV_CQSQ], double HvvLambda_qsq[SIZE_HVV_LAMBDAQSQ][SIZE_HVV_CQSQ], bool useWWcoupl)
void SetJHUGenDistinguishWWCouplings(bool doAllow)
void SetJHUGenAZffCouplings(bool needAZff, double AZffcoupl[SIZE_AZff][2])
double GetDecayWidth(int ipart)
double GetCKMElement(int iquark, int jquark)
TString GetMCFMParticleLabel(const int &pid, bool useQJ, bool useExtendedConventions)
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)
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)
std::pair< TLorentzVector, TLorentzVector > ComplexBoost(TVector3 const &beta, TLorentzVector const &p4)
void SetEwkCouplingParameters(double ext_Gf, double ext_aemmz, double ext_mW, double ext_mZ, double ext_xW, int ext_ewscheme)
void computeAnglesCS(float const &pbeam, 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)
std::pair< TLorentzVector, TLorentzVector > removeMassFromPair(TLorentzVector const &jet1, int const &jet1Id, TLorentzVector const &jet2, int const &jet2Id, double m1=0, double m2=0)
TVar::FermionMassRemoval LeptonMassScheme
void GetMassWidth(int ipart, double &m, double &ga)
MELAThreeBodyDecayCandidate * ConvertThreeBodyDecayCandidate(SimpleParticleCollection_t *tbdDaughters, std::vector< MELAParticle * > *particleList, std::vector< MELAThreeBodyDecayCandidate * > *tbdCandList)
void ComputePDF(const TLorentzVector &p0, const TLorentzVector &p1, double fx1[nmsq], double fx2[nmsq], const double &EBEAM, const TVar::VerbosityLevel &verbosity)
double GetMass(int ipart)
ResonancePropagatorScheme
double alphas_(double *q, double *amz, int *nloop)
bool MCFM_SetupParticleCouplings(const TVar::Process &process, const TVar::Production &production, const TVar::VerbosityLevel &verbosity, const TVar::simple_event_record &mela_event, std::vector< int > *partOrder, std::vector< int > *apartOrder)
void SetMCFMaTQGCCouplings(bool useBSM, aTQGCCouplings const *couplings)
void computeFakeJet(TLorentzVector const &realJet, TLorentzVector const &others, TLorentzVector &fakeJet)
double TTHiggsMatEl(const TVar::Process &process, const TVar::Production &production, const TVar::MatrixElement &matrixElement, TVar::event_scales_type *event_scales, MelaIO *RcdME, const double &EBEAM, int topDecay, int topProcess, TVar::VerbosityLevel verbosity)
void constrainedRemovePairMass(TLorentzVector &p1, TLorentzVector &p2, double m1=0, double m2=0)
void InitJHUGenMELA(const char *pathtoPDFSet, int PDFMember, double collider_sqrts)
double SumMatrixElementPDF(const TVar::Process &process, const TVar::Production &production, const TVar::MatrixElement &matrixElement, const TVar::LeptonInterference &leptonInterf, TVar::event_scales_type *event_scales, MelaIO *RcdME, const double &EBEAM, TVar::VerbosityLevel verbosity)
TVar::FermionMassRemoval JetMassScheme
void adjustTopDaughters(SimpleParticleCollection_t &daughters)
bool CheckPartonMomFraction(const TLorentzVector &p0, const TLorentzVector &p1, double xx[2], const double &EBEAM, const TVar::VerbosityLevel &verbosity)
void SetJHUGenSpinZeroQQCouplings(double Hqqcoupl[SIZE_HQQ][2])
std::vector< SimpleParticle_t > SimpleParticleCollection_t
void setJetMassScheme(TVar::FermionMassRemoval scheme=TVar::ConserveDifermionMass)
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)
void GetBoostedParticleVectors(MELACandidate *melaCand, TVar::simple_event_record &mela_event, TVar::VerbosityLevel verbosity=TVar::DEBUG)
double JHUGenMatEl(const TVar::Process &process, const TVar::Production &production, const TVar::MatrixElement &matrixElement, TVar::event_scales_type *event_scales, MelaIO *RcdME, const double &EBEAM, TVar::VerbosityLevel verbosity)
double VHiggsMatEl(const TVar::Process &process, const TVar::Production &production, const TVar::MatrixElement &matrixElement, TVar::event_scales_type *event_scales, MelaIO *RcdME, const double &EBEAM, bool includeHiggsDecay, TVar::VerbosityLevel verbosity)
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)
void ResetAmplitudeIncludes()
void applyJetMassCorrection(bool flag=true)