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.
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ZZMatrixElement Class Reference

#include <ZZMatrixElement.h>

Collaboration diagram for ZZMatrixElement:
Collaboration graph
[legend]

Public Member Functions

 ZZMatrixElement (const char *pathtoPDFSet, int PDFMember, const char *pathtoHiggsCSandWidth, double ebeam, TVar::VerbosityLevel verbosity)
 
 ZZMatrixElement (const ZZMatrixElement &other)
 
 ~ZZMatrixElement ()
 
void computeXS (float &mevalue)
 
void computeProdXS_VVHVV (float &mevalue)
 
void computeProdXS_JJH (float &mevalue)
 
void computeProdXS_JH (float &mevalue)
 
void computeProdXS_VH (float &mevalue, bool includeHiggsDecay=false)
 
void computeProdXS_ttH (float &mevalue, int topProcess, int topDecay=0)
 
void get_XPropagator (TVar::ResonancePropagatorScheme scheme, float &prop)
 
void set_Process (TVar::Process process_, TVar::MatrixElement me_, TVar::Production production_)
 
void set_Verbosity (TVar::VerbosityLevel verbosity_)
 
void set_LeptonInterference (TVar::LeptonInterference myLepInterf)
 
void set_TempCandidate (SimpleParticleCollection_t *pDaughters, SimpleParticleCollection_t *pAssociated=0, SimpleParticleCollection_t *pMothers=0, bool isGen=false)
 
void set_RenFacScaleMode (TVar::EventScaleScheme renormalizationSch, TVar::EventScaleScheme factorizationSch, double ren_sf, double fac_sf)
 
void set_LHAgrid (const char *path, int pdfmember=0)
 
void set_PrimaryHiggsMass (double mh)
 
void set_CandidateDecayMode (TVar::CandidateDecayMode mode)
 
void set_CurrentCandidateFromIndex (unsigned int icand)
 
void set_CurrentCandidate (MELACandidate *cand)
 
void set_InputEvent (SimpleParticleCollection_t *pDaughters, SimpleParticleCollection_t *pAssociated=0, SimpleParticleCollection_t *pMothers=0, bool isGen=false)
 
void append_TopCandidate (SimpleParticleCollection_t *TopDaughters)
 
void set_mHiggs (double mh_, int index)
 
void set_wHiggs (double gah_, int index)
 
void set_mHiggs_wHiggs (double mh_, double gah_, int index)
 
void reset_Mass (double inmass, int ipart)
 
void reset_Width (double inmass, int ipart)
 
void reset_QuarkMasses ()
 
void reset_MCFM_EWKParameters (double ext_Gf, double ext_aemmz, double ext_mW, double ext_mZ, double ext_xW, int ext_ewscheme=3)
 
void resetPerEvent ()
 
void reset_InputEvent ()
 
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], bool diffHWW=false)
 
void set_SpinZeroContact (double selfDHzzpcoupl[SIZE_HVV][2], double selfDHzpzpcoupl[SIZE_HVV][2], double selfDHwwpcoupl[SIZE_HVV][2], double selfDHwpwpcoupl[SIZE_HVV][2])
 
void set_SpinOneCouplings (double selfDZqqcoupl[SIZE_ZQQ][2], double selfDZvvcoupl[SIZE_ZVV][2])
 
void set_SpinTwoCouplings (double selfDGqqcoupl[SIZE_GQQ][2], double selfDGggcoupl[SIZE_GGG][2], double selfDGvvcoupl[SIZE_GVV][2])
 
void set_SpinTwoContact (double selfDGvvpcoupl[SIZE_GVV][2], double selfDGvpvpcoupl[SIZE_GVV][2])
 
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 set_aTQGCCouplings (double selfDaTQGCcoupl[SIZE_ATQGC][2])
 
void set_AZffCouplings (double selfDAZffcoupl[SIZE_AZff][2])
 
std::vector< TLorentzVector > Calculate4Momentum (double Mx, double M1, double M2, double theta, double theta1, double theta2, double Phi1, double Phi)
 
MelaIOget_IORecord ()
 
double get_PrimaryMass (int ipart)
 
double get_PrimaryHiggsMass ()
 
double get_PrimaryWidth (int ipart)
 
double get_HiggsWidthAtPoleMass (double mass)
 
MELACandidateget_CurrentCandidate ()
 
int get_CurrentCandidateIndex ()
 
int get_NCandidates ()
 
std::vector< MELATopCandidate_t * > * get_TopCandidateCollection ()
 

Protected Member Functions

void build ()
 

Protected Attributes

TVar::Process processModel
 
TVar::MatrixElement processME
 
TVar::Production processProduction
 
TVar::VerbosityLevel processVerbosity
 
TVar::LeptonInterference processLeptonInterference
 
double EBEAM
 
double mHiggs [nSupportedHiggses]
 
double wHiggs [nSupportedHiggses]
 
TEvtProb Xcal2
 
SpinZeroCouplingsselfD_SpinZeroCouplings
 
SpinOneCouplingsselfD_SpinOneCouplings
 
SpinTwoCouplingsselfD_SpinTwoCouplings
 
VprimeCouplingsselfD_VprimeCouplings
 
aTQGCCouplingsselfD_aTQGCCouplings
 
AZffCouplingsselfD_AZffCouplings
 
MELACandidatemelaCand
 
std::vector< MELAParticle * > tmpPartList
 
std::vector< MELACandidate * > tmpCandList
 

Detailed Description

Definition at line 9 of file ZZMatrixElement.h.

Constructor & Destructor Documentation

◆ ZZMatrixElement() [1/2]

ZZMatrixElement::ZZMatrixElement ( const char *  pathtoPDFSet,
int  PDFMember,
const char *  pathtoHiggsCSandWidth,
double  ebeam,
TVar::VerbosityLevel  verbosity 
)

Definition at line 12 of file ZZMatrixElement.cc.

18  :
19  processVerbosity(verbosity),
21  EBEAM(ebeam),
22  Xcal2(pathtoHiggsCSandWidth, EBEAM, pathtoPDFSet, PDFMember, verbosity),
23  melaCand(0)
24 {
25  if (processVerbosity>=TVar::DEBUG) MELAout << "Begin ZZMatrixElement constructor" << endl;
26  build();
27  if (processVerbosity>=TVar::DEBUG) MELAout << "End ZZMatrixElement constructor" << endl;
28 }

◆ ZZMatrixElement() [2/2]

ZZMatrixElement::ZZMatrixElement ( const ZZMatrixElement other)

Definition at line 29 of file ZZMatrixElement.cc.

29  :
32 EBEAM(other.EBEAM),
33 Xcal2(other.Xcal2),
34 melaCand(0) // 0 is correct in the copy constructor
35 {
36  if (processVerbosity>=TVar::DEBUG) MELAout << "Begin ZZMatrixElement copy constructor" << endl;
37  build();
38  if (processVerbosity>=TVar::DEBUG) MELAout << "End ZZMatrixElement copy constructor" << endl;
39 }

◆ ~ZZMatrixElement()

