CMS 3D CMS Logo

TauSpinnerTableProducer.cc
Go to the documentation of this file.
1 
22 
23 #include "Tauola/Tauola.h"
24 #include "TauSpinner/SimpleParticle.h"
25 #include "TauSpinner/tau_reweight_lib.h"
26 
27 #include <atomic>
28 
29 class TauSpinnerTableProducer : public edm::one::EDProducer<edm::one::SharedResources> {
30 public:
32 
33  void produce(edm::Event &, const edm::EventSetup &) final;
34  void beginJob() final;
35  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
37  desc.add<edm::InputTag>("src")->setComment("input genParticle collection");
38  desc.add<std::string>("name")->setComment("name of the TauSpinner weights table");
39  desc.add<std::vector<double>>("theta")->setComment("values of CP-even and CP-odd tau Yukawa mixing angle");
40  desc.ifValue(edm::ParameterDescription<int>("bosonPdgId", 25, true), edm::allowedValues<int>(25, 35, 36))
41  ->setComment("boson pdgId, default: 25"); // Allow only neutral Higgs bosons
42  desc.add<double>("defaultWeight", 1)
43  ->setComment("default weight stored in case of presence of a tau decay unsupported by TauSpinner");
44  descriptions.addWithDefaultLabel(desc);
45  }
46 
47 private:
50  static void getTaus(reco::GenParticleRefVector &taus, const reco::GenParticle &boson);
51  static bool getTauDaughters(reco::GenParticleRefVector &tau_daughters, const reco::GenParticle &tau);
52  TauSpinner::SimpleParticle convertToSimplePart(const reco::GenParticle &input_part) const {
53  return TauSpinner::SimpleParticle(
54  input_part.px(), input_part.py(), input_part.pz(), input_part.energy(), input_part.pdgId());
55  }
56 
57  static std::vector<std::pair<std::string, double>> nameAndValue(const std::vector<double> &val_vec) {
58  std::vector<std::pair<std::string, double>> out;
59  for (auto val : val_vec) {
61  name.erase(name.find_last_not_of('0') + 1, std::string::npos);
62  name.erase(name.find_last_not_of('.') + 1, std::string::npos);
63  size_t pos = name.find(".");
64  if (pos != std::string::npos)
65  name.replace(pos, 1, "p");
66  pos = name.find("-");
67  if (pos != std::string::npos)
68  name.replace(pos, 1, "minus");
69  out.push_back(std::make_pair(name, val));
70  }
71  return out;
72  }
73 
75  std::cout << std::string(78, '-') << "\n";
76  std::cout << config.getParameter<std::string>("@module_type") << '/'
77  << config.getParameter<std::string>("@module_label") << "\n";
78  std::cout << "Input collection: " << config.getParameter<edm::InputTag>("src").encode() << '\n';
79  std::cout << "Table name: " << config.getParameter<std::string>("name") << '\n';
80  std::string thetaStr;
81  for (const auto &theta : theta_vec_)
82  thetaStr += theta.first + ",";
83  std::cout << "Theta: " << thetaStr << std::endl;
84  }
85 
88  const std::vector<std::pair<std::string, double>> theta_vec_;
89  const int bosonPdgId_;
91  const bool ipp_;
92  const int ipol_;
93  const int nonSM2_;
94  const int nonSMN_;
95  const double cmsE_;
96  const double default_weight_;
97 
98  std::atomic<unsigned int> nWarnings{0};
99  static const unsigned int nMaxWarnings = 10;
100 };
101 
103  : genPartsToken_(consumes(config.getParameter<edm::InputTag>("src"))),
104  name_(config.getParameter<std::string>("name")),
105  theta_vec_(nameAndValue(config.getParameter<std::vector<double>>("theta"))),
106  bosonPdgId_(config.getParameter<int>("bosonPdgId")),
107  tauSpinnerPDF_(
108  "NNPDF31_nnlo_hessian_pdfas"), // PDF set for TauSpinner, relevant only in case of Z/gamma* polarization weights (set "sensible" default)
109  ipp_(true), // pp collisions
110  ipol_(0),
111  nonSM2_(0),
112  nonSMN_(0),
113  cmsE_(
114  13600), // collision energy in GeV, relevant only in case of Z/gamma* polarization weights (set "sensible" default)
115  default_weight_(config.getParameter<double>(
116  "defaultWeight")) // default weight stored in case of presence of a tau decay unsupported by TauSpinner
117 {
119 
120  // State that we use tauola/tauspinner resource
121  usesResource(edm::SharedResourceNames::kTauola);
122 
123  produces<nanoaod::FlatTable>();
124 }
125 
127  const edm::View<reco::GenParticle> &parts) const {
128  unsigned idx = 0;
129  for (const auto &part : parts) {
130  if (std::abs(part.pdgId()) == bosonPdgId_ && part.isLastCopy()) {
132  bosons.push_back(partRef);
133  }
134  ++idx;
135  }
136 }
137 
139  if (part->statusFlags().isLastCopy())
140  return part;
141  for (const auto &daughter : part->daughterRefVector()) {
142  if (daughter->pdgId() == part->pdgId() && daughter->statusFlags().fromHardProcess()) {
143  return getLastCopy(daughter);
144  }
145  }
146  throw std::runtime_error("getLastCopy: no last copy found");
147 }
148 
150  for (const auto &daughterRef : boson.daughterRefVector()) {
151  if (std::abs(daughterRef->pdgId()) == 15)
152  taus.push_back(getLastCopy(daughterRef));
153  }
154 }
155 
157  static const std::set<int> directTauProducts = {11, 12, 13, 14, 16, 22};
158  static const std::set<int> finalHadrons = {111, 130, 211, 310, 311, 321};
159  static const std::set<int> intermediateHadrons = {221, 223, 323};
160  for (auto daughterRef : tau.daughterRefVector()) {
161  const int daughter_pdgId = std::abs(daughterRef->pdgId());
162  if ((std::abs(tau.pdgId()) == 15 && directTauProducts.count(daughter_pdgId)) ||
163  finalHadrons.count(daughter_pdgId)) {
164  tau_daughters.push_back(daughterRef);
165  } else if (intermediateHadrons.count(daughter_pdgId)) {
166  if (!getTauDaughters(tau_daughters, *daughterRef))
167  return false;
168  } else {
169  edm::LogWarning("TauSpinnerTableProducer::getTauDaughters")
170  << "Unsupported decay with " << daughter_pdgId << " being daughter of " << std::abs(tau.pdgId()) << "\n";
171  return false;
172  }
173  }
174  return true;
175 }
176 
178  // Initialize TauSpinner
179  Tauolapp::Tauola::setNewCurrents(0);
181  LHAPDF::initPDFSetByName(tauSpinnerPDF_);
182  TauSpinner::initialize_spinner(ipp_, ipol_, nonSM2_, nonSMN_, cmsE_);
183 }
184 
186  // Input gen-particles collection
187  auto const &genParts = event.get(genPartsToken_);
188 
189  // Output table
190  auto wtTable = std::make_unique<nanoaod::FlatTable>(1, name_, true);
191  wtTable->setDoc("TauSpinner weights");
192 
193  // Search for boson
195  getBosons(bosons, genParts);
196  if (bosons.size() !=
197  1) { // no boson found or more than one found, produce table with default weights (expected for non HTT samples)
198  if (++nWarnings < nMaxWarnings)
199  edm::LogWarning("TauSpinnerTableProducer::produce")
200  << "Current event has " << bosons.size()
201  << " Higgs bosons while there must be exactly one; table with default weights is produced.\n";
202  // Fill table with default values
203  for (const auto &theta : theta_vec_) {
204  wtTable->addColumnValue<double>(
205  "weight_cp_" + theta.first, default_weight_, "TauSpinner weight for theta_CP = " + theta.first);
206  wtTable->addColumnValue<double>(
207  "weight_cp_" + theta.first + "_alt",
209  "TauSpinner weight for theta_CP = " + theta.first + " (alternative hadronic currents)");
210  }
211  event.put(std::move(wtTable));
212  return;
213  }
214 
215  // Search for taus from boson decay
217  getTaus(taus, *bosons[0]);
218  if (taus.size() !=
219  2) { // boson does not decay to tau pair, produce table with default weights (expected for non HTT samples)
220  if (++nWarnings < nMaxWarnings)
221  edm::LogWarning("TauSpinnerTableProducer::produce")
222  << "Current event has " << taus.size()
223  << " taus from boson decay while there must be exactly one pair; table with default weights is produced.\n";
224  // Fill table with default values
225  for (const auto &theta : theta_vec_) {
226  wtTable->addColumnValue<double>(
227  "weight_cp_" + theta.first, default_weight_, "TauSpinner weight for theta_CP = " + theta.first);
228  wtTable->addColumnValue<double>(
229  "weight_cp_" + theta.first + "_alt",
231  "TauSpinner weight for theta_CP = " + theta.first + " (alternative hadronic currents)");
232  }
233  event.put(std::move(wtTable));
234  return;
235  }
236 
237  // Get tau daughters and convert all particles to TauSpinner format
238  TauSpinner::SimpleParticle simple_boson = convertToSimplePart(*bosons[0]);
239  std::array<TauSpinner::SimpleParticle, 2> simple_taus;
240  std::array<std::vector<TauSpinner::SimpleParticle>, 2> simple_tau_daughters;
241  bool supportedDecays = true;
242  for (size_t tau_idx = 0; tau_idx < 2; ++tau_idx) {
243  simple_taus[tau_idx] = convertToSimplePart(*taus[tau_idx]);
244  reco::GenParticleRefVector tau_daughters;
245  supportedDecays &= getTauDaughters(tau_daughters, *taus[tau_idx]);
246  for (const auto &daughterRef : tau_daughters)
247  simple_tau_daughters[tau_idx].push_back(convertToSimplePart(*daughterRef));
248  }
249 
250  // Compute TauSpinner weights and fill table
251  std::array<double, 2> weights;
252  for (const auto &theta : theta_vec_) {
253  // Can make this more general by having boson pdgid as input or have option for set boson type
254  TauSpinner::setHiggsParametersTR(-cos(2 * M_PI * theta.second),
255  cos(2 * M_PI * theta.second),
256  -sin(2 * M_PI * theta.second),
257  -sin(2 * M_PI * theta.second));
258  for (size_t i = 0; i < weights.size(); ++i) {
259  Tauolapp::Tauola::setNewCurrents(i);
260  weights[i] =
261  supportedDecays
262  ? TauSpinner::calculateWeightFromParticlesH(
263  simple_boson, simple_taus[0], simple_taus[1], simple_tau_daughters[0], simple_tau_daughters[1])
264  : default_weight_;
265  }
266  // Nominal weights for setNewCurrents(0)
267  wtTable->addColumnValue<double>(
268  "weight_cp_" + theta.first, weights[0], "TauSpinner weight for theta_CP=" + theta.first);
269  // Weights for alternative hadronic currents (can be used for uncertainty estimates)
270  wtTable->addColumnValue<double>(
271  "weight_cp_" + theta.first + "_alt",
272  weights[1],
273  "TauSpinner weight for theta_CP=" + theta.first + " (alternative hadronic currents)");
274  }
275 
276  event.put(std::move(wtTable));
277 }
278 
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)
Definition: Sin.h:22
static const std::string kTauola
Definition: config.py:1
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)
Definition: Cos.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static bool getTauDaughters(reco::GenParticleRefVector &tau_daughters, const reco::GenParticle &tau)
TauSpinnerTableProducer(const edm::ParameterSet &)
double py() const final
y coordinate of momentum vector
#define M_PI
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.
Definition: RefVector.h:102
const std::vector< std::pair< std::string, double > > theta_vec_
std::atomic< unsigned int > nWarnings
static reco::GenParticleRef getLastCopy(const reco::GenParticleRef &part)
part
Definition: HCALResponse.h:20
def encode(args, files)
HLT enums.
static const unsigned int nMaxWarnings
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
Log< level::Warning, false > LogWarning
void produce(edm::Event &, const edm::EventSetup &) final
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1
double energy() const final
energy