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.
TUtil.hh
Go to the documentation of this file.
1 // March 28 2011
2 // S. Jindariani ([email protected])
3 // Y. Gao ([email protected])
4 // K. Burkett ([email protected])
5 
6 
7 #ifndef ZZ_COMMON
8 #define ZZ_COMMON
9 #include <string>
10 #include <vector>
11 // MelaIO class
12 #include "MelaIO.h"
13 // Couplings classes
14 #include "TCouplings.hh"
15 // MCFM utilities
16 #include "TMCFMUtils.hh"
17 // Mod_Parameters
18 #include "TModParameters.hh"
19 // NNPDF Driver for JHUGen
20 #include "TNNPDFDriver.hh"
21 // Mod Kinematics
22 #include "TModKinematics.hh"
23 // JHUGenMELA
24 #include "TModJHUGen.hh"
25 #include "TModJHUGenMELA.hh"
26 // Higgs + 0 jet
27 #include "TModHiggsMatEl.hh"
28 #include "TModGravitonMatEl.hh"
29 #include "TModZprimeMatEl.hh"
30 // Higgs + 1/2 jets
31 #include "TModHiggsJJMatEl.hh"
32 #include "TModHiggsJMatEl.hh"
33 // VH
34 #include "TModVHiggsMatEl.hh"
35 // ttH
36 #include "TModTTBHMatEl.hh"
37 // ROOT includes
38 #include "TLorentzVector.h"
39 #include "TTree.h"
40 #include "TH1F.h"
41 #include "TH2F.h"
42 #include "TH1D.h"
43 #include "TH2D.h"
44 #include "TFile.h"
45 #include "TF1.h"
46 
47 
48 namespace TUtil{
50  extern bool forbidMassiveLeptons;
51  extern bool forbidMassiveJets;
54 
55  void applyLeptonMassCorrection(bool flag=true);
56  void applyJetMassCorrection(bool flag=true);
59  // This version makes the masses of p1 and p2 to be m1 and m2, leaving p1+p2 unchanged.
60  void constrainedRemovePairMass(TLorentzVector& p1, TLorentzVector& p2, double m1=0, double m2=0);
61  // This version simply scales momentum to match energy for the desired mass.
62  void scaleMomentumToEnergy(const TLorentzVector& massiveJet, TLorentzVector& masslessJet, double mass=0);
63  // Function that has generic removal features
64  std::pair<TLorentzVector, TLorentzVector> removeMassFromPair(
65  TLorentzVector const& jet1, int const& jet1Id,
66  TLorentzVector const& jet2, int const& jet2Id,
67  double m1=0, double m2=0
68  );
69  // Function that adjusts top daughter kinematics
70  void adjustTopDaughters(SimpleParticleCollection_t& daughters); // Daughters are arranged as b, Wf, Wfb
71  // Compute a fake jet from the massless jets
72  void computeFakeJet(TLorentzVector const& realJet, TLorentzVector const& others, TLorentzVector& fakeJet); // Input massive + higgs -> output massless fake jet
73 
74  // TLorentzVector::Boost in complex plane
75  std::pair<TLorentzVector, TLorentzVector> ComplexBoost(TVector3 const& beta, TLorentzVector const& p4);
76 
82  void computeAngles(
83  float& costhetastar,
84  float& costheta1,
85  float& costheta2,
86  float& Phi,
87  float& Phi1,
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
92  );
93  void computeAnglesCS(
94  float const& pbeam,
95  float& costhetastar,
96  float& costheta1,
97  float& costheta2,
98  float& Phi,
99  float& Phi1,
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
104  );
105  // Angles of associated production
106  void computeVBFAngles(
107  float& costhetastar,
108  float& costheta1,
109  float& costheta2,
110  float& Phi,
111  float& Phi1,
112  float& Q2V1,
113  float& Q2V2,
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, // Gen. partons in lab frame
121  TLorentzVector* injet2=0, int injet2Id=0
122  );
124  float& costhetastar,
125  float& costheta1_real, float& costheta1_imag,
126  float& costheta2_real, float& costheta2_imag,
127  float& Phi,
128  float& Phi1,
129  float& Q2V1,
130  float& Q2V2,
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, // Gen. partons in lab frame
138  TLorentzVector* injet2=0, int injet2Id=0
139  );
140  void computeVHAngles(
141  float& costhetastar,
142  float& costheta1,
143  float& costheta2,
144  float& Phi,
145  float& Phi1,
146  float& m1,
147  float& m2,
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, // Gen. partons in lab frame
155  TLorentzVector* injet2=0, int injet2Id=0
156  );
157  // hs: cos theta*
158  // h[x]: cos theta_X
159  void computeTTHAngles(
160  // ttH system
161  float& hs,
162  float& hincoming,
163  float& hTT,
164  float& PhiTT,
165  float& Phi1,
166 
167  // tt system
168  float& hbb,
169  float& hWW,
170  float& Phibb,
171  float& Phi1bb,
172 
173  // Wplus system
174  float& hWplusf,
175  float& PhiWplusf,
176 
177  // Wminus system
178  float& hWminusf,
179  float& PhiWminusf,
180 
181  // H->ZZ/ZA/AA vectors. For A decays, leave second vector (0,0,0,0) with id=-9000
182  // Z/A1
183  TLorentzVector p4M11, int Z1_lept1Id,
184  TLorentzVector p4M12, int Z1_lept2Id,
185  // Z/A2
186  TLorentzVector p4M21, int Z2_lept1Id,
187  TLorentzVector p4M22, int Z2_lept2Id,
188  // TT system
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, // Gen. partons in lab frame
196  TLorentzVector* injet2=0, int injet2Id=0
197  );
198 
199  // Parameter settings
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);
202  void SetDecayWidth(double inwidth, int ipart);
203  double GetMass(int ipart);
204  double GetDecayWidth(int ipart);
205  double GetMass(const MELAParticle* part);
206  double GetDecayWidth(const MELAParticle* part);
207  void GetMassWidth(int ipart, double& m, double& ga);
208  void GetMassWidth(const MELAParticle* part, double& m, double& ga);
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);
210  double GetCKMElement(int iquark, int jquark);
211  double InterpretScaleScheme(const TVar::Production& production, const TVar::MatrixElement& matrixElement, const TVar::EventScaleScheme& scheme, TLorentzVector p[mxpart]);
212  void SetAlphaS(double& Q_ren, double& Q_fac, double multiplier_ren, double multiplier_fac, int mynloop, int mynflav, std::string mypartons); // Q_ren/fac -> Q_ren/fac * multiplier_ren/fac
213  void GetAlphaS(double* alphas_, double* alphasmz_); // Get last alpha_s value set
214 
215  // chooser.f split into 3 different functions
216  bool MCFM_chooser(
217  const TVar::Process& process, const TVar::Production& production, const TVar::LeptonInterference& leptonInterf,
218  const TVar::VerbosityLevel& verbosity,
219  const TVar::simple_event_record& mela_event
220  );
223  const TVar::VerbosityLevel& verbosity,
224  const TVar::simple_event_record& mela_event,
225  std::vector<int>* partOrder, std::vector<int>* apartOrder
226  );
227  TString GetMCFMParticleLabel(const int& pid, bool useQJ, bool useExtendedConventions);
228 
229  // JHUGen-specific wrappers
230  void InitJHUGenMELA(const char* pathtoPDFSet, int PDFMember, double collider_sqrts);
231  void SetJHUGenHiggsMassWidth(double MReso, double GaReso);
232  void SetJHUGenDistinguishWWCouplings(bool doAllow);
233  void ResetAmplitudeIncludes();
234 
235  // Spin-0 couplings
236  void SetMCFMSpinZeroCouplings(bool useBSM, SpinZeroCouplings const* Hcouplings, bool forceZZ);
237  void SetMCFMaTQGCCouplings(bool useBSM, aTQGCCouplings const* couplings);
238  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);
239  void SetJHUGenSpinZeroGGCouplings(double Hggcoupl[SIZE_HGG][2]);
240  void SetJHUGenSpinZeroQQCouplings(double Hqqcoupl[SIZE_HQQ][2]);
241  // Spin-1 couplings
242  void SetJHUGenSpinOneCouplings(double Zqqcoupl[SIZE_ZQQ][2], double Zvvcoupl[SIZE_ZVV][2]);
243  // Spin-2 couplings
244  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]);
245  // Vprime / contact couplings
246  void SetJHUGenVprimeContactCouplings(double Zpffcoupl[SIZE_Vpff][2], double Wpffcoupl[SIZE_Vpff][2]);
247  // AZff
248  void SetMCFMAZffCouplings(bool useBSM, AZffCouplings const* Zcouplings);
249  void SetJHUGenAZffCouplings(bool needAZff, double AZffcoupl[SIZE_AZff][2]);
250 
251  // PS cuts, unused
252  bool MCFM_masscuts(double s[][mxpart], const TVar::Process& process);
253  bool MCFM_smalls(double s[][mxpart], int npart);
254 
255  // ME computations
256  double SumMatrixElementPDF(
257  const TVar::Process& process, const TVar::Production& production, const TVar::MatrixElement& matrixElement, const TVar::LeptonInterference& leptonInterf,
258  TVar::event_scales_type* event_scales, MelaIO* RcdME,
259  const double& EBEAM,
260  TVar::VerbosityLevel verbosity
261  );
262  double JHUGenMatEl(
263  const TVar::Process& process, const TVar::Production& production, const TVar::MatrixElement& matrixElement,
264  TVar::event_scales_type* event_scales, MelaIO* RcdME,
265  const double& EBEAM,
266  TVar::VerbosityLevel verbosity
267  );
268  double HJJMatEl(
269  const TVar::Process& process, const TVar::Production& production, const TVar::MatrixElement& matrixElement,
270  TVar::event_scales_type* event_scales, MelaIO* RcdME,
271  const double& EBEAM,
272  TVar::VerbosityLevel verbosity
273  );
274  double VHiggsMatEl(
275  const TVar::Process& process, const TVar::Production& production, const TVar::MatrixElement& matrixElement,
276  TVar::event_scales_type* event_scales, MelaIO* RcdME,
277  const double& EBEAM,
278  bool includeHiggsDecay,
279  TVar::VerbosityLevel verbosity
280  );
281  double TTHiggsMatEl(
282  const TVar::Process& process, const TVar::Production& production, const TVar::MatrixElement& matrixElement,
283  TVar::event_scales_type* event_scales, MelaIO* RcdME,
284  const double& EBEAM,
285  int topDecay, int topProcess,
286  TVar::VerbosityLevel verbosity
287  );
288  double BBHiggsMatEl(
289  const TVar::Process& process, const TVar::Production& production, const TVar::MatrixElement& matrixElement,
290  TVar::event_scales_type* event_scales, MelaIO* RcdME,
291  const double& EBEAM,
292  int botProcess,
293  TVar::VerbosityLevel verbosity
294  );
295 
296  int WipeMEArray(const TVar::Process& process, const TVar::Production& production, const int id[mxpart], double msq[nmsq][nmsq], const TVar::VerbosityLevel& verbosity);
297  bool CheckPartonMomFraction(const TLorentzVector& p0, const TLorentzVector& p1, double xx[2], const double& EBEAM, const TVar::VerbosityLevel& verbosity);
298  void ComputePDF(const TLorentzVector& p0, const TLorentzVector& p1, double fx1[nmsq], double fx2[nmsq], const double& EBEAM, const TVar::VerbosityLevel& verbosity);
299  double SumMEPDF(const TLorentzVector& p0, const TLorentzVector& p1, double msq[nmsq][nmsq], MelaIO* RcdME, const double& EBEAM, const TVar::VerbosityLevel& verbosity);
300 
301  // Propagator reweighting
302  double ResonancePropagator(double const& sqrts, TVar::ResonancePropagatorScheme scheme);
303 
304  // Boost the particles with or without associated ones to pT=0 frame and return std::vectors filled with (id, momentum) pairs
306  MELACandidate* melaCand,
307  TVar::simple_event_record& mela_event,
309  );
310 
311  // Convert vectors of simple particles to MELAParticles and create a MELACandidate
312  // The output lists could be members of TEvtProb directly.
314  // Inputs
315  SimpleParticleCollection_t* pDaughters, // Cannot be 0
316  SimpleParticleCollection_t* pAssociated, // Allowed to be 0
317  SimpleParticleCollection_t* pMothers, // Allowed to be 0
318  bool isGen,
319  // Outputs
320  std::vector<MELAParticle*>* particleList,
321  std::vector<MELACandidate*>* candList
322  );
323  // Convert the vector of three body decay daughters (as simple particles) to MELAParticles and create a MELAThreeBodyDecayCandidate
324  // The output lists could be members of TEvtProb directly.
326  // Input
327  SimpleParticleCollection_t* tbdDaughters,
328  // Outputs
329  std::vector<MELAParticle*>* particleList,
330  std::vector<MELAThreeBodyDecayCandidate*>* tbdCandList
331  );
334 
335 }
336 
337 #endif
338 
TUtil::scaleMomentumToEnergy
void scaleMomentumToEnergy(const TLorentzVector &massiveJet, TLorentzVector &masslessJet, double mass=0)
Definition: TUtil.cc:66
TUtil::MCFM_chooser
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)
Definition: TUtil.cc:1802
TUtil::ResonancePropagator
double ResonancePropagator(double const &sqrts, TVar::ResonancePropagatorScheme scheme)
Definition: TUtil.cc:7929
aTQGCCouplings
Definition: TCouplings.hh:128
TUtil::SetJHUGenSpinZeroGGCouplings
void SetJHUGenSpinZeroGGCouplings(double Hggcoupl[SIZE_HGG][2])
Definition: TUtil.cc:4084
TUtil::WipeMEArray
int WipeMEArray(const TVar::Process &process, const TVar::Production &production, const int id[mxpart], double msq[nmsq][nmsq], const TVar::VerbosityLevel &verbosity)
Definition: TUtil.cc:7729
mela.daughters
string daughters
Definition: mela.py:767
TVar::FermionMassRemoval
FermionMassRemoval
Definition: TVar.hh:111
MelaIO
Definition: MelaIO.h:8
TVar::LeptonInterference
LeptonInterference
Definition: TVar.hh:106
SIZE_GVV
@ SIZE_GVV
Definition: raw_couplings.txt:158
TUtil::MCFM_smalls
bool MCFM_smalls(double s[][mxpart], int npart)
Definition: TUtil.cc:3298
TVar::VerbosityLevel
VerbosityLevel
Definition: TVar.hh:47
TUtil::SumMEPDF
double SumMEPDF(const TLorentzVector &p0, const TLorentzVector &p1, double msq[nmsq][nmsq], MelaIO *RcdME, const double &EBEAM, const TVar::VerbosityLevel &verbosity)
Definition: TUtil.cc:7911
TModParameters.hh
TUtil::BBHiggsMatEl
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)
Definition: TUtil.cc:7562
TModZprimeMatEl.hh
SIZE_HVV_CQSQ
@ SIZE_HVV_CQSQ
Definition: raw_couplings.txt:74
TModKinematics.hh
TUtil::InterpretScaleScheme
double InterpretScaleScheme(const TVar::Production &production, const TVar::MatrixElement &matrixElement, const TVar::EventScaleScheme &scheme, TLorentzVector p[mxpart])
Definition: TUtil.cc:1601
MelaIO.h
TUtil::SetJHUGenHiggsMassWidth
void SetJHUGenHiggsMassWidth(double MReso, double GaReso)
Definition: TUtil.cc:3347
TUtil::HJJMatEl
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)
Definition: TUtil.cc:5118
TUtil::SetCKMElements
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)
Definition: TUtil.cc:1507
TUtil::SetMCFMAZffCouplings
void SetMCFMAZffCouplings(bool useBSM, AZffCouplings const *Zcouplings)
Definition: TUtil.cc:4044
TModJHUGenMELA.hh
TVar::event_scales_type
Definition: TVar.hh:260
SIZE_GGG
@ SIZE_GGG
Definition: raw_couplings.txt:131
TModHiggsJJMatEl.hh
AZffCouplings
Definition: TCouplings.hh:143
TNNPDFDriver.hh
TUtil::SetAlphaS
void SetAlphaS(double &Q_ren, double &Q_fac, double multiplier_ren, double multiplier_fac, int mynloop, int mynflav, std::string mypartons)
testME_all.production
def production
Definition: testME_all.py:115
TUtil::computeTTHAngles
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)
Definition: TUtil.cc:989
SIZE_ZQQ
@ SIZE_ZQQ
Definition: raw_couplings.txt:107
TUtil::applyLeptonMassCorrection
void applyLeptonMassCorrection(bool flag=true)
Definition: TUtil.cc:35
TUtil::SetMass
void SetMass(double inmass, int ipart)
Definition: TUtil.cc:1445
TVar::Process
Process
Definition: TVar.hh:125
TVar::EventScaleScheme
EventScaleScheme
Definition: TVar.hh:196
TUtil::SetJHUGenVprimeContactCouplings
void SetJHUGenVprimeContactCouplings(double Zpffcoupl[SIZE_Vpff][2], double Wpffcoupl[SIZE_Vpff][2])
Definition: TUtil.cc:4090
testME_all.p
p
Definition: testME_all.py:11
mela.m
m
Definition: mela.py:715
TUtil::PrintCandidateSummary
void PrintCandidateSummary(MELACandidate *cand)
Definition: TUtil.cc:8590
TUtil::SetMCFMSpinZeroCouplings
void SetMCFMSpinZeroCouplings(bool useBSM, SpinZeroCouplings const *Hcouplings, bool forceZZ)
Definition: TUtil.cc:3360
hto_masses::msq
real *8, parameter msq
Definition: CALLING_cpHTO.f:83
TUtil::SetJHUGenSpinOneCouplings
void SetJHUGenSpinOneCouplings(double Zqqcoupl[SIZE_ZQQ][2], double Zvvcoupl[SIZE_ZVV][2])
Definition: TUtil.cc:4086
SIZE_ZVV
@ SIZE_ZVV
Definition: raw_couplings.txt:114
TUtil::forbidMassiveLeptons
bool forbidMassiveLeptons
Remove fermion mass if the flag is set to true.
Definition: TUtil.cc:25
TUtil::GetAlphaS
void GetAlphaS(double *alphas_, double *alphasmz_)
Definition: TUtil.cc:1794
TModJHUGen.hh
TUtil::SetDecayWidth
void SetDecayWidth(double inwidth, int ipart)
Definition: TUtil.cc:1491
TUtil::MCFM_masscuts
bool MCFM_masscuts(double s[][mxpart], const TVar::Process &process)
Definition: TUtil.cc:3289
TUtil::setLeptonMassScheme
void setLeptonMassScheme(TVar::FermionMassRemoval scheme=TVar::ConserveDifermionMass)
Definition: TUtil.cc:37
TUtil::SetJHUGenSpinTwoCouplings
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])
Definition: TUtil.cc:4087
SIZE_HQQ
@ SIZE_HQQ
Definition: raw_couplings.txt:5
TUtil::ConvertVectorFormat
MELACandidate * ConvertVectorFormat(SimpleParticleCollection_t *pDaughters, SimpleParticleCollection_t *pAssociated, SimpleParticleCollection_t *pMothers, bool isGen, std::vector< MELAParticle * > *particleList, std::vector< MELACandidate * > *candList)
Definition: TUtil.cc:8410
anonymous_namespace{TCouplingsBase.hh}::SIZE_HGG
@ SIZE_HGG
Definition: TCouplingsBase.hh:40
TUtil::SetJHUGenSpinZeroVVCouplings
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)
Definition: TUtil.cc:4078
TUtil::SetJHUGenDistinguishWWCouplings
void SetJHUGenDistinguishWWCouplings(bool doAllow)
Definition: TUtil.cc:3353
TUtil::SetJHUGenAZffCouplings
void SetJHUGenAZffCouplings(bool needAZff, double AZffcoupl[SIZE_AZff][2])
Definition: TUtil.cc:4093
TModHiggsMatEl.hh
TUtil::GetDecayWidth
double GetDecayWidth(int ipart)
Definition: TUtil.cc:1566
TUtil::GetCKMElement
double GetCKMElement(int iquark, int jquark)
Definition: TUtil.cc:1525
TModHiggsJMatEl.hh
TModGravitonMatEl.hh
TUtil::GetMCFMParticleLabel
TString GetMCFMParticleLabel(const int &pid, bool useQJ, bool useExtendedConventions)
Definition: TUtil.cc:3249
TUtil::computeVHAngles
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)
Definition: TUtil.cc:873
TModTTBHMatEl.hh
SIZE_Vpff
@ SIZE_Vpff
Definition: raw_couplings.txt:100
MELAParticle
Definition: MELAParticle.h:13
TUtil::computeAngles
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)
Definition: TUtil.cc:208
TUtil::ComplexBoost
std::pair< TLorentzVector, TLorentzVector > ComplexBoost(TVector3 const &beta, TLorentzVector const &p4)
Definition: TUtil.cc:166
TUtil::SetEwkCouplingParameters
void SetEwkCouplingParameters(double ext_Gf, double ext_aemmz, double ext_mW, double ext_mZ, double ext_xW, int ext_ewscheme)
Definition: TUtil.cc:1386
TVar::ConserveDifermionMass
@ ConserveDifermionMass
Definition: TVar.hh:113
TUtil::computeAnglesCS
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)
Definition: TUtil.cc:374
TVar::MatrixElement
MatrixElement
Definition: TVar.hh:55
nmsq
@ nmsq
Definition: TMCFM.hh:24
TUtil::removeMassFromPair
std::pair< TLorentzVector, TLorentzVector > removeMassFromPair(TLorentzVector const &jet1, int const &jet1Id, TLorentzVector const &jet2, int const &jet2Id, double m1=0, double m2=0)
Definition: TUtil.cc:74
TUtil::LeptonMassScheme
TVar::FermionMassRemoval LeptonMassScheme
Definition: TUtil.cc:27
TUtil::GetMassWidth
void GetMassWidth(int ipart, double &m, double &ga)
Definition: TUtil.cc:1592
SIZE_HVV_LAMBDAQSQ
@ SIZE_HVV_LAMBDAQSQ
Definition: raw_couplings.txt:66
TUtil::ConvertThreeBodyDecayCandidate
MELAThreeBodyDecayCandidate * ConvertThreeBodyDecayCandidate(SimpleParticleCollection_t *tbdDaughters, std::vector< MELAParticle * > *particleList, std::vector< MELAThreeBodyDecayCandidate * > *tbdCandList)
Definition: TUtil.cc:8554
TUtil::ComputePDF
void ComputePDF(const TLorentzVector &p0, const TLorentzVector &p1, double fx1[nmsq], double fx2[nmsq], const double &EBEAM, const TVar::VerbosityLevel &verbosity)
Definition: TUtil.cc:7877
TUtil::GetMass
double GetMass(int ipart)
Definition: TUtil.cc:1529
TVar::ResonancePropagatorScheme
ResonancePropagatorScheme
Definition: TVar.hh:117
alphas_
double alphas_(double *q, double *amz, int *nloop)
TUtil::MCFM_SetupParticleCouplings
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)
Definition: TUtil.cc:2227
mxpart
@ mxpart
Definition: TMCFM.hh:19
TUtil::forbidMassiveJets
bool forbidMassiveJets
Definition: TUtil.cc:26
TUtil::SetMCFMaTQGCCouplings
void SetMCFMaTQGCCouplings(bool useBSM, aTQGCCouplings const *couplings)
Definition: TUtil.cc:3993
TUtil::computeFakeJet
void computeFakeJet(TLorentzVector const &realJet, TLorentzVector const &others, TLorentzVector &fakeJet)
Definition: TUtil.cc:156
SIZE_AZff
@ SIZE_AZff
Definition: raw_couplings.txt:194
TUtil::TTHiggsMatEl
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)
Definition: TUtil.cc:7187
TUtil::constrainedRemovePairMass
void constrainedRemovePairMass(TLorentzVector &p1, TLorentzVector &p2, double m1=0, double m2=0)
Definition: TUtil.cc:39
SpinZeroCouplings
Definition: TCouplings.hh:7
TUtil::InitJHUGenMELA
void InitJHUGenMELA(const char *pathtoPDFSet, int PDFMember, double collider_sqrts)
Definition: TUtil.cc:3337
TUtil::SumMatrixElementPDF
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)
Definition: TUtil.cc:4103
npart
int npart
Definition: TMCFM.hh:77
TUtil::JetMassScheme
TVar::FermionMassRemoval JetMassScheme
Definition: TUtil.cc:28
TUtil::adjustTopDaughters
void adjustTopDaughters(SimpleParticleCollection_t &daughters)
Definition: TUtil.cc:132
MELACandidate
Definition: MELACandidate.h:7
TUtil::CheckPartonMomFraction
bool CheckPartonMomFraction(const TLorentzVector &p0, const TLorentzVector &p1, double xx[2], const double &EBEAM, const TVar::VerbosityLevel &verbosity)
Definition: TUtil.cc:7852
TUtil::SetJHUGenSpinZeroQQCouplings
void SetJHUGenSpinZeroQQCouplings(double Hqqcoupl[SIZE_HQQ][2])
Definition: TUtil.cc:4085
SIZE_HVV
@ SIZE_HVV
Definition: raw_couplings.txt:57
TVar::DEBUG
@ DEBUG
Definition: TVar.hh:51
mela.couplings
tuple couplings
Definition: mela.py:791
TMCFMUtils.hh
MELAThreeBodyDecayCandidate
Definition: MELAThreeBodyDecayCandidate.h:7
SimpleParticleCollection_t
std::vector< SimpleParticle_t > SimpleParticleCollection_t
Definition: TVar.hh:25
TUtil
Definition: TUtil.hh:48
sqrts
double sqrts
Definition: TMCFM.hh:290
TVar::Production
Production
Definition: TVar.hh:60
TUtil::setJetMassScheme
void setJetMassScheme(TVar::FermionMassRemoval scheme=TVar::ConserveDifermionMass)
Definition: TUtil.cc:38
TCouplings.hh
TUtil::computeVBFAngles
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)
Definition: TUtil.cc:596
TVar::simple_event_record
Definition: TVar.hh:227
TUtil::GetBoostedParticleVectors
void GetBoostedParticleVectors(MELACandidate *melaCand, TVar::simple_event_record &mela_event, TVar::VerbosityLevel verbosity=TVar::DEBUG)
Definition: TUtil.cc:7947
modparameters::process
integer, public process
Definition: mod_Parameters.F90:17
TUtil::JHUGenMatEl
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)
Definition: TUtil.cc:4746
TModVHiggsMatEl.hh
TUtil::VHiggsMatEl
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)
Definition: TUtil.cc:6578
SIZE_GQQ
@ SIZE_GQQ
Definition: raw_couplings.txt:121
TUtil::computeVBFAngles_ComplexBoost
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)
Definition: TUtil.cc:730
TUtil::ResetAmplitudeIncludes
void ResetAmplitudeIncludes()
Definition: TUtil.cc:3357
TUtil::applyJetMassCorrection
void applyJetMassCorrection(bool flag=true)
Definition: TUtil.cc:36