TauolappInterface

Macros

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97
#ifndef gen_TauolaInterface_TauolappInterface_h
#define gen_TauolaInterface_TauolappInterface_h

#include "HepPDT/ParticleDataTable.hh"

#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/Framework/interface/ConsumesCollector.h"
#include "GeneratorInterface/TauolaInterface/interface/TauolaInterfaceBase.h"
#include "TLorentzVector.h"
#include "TVector.h"

namespace HepMC {
  class GenEvent;
}

namespace CLHEP {
  class HepRandomEngine;
}

namespace gen {
  extern "C" {
  void ranmar_(float* rvec, int* lenv);
  void rmarin_(int*, int*, int*);
  }

  class TauolappInterface : public TauolaInterfaceBase {
  public:
    // ctor & dtor
    TauolappInterface(const edm::ParameterSet&, edm::ConsumesCollector);
    ~TauolappInterface() override;

    void enablePolarization() override {
      fPolarization = true;
      return;
    }
    void disablePolarization() override {
      fPolarization = false;
      return;
    }
    void init(const edm::EventSetup&) override;
    const std::vector<int>& operatesOnParticles() override { return fPDGs; }
    HepMC::GenEvent* decay(HepMC::GenEvent*) override;
    void statistics() override;
    void SetLHE(lhef::LHEEvent* l) override { lhe = l; }
    void setRandomEngine(CLHEP::HepRandomEngine* v) override { fRandomEngine = v; }
    static double flat();

  private:
    // member function(s)
    void decodeMDTAU(int);
    void selectDecayByMDTAU();
    int selectLeptonic();
    int selectHadronic();

    HepMC::GenEvent* make_simple_tau_event(const TLorentzVector& l, int pdgid, int status);
    void update_particles(HepMC::GenParticle* partHep,
                          HepMC::GenEvent* theEvent,
                          HepMC::GenParticle* p,
                          TVector3& boost);
    bool isLastTauInChain(const HepMC::GenParticle* tau);
    HepMC::GenParticle* GetMother(HepMC::GenParticle* tau);
    double MatchedLHESpinUp(HepMC::GenParticle* tau,
                            std::vector<HepMC::GenParticle>& p,
                            std::vector<double>& spinup,
                            std::vector<int>& m_idx);
    HepMC::GenParticle* FirstTauInChain(HepMC::GenParticle* tau);
    void BoostProdToLabLifeTimeInDecays(HepMC::GenParticle* p, TLorentzVector& lab, TLorentzVector& prod);

    //
    static CLHEP::HepRandomEngine* fRandomEngine;
    std::vector<int> fPDGs;
    bool fPolarization;
    edm::ESHandle<HepPDT::ParticleDataTable> fPDGTable;
    edm::ESGetToken<HepPDT::ParticleDataTable, edm::DefaultRecord> fPDGTableToken;
    edm::ParameterSet* fPSet;
    bool fIsInitialized;

    int fMDTAU;
    bool fSelectDecayByEvent;
    std::vector<int> fLeptonModes;
    std::vector<int> fHadronModes;
    std::vector<double> fScaledLeptonBrRatios;
    std::vector<double> fScaledHadronBrRatios;
    lhef::LHEEvent* lhe;

    double dmMatch;
    bool dolhe;
    bool dolheBosonCorr;
    int ntries;
    double lifetime;
  };

}  // namespace gen

#endif