37 #include <Pythia8/Pythia.h>
42 #include "TLorentzVector.h"
59 virtual void endJob()
override;
69 std::map<int,TH1D*>
h_p;
70 std::map<int,TH1D*>
h_v;
105 pythia =
new Pythia8::Pythia();
120 for(
size_t i = 0;
i<
pids.size();++
i){
125 if(!pdt.isParticle(pid)){
126 std::cout <<
"ERROR: BAD PARTICLE, pythia is not aware of pid " << pid << std::endl;
129 Pythia8::ParticleDataEntry * pd = pdt.particleDataEntryPtr(pid);
132 double m0 = pd->m0();
133 double w = pd->mWidth();
136 mmin = m0 - m0/1000.;
137 mmax = m0 + m0/1000.;
143 std::stringstream strstr;
144 strstr <<
"mass_" <<
pid;
145 h_mass[
pid] =
new TH1D(strstr.str().c_str(),strstr.str().c_str(),100,mmin,mmax);
147 h_mass_ref[
pid]->SetTitle(h_mass_ref[pid]->GetName());
149 h_mass_ref[
pid]->Fill(m0);
151 for(
int b =1;
b<=h_mass_ref[
pid]->GetNbinsX();++
b){
152 double _val = h_mass_ref[
pid]->GetBinCenter(
b);
153 h_mass_ref[
pid]->SetBinContent(
b,TMath::BreitWigner(_val,m0,w));
159 strstr <<
"p_" <<
pid;
160 h_p[
pid] =
new TH1D(strstr.str().c_str(),strstr.str().c_str(),100,0,20);
164 strstr <<
"v_" <<
pid;
165 h_v[
pid] =
new TH1D(strstr.str().c_str(),strstr.str().c_str(),100,0,1.);
168 double ctau0 = pd->tau0()/10.;
170 strstr <<
"plt_" <<
pid;
171 h_plt[
pid] =
new TH1D(strstr.str().c_str(),strstr.str().c_str(),100,0,
std::min(5.*ctau0,500.));
173 h_plt_ref[
pid]->SetTitle(h_plt_ref[pid]->GetName());
174 for(
int b =1;
b<=h_plt_ref[
pid]->GetNbinsX();++
b){
175 double _val = h_plt_ref[
pid]->GetBinCenter(
b);
176 h_plt_ref[
pid]->SetBinContent(
b,TMath::Exp(-_val/ctau0));
182 strstr <<
"br_" <<
pid;
183 h_br[
pid] =
new TH1D(strstr.str().c_str(),strstr.str().c_str(),0,0,0);
184 h_br[
pid]->SetBit(TH1::kCanRebin);
188 for(
int d = 0;d<pd->sizeChannels();++d){
189 Pythia8::DecayChannel & channel = pd->channel(d);
190 std::vector<int>
prod;
191 for(
int p = 0;
p<channel.multiplicity();++
p){
192 int pId =
abs(channel.product(
p));
194 bool particleCut = ( pId > 10 && pId != 12 && pId != 14 &&
195 pId != 16 && pId != 18 && pId != 21 &&
196 (pId < 23 || pId > 40 ) &&
197 (pId < 81 || pId > 100 ) && pId != 2101 &&
198 pId != 3101 && pId != 3201 && pId != 1103 &&
199 pId != 2103 && pId != 2203 && pId != 3103 &&
200 pId != 3203 && pId != 3303 );
202 prod.push_back(
abs(channel.product(
p)));
206 for(
size_t p = 0;
p<prod.size();++
p){
207 strstr <<
"_" << prod[
p];
210 h_br[
pid]->Fill(label.c_str(),0.);
211 h_br_ref[
pid]->Fill(label.c_str(),channel.bRatio());
218 strstr <<
"originVertexRho_" <<
pid;
221 strstr <<
"originVertexZ_" <<
pid;
224 strstr <<
"decayVertexRho_" <<
pid;
227 strstr <<
"decayVertexZ_" <<
pid;
228 h_decayVertexZ[
pid] =
new TH1D(strstr.str().c_str(),strstr.str().c_str(),100,0,400);
239 TFile *
f = TFile::Open(
outputFile.c_str(),
"RECREATE");
241 f->mkdir(
"observed");
242 f->mkdir(
"prediction");
243 for(
size_t i = 0;
i<
pids.size();++
i){
278 iEvent.
getByLabel(
"famosSimHits",simvertices);
283 std::map<size_t,std::vector<size_t> > childMap;
284 std::map<size_t,int> parentMap;
285 for(
size_t j = 0;
j<simtracks->size();
j++){
286 childMap[
j] = std::vector<size_t>();
291 for(
size_t j = 0;
j<simtracks->size();
j++){
292 size_t childIndex =
j;
300 childMap[parentIndex].push_back(childIndex);
301 parentMap[childIndex] = int(parentIndex);
305 for(
size_t j = 0;
j<simtracks->size();
j++){
312 double mass = parent.
momentum().M();
324 if(childMap[
j].
size() > 0){
333 for(
std::map<
size_t,std::vector<size_t> >::iterator it = childMap.begin();it != childMap.end();++it){
336 size_t parentIndex = it->first;
339 vector<size_t> & childIndices = it->second;
340 if(childIndices.size() == 0)
348 const SimTrack & child0 = simtracks->at(childIndices[0]);
353 TLorentzVector lv_dist = lv_decay_vertex - lv_origin_vertex;
355 TVector3 boost = lv_parent.BoostVector();
356 lv_dist.Boost(-boost);
357 h_v[
pid]->Fill(boost.Mag());
358 double plt = lv_dist.T();
362 std::vector<int>
prod;
363 for(
size_t d = 0;d<childIndices.size();++d){
364 prod.push_back(
abs(simtracks->at(childIndices[d]).type()));
367 std::stringstream strstr;
368 for(
size_t p = 0;
p<prod.size();++
p){
369 strstr <<
"_" << prod[
p];
374 h_br[
pid]->Fill(label.c_str(),1.);
T getParameter(std::string const &) const
std::map< int, std::vector< string > > knownDecayModes
#define DEFINE_FWK_MODULE(type)
virtual void endJob() override
std::map< int, TH1D * > h_v
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
virtual void analyze(const edm::Event &, const edm::EventSetup &) override
std::map< int, TH1D * > h_decayVertexZ
void addDefault(ParameterSetDescription const &psetDescription)
std::map< int, TH1D * > h_br_ref
std::map< int, TH1D * > h_plt_ref
std::map< int, TH1D * > h_decayVertexRho
std::map< int, TH1D * > h_originVertexZ
Abs< T >::type abs(const T &t)
const math::XYZTLorentzVectorD & position() const
HepPDT::ParticleData ParticleData
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
std::map< int, TH1D * > h_plt
TestPythiaDecays(const edm::ParameterSet &)
int type() const
particle type (HEP PDT convension)
std::map< int, TH1D * > h_br
const math::XYZTLorentzVectorD & momentum() const
std::map< int, TH1D * > h_mass
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::map< int, TH1D * > h_p
std::map< int, TH1D * > h_mass_ref
virtual void beginJob() override
std::map< int, TH1D * > h_originVertexRho
tuple size
Write out results.