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.
Mela.h
Go to the documentation of this file.
1 
11 /*
12 ************* HEADER: CMS MELA interface to MCFM/JHUGen-MELA *************
13 Please see the ../src/Mela.cc file for the instructions.
14 */
15 
16 #ifndef MELA_Mela_h
17 #define MELA_Mela_h
18 
19 #include <vector>
20 #include "TLorentzVector.h"
21 #include "TRandom3.h"
22 
23 
24 class TFile;
25 class TGraph;
26 class TH1F;
27 class TH2F;
28 class TH3F;
29 class RooRealVar;
30 class RooAbsPdf;
31 class RooArgSet;
33 class VectorPdfFactory;
34 class TensorPdfFactory;
36 class ZZMatrixElement;
37 class SuperMELA;
38 
39 #include "TVar.hh"
40 #include "TEvtProb.hh"
41 #include "MelaPConstant.h"
42 #include "SuperDijetMela.h"
43 #include "ScalarPdfFactory_HVV.h"
44 #include "VectorPdfFactory.h"
45 #include "TensorPdfFactory_ppHVV.h"
47 
48 class Mela{
49 
50 public:
58  Mela(double LHCsqrts_=13., double mh_=125., TVar::VerbosityLevel verbosity_=TVar::ERROR); // Higgs mass for supermela
59 
65  Mela(const Mela& other);
66 
70  ~Mela();
71 
77  void build(double mh_);
78 
87  void setProcess(TVar::Process myModel, TVar::MatrixElement myME, TVar::Production myProduction);
88 
95 
102 
109 
116  void setRemoveLeptonMasses(bool MasslessLeptonSwitch=true);
117 
124  void setRemoveJetMasses(bool MasslessLeptonSwitch=true);
125 
135  void setMelaPrimaryHiggsMass(double myHiggsMass);
136 
144  void setMelaHiggsMass(double myHiggsMass, int index=0);
145 
153  void setMelaHiggsWidth(double myHiggsWidth=-1, int index=0);
154 
163  void setMelaHiggsMassWidth(double myHiggsMass, double myHiggsWidth, int index);
164 
174  void setRenFacScaleMode(TVar::EventScaleScheme renormalizationSch, TVar::EventScaleScheme factorizationSch, double ren_sf, double fac_sf);
175 
183 
191  void setCurrentCandidateFromIndex(unsigned int icand); // Switches to another candidate
192 
199  void setCurrentCandidate(MELACandidate* cand); // Switches to another candidate
200 
211  void setInputEvent(
212  SimpleParticleCollection_t* pDaughters,
213  SimpleParticleCollection_t* pAssociated=0,
214  SimpleParticleCollection_t* pMothers=0,
215  bool isGen=false
216  ); // Adds another candidate
217 
224  void resetInputEvent(); // Reset the input candidates. Important to call in order to clean up TEvtProb!
225 
235  void setTempCandidate(
236  SimpleParticleCollection_t* pDaughters,
237  SimpleParticleCollection_t* pAssociated=0,
238  SimpleParticleCollection_t* pMothers=0,
239  bool isGen=false
240  ); // Adds a temp. candidate
241 
248  void appendTopCandidate(SimpleParticleCollection_t* TopDaughters); // Adds a top
249 
250  // Function to set EW parameters in MCFM/JHUGen
309  void resetMass(double inmass, int ipart);
310 
359  void resetWidth(double inwidth, int ipart);
360 
378  void resetQuarkMasses();
379 
384  void resetMCFM_EWKParameters(double ext_Gf, double ext_aemmz, double ext_mW, double ext_mZ, double ext_xW, int ext_ewscheme=3);
385 
386  // Function to get current primary EW/QCD parameters from MCFM/JHUGen (notice Higgs mass/width used in the ME could be different)
438  double getPrimaryMass(int ipart);
439 
482  double getPrimaryWidth(int ipart);
483 
492  double getHiggsWidthAtPoleMass(double mass);
493 
498  MelaIO* getIORecord(); // Full parton-by-parton ME record
499 
506 
512 
517  int getNCandidates();
518 
523  std::vector<MELATopCandidate_t*>* getTopCandidateCollection();
524 
545  void getConstant(float& prob); // <ME> constants
546 
547 
548  void getPAux(float& prob); // SuperProb
549 
567 
585  void computeDecayAngles(
586  float& qH,
587  float& m1,
588  float& m2,
589  float& costheta1,
590  float& costheta2,
591  float& Phi,
592  float& costhetastar,
593  float& Phi1
594  );
595 
612  void computeVBFAngles(
613  float& Q2V1,
614  float& Q2V2,
615  float& costheta1,
616  float& costheta2,
617  float& Phi,
618  float& costhetastar,
619  float& Phi1
620  );
622  float& Q2V1,
623  float& Q2V2,
624  float& costheta1_real, float& costheta1_imag,
625  float& costheta2_real, float& costheta2_imag,
626  float& Phi,
627  float& costhetastar,
628  float& Phi1
629  );
630 
647  void computeVHAngles(
648  float& mVstar,
649  float& mV,
650  float& costheta1,
651  float& costheta2,
652  float& Phi,
653  float& costhetastar,
654  float& Phi1
655  );
656  void computeTTHAngles(
657  int topDecay,
658 
659  float& mT1,
660  float& mW1,
661  float& mT2,
662  float& mW2,
663 
664  // TTH system
665  float& costheta1,
666  float& costheta2,
667  float& Phi,
668  float& costhetastar,
669  float& Phi1,
670 
671  // TT system
672  float& hbb,
673  float& hWW,
674  float& Phibb,
675  float& Phi1bb,
676 
677  // Wplus system
678  float& hWplusf,
679  float& PhiWplusf,
680 
681  // Wminus system
682  float& hWminusf,
683  float& PhiWminusf
684  );
685 
695  void computeP_selfDspin0(
696  double selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2],
697  float& prob,
698  bool useConstant=true
699  );
700 
711  void computeP_selfDspin1(
712  double selfDZqqcoupl_input[SIZE_ZQQ][2],
713  double selfDZvvcoupl_input[SIZE_ZVV][2],
714  float& prob,
715  bool useConstant=true
716  );
717 
727  void computeP_selfDspin1(
728  double selfDZvvcoupl_input[SIZE_ZVV][2],
729  float& prob,
730  bool useConstant=true
731  );
732 
744  void computeP_selfDspin2(
745  double selfDGggcoupl_input[SIZE_GGG][2],
746  double selfDGqqcoupl_input[SIZE_GQQ][2],
747  double selfDGvvcoupl_input[SIZE_GVV][2],
748  float& prob,
749  bool useConstant=true
750  );
751 
762  void computeP_selfDspin2(
763  double selfDGggcoupl_input[SIZE_GGG][2],
764  double selfDGvvcoupl_input[SIZE_GVV][2],
765  float& prob,
766  bool useConstant=true
767  );
768 
777  void computeP(
778  float& prob,
779  bool useConstant=true
780  );
781 
789  void computeD_CP(
790  TVar::MatrixElement myME,
791  TVar::Process myType,
792  float& prob
793  );
794 
795  //****VVH Spin-0****//
809  void computeProdDecP(
810  double selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2],
811  double selfDHwwcoupl_input[nSupportedHiggses][SIZE_HVV][2],
812  double selfDaTQGCcoupl_input[SIZE_ATQGC][2],
813  double selfDAZffcoupl_input[SIZE_AZff][2],
814  float& prob,
815  bool useConstant=true
816  );
817 
849  void computeProdDecP(
850  float& prob,
851  bool useConstant=true
852  );
853 
854  //****HJ/HJJ/VBF Spin-0****//
866  void computeProdP(
867  double selfDHggcoupl_input[SIZE_HGG][2],
868  double selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2],
869  double selfDHwwcoupl_input[nSupportedHiggses][SIZE_HVV][2],
870  float& prob,
871  bool useConstant=true
872  );
873 
880  void computeProdP(
881  float& prob,
882  bool useConstant=true
883  );
884 
885  //****VH Spin-0****//
886  void computeProdP_VH(
887  double selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2],
888  float& prob,
889  bool includeHiggsDecay=false,
890  bool useConstant=true
891  );
892  void computeProdP_VH(
893  float& prob,
894  bool includeHiggsDecay=false,
895  bool useConstant=true
896  );
897 
898  //***ttH Spin-0****//
899  void computeProdP_ttH(
900  float& prob,
901  int topProcess=2,
902  int topDecay=0,
903  bool useConstant=true
904  );
905 
906  // Calculation weight to correct for fermion interference
907  void compute4FermionWeight(float& w);
908 
909  // Calculation of X propagator
910  void getXPropagator(TVar::ResonancePropagatorScheme scheme, float& prop);
911 
912  //*** SuperMela ***//
913  void computePM4l(
914  TVar::SuperMelaSyst syst,
915  float& prob
916  );
917 
918  //*** SuperJJMela ***//
919  void computeDijetConvBW(float& prob, bool useTrueBW=false);
920 
921  //*** Dgg10 ***//
922  void computeD_gg(
923  TVar::MatrixElement myME,
924  TVar::Process myType,
925  float& prob
926  );
927 
928  // Access ZZMEs Calculate4Momentum
929  std::vector<TLorentzVector> calculate4Momentum(double Mx, double M1, double M2, double theta, double theta1, double theta2, double Phi1, double Phi);
930 
931  /********************/
932  /*** Data members ***/
933  /********************/
934  TRandom3 melaRandomNumber; // Used in SuperMELA smearing
935  RooRealVar* mzz_rrv;
936  RooRealVar* z1mass_rrv;
937  RooRealVar* z2mass_rrv;
938  RooRealVar* costhetastar_rrv;
939  RooRealVar* costheta1_rrv;
940  RooRealVar* costheta2_rrv;
941  RooRealVar* phi_rrv;
942  RooRealVar* phi1_rrv;
943  RooRealVar* Y_rrv;
944  RooRealVar* upFrac_rrv;
945 
946  RooAbsPdf* pdf;
951 
953 
954  // Self-define arrays are now members of MELA.
955  // There are a lot of them!
956  //****Spin-0****//
957  // The first dimension (of size [nSupportedHiggses=2]) supports a second resonance present in MCFM
958 
991  //****Spin-1****//
994  //****Spin-2****//
1000  //****aTQGC****//
1002  //****Anomnalous Zff**//
1004  // That is a lot of them!
1005 
1006  static void cleanLinkedFiles();
1007 
1008 protected:
1009  /********************/
1010  /*** Data members ***/
1011  /********************/
1012  double LHCsqrts;
1018 
1021 
1022 
1024 
1025  MELACandidate* melaCand; // Pointer to persistent TEvtProb object
1026 
1027  /***** ME CONSTANT HANDLES *****/
1028  // Constants that vary with sqrts due to application of PDFs
1029  //
1031  //
1033  //
1035  //
1037  //
1039  // Decay ME constants that do not use PDFs
1040  //
1044  //
1048  //
1052  //
1056  //
1060  //
1064  //
1068  //
1072  //
1076  //
1080  //
1084  //
1086 
1087 
1088  /*****************/
1089  /*** Functions ***/
1090  /*****************/
1091  void printLogo() const;
1092 
1093  void setSpinZeroCouplings();
1094  void setSpinOneCouplings();
1095  void setSpinTwoCouplings();
1096  void setATQGCCouplings();
1097  void setAZffCouplings();
1098 
1099  bool configureAnalyticalPDFs();
1100  void reset_SelfDCouplings();
1101  void reset_PAux(); // SuperProb reset
1102  void reset_CandRef();
1103 
1104  void constructDggr(
1105  float bkg_VAMCFM_noscale,
1106  float ggzz_VAMCFM_noscale,
1107  float ggHZZ_prob_pure_noscale,
1108  float ggHZZ_prob_int_noscale,
1109  float widthScale,
1110  float& myDggr
1111  );
1112 
1113  void getPConstantHandles();
1114  void deletePConstantHandles();
1116  TVar::MatrixElement me_,
1117  TVar::Production prod_,
1118  TVar::Process proc_,
1119  TString relpath,
1120  TString spname,
1121  const bool useSqrts=false
1122  );
1123  void deletePConstantHandle(MelaPConstant*& handle);
1124  void computeConstant(float& prob);
1125  void setConstant();
1127  float getConstant_4l();
1128  float getConstant_2l2q();
1129  float getConstant_4q();
1130  float getConstant_FourFermionDecay(const int& decid);
1131 
1132 };
1133 
1134 #endif
1135 
TVar::ERROR
@ ERROR
Definition: TVar.hh:49
Mela::pAvgSmooth_MCFM_Had_ZH_S_HSMHiggs_4e
MelaPConstant * pAvgSmooth_MCFM_Had_ZH_S_HSMHiggs_4e
Definition: Mela.h:1054
Mela::pAvgSmooth_MCFM_Had_ZH_bkgZZ_4e
MelaPConstant * pAvgSmooth_MCFM_Had_ZH_bkgZZ_4e
Definition: Mela.h:1074
Mela::getConstant
void getConstant(float &prob)
This function returns a multiplicative constant to the matrix element calculation.
Definition: Mela.cc:2168
Mela::computeP_selfDspin2
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...
Definition: Mela.cc:1146
Mela::getConstant_FourFermionDecay
float getConstant_FourFermionDecay(const int &decid)
Definition: Mela.cc:2313
Mela::pAvgSmooth_MCFM_ZZQQB_bkgZZ_4mu
MelaPConstant * pAvgSmooth_MCFM_ZZQQB_bkgZZ_4mu
Definition: Mela.h:1065
Mela::setAZffCouplings
void setAZffCouplings()
Definition: Mela.cc:418
MelaIO
Definition: MelaIO.h:8
Mela::pAvgSmooth_JHUGen_JQCD_HSMHiggs
MelaPConstant * pAvgSmooth_JHUGen_JQCD_HSMHiggs[TVar::nFermionMassRemovalSchemes-1]
Definition: Mela.h:1030
TVar::LeptonInterference
LeptonInterference
Definition: TVar.hh:106
SIZE_GVV
@ SIZE_GVV
Definition: raw_couplings.txt:158
Mela::pAvgSmooth_MCFM_JJQCD_bkgZZ_4e
MelaPConstant * pAvgSmooth_MCFM_JJQCD_bkgZZ_4e
Definition: Mela.h:1082
Mela::setInputEvent
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.
Definition: Mela.cc:343
Mela::selfDGqqcoupl
double selfDGqqcoupl[SIZE_GQQ][2]
Definition: Mela.h:995
Mela::selfDGvpvpcoupl
double selfDGvpvpcoupl[SIZE_GVV][2]
Definition: Mela.h:999
Mela::Y_rrv
RooRealVar * Y_rrv
Definition: Mela.h:943
TVar::VerbosityLevel
VerbosityLevel
Definition: TVar.hh:47
Mela::getPConstantHandles
void getPConstantHandles()
Definition: Mela.cc:2470
RooqqZZ_JHU_ZgammaZZ_fast
Definition: RooqqZZ_JHU_ZgammaZZ_fast.h:18
Mela::build
void build(double mh_)
This is the actual building of the tool that occurs in each instance of the Mela::Mela constructor.
Definition: Mela.cc:175
Mela::resetMCFM_EWKParameters
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.
Definition: Mela.cc:511
Mela::getPrimaryWidth
double getPrimaryWidth(int ipart)
A function to get the current primary EW/QCD parameters from MELA.
Definition: Mela.cc:516
Mela::getNCandidates
int getNCandidates()
Returns the size of the candidate list TEvtProb::candList.
Definition: Mela.cc:548
SIZE_HVV_CQSQ
@ SIZE_HVV_CQSQ
Definition: raw_couplings.txt:74
Mela::computeDecayAngles
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
Definition: Mela.cc:558
ScalarPdfFactory_HVV.h
Mela::selfDM_Zprime
double selfDM_Zprime
Definition: Mela.h:987
Mela::getPConstantHandle
MelaPConstant * getPConstantHandle(TVar::MatrixElement me_, TVar::Production prod_, TVar::Process proc_, TString relpath, TString spname, const bool useSqrts=false)
Definition: Mela.cc:2665
Mela::computeProdP_VH
void computeProdP_VH(double selfDHvvcoupl_input[nSupportedHiggses][SIZE_HVV][2], float &prob, bool includeHiggsDecay=false, bool useConstant=true)
Definition: Mela.cc:1730
MelaPConstant.h
Mela::pAvgSmooth_MCFM_Had_ZH_bkgZZ_2mu2e
MelaPConstant * pAvgSmooth_MCFM_Had_ZH_bkgZZ_2mu2e
Definition: Mela.h:1075
Mela::compute4FermionWeight
void compute4FermionWeight(float &w)
Definition: Mela.cc:1807
Mela::pAvgSmooth_MCFM_ZZGG_HSMHiggs_4e
MelaPConstant * pAvgSmooth_MCFM_ZZGG_HSMHiggs_4e
Definition: Mela.h:1046
Mela::selfDHwwCLambda_qsq
int selfDHwwCLambda_qsq[nSupportedHiggses][SIZE_HVV_CQSQ]
Definition: Mela.h:979
Mela::pAvgSmooth_MCFM_JJVBF_S_HSMHiggs_4mu
MelaPConstant * pAvgSmooth_MCFM_JJVBF_S_HSMHiggs_4mu
Definition: Mela.h:1049
TVar::nFermionMassRemovalSchemes
@ nFermionMassRemovalSchemes
Definition: TVar.hh:115
Mela::auxiliaryProb
float auxiliaryProb
Definition: Mela.h:1023
SIZE_GGG
@ SIZE_GGG
Definition: raw_couplings.txt:131
Mela::getIORecord
MelaIO * getIORecord()
Returns the MELAIO object, and by consequence, the entire parton-by-parton matrix element record.
Definition: Mela.cc:544
Mela::pAvgSmooth_MCFM_Had_WH_bkgZZ_4mu
MelaPConstant * pAvgSmooth_MCFM_Had_WH_bkgZZ_4mu
Definition: Mela.h:1077
SuperDijetMela
Definition: SuperDijetMela.h:8
ScalarPdfFactory_HVV
Definition: ScalarPdfFactory_HVV.h:8
Mela::pAvgSmooth_MCFM_Had_WH_S_HSMHiggs_2mu2e
MelaPConstant * pAvgSmooth_MCFM_Had_WH_S_HSMHiggs_2mu2e
Definition: Mela.h:1059
Mela::setCurrentCandidateFromIndex
void setCurrentCandidateFromIndex(unsigned int icand)
Switches the candidate that you are working on to another candidate based off of an index.
Definition: Mela.cc:341
SuperMELA
Definition: SuperMELA.h:17
Mela::pAvgSmooth_MCFM_ZZGG_bkgZZ_4mu
MelaPConstant * pAvgSmooth_MCFM_ZZGG_bkgZZ_4mu
Definition: Mela.h:1061
Mela::z1mass_rrv
RooRealVar * z1mass_rrv
Definition: Mela.h:936
Mela::selfDHttcoupl
double selfDHttcoupl[nSupportedHiggses][SIZE_HQQ][2]
Definition: Mela.h:971
Mela::pAvgSmooth_MCFM_JJQCD_bkgZZ_2mu2e
MelaPConstant * pAvgSmooth_MCFM_JJQCD_bkgZZ_2mu2e
Definition: Mela.h:1083
Mela::setRemoveLeptonMasses
void setRemoveLeptonMasses(bool MasslessLeptonSwitch=true)
either permits or forbids massive leptons.
Definition: Mela.cc:519
Mela::getPAux
void getPAux(float &prob)
Definition: Mela.cc:554
Mela::getMeasurablesRRV
RooSpin::modelMeasurables getMeasurablesRRV()
Returns a RooSpin::modelMeasureables object containing many observable quantities.
Definition: Mela.cc:528
Mela::setVerbosity
void setVerbosity(TVar::VerbosityLevel verbosity_=TVar::ERROR)
Sets the verbosity for MELA outside of the initial constructor.
Definition: Mela.cc:325
Mela::setCurrentCandidate
void setCurrentCandidate(MELACandidate *cand)
Switches the candidate that you are working on to another candidate object specified.
Definition: Mela.cc:342
TVar::CandidateDecayMode
CandidateDecayMode
Definition: TVar.hh:37
Mela::resetWidth
void resetWidth(double inwidth, int ipart)
Resets the width for a particle that is an electroweak parameter according to its id.
Definition: Mela.cc:509
Mela::mzz_rrv
RooRealVar * mzz_rrv
Definition: Mela.h:935
Mela::setMelaHiggsWidth
void setMelaHiggsWidth(double myHiggsWidth=-1, int index=0)
Sets the width of your chosen Higgs.
Definition: Mela.cc:337
Mela::setSpinZeroCouplings
void setSpinZeroCouplings()
Definition: Mela.cc:366
Mela::getTopCandidateCollection
std::vector< MELATopCandidate_t * > * getTopCandidateCollection()
Same as getNCandidates but specifically for Top Quark Candidates.
Definition: Mela.cc:549
Mela::computeD_CP
void computeD_CP(TVar::MatrixElement myME, TVar::Process myType, float &prob)
computes the value of D_CP
Definition: Mela.cc:1339
Mela::pAvgSmooth_MCFM_Had_WH_S_HSMHiggs_4mu
MelaPConstant * pAvgSmooth_MCFM_Had_WH_S_HSMHiggs_4mu
Definition: Mela.h:1057
Mela::pAvgSmooth_JHUGen_Had_WH_HSMHiggs
MelaPConstant * pAvgSmooth_JHUGen_Had_WH_HSMHiggs[TVar::nFermionMassRemovalSchemes-1]
Definition: Mela.h:1038
Mela::computeProdP
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
Definition: Mela.cc:1501
Mela::pAvgSmooth_MCFM_ZZQQB_bkgZZ_4e
MelaPConstant * pAvgSmooth_MCFM_ZZQQB_bkgZZ_4e
Definition: Mela.h:1066
Mela::selfDHzpzpcoupl
double selfDHzpzpcoupl[SIZE_HVV][2]
Definition: Mela.h:982
SIZE_ZQQ
@ SIZE_ZQQ
Definition: raw_couplings.txt:107
Mela::setRenFacScaleMode
void setRenFacScaleMode(TVar::EventScaleScheme renormalizationSch, TVar::EventScaleScheme factorizationSch, double ren_sf, double fac_sf)
Sets the renormalization and the factorization schemes.
Definition: Mela.cc:521
TVar::Process
Process
Definition: TVar.hh:125
Mela::selfDGa_Zprime
double selfDGa_Zprime
Definition: Mela.h:988
TVar::EventScaleScheme
EventScaleScheme
Definition: TVar.hh:196
Mela::selfDWpffcoupl
double selfDWpffcoupl[SIZE_Vpff][2]
Definition: Mela.h:986
Mela::getConstant_4l
float getConstant_4l()
Definition: Mela.cc:2286
Mela::reset_SelfDCouplings
void reset_SelfDCouplings()
Definition: Mela.cc:426
Mela::selfDZpffcoupl
double selfDZpffcoupl[SIZE_Vpff][2]
Definition: Mela.h:983
Mela::reset_CandRef
void reset_CandRef()
Definition: Mela.cc:550
Mela::selfDaTQGCcoupl
double selfDaTQGCcoupl[SIZE_ATQGC][2]
Definition: Mela.h:1001
Mela::getConstant_JHUGenUndecayed
float getConstant_JHUGenUndecayed()
Definition: Mela.cc:2251
VectorPdfFactory
Definition: VectorPdfFactory.h:12
Mela::pAvgSmooth_JHUGen_ZZGG_HSMHiggs_2mu2e
MelaPConstant * pAvgSmooth_JHUGen_ZZGG_HSMHiggs_2mu2e
Definition: Mela.h:1043
Mela::computeVBFAngles
void computeVBFAngles(float &Q2V1, float &Q2V2, float &costheta1, float &costheta2, float &Phi, float &costhetastar, float &Phi1)
computes the decay angles for Vector Boson Fusion (VBF) -> H -> ZZ as defined in Figure 1 of arXiv:12...
Definition: Mela.cc:620
Mela::setProcess
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....
Definition: Mela.cc:310
Mela::Mela
Mela(double LHCsqrts_=13., double mh_=125., TVar::VerbosityLevel verbosity_=TVar::ERROR)
the MELA constructor
Definition: Mela.cc:94
Mela::qqZZmodel
RooqqZZ_JHU_ZgammaZZ_fast * qqZZmodel
Definition: Mela.h:950
Mela::resetQuarkMasses
void resetQuarkMasses()
Resets the masses of each quark to their original values.
Definition: Mela.cc:510
Mela::pAvgSmooth_MCFM_Had_ZH_S_HSMHiggs_4mu
MelaPConstant * pAvgSmooth_MCFM_Had_ZH_S_HSMHiggs_4mu
Definition: Mela.h:1053
Mela::selfDHqqcoupl
double selfDHqqcoupl[nSupportedHiggses][SIZE_HQQ][2]
Definition: Mela.h:969
Mela::setConstant
void setConstant()
Definition: Mela.cc:2175
Mela::setTempCandidate
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 pushin...
Definition: Mela.cc:357
SIZE_ZVV
@ SIZE_ZVV
Definition: raw_couplings.txt:114
Mela::deletePConstantHandle
void deletePConstantHandle(MelaPConstant *&handle)
Definition: Mela.cc:2785
Mela::ggSpin0Model
ScalarPdfFactory_HVV * ggSpin0Model
Definition: Mela.h:947
Mela::costheta1_rrv
RooRealVar * costheta1_rrv
Definition: Mela.h:939
Mela::pAvgSmooth_MCFM_ZZGG_bkgZZ_2mu2e
MelaPConstant * pAvgSmooth_MCFM_ZZGG_bkgZZ_2mu2e
Definition: Mela.h:1063
TEvtProb.hh
Mela::deletePConstantHandles
void deletePConstantHandles()
Definition: Mela.cc:2729
Mela::selfDGvvcoupl
double selfDGvvcoupl[SIZE_GVV][2]
Definition: Mela.h:997
Mela::computeP_selfDspin0
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...
Definition: Mela.cc:1077
nSupportedHiggses
@ nSupportedHiggses
Definition: TMCFM.hh:15
Mela::pAvgSmooth_MCFM_ZZQQB_bkgZZ_2mu2e
MelaPConstant * pAvgSmooth_MCFM_ZZQQB_bkgZZ_2mu2e
Definition: Mela.h:1067
TVar::DefaultLeptonInterf
@ DefaultLeptonInterf
Definition: TVar.hh:107
Mela::pAvgSmooth_MCFM_Had_WH_bkgZZ_2mu2e
MelaPConstant * pAvgSmooth_MCFM_Had_WH_bkgZZ_2mu2e
Definition: Mela.h:1079
TensorPdfFactory_ppHVV.h
Mela::differentiate_HWW_HZZ
bool differentiate_HWW_HZZ
Definition: Mela.h:980
Mela::melaCand
MELACandidate * melaCand
Definition: Mela.h:1025
SIZE_HQQ
@ SIZE_HQQ
Definition: raw_couplings.txt:5
Mela::getHiggsWidthAtPoleMass
double getHiggsWidthAtPoleMass(double mass)
Returns the width of the Higgs at a given pole mass as a calculation.
Definition: Mela.cc:517
Mela::getConstant_2l2q
float getConstant_2l2q()
Definition: Mela.cc:2297
Mela::computeTTHAngles
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)
Definition: Mela.cc:895
anonymous_namespace{TCouplingsBase.hh}::SIZE_HGG
@ SIZE_HGG
Definition: TCouplingsBase.hh:40
Mela::getConstant_4q
float getConstant_4q()
Definition: Mela.cc:2306
Mela::pAvgSmooth_JHUGen_ZZGG_HSMHiggs_4e
MelaPConstant * pAvgSmooth_JHUGen_ZZGG_HSMHiggs_4e
Definition: Mela.h:1042
Mela::myModel_
TVar::Process myModel_
Definition: Mela.h:1013
Mela::cleanLinkedFiles
static void cleanLinkedFiles()
Definition: Mela.cc:162
Mela::myME_
TVar::MatrixElement myME_
Definition: Mela.h:1014
Mela::computeVHAngles
void computeVHAngles(float &mVstar, float &mV, float &costheta1, float &costheta2, float &Phi, float &costhetastar, float &Phi1)
computes the decay angles for Vector Boson Fusion (VBF) -> H -> ZZ as defined in Figure 1 of arXiv:12...
Definition: Mela.cc:786
TensorPdfFactory_ppHVV
Definition: TensorPdfFactory_ppHVV.h:8
Mela::getVerbosity
TVar::VerbosityLevel getVerbosity()
Gets the current verbosity level for MELA.
Definition: Mela.cc:333
Mela::selfDHg4g4coupl
double selfDHg4g4coupl[nSupportedHiggses][SIZE_HGG][2]
Definition: Mela.h:967
Mela::computeP
void computeP(float &prob, bool useConstant=true)
Computes the probability for the probabilities on the decay side of things using the constituent daug...
Definition: Mela.cc:1183
SuperDijetMela.h
Mela::pAvgSmooth_MCFM_Had_WH_bkgZZ_4e
MelaPConstant * pAvgSmooth_MCFM_Had_WH_bkgZZ_4e
Definition: Mela.h:1078
Mela::ZZME
ZZMatrixElement * ZZME
Definition: Mela.h:1019
Mela::computePM4l
void computePM4l(TVar::SuperMelaSyst syst, float &prob)
Definition: Mela.cc:1864
Mela::setCandidateDecayMode
void setCandidateDecayMode(TVar::CandidateDecayMode mode)
Sets the decay mode for your event.
Definition: Mela.cc:340
VectorPdfFactory.h
Mela::pAvgSmooth_MCFM_Had_ZH_S_HSMHiggs_2mu2e
MelaPConstant * pAvgSmooth_MCFM_Had_ZH_S_HSMHiggs_2mu2e
Definition: Mela.h:1055
Mela::setMelaHiggsMass
void setMelaHiggsMass(double myHiggsMass, int index=0)
Sets the mass of your chosen Higgs.
Definition: Mela.cc:336
Mela::LHCsqrts
double LHCsqrts
Definition: Mela.h:1012
Mela::spin2Model
TensorPdfFactory_ppHVV * spin2Model
Definition: Mela.h:949
Mela::getCurrentCandidateIndex
int getCurrentCandidateIndex()
Returns the index of the current MELA candidate - returns -1 if there is no candidate to be found.
Definition: Mela.cc:547
Mela::pAvgSmooth_JHUGen_ZZGG_HSMHiggs_4mu
MelaPConstant * pAvgSmooth_JHUGen_ZZGG_HSMHiggs_4mu
Definition: Mela.h:1041
Mela::selfDZqqcoupl
double selfDZqqcoupl[SIZE_ZQQ][2]
Definition: Mela.h:992
Mela::selfDHb4b4coupl
double selfDHb4b4coupl[nSupportedHiggses][SIZE_HQQ][2]
Definition: Mela.h:972
Mela::selfDHwpwpcoupl
double selfDHwpwpcoupl[SIZE_HVV][2]
Definition: Mela.h:985
Mela::pAvgSmooth_JHUGen_JJQCD_HSMHiggs
MelaPConstant * pAvgSmooth_JHUGen_JJQCD_HSMHiggs[TVar::nFermionMassRemovalSchemes-1]
Definition: Mela.h:1032
Mela::printLogo
void printLogo() const
Definition: Mela.cc:278
Mela::pAvgSmooth_MCFM_JJVBF_S_HSMHiggs_4e
MelaPConstant * pAvgSmooth_MCFM_JJVBF_S_HSMHiggs_4e
Definition: Mela.h:1050
SIZE_Vpff
@ SIZE_Vpff
Definition: raw_couplings.txt:100
Mela::selfDHzzLambda_qsq
double selfDHzzLambda_qsq[nSupportedHiggses][SIZE_HVV_LAMBDAQSQ][SIZE_HVV_CQSQ]
Definition: Mela.h:976
Mela::phi_rrv
RooRealVar * phi_rrv
Definition: Mela.h:941
Mela::getXPropagator
void getXPropagator(TVar::ResonancePropagatorScheme scheme, float &prop)
Definition: Mela.cc:1799
Mela::~Mela
~Mela()
MELA destructor.
Definition: Mela.cc:122
Mela::computeProdDecP
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
Definition: Mela.cc:1430
Mela::selfDHwwcoupl
double selfDHwwcoupl[nSupportedHiggses][SIZE_HVV][2]
Definition: Mela.h:975
Mela::pAvgSmooth_MCFM_JJVBF_bkgZZ_4mu
MelaPConstant * pAvgSmooth_MCFM_JJVBF_bkgZZ_4mu
Definition: Mela.h:1069
Mela::constructDggr
void constructDggr(float bkg_VAMCFM_noscale, float ggzz_VAMCFM_noscale, float ggHZZ_prob_pure_noscale, float ggHZZ_prob_int_noscale, float widthScale, float &myDggr)
Definition: Mela.cc:1920
Mela::pAvgSmooth_MCFM_ZZGG_bkgZZ_4e
MelaPConstant * pAvgSmooth_MCFM_ZZGG_bkgZZ_4e
Definition: Mela.h:1062
Mela::setRemoveJetMasses
void setRemoveJetMasses(bool MasslessLeptonSwitch=true)
either permits or forbids massive jets.
Definition: Mela.cc:520
Mela::superDijet
SuperDijetMela * superDijet
Definition: Mela.h:1020
Mela::selfDGvvpcoupl
double selfDGvvpcoupl[SIZE_GVV][2]
Definition: Mela.h:998
Mela::resetMass
void resetMass(double inmass, int ipart)
Resets the mass for a particle that is an electroweak parameter according to its id.
Definition: Mela.cc:508
Mela::resetInputEvent
void resetInputEvent()
Resets the event in preparation for the next iteration of the event loop.
Definition: Mela.cc:356
TVar::MatrixElement
MatrixElement
Definition: TVar.hh:55
Mela::pdf
RooAbsPdf * pdf
Definition: Mela.h:946
SIZE_ATQGC
@ SIZE_ATQGC
Definition: raw_couplings.txt:176
Mela::selfDAZffcoupl
double selfDAZffcoupl[SIZE_AZff][2]
Definition: Mela.h:1003
Mela::setATQGCCouplings
void setATQGCCouplings()
Definition: Mela.cc:413
Mela::selfDHzzCLambda_qsq
int selfDHzzCLambda_qsq[nSupportedHiggses][SIZE_HVV_CQSQ]
Definition: Mela.h:978
Mela::configureAnalyticalPDFs
bool configureAnalyticalPDFs()
Definition: Mela.cc:1969
Mela::pAvgSmooth_MCFM_JJVBF_bkgZZ_4e
MelaPConstant * pAvgSmooth_MCFM_JJVBF_bkgZZ_4e
Definition: Mela.h:1070
Mela::setSpinOneCouplings
void setSpinOneCouplings()
Definition: Mela.cc:398
Mela::myProduction_
TVar::Production myProduction_
Definition: Mela.h:1015
Mela::setMelaHiggsMassWidth
void setMelaHiggsMassWidth(double myHiggsMass, double myHiggsWidth, int index)
a combination of setMelaHiggsMass and setMelaHiggsWidth.
Definition: Mela.cc:338
SIZE_HVV_LAMBDAQSQ
@ SIZE_HVV_LAMBDAQSQ
Definition: raw_couplings.txt:66
Mela::selfDZvvcoupl
double selfDZvvcoupl[SIZE_ZVV][2]
Definition: Mela.h:993
Mela::pAvgSmooth_JHUGen_JJVBF_HSMHiggs
MelaPConstant * pAvgSmooth_JHUGen_JJVBF_HSMHiggs[TVar::nFermionMassRemovalSchemes-1]
Definition: Mela.h:1034
Mela::computeVBFAngles_ComplexBoost
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)
Definition: Mela.cc:702
Mela
Definition: Mela.h:48
Mela::selfDM_Wprime
double selfDM_Wprime
Definition: Mela.h:989
TVar::ResonancePropagatorScheme
ResonancePropagatorScheme
Definition: TVar.hh:117
Mela::reset_PAux
void reset_PAux()
Definition: Mela.cc:555
Mela::melaRandomNumber
TRandom3 melaRandomNumber
Definition: Mela.h:934
Mela::computeProdP_ttH
void computeProdP_ttH(float &prob, int topProcess=2, int topDecay=0, bool useConstant=true)
Definition: Mela.cc:1778
Mela::computeD_gg
void computeD_gg(TVar::MatrixElement myME, TVar::Process myType, float &prob)
Definition: Mela.cc:1934
Mela::computeConstant
void computeConstant(float &prob)
Definition: Mela.cc:2169
Mela::super
SuperMELA * super
Definition: Mela.h:952
Mela::selfDGa_Wprime
double selfDGa_Wprime
Definition: Mela.h:990
Mela::selfDHt4t4coupl
double selfDHt4t4coupl[nSupportedHiggses][SIZE_HQQ][2]
Definition: Mela.h:973
Mela::selfDHwwLambda_qsq
double selfDHwwLambda_qsq[nSupportedHiggses][SIZE_HVV_LAMBDAQSQ][SIZE_HVV_CQSQ]
Definition: Mela.h:977
Mela::setMelaPrimaryHiggsMass
void setMelaPrimaryHiggsMass(double myHiggsMass)
Sets the mass of the "primary" higgs.
Definition: Mela.cc:335
Mela::z2mass_rrv
RooRealVar * z2mass_rrv
Definition: Mela.h:937
Mela::pAvgSmooth_MCFM_ZZGG_HSMHiggs_2mu2e
MelaPConstant * pAvgSmooth_MCFM_ZZGG_HSMHiggs_2mu2e
Definition: Mela.h:1047
SIZE_AZff
@ SIZE_AZff
Definition: raw_couplings.txt:194
Mela::computeDijetConvBW
void computeDijetConvBW(float &prob, bool useTrueBW=false)
Definition: Mela.cc:2792
Mela::costhetastar_rrv
RooRealVar * costhetastar_rrv
Definition: Mela.h:938
Mela::pAvgSmooth_MCFM_JJQCD_bkgZZ_4mu
MelaPConstant * pAvgSmooth_MCFM_JJQCD_bkgZZ_4mu
Definition: Mela.h:1081
TVar::SuperMelaSyst
SuperMelaSyst
Definition: TVar.hh:186
Mela::costheta2_rrv
RooRealVar * costheta2_rrv
Definition: Mela.h:940
Mela::selfDHwwpcoupl
double selfDHwwpcoupl[SIZE_HVV][2]
Definition: Mela.h:984
Mela::myVerbosity_
TVar::VerbosityLevel myVerbosity_
Definition: Mela.h:1017
RooqqZZ_JHU_ZgammaZZ_fast.h
Mela::pAvgSmooth_MCFM_ZZGG_HSMHiggs_4mu
MelaPConstant * pAvgSmooth_MCFM_ZZGG_HSMHiggs_4mu
Definition: Mela.h:1045
ZZMatrixElement
Definition: ZZMatrixElement.h:9
TVar.hh
Mela::selfDHzzcoupl
double selfDHzzcoupl[nSupportedHiggses][SIZE_HVV][2]
Definition: Mela.h:974
Mela::calculate4Momentum
std::vector< TLorentzVector > calculate4Momentum(double Mx, double M1, double M2, double theta, double theta1, double theta2, double Phi1, double Phi)
Definition: Mela.cc:524
Mela::computeP_selfDspin1
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...
Definition: Mela.cc:1112
TensorPdfFactory
Definition: TensorPdfFactory.h:10
Mela::pAvgSmooth_MCFM_JJQCD_bkgZJets_2l2q
MelaPConstant * pAvgSmooth_MCFM_JJQCD_bkgZJets_2l2q
Definition: Mela.h:1085
Mela::pAvgSmooth_MCFM_Had_ZH_bkgZZ_4mu
MelaPConstant * pAvgSmooth_MCFM_Had_ZH_bkgZZ_4mu
Definition: Mela.h:1073
Mela::phi1_rrv
RooRealVar * phi1_rrv
Definition: Mela.h:942
Mela::spin1Model
VectorPdfFactory * spin1Model
Definition: Mela.h:948
Mela::selfDHzzpcoupl
double selfDHzzpcoupl[SIZE_HVV][2]
Definition: Mela.h:981
MELACandidate
Definition: MELACandidate.h:7
SIZE_HVV
@ SIZE_HVV
Definition: raw_couplings.txt:57
Mela::pAvgSmooth_MCFM_Had_WH_S_HSMHiggs_4e
MelaPConstant * pAvgSmooth_MCFM_Had_WH_S_HSMHiggs_4e
Definition: Mela.h:1058
Mela::setSpinTwoCouplings
void setSpinTwoCouplings()
Definition: Mela.cc:401
SimpleParticleCollection_t
std::vector< SimpleParticle_t > SimpleParticleCollection_t
Definition: TVar.hh:25
MelaPConstant
Definition: MelaPConstant.h:9
Mela::selfDHggcoupl
double selfDHggcoupl[nSupportedHiggses][SIZE_HGG][2]
Definition: Mela.h:965
Mela::setMelaLeptonInterference
void setMelaLeptonInterference(TVar::LeptonInterference myLepInterf=TVar::DefaultLeptonInterf)
Sets the MELA Lepton Interference.
Definition: Mela.cc:339
Mela::pAvgSmooth_MCFM_JJVBF_S_HSMHiggs_2mu2e
MelaPConstant * pAvgSmooth_MCFM_JJVBF_S_HSMHiggs_2mu2e
Definition: Mela.h:1051
Mela::appendTopCandidate
void appendTopCandidate(SimpleParticleCollection_t *TopDaughters)
Adds a top quark as a MELA candidate.
Definition: Mela.cc:363
RooSpin::modelMeasurables
Definition: RooSpin.h:50
Mela::upFrac_rrv
RooRealVar * upFrac_rrv
Definition: Mela.h:944
Mela::selfDHbbcoupl
double selfDHbbcoupl[nSupportedHiggses][SIZE_HQQ][2]
Definition: Mela.h:970
Mela::myLepInterf_
TVar::LeptonInterference myLepInterf_
Definition: Mela.h:1016
TVar::Production
Production
Definition: TVar.hh:60
Mela::selfDGggcoupl
double selfDGggcoupl[SIZE_GGG][2]
Definition: Mela.h:996
Mela::pAvgSmooth_MCFM_JJVBF_bkgZZ_2mu2e
MelaPConstant * pAvgSmooth_MCFM_JJVBF_bkgZZ_2mu2e
Definition: Mela.h:1071
Mela::getCurrentCandidate
MELACandidate * getCurrentCandidate()
Gets the current MELA top-level (input) candList object.
Definition: Mela.cc:546
SIZE_GQQ
@ SIZE_GQQ
Definition: raw_couplings.txt:121
Mela::getPrimaryMass
double getPrimaryMass(int ipart)
A function to get the current primary EW/QCD parameters from MELA.
Definition: Mela.cc:515
Mela::pAvgSmooth_JHUGen_Had_ZH_HSMHiggs
MelaPConstant * pAvgSmooth_JHUGen_Had_ZH_HSMHiggs[TVar::nFermionMassRemovalSchemes-1]
Definition: Mela.h:1036