ZZMatrixElement::~ZZMatrixElement ( )

Definition at line 59 of file ZZMatrixElement.cc.

59  {
60  if (processVerbosity>=TVar::DEBUG) MELAout << "Begin ZZMatrixElement destructor" << endl;
61  resetPerEvent();
63  if (processVerbosity>=TVar::DEBUG) MELAout << "End ZZMatrixElement destructor" << endl;
64 }

Member Function Documentation

◆ append_TopCandidate()

void ZZMatrixElement::append_TopCandidate ( SimpleParticleCollection_t TopDaughters)

Definition at line 178 of file ZZMatrixElement.cc.

178 { Xcal2.AppendTopCandidate(TopDaughters); }

◆ build()

void ZZMatrixElement::build ( )
protected

Definition at line 41 of file ZZMatrixElement.cc.

41  {
42  if (processVerbosity>=TVar::DEBUG) MELAout << "Begin ZZMatrixElement::build" << endl;
43 
44  // Set default parameters explicitly
46  set_mHiggs(125., 0); set_wHiggs(-1., 0);
47  set_mHiggs(-1., 1); set_wHiggs(-1, 1);
48 
55 
56  if (processVerbosity>=TVar::DEBUG) MELAout << "End ZZMatrixElement::build" << endl;
57 }

◆ Calculate4Momentum()

std::vector< TLorentzVector > ZZMatrixElement::Calculate4Momentum ( double  Mx,
double  M1,
double  M2,
double  theta,
double  theta1,
double  theta2,
double  Phi1,
double  Phi 
)

Definition at line 67 of file ZZMatrixElement.cc.

67  {
68  double phi1, phi2;
69  phi1=TMath::Pi()-Phi1;
70  phi2=Phi1+Phi;
71 
72  double gamma1=1, gamma2=1, beta1=0, beta2=0;
73 
74  if (M1>0. && Mx>0.){
75  gamma1=(Mx*Mx+M1*M1-M2*M2)/(2*Mx*M1);
76  beta1=sqrt(1.-1./(gamma1*gamma1));
77  }
78  if (M2>0. && Mx>0.){
79  gamma2=(Mx*Mx-M1*M1+M2*M2)/(2*Mx*M2);
80  beta2=sqrt(1.-1./(gamma2*gamma2));
81  }
82 
83  //gluon 4 vectors
84  TLorentzVector p1CM(0, 0, Mx/2, Mx/2);
85  TLorentzVector p2CM(0, 0, -Mx/2, Mx/2);
86 
87  //vector boson 4 vectors
88  TLorentzVector kZ1(gamma1*M1*sin(theta)*beta1, 0, gamma1*M1*cos(theta)*beta1, gamma1*M1);
89  TLorentzVector kZ2(-gamma2*M2*sin(theta)*beta2, 0, -gamma2*M2*cos(theta)*beta2, gamma2*M2);
90 
91  //Rotation and Boost matrices. Note gamma1*beta1*M1=gamma2*beta2*M2.
92 
93  TLorentzRotation Z1ToZ, Z2ToZ;
94 
95  Z1ToZ.Boost(0, 0, beta1);
96  Z2ToZ.Boost(0, 0, beta2);
97  Z1ToZ.RotateY(theta);
98  Z2ToZ.RotateY(TMath::Pi()+theta);
99 
100 
101  //fermion 4 vectors in vector boson rest frame
102 
103  TLorentzVector p3Z1((M1/2)*sin(theta1)*cos(phi1), (M1/2)*sin(theta1)*sin(phi1), (M1/2)*cos(theta1), (M1/2)*1);
104  TLorentzVector p4Z1(-(M1/2)*sin(theta1)*cos(phi1), -(M1/2)*sin(theta1)*sin(phi1), -(M1/2)*cos(theta1), (M1/2)*1);
105  TLorentzVector p5Z2((M2/2)*sin(theta2)*cos(phi2), (M2/2)*sin(theta2)*sin(phi2), (M2/2)*cos(theta2), (M2/2)*1);
106  TLorentzVector p6Z2(-(M2/2)*sin(theta2)*cos(phi2), -(M2/2)*sin(theta2)*sin(phi2), -(M2/2)*cos(theta2), (M2/2)*1);
107 
108  // fermions 4 vectors in CM frame
109 
110  TLorentzVector p3CM, p4CM, p5CM, p6CM;
111 
112  p3CM=Z1ToZ*p3Z1;
113  p4CM=Z1ToZ*p4Z1;
114  p5CM=Z2ToZ*p5Z2;
115  p6CM=Z2ToZ*p6Z2;
116 
117  vector<TLorentzVector> p;
118 
119  p.push_back(p3CM);
120  p.push_back(p4CM);
121  p.push_back(p5CM);
122  p.push_back(p6CM);
123 
124  return p;
125 }

◆ computeProdXS_JH()

void ZZMatrixElement::computeProdXS_JH ( float &  mevalue)

Definition at line 398 of file ZZMatrixElement.cc.

400  {
402 
403  if (melaCand){
404  mevalue = Xcal2.XsecCalcXJ();
405  }
406 
407  resetPerEvent();
408  return;
409 }

◆ computeProdXS_JJH()

void ZZMatrixElement::computeProdXS_JJH ( float &  mevalue)

Definition at line 384 of file ZZMatrixElement.cc.

386  {
388 
389  if (melaCand){
390  mevalue = Xcal2.XsecCalcXJJ();
391  }
392 
393  resetPerEvent();
394  return;
395 }

◆ computeProdXS_ttH()

void ZZMatrixElement::computeProdXS_ttH ( float &  mevalue,
int  topProcess,
int  topDecay = 0 
)

Definition at line 435 of file ZZMatrixElement.cc.

439  {
441 
442  if (melaCand){
443  double zzmass = melaCand->m();
444  if (processME == TVar::MCFM){ for (int jh=0; jh<(int)nSupportedHiggses; jh++) Xcal2.SetHiggsMass(mHiggs[jh], wHiggs[jh], jh+1); }
445  else Xcal2.SetHiggsMass(zzmass, wHiggs[0], -1);
446 
448  mevalue = Xcal2.XsecCalc_TTX(
449  topProcess, topDecay
450  );
451  }
452  }
453 
454  resetPerEvent();
455  return;
456 }

◆ computeProdXS_VH()

void ZZMatrixElement::computeProdXS_VH ( float &  mevalue,
bool  includeHiggsDecay = false 
)

Definition at line 412 of file ZZMatrixElement.cc.

415  {
417 
418  if (melaCand){
419  double zzmass = melaCand->m();
420  if (processME == TVar::MCFM){ for (int jh=0; jh<(int)nSupportedHiggses; jh++) Xcal2.SetHiggsMass(mHiggs[jh], wHiggs[jh], jh+1); }
421  else Xcal2.SetHiggsMass(zzmass, wHiggs[0], -1);
422 
423  mevalue = Xcal2.XsecCalc_VX(
424  includeHiggsDecay
425  );
426  }
427 
428  resetPerEvent();
429  return;
430 }

