|
| 1 | +/**\class TauSpinnerTableProducer |
| 2 | +
|
| 3 | + Description: Produces FlatTable with TauSpinner weights for H->tau,tau events |
| 4 | +
|
| 5 | + Original Author: D. Winterbottom (IC) |
| 6 | + Update and adaptation to NanoAOD: M. Bluj (NCBJ) |
| 7 | +
|
| 8 | +*/ |
| 9 | + |
| 10 | +#include "FWCore/Framework/interface/one/EDProducer.h" |
| 11 | +#include "FWCore/Framework/interface/Event.h" |
| 12 | +#include "FWCore/Framework/interface/EventSetup.h" |
| 13 | +#include "FWCore/Concurrency/interface/SharedResourceNames.h" |
| 14 | +#include "FWCore/MessageLogger/interface/MessageLogger.h" |
| 15 | +#include "FWCore/ParameterSet/interface/ParameterSet.h" |
| 16 | +#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h" |
| 17 | +#include "FWCore/ParameterSet/interface/ParameterSetDescription.h" |
| 18 | +#include "FWCore/ParameterSet/interface/allowedValues.h" |
| 19 | +#include "DataFormats/Common/interface/View.h" |
| 20 | +#include "DataFormats/HepMCCandidate/interface/GenParticle.h" |
| 21 | +#include "DataFormats/NanoAOD/interface/FlatTable.h" |
| 22 | + |
| 23 | +#include "Tauola/Tauola.h" |
| 24 | +#include "TauSpinner/SimpleParticle.h" |
| 25 | +#include "TauSpinner/tau_reweight_lib.h" |
| 26 | + |
| 27 | +class TauSpinnerTableProducer : public edm::one::EDProducer<edm::one::SharedResources> { |
| 28 | +public: |
| 29 | + explicit TauSpinnerTableProducer(const edm::ParameterSet &); |
| 30 | + |
| 31 | + void produce(edm::Event &, const edm::EventSetup &) final; |
| 32 | + void beginJob() final; |
| 33 | + static void fillDescriptions(edm::ConfigurationDescriptions &descriptions) { |
| 34 | + edm::ParameterSetDescription desc; |
| 35 | + desc.add<edm::InputTag>("src")->setComment("input genParticle collection"); |
| 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"); |
| 38 | + desc.ifValue(edm::ParameterDescription<int>("bosonPdgId", 25, true), edm::allowedValues<int>(25, 35, 36)) |
| 39 | + ->setComment("boson pdgId, default: 25"); // Allow only neutral Higgs bosons |
| 40 | + desc.add<double>("defaultWeight", 1) |
| 41 | + ->setComment("default weight stored in case of presence of a tau decay unsupported by TauSpinner"); |
| 42 | + descriptions.addWithDefaultLabel(desc); |
| 43 | + } |
| 44 | + |
| 45 | +private: |
| 46 | + void getBosons(edm::RefVector<edm::View<reco::GenParticle>> &bosons, const edm::View<reco::GenParticle> &parts) const; |
| 47 | + static reco::GenParticleRef getLastCopy(const reco::GenParticleRef &part); |
| 48 | + static void getTaus(reco::GenParticleRefVector &taus, const reco::GenParticle &boson); |
| 49 | + static bool getTauDaughters(reco::GenParticleRefVector &tau_daughters, const reco::GenParticle &tau); |
| 50 | + TauSpinner::SimpleParticle convertToSimplePart(const reco::GenParticle &input_part) const { |
| 51 | + return TauSpinner::SimpleParticle( |
| 52 | + input_part.px(), input_part.py(), input_part.pz(), input_part.energy(), input_part.pdgId()); |
| 53 | + } |
| 54 | + |
| 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) { |
| 58 | + std::string name = std::to_string(val); |
| 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); |
| 61 | + size_t pos = name.find("."); |
| 62 | + if (pos != std::string::npos) |
| 63 | + name.replace(pos, 1, "p"); |
| 64 | + pos = name.find("-"); |
| 65 | + if (pos != std::string::npos) |
| 66 | + name.replace(pos, 1, "minus"); |
| 67 | + out.push_back(std::make_pair(name, val)); |
| 68 | + } |
| 69 | + return out; |
| 70 | + } |
| 71 | + |
| 72 | + void printModuleInfo(edm::ParameterSet const &config) const { |
| 73 | + std::cout << std::string(78, '-') << "\n"; |
| 74 | + std::cout << config.getParameter<std::string>("@module_type") << '/' |
| 75 | + << config.getParameter<std::string>("@module_label") << "\n"; |
| 76 | + std::cout << "Input collection: " << config.getParameter<edm::InputTag>("src").encode() << '\n'; |
| 77 | + std::cout << "Table name: " << config.getParameter<std::string>("name") << '\n'; |
| 78 | + std::string thetaStr; |
| 79 | + for (const auto &theta : theta_vec_) |
| 80 | + thetaStr += theta.first + ","; |
| 81 | + std::cout << "Theta: " << thetaStr << std::endl; |
| 82 | + } |
| 83 | + |
| 84 | + const edm::EDGetTokenT<edm::View<reco::GenParticle>> genPartsToken_; |
| 85 | + const std::string name_; |
| 86 | + const std::vector<std::pair<std::string, double>> theta_vec_; |
| 87 | + const int bosonPdgId_; |
| 88 | + const std::string tauSpinnerPDF_; |
| 89 | + const bool ipp_; |
| 90 | + const int ipol_; |
| 91 | + const int nonSM2_; |
| 92 | + const int nonSMN_; |
| 93 | + const double cmsE_; |
| 94 | + const double default_weight_; |
| 95 | +}; |
| 96 | + |
| 97 | +TauSpinnerTableProducer::TauSpinnerTableProducer(const edm::ParameterSet &config) |
| 98 | + : genPartsToken_(consumes(config.getParameter<edm::InputTag>("src"))), |
| 99 | + name_(config.getParameter<std::string>("name")), |
| 100 | + theta_vec_(nameAndValue(config.getParameter<std::vector<double>>("theta"))), |
| 101 | + bosonPdgId_(config.getParameter<int>("bosonPdgId")), |
| 102 | + tauSpinnerPDF_( |
| 103 | + "NNPDF31_nnlo_hessian_pdfas"), // PDF set for TauSpinner, relevant only in case of Z/gamma* polarization weights (set "sensible" default) |
| 104 | + ipp_(true), // pp collisions |
| 105 | + ipol_(0), |
| 106 | + nonSM2_(0), |
| 107 | + nonSMN_(0), |
| 108 | + cmsE_( |
| 109 | + 13600), // collision energy in GeV, relevant only in case of Z/gamma* polarization weights (set "sensible" default) |
| 110 | + default_weight_(config.getParameter<double>( |
| 111 | + "defaultWeight")) // default weight stored in case of presence of a tau decay unsupported by TauSpinner |
| 112 | +{ |
| 113 | + printModuleInfo(config); |
| 114 | + |
| 115 | + // State that we use tauola/tauspinner resource |
| 116 | + usesResource(edm::SharedResourceNames::kTauola); |
| 117 | + |
| 118 | + produces<nanoaod::FlatTable>(); |
| 119 | +} |
| 120 | + |
| 121 | +void TauSpinnerTableProducer::getBosons(edm::RefVector<edm::View<reco::GenParticle>> &bosons, |
| 122 | + const edm::View<reco::GenParticle> &parts) const { |
| 123 | + unsigned idx = 0; |
| 124 | + for (const auto &part : parts) { |
| 125 | + if (std::abs(part.pdgId()) == bosonPdgId_ && part.isLastCopy()) { |
| 126 | + edm::Ref<edm::View<reco::GenParticle>> partRef(&parts, idx); |
| 127 | + bosons.push_back(partRef); |
| 128 | + } |
| 129 | + ++idx; |
| 130 | + } |
| 131 | +} |
| 132 | + |
| 133 | +reco::GenParticleRef TauSpinnerTableProducer::getLastCopy(const reco::GenParticleRef &part) { |
| 134 | + if (part->statusFlags().isLastCopy()) |
| 135 | + return part; |
| 136 | + for (const auto &daughter : part->daughterRefVector()) { |
| 137 | + if (daughter->pdgId() == part->pdgId() && daughter->statusFlags().fromHardProcess()) { |
| 138 | + return getLastCopy(daughter); |
| 139 | + } |
| 140 | + } |
| 141 | + throw std::runtime_error("getLastCopy: no last copy found"); |
| 142 | +} |
| 143 | + |
| 144 | +void TauSpinnerTableProducer::getTaus(reco::GenParticleRefVector &taus, const reco::GenParticle &boson) { |
| 145 | + for (const auto &daughterRef : boson.daughterRefVector()) { |
| 146 | + if (std::abs(daughterRef->pdgId()) == 15) |
| 147 | + taus.push_back(getLastCopy(daughterRef)); |
| 148 | + } |
| 149 | +} |
| 150 | + |
| 151 | +bool TauSpinnerTableProducer::getTauDaughters(reco::GenParticleRefVector &tau_daughters, const reco::GenParticle &tau) { |
| 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)) { |
| 159 | + tau_daughters.push_back(daughterRef); |
| 160 | + } else if (intermediateHadrons.count(daughter_pdgId)) { |
| 161 | + if (!getTauDaughters(tau_daughters, *daughterRef)) |
| 162 | + return false; |
| 163 | + } else { |
| 164 | + edm::LogWarning("TauSpinnerTableProducer::getTauDaughters") |
| 165 | + << "Unsupported decay with " << daughter_pdgId << " being daughter of " << std::abs(tau.pdgId()) << "\n"; |
| 166 | + return false; |
| 167 | + } |
| 168 | + } |
| 169 | + return true; |
| 170 | +} |
| 171 | + |
| 172 | +void TauSpinnerTableProducer::beginJob() { |
| 173 | + // Initialize TauSpinner |
| 174 | + Tauolapp::Tauola::setNewCurrents(0); |
| 175 | + Tauolapp::Tauola::initialize(); |
| 176 | + LHAPDF::initPDFSetByName(tauSpinnerPDF_); |
| 177 | + TauSpinner::initialize_spinner(ipp_, ipol_, nonSM2_, nonSMN_, cmsE_); |
| 178 | +} |
| 179 | + |
| 180 | +void TauSpinnerTableProducer::produce(edm::Event &event, const edm::EventSetup &setup) { |
| 181 | + // Input gen-particles collection |
| 182 | + auto const &genParts = event.get(genPartsToken_); |
| 183 | + |
| 184 | + // Output table |
| 185 | + auto wtTable = std::make_unique<nanoaod::FlatTable>(1, name_, true); |
| 186 | + wtTable->setDoc("TauSpinner weights"); |
| 187 | + |
| 188 | + // Search for boson |
| 189 | + edm::RefVector<edm::View<reco::GenParticle>> bosons; |
| 190 | + getBosons(bosons, genParts); |
| 191 | + if (bosons.size() != |
| 192 | + 1) { // no boson found or more than one found, produce empty table (expected for non HTT samples) |
| 193 | + event.put(std::move(wtTable)); |
| 194 | + return; |
| 195 | + } |
| 196 | + |
| 197 | + // Search for taus from boson decay |
| 198 | + reco::GenParticleRefVector taus; |
| 199 | + getTaus(taus, *bosons[0]); |
| 200 | + if (taus.size() != 2) { // boson does not decay to tau pair, produce empty table (expected for non HTT samples) |
| 201 | + event.put(std::move(wtTable)); |
| 202 | + return; |
| 203 | + } |
| 204 | + |
| 205 | + // Get tau daughters and convert all particles to TauSpinner format |
| 206 | + TauSpinner::SimpleParticle simple_boson = convertToSimplePart(*bosons[0]); |
| 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) { |
| 211 | + simple_taus[tau_idx] = convertToSimplePart(*taus[tau_idx]); |
| 212 | + reco::GenParticleRefVector tau_daughters; |
| 213 | + supportedDecays &= getTauDaughters(tau_daughters, *taus[tau_idx]); |
| 214 | + for (const auto &daughterRef : tau_daughters) |
| 215 | + simple_tau_daughters[tau_idx].push_back(convertToSimplePart(*daughterRef)); |
| 216 | + } |
| 217 | + |
| 218 | + // Compute TauSpinner weights and fill table |
| 219 | + std::array<double, 2> weights; |
| 220 | + for (const auto &theta : theta_vec_) { |
| 221 | + // Can make this more general by having boson pdgid as input or have option for set boson type |
| 222 | + TauSpinner::setHiggsParametersTR(-cos(2 * M_PI * theta.second), |
| 223 | + cos(2 * M_PI * theta.second), |
| 224 | + -sin(2 * M_PI * theta.second), |
| 225 | + -sin(2 * M_PI * theta.second)); |
| 226 | + for (size_t i = 0; i < weights.size(); ++i) { |
| 227 | + Tauolapp::Tauola::setNewCurrents(i); |
| 228 | + weights[i] = |
| 229 | + supportedDecays |
| 230 | + ? TauSpinner::calculateWeightFromParticlesH( |
| 231 | + simple_boson, simple_taus[0], simple_taus[1], simple_tau_daughters[0], simple_tau_daughters[1]) |
| 232 | + : default_weight_; |
| 233 | + } |
| 234 | + // Nominal weights for setNewCurrents(0) |
| 235 | + wtTable->addColumnValue<double>( |
| 236 | + "weight_cp_" + theta.first, weights[0], "TauSpinner weight for theta_CP=" + theta.first); |
| 237 | + // Weights for alternative hadronic currents (can be used for uncertainty estimates) |
| 238 | + wtTable->addColumnValue<double>( |
| 239 | + "weight_cp_" + theta.first + "_alt", |
| 240 | + weights[1], |
| 241 | + "TauSpinner weight for theta_CP=" + theta.first + " (alternative hadronic currents)"); |
| 242 | + } |
| 243 | + |
| 244 | + event.put(std::move(wtTable)); |
| 245 | +} |
| 246 | + |
| 247 | +#include "FWCore/Framework/interface/MakerMacros.h" |
| 248 | +DEFINE_FWK_MODULE(TauSpinnerTableProducer); |
0 commit comments