23 #include "Tauola/Tauola.h" 24 #include "TauSpinner/SimpleParticle.h" 25 #include "TauSpinner/tau_reweight_lib.h" 36 desc.add<
std::string>(
"name")->setComment(
"name of the TauSpinner weights table");
37 desc.add<std::vector<double>>(
"theta")->setComment(
"values of CP-even and CP-odd tau Yukawa mixing angle");
39 ->setComment(
"boson pdgId, default: 25");
40 desc.add<
double>(
"defaultWeight", 1)
41 ->setComment(
"default weight stored in case of presence of a tau decay unsupported by TauSpinner");
51 return TauSpinner::SimpleParticle(
52 input_part.
px(), input_part.
py(), input_part.
pz(), input_part.
energy(), input_part.
pdgId());
55 static std::vector<std::pair<std::string, double>>
nameAndValue(
const std::vector<double> &val_vec) {
56 std::vector<std::pair<std::string, double>>
out;
57 for (
auto val : val_vec) {
59 name.erase(
name.find_last_not_of(
'0') + 1, std::string::npos);
60 name.erase(
name.find_last_not_of(
'.') + 1, std::string::npos);
62 if (
pos != std::string::npos)
65 if (
pos != std::string::npos)
80 thetaStr +=
theta.first +
",";
81 std::cout <<
"Theta: " << thetaStr << std::endl;
86 const std::vector<std::pair<std::string, double>>
theta_vec_;
100 theta_vec_(nameAndValue(
config.getParameter<
std::
vector<double>>(
"theta"))),
101 bosonPdgId_(
config.getParameter<
int>(
"bosonPdgId")),
103 "NNPDF31_nnlo_hessian_pdfas"),
110 default_weight_(
config.getParameter<double>(
118 produces<nanoaod::FlatTable>();
127 bosons.push_back(partRef);
134 if (
part->statusFlags().isLastCopy())
136 for (
const auto &daughter :
part->daughterRefVector()) {
137 if (daughter->pdgId() ==
part->pdgId() && daughter->statusFlags().fromHardProcess()) {
141 throw std::runtime_error(
"getLastCopy: no last copy found");
146 if (
std::abs(daughterRef->pdgId()) == 15)
152 static const std::set<int> directTauProducts = {11, 12, 13, 14, 16, 22};
153 static const std::set<int> finalHadrons = {111, 130, 211, 310, 311, 321};
154 static const std::set<int> intermediateHadrons = {221, 223, 323};
155 for (
auto daughterRef :
tau.daughterRefVector()) {
156 const int daughter_pdgId =
std::abs(daughterRef->pdgId());
157 if ((
std::abs(
tau.pdgId()) == 15 && directTauProducts.count(daughter_pdgId)) ||
158 finalHadrons.count(daughter_pdgId)) {
160 }
else if (intermediateHadrons.count(daughter_pdgId)) {
165 <<
"Unsupported decay with " << daughter_pdgId <<
" being daughter of " <<
std::abs(
tau.pdgId()) <<
"\n";
174 Tauolapp::Tauola::setNewCurrents(0);
185 auto wtTable = std::make_unique<nanoaod::FlatTable>(1,
name_,
true);
186 wtTable->setDoc(
"TauSpinner weights");
200 if (
taus.size() != 2) {
207 std::array<TauSpinner::SimpleParticle, 2> simple_taus;
208 std::array<std::vector<TauSpinner::SimpleParticle>, 2> simple_tau_daughters;
209 bool supportedDecays =
true;
210 for (
size_t tau_idx = 0; tau_idx < 2; ++tau_idx) {
214 for (
const auto &daughterRef : tau_daughters)
222 TauSpinner::setHiggsParametersTR(-
cos(2 *
M_PI *
theta.second),
227 Tauolapp::Tauola::setNewCurrents(
i);
230 ? TauSpinner::calculateWeightFromParticlesH(
231 simple_boson, simple_taus[0], simple_taus[1], simple_tau_daughters[0], simple_tau_daughters[1])
235 wtTable->addColumnValue<
double>(
236 "weight_cp_" +
theta.first,
weights[0],
"TauSpinner weight for theta_CP=" +
theta.first);
238 wtTable->addColumnValue<
double>(
239 "weight_cp_" +
theta.first +
"_alt",
241 "TauSpinner weight for theta_CP=" +
theta.first +
" (alternative hadronic currents)");
TauSpinner::SimpleParticle convertToSimplePart(const reco::GenParticle &input_part) const
static AlgebraicMatrix initialize()
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
void printModuleInfo(edm::ParameterSet const &config) const
double pz() const final
z coordinate of momentum vector
static void getTaus(reco::GenParticleRefVector &taus, const reco::GenParticle &boson)
Sin< T >::type sin(const T &t)
static const std::string kTauola
const daughters & daughterRefVector() const
references to daughtes
static std::vector< std::pair< std::string, double > > nameAndValue(const std::vector< double > &val_vec)
static std::string to_string(const XMLCh *ch)
int pdgId() const final
PDG identifier.
double px() const final
x coordinate of momentum vector
Cos< T >::type cos(const T &t)
Abs< T >::type abs(const T &t)
#define DEFINE_FWK_MODULE(type)
static bool getTauDaughters(reco::GenParticleRefVector &tau_daughters, const reco::GenParticle &tau)
const double default_weight_
TauSpinnerTableProducer(const edm::ParameterSet &)
double py() const final
y coordinate of momentum vector
const std::string tauSpinnerPDF_
const edm::EDGetTokenT< edm::View< reco::GenParticle > > genPartsToken_
void getBosons(edm::RefVector< edm::View< reco::GenParticle >> &bosons, const edm::View< reco::GenParticle > &parts) const
size_type size() const
Size of the RefVector.
const std::vector< std::pair< std::string, double > > theta_vec_
static reco::GenParticleRef getLastCopy(const reco::GenParticleRef &part)
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Log< level::Warning, false > LogWarning
void produce(edm::Event &, const edm::EventSetup &) final
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double energy() const final
energy