◆ computeProdXS_VVHVV()

void ZZMatrixElement::computeProdXS_VVHVV ( float &  mevalue)

Definition at line 364 of file ZZMatrixElement.cc.

366  {
368 
369  if (melaCand){
370  double zzmass = melaCand->m();
371  if (processME == TVar::MCFM){
372  for (int jh=0; jh<(int)nSupportedHiggses; jh++) Xcal2.SetHiggsMass(mHiggs[jh], wHiggs[jh], jh+1);
373  }
374  else Xcal2.SetHiggsMass(zzmass, wHiggs[0], -1);
375 
376  mevalue = Xcal2.XsecCalc_VVXVV();
377  }
378 
379  resetPerEvent();
380  return;
381 }

◆ computeXS()

void ZZMatrixElement::computeXS ( float &  mevalue)

Definition at line 344 of file ZZMatrixElement.cc.

346  {
348 
349  if (melaCand){
350  double zzmass = melaCand->m();
351  if (processME == TVar::MCFM){
352  for (int jh=0; jh<(int)nSupportedHiggses; jh++) Xcal2.SetHiggsMass(mHiggs[jh], wHiggs[jh], jh+1);
353  }
354  else Xcal2.SetHiggsMass(zzmass, wHiggs[0], -1);
355 
356  mevalue = Xcal2.XsecCalc_XVV();
357  }
358 
359  resetPerEvent();
360  return;
361 }

◆ get_CurrentCandidate()

MELACandidate * ZZMatrixElement::get_CurrentCandidate ( )

Definition at line 232 of file ZZMatrixElement.cc.

232 { return Xcal2.GetCurrentCandidate(); }

◆ get_CurrentCandidateIndex()

int ZZMatrixElement::get_CurrentCandidateIndex ( )

Definition at line 233 of file ZZMatrixElement.cc.

233 { return Xcal2.GetCurrentCandidateIndex(); }

◆ get_HiggsWidthAtPoleMass()

double ZZMatrixElement::get_HiggsWidthAtPoleMass ( double  mass)

Definition at line 231 of file ZZMatrixElement.cc.

231 { return Xcal2.GetHiggsWidthAtPoleMass(mass); }

◆ get_IORecord()

MelaIO * ZZMatrixElement::get_IORecord ( )

Definition at line 228 of file ZZMatrixElement.cc.

228 { return Xcal2.GetIORecord(); }

◆ get_NCandidates()

int ZZMatrixElement::get_NCandidates ( )

Definition at line 234 of file ZZMatrixElement.cc.

234 { return Xcal2.GetNCandidates(); }

◆ get_PrimaryHiggsMass()

double ZZMatrixElement::get_PrimaryHiggsMass ( )
inline

Definition at line 141 of file ZZMatrixElement.h.

141 { return get_PrimaryMass(25); }

◆ get_PrimaryMass()

double ZZMatrixElement::get_PrimaryMass ( int  ipart)

Definition at line 229 of file ZZMatrixElement.cc.

229 { return Xcal2.GetPrimaryMass(ipart); }

◆ get_PrimaryWidth()

double ZZMatrixElement::get_PrimaryWidth ( int  ipart)

Definition at line 230 of file ZZMatrixElement.cc.

230 { return Xcal2.GetPrimaryWidth(ipart); }

◆ get_TopCandidateCollection()

std::vector< MELATopCandidate_t * > * ZZMatrixElement::get_TopCandidateCollection ( )

Definition at line 235 of file ZZMatrixElement.cc.

235 { return Xcal2.GetTopCandidates(); }

◆ get_XPropagator()

void ZZMatrixElement::get_XPropagator ( TVar::ResonancePropagatorScheme  scheme,
float &  prop 
)

Definition at line 459 of file ZZMatrixElement.cc.

459  {
460  prop=0.;
462 
463  if (melaCand){
464  Xcal2.SetHiggsMass(mHiggs[0], wHiggs[0], -1);
465  prop=Xcal2.GetXPropagator(scheme);
466  }
467 
468  resetPerEvent();
469  return;
470 }

◆ reset_InputEvent()

void ZZMatrixElement::reset_InputEvent ( )

Definition at line 225 of file ZZMatrixElement.cc.

225 { Xcal2.ResetInputEvent(); }

◆ reset_Mass()

void ZZMatrixElement::reset_Mass ( double  inmass,
int  ipart 
)

Definition at line 197 of file ZZMatrixElement.cc.

197 { Xcal2.ResetMass(inmass, ipart); }

◆ reset_MCFM_EWKParameters()

void ZZMatrixElement::reset_MCFM_EWKParameters ( double  ext_Gf,
double  ext_aemmz,
double  ext_mW,
double  ext_mZ,
double  ext_xW,
int  ext_ewscheme = 3 
)

Definition at line 200 of file ZZMatrixElement.cc.

200  {
201  Xcal2.ResetMCFM_EWKParameters(ext_Gf, ext_aemmz, ext_mW, ext_mZ, ext_xW, ext_ewscheme);
202 }

◆ reset_QuarkMasses()

void ZZMatrixElement::reset_QuarkMasses ( )

Definition at line 199 of file ZZMatrixElement.cc.

199 { Xcal2.ResetQuarkMasses(); }

◆ reset_Width()

void ZZMatrixElement::reset_Width ( double  inmass,
int  ipart 
)

Definition at line 198 of file ZZMatrixElement.cc.

198 { Xcal2.ResetWidth(inwidth, ipart); }

◆ resetPerEvent()

void ZZMatrixElement::resetPerEvent ( )

Definition at line 205 of file ZZMatrixElement.cc.

205  {
206  // Protection against forgetfulness; custom width has to be set per-computation
207  set_mHiggs(Xcal2.GetPrimaryHiggsMass(), 0); // Sets mHiggs[0]
208  if (wHiggs[0]>=0.) set_wHiggs(-1., 0);
209 
210  if (mHiggs[1]>=0.) set_mHiggs(-1., 1);
211  if (wHiggs[1]>=0.) set_wHiggs(-1., 1);
212 
213  Xcal2.SetHiggsMass(mHiggs[0], -1, -1);
214 
215  // Return back to default lepton interference settings after each calculation
217 
218  // Delete the temporary input objects owned
219  for (unsigned int ic=0; ic<tmpCandList.size(); ic++){ if (tmpCandList.at(ic)) delete tmpCandList.at(ic); }
220  //for (unsigned int itc=0; itc<tmpTopCandList.size(); itc++){ if (tmpTopCandList.at(itc)) delete tmpTopCandList.at(itc); }
221  for (unsigned int ip=0; ip<tmpPartList.size(); ip++){ if (tmpPartList.at(ip)) delete tmpPartList.at(ip); }
222  melaCand=nullptr;
223 }

◆ set_aTQGCCouplings()

void ZZMatrixElement::set_aTQGCCouplings ( double  selfDaTQGCcoupl[SIZE_ATQGC][2])

Definition at line 330 of file ZZMatrixElement.cc.

