CharmHadronization.cxx
Go to the documentation of this file.
1 //____________________________________________________________________________
2 /*
3  Copyright (c) 2003-2019, The GENIE Collaboration
4  For the full text of the license visit http://copyright.genie-mc.org
5  or see $GENIE/LICENSE
6 
7  Author: Costas Andreopoulos <costas.andreopoulos \at stfc.ac.uk>
8  University of Liverpool & STFC Rutherford Appleton Lab
9 
10  Changes required to implement the GENIE Boosted Dark Matter module
11  were installed by Josh Berger (Univ. of Wisconsin)
12 
13 */
14 //____________________________________________________________________________
15 
16 #include <RVersion.h>
17 #include <TPythia6.h>
18 #include <TVector3.h>
19 #include <TF1.h>
20 #include <TROOT.h>
21 
22 #if ROOT_VERSION_CODE >= ROOT_VERSION(5,15,6)
23 #include <TMCParticle.h>
24 #else
25 #include <TMCParticle6.h>
26 #endif
27 
49 
50 using namespace genie;
51 using namespace genie::constants;
52 using namespace genie::controls;
53 
54 extern "C" void py1ent_(int *, int *, double *, double *, double *);
55 extern "C" void py2ent_(int *, int *, int *, double *);
56 
57 //____________________________________________________________________________
59 EventRecordVisitorI("genie::CharmHadronization")
60 {
61  this->Initialize();
62 }
63 //____________________________________________________________________________
65 EventRecordVisitorI("genie::CharmHadronization", config)
66 {
67  this->Initialize();
68 }
69 //____________________________________________________________________________
71 {
72  delete fCharmPT2pdf;
73  fCharmPT2pdf = 0;
74 
75  delete fD0FracSpl;
76  fD0FracSpl = 0;
77 
78  delete fDpFracSpl;
79  fDpFracSpl = 0;
80 
81  delete fDsFracSpl;
82  fDsFracSpl = 0;
83 }
84 //____________________________________________________________________________
86 {
87  fPythia = TPythia6::Instance();
88 }
89 //____________________________________________________________________________
91 {
92  Interaction * interaction = event->Summary();
93  TClonesArray * particle_list = this->Hadronize(interaction);
94 
95  if(! particle_list ) {
96  LOG("CharmHadronization", pWARN) << "Got an empty particle list. Hadronizer failed!";
97  LOG("CharmHadronization", pWARN) << "Quitting the current event generation thread";
98 
99  event->EventFlags()->SetBitNumber(kHadroSysGenErr, true);
100 
102  exception.SetReason("Could not simulate the hadronic system");
103  exception.SwitchOnFastForward();
104  throw exception;
105 
106  return;
107  }
108 
109  int mom = event->FinalStateHadronicSystemPosition();
110  assert(mom!=-1);
111 
112  // find the proper status for the particles we are going to put in event record
113  bool is_nucleus = interaction->InitState().Tgt().IsNucleus();
114  GHepStatus_t istfin = (is_nucleus) ?
116 
117  // retrieve the hadronic blob lorentz boost
118  // Because Hadronize() returned particles not in the LAB reference frame
119  const TLorentzVector * had_syst = event -> Particle(mom) -> P4() ;
120  TVector3 beta( 0., 0., had_syst -> P()/ had_syst -> Energy() ) ;
121 
122  // Vector defining rotation from LAB to LAB' (z:= \vec{phad})
123  TVector3 unitvq = had_syst -> Vect().Unit();
124 
125  GHepParticle * neutrino = event->Probe();
126  const TLorentzVector & vtx = *(neutrino->X4());
127 
128  GHepParticle * particle = 0;
129  TIter particle_iter(particle_list);
130  while ((particle = (GHepParticle *) particle_iter.Next())) {
131 
132  int pdgc = particle -> Pdg() ;
133 
134  // bring the particle in the LAB reference frame
135  particle -> P4() -> Boost( beta ) ;
136  particle -> P4() -> RotateUz( unitvq ) ;
137 
138  // set the proper status according to a number of things:
139  // interaction on a nucleaus or nucleon, particle type
140  GHepStatus_t ist = ( particle -> Status() ==1 ) ? istfin : kIStDISPreFragmHadronicState;
141 
142  // handle gammas, and leptons that might come from internal pythia decays
143  // mark them as final state particles
144  bool not_hadron = ( pdgc == kPdgGamma ||
145  pdg::IsNeutralLepton(pdgc) ||
146  pdg::IsChargedLepton(pdgc) ) ;
147 
148  if(not_hadron) { ist = kIStStableFinalState; }
149  particle -> SetStatus( ist ) ;
150 
151  int im = mom + 1 + particle -> FirstMother() ;
152  int ifc = ( particle -> FirstDaughter() == -1) ? -1 : mom + 1 + particle -> FirstDaughter();
153  int ilc = ( particle -> LastDaughter() == -1) ? -1 : mom + 1 + particle -> LastDaughter();
154 
155  particle -> SetFirstMother( im ) ;
156  if ( ifc > -1 ) {
157  particle -> SetFirstDaughter( ifc ) ;
158  particle -> SetLastDaughter( ilc ) ;
159  }
160 
161  // the Pythia particle position is overridden
162  particle -> SetPosition( vtx ) ;
163 
164  event->AddParticle(*particle);
165  }
166 
167  particle_list -> Delete() ;
168  delete particle_list ;
169 
170  // update the weight of the event
171  event -> SetWeight ( Weight() * event->Weight() );
172 
173 }
174 //____________________________________________________________________________
176  const Interaction * interaction) const
177 {
178  LOG("CharmHad", pNOTICE) << "** Running CHARM hadronizer";
179 
180  PDGLibrary * pdglib = PDGLibrary::Instance();
181  RandomGen * rnd = RandomGen::Instance();
182 
183  // ....................................................................
184  // Get information on the input event
185  //
186  const InitialState & init_state = interaction -> InitState();
187  const Kinematics & kinematics = interaction -> Kine();
188  const Target & target = init_state.Tgt();
189 
190  const TLorentzVector & p4Had = kinematics.HadSystP4();
191 
192  double Ev = init_state.ProbeE(kRfLab);
193  double W = kinematics.W(true);
194 
195  TVector3 beta = -1 * p4Had.BoostVector(); // boost vector for LAB' -> HCM'
196  TLorentzVector p4H(0,0,0,W); // hadronic system 4p @ HCM'
197 
198  double Eh = p4Had.Energy();
199 
200  LOG("CharmHad", pNOTICE) << "Ehad (LAB) = " << Eh << ", W = " << W;
201 
202  int nu_pdg = init_state.ProbePdg();
203  int nuc_pdg = target.HitNucPdg();
204 //int qpdg = target.HitQrkPdg();
205 //bool sea = target.HitSeaQrk();
206  bool isp = pdg::IsProton (nuc_pdg);
207  bool isn = pdg::IsNeutron(nuc_pdg);
208  bool isnu = pdg::IsNeutrino(nu_pdg);
209  bool isnub = pdg::IsAntiNeutrino(nu_pdg);
210  bool isdm = pdg::IsDarkMatter(nu_pdg);
211 
212  // ....................................................................
213  // Attempt to generate a charmed hadron & its 4-momentum
214  //
215 
216  TLorentzVector p4C(0,0,0,0);
217  int ch_pdg = -1;
218 
219  bool got_charmed_hadron = false;
220  unsigned int itry=0;
221 
222  while(itry++ < kRjMaxIterations && !got_charmed_hadron) {
223 
224  // Generate a charmed hadron PDG code
225  int pdg = this->GenerateCharmHadron(nu_pdg,Ev); // generate hadron
226  double mc = pdglib->Find(pdg)->Mass(); // lookup mass
227 
228  LOG("CharmHad", pNOTICE)
229  << "Trying charm hadron = " << pdg << "(m = " << mc << ")";
230 
231  if(mc>=W) continue; // dont' accept
232 
233  // Generate the charmed hadron energy based on the input
234  // fragmentation function
235  double z = fFragmFunc->GenerateZ(); // generate z(=Eh/Ev)
236  double Ec = z*Eh; // @ LAB'
237  double mc2 = TMath::Power(mc,2);
238  double Ec2 = TMath::Power(Ec,2);
239  double pc2 = Ec2-mc2;
240 
241  LOG("CharmHad", pINFO)
242  << "Trying charm hadron z = " << z << ", E = " << Ec;
243 
244  if(pc2<=0) continue;
245 
246  // Generate the charm hadron pT^2 and pL^2 (with respect to the
247  // hadronic system direction @ the LAB)
248  double ptc2 = fCharmPT2pdf->GetRandom();
249  double plc2 = Ec2 - ptc2 - mc2;
250  LOG("CharmHad", pINFO)
251  << "Trying charm hadron pT^2 (tranv to pHad) = " << ptc2;
252  if(plc2<0) continue;
253 
254  // Generate the charm hadron momentum components (@ LAB', z:\vec{pHad})
255  double ptc = TMath::Sqrt(ptc2);
256  double plc = TMath::Sqrt(plc2);
257  double phi = (2*kPi) * rnd->RndHadro().Rndm();
258  double pxc = ptc * TMath::Cos(phi);
259  double pyc = ptc * TMath::Sin(phi);
260  double pzc = plc;
261 
262  p4C.SetPxPyPzE(pxc,pyc,pzc,Ec); // @ LAB'
263 
264  // Boost charm hadron 4-momentum from the LAB' to the HCM' frame
265  //
266  LOG("CharmHad", pDEBUG)
267  << "Charm hadron p4 (@LAB') = " << utils::print::P4AsString(&p4C);
268 
269  p4C.Boost(beta);
270 
271  LOG("CharmHad", pDEBUG)
272  << "Charm hadron p4 (@HCM') = " << utils::print::P4AsString(&p4C);
273 
274  // Hadronic non-charm remnant 4p at HCM'
275  TLorentzVector p4 = p4H - p4C;
276  double wr = p4.M();
277  LOG("CharmHad", pINFO)
278  << "Invariant mass of remnant hadronic system= " << wr;
279  if(wr < kNucleonMass + kPionMass + kASmallNum) {
280  LOG("CharmHad", pINFO) << "Too small hadronic remnant mass!";
281  continue;
282  }
283 
284  ch_pdg = pdg;
285  got_charmed_hadron = true;
286 
287  LOG("CharmHad", pNOTICE)
288  << "Generated charm hadron = " << pdg << "(m = " << mc << ")";
289  LOG("CharmHad", pNOTICE)
290  << "Generated charm hadron z = " << z << ", E = " << Ec;
291  }
292 
293  // ....................................................................
294  // Check whether the code above had difficulty generating the charmed
295  // hadron near the W - threshold.
296  // If yes, attempt a phase space decay of a low mass charm hadron + nucleon
297  // pair that maintains the charge.
298  // That's a desperate solution but don't want to quit too early as that
299  // would distort the generated dsigma/dW distribution near threshold.
300  //
301  bool used_lowW_strategy = false;
302  int fs_nucleon_pdg = -1;
303  if(ch_pdg==-1 && W < 3.){
304  LOG("CharmHad", pNOTICE)
305  << "Had difficulty generating charm hadronic system near W threshold";
306  LOG("CharmHad", pNOTICE)
307  << "Trying an alternative strategy";
308 
309  double qfsl = interaction->FSPrimLepton()->Charge() / 3.;
310  double qinit = pdglib->Find(nuc_pdg)->Charge() / 3.;
311  int qhad = (int) (qinit - qfsl);
312 
313  int remn_pdg = -1;
314  int chrm_pdg = -1;
315 
316  //cc-only: qhad(nu) = +1,+2, qhad(nubar)= -1,0
317  //
318  if(qhad == 2) {
319  chrm_pdg = kPdgDP; remn_pdg = kPdgProton;
320  } else if(qhad == 1) {
321  if(rnd->RndHadro().Rndm() > 0.5) {
322  chrm_pdg = kPdgD0; remn_pdg = kPdgProton;
323  } else {
324  chrm_pdg = kPdgDP; remn_pdg = kPdgNeutron;
325  }
326  } else if(qhad == 0) {
327  chrm_pdg = kPdgAntiD0; remn_pdg = kPdgNeutron;
328  } else if(qhad == -1) {
329  chrm_pdg = kPdgDM; remn_pdg = kPdgNeutron;
330  }
331 
332  double mc = pdglib->Find(chrm_pdg)->Mass();
333  double mn = pdglib->Find(remn_pdg)->Mass();
334 
335  if(mc+mn < W) {
336  // Set decay
337  double mass[2] = {mc, mn};
338  bool permitted = fPhaseSpaceGenerator.SetDecay(p4H, 2, mass);
339  assert(permitted);
340 
341  // Get the maximum weight
342  double wmax = -1;
343  for(int i=0; i<200; i++) {
344  double w = fPhaseSpaceGenerator.Generate();
345  wmax = TMath::Max(wmax,w);
346  }
347 
348  if(wmax>0) {
349  wmax *= 2;
350 
351  // Generate unweighted decay
352  bool accept_decay=false;
353  unsigned int idecay_try=0;
354  while(!accept_decay)
355  {
356  idecay_try++;
357 
358  if(idecay_try>kMaxUnweightDecayIterations) {
359  LOG("CharmHad", pWARN)
360  << "Couldn't generate an unweighted phase space decay after "
361  << idecay_try << " attempts";
362  }
363  double w = fPhaseSpaceGenerator.Generate();
364  if(w > wmax) {
365  LOG("CharmHad", pWARN)
366  << "Decay weight = " << w << " > max decay weight = " << wmax;
367  }
368  double gw = wmax * rnd->RndHadro().Rndm();
369  accept_decay = (gw<=w);
370 
371  if(accept_decay) {
372  used_lowW_strategy = true;
373  TLorentzVector * p4 = fPhaseSpaceGenerator.GetDecay(0);
374  p4C = *p4;
375  ch_pdg = chrm_pdg;
376  fs_nucleon_pdg = remn_pdg;
377  }
378  } // decay loop
379  }//wmax>0
380 
381  }// allowed decay
382  } // alt low-W strategy
383 
384  // ....................................................................
385  // Check success in generating the charm hadron & compute 4p for
386  // remnant system
387  //
388  if(ch_pdg==-1){
389  LOG("CharmHad", pWARN)
390  << "Couldn't generate charm hadron for: " << *interaction;
391  return 0;
392  }
393 
394  TLorentzVector p4R = p4H - p4C;
395  double WR = p4R.M();
396  //double MC = pdglib->Find(ch_pdg)->Mass();
397 
398  LOG("CharmHad", pNOTICE) << "Remnant hadronic system mass = " << WR;
399 
400  // ....................................................................
401  // Handle case where the user doesn't want to remnant system to be
402  // hadronized (add as 'hadronic blob')
403  //
404  if(fCharmOnly) {
405  // Create particle list (fragmentation record)
406  TClonesArray * particle_list = new TClonesArray("genie::GHepParticle", 2);
407  particle_list->SetOwner(true);
408 
409  // insert the generated particles
410  new ((*particle_list)[0]) GHepParticle (ch_pdg,kIStStableFinalState,
411  -1,-1,-1,-1, p4C.Px(),p4C.Py(),p4C.Pz(),p4C.E(), 0,0,0,0);
412  new ((*particle_list)[1]) GHepParticle (kPdgHadronicBlob,kIStStableFinalState,
413  -1,-1,-1,-1, p4R.Px(),p4R.Py(),p4R.Pz(),p4R.E(), 0,0,0,0);
414 
415  return particle_list;
416  }
417 
418  // ....................................................................
419  // Handle case where the remnant system is already known and doesn't
420  // have to be hadronized. That happens when (close to the W threshold)
421  // the hadronic system was generated by a simple 2-body decay
422  //
423  if(used_lowW_strategy) {
424  // Create particle list (fragmentation record)
425  TClonesArray * particle_list = new TClonesArray("genie::GHepParticle", 3);
426  particle_list->SetOwner(true);
427 
428  // insert the generated particles
429  new ((*particle_list)[0]) GHepParticle (ch_pdg,kIStStableFinalState,
430  -1,-1,-1,-1, p4C.Px(),p4C.Py(),p4C.Pz(),p4C.E(), 0,0,0,0);
431  new ((*particle_list)[1]) GHepParticle (kPdgHadronicBlob,kIStNucleonTarget,
432  -1,-1,2,2, p4R.Px(),p4R.Py(),p4R.Pz(),p4R.E(), 0,0,0,0);
433  new ((*particle_list)[2]) GHepParticle (fs_nucleon_pdg,kIStStableFinalState,
434  1,1,-1,-1, p4R.Px(),p4R.Py(),p4R.Pz(),p4R.E(), 0,0,0,0);
435 
436  return particle_list;
437  }
438 
439  // ....................................................................
440  // --------------------------------------------------------------------
441  // Hadronize non-charm hadronic blob using PYTHIA/JETSET
442  // --------------------------------------------------------------------
443  // ....................................................................
444 
445  // Create output event record
446  // Insert the generated charm hadron & the hadronic (non-charm) blob.
447  // In this case the hadronic blob is entered as a pre-fragm. state.
448 
449  TClonesArray * particle_list = new TClonesArray("genie::GHepParticle");
450  particle_list->SetOwner(true);
451 
452  new ((*particle_list)[0]) GHepParticle (ch_pdg,kIStStableFinalState,
453  -1,-1,-1,-1, p4C.Px(),p4C.Py(),p4C.Pz(),p4C.E(), 0,0,0,0);
454  new ((*particle_list)[1]) GHepParticle (kPdgHadronicBlob,kIStNucleonTarget,
455  -1,-1,2,3, p4R.Px(),p4R.Py(),p4R.Pz(),p4R.E(), 0,0,0,0);
456 
457  unsigned int rpos =2; // offset in event record
458 
459  bool use_pythia = (WR>1.5);
460 
461 /*
462  // Determining quark systems to input to PYTHIA based on simple quark model
463  // arguments
464  //
465  // Neutrinos
466  // ------------------------------------------------------------------
467  // Scattering off valence q
468  // ..................................................................
469  // p: [uu]+d
470  // |--> c --> D0 <c+\bar(u)> : [u]
471  // --> D+ <c+\bar(d)> : [d]
472  // --> Ds+ <c+\bar(s)> : [s]
473  // --> Lamda_c+ <c+ud > : [\bar(ud)]
474  //
475  // (for n: [uu] -> 50%[ud]_{0} + 50%[ud]_{1})
476  //
477  // Scattering off sea q
478  // ..................................................................
479  // p: [uud] + [\bar(d)]d (or)
480  // [\bar(s)]s
481  // |--> c --> D0 <c+\bar(u)> : [u]
482  // --> D+ <c+\bar(d)> : [d]
483  // --> Ds+ <c+\bar(s)> : [s]
484  // --> Lamda_c+ <c+ud > : [\bar(ud)]
485  // Anti-Neutrinos
486  // ------------------------------------------------------------------
487  // Scattering off sea q
488  // ..................................................................
489  // p: [uud] + [d] \bar(d) (or)
490  // [s] \bar(s)
491  // |----> \bar(c) --> \bar(D0) <\bar(c)+u> : [\bar(u)]
492  // --> D- <\bar(c)+d> : [\bar(d)]
493  // --> Ds- <\bar(c)+s> : [\bar(s)]
494  // [Summary]
495  // Qq
496  // | v + p [val/d] --> D0 + { u uu }(+2) / u,uu
497  // | v + p [val/d] --> D+ + { d uu }(+1) / d,uu
498  // | v + p [val/d] --> Ds+ + { s uu }(+1) / s,uu
499  // | v + p [val/d] --> Lc+ + { \bar(ud) uu }(+1) / \bar(d),u
500  // | v + n [val/d] --> D0 + { u ud }(+1) / u,ud
501  // | v + n [val/d] --> D+ + { d ud }( 0) / d,ud
502  // | v + n [val/d] --> Ds+ + { s ud }( 0) / s,ud
503  // | v + n [val/d] --> Lc+ + { \bar(ud) ud }( 0) / \bar(d),d
504  // | v + p [sea/d] --> D0 + { uud \bar(d) u }(+2) / u,uu
505  // | v + p [sea/d] --> D+ + { uud \bar(d) d }(+1) / d,uu
506  // | v + p [sea/d] --> Ds+ + { uud \bar(d) s }(+1) / s,uu
507  // | v + p [sea/d] --> Lc+ + { uud \bar(d) \bar(ud) }(+1) / \bar(d),u
508  // | v + n [sea/d] --> D0 + { udd \bar(d) u }(+1) / u,ud
509  // | v + n [sea/d] --> D+ + { udd \bar(d) d }( 0) / d,ud
510  // | v + n [sea/d] --> Ds+ + { udd \bar(d) s }( 0) / s,ud
511  // | v + n [sea/d] --> Lc+ + { udd \bar(d) \bar(ud) }( 0) / \bar(d),d
512  // | v + p [sea/s] --> D0 + { uud \bar(s) u }(+2) / u,uu
513  // | v + p [sea/s] --> D+ + { uud \bar(s) d }(+1) / d,uu
514  // | v + p [sea/s] --> Ds+ + { uud \bar(s) s }(+1) / s,uu
515  // | v + p [sea/s] --> Lc+ + { uud \bar(s) \bar(ud) }(+1) / \bar(d),u
516  // | v + n [sea/s] --> D0 + { udd \bar(s) u }(+1) / u,ud
517  // | v + n [sea/s] --> D+ + { udd \bar(s) d }( 0) / d,ud
518  // | v + n [sea/s] --> Ds+ + { udd \bar(s) s }( 0) / s,ud
519  // | v + n [sea/s] --> Lc+ + { udd \bar(s) \bar(ud) }( 0) / \bar(d),d
520 
521  // | \bar(v) + p [sea/\bar(d)] --> \bar(D0) + { uud d \bar(u) }( 0) / d,ud
522  // | \bar(v) + p [sea/\bar(d)] --> D- + { uud d \bar(d) }(+1) / d,uu
523  // | \bar(v) + p [sea/\bar(d)] --> Ds- + { uud d \bar(s) }(+1) / d,uu
524  // | \bar(v) + n [sea/\bar(d)] --> \bar(D0) + { udd d \bar(u) }(-1) / d,dd
525  // | \bar(v) + n [sea/\bar(d)] --> D- + { udd d \bar(d) }( 0) / d,ud
526  // | \bar(v) + n [sea/\bar(d)] --> Ds- + { udd d \bar(s) }( 0) / d,ud
527  // | \bar(v) + p [sea/\bar(s)] --> \bar(D0) + { uud s \bar(u) }( 0) / d,ud
528  // | \bar(v) + p [sea/\bar(s)] --> D- + { uud s \bar(d) }(+1) / d,uu
529  // | \bar(v) + p [sea/\bar(s)] --> Ds- + { uud s \bar(s) }(+1) / d,uu
530  // | \bar(v) + n [sea/\bar(s)] --> \bar(D0) + { udd s \bar(u) }(-1) / d,dd
531  // | \bar(v) + n [sea/\bar(s)] --> D- + { udd s \bar(d) }( 0) / d,ud
532  // | \bar(v) + n [sea/\bar(s)] --> Ds- + { udd s \bar(s) }( 0) / d,ud
533  //
534  //
535  // Taking some short-cuts below :
536  // In the process of obtaining 2 q systems (while conserving the charge) I might tread
537  // d,s or \bar(d),\bar(s) as the same
538  // In the future I should perform the first steps of the multi-q hadronization manualy
539  // (to reduce the number of q's input to PYTHIA) or use py3ent_(), py4ent_() ...
540  //
541 */
542 
543  if(use_pythia) {
544  int qrkSyst1 = 0;
545  int qrkSyst2 = 0;
546  if(isnu||isdm) { // neutrinos
547  if(ch_pdg==kPdgLambdaPc) {
548  if(isp) { qrkSyst1 = kPdgAntiDQuark; qrkSyst2 = kPdgUQuark; }
549  if(isn) { qrkSyst1 = kPdgAntiDQuark; qrkSyst2 = kPdgDQuark; }
550  } else {
551  if(isp) { qrkSyst2 = kPdgUUDiquarkS1; }
552  if(isn) { qrkSyst2 = (rnd->RndHadro().Rndm()<0.5) ? kPdgUDDiquarkS0 : kPdgUDDiquarkS1; }
553  if (ch_pdg==kPdgD0 ) { qrkSyst1 = kPdgUQuark; }
554  if (ch_pdg==kPdgDP ) { qrkSyst1 = kPdgDQuark; }
555  if (ch_pdg==kPdgDPs ) { qrkSyst1 = kPdgSQuark; }
556  }
557  }
558  if(isnub) { // antineutrinos
559  qrkSyst1 = kPdgDQuark;
560  if (isp && ch_pdg==kPdgAntiD0) { qrkSyst2 = (rnd->RndHadro().Rndm()<0.5) ? kPdgUDDiquarkS0 : kPdgUDDiquarkS1; }
561  if (isp && ch_pdg==kPdgDM ) { qrkSyst2 = kPdgUUDiquarkS1; }
562  if (isp && ch_pdg==kPdgDMs ) { qrkSyst2 = kPdgUUDiquarkS1; }
563  if (isn && ch_pdg==kPdgAntiD0) { qrkSyst2 = kPdgDDDiquarkS1; }
564  if (isn && ch_pdg==kPdgDM ) { qrkSyst2 = (rnd->RndHadro().Rndm()<0.5) ? kPdgUDDiquarkS0 : kPdgUDDiquarkS1; }
565  if (isn && ch_pdg==kPdgDMs ) { qrkSyst2 = (rnd->RndHadro().Rndm()<0.5) ? kPdgUDDiquarkS0 : kPdgUDDiquarkS1; }
566  }
567  if(qrkSyst1 == 0 && qrkSyst2 == 0) {
568  LOG("CharmHad", pWARN)
569  << "Couldn't generate quark systems for PYTHIA in: " << *interaction;
570  return 0;
571  }
572 
573  //
574  // Run PYTHIA for the hadronization of remnant system
575  //
576  fPythia->SetMDCY(fPythia->Pycomp(kPdgPi0), 1,0); // don't decay pi0
577  fPythia->SetMDCY(fPythia->Pycomp(kPdgP33m1232_DeltaM), 1,1); // decay Delta+
578  fPythia->SetMDCY(fPythia->Pycomp(kPdgP33m1232_Delta0), 1,1); // decay Delta++
579  fPythia->SetMDCY(fPythia->Pycomp(kPdgP33m1232_DeltaP), 1,1); // decay Delta++
580  fPythia->SetMDCY(fPythia->Pycomp(kPdgP33m1232_DeltaPP), 1,1); // decay Delta++
581 // fPythia->SetMDCY(fPythia->Pycomp(kPdgDeltaP), 1,1); // decay Delta+
582 // fPythia->SetMDCY(fPythia->Pycomp(kPdgDeltaPP), 1,1); // decay Delta++
583 
584  int ip = 0;
585  py2ent_(&ip, &qrkSyst1, &qrkSyst2, &WR); // hadronize
586 
587  fPythia->SetMDCY(fPythia->Pycomp(kPdgPi0),1,1); // restore
588 
589  //-- Get PYTHIA's LUJETS event record
590  TClonesArray * pythia_remnants = 0;
591  fPythia->GetPrimaries();
592  pythia_remnants = dynamic_cast<TClonesArray *>(fPythia->ImportParticles("All"));
593  if(!pythia_remnants) {
594  LOG("CharmHad", pWARN) << "Couldn't hadronize (non-charm) remnants!";
595  return 0;
596  }
597 
598  int np = pythia_remnants->GetEntries();
599  assert(np>0);
600 
601  // PYTHIA performs the hadronization at the *remnant hadrons* centre of mass
602  // frame (not the hadronic centre of mass frame).
603  // Boost all hadronic blob fragments to the HCM', fix their mother/daughter
604  // assignments and add them to the fragmentation record.
605 
606  TVector3 rmnbeta = +1 * p4R.BoostVector(); // boost velocity
607 
608  TMCParticle * pythia_remn = 0; // remnant
609  GHepParticle * bremn = 0; // boosted remnant
610  TIter remn_iter(pythia_remnants);
611  while( (pythia_remn = (TMCParticle *) remn_iter.Next()) ) {
612 
613  // insert and get a pointer to inserted object for mods
614  bremn = new ((*particle_list)[rpos++]) GHepParticle ( pythia_remn->GetKF(), // pdg
615  GHepStatus_t(pythia_remn->GetKS()), // status
616  pythia_remn->GetParent(), // first parent
617  -1, // second parent
618  pythia_remn->GetFirstChild(), // first daughter
619  pythia_remn->GetLastChild(), // second daughter
620  pythia_remn -> GetPx(), // px
621  pythia_remn -> GetPy(), // py
622  pythia_remn -> GetPz(), // pz
623  pythia_remn -> GetEnergy(), // e
624  pythia_remn->GetVx(), // x
625  pythia_remn->GetVy(), // y
626  pythia_remn->GetVz(), // z
627  pythia_remn->GetTime() // t
628  );
629 
630  // boost
631  bremn -> P4() -> Boost( rmnbeta ) ;
632 
633  // handle insertion of charmed hadron
634  int jp = bremn->FirstMother();
635  int ifc = bremn->FirstDaughter();
636  int ilc = bremn->LastDaughter();
637 
638  bremn -> SetFirstMother( (jp == 0 ? 1 : jp +1) );
639  bremn -> SetFirstDaughter ( (ifc == 0 ? -1 : ifc+1) );
640  bremn -> SetLastDaughter ( (ilc == 0 ? -1 : ilc+1) );
641  }
642  } // use_pythia
643 
644  // ....................................................................
645  // Hadronizing low-W non-charm hadronic blob using a phase space decay
646  // ....................................................................
647 
648  else {
649  // Just a small fraction of events (low-W remnant syste) causing trouble
650  // to PYTHIA/JETSET
651  // Set a 2-body N+pi system that matches the remnant system charge and
652  // do a simple phase space decay
653  //
654  // q(remn) remn/syst
655  // +2 : (p pi+)
656  // +1 : 50%(p pi0) + 50% n pi+
657  // 0 : 50%(p pi-) + 50% n pi0
658  // -1 : (n pi-)
659  //
660  double qfsl = interaction->FSPrimLepton()->Charge() / 3.;
661  double qinit = pdglib->Find(nuc_pdg)->Charge() / 3.;
662  double qch = pdglib->Find(ch_pdg)->Charge() / 3.;
663  int Q = (int) (qinit - qfsl - qch); // remnant hadronic system charge
664 
665  bool allowdup=true;
666  PDGCodeList pd(allowdup);
667  if(Q==2) {
669  else if (Q==1) {
670  if(rnd->RndHadro().Rndm()<0.5) {
672  else {
674  }
675  else if (Q==0) {
676  if(rnd->RndHadro().Rndm()<0.5) {
678  else {
680  }
681  else if (Q==-1) {
683 
684  double mass[2] = {
685  pdglib->Find(pd[0])->Mass(), pdglib->Find(pd[1])->Mass()
686  };
687 
688  // Set the decay
689  bool permitted = fPhaseSpaceGenerator.SetDecay(p4R, 2, mass);
690  if(!permitted) {
691  LOG("CharmHad", pERROR) << " *** Phase space decay is not permitted";
692  return 0;
693  }
694  // Get the maximum weight
695  double wmax = -1;
696  for(int i=0; i<200; i++) {
697  double w = fPhaseSpaceGenerator.Generate();
698  wmax = TMath::Max(wmax,w);
699  }
700  if(wmax<=0) {
701  LOG("CharmHad", pERROR) << " *** Non-positive maximum weight";
702  LOG("CharmHad", pERROR) << " *** Can not generate an unweighted phase space decay";
703  return 0;
704  }
705 
706  LOG("CharmHad", pINFO)
707  << "Max phase space gen. weight @ current hadronic system: " << wmax;
708 
709  // *** generating an un-weighted decay ***
710  wmax *= 1.3;
711  bool accept_decay=false;
712  unsigned int idectry=0;
713  while(!accept_decay)
714  {
715  idectry++;
716  if(idectry>kMaxUnweightDecayIterations) {
717  // report, clean-up and return
718  LOG("Char,Had", pWARN)
719  << "Couldn't generate an unweighted phase space decay after "
720  << itry << " attempts";
721  return 0;
722  }
723  double w = fPhaseSpaceGenerator.Generate();
724  if(w > wmax) {
725  LOG("CharmHad", pWARN)
726  << "Decay weight = " << w << " > max decay weight = " << wmax;
727  }
728  double gw = wmax * rnd->RndHadro().Rndm();
729  accept_decay = (gw<=w);
730  LOG("CharmHad", pDEBUG)
731  << "Decay weight = " << w << " / R = " << gw << " - accepted: " << accept_decay;
732  }
733  for(unsigned int i=0; i<2; i++) {
734  int pdgc = pd[i];
735  TLorentzVector * p4d = fPhaseSpaceGenerator.GetDecay(i);
736  new ( (*particle_list)[rpos+i] ) GHepParticle(
737  pdgc,kIStStableFinalState,1,1,-1,-1,p4d->Px(),p4d->Py(),p4d->Pz(),p4d->Energy(),
738  0,0,0,0);
739  }
740  }
741 
742  //-- Print & return the fragmentation record
743  //utils::fragmrec::Print(particle_list);
744  return particle_list;
745 }
746 //____________________________________________________________________________
747 int CharmHadronization::GenerateCharmHadron(int nu_pdg, double EvLab) const
748 {
749  // generate a charmed hadron pdg code using a charm fraction table
750 
751  RandomGen * rnd = RandomGen::Instance();
752  double r = rnd->RndHadro().Rndm();
753 
754  if(pdg::IsNeutrino(nu_pdg)) {
755  double tf = 0;
756  if (r < (tf+=fD0FracSpl->Evaluate(EvLab))) return kPdgD0; // D^0
757  else if (r < (tf+=fDpFracSpl->Evaluate(EvLab))) return kPdgDP; // D^+
758  else if (r < (tf+=fDsFracSpl->Evaluate(EvLab))) return kPdgDPs; // Ds^+
759  else return kPdgLambdaPc; // Lamda_c^+
760 
761  } else if(pdg::IsAntiNeutrino(nu_pdg)) {
762  if (r < fD0BarFrac) return kPdgAntiD0;
763  else if (r < fD0BarFrac+fDmFrac) return kPdgDM;
764  else return kPdgDMs;
765  }
766 
767  LOG("CharmHad", pERROR) << "Could not generate a charm hadron!";
768  return 0;
769 }
770 //____________________________________________________________________________
771 double CharmHadronization::Weight(void) const
772 {
773  return 1. ;
774 }
775 
776 //____________________________________________________________________________
778 {
779  Algorithm::Configure(config);
780  this->LoadConfig();
781 }
782 //____________________________________________________________________________
784 {
785  Algorithm::Configure(config);
786  this->LoadConfig();
787 }
788 //____________________________________________________________________________
790 {
791 
792  bool hadronize_remnants ;
793  GetParamDef( "HadronizeRemnants", hadronize_remnants, true ) ;
794 
795  fCharmOnly = ! hadronize_remnants ;
796 
797  //-- Get a fragmentation function
798  fFragmFunc = dynamic_cast<const FragmentationFunctionI *> (
799  this->SubAlg("FragmentationFunc"));
800  assert(fFragmFunc);
801 
802  string pt_function ;
803  this -> GetParam( "PTFunction", pt_function ) ;
804 
805  fCharmPT2pdf = new TF1("fCharmPT2pdf", pt_function.c_str(),0,0.6);
806 
807  // stop ROOT from deleting this object of its own volition
808  gROOT->GetListOfFunctions()->Remove(fCharmPT2pdf);
809 
810  // neutrino charm fractions: D^0, D^+, Ds^+ (remainder: Lamda_c^+)
811  std::vector<double> ec, d0frac, dpfrac, dsfrac ;
812 
813  std::string raw ;
814  std::vector<std::string> bits ;
815 
816  bool invalid_configuration = false ;
817 
818  // load energy points
819  this -> GetParam( "CharmFrac-E", raw ) ;
820  bits = utils::str::Split( raw, ";" ) ;
821 
822  if ( ! utils::str::Convert(bits, ec) ) {
823  LOG("CharmHadronization", pFATAL) <<
824  "Failed to decode CharmFrac-E string: ";
825  LOG("CharmHadronization", pFATAL) << "string: "<< raw ;
826  invalid_configuration = true ;
827  }
828 
829  // load D0 fractions
830  this -> GetParam( "CharmFrac-D0", raw ) ;
831  bits = utils::str::Split( raw, ";" ) ;
832 
833  if ( ! utils::str::Convert(bits, d0frac) ) {
834  LOG("CharmHadronization", pFATAL) <<
835  "Failed to decode CharmFrac-D0 string: ";
836  LOG("CharmHadronization", pFATAL) << "string: "<< raw ;
837  invalid_configuration = true ;
838  }
839 
840  // check the size
841  if ( d0frac.size() != ec.size() ) {
842  LOG("CharmHadronization", pFATAL) << "E entries don't match D0 fraction entries";
843  LOG("CharmHadronization", pFATAL) << "E: " << ec.size() ;
844  LOG("CharmHadronization", pFATAL) << "D0: " << d0frac.size() ;
845  invalid_configuration = true ;
846  }
847 
848  // load D+ fractions
849  this -> GetParam( "CharmFrac-D+", raw ) ;
850  bits = utils::str::Split( raw, ";" ) ;
851 
852  if ( ! utils::str::Convert(bits, dpfrac) ) {
853  LOG("CharmHadronization", pFATAL) <<
854  "Failed to decode CharmFrac-D+ string: ";
855  LOG("CharmHadronization", pFATAL) << "string: "<< raw ;
856  invalid_configuration = true ;
857  }
858 
859  // check the size
860  if ( dpfrac.size() != ec.size() ) {
861  LOG("CharmHadronization", pFATAL) << "E entries don't match D+ fraction entries";
862  LOG("CharmHadronization", pFATAL) << "E: " << ec.size() ;
863  LOG("CharmHadronization", pFATAL) << "D+: " << dpfrac.size() ;
864  invalid_configuration = true ;
865  }
866 
867  // load D_s fractions
868  this -> GetParam( "CharmFrac-Ds", raw ) ;
869  bits = utils::str::Split( raw, ";" ) ;
870 
871  if ( ! utils::str::Convert(bits, dsfrac) ) {
872  LOG("CharmHadronization", pFATAL) <<
873  "Failed to decode CharmFrac-Ds string: ";
874  LOG("CharmHadronization", pFATAL) << "string: "<< raw ;
875  invalid_configuration = true ;
876  }
877 
878  // check the size
879  if ( dsfrac.size() != ec.size() ) {
880  LOG("CharmHadronization", pFATAL) << "E entries don't match Ds fraction entries";
881  LOG("CharmHadronization", pFATAL) << "E: " << ec.size() ;
882  LOG("CharmHadronization", pFATAL) << "Ds: " << dsfrac.size() ;
883  invalid_configuration = true ;
884  }
885 
886  fD0FracSpl = new Spline( ec.size(), & ec[0], & d0frac[0] );
887  fDpFracSpl = new Spline( ec.size(), & ec[0], & dpfrac[0] );
888  fDsFracSpl = new Spline( ec.size(), & ec[0], & dsfrac[0] );
889 
890  // anti-neutrino charm fractions: bar(D^0), D^-, (remainder: Ds^-)
891 
892  this -> GetParam( "CharmFrac-D0bar", fD0BarFrac ) ;
893  this -> GetParam( "CharmFrac-D-", fDmFrac ) ;
894 
895  if ( invalid_configuration ) {
896 
897  LOG("CharmHadronization", pFATAL)
898  << "Invalid configuration: Exiting" ;
899 
900  // From the FreeBSD Library Functions Manual
901  //
902  // EX_CONFIG (78) Something was found in an unconfigured or miscon-
903  // figured state.
904 
905  exit( 78 ) ;
906  }
907 }
908 //____________________________________________________________________________
const int kPdgAntiD0
Definition: PDGCodes.h:165
const int kPdgP33m1232_DeltaPP
Definition: PDGCodes.h:91
const int kPdgDPs
Definition: PDGCodes.h:166
const int kPdgUUDiquarkS1
Definition: PDGCodes.h:58
double W(bool selected=false) const
Definition: Kinematics.cxx:167
double fDmFrac
nubar D- charm fraction
Basic constants.
TPythia6 * fPythia
remnant (non-charm) hadronizer
bool IsNeutrino(int pdgc)
Definition: PDGUtils.cxx:108
THE MAIN GENIE PROJECT NAMESPACE
Definition: AlgCmp.h:26
Spline * fDpFracSpl
nu charm fraction vs Ev: D+
#define pERROR
Definition: Messenger.h:60
static const double kNucleonMass
Definition: Constants.h:78
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:79
int HitNucPdg(void) const
Definition: Target.cxx:321
Defines the EventRecordVisitorI interface. Concrete implementations of this interface use the &#39;Visito...
virtual double GenerateZ(void) const =0
std::string string
Definition: nybbler.cc:12
int FirstDaughter(void) const
Definition: GHepParticle.h:69
const int kPdgHadronicBlob
Definition: PDGCodes.h:192
enum genie::EGHepStatus GHepStatus_t
A numeric analysis tool class for interpolating 1-D functions.
Definition: Spline.h:47
#define pFATAL
Definition: Messenger.h:57
bool IsNucleus(void) const
Definition: Target.cxx:289
const int kPdgUQuark
Definition: PDGCodes.h:42
void py1ent_(int *, int *, double *, double *, double *)
Generated/set kinematical variables for an event.
Definition: Kinematics.h:40
string P4AsString(const TLorentzVector *p)
Definition: PrintUtils.cxx:34
static const unsigned int kMaxUnweightDecayIterations
Definition: Controls.h:62
bool IsDarkMatter(int pdgc)
Definition: PDGUtils.cxx:125
bool IsChargedLepton(int pdgc)
Definition: PDGUtils.cxx:99
virtual double Weight(void) const
Definition: GHepRecord.h:125
const TLorentzVector & HadSystP4(void) const
Definition: Kinematics.h:67
Raw data description.
Definition: tf_graph.h:23
double Evaluate(double x) const
Definition: Spline.cxx:362
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition: RandomGen.h:30
void py2ent_(int *, int *, int *, double *)
const int kPdgSQuark
Definition: PDGCodes.h:46
A list of PDG codes.
Definition: PDGCodeList.h:33
const int kPdgP33m1232_DeltaP
Definition: PDGCodes.h:90
const int kPdgP33m1232_DeltaM
Definition: PDGCodes.h:88
TClonesArray * Hadronize(const Interaction *) const
int FirstMother(void) const
Definition: GHepParticle.h:67
bool IsNeutron(int pdgc)
Definition: PDGUtils.cxx:310
Summary information for an interaction.
Definition: Interaction.h:55
int LastDaughter(void) const
Definition: GHepParticle.h:70
An exception thrown by EventRecordVisitorI when the normal processing sequence has to be disrupted (f...
bool IsProton(int pdgc)
Definition: PDGUtils.cxx:305
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:97
const int kPdgUDDiquarkS1
Definition: PDGCodes.h:57
static Config * config
Definition: config.cpp:1054
const int kPdgGamma
Definition: PDGCodes.h:170
bool IsAntiNeutrino(int pdgc)
Definition: PDGUtils.cxx:116
const int kPdgAntiDQuark
Definition: PDGCodes.h:45
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:41
const FragmentationFunctionI * fFragmFunc
charm hadron fragmentation func
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:71
int ProbePdg(void) const
Definition: InitialState.h:65
double fD0BarFrac
nubar {D0} charm fraction
const int kPdgLambdaPc
Definition: PDGCodes.h:83
double z
const int kPdgPiP
Definition: PDGCodes.h:139
const int kPdgPi0
Definition: PDGCodes.h:141
void Configure(const Registry &config)
const int kPdgDQuark
Definition: PDGCodes.h:44
static const double kASmallNum
Definition: Controls.h:41
#define pINFO
Definition: Messenger.h:63
const int kPdgUDDiquarkS0
Definition: PDGCodes.h:56
Misc GENIE control constants.
#define pWARN
Definition: Messenger.h:61
const int kPdgP33m1232_Delta0
Definition: PDGCodes.h:89
bool Convert(const vector< std::string > &input, std::vector< T > &v)
bool fCharmOnly
don&#39;t hadronize non-charm blob
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
Pure abstract base class. Defines the FragmentationFunctionI interface to be implemented by any algor...
bool IsNeutralLepton(int pdgc)
Definition: PDGUtils.cxx:93
const int kPdgDP
Definition: PDGCodes.h:162
TRandom3 & RndHadro(void) const
rnd number generator used by hadronization models
Definition: RandomGen.h:54
void ProcessEventRecord(GHepRecord *event) const
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:43
vector< string > Split(string input, string delim)
Definition: StringUtils.cxx:42
static const double kPionMass
Definition: Constants.h:74
Singleton class to load & serve a TDatabasePDG.
Definition: PDGLibrary.h:30
const TLorentzVector * X4(void) const
Definition: GHepParticle.h:80
int GenerateCharmHadron(int nupdg, double EvLab) const
A registry. Provides the container for algorithm configuration parameters.
Definition: Registry.h:66
const int kPdgDMs
Definition: PDGCodes.h:167
const int kPdgDM
Definition: PDGCodes.h:163
TParticlePDG * Find(int pdgc)
Definition: PDGLibrary.cxx:61
const int kPdgPiM
Definition: PDGCodes.h:140
Spline * fDsFracSpl
nu charm fraction vs Ev: Ds+
const InitialState & InitState(void) const
Definition: Interaction.h:68
const int kPdgProton
Definition: PDGCodes.h:65
#define pNOTICE
Definition: Messenger.h:62
bool GetParamDef(const RgKey &name, T &p, const T &def) const
bool GetParam(const RgKey &name, T &p, bool is_top_call=true) const
const Target & Tgt(void) const
Definition: InitialState.h:67
static const unsigned int kRjMaxIterations
Definition: Controls.h:27
double ProbeE(RefFrame_t rf) const
const int kPdgNeutron
Definition: PDGCodes.h:67
GENIE&#39;s GHEP MC event record.
Definition: GHepRecord.h:46
static const double kPi
Definition: Constants.h:38
const int kPdgD0
Definition: PDGCodes.h:164
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils...
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:40
Spline * fD0FracSpl
nu charm fraction vs Ev: D0
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
TF1 * fCharmPT2pdf
charm hadron pT^2 pdf
void push_back(int pdg_code)
Definition: PDGCodeList.cxx:67
Event finding and building.
TGenPhaseSpace fPhaseSpaceGenerator
a phase space generator
Initial State information.
Definition: InitialState.h:49
#define pDEBUG
Definition: Messenger.h:64
const int kPdgDDDiquarkS1
Definition: PDGCodes.h:55
const Algorithm * SubAlg(const RgKey &registry_key) const
Definition: Algorithm.cxx:354