MECGenerator.cxx
Go to the documentation of this file.
1 //____________________________________________________________________________
2 /*
3  Copyright (c) 2003-2020, The GENIE Collaboration
4  For the full text of the license visit http://copyright.genie-mc.org
5 
6  Costas Andreopoulos <constantinos.andreopoulos \at cern.ch>
7  University of Liverpool & STFC Rutherford Appleton Laboratory
8 
9  Steve Dytman <dytman+ \at pitt.edu>
10  Pittsburgh University
11 */
12 //____________________________________________________________________________
13 
14 #include <TMath.h>
15 #include "Math/Minimizer.h"
16 #include "Math/Factory.h"
17 
34 
36 //#include "Physics/Multinucleon/XSection/MECHadronTensor.h"
44 
45 using namespace genie;
46 using namespace genie::utils;
47 using namespace genie::constants;
48 using namespace genie::controls;
49 
50 //___________________________________________________________________________
52 EventRecordVisitorI("genie::MECGenerator")
53 {
54 
55 }
56 //___________________________________________________________________________
58 EventRecordVisitorI("genie::MECGenerator", config)
59 {
60 
61 }
62 //___________________________________________________________________________
64 {
65 
66 }
67 //___________________________________________________________________________
69 {
70  //-- Access cross section algorithm for running thread
72  const EventGeneratorI * evg = rtinfo->RunningThread();
73  fXSecModel = evg->CrossSectionAlg();
74  if (fXSecModel->Id().Name() == "genie::EmpiricalMECPXSec2015") {
75  this -> AddTargetRemnant (event); /// shortly, this will be handled by the InitialStateAppender module
76  this -> GenerateFermiMomentum(event);
77  this -> SelectEmpiricalKinematics(event);
78  // TODO: test removing `AddFinalStateLepton` and replacing it with
79  // PrimaryLeptonGenerator::ProcessEventRecord(evrec);
80  this -> AddFinalStateLepton(event);
81  // TODO: maybe put `RecoilNucleonCluster` in a `MECHadronicSystemGenerator` class,
82  // name it `GenerateEmpiricalDiNucleonCluster` or something...
83  this -> RecoilNucleonCluster(event);
84  // TODO: `DecayNucleonCluster` should probably be in `MECHadronicSystemGenerator`,
85  // if we make that...
86  this -> DecayNucleonCluster(event);
87  } else if (fXSecModel->Id().Name() == "genie::NievesSimoVacasMECPXSec2016") {
88  this -> SelectNSVLeptonKinematics(event);
89  this -> AddTargetRemnant(event);
90  this -> GenerateNSVInitialHadrons(event);
91  // Note: this method in `MECTensor/MECTensorGenerator.cxx` appeared to be a straight
92  // copy of an earlier version of the `DecayNucleonCluster` method here - but, watch
93  // for this...
94  this -> DecayNucleonCluster(event);
95  } else if (fXSecModel->Id().Name() == "genie::SuSAv2MECPXSec") {
96  event->Print();
97  this -> SelectSuSALeptonKinematics(event);
98  event->Print();
99  this -> AddTargetRemnant(event);
100  event->Print();
101  this -> GenerateNSVInitialHadrons(event);
102  event->Print();
103  // Note: this method in `MECTensor/MECTensorGenerator.cxx` appeared to be a straight
104  // copy of an earlier version of the `DecayNucleonCluster` method here - but, watch
105  // for this...
106  this -> DecayNucleonCluster(event);
107  }
108  else {
109  LOG("MECGenerator",pFATAL) <<
110  "ProcessEventRecord >> Cannot calculate kinematics for " <<
111  fXSecModel->Id().Name();
112  }
113 
114 
115 }
116 //___________________________________________________________________________
118 {
119 // Add the remnant nucleus (= initial nucleus - nucleon cluster) in the
120 // event record.
121 
122  GHepParticle * target = event->TargetNucleus();
123  GHepParticle * cluster = event->HitNucleon();
124 
125  int Z = target->Z();
126  int A = target->A();
127 
128  if(cluster->Pdg() == kPdgClusterNN) { A-=2; ; }
129  if(cluster->Pdg() == kPdgClusterNP) { A-=2; Z-=1; }
130  if(cluster->Pdg() == kPdgClusterPP) { A-=2; Z-=2; }
131 
132  int ipdgc = pdg::IonPdgCode(A, Z);
133 
134  const TLorentzVector & p4cluster = *(cluster->P4());
135  const TLorentzVector & p4tgt = *(target ->P4());
136 
137  const TLorentzVector p4 = p4tgt - p4cluster;
138  const TLorentzVector v4(0.,0.,0., 0.);
139 
140  int momidx = event->TargetNucleusPosition();
141 
142  /*
143  if( A == 2 && Z == 2){
144  // residual nucleus was just two protons, not a nucleus we know.
145  // this might not conserve energy...
146  event->AddParticle(kPdgProton,kIStStableFinalState, momidx,-1,-1,-1, p4,v4);
147  // v4,v4 because I'm lazy, give the four momentum to one of the protons, not the other
148  event->AddParticle(kPdgProton,kIStStableFinalState, momidx,-1,-1,-1, v4,v4);
149  } else if ( A == 2 && Z == 0){
150  // residual nucleus was just two neutrons, not a nucleus we know.
151  // this might not conserve energy...
152  event->AddParticle(kPdgNeutron,kIStStableFinalState, momidx,-1,-1,-1, p4,v4);
153  // v4,v4 because I'm lazy, give the four momentum to one of the protons, not the other
154  event->AddParticle(kPdgNeutron,kIStStableFinalState, momidx,-1,-1,-1, v4,v4);
155  } else {
156  // regular nucleus, including deuterium
157  event->AddParticle(ipdgc,kIStStableFinalState, momidx,-1,-1,-1, p4,v4);
158  }
159  */
160 
161  event->AddParticle(ipdgc,kIStStableFinalState, momidx,-1,-1,-1, p4,v4);
162 
163 }
164 //___________________________________________________________________________
166 {
167 // Generate the initial state di-nucleon cluster 4-momentum.
168 // Draw Fermi momenta for each of the two nucleons.
169 // Sum the two Fermi momenta to calculate the di-nucleon momentum.
170 // For simplicity, keep the di-nucleon cluster on the mass shell.
171 //
172  GHepParticle * target_nucleus = event->TargetNucleus();
173  assert(target_nucleus);
174  GHepParticle * nucleon_cluster = event->HitNucleon();
175  assert(nucleon_cluster);
176  GHepParticle * remnant_nucleus = event->RemnantNucleus();
177  assert(remnant_nucleus);
178 
179  // generate a Fermi momentum for each nucleon
180 
181  Target tgt(target_nucleus->Pdg());
182  PDGCodeList pdgv = this->NucleonClusterConstituents(nucleon_cluster->Pdg());
183  assert(pdgv.size()==2);
184  tgt.SetHitNucPdg(pdgv[0]);
186  TVector3 p3a = fNuclModel->Momentum3();
187  tgt.SetHitNucPdg(pdgv[1]);
189  TVector3 p3b = fNuclModel->Momentum3();
190 
191  LOG("FermiMover", pINFO)
192  << "1st nucleon (code = " << pdgv[0] << ") generated momentum: ("
193  << p3a.Px() << ", " << p3a.Py() << ", " << p3a.Pz() << "), "
194  << "|p| = " << p3a.Mag();
195  LOG("FermiMover", pINFO)
196  << "2nd nucleon (code = " << pdgv[1] << ") generated momentum: ("
197  << p3b.Px() << ", " << p3b.Py() << ", " << p3b.Pz() << "), "
198  << "|p| = " << p3b.Mag();
199 
200  // calcute nucleon cluster momentum
201 
202  TVector3 p3 = p3a + p3b;
203 
204  LOG("FermiMover", pINFO)
205  << "di-nucleon cluster momentum: ("
206  << p3.Px() << ", " << p3.Py() << ", " << p3.Pz() << "), "
207  << "|p| = " << p3.Mag();
208 
209  // target (initial) nucleus and nucleon cluster mass
210 
211  double Mi = PDGLibrary::Instance()->Find(target_nucleus->Pdg() )-> Mass(); // initial nucleus mass
212  double M2n = PDGLibrary::Instance()->Find(nucleon_cluster->Pdg())-> Mass(); // nucleon cluster mass
213 
214  // nucleon cluster energy
215 
216  double EN = TMath::Sqrt(p3.Mag2() + M2n*M2n);
217 
218  // set the remnant nucleus and nucleon cluster 4-momenta at GHEP record
219 
220  TLorentzVector p4nclust ( p3.Px(), p3.Py(), p3.Pz(), EN );
221  TLorentzVector p4remnant (-1*p3.Px(), -1*p3.Py(), -1*p3.Pz(), Mi-EN);
222 
223  nucleon_cluster->SetMomentum(p4nclust);
224  remnant_nucleus->SetMomentum(p4remnant);
225 
226  // set the nucleon cluster 4-momentum at the interaction summary
227 
228  event->Summary()->InitStatePtr()->TgtPtr()->SetHitNucP4(p4nclust);
229 }
230 //___________________________________________________________________________
232 {
233 // Select interaction kinematics using the rejection method.
234 //
235 
236  // Access cross section algorithm for running thread
238  const EventGeneratorI * evg = rtinfo->RunningThread();
239  fXSecModel = evg->CrossSectionAlg();
240 
241  Interaction * interaction = event->Summary();
242  double Ev = interaction->InitState().ProbeE(kRfHitNucRest);
243 
244  // **** NOTE / TODO:
245  // **** Hardcode bogus limits for the time-being
246  // **** Should be able to get limits via Interaction::KPhaseSpace
247  double Q2min = 0.01;
248  double Q2max = 8.00;
249  double Wmin = 1.88;
250  double Wmax = 3.00;
251 
252  // Scan phase-space for the maximum differential cross section
253  // at the current neutrino energy
254  const int nq=30;
255  const int nw=20;
256  double dQ2 = (Q2max-Q2min) / (nq-1);
257  double dW = (Wmax-Wmin ) / (nw-1);
258  double xsec_max = 0;
259  for(int iw=0; iw<nw; iw++) {
260  for(int iq=0; iq<nq; iq++) {
261  double Q2 = Q2min + iq*dQ2;
262  double W = Wmin + iw*dW;
263  interaction->KinePtr()->SetQ2(Q2);
264  interaction->KinePtr()->SetW (W);
265  double xsec = fXSecModel->XSec(interaction, kPSWQ2fE);
266  xsec_max = TMath::Max(xsec, xsec_max);
267  }
268  }
269  LOG("MEC", pNOTICE) << "xsec_max (E = " << Ev << " GeV) = " << xsec_max;
270 
271  // Select kinematics
272  RandomGen * rnd = RandomGen::Instance();
273  unsigned int iter = 0;
274  bool accept = false;
275  while(1) {
276  iter++;
277  if(iter > kRjMaxIterations) {
278  LOG("MEC", pWARN)
279  << "Couldn't select a valid W, Q^2 pair after "
280  << iter << " iterations";
281  event->EventFlags()->SetBitNumber(kKineGenErr, true);
283  exception.SetReason("Couldn't select kinematics");
284  exception.SwitchOnFastForward();
285  throw exception;
286  }
287 
288  // Generate next pair
289  double gQ2 = Q2min + (Q2max-Q2min) * rnd->RndKine().Rndm();
290  double gW = Wmin + (Wmax -Wmin ) * rnd->RndKine().Rndm();
291 
292  // Calculate d2sigma/dQ2dW
293  interaction->KinePtr()->SetQ2(gQ2);
294  interaction->KinePtr()->SetW (gW);
295  double xsec = fXSecModel->XSec(interaction, kPSWQ2fE);
296 
297  // Decide whether to accept the current kinematics
298  double t = xsec_max * rnd->RndKine().Rndm();
299  double J = 1; // jacobean
300  accept = (t < J*xsec);
301 
302  // If the generated kinematics are accepted, finish-up module's job
303  if(accept) {
304  LOG("MEC", pINFO) << "Selected: Q^2 = " << gQ2 << ", W = " << gW;
305  double gx = 0;
306  double gy = 0;
307  // More accurate calculation of the mass of the cluster than 2*Mnucl
308  int nucleon_cluster_pdg = interaction->InitState().Tgt().HitNucPdg();
309  double M2n = PDGLibrary::Instance()->Find(nucleon_cluster_pdg)->Mass();
310  //bool is_em = interaction->ProcInfo().IsEM();
311  kinematics::WQ2toXY(Ev,M2n,gW,gQ2,gx,gy);
312 
313  LOG("MEC", pINFO) << "x = " << gx << ", y = " << gy;
314  // lock selected kinematics & clear running values
315  interaction->KinePtr()->SetQ2(gQ2, true);
316  interaction->KinePtr()->SetW (gW, true);
317  interaction->KinePtr()->Setx (gx, true);
318  interaction->KinePtr()->Sety (gy, true);
319  interaction->KinePtr()->ClearRunningValues();
320 
321  return;
322  }//accepted?
323  }//iter
324 }
325 //___________________________________________________________________________
327 {
328 // Add the final-state primary lepton in the event record.
329 // Compute its 4-momentum based on the selected interaction kinematics.
330 //
331  Interaction * interaction = event->Summary();
332 
333  // The boost back to the lab frame was missing, that is included now with the introduction of the beta factor
334  const InitialState & init_state = interaction->InitState();
335  const TLorentzVector & pnuc4 = init_state.Tgt().HitNucP4(); //[@LAB]
336  TVector3 beta = pnuc4.BoostVector();
337 
338  // Boosting the incoming neutrino to the NN-cluster rest frame
339  // Neutrino 4p
340  TLorentzVector * p4v = event->Probe()->GetP4(); // v 4p @ LAB
341  p4v->Boost(-1.*beta); // v 4p @ NN-cluster rest frame
342 
343  // Look-up selected kinematics
344  double Q2 = interaction->Kine().Q2(true);
345  double y = interaction->Kine().y(true);
346 
347  // Auxiliary params
348  double Ev = interaction->InitState().ProbeE(kRfHitNucRest);
349  LOG("MEC", pNOTICE) << "neutrino energy = " << Ev;
350  double ml = interaction->FSPrimLepton()->Mass();
351  double ml2 = TMath::Power(ml,2);
352 
353  // Compute the final state primary lepton energy and momentum components
354  // along and perpendicular the neutrino direction
355  double El = (1-y)*Ev;
356  double plp = El - 0.5*(Q2+ml2)/Ev; // p(//)
357  double plt = TMath::Sqrt(TMath::Max(0.,El*El-plp*plp-ml2)); // p(-|)
358 
359  LOG("MEC", pNOTICE)
360  << "fsl: E = " << El << ", |p//| = " << plp << ", |pT| = " << plt;
361 
362  // Randomize transverse components
363  RandomGen * rnd = RandomGen::Instance();
364  double phi = 2*kPi * rnd->RndLep().Rndm();
365  double pltx = plt * TMath::Cos(phi);
366  double plty = plt * TMath::Sin(phi);
367 
368  // Take a unit vector along the neutrino direction
369  // WE NEED THE UNIT VECTOR ALONG THE NEUTRINO DIRECTION IN THE NN-CLUSTER REST FRAME, NOT IN THE LAB FRAME
370  TVector3 unit_nudir = p4v->Vect().Unit(); //We use this one, which is in the NN-cluster rest frame
371  // Rotate lepton momentum vector from the reference frame (x'y'z') where
372  // {z':(neutrino direction), z'x':(theta plane)} to the LAB
373  TVector3 p3l(pltx,plty,plp);
374  p3l.RotateUz(unit_nudir);
375 
376  // Lepton 4-momentum in LAB
377  TLorentzVector p4l(p3l,El);
378 
379  // Boost final state primary lepton to the lab frame
380  p4l.Boost(beta); // active Lorentz transform
381 
382  // Figure out the final-state primary lepton PDG code
383  int pdgc = interaction->FSPrimLepton()->PdgCode();
384 
385  // Lepton 4-position (= interacton vtx)
386  TLorentzVector v4(*event->Probe()->X4());
387 
388  // Add the final-state lepton to the event record
389  int momidx = event->ProbePosition();
390  event->AddParticle(
391  pdgc, kIStStableFinalState, momidx, -1, -1, -1, p4l, v4);
392 
393  // Set its polarization
395 }
396 //___________________________________________________________________________
398 {
399  // get di-nucleon cluster & its 4-momentum
400  GHepParticle * nucleon_cluster = event->HitNucleon();
401  assert(nucleon_cluster);
402  TLorentzVector * tmp=nucleon_cluster->GetP4();
403  TLorentzVector p4cluster(*tmp);
404  delete tmp;
405 
406  // get neutrino & its 4-momentum
407  GHepParticle * neutrino = event->Probe();
408  assert(neutrino);
409  TLorentzVector p4v(*neutrino->P4());
410 
411  // get final state primary lepton & its 4-momentum
412  GHepParticle * fsl = event->FinalStatePrimaryLepton();
413  assert(fsl);
414  TLorentzVector p4l(*fsl->P4());
415 
416  // calculate 4-momentum transfer
417  TLorentzVector q = p4v - p4l;
418 
419  // calculate recoil nucleon cluster 4-momentum
420  TLorentzVector p4cluster_recoil = p4cluster + q;
421 
422  // work-out recoil nucleon cluster code
423  LOG("MEC", pINFO) << "Interaction summary";
424  LOG("MEC", pINFO) << *event->Summary();
425  int recoil_nucleon_cluster_pdg = event->Summary()->RecoilNucleonPdg();
426 
427  // vtx position in nucleus coord system
428  TLorentzVector v4(*neutrino->X4());
429 
430  // add to the event record
431  event->AddParticle(
432  recoil_nucleon_cluster_pdg, kIStDecayedState,
433  2, -1, -1, -1, p4cluster_recoil, v4);
434 }
435 //___________________________________________________________________________
437 {
438 // Perform a phase-space decay of the nucleon cluster and add its decay
439 // products in the event record
440 //
441  LOG("MEC", pINFO) << "Decaying nucleon cluster...";
442 
443  // get di-nucleon cluster
444  int nucleon_cluster_id = 5;
445  GHepParticle * nucleon_cluster = event->Particle(nucleon_cluster_id);
446  assert(nucleon_cluster);
447 
448  // get decay products
449  PDGCodeList pdgv = this->NucleonClusterConstituents(nucleon_cluster->Pdg());
450  LOG("MEC", pINFO) << "Decay product IDs: " << pdgv;
451 
452  // Get the decay product masses
454  int i = 0;
455  double * mass = new double[pdgv.size()];
456  double sum = 0;
457  for(pdg_iter = pdgv.begin(); pdg_iter != pdgv.end(); ++pdg_iter) {
458  int pdgc = *pdg_iter;
459  double m = PDGLibrary::Instance()->Find(pdgc)->Mass();
460  mass[i++] = m;
461  sum += m;
462  }
463 
464  LOG("MEC", pINFO)
465  << "Performing a phase space decay to "
466  << pdgv.size() << " particles / total mass = " << sum;
467 
468  TLorentzVector * p4d = nucleon_cluster->GetP4();
469  TLorentzVector * v4d = nucleon_cluster->GetX4();
470 
471  LOG("MEC", pINFO)
472  << "Decaying system p4 = " << utils::print::P4AsString(p4d);
473 
474  // Set the decay
475  bool permitted = fPhaseSpaceGenerator.SetDecay(*p4d, pdgv.size(), mass);
476  if(!permitted) {
477  LOG("MEC", pERROR)
478  << " *** Phase space decay is not permitted \n"
479  << " Total particle mass = " << sum << "\n"
480  << " Decaying system p4 = " << utils::print::P4AsString(p4d);
481  // clean-up
482  delete [] mass;
483  delete p4d;
484  delete v4d;
485  // throw exception
486  event->EventFlags()->SetBitNumber(kHadroSysGenErr, true);
488  exception.SetReason("Decay not permitted kinematically");
489  exception.SwitchOnStepBack();
490  exception.SetReturnStep(0);
491  throw exception;
492  }
493 
494  // Get the maximum weight
495  double wmax = -1;
496  for(int idec=0; idec<200; idec++) {
497  double w = fPhaseSpaceGenerator.Generate();
498  wmax = TMath::Max(wmax,w);
499  }
500  assert(wmax>0);
501  wmax *= 2;
502 
503  LOG("MEC", pNOTICE)
504  << "Max phase space gen. weight = " << wmax;
505 
506  RandomGen * rnd = RandomGen::Instance();
507  bool accept_decay=false;
508  unsigned int itry=0;
509  while(!accept_decay)
510  {
511  itry++;
512 
514  // report, clean-up and return
515  LOG("MEC", pWARN)
516  << "Couldn't generate an unweighted phase space decay after "
517  << itry << " attempts";
518  // clean up
519  delete [] mass;
520  delete p4d;
521  delete v4d;
522  // throw exception
523  event->EventFlags()->SetBitNumber(kHadroSysGenErr, true);
525  exception.SetReason("Couldn't select decay after N attempts");
526  exception.SwitchOnStepBack();
527  exception.SetReturnStep(0);
528  throw exception;
529  }
530  double w = fPhaseSpaceGenerator.Generate();
531  if(w > wmax) {
532  LOG("MEC", pWARN)
533  << "Decay weight = " << w << " > max decay weight = " << wmax;
534  }
535  double gw = wmax * rnd->RndDec().Rndm();
536  accept_decay = (gw<=w);
537 
538  LOG("MEC", pINFO)
539  << "Decay weight = " << w << " / R = " << gw
540  << " - accepted: " << accept_decay;
541 
542  } //!accept_decay
543 
544  // Insert the decay products in the event record
545  TLorentzVector v4(*v4d);
547  int idp = 0;
548  for(pdg_iter = pdgv.begin(); pdg_iter != pdgv.end(); ++pdg_iter) {
549  int pdgc = *pdg_iter;
550  TLorentzVector * p4fin = fPhaseSpaceGenerator.GetDecay(idp);
551  event->AddParticle(pdgc, ist, nucleon_cluster_id,-1,-1,-1, *p4fin, v4);
552  idp++;
553  }
554 
555  // Clean-up
556  delete [] mass;
557  delete p4d;
558  delete v4d;
559 }
560 //___________________________________________________________________________
562 {
563  bool allowdup = true;
564  PDGCodeList pdgv(allowdup);
565 
566  if(pdgc == kPdgClusterNN) {
567  pdgv.push_back(kPdgNeutron);
568  pdgv.push_back(kPdgNeutron);
569  }
570  else
571  if(pdgc == kPdgClusterNP) {
572  pdgv.push_back(kPdgNeutron);
573  pdgv.push_back(kPdgProton);
574  }
575  else
576  if(pdgc == kPdgClusterPP) {
577  pdgv.push_back(kPdgProton);
578  pdgv.push_back(kPdgProton);
579  }
580  else
581  {
582  LOG("MEC", pERROR)
583  << "Unknown di-nucleon cluster PDG code (" << pdgc << ")";
584  }
585 
586  return pdgv;
587 }
588 //___________________________________________________________________________
590 {
591  // -- implementation -- //
592  // The IFIC Valencia model can provide three different hadron tensors.
593  // The user probably wants all CC QE-like 2p2h events
594  // But could in principle get the no-delta component if they want (deactivated incode)
595  int FullDeltaNodelta = 1; // 1: full, 2: only delta, 3: zero delta
596 
597  // -- Event Properties -----------------------------//
598  Interaction * interaction = event->Summary();
599  Kinematics * kinematics = interaction->KinePtr();
600 
601  double Enu = interaction->InitState().ProbeE(kRfHitNucRest);
602 
603  int NuPDG = interaction->InitState().ProbePdg();
604  int TgtPDG = interaction->InitState().TgtPdg();
605  // interacton vtx
606  TLorentzVector v4(*event->Probe()->X4());
607  TLorentzVector tempp4(0.,0.,0.,0.);
608 
609  // -- Lepton Kinematic Limits ----------------------------------------- //
610  double Costh = 0.0; // lepton angle
611  double CosthMax = 1.0;
612  double CosthMin = -1.0;
613 
614  double T = 0.0; // lepton kinetic energy
615  double TMax = std::numeric_limits<double>::max();
616  double TMin = 0.0;
617 
618  double Plep = 0.0; // lepton 3 momentum
619  double Elep = 0.0; // lepton energy
620  double LepMass = interaction->FSPrimLepton()->Mass();
621 
622  double Q0 = 0.0; // energy component of q four vector
623  double Q3 = 0.0; // magnitude of transfered 3 momentum
624  double Q2 = 0.0; // properly Q^2 (Q squared) - transfered 4 momentum.
625 
626  // Set lepton KE TMax for for throwing rndm in the accept/reject loop.
627  // We can accidentally set it too high, because the xsec will return zero.
628  // This way if someone reuses this code, they are not tripped up by it.
629  TMax = Enu - LepMass;
630 
631  // Set Tmin for throwing rndm in the accept/reject loop
632  // the hadron tensors we expect will be limited in q3
633  // therefore also the outgoing lepton KE can't be too low or costheta too backward
634  // make the accept/reject loop more efficient by using Min values.
635  if(Enu < fQ3Max){
636  TMin = 0 ;
637  CosthMin = -1 ;
638  } else {
639  TMin = TMath::Sqrt(TMath::Power(LepMass, 2) + TMath::Power((Enu - fQ3Max), 2)) - LepMass;
640  CosthMin = TMath::Sqrt(1 - TMath::Power((fQ3Max / Enu ), 2));
641  }
642 
643  // Compute the maximum xsec value:
644  Range1D_t Tl_range ( TMin, TMax ) ;
645  Range1D_t ctl_range ( CosthMin, CosthMax ) ;
646  double XSecMax = GetXSecMaxTlctl( *interaction, Tl_range, ctl_range ) ;
647 
648  // -- Generate and Test the Kinematics----------------------------------//
649 
650  RandomGen * rnd = RandomGen::Instance();
651  bool accept = false;
652  unsigned int iter = 0;
653 
654  // loop over different (randomly) selected T and Costh
655  while (!accept) {
656  iter++;
657  if(iter > kRjMaxIterations) {
658  // error if try too many times
659  LOG("MEC", pWARN)
660  << "Couldn't select a valid Tmu, CosTheta pair after "
661  << iter << " iterations";
662  event->EventFlags()->SetBitNumber(kKineGenErr, true);
664  exception.SetReason("Couldn't select lepton kinematics");
665  exception.SwitchOnFastForward();
666  throw exception;
667  }
668 
669  // generate random kinetic energy T and Costh
670  T = TMin + (TMax-TMin)*rnd->RndKine().Rndm();
671  Costh = CosthMin + (CosthMax-CosthMin)*rnd->RndKine().Rndm();
672 
673  // Calculate useful values for judging this choice
674  genie::utils::mec::Getq0q3FromTlCostl(T, Costh, Enu, LepMass, Q0, Q3);
675 
676  // Don't bother doing hard work if the selected Q3 is greater than Q3Max
677  if (Q3 > fQ3Max) continue ;
678 
679  Plep = TMath::Sqrt( T * (T + (2.0 * LepMass))); // ok is sqrt(E2 - m2)
680  kinematics->SetKV(kKVTl, T);
681  kinematics->SetKV(kKVctl, Costh);
682  kinematics->SetKV( kKVQ0, Q0 ) ;
683  kinematics->SetKV( kKVQ3, Q3 ) ;
684 
685  // decide whether to accept or reject these kinematics
686  // AND set the chosen two-nucleon system
687 
688  if (FullDeltaNodelta == 1){
689  // this block for the user who wants all CC QE-like 2p2h events
690  // We need four different cross sections. Right now, pursue the
691  // inelegant method of calling XSec four times - there is
692  // definitely some runtime inefficiency here, but it is not awful
693 
694  // first, get delta-less all
695  if (NuPDG > 0) {
696  interaction->InitStatePtr()->TgtPtr()->SetHitNucPdg(kPdgClusterNN);
697  }
698  else {
699  interaction->InitStatePtr()->TgtPtr()->SetHitNucPdg(kPdgClusterPP);
700  }
701  double XSec = fXSecModel->XSec(interaction, kPSTlctl);
702  // now get all with delta
703  interaction->ExclTagPtr()->SetResonance(genie::kP33_1232);
704  double XSecDelta = fXSecModel->XSec(interaction, kPSTlctl);
705  // get PN with delta
706  interaction->InitStatePtr()->TgtPtr()->SetHitNucPdg(kPdgClusterNP);
707  double XSecDeltaPN = fXSecModel->XSec(interaction, kPSTlctl);
708  // now get delta-less PN
710  double XSecPN = fXSecModel->XSec(interaction, kPSTlctl);
711 
712  if (XSec > XSecMax) {
713  LOG("MEC", pERROR) << "XSec is > XSecMax for nucleus " << TgtPDG << " "
714  << XSec << " > " << XSecMax
715  << " don't let this happen.";
716  }
717  assert(XSec <= XSecMax);
718  accept = XSec > XSecMax*rnd->RndKine().Rndm();
719  LOG("MEC", pINFO) << "Xsec, Max, Accept: " << XSec << ", "
720  << XSecMax << ", " << accept;
721 
722  if(accept){
723  // If it passes the All cross section we still need to do two things:
724  // * Was the initial state pn or not?
725  // * Do we assign the reaction to have had a Delta on the inside?
726 
727  // PDD means from the part of the XSec with an internal Delta line
728  // that (at the diagram level) did not produce a pion in the final state.
729 
730  bool isPDD = false;
731 
732  // Find out if we should use a pn initial state
733  double myrand = rnd->RndKine().Rndm();
734  double pnFraction = XSecPN / XSec;
735  LOG("MEC", pDEBUG) << "Test for pn: xsec_pn = " << XSecPN
736  << "; xsec = " << XSec
737  << "; pn_fraction = " << pnFraction
738  << "; random number val = " << myrand;
739 
740  if (myrand <= pnFraction) {
741  // yes it is, add a PN initial state to event record
742  event->AddParticle(kPdgClusterNP, kIStNucleonTarget,
743  1, -1, -1, -1, tempp4, v4);
744  interaction->InitStatePtr()->TgtPtr()->SetHitNucPdg(kPdgClusterNP);
745 
746  // Its a pn, so test for Delta by comparing DeltaPN/PN
747  if (rnd->RndKine().Rndm() <= XSecDeltaPN / XSecPN) {
748  isPDD = true;
749  }
750  }
751  else {
752  // no it is not a PN, add either NN or PP initial state to event record.
753  if (NuPDG > 0) {
754  event->AddParticle(kPdgClusterNN, kIStNucleonTarget,
755  1, -1, -1, -1, tempp4, v4);
756  interaction->InitStatePtr()->TgtPtr()->SetHitNucPdg(kPdgClusterNN);
757  }
758  else {
759  event->AddParticle(kPdgClusterPP, kIStNucleonTarget,
760  1, -1, -1, -1, tempp4, v4);
761  interaction->InitStatePtr()->TgtPtr()->SetHitNucPdg(kPdgClusterPP);
762  }
763  // its not pn, so test for Delta (XSecDelta-XSecDeltaPN)/(XSec-XSecPN)
764  // right, both numerator and denominator are total not pn.
765  if (rnd->RndKine().Rndm() <=
766  (XSecDelta - XSecDeltaPN) / (XSec - XSecPN)) {
767  isPDD = true;
768  }
769  }
770 
771  // now test whether we tagged this as a pion event
772  // and assign that fact to the Exclusive State tag
773  // later, we can query const XclsTag & xcls = interaction->ExclTag()
774  if (isPDD){
775  interaction->ExclTagPtr()->SetResonance(genie::kP33_1232);
776  }
777 
778 
779  } // end if accept
780  } // end if delta == 1
781 
782  /* One can make simpler versions of the above for the
783  FullDeltaNodelta == 2 (only delta)
784  or
785  FullDeltaNodelta == 3 (set Delta FF = 1, lose interference effect).
786  but I (Rik) don't see what the use-case is for these, genratorly speaking.
787  */
788 
789  } // end while
790 
791  // -- finish lepton kinematics
792  // If the code got here, then we accepted some kinematics
793  // and we can proceed to generate the final state.
794 
795  // define coordinate system wrt neutrino: z along neutrino, xy perp
796 
797  // Cos theta gives us z, the rest in xy:
798  double PlepZ = Plep * Costh;
799  double PlepXY = Plep * TMath::Sqrt(1. - TMath::Power(Costh,2));
800 
801  // random rotation about unit vector for phi direction
802  double phi= 2 * kPi * rnd->RndLep().Rndm();
803  // now fill x and y from PlepXY
804  double PlepX = PlepXY * TMath::Cos(phi);
805  double PlepY = PlepXY * TMath::Sin(phi);
806 
807  // Rotate lepton momentum vector from the reference frame (x'y'z') where
808  // {z':(neutrino direction), z'x':(theta plane)} to the LAB
809  TVector3 unit_nudir = event->Probe()->P4()->Vect().Unit();
810  TVector3 p3l(PlepX, PlepY, PlepZ);
811  p3l.RotateUz(unit_nudir);
812 
813  // Lepton 4-momentum in LAB
814  Elep = TMath::Sqrt(LepMass*LepMass + PlepX*PlepX + PlepY*PlepY + PlepZ*PlepZ);
815  TLorentzVector p4l(p3l,Elep);
816 
817  // Figure out the final-state primary lepton PDG code
818  int pdgc = interaction->FSPrimLepton()->PdgCode();
819  int momidx = event->ProbePosition();
820 
821  // -- Store Values ------------------------------------------//
822  // -- Interaction: Q2
823  Q2 = Q3*Q3 - Q0*Q0;
824  double gy = Q0 / Enu;
825  double gx = kinematics::Q2YtoX(Enu, 2 * kNucleonMass, Q2, gy);
826  double gW = kinematics::XYtoW(Enu, 2 * kNucleonMass, gx, gy);
827 
828  interaction->KinePtr()->SetQ2(Q2, true);
829  interaction->KinePtr()->Sety(gy, true);
830  interaction->KinePtr()->Setx(gx, true);
831  interaction->KinePtr()->SetW(gW, true);
832  interaction->KinePtr()->SetFSLeptonP4(p4l);
833  // in later methods
834  // will also set the four-momentum and W^2 of the hadron system.
835 
836  // -- Lepton
837  event->AddParticle( pdgc, kIStStableFinalState, momidx, -1, -1, -1, p4l, v4);
838 
839  // Set the final-state lepton polarization
841 
842  LOG("MEC",pDEBUG) << "~~~ LEPTON DONE ~~~";
843 }
844 //___________________________________________________________________________
846 {
847  // Event Properties
848  Interaction* interaction = event->Summary();
849  Kinematics* kinematics = interaction->KinePtr();
850 
851  // Choose the appropriate minimum Q^2 value based on the interaction
852  // mode (this is important for EM interactions since the differential
853  // cross section blows up as Q^2 --> 0)
854  double Q2min = genie::controls::kMinQ2Limit; // CC/NC limit
855  if ( interaction->ProcInfo().IsEM() ) Q2min = genie::utils::kinematics
856  ::electromagnetic::kMinQ2Limit; // EM limit
857 
858  LOG("MEC", pDEBUG) << "Q2min = " << Q2min;
859 
860  double Enu = interaction->InitState().ProbeE( kRfLab );
861 
862  int NuPDG = interaction->InitState().ProbePdg();
863  int TgtPDG = interaction->InitState().TgtPdg();
864 
865  // Interacton vtx
866  TLorentzVector v4( *event->Probe()->X4() );
867  TLorentzVector tempp4( 0., 0., 0., 0. );
868 
869  // Lepton Kinematic Limits
870  double Costh = 0.0; // lepton angle
871  double CosthMax = 1.0;
872  double CosthMin = -1.0;
873 
874  double T = 0.0; // lepton kinetic energy
875  double TMax = std::numeric_limits<double>::max();
876  double TMin = 0.0;
877 
878  double Plep = 0.0; // lepton 3 momentum
879  double Elep = 0.0; // lepton energy
880  double LepMass = interaction->FSPrimLepton()->Mass();
881 
882  double Q0 = 0.0; // energy component of q four vector
883  double Q3 = 0.0; // magnitude of transfered 3 momentum
884  double Q2 = 0.0; // properly Q^2 (Q squared) - transfered 4 momentum.
885 
886  // Set lepton KE TMax for for throwing rndm in the accept/reject loop.
887  // We can accidentally set it too high, because the xsec will return zero.
888  // This way if someone reuses this code, they are not tripped up by it.
889  TMax = Enu - LepMass;
890 
891  // TODO: double-check the limits below
892 
893  // Set Tmin for throwing rndm in the accept/reject loop
894  // the hadron tensors we expect will be limited in q3
895  // therefore also the outgoing lepton KE can't be too low or costheta too backward
896  // make the accept/reject loop more efficient by using Min values.
897  if ( Enu < fQ3Max ) {
898  TMin = 0;
899  CosthMin = -1;
900  } else {
901  TMin = TMath::Sqrt( TMath::Power(LepMass, 2) + TMath::Power(Enu - fQ3Max, 2) ) - LepMass;
902  CosthMin = TMath::Sqrt( 1. - TMath::Power(fQ3Max / Enu, 2) );
903  }
904 
905  // Generate and Test the Kinematics
906 
908  bool accept = false;
909  unsigned int iter = 0;
910  unsigned int maxIter = kRjMaxIterations;
911 
912  // TODO: revisit this
913  // e-scat xsecs blow up close to theta=0, MC methods won't work ...
914  if ( NuPDG == 11 ) maxIter *= 100000;
915 
916  // Scan the accessible phase space to find the maximum differential cross
917  // section to throw against
918  double XSecMax = utils::mec::GetMaxXSecTlctl( *fXSecModel, *interaction );
919 
920  // loop over different (randomly) selected T and Costh
921  while ( !accept ) {
922  ++iter;
923  if ( iter > maxIter ) {
924  // error if try too many times
925  LOG("MEC", pWARN)
926  << "Couldn't select a valid Tmu, CosTheta pair after "
927  << iter << " iterations";
928  event->EventFlags()->SetBitNumber( kKineGenErr, true );
930  exception.SetReason( "Couldn't select lepton kinematics" );
931  exception.SwitchOnFastForward();
932  throw exception;
933  }
934 
935  // generate random kinetic energy T and Costh
936  T = TMin + (TMax-TMin)*rnd->RndKine().Rndm();
937  Costh = CosthMin + (CosthMax-CosthMin)*rnd->RndKine().Rndm();
938 
939  // Calculate useful values for judging this choice
940  Plep = TMath::Sqrt( T * (T + (2.0 * LepMass))); // ok is sqrt(E2 - m2)
941  Q3 = TMath::Sqrt(Plep*Plep + Enu*Enu - 2.0 * Plep * Enu * Costh);
942 
943  // TODO: implement this more cleanly (throw Costh from restricted range)
944  Q0 = Enu - (T + LepMass);
945  Q2 = Q3*Q3 - Q0*Q0;
946 
947  LOG("MEC", pDEBUG) << "T = " << T << ", Costh = " << Costh
948  << ", Q2 = " << Q2;
949 
950  // Don't bother doing hard work if the selected Q3 is greater than Q3Max
951  // or if Q2 falls below the minimum allowed Q^2 value
952  if ( Q3 < fQ3Max && Q2 >= Q2min ) {
953 
954  kinematics->SetKV(kKVTl, T);
955  kinematics->SetKV(kKVctl, Costh);
956 
957  // decide whether to accept or reject these kinematics
958  // AND set the chosen two-nucleon system
959 
960  LOG("MEC", pDEBUG) << " T, Costh: " << T << ", " << Costh ;
961 
962  // Get total xsec (nn+np)
963  double XSec = fXSecModel->XSec( interaction, kPSTlctl );
964 
965  if ( XSec > XSecMax ) {
966  LOG("MEC", pERROR) << "XSec is > XSecMax for nucleus " << TgtPDG << " "
967  << XSec << " > " << XSecMax << " don't let this happen.";
968 
969  double percent_deviation = 200. * ( XSec - XSecMax ) / ( XSecMax + XSec );
970 
971  if ( percent_deviation > fSuSAMaxXSecDiffTolerance ) {
972  LOG( "Kinematics", pFATAL ) << "xsec: (curr) = " << XSec
973  << " > (max) = " << XSecMax << "\n for " << *interaction;
974  LOG( "Kinematics", pFATAL )
975  << "*** Exceeding estimated maximum differential cross section";
976  std::terminate();
977  }
978  else {
979  LOG( "Kinematics", pWARN ) << "xsec: (curr) = " << XSec
980  << " > (max) = " << XSecMax << "\n for " << *interaction;
981  LOG("Kinematics", pWARN) << "*** The fractional deviation of "
982  << percent_deviation << " % was allowed";
983  }
984  }
985 
986  accept = XSec > XSecMax*rnd->RndKine().Rndm();
987  LOG("MEC", pINFO) << "Xsec, Max, Accept: " << XSec << ", "
988  << XSecMax << ", " << accept;
989 
990  if ( accept ) {
991  // Now that we've selected kinematics, we also need to choose the
992  // isospin of the initial hit nucleon pair
993 
994  // Find out if we should use a pn initial state
995  double myrand = rnd->RndKine().Rndm();
996  double pnFraction = dynamic_cast< const SuSAv2MECPXSec* >( fXSecModel )
997  ->PairRatio( interaction );
998 
999  LOG("MEC", pINFO) << "Test for pn: "
1000  << "; xsec = " << XSec << "; pn_fraction = " << pnFraction
1001  << "; random number val = " << myrand;
1002 
1003  if ( myrand <= pnFraction ) {
1004  // yes it is, add a PN initial state to event record
1005  event->AddParticle(kPdgClusterNP, kIStNucleonTarget,
1006  1, -1, -1, -1, tempp4, v4);
1007  interaction->InitStatePtr()->TgtPtr()->SetHitNucPdg( kPdgClusterNP );
1008  }
1009  else {
1010  // no it is not a PN, add either NN or PP initial state to event record.
1011  if ( NuPDG > 0 ) {
1012  event->AddParticle(kPdgClusterNN, kIStNucleonTarget,
1013  1, -1, -1, -1, tempp4, v4);
1014  interaction->InitStatePtr()->TgtPtr()->SetHitNucPdg( kPdgClusterNN );
1015  }
1016  else {
1017  event->AddParticle(kPdgClusterPP, kIStNucleonTarget,
1018  1, -1, -1, -1, tempp4, v4);
1019  interaction->InitStatePtr()->TgtPtr()->SetHitNucPdg( kPdgClusterPP );
1020  }
1021  }
1022  } // end if accept
1023  } // end if passes q3 test
1024  } // end while
1025 
1026  // -- finish lepton kinematics
1027  // If the code got here, then we accepted some kinematics
1028  // and we can proceed to generate the final state.
1029 
1030  // define coordinate system wrt neutrino: z along neutrino, xy perp
1031 
1032  // Cos theta gives us z, the rest in xy:
1033  double PlepZ = Plep * Costh;
1034  double PlepXY = Plep * TMath::Sqrt( 1. - TMath::Power(Costh,2) );
1035 
1036  // random rotation about unit vector for phi direction
1037  double phi = 2. * kPi * rnd->RndLep().Rndm();
1038  // now fill x and y from PlepXY
1039  double PlepX = PlepXY * TMath::Cos(phi);
1040  double PlepY = PlepXY * TMath::Sin(phi);
1041 
1042  // Rotate lepton momentum vector from the reference frame (x'y'z') where
1043  // {z':(neutrino direction), z'x':(theta plane)} to the LAB
1044  TVector3 unit_nudir = event->Probe()->P4()->Vect().Unit();
1045  TVector3 p3l( PlepX, PlepY, PlepZ );
1046  p3l.RotateUz( unit_nudir );
1047 
1048  // Lepton 4-momentum in LAB
1049  Elep = TMath::Sqrt( LepMass*LepMass + PlepX*PlepX + PlepY*PlepY + PlepZ*PlepZ );
1050  TLorentzVector p4l( p3l, Elep );
1051 
1052  // Figure out the final-state primary lepton PDG code
1053  int pdgc = interaction->FSPrimLepton()->PdgCode();
1054  int momidx = event->ProbePosition();
1055 
1056  // -- Store Values ------------------------------------------//
1057  // -- Interaction: Q2
1058  Q0 = Enu - Elep;
1059  Q2 = Q3*Q3 - Q0*Q0;
1060  double gy = Q0 / Enu;
1061  double gx = kinematics::Q2YtoX(Enu, 2 * kNucleonMass, Q2, gy);
1062  double gW = kinematics::XYtoW(Enu, 2 * kNucleonMass, gx, gy);
1063 
1064  interaction->KinePtr()->SetQ2(Q2, true);
1065  interaction->KinePtr()->Sety(gy, true);
1066  interaction->KinePtr()->Setx(gx, true);
1067  interaction->KinePtr()->SetW(gW, true);
1068  interaction->KinePtr()->SetFSLeptonP4(p4l);
1069  // in later methods
1070  // will also set the four-momentum and W^2 of the hadron system.
1071 
1072  // -- Lepton
1073  event->AddParticle( pdgc, kIStStableFinalState, momidx, -1, -1, -1, p4l, v4 );
1074 
1075  LOG("MEC", pDEBUG) << "~~~ LEPTON DONE ~~~";
1076 }
1077 //___________________________________________________________________________
1079 {
1080  // We need a kinematic limits accept/reject loop here, so generating the
1081  // initial hadrons is combined with generating the recoil hadrons...
1082 
1083  LOG("MEC",pDEBUG) << "Generate Initial Hadrons - Start";
1084 
1085  Interaction * interaction = event->Summary();
1086  GHepParticle * neutrino = event->Probe();
1087  assert(neutrino);
1088  TLorentzVector p4nu(*neutrino->P4());
1089 
1090  // get final state primary lepton & its 4-momentum
1091  GHepParticle * fsl = event->FinalStatePrimaryLepton();
1092  assert(fsl);
1093  TLorentzVector p4l(*fsl->P4());
1094 
1095  // calculate 4-momentum transfer from these
1096  TLorentzVector Q4 = p4nu - p4l;
1097 
1098  // get the target two-nucleon cluster and nucleus.
1099  // the remnant nucleus is apparently set, except for its momentum.
1100  GHepParticle * target_nucleus = event->TargetNucleus();
1101  assert(target_nucleus);
1102  GHepParticle * initial_nucleon_cluster = event->HitNucleon();
1103  assert(initial_nucleon_cluster);
1104  GHepParticle * remnant_nucleus = event->RemnantNucleus();
1105  assert(remnant_nucleus);
1106 
1107  // -- make a two-nucleon system, then give them some momenta.
1108 
1109  // instantiate an empty local target nucleus, so I can use existing methods
1110  // to get a momentum from the prevailing Fermi-motion distribution
1111  Target tgt(target_nucleus->Pdg());
1112 
1113  // NucleonClusterConstituents is an implementation within this class, called with this
1114  // It using the nucleon cluster from the earlier tests for a pn state,
1115  // the method returns a vector of pdgs, which hopefully will be of size two.
1116 
1117  PDGCodeList pdgv = this->NucleonClusterConstituents(initial_nucleon_cluster->Pdg());
1118  assert(pdgv.size()==2);
1119 
1120  // These things need to be saved through to the end of the accept loop.
1121  bool accept = false;
1122  TVector3 p31i;
1123  TVector3 p32i;
1124  unsigned int iter = 0;
1125 
1126  int initial_nucleon_cluster_pdg = initial_nucleon_cluster->Pdg();
1127  int final_nucleon_cluster_pdg = 0;
1128 
1129  //e-scat case:
1130  if (neutrino->Pdg() == 11) {
1131  final_nucleon_cluster_pdg = initial_nucleon_cluster->Pdg();
1132  }
1133  //neutrino case
1134  else if (neutrino->Pdg() > 0) {
1135  if (initial_nucleon_cluster->Pdg() == kPdgClusterNP) {
1136  final_nucleon_cluster_pdg = kPdgClusterPP;
1137  }
1138  else if (initial_nucleon_cluster->Pdg() == kPdgClusterNN) {
1139  final_nucleon_cluster_pdg = kPdgClusterNP;
1140  }
1141  else {
1142  LOG("MEC", pERROR) << "Wrong pdg for a CC neutrino MEC interaction"
1143  << initial_nucleon_cluster->Pdg();
1144  }
1145  }
1146  else if (neutrino->Pdg() < 0) {
1147  if (initial_nucleon_cluster->Pdg() == kPdgClusterNP) {
1148  final_nucleon_cluster_pdg = kPdgClusterNN;
1149  }
1150  else if (initial_nucleon_cluster->Pdg() == kPdgClusterPP) {
1151  final_nucleon_cluster_pdg = kPdgClusterNP;
1152  }
1153  else {
1154  LOG("MEC", pERROR) << "Wrong pdg for a CC anti-neutrino MEC interaction"
1155  << initial_nucleon_cluster->Pdg();
1156  }
1157  }
1158 
1159  TLorentzVector p4initial_cluster;
1160  TLorentzVector p4final_cluster;
1161  TLorentzVector p4remnant_nucleus;
1162  double removalenergy1;
1163  double removalenergy2;
1164 
1165  //===========================================================================
1166  // Choose two nucleons from the prevailing fermi-motion distribution.
1167  // Some produce kinematically unallowed combinations initial cluster and Q2
1168  // Find out, and if so choose them again with this accept/reject loop.
1169  // Some kinematics are especially tough
1170  while(!accept){
1171  iter++;
1172  if(iter > kRjMaxIterations*1000) {
1173  // error if try too many times
1174  LOG("MEC", pWARN)
1175  << "Couldn't select a valid W, Q^2 pair after "
1176  << iter << " iterations";
1177  event->EventFlags()->SetBitNumber(kKineGenErr, true);
1179  exception.SetReason("Couldn't select initial hadron kinematics");
1180  exception.SwitchOnFastForward();
1181  throw exception;
1182  }
1183 
1184  // generate nucleons
1185  // tgt is a Target object for local use, just waiting to be filled.
1186  // this sets the pdg of each nucleon and its momentum from a Fermi-gas
1187  // Nieves et al. would use a local Fermi gas here, not this, but ok.
1188  // so momentum from global Fermi gas, local Fermi gas, or spectral function
1189  // and removal energy ~0.025 GeV, correlated with density, or from SF distribution
1190  tgt.SetHitNucPdg(pdgv[0]);
1192  p31i = fNuclModel->Momentum3();
1193  removalenergy1 = fNuclModel->RemovalEnergy();
1194  tgt.SetHitNucPdg(pdgv[1]);
1196  p32i = fNuclModel->Momentum3();
1197  removalenergy2 = fNuclModel->RemovalEnergy();
1198 
1199  // not sure -- could give option to use Nieves q-value here.
1200 
1201  // Now write down the initial cluster four-vector for this choice
1202  TVector3 p3i = p31i + p32i;
1203  double mass2 = PDGLibrary::Instance()->Find( initial_nucleon_cluster_pdg )->Mass();
1204  mass2 *= mass2;
1205  double energy = TMath::Sqrt(p3i.Mag2() + mass2);
1206  p4initial_cluster.SetPxPyPzE(p3i.Px(),p3i.Py(),p3i.Pz(),energy);
1207 
1208  // cast the removal energy as the energy component of a 4-vector for later.
1209  TLorentzVector tLVebind(0., 0., 0., -1.0 * (removalenergy1 + removalenergy2));
1210 
1211  // RIK: You might ask why is this the right place to subtract ebind?
1212  // It is okay. Physically, I'm subtracting it from q. The energy
1213  // transfer to the nucleon is 50 MeV less. the energy cost to put this
1214  // cluster on-shell. What Jan says he does in PRC.86.015504 is this:
1215  // The nucleons are assumed to be in a potential well
1216  // V = Efermi + 8 MeV.
1217  // The Fermi energy is subtracted from each initial-state nucleon
1218  // (I guess he does it dynamically based on Ef = p2/2M or so which
1219  // is what we are doing above, on average). Then after the reaction,
1220  // another 8 MeV is subtracted at that point; a small adjustment to the
1221  // momentum is needed to keep the nucleon on shell.
1222 
1223  // calculate recoil nucleon cluster 4-momentum (tLVebind is negative)
1224  p4final_cluster = p4initial_cluster + Q4 + tLVebind;
1225 
1226  // Test if the resulting four-vector corresponds to a high-enough invariant mass.
1227  // Fail the accept if we couldn't put this thing on-shell.
1228  if (p4final_cluster.M() <
1229  PDGLibrary::Instance()->Find(final_nucleon_cluster_pdg )->Mass()) {
1230  accept = false;
1231  } else {
1232  accept = true;
1233  }
1234 
1235  } // end accept loop
1236 
1237  // we got here if we accepted the final state two-nucleon system
1238  // so now we need to write everything to ghep
1239 
1240  // First the initial state nucleons.
1241  initial_nucleon_cluster->SetMomentum(p4initial_cluster);
1242 
1243  // and the remnant nucleus
1244  double Mi = PDGLibrary::Instance()->Find(target_nucleus->Pdg() )-> Mass();
1245  remnant_nucleus->SetMomentum(-1.0*p4initial_cluster.Px(),
1246  -1.0*p4initial_cluster.Py(),
1247  -1.0*p4initial_cluster.Pz(),
1248  Mi - p4initial_cluster.E() + removalenergy1 + removalenergy2);
1249 
1250  // Now the final nucleon cluster.
1251 
1252  // Getting this v4 is a position, i.e. a position within the nucleus (?)
1253  // possibly it takes on a non-trivial value only for local Fermi gas
1254  // or for sophisticated treatments of intranuclear rescattering.
1255  TLorentzVector v4(*neutrino->X4());
1256 
1257  // Now write the (undecayed) final two-nucleon system
1258  GHepParticle p1(final_nucleon_cluster_pdg, kIStDecayedState,
1259  2, -1, -1, -1, p4final_cluster, v4);
1260 
1261  //p1.SetRemovalEnergy(removalenergy1 + removalenergy2);
1262  // The "bound particle" concept applies only to p or n.
1263  // Instead, add this directly to the remnant nucleon a few lines above.
1264 
1265  // actually, this is not an status1 particle, so it is not picked up
1266  // by the aggregator. anyway, the aggregator does not run until the very end.
1267 
1268  event->AddParticle(p1);
1269 
1270  interaction->KinePtr()->SetHadSystP4(p4final_cluster);
1271 }
1272 //___________________________________________________________________________
1274 {
1275  Algorithm::Configure(config);
1276  this->LoadConfig();
1277 }
1278 //___________________________________________________________________________
1280 {
1281  Algorithm::Configure(config);
1282  this->LoadConfig();
1283 }
1284 //___________________________________________________________________________
1286 {
1287  fNuclModel = 0;
1288  RgKey nuclkey = "NuclearModel";
1289  fNuclModel = dynamic_cast<const NuclearModelI *> (this->SubAlg(nuclkey));
1290  assert(fNuclModel);
1291 
1292  GetParamDef( "MaxXSec-SafetyFactor", fSafetyFactor, 1.6 ) ;
1293  GetParam( "MaxXSec-FunctionCalls", fFunctionCalls ) ;
1294  GetParam( "MaxXSec-RelativeTolerance", fRelTolerance ) ;
1295  GetParam( "MaxXSec-MinScanPointsTmu", fMinScanPointsTmu ) ;
1296  GetParam( "MaxXSec-MinScanPointsCosth", fMinScanPointsCosth ) ;
1297 
1298  GetParam( "NSV-Q3Max", fQ3Max ) ;
1299 
1300  // Maximum allowed percentage deviation from the maximum cross section used
1301  // in the accept/reject loop for selecting lepton kinematics for SuSAv2.
1302  // Similar to the tolerance used by QELEventGenerator.
1303  GetParamDef( "SuSA-MaxXSec-DiffTolerance", fSuSAMaxXSecDiffTolerance, 999999. );
1304 }
1305 //___________________________________________________________________________
1307  const Range1D_t & Tl_range,
1308  const Range1D_t & ctl_range ) const {
1309 
1310  ROOT::Math::Minimizer * min = ROOT::Math::Factory::CreateMinimizer("Minuit2");
1311 
1312  double Enu = in.InitState().ProbeE(kRfHitNucRest);
1313  double LepMass = in.FSPrimLepton()->Mass();
1314 
1315  genie::utils::mec::gsl::d2Xsec_dTCosth f(fXSecModel,in, Enu, LepMass, -1.) ;
1316 
1317  std::array<string,2> names = { "Tl", "CosThetal" } ;
1318  std::array<Range1D_t,2> ranges = { Tl_range, ctl_range } ;
1319 
1320  std::array<double,2> start, steps, temp_point ;
1321  steps[0] = ( ranges[0].max - ranges[0].min ) / ( fMinScanPointsTmu + 1 ) ;
1322  steps[1] = ( ranges[1].max - ranges[1].min ) / ( fMinScanPointsCosth + 1 ) ;
1323 
1324  double xsec = 0 ;
1325 
1326  // preliimnary scan
1327  for ( unsigned int i = 1 ; i <= (unsigned int) fMinScanPointsTmu ; ++i ) {
1328  temp_point[0] = ranges[0].min + steps[0]*i ;
1329 
1330  for ( unsigned int j = 1 ; j <= (unsigned int) fMinScanPointsCosth ; ++j ) {
1331  temp_point[1] = ranges[1].min + steps[1]*j ;
1332 
1333  double temp_xsec = - f( temp_point.data() ) ;
1334  if ( temp_xsec > xsec ) {
1335  start = temp_point ;
1336  xsec = temp_xsec ;
1337  }
1338  }
1339  }
1340 
1341  // Set minimizer function and absolute tolerance :
1342  min->SetFunction( f );
1343  min->SetMaxFunctionCalls(fFunctionCalls);
1344  // min->SetTolerance( fRelTolerance * xsec );
1345 
1346  for ( unsigned int i = 0 ; i < ranges.size() ; ++i ) {
1347  min -> SetLimitedVariable( i, names[i], start[i], steps[i], ranges[i].min, ranges[i].max ) ;
1348  }
1349 
1350  min->Minimize();
1351 
1352  double max_xsec = -min->MinValue(); //back to positive xsec
1353 
1354  // Apply safety factor, since value retrieved from the cache might
1355  // correspond to a slightly different energy.
1356  max_xsec *= fSafetyFactor;
1357 
1358  return max_xsec ;
1359 }
1360 
1361 //___________________________________________________________________________
int Z(void) const
void ProcessEventRecord(GHepRecord *event) const
Basic constants.
TLorentzVector * GetX4(void) const
double beta(double KE, const simb::MCParticle *part)
void SelectNSVLeptonKinematics(GHepRecord *event) const
double RemovalEnergy(void) const
Definition: NuclearModelI.h:65
const NuclearModelI * fNuclModel
Definition: MECGenerator.h:72
double J(double q0, double q3, double Enu, double ml)
Definition: MECUtils.cxx:147
THE MAIN GENIE PROJECT NAMESPACE
Definition: AlgCmp.h:25
#define pERROR
Definition: Messenger.h:59
TRandom3 & RndLep(void) const
rnd number generator used by final state primary lepton generators
Definition: RandomGen.h:62
static const double kNucleonMass
Definition: Constants.h:77
double Q2(const Interaction *const i)
Definition: KineUtils.cxx:1064
static RandomGen * Instance()
Access instance.
Definition: RandomGen.cxx:71
int HitNucPdg(void) const
Definition: Target.cxx:304
void SetQ2(double Q2, bool selected=false)
Definition: Kinematics.cxx:255
const TLorentzVector * P4(void) const
Definition: GHepParticle.h:78
Kinematics * KinePtr(void) const
Definition: Interaction.h:76
Defines the EventRecordVisitorI interface. Concrete implementations of this interface use the &#39;Visito...
virtual double XSec(const Interaction *i, KinePhaseSpace_t k=kPSfE) const =0
Compute the cross section for the input interaction.
A simple [min,max] interval for doubles.
Definition: Range1.h:42
const int kPdgClusterNP
Definition: PDGCodes.h:215
void Configure(const Registry &config)
#define pFATAL
Definition: Messenger.h:56
double GetMaxXSecTlctl(const XSecAlgorithmI &xsec_model, const Interaction &inter, const double tolerance=0.01, const double safety_factor=1.2, const int max_n_layers=100)
Definition: MECUtils.cxx:334
Defines the EventGeneratorI interface.
Generated/set kinematical variables for an event.
Definition: Kinematics.h:39
double GetXSecMaxTlctl(const Interaction &inter, const Range1D_t &Tl_range, const Range1D_t &ctl_range) const
string P4AsString(const TLorentzVector *p)
Definition: PrintUtils.cxx:27
static const unsigned int kMaxUnweightDecayIterations
Definition: Controls.h:61
Pure abstract base class. Defines the NuclearModelI interface to be implemented by any physics model ...
Definition: NuclearModelI.h:46
double Mass(Resonance_t res)
resonance mass (GeV)
void SetMomentum(const TLorentzVector &p4)
intermediate_table::const_iterator const_iterator
Cluster finding and building.
A singleton holding random number generator classes. All random number generation in GENIE should tak...
Definition: RandomGen.h:29
void SetResonance(Resonance_t res)
Definition: XclsTag.cxx:128
double XYtoW(double Ev, double M, double x, double y)
Definition: KineUtils.cxx:1179
Computes the SuSAv2-MEC model differential cross section. Uses precomputed hadron tensor tables...
const TVector3 & Momentum3(void) const
Definition: NuclearModelI.h:75
A list of PDG codes.
Definition: PDGCodeList.h:32
const int kPdgClusterNN
Definition: PDGCodes.h:214
static const double kMinQ2Limit
Definition: Controls.h:41
virtual GHepParticle * Probe(void) const
Definition: GHepRecord.cxx:277
double y(bool selected=false) const
Definition: Kinematics.cxx:112
int Pdg(void) const
Definition: GHepParticle.h:63
void GenerateNSVInitialHadrons(GHepRecord *event) const
Summary information for an interaction.
Definition: Interaction.h:56
const TLorentzVector & HitNucP4(void) const
Definition: Target.h:91
An exception thrown by EventRecordVisitorI when the normal processing sequence has to be disrupted (f...
#define LOG(stream, priority)
A macro that returns the requested log4cpp::Category appending a string (using the FILE...
Definition: Messenger.h:96
void SetFSLeptonP4(const TLorentzVector &p4)
Definition: Kinematics.cxx:297
string Name(void) const
Definition: AlgId.h:44
static Config * config
Definition: config.cpp:1054
void RecoilNucleonCluster(GHepRecord *event) const
void DecayNucleonCluster(GHepRecord *event) const
TLorentzVector * GetP4(void) const
void WQ2toXY(double Ev, double M, double W, double Q2, double &x, double &y)
Definition: KineUtils.cxx:1119
const Kinematics & Kine(void) const
Definition: Interaction.h:71
A Neutrino Interaction Target. Is a transparent encapsulation of quite different physical systems suc...
Definition: Target.h:40
virtual void Configure(const Registry &config)
Definition: Algorithm.cxx:62
int ProbePdg(void) const
Definition: InitialState.h:64
bool Getq0q3FromTlCostl(double Tl, double costl, double Enu, double ml, double &q0, double &q3)
Definition: MECUtils.cxx:121
#define pINFO
Definition: Messenger.h:62
string tmp
Definition: languages.py:63
static int max(int a, int b)
TRandom3 & RndKine(void) const
rnd number generator used by kinematics generators
Definition: RandomGen.h:50
void SelectSuSALeptonKinematics(GHepRecord *event) const
Misc GENIE control constants.
PDGCodeList NucleonClusterConstituents(int pdgc) const
void SelectEmpiricalKinematics(GHepRecord *event) const
#define pWARN
Definition: Messenger.h:60
virtual const XSecAlgorithmI * CrossSectionAlg(void) const =0
TParticlePDG * FSPrimLepton(void) const
final state primary lepton
bool IsEM(void) const
void Setx(double x, bool selected=false)
Definition: Kinematics.cxx:231
XclsTag * ExclTagPtr(void) const
Definition: Interaction.h:77
void SetKV(KineVar_t kv, double value)
Definition: Kinematics.cxx:335
static RunningThreadInfo * Instance(void)
Kinematical utilities.
virtual const AlgId & Id(void) const
Get algorithm ID.
Definition: Algorithm.h:97
void SetW(double W, bool selected=false)
Definition: Kinematics.cxx:279
static PDGLibrary * Instance(void)
Definition: PDGLibrary.cxx:57
T min(sqlite3 *const db, std::string const &table_name, std::string const &column_name)
Definition: statistics.h:55
double Q2YtoX(double Ev, double M, double Q2, double y)
Definition: KineUtils.cxx:1209
string RgKey
const TLorentzVector * X4(void) const
Definition: GHepParticle.h:79
A registry. Provides the container for algorithm configuration parameters.
Definition: Registry.h:65
void SetHitNucPdg(int pdgc)
Definition: Target.cxx:171
void Sety(double y, bool selected=false)
Definition: Kinematics.cxx:243
int TgtPdg(void) const
Target * TgtPtr(void) const
Definition: InitialState.h:67
int IonPdgCode(int A, int Z)
Definition: PDGUtils.cxx:68
const XSecAlgorithmI * fXSecModel
Definition: MECGenerator.h:70
#define A
Definition: memgrp.cpp:38
void SetHadSystP4(const TLorentzVector &p4)
Definition: Kinematics.cxx:307
double gQ2
Definition: gtestDISSF.cxx:55
InitialState * InitStatePtr(void) const
Definition: Interaction.h:74
virtual bool GenerateNucleon(const Target &) const =0
const InitialState & InitState(void) const
Definition: Interaction.h:69
void AddTargetRemnant(GHepRecord *event) const
const ProcessInfo & ProcInfo(void) const
Definition: Interaction.h:70
void GenerateFermiMomentum(GHepRecord *event) const
const int kPdgProton
Definition: PDGCodes.h:81
TParticlePDG * Find(int pdgc, bool must_exist=true)
Definition: PDGLibrary.cxx:75
int A(void) const
double Q2(bool selected=false) const
Definition: Kinematics.cxx:125
#define pNOTICE
Definition: Messenger.h:61
void ClearRunningValues(void)
Definition: Kinematics.cxx:347
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:66
static std::vector< std::string > const names
Definition: FragmentType.hh:8
static const unsigned int kRjMaxIterations
Definition: Controls.h:26
const EventGeneratorI * RunningThread(void)
void SetPrimaryLeptonPolarization(GHepRecord *ev)
double ProbeE(RefFrame_t rf) const
const int kPdgNeutron
Definition: PDGCodes.h:83
GENIE&#39;s GHEP MC event record.
Definition: GHepRecord.h:45
static const double kPi
Definition: Constants.h:37
Most commonly used PDG codes. A set of utility functions to handle PDG codes is provided in PDGUtils...
TRandom3 & RndDec(void) const
rnd number generator used by decay models
Definition: RandomGen.h:56
STDHEP-like event record entry that can fit a particle or a nucleus.
Definition: GHepParticle.h:39
Keep info on the event generation thread currently on charge. This is used so that event generation m...
void AddFinalStateLepton(GHepRecord *event) const
Root of GENIE utility namespaces.
cet::coded_exception< error, detail::translate > exception
Definition: exception.h:33
void push_back(int pdg_code)
Definition: PDGCodeList.cxx:58
Event finding and building.
TGenPhaseSpace fPhaseSpaceGenerator
Definition: MECGenerator.h:71
enum genie::EGHepStatus GHepStatus_t
double fSuSAMaxXSecDiffTolerance
Definition: MECGenerator.h:84
Initial State information.
Definition: InitialState.h:48
#define pDEBUG
Definition: Messenger.h:63
const int kPdgClusterPP
Definition: PDGCodes.h:216
const Algorithm * SubAlg(const RgKey &registry_key) const
Definition: Algorithm.cxx:345