332  {
333  for (int ic=0; ic<SIZE_ATQGC; ic++) selfD_aTQGCCouplings->SetATQGCCouplings(ic, selfDaTQGCcoupl[ic][0], selfDaTQGCcoupl[ic][1]);
334 }

◆ set_AZffCouplings()

void ZZMatrixElement::set_AZffCouplings ( double  selfDAZffcoupl[SIZE_AZff][2])

Definition at line 335 of file ZZMatrixElement.cc.

337  {
338  for (int ic=0; ic<SIZE_AZff; ic++) selfD_AZffCouplings->SetAZffCouplings(ic, selfDAZffcoupl[ic][0], selfDAZffcoupl[ic][1]);
339 }

◆ set_CandidateDecayMode()

void ZZMatrixElement::set_CandidateDecayMode ( TVar::CandidateDecayMode  mode)

Definition at line 141 of file ZZMatrixElement.cc.

141 { Xcal2.SetCandidateDecayMode(mode); }

◆ set_CurrentCandidate()

void ZZMatrixElement::set_CurrentCandidate ( MELACandidate cand)

Definition at line 144 of file ZZMatrixElement.cc.

144 { Xcal2.SetCurrentCandidate(cand); }

◆ set_CurrentCandidateFromIndex()

void ZZMatrixElement::set_CurrentCandidateFromIndex ( unsigned int  icand)

Definition at line 143 of file ZZMatrixElement.cc.

◆ set_InputEvent()

void ZZMatrixElement::set_InputEvent ( SimpleParticleCollection_t pDaughters,
SimpleParticleCollection_t pAssociated = 0,
SimpleParticleCollection_t pMothers = 0,
bool  isGen = false 
)

Definition at line 145 of file ZZMatrixElement.cc.

150  { // Adds a new candidate to Xcal2
152  pDaughters,
153  pAssociated,
154  pMothers,
155  isGen
156  );
157 }

◆ set_LeptonInterference()

void ZZMatrixElement::set_LeptonInterference ( TVar::LeptonInterference  myLepInterf)

◆ set_LHAgrid()

void ZZMatrixElement::set_LHAgrid ( const char *  path,
int  pdfmember = 0 
)

Definition at line 137 of file ZZMatrixElement.cc.

137 { Xcal2.Set_LHAgrid(path, pdfmember); }

◆ set_mHiggs()

void ZZMatrixElement::set_mHiggs ( double  mh_,
int  index 
)

Definition at line 180 of file ZZMatrixElement.cc.

180  {
181  if (index<nSupportedHiggses && index>=0) mHiggs[index] = mh_;
182  else MELAerr << "ZZMatrixElement::set_mHiggs: Only resonances 0 (regular) and 1 (additional, possibly high-mass) are supported" << endl;
183 }

◆ set_mHiggs_wHiggs()

void ZZMatrixElement::set_mHiggs_wHiggs ( double  mh_,
double  gah_,
int  index 
)

Definition at line 188 of file ZZMatrixElement.cc.

188  {
189  if (index<nSupportedHiggses && index>=0){
190  mHiggs[index] = mh_;
191  wHiggs[index] = gah_;
192  }
193  else MELAerr << "ZZMatrixElement::set_mHiggs_wHiggs: Only resonances 0 (regular) and 1 (additional, possibly high-mass) are supported" << endl;
194 }

◆ set_PrimaryHiggsMass()

void ZZMatrixElement::set_PrimaryHiggsMass ( double  mh)

Definition at line 142 of file ZZMatrixElement.cc.

142 { Xcal2.SetPrimaryHiggsMass(mh); }

◆ set_Process()

void ZZMatrixElement::set_Process ( TVar::Process  process_,
TVar::MatrixElement  me_,
TVar::Production  production_ 
)

Definition at line 128 of file ZZMatrixElement.cc.

128  {
129  processModel = process_;
130  processME = me_;
131  processProduction = production_;
133 }

◆ set_RenFacScaleMode()

void ZZMatrixElement::set_RenFacScaleMode ( TVar::EventScaleScheme  renormalizationSch,
TVar::EventScaleScheme  factorizationSch,
double  ren_sf,
double  fac_sf 
)

Definition at line 138 of file ZZMatrixElement.cc.

138  {
139  Xcal2.SetRenFacScaleMode(renormalizationSch, factorizationSch, ren_sf, fac_sf);
140 }

◆ set_SpinOneCouplings()

void ZZMatrixElement::set_SpinOneCouplings ( double  selfDZqqcoupl[SIZE_ZQQ][2],
double  selfDZvvcoupl[SIZE_ZVV][2] 
)

Definition at line 290 of file ZZMatrixElement.cc.

293  {
294  for (int ic=0; ic<SIZE_ZQQ; ic++) selfD_SpinOneCouplings->SetZQQCouplings(ic, selfDZqqcoupl[ic][0], selfDZqqcoupl[ic][1]);
295  for (int ic=0; ic<SIZE_ZVV; ic++) selfD_SpinOneCouplings->SetZVVCouplings(ic, selfDZvvcoupl[ic][0], selfDZvvcoupl[ic][1]);
296 }

◆ set_SpinTwoContact()

void ZZMatrixElement::set_SpinTwoContact ( double  selfDGvvpcoupl[SIZE_GVV][2],
double  selfDGvpvpcoupl[SIZE_GVV][2] 
)

Definition at line 306 of file ZZMatrixElement.cc.

309  {
310  for (int ic=0; ic<SIZE_GVV; ic++){
313  }
314 }

◆ set_SpinTwoCouplings()

void ZZMatrixElement::set_SpinTwoCouplings ( double  selfDGqqcoupl[SIZE_GQQ][2],
double  selfDGggcoupl[SIZE_GGG][2],
double  selfDGvvcoupl[SIZE_GVV][2] 
)

Definition at line 297 of file ZZMatrixElement.cc.

301  {
302  for (int ic=0; ic<SIZE_GQQ; ic++) selfD_SpinTwoCouplings->SetGQQCouplings(ic, selfDGqqcoupl[ic][0], selfDGqqcoupl[ic][1]);
303  for (int ic=0; ic<SIZE_GGG; ic++) selfD_SpinTwoCouplings->SetGGGCouplings(ic, selfDGggcoupl[ic][0], selfDGggcoupl[ic][1]);
304  for (int ic=0; ic<SIZE_GVV; ic++) selfD_SpinTwoCouplings->SetGVVCouplings(ic, selfDGvvcoupl[ic][0], selfDGvvcoupl[ic][1]);
305 }

◆ set_SpinZeroContact()

void ZZMatrixElement::set_SpinZeroContact ( double  selfDHzzpcoupl[SIZE_HVV][2],
double  selfDHzpzpcoupl[SIZE_HVV][2],
double  selfDHwwpcoupl[SIZE_HVV][2],
double  selfDHwpwpcoupl[SIZE_HVV][2] 
)

Definition at line 277 of file ZZMatrixElement.cc.

282  {
283  for (int ic=0; ic<SIZE_HVV; ic++){
288  }
289 }

◆ set_SpinZeroCouplings()

