4 #include "Tauola/Tauola.h" 5 #include "TauSpinner/tau_reweight_lib.h" 6 #include "TauSpinner/Tauola_wrapper.h" 8 #include "TLorentzVector.h" 10 #include "CLHEP/Random/RandomEngine.h" 25 ,isReco_(pset.getParameter<
bool>(
"isReco"))
26 ,isTauolaConfigured_(pset.getParameter<
bool>(
"isTauolaConfigured" ))
27 ,isLHPDFConfigured_(pset.getParameter<
bool>(
"isLHPDFConfigured" ))
28 ,LHAPDFname_(pset.getUntrackedParameter(
"LHAPDFname",(
string)(
"MSTW2008nnlo90cl.LHgrid")))
29 ,CMSEnergy_(pset.getParameter<double>(
"CMSEnergy"))
30 ,gensrc_(pset.getParameter<
edm::
InputTag>(
"gensrc"))
31 ,MotherPDGID_(pset.getUntrackedParameter(
"MotherPDGID",(
int)(-1)))
32 ,Ipol_(pset.getUntrackedParameter(
"Ipol",(
int)(0)))
33 ,nonSM2_(pset.getUntrackedParameter(
"nonSM2",(
int)(0)))
34 ,nonSMN_(pset.getUntrackedParameter(
"nonSMN",(
int)(0)))
35 ,roundOff_(pset.getUntrackedParameter(
"roundOff",(double)(0.01)))
40 produces<bool>(
"TauSpinnerWTisValid").setBranchAlias(
"TauSpinnerWTisValid");
41 produces<double>(
"TauSpinnerWT").setBranchAlias(
"TauSpinnerWT");
42 produces<double>(
"TauSpinnerWTFlip").setBranchAlias(
"TauSpinnerWTFlip");
43 produces<double>(
"TauSpinnerWThplus").setBranchAlias(
"TauSpinnerWThplus");
44 produces<double>(
"TauSpinnerWThminus").setBranchAlias(
"TauSpinnerWThminus");
84 Tauolapp::jaki_.ktom = 1;
89 std::vector<SimpleParticle> tau_daughters, tau_daughters2;
98 HepMC::GenEvent *Evt =
new HepMC::GenEvent(*(evt->
GetEvent()));
104 if(
abs(X.pdgid())==24 ||
abs(X.pdgid())==37 ){
105 TLorentzVector tau_1r(0,0,0,0);
106 TLorentzVector tau_1(tau.px(),tau.py(),tau.pz(),tau.e());
107 for(
unsigned int i=0;
i<tau_daughters.size();
i++){
108 tau_1r+=TLorentzVector(tau_daughters.at(
i).px(),tau_daughters.at(
i).py(),tau_daughters.at(
i).pz(),tau_daughters.at(
i).e());
110 if(fabs(tau_1r.M()-tau_1.M())<
roundOff_){
111 WT = TauSpinner::calculateWeightFromParticlesWorHpn(X, tau, tau2, tau_daughters);
116 else if( X.pdgid()==25 || X.pdgid()==36 || X.pdgid()==22 || X.pdgid()==23 ){
117 TLorentzVector tau_1r(0,0,0,0), tau_2r(0,0,0,0);
118 TLorentzVector tau_1(tau.px(),tau.py(),tau.pz(),tau.e()), tau_2(tau2.px(),tau2.py(),tau2.pz(),tau2.e());
119 for(
unsigned int i=0;
i<tau_daughters.size();
i++){
120 tau_1r+=TLorentzVector(tau_daughters.at(
i).px(),tau_daughters.at(
i).py(),tau_daughters.at(
i).pz(),tau_daughters.at(
i).e());
122 for(
unsigned int i=0;
i<tau_daughters2.size();
i++){
123 tau_2r+=TLorentzVector(tau_daughters2.at(
i).px(),tau_daughters2.at(
i).py(),tau_daughters2.at(
i).pz(),tau_daughters2.at(
i).e());
127 WT = TauSpinner::calculateWeightFromParticlesH(X, tau, tau2, tau_daughters,tau_daughters2);
130 if(X.pdgid()==25 || X.pdgid()==22 || X.pdgid()==23 ){
131 if(X.pdgid()==25) X.setPdgid(23);
132 if( X.pdgid()==22 || X.pdgid()==23 ) X.setPdgid(25);
134 double WTother=TauSpinner::calculateWeightFromParticlesH(X, tau, tau2, tau_daughters,tau_daughters2);
140 cout<<
"TauSpinner: WARNING: Unexpected PDG for tau mother: "<<X.pdgid()<<endl;
145 if(!(0<=WT && WT<10)){isValid=
false; WT=1.0; WTFlip=1.0;}
146 std::unique_ptr<bool> TauSpinnerWeightisValid(
new bool);
147 *TauSpinnerWeightisValid =isValid;
148 e.
put(
std::move(TauSpinnerWeightisValid),
"TauSpinnerWTisValid");
151 std::unique_ptr<double> TauSpinnerWeight(
new double);
152 *TauSpinnerWeight =WT;
156 std::unique_ptr<double> TauSpinnerWeightFlip(
new double);
157 *TauSpinnerWeightFlip =WTFlip;
158 e.
put(
std::move(TauSpinnerWeightFlip),
"TauSpinnerWTFlip");
162 if(polSM<0.0 && polSM!=-999 && isValid) WThplus=0;
163 std::unique_ptr<double> TauSpinnerWeighthplus(
new double);
164 *TauSpinnerWeighthplus = WThplus;
165 e.
put(
std::move(TauSpinnerWeighthplus),
"TauSpinnerWThplus");
169 if(polSM>0.0&& polSM!=-999 && isValid) WThminus=0;
170 std::unique_ptr<double> TauSpinnerWeighthminus(
new double);
171 *TauSpinnerWeighthminus = WThminus;
172 e.
put(
std::move(TauSpinnerWeighthminus),
"TauSpinnerWThminus");
177 std::vector<SimpleParticle> &tau_daughters,std::vector<SimpleParticle> &tau2_daughters){
180 for(reco::GenParticleCollection::const_iterator itr = genParticles->begin(); itr!= genParticles->end(); ++itr){
182 if(pdgid==24 || pdgid==37 || pdgid ==25 || pdgid==36 || pdgid==22 || pdgid==23 ){
188 unsigned int ntau(0),ntauornu(0);
189 for(
unsigned int i=0;
i<itr->numberOfDaughters();
i++){
193 if(ntau==0 &&
abs(dau->
pdgId())==15){
206 if((ntau==2 && ntauornu==2) || (ntau==1 && ntauornu==2)){
207 X.setPx(itr->p4().Px());
208 X.setPy(itr->p4().Py());
209 X.setPz(itr->p4().Pz());
210 X.setE (itr->p4().E());
211 X.setPdgid(itr->pdgId());
212 tau.setPx(recotau1->
p4().Px());
213 tau.setPy(recotau1->
p4().Py());
214 tau.setPz(recotau1->
p4().Pz());
215 tau.setE (recotau1->
p4().E());
216 tau.setPdgid(recotau1->
pdgId());
218 tau2.setPx(recotau2->
p4().Px());
219 tau2.setPy(recotau2->
p4().Py());
220 tau2.setPz(recotau2->
p4().Pz());
221 tau2.setE (recotau2->
p4().E());
222 tau2.setPdgid(recotau2->
pdgId());
253 SimpleParticle tp(Particle->
p4().Px(), Particle->
p4().Py(), Particle->
p4().Pz(), Particle->
p4().E(), Particle->
pdgId());
254 daughters.push_back(tp);
267 <<
"TauSpinnerCMS::flat: Attempt to generate random number when engine pointer is null\n" 268 <<
"This might mean that the code was modified to generate a random number outside the\n" 269 <<
"event and beginLuminosityBlock methods, which is not allowed.\n";
static AlgebraicMatrix initialize()
int pdgId() const final
PDG identifier.
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
void GetLastSelf(const reco::GenParticle *Particle)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
void produce(edm::Event &, const edm::EventSetup &) final
bool getByToken(EDGetToken token, Handle< PROD > &result) const
#define DEFINE_FWK_MODULE(type)
virtual void initialize()
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
size_t numberOfMothers() const override
number of mothers
bool isFirst(const reco::GenParticle *Particle)
static const std::string kTauola
static CLHEP::HepRandomEngine * fRandomEngine
edm::EDGetTokenT< reco::GenParticleCollection > GenParticleCollectionToken_
void GetRecoDaughters(const reco::GenParticle *Particle, std::vector< TauSpinner::SimpleParticle > &daughters, int parentpdgid)
static bool isTauSpinnerConfigure
virtual int pdgId() const =0
PDG identifier.
Abs< T >::type abs(const T &t)
const LorentzVector & p4() const final
four-momentum Lorentz vector
int readParticlesfromReco(edm::Event &e, TauSpinner::SimpleParticle &X, TauSpinner::SimpleParticle &tau, TauSpinner::SimpleParticle &tau2, std::vector< TauSpinner::SimpleParticle > &tau_daughters, std::vector< TauSpinner::SimpleParticle > &tau2_daughters)
int readParticlesFromHepMC(const HepMC::GenEvent *event, TauSpinner::SimpleParticle &X, TauSpinner::SimpleParticle &tau, TauSpinner::SimpleParticle &tau2, std::vector< TauSpinner::SimpleParticle > &tau_daughters, std::vector< TauSpinner::SimpleParticle > &tau2_daughters)
const HepMC::GenEvent * GetEvent() const
size_t numberOfDaughters() const override
number of daughters
TauSpinnerCMS(const edm::ParameterSet &)
StreamID streamID() const
const Candidate * mother(size_type=0) const override
return mother at a given position, i = 0, ... numberOfMothers() - 1 (read only mode) ...