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