void ZZMatrixElement::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],
bool  diffHWW = false 
)

Definition at line 238 of file ZZMatrixElement.cc.

253  {
255  for (int jh=1; jh<=(int)nSupportedHiggses; jh++){
256  for (int ic=0; ic<SIZE_HGG; ic++) selfD_SpinZeroCouplings->SetHGGCouplings(ic, selfDHggcoupl[jh-1][ic][0], selfDHggcoupl[jh-1][ic][1], 1, jh);
257  for (int ic=0; ic<SIZE_HGG; ic++) selfD_SpinZeroCouplings->SetHGGCouplings(ic, selfDHg4g4coupl[jh-1][ic][0], selfDHg4g4coupl[jh-1][ic][1], 2, jh);
258 
259  for (int ic=0; ic<SIZE_HQQ; ic++) selfD_SpinZeroCouplings->SetHQQCouplings(ic, selfDHqqcoupl[jh-1][ic][0], selfDHqqcoupl[jh-1][ic][1], 0, jh);
260  for (int ic=0; ic<SIZE_HQQ; ic++) selfD_SpinZeroCouplings->SetHQQCouplings(ic, selfDHbbcoupl[jh-1][ic][0], selfDHbbcoupl[jh-1][ic][1], 5, jh);
261  for (int ic=0; ic<SIZE_HQQ; ic++) selfD_SpinZeroCouplings->SetHQQCouplings(ic, selfDHttcoupl[jh-1][ic][0], selfDHttcoupl[jh-1][ic][1], 6, jh);
262  for (int ic=0; ic<SIZE_HQQ; ic++) selfD_SpinZeroCouplings->SetHQQCouplings(ic, selfDHb4b4coupl[jh-1][ic][0], selfDHb4b4coupl[jh-1][ic][1], 7, jh);
263  for (int ic=0; ic<SIZE_HQQ; ic++) selfD_SpinZeroCouplings->SetHQQCouplings(ic, selfDHt4t4coupl[jh-1][ic][0], selfDHt4t4coupl[jh-1][ic][1], 8, jh);
264 
265  for (int ic=0; ic<SIZE_HVV; ic++) selfD_SpinZeroCouplings->SetHVVCouplings(ic, selfDHzzcoupl[jh-1][ic][0], selfDHzzcoupl[jh-1][ic][1], false, jh);
266  for (int ic=0; ic<SIZE_HVV; ic++) selfD_SpinZeroCouplings->SetHVVCouplings(ic, selfDHwwcoupl[jh-1][ic][0], selfDHwwcoupl[jh-1][ic][1], true, jh);
267  for (int ik=0; ik<SIZE_HVV_CQSQ; ik++){
268  for (int ig=0; ig<SIZE_HVV_LAMBDAQSQ; ig++){
269  selfD_SpinZeroCouplings->SetHVVLambdaQ2(ig, ik, selfDHzzLambda_qsq[jh-1][ig][ik], false, jh);
270  selfD_SpinZeroCouplings->SetHVVLambdaQ2(ig, ik, selfDHwwLambda_qsq[jh-1][ig][ik], true, jh);
271  }
272  selfD_SpinZeroCouplings->SetHVVSignCQ2(ik, selfDHzzCLambda_qsq[jh-1][ik], false, jh);
274  }
275  }
276 }

◆ set_TempCandidate()

void ZZMatrixElement::set_TempCandidate ( SimpleParticleCollection_t pDaughters,
SimpleParticleCollection_t pAssociated = 0,
SimpleParticleCollection_t pMothers = 0,
bool  isGen = false 
)

Definition at line 159 of file ZZMatrixElement.cc.

164  {
166  pDaughters,
167  pAssociated,
168  pMothers,
169  isGen,
170  &tmpPartList, &tmpCandList // push_back is done automatically
171  );
172  if (cand){
173  melaCand=cand;
175  }
176 }

◆ set_Verbosity()

void ZZMatrixElement::set_Verbosity ( TVar::VerbosityLevel  verbosity_)

Definition at line 134 of file ZZMatrixElement.cc.

134 { processVerbosity = verbosity_; Xcal2.SetVerbosity(verbosity_); }

◆ set_VprimeContactCouplings()

void ZZMatrixElement::set_VprimeContactCouplings ( double  selfDZpffcoupl[SIZE_Vpff][2],
double  selfDWpffcoupl[SIZE_Vpff][2],
double  M_Zprime,
double  Ga_Zprime,
double  M_Wprime,
double  Ga_Wprime 
)

Definition at line 315 of file ZZMatrixElement.cc.

322  {
323  for (int ic=0; ic<SIZE_Vpff; ic++){
326  }
327  selfD_VprimeCouplings->SetZPrimeMassWidth(M_Zprime, Ga_Zprime);
328  selfD_VprimeCouplings->SetWPrimeMassWidth(M_Wprime, Ga_Wprime);
329 }

◆ set_wHiggs()

void ZZMatrixElement::set_wHiggs ( double  gah_,
int  index 
)

Definition at line 184 of file ZZMatrixElement.cc.

184  {
185  if (index<nSupportedHiggses && index>=0) wHiggs[index] = (double)gah_;
186  else MELAerr << "ZZMatrixElement::set_wHiggs: Only resonances 0 (regular) and 1 (additional, possibly high-mass) are supported" << endl;
187 }

Member Data Documentation

◆ EBEAM

double ZZMatrixElement::EBEAM
protected

Definition at line 157 of file ZZMatrixElement.h.

◆ melaCand

MELACandidate* ZZMatrixElement::melaCand
protected

Definition at line 169 of file ZZMatrixElement.h.

◆ mHiggs

double ZZMatrixElement::mHiggs[nSupportedHiggses]
protected

Definition at line 158 of file ZZMatrixElement.h.

◆ processLeptonInterference

TVar::LeptonInterference ZZMatrixElement::processLeptonInterference
protected

Definition at line 155 of file ZZMatrixElement.h.

◆ processME

TVar::MatrixElement ZZMatrixElement::processME
protected

Definition at line 152 of file ZZMatrixElement.h.

◆ processModel

TVar::Process ZZMatrixElement::processModel
protected

Definition at line 151 of file ZZMatrixElement.h.

◆ processProduction

TVar::Production ZZMatrixElement::processProduction
protected

Definition at line 153 of file ZZMatrixElement.h.

◆ processVerbosity

TVar::VerbosityLevel ZZMatrixElement::processVerbosity
protected

Definition at line 154 of file ZZMatrixElement.h.

◆ selfD_aTQGCCouplings

aTQGCCouplings* ZZMatrixElement::selfD_aTQGCCouplings
protected

Definition at line 166 of file ZZMatrixElement.h.

◆ selfD_AZffCouplings

AZffCouplings* ZZMatrixElement::selfD_AZffCouplings
protected

Definition at line 167 of file ZZMatrixElement.h.

◆ selfD_SpinOneCouplings

SpinOneCouplings* ZZMatrixElement::selfD_SpinOneCouplings
protected

Definition at line 163 of file ZZMatrixElement.h.

◆ selfD_SpinTwoCouplings

SpinTwoCouplings* ZZMatrixElement::selfD_SpinTwoCouplings
protected

Definition at line 164 of file ZZMatrixElement.h.

◆ selfD_SpinZeroCouplings

SpinZeroCouplings* ZZMatrixElement::selfD_SpinZeroCouplings
protected

Definition at line 162 of file ZZMatrixElement.h.

◆ selfD_VprimeCouplings

VprimeCouplings* ZZMatrixElement::selfD_VprimeCouplings
protected

Definition at line 165 of file ZZMatrixElement.h.

◆ tmpCandList

std::vector<MELACandidate*> ZZMatrixElement::tmpCandList
protected

Definition at line 173 of file ZZMatrixElement.h.

◆ tmpPartList

std::vector<MELAParticle*> ZZMatrixElement::tmpPartList
protected

Definition at line 170 of file ZZMatrixElement.h.

◆ wHiggs

double ZZMatrixElement::wHiggs[nSupportedHiggses]
protected

Definition at line 159 of file ZZMatrixElement.h.

◆ Xcal2

TEvtProb ZZMatrixElement::Xcal2
protected

Definition at line 160 of file ZZMatrixElement.h.


The documentation for this class was generated from the following files:
TEvtProb::XsecCalc_TTX
double XsecCalc_TTX(int topProcess, int topDecay)
Definition: TEvtProb.cc:975
VprimeCouplings::SetWPrimeMassWidth
void SetWPrimeMassWidth(double inmass, double inwidth)
Definition: TCouplings.cc:435
selfDHg4g4coupl
double selfDHg4g4coupl[nSupportedHiggses][SIZE_HGG][2]
Definition: raw_names.txt:2
SpinTwoCouplings::SetGQQCouplings
void SetGQQCouplings(unsigned int index, double c_real, double c_imag)
Definition: TCouplings.cc:364
SpinTwoCouplings::SetGVVpCouplings
void SetGVVpCouplings(unsigned int index, double c_real, double c_imag)
Definition: TCouplings.cc:350
TEvtProb::SetCandidateDecayMode
void SetCandidateDecayMode(TVar::CandidateDecayMode mode)
Definition: TEvtProb.cc:208
TEvtProb::GetHiggsWidthAtPoleMass
double GetHiggsWidthAtPoleMass(double mass)
Definition: TEvtProb.cc:371
SIZE_GVV
@ SIZE_GVV
Definition: raw_couplings.txt:158
selfDHzzpcoupl
double selfDHzzpcoupl[SIZE_HVV][2]
Definition: raw_names.txt:15
TEvtProb::GetNCandidates
int GetNCandidates()
Definition: TEvtProb.cc:384
TEvtProb::GetSelfDaTQGCCouplings
aTQGCCouplings * GetSelfDaTQGCCouplings()
Definition: TEvtProb.cc:360
TEvtProb::GetSelfDVprimeCouplings
VprimeCouplings * GetSelfDVprimeCouplings()
Definition: TEvtProb.cc:359
selfDaTQGCcoupl
double selfDaTQGCcoupl[SIZE_ATQGC][2]
Definition: raw_names.txt:32
ZZMatrixElement::tmpPartList
std::vector< MELAParticle * > tmpPartList
Definition: ZZMatrixElement.h:170
selfDHzzLambda_qsq
double selfDHzzLambda_qsq[nSupportedHiggses][SIZE_HVV_LAMBDAQSQ][SIZE_HVV_CQSQ]
Definition: raw_names.txt:10
selfDGggcoupl
double selfDGggcoupl[SIZE_GGG][2]
Definition: raw_names.txt:28
TEvtProb::XsecCalcXJ
double XsecCalcXJ()
Definition: TEvtProb.cc:868
SIZE_HVV_CQSQ
@ SIZE_HVV_CQSQ
Definition: raw_couplings.txt:74
selfDZpffcoupl
double selfDZpffcoupl[SIZE_Vpff][2]
Definition: raw_names.txt:17
selfDHb4b4coupl
double selfDHb4b4coupl[nSupportedHiggses][SIZE_HQQ][2]
Definition: raw_names.txt:6
SpinTwoCouplings::SetGVVCouplings
void SetGVVCouplings(unsigned int index, double c_real, double c_imag)
Definition: TCouplings.cc:343
selfDZqqcoupl
double selfDZqqcoupl[SIZE_ZQQ][2]
Definition: raw_names.txt:25
TEvtProb::SetHiggsMass
void SetHiggsMass(double mass, double wHiggs=-1., int whichResonance=-1)
Definition: TEvtProb.cc:217
SpinZeroCouplings::SetHQQCouplings
void SetHQQCouplings(unsigned int index, double c_real, double c_imag, int qid=0, int whichResonance=1)
Definition: TCouplings.cc:190
SIZE_GGG
@ SIZE_GGG
Definition: raw_couplings.txt:131
ZZMatrixElement::processLeptonInterference
TVar::LeptonInterference processLeptonInterference
Definition: ZZMatrixElement.h:155
selfDGqqcoupl
double selfDGqqcoupl[SIZE_GQQ][2]
Definition: raw_names.txt:27
TEvtProb::GetTopCandidates
std::vector< MELATopCandidate_t * > * GetTopCandidates()
Definition: TEvtProb.cc:385
TEvtProb::GetSelfDSpinOneCouplings
SpinOneCouplings * GetSelfDSpinOneCouplings()
Definition: TEvtProb.cc:357
hto_betacom::beta2
real *8, dimension(3:6) beta2
Definition: CALLING_cpHTO.f:2080
TEvtProb::GetPrimaryWidth
double GetPrimaryWidth(int ipart)
Definition: TEvtProb.cc:367
TEvtProb::XsecCalc_XVV
double XsecCalc_XVV()
Definition: TEvtProb.cc:407
SpinOneCouplings::SetZVVCouplings
void SetZVVCouplings(unsigned int index, double c_real, double c_imag)
Definition: TCouplings.cc:296
selfDHt4t4coupl
double selfDHt4t4coupl[nSupportedHiggses][SIZE_HQQ][2]
Definition: raw_names.txt:7
AZffCouplings::SetAZffCouplings
void SetAZffCouplings(unsigned int index, double c_real, double c_imag)
Definition: TCouplings.cc:480
TEvtProb::GetPrimaryHiggsMass
double GetPrimaryHiggsMass()
Definition: TEvtProb.cc:362
SIZE_ZQQ
@ SIZE_ZQQ
Definition: raw_couplings.txt:107
TEvtProb::ResetMCFM_EWKParameters
void ResetMCFM_EWKParameters(double ext_Gf, double ext_aemmz, double ext_mW, double ext_mZ, double ext_xW, int ext_ewscheme=3)
Definition: TEvtProb.cc:327
TEvtProb::ResetInputEvent
void ResetInputEvent()
Definition: TEvtProb.cc:338
testME_all.p
p
Definition: testME_all.py:11
selfDHbbcoupl
double selfDHbbcoupl[nSupportedHiggses][SIZE_HQQ][2]
Definition: raw_names.txt:4
ZZMatrixElement::wHiggs
double wHiggs[nSupportedHiggses]
Definition: ZZMatrixElement.h:159
SpinZeroCouplings::SetHGGCouplings
void SetHGGCouplings(unsigned int index, double c_real, double c_imag, int whichLoop=1, int whichResonance=1)
Definition: TCouplings.cc:163
ZZMatrixElement::processModel
TVar::Process processModel
Definition: ZZMatrixElement.h:151
ZZMatrixElement::EBEAM
double EBEAM
Definition: ZZMatrixElement.h:157
selfDHzpzpcoupl
double selfDHzpzpcoupl[SIZE_HVV][2]
Definition: raw_names.txt:16
ZZMatrixElement::mHiggs
double mHiggs[nSupportedHiggses]
Definition: ZZMatrixElement.h:158
TEvtProb::AppendTopCandidate
void AppendTopCandidate(SimpleParticleCollection_t *TopDaughters)
Definition: TEvtProb.cc:287
SIZE_ZVV
@ SIZE_ZVV
Definition: raw_couplings.txt:114
selfDHwwCLambda_qsq
int selfDHwwCLambda_qsq[nSupportedHiggses][SIZE_HVV_CQSQ]
Definition: raw_names.txt:13
ZZMatrixElement::Xcal2
TEvtProb Xcal2
Definition: ZZMatrixElement.h:160
SpinZeroCouplings::SetHVVCouplings
void SetHVVCouplings(unsigned int index, double c_real, double c_imag, bool setWW=false, int whichResonance=1)
Definition: TCouplings.cc:105
ZZMatrixElement::processME
TVar::MatrixElement processME
Definition: ZZMatrixElement.h:152
ZZMatrixElement::tmpCandList
std::vector< MELACandidate * > tmpCandList
Definition: ZZMatrixElement.h:173
ZZMatrixElement::get_PrimaryMass
double get_PrimaryMass(int ipart)
Definition: ZZMatrixElement.cc:229
ZZMatrixElement::get_CurrentCandidate
MELACandidate * get_CurrentCandidate()
Definition: ZZMatrixElement.cc:232
TEvtProb::SetRenFacScaleMode
void SetRenFacScaleMode(TVar::EventScaleScheme renormalizationSch, TVar::EventScaleScheme factorizationSch, double ren_sf, double fac_sf)
Definition: TEvtProb.cc:209
SpinTwoCouplings::SetGVpVpCouplings
void SetGVpVpCouplings(unsigned int index, double c_real, double c_imag)
Definition: TCouplings.cc:357
nSupportedHiggses
@ nSupportedHiggses
Definition: TMCFM.hh:15
VprimeCouplings::SetVpffCouplings
void SetVpffCouplings(unsigned int index, double c_real, double c_imag, bool setWpff=false, int whichResonance=1)
Definition: TCouplings.cc:411
MELAStreamHelpers::MELAout
MELAOutputStreamer MELAout
selfDHzzcoupl
double selfDHzzcoupl[nSupportedHiggses][SIZE_HVV][2]
Definition: raw_names.txt:8
TVar::DefaultLeptonInterf
@ DefaultLeptonInterf
Definition: TVar.hh:107
TEvtProb::SetLeptonInterf
void SetLeptonInterf(TVar::LeptonInterference tmp)
Definition: TEvtProb.cc:207
MELAParticle::m
double m() const
Definition: MELAParticle.h:66
ZZMatrixElement::set_mHiggs
void set_mHiggs(double mh_, int index)
Definition: ZZMatrixElement.cc:180
SIZE_HQQ
@ SIZE_HQQ
Definition: raw_couplings.txt:5
ZZMatrixElement::build
void build()
Definition: ZZMatrixElement.cc:41
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
TEvtProb::SetCurrentCandidateFromIndex
void SetCurrentCandidateFromIndex(unsigned int icand)
Definition: TEvtProb.cc:299
anonymous_namespace{TCouplingsBase.hh}::SIZE_HGG
@ SIZE_HGG
Definition: TCouplingsBase.hh:40
ZZMatrixElement::melaCand
MELACandidate * melaCand
Definition: ZZMatrixElement.h:169
testME_all.int
int
Definition: testME_all.py:13
selfDHqqcoupl
double selfDHqqcoupl[nSupportedHiggses][SIZE_HQQ][2]
Definition: raw_names.txt:3
SpinZeroCouplings::SetHVVSignCQ2
void SetHVVSignCQ2(unsigned int index, int csign, bool setWW=false, int whichResonance=1)
Definition: TCouplings.cc:147
ZZMatrixElement::selfD_AZffCouplings
AZffCouplings * selfD_AZffCouplings
Definition: ZZMatrixElement.h:167
ZZMatrixElement::processVerbosity
TVar::VerbosityLevel processVerbosity
Definition: ZZMatrixElement.h:154
VprimeCouplings::SetZPrimeMassWidth
void SetZPrimeMassWidth(double inmass, double inwidth)
Definition: TCouplings.cc:434
TEvtProb::SetInputEvent
void SetInputEvent(SimpleParticleCollection_t *pDaughters, SimpleParticleCollection_t *pAssociated=0, SimpleParticleCollection_t *pMothers=0, bool isGen=false)
Definition: TEvtProb.cc:272
ZZMatrixElement::selfD_VprimeCouplings
VprimeCouplings * selfD_VprimeCouplings
Definition: ZZMatrixElement.h:165
TEvtProb::XsecCalc_VVXVV
double XsecCalc_VVXVV()
Definition: TEvtProb.cc:729
TEvtProb::GetSelfDSpinTwoCouplings
SpinTwoCouplings * GetSelfDSpinTwoCouplings()
Definition: TEvtProb.cc:358
selfDHzzCLambda_qsq
int selfDHzzCLambda_qsq[nSupportedHiggses][SIZE_HVV_CQSQ]
Definition: raw_names.txt:12
TEvtProb::XsecCalcXJJ
double XsecCalcXJJ()
Definition: TEvtProb.cc:779
TVar::MCFM
@ MCFM
Definition: TVar.hh:56
SIZE_Vpff
@ SIZE_Vpff
Definition: raw_couplings.txt:100
selfDHggcoupl
double selfDHggcoupl[nSupportedHiggses][SIZE_HGG][2]
Definition: raw_names.txt:1
aTQGCCouplings::SetATQGCCouplings
void SetATQGCCouplings(unsigned int index, double c_real, double c_imag)
Definition: TCouplings.cc:455
TEvtProb::GetIORecord
MelaIO * GetIORecord()
Definition: TEvtProb.cc:375
TEvtProb::ResetMass
void ResetMass(double inmass, int ipart)
Definition: TEvtProb.cc:312
TEvtProb::SetProcess
void SetProcess(TVar::Process proc, TVar::MatrixElement me, TVar::Production prod)
Definition: TEvtProb.cc:192
SIZE_ATQGC
@ SIZE_ATQGC
Definition: raw_couplings.txt:176
ZZMatrixElement::selfD_aTQGCCouplings
aTQGCCouplings * selfD_aTQGCCouplings
Definition: ZZMatrixElement.h:166
TEvtProb::GetPrimaryMass
double GetPrimaryMass(int ipart)
Definition: TEvtProb.cc:363
TEvtProb::GetCurrentCandidateIndex
int GetCurrentCandidateIndex()
Definition: TEvtProb.cc:377
TEvtProb::SetVerbosity
void SetVerbosity(TVar::VerbosityLevel tmp)
Definition: TEvtProb.cc:206
selfDHwwpcoupl
double selfDHwwpcoupl[SIZE_HVV][2]
Definition: raw_names.txt:18
TEvtProb::Set_LHAgrid
void Set_LHAgrid(const char *path, int pdfmember=0)
Definition: TEvtProb.cc:185
SIZE_HVV_LAMBDAQSQ
@ SIZE_HVV_LAMBDAQSQ
Definition: raw_couplings.txt:66
SpinZeroCouplings::SetHVVpCouplings
void SetHVVpCouplings(unsigned int index, double c_real, double c_imag, bool setWWp=false, int whichResonance=1)
Definition: TCouplings.cc:242
MELAStreamHelpers::MELAerr
MELAOutputStreamer MELAerr
selfDHttcoupl
double selfDHttcoupl[nSupportedHiggses][SIZE_HQQ][2]
Definition: raw_names.txt:5
ZZMatrixElement::resetPerEvent
void resetPerEvent()
Definition: ZZMatrixElement.cc:205
SpinZeroCouplings::SetHVpVpCouplings
void SetHVpVpCouplings(unsigned int index, double c_real, double c_imag, bool setWpWp=false, int whichResonance=1)
Definition: TCouplings.cc:256
TEvtProb::GetSelfDSpinZeroCouplings
SpinZeroCouplings * GetSelfDSpinZeroCouplings()
Definition: TEvtProb.cc:356
TEvtProb::GetCurrentCandidate
MELACandidate * GetCurrentCandidate()
Definition: TEvtProb.cc:376
SIZE_AZff
@ SIZE_AZff
Definition: raw_couplings.txt:194
TEvtProb::ResetWidth
void ResetWidth(double inwidth, int ipart)
Definition: TEvtProb.cc:315
ZZMatrixElement::selfD_SpinZeroCouplings
SpinZeroCouplings * selfD_SpinZeroCouplings
Definition: ZZMatrixElement.h:162
TEvtProb::AllowSeparateWWCouplings
void AllowSeparateWWCouplings(bool doAllow=false)
Definition: TEvtProb.cc:215
TVar::ttH
@ ttH
Definition: TVar.hh:66
TEvtProb::SetPrimaryHiggsMass
void SetPrimaryHiggsMass(double mass)
Definition: TEvtProb.cc:216
TEvtProb::XsecCalc_VX
double XsecCalc_VX(bool includeHiggsDecay)
Definition: TEvtProb.cc:899
selfDWpffcoupl
double selfDWpffcoupl[SIZE_Vpff][2]
Definition: raw_names.txt:20
TEvtProb::GetSelfDAZffCouplings
AZffCouplings * GetSelfDAZffCouplings()
Definition: TEvtProb.cc:361
ZZMatrixElement::set_CurrentCandidate
void set_CurrentCandidate(MELACandidate *cand)
Definition: ZZMatrixElement.cc:144
selfDGvpvpcoupl
double selfDGvpvpcoupl[SIZE_GVV][2]
Definition: raw_names.txt:31
selfDAZffcoupl
double selfDAZffcoupl[SIZE_AZff][2]
Definition: raw_names.txt:33
MELACandidate
Definition: MELACandidate.h:7
ZZMatrixElement::set_LeptonInterference
void set_LeptonInterference(TVar::LeptonInterference myLepInterf)
Definition: ZZMatrixElement.cc:135
SIZE_HVV
@ SIZE_HVV
Definition: raw_couplings.txt:57
TVar::DEBUG
@ DEBUG
Definition: TVar.hh:51
ZZMatrixElement::selfD_SpinOneCouplings
SpinOneCouplings * selfD_SpinOneCouplings
Definition: ZZMatrixElement.h:163
selfDHwpwpcoupl
double selfDHwpwpcoupl[SIZE_HVV][2]
Definition: raw_names.txt:19
ZZMatrixElement::reset_InputEvent
void reset_InputEvent()
Definition: ZZMatrixElement.cc:225
TEvtProb::SetCurrentCandidate
void SetCurrentCandidate(MELACandidate *cand)
Definition: TEvtProb.cc:303
selfDHwwcoupl
double selfDHwwcoupl[nSupportedHiggses][SIZE_HVV][2]
Definition: raw_names.txt:9
ZZMatrixElement::processProduction
TVar::Production processProduction
Definition: ZZMatrixElement.h:153
TEvtProb::GetXPropagator
double GetXPropagator(TVar::ResonancePropagatorScheme scheme)
Definition: TEvtProb.cc:1011
SpinZeroCouplings::SetHVVLambdaQ2
void SetHVVLambdaQ2(unsigned int gType, unsigned int index, double lambda, bool setWW=false, int whichResonance=1)
Definition: TCouplings.cc:132
selfDGvvpcoupl
double selfDGvvpcoupl[SIZE_GVV][2]
Definition: raw_names.txt:30
SpinTwoCouplings::SetGGGCouplings
void SetGGGCouplings(unsigned int index, double c_real, double c_imag)
Definition: TCouplings.cc:371
selfDHwwLambda_qsq
double selfDHwwLambda_qsq[nSupportedHiggses][SIZE_HVV_LAMBDAQSQ][SIZE_HVV_CQSQ]
Definition: raw_names.txt:11
TVar::bbH
@ bbH
Definition: TVar.hh:67
selfDZvvcoupl
double selfDZvvcoupl[SIZE_ZVV][2]
Definition: raw_names.txt:26
SIZE_GQQ
@ SIZE_GQQ
Definition: raw_couplings.txt:121
TEvtProb::ResetQuarkMasses
void ResetQuarkMasses()
Definition: TEvtProb.cc:316
ZZMatrixElement::selfD_SpinTwoCouplings
SpinTwoCouplings * selfD_SpinTwoCouplings
Definition: ZZMatrixElement.h:164
hto_betacom::beta1
real *8, dimension(3:6) beta1
Definition: CALLING_cpHTO.f:2080
selfDGvvcoupl
double selfDGvvcoupl[SIZE_GVV][2]
Definition: raw_names.txt:29
SpinOneCouplings::SetZQQCouplings
void SetZQQCouplings(unsigned int index, double c_real, double c_imag)
Definition: TCouplings.cc:303
ZZMatrixElement::set_PrimaryHiggsMass
void set_PrimaryHiggsMass(double mh)
Definition: ZZMatrixElement.cc:142
ZZMatrixElement::set_wHiggs
void set_wHiggs(double gah_, int index)
Definition: ZZMatrixElement.cc:184