Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutrinoNucleusModel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id: G4NeutrinoNucleusModel.cc 91806 2015-08-06 12:20:45Z gcosmo $
27//
28// Geant4 Header : G4NeutrinoNucleusModel
29//
30// Author : V.Grichine 12.2.19
31//
32
34
35#include "G4SystemOfUnits.hh"
36#include "G4ParticleTable.hh"
38#include "G4IonTable.hh"
39#include "Randomize.hh"
40#include "G4RandomDirection.hh"
41
42//#include "G4Integrator.hh"
43#include "G4DataVector.hh"
44#include "G4PhysicsTable.hh"
45
46#include "G4CascadeInterface.hh"
47// #include "G4BinaryCascade.hh"
48#include "G4TheoFSGenerator.hh"
51#include "G4PreCompoundModel.hh"
54#include "G4FTFModel.hh"
55// #include "G4BinaryCascade.hh"
56#include "G4HadFinalState.hh"
57#include "G4HadSecondary.hh"
59// #include "G4INCLXXInterface.hh"
60#include "G4KineticTrack.hh"
63#include "G4Fragment.hh"
65
66// #include "G4QGSModel.hh"
67// #include "G4QGSMFragmentation.hh"
68// #include "G4QGSParticipants.hh"
69
70
71#include "G4MuonMinus.hh"
72#include "G4MuonPlus.hh"
73#include "G4Nucleus.hh"
74#include "G4LorentzVector.hh"
75
76using namespace std;
77using namespace CLHEP;
78
80
81const G4double G4NeutrinoNucleusModel::fResMass[6] = // [fResNumber] =
82 {2190., 1920., 1700., 1600., 1440., 1232. };
83
85
86const G4double G4NeutrinoNucleusModel::fMesMass[4] = {1260., 980., 770., 139.57};
87const G4int G4NeutrinoNucleusModel::fMesPDG[4] = {20213, 9000211, 213, 211};
88
89// const G4double G4NeutrinoNucleusModel::fBarMass[4] = {1905., 1600., 1232., 939.57};
90// const G4int G4NeutrinoNucleusModel::fBarPDG[4] = {2226, 32224, 2224, 2212};
91
92const G4double G4NeutrinoNucleusModel::fBarMass[4] = {1700., 1600., 1232., 939.57};
93const G4int G4NeutrinoNucleusModel::fBarPDG[4] = {12224, 32224, 2224, 2212};
94
96115.603, 133.424, 153.991, 177.729, 205.126, 236.746, 273.24, 315.361, 363.973, 420.08, 484.836, 559.573, 645.832,
97745.387, 860.289, 992.903, 1145.96, 1322.61, 1526.49, 1761.8, 2033.38, 2346.83, 2708.59, 3126.12, 3608.02, 4164.19,
984806.1, 5546.97, 6402.04, 7388.91, 8527.92, 9842.5, 11359.7, 13110.8, 15131.9, 17464.5, 20156.6, 23263.8, 26849.9,
9930988.8, 35765.7, 41279, 47642.2, 54986.3, 63462.4, 73245.2, 84536, 97567.2, 112607, 129966 };
100
101
104G4double G4NeutrinoNucleusModel::fNuMuQarrayKR[50][51][51] = {{{1.0}}};
105G4double G4NeutrinoNucleusModel::fNuMuQdistrKR[50][51][50] = {{{1.0}}};
106
107///////////////////////////////////////////
108
111{
112 SetMinEnergy( 0.0*GeV );
113 SetMaxEnergy( 100.*TeV );
114 SetMinEnergy(1.e-6*eV);
115
116 fNbin = 50;
117 fEindex = fXindex = fQindex = 0;
118 fOnePionIndex = 58;
119 fIndex = 50;
120 fCascade = fString = fProton = f2p2h = fBreak = false;
121
122 fNuEnergy = fQ2 = fQtransfer = fXsample = fDp = fTr = 0.;
123 fCosTheta = fCosThetaPi = 1.;
124 fEmuPi = fW2 = fW2pi = 0.;
125
126 fMu = 105.6583745*MeV;
127 fMpi = 139.57018*MeV;
128 fM1 = 939.5654133*MeV; // for nu_mu -> mu-, and n -> p
129 fM2 = 938.2720813*MeV;
130
131 fEmu = fMu;
132 fEx = fM1;
133 fMr = 1232.*MeV;
134 fMt = fM2; // threshold for N*-diffraction
135
137
138 fLVh = G4LorentzVector(0.,0.,0.,0.);
139 fLVl = G4LorentzVector(0.,0.,0.,0.);
140 fLVt = G4LorentzVector(0.,0.,0.,0.);
141 fLVcpi = G4LorentzVector(0.,0.,0.,0.);
142
145
146 // PDG2016: sin^2 theta Weinberg
147
148 fSin2tW = 0.23129; // 0.2312;
149
150 fCutEnergy = 0.; // default value
151
152
153 /*
154 // G4VPreCompoundModel* ptr;
155 // reuse existing pre-compound model as in binary cascade
156
157 fPrecoInterface = new G4GeneratorPrecompoundInterface ;
158
159 if( !fPreCompound )
160 {
161 G4HadronicInteraction* p =
162 G4HadronicInteractionRegistry::Instance()->FindModel("PRECO");
163 G4VPreCompoundModel* fPreCompound = static_cast<G4VPreCompoundModel*>(p);
164
165 if(!fPreCompound)
166 {
167 fPreCompound = new G4PreCompoundModel();
168 }
169 fPrecoInterface->SetDeExcitation(fPreCompound);
170 }
171 fDeExcitation = GetDeExcitation()->GetExcitationHandler();
172 */
173
178
179 fPDGencoding = 0; // unphysical as default
180 fRecoil = nullptr;
181
182}
183
184
186{
188}
189
190
191void G4NeutrinoNucleusModel::ModelDescription(std::ostream& outFile) const
192{
193
194 outFile << "G4NeutrinoNucleusModel is a neutrino-nucleus general\n"
195 << "model which uses the standard model \n"
196 << "transfer parameterization. The model is fully relativistic\n";
197
198}
199
200/////////////////////////////////////////////////////////
201
203 G4Nucleus & targetNucleus)
204{
205 G4bool result = false;
206 G4String pName = aPart.GetDefinition()->GetParticleName();
207 G4double energy = aPart.GetTotalEnergy();
208
209 if( pName == "nu_mu" // || pName == "anti_nu_mu" )
210 &&
211 energy > fMinNuEnergy )
212 {
213 result = true;
214 }
215 G4int Z = targetNucleus.GetZ_asInt();
216 Z *= 1;
217
218 return result;
219}
220
221//////////////////////////////////////
222
224{
225 G4int i(0), nBin(50);
226 G4double xx(0.), prob = G4UniformRand();
227
228 for( i = 0; i < nBin; ++i )
229 {
230 if( energy <= fNuMuEnergyLogVector[i] ) break;
231 }
232 if( i <= 0) // E-edge
233 {
234 fEindex = 0;
235 xx = GetXkr( 0, prob);
236 }
237 else if ( i >= nBin)
238 {
239 fEindex = nBin-1;
240 xx = GetXkr( nBin-1, prob);
241 }
242 else
243 {
244 fEindex = i;
245 G4double x1 = GetXkr(i-1,prob);
246 G4double x2 = GetXkr(i,prob);
247
250 G4double e = G4Log(energy);
251
252 if( e2 <= e1) xx = x1 + G4UniformRand()*(x2-x1);
253 else xx = x1 + (e-e1)*(x2-x1)/(e2-e1); // lin in energy log-scale
254 }
255 return xx;
256}
257
258//////////////////////////////////////////////
259//
260// sample X according to prob (xmin,1) at a given energy index iEnergy
261
263{
264 G4int i(0), nBin=50;
265 G4double xx(0.);
266
267 for( i = 0; i < nBin; ++i )
268 {
269 if( prob <= fNuMuXdistrKR[iEnergy][i] )
270 break;
271 }
272 if(i <= 0 ) // X-edge
273 {
274 fXindex = 0;
275 xx = fNuMuXarrayKR[iEnergy][0];
276 }
277 if ( i >= nBin )
278 {
279 fXindex = nBin;
280 xx = fNuMuXarrayKR[iEnergy][nBin];
281 }
282 else
283 {
284 fXindex = i;
285 G4double x1 = fNuMuXarrayKR[iEnergy][i];
286 G4double x2 = fNuMuXarrayKR[iEnergy][i+1];
287
288 G4double p1 = 0.;
289
290 if( i > 0 ) p1 = fNuMuXdistrKR[iEnergy][i-1];
291
292 G4double p2 = fNuMuXdistrKR[iEnergy][i];
293
294 if( p2 <= p1 ) xx = x1 + G4UniformRand()*(x2-x1);
295 else xx = x1 + (prob-p1)*(x2-x1)/(p2-p1);
296 }
297 return xx;
298}
299
300//////////////////////////////////////
301//
302// Sample fQtransfer at a given Enu and fX
303
305{
306 G4int nBin(50), iE=fEindex, jX=fXindex;
307 G4double qq(0.), qq1(0.), qq2(0.);
308 G4double prob = G4UniformRand();
309
310 // first E
311
312 if( iE <= 0 )
313 {
314 qq1 = GetQkr( 0, jX, prob);
315 }
316 else if ( iE >= nBin-1)
317 {
318 qq1 = GetQkr( nBin-1, jX, prob);
319 }
320 else
321 {
322 G4double q1 = GetQkr(iE-1,jX, prob);
323 G4double q2 = GetQkr(iE,jX, prob);
324
327 G4double e = G4Log(energy);
328
329 if( e2 <= e1) qq1 = q1 + G4UniformRand()*(q2-q1);
330 else qq1 = q1 + (e-e1)*(q2-q1)/(e2-e1); // lin in energy log-scale
331 }
332
333 // then X
334
335 if( jX <= 0 )
336 {
337 qq2 = GetQkr( iE, 0, prob);
338 }
339 else if ( jX >= nBin)
340 {
341 qq2 = GetQkr( iE, nBin, prob);
342 }
343 else
344 {
345 G4double q1 = GetQkr(iE,jX-1, prob);
346 G4double q2 = GetQkr(iE,jX, prob);
347
348 G4double e1 = G4Log(fNuMuXarrayKR[iE][jX-1]);
349 G4double e2 = G4Log(fNuMuXarrayKR[iE][jX]);
350 G4double e = G4Log(xx);
351
352 if( e2 <= e1) qq2 = q1 + G4UniformRand()*(q2-q1);
353 else qq2 = q1 + (e-e1)*(q2-q1)/(e2-e1); // lin in energy log-scale
354 }
355 qq = 0.5*(qq1+qq2);
356
357 return qq;
358}
359
360//////////////////////////////////////////////
361//
362// sample Q according to prob (qmin,qmax) at a given energy index iE and X index jX
363
365{
366 G4int i(0), nBin=50;
367 G4double qq(0.);
368
369 for( i = 0; i < nBin; ++i )
370 {
371 if( prob <= fNuMuQdistrKR[iE][jX][i] )
372 break;
373 }
374 if(i <= 0 ) // Q-edge
375 {
376 fQindex = 0;
377 qq = fNuMuQarrayKR[iE][jX][0];
378 }
379 if ( i >= nBin )
380 {
381 fQindex = nBin;
382 qq = fNuMuQarrayKR[iE][jX][nBin];
383 }
384 else
385 {
386 fQindex = i;
387 G4double q1 = fNuMuQarrayKR[iE][jX][i];
388 G4double q2 = fNuMuQarrayKR[iE][jX][i+1];
389
390 G4double p1 = 0.;
391
392 if( i > 0 ) p1 = fNuMuQdistrKR[iE][jX][i-1];
393
394 G4double p2 = fNuMuQdistrKR[iE][jX][i];
395
396 if( p2 <= p1 ) qq = q1 + G4UniformRand()*(q2-q1);
397 else qq = q1 + (prob-p1)*(q2-q1)/(p2-p1);
398 }
399 return qq;
400}
401
402
403///////////////////////////////////////////////////////////
404//
405// Final meson to theParticleChange
406
408{
409 G4int pdg = pdgM;
410 // if ( qM == 0 ) pdg = pdgM - 100;
411 // else if ( qM == -1 ) pdg = -pdgM;
412
413 if( pdg == 211 || pdg == -211 || pdg == 111) // pions
414 {
416 G4DynamicParticle* dp2 = new G4DynamicParticle( pd2, lvM);
418 }
419 else // meson resonances
420 {
422 FindParticle(pdg);
423 G4KineticTrack ddkt( rePart, 0., G4ThreeVector(0.,0.,0.), lvM);
424 G4KineticTrackVector* ddktv = ddkt.Decay();
425
426 G4DecayKineticTracks decay( ddktv );
427
428 for( unsigned int i = 0; i < ddktv->size(); i++ ) // add products to partchange
429 {
430 G4DynamicParticle * aNew =
431 new G4DynamicParticle( ddktv->operator[](i)->GetDefinition(),
432 ddktv->operator[](i)->Get4Momentum());
433
434 // G4cout<<" "<<i<<", "<<aNew->GetDefinition()->GetParticleName()<<", "<<aNew->Get4Momentum()<<G4endl;
435
437 delete ddktv->operator[](i);
438 }
439 delete ddktv;
440 }
441}
442
443////////////////////////////////////////////////////////
444//
445// Final barion to theParticleChange, and recoil nucleus treatment
446
448{
449 G4int A(0), Z(0), pdg = pdgB;
450 // G4bool FiNucleon(false);
451
452 // if ( qB == 1 ) pdg = pdgB - 10;
453 // else if ( qB == 0 ) pdg = pdgB - 110;
454 // else if ( qB == -1 ) pdg = pdgB - 1110;
455
456 if( pdg == 2212 || pdg == 2112)
457 {
459 // FiNucleon = true;
460 }
461 else fMr = lvB.m();
462
464 lvB.boost(-bst); // in fLVt rest system
465
466 G4double eX = lvB.e();
467 G4double det(0.), det2(0.), rM(0.), mX = lvB.m();
468 G4ThreeVector dX = (lvB.vect()).unit();
469 G4double pX = sqrt(eX*eX-mX*mX);
470
471 if( fRecoil )
472 {
473 Z = fRecoil->GetZ_asInt();
474 A = fRecoil->GetA_asInt();
475 rM = fRecoil->AtomicMass(A,Z); //->AtomicMass(); //
476 rM = fLVt.m();
477 }
478 else // A=0 nu+p
479 {
480 A = 0;
481 Z = 1;
482 rM = electron_mass_c2;
483 }
484 // G4cout<<A<<", ";
485
486 G4double sumE = eX + rM;
487 G4double B = sumE*sumE + rM*rM - fMr*fMr - pX*pX;
488 G4double a = 4.*(sumE*sumE - pX*pX);
489 G4double b = -4.*B*pX;
490 G4double c = 4.*sumE*sumE*rM*rM - B*B;
491 det2 = b*b-4.*a*c;
492 if( det2 <= 0. ) det = 0.;
493 else det = sqrt(det2);
494 G4double dP = 0.5*(-b - det )/a;
495
496 fDp = dP;
497
498 pX -= dP;
499
500 if(pX < 0.) pX = 0.;
501
502 // if( A == 0 ) G4cout<<pX/MeV<<", ";
503
504 eX = sqrt( pX*pX + fMr*fMr );
505 G4LorentzVector lvN( pX*dX, eX );
506 lvN.boost(bst); // back to lab
507
508 if( pdg == 2212 || pdg == 2112) // nucleons mX >= fMr, dP >= 0
509 {
511 G4DynamicParticle* dp2 = new G4DynamicParticle( pd2, lvN);
513
514 }
515 else // delta resonances
516 {
518 G4KineticTrack ddkt( rePart, 0., G4ThreeVector( 0., 0., 0.), lvN);
519 G4KineticTrackVector* ddktv = ddkt.Decay();
520
521 G4DecayKineticTracks decay( ddktv );
522
523 for( unsigned int i = 0; i < ddktv->size(); i++ ) // add products to partchange
524 {
525 G4DynamicParticle * aNew =
526 new G4DynamicParticle( ddktv->operator[](i)->GetDefinition(),
527 ddktv->operator[](i)->Get4Momentum() );
528
529 // G4cout<<" "<<i<<", "<<aNew->GetDefinition()->GetParticleName()<<", "<<aNew->Get4Momentum()<<G4endl;
530
532 delete ddktv->operator[](i);
533 }
534 delete ddktv;
535 }
536 // recoil nucleus
537
538 G4double eRecoil = sqrt( rM*rM + dP*dP );
539 fTr = eRecoil - rM;
540 G4ThreeVector vRecoil(dP*dX);
541 // dP += G4UniformRand()*10.*MeV;
542 G4LorentzVector rec4v(vRecoil, 0.);
543 rec4v.boost(bst); // back to lab
544 fLVt += rec4v;
545 const G4LorentzVector lvTarg = fLVt; // (vRecoil, eRecoil);
546
547
548 if( fRecoil ) // proton*?
549 {
551 G4double exE = fLVt.m() - grM;
552 if( exE < 5.*MeV ) exE = 5.*MeV + G4UniformRand()*10.*MeV;
553
554 const G4LorentzVector in4v( G4ThreeVector( 0., 0., 0.), grM );
555 G4Fragment fragment( A, Z, in4v); // lvTarg );
556 fragment.SetNumberOfHoles(1);
557 fragment.SetExcEnergyAndMomentum( exE, lvTarg );
558
559 RecoilDeexcitation(fragment);
560 }
561 else // momentum?
562 {
564 }
565}
566
567
568///////////////////////////////////////////////////////
569//
570// Get final particles from excited recoil nucleus and write them to theParticleChange, delete the particle vector
571
573{
574 G4ReactionProductVector* products = fPreCompound->DeExcite(fragment);
575
576 if( products != nullptr )
577 {
578 for( auto iter = products->cbegin(); iter != products->cend(); ++iter )
579 // for( auto & prod : products ) // prod = (*iter) is the pointer to final hadronic particle
580 {
581 theParticleChange.AddSecondary(new G4DynamicParticle( (*iter)->GetDefinition(),
582 (*iter)->GetTotalEnergy(),
583 (*iter)->GetMomentum() ) );
584 // delete prod;
585 }
586 // delete products;
587 products->clear();
588 }
589 return;
590}
591
592
593///////////////////////////////////////////
594//
595// Fragmentation of lvX directly to pion and recoil nucleus (A,Z): fLVh + fLVt -> pi + A*
596
598{
599 G4int A(0), Z(0), pdg = pdgP;
600 fLVcpi = G4LorentzVector(0.,0.,0.,0.);
601
602 G4double rM(0.), mN(938.), det(0.), det2(0.);
603 G4double mI(0.);
604 mN = G4ParticleTable::GetParticleTable()->FindParticle(2212)->GetPDGMass(); // *0.85; // *0.9; //
605
606 // mN = 1.*139.57 + G4UniformRand()*(938. - 1.*139.57);
607
608 G4ThreeVector vN = lvP.boostVector(), bst(0.,0.,0.);
609 // G4double gN = lvP.e()/lvP.m();
610 // G4LorentzVector lvNu(vN*gN*mN, mN*gN);
611 G4LorentzVector lvNu(0.,0.,0., mN); // lvNu(bst, mN);
612 lvP.boost(-vN); // 9-3-20
613 lvP = lvP - lvNu; // 9-3-20 already 1pi
614 lvP.boost(vN); // 9-3-20
615 lvNu.boost(vN); // 9-3-20
616
617 // G4cout<<vN-lvP.boostVector()<<", ";
618
619 Z = targetNucleus.GetZ_asInt();
620 A = targetNucleus.GetA_asInt();
621 rM = targetNucleus.AtomicMass(A,Z); //->AtomicMass(); //
622
623 // G4cout<<rM<<", ";
624 // G4cout<<A<<", ";
625
626 if( A == 1 )
627 {
628 bst = vN; // lvNu.boostVector(); // 9-3-20
629 // mI = 0.; // 9-3-20
630 rM = mN;
631 }
632 else
633 {
634 G4Nucleus targ(A-1,Z);
635 mI = targ.AtomicMass(A-1,Z);
636 G4LorentzVector lvTar(0.,0.,0.,mI);
637 lvNu = lvNu + lvTar;
638 bst = lvNu.boostVector();
639 // bst = fLVt.boostVector(); // to recoil rest frame
640 // G4cout<<fLVt<<" "<<bst<<G4endl;
641 }
642 lvP.boost(-bst); // 9-3-20
644 G4double eX = lvP.e();
645 G4double mX = lvP.m();
646 // G4cout<<mX-fMr<<", ";
647 G4ThreeVector dX = (lvP.vect()).unit();
648 // G4cout<<dX<<", ";
649 G4double pX = sqrt(eX*eX-mX*mX);
650 // G4cout<<pX<<", ";
651 G4double sumE = eX + rM;
652 G4double B = sumE*sumE + rM*rM - fMr*fMr - pX*pX;
653 G4double a = 4.*(sumE*sumE - pX*pX);
654 G4double b = -4.*B*pX;
655 G4double c = 4.*sumE*sumE*rM*rM - B*B;
656 det2 = b*b-4.*a*c;
657 if(det2 > 0.) det = sqrt(det2);
658 G4double dP = 0.5*(-b - det )/a;
659
660 // dP = FinalMomentum( mI, rM, fMr, lvP);
661 dP = FinalMomentum( rM, rM, fMr, lvP); // 9-3-20
662
663 // G4cout<<dP<<", ";
664 pX -= dP;
665 if( pX < 0. ) pX = 0.;
666
667 eX = sqrt( dP*dP + fMr*fMr );
668 G4LorentzVector lvN( dP*dX, eX );
669
670 if( A >= 1 ) lvN.boost(bst); // 9-3-20 back to lab
671
672 fLVcpi = lvN;
673
675 G4DynamicParticle* dp2 = new G4DynamicParticle( pd2, lvN);
676 theParticleChange.AddSecondary( dp2 ); // coherent pion
677
678 // recoil nucleus
679
680 G4double eRecoil = sqrt( rM*rM + pX*pX );
681 G4ThreeVector vRecoil(pX*dX);
682 G4LorentzVector lvTarg1( vRecoil, eRecoil);
683 lvTarg1.boost(bst);
684
685 const G4LorentzVector lvTarg = lvTarg1;
686
687 if( A > 1 ) // recoil target nucleus*
688 {
690 G4double exE = fLVt.m() - grM;
691
692 if( exE < 5.*MeV ) exE = 5.*MeV + G4UniformRand()*10.*MeV; // vmg???
693
694 const G4LorentzVector in4v( G4ThreeVector( 0., 0., 0.), grM );
695 G4Fragment fragment( A, Z, in4v); // lvTarg );
696 fragment.SetNumberOfHoles(1);
697 fragment.SetExcEnergyAndMomentum( exE, lvTarg );
698
699 RecoilDeexcitation(fragment);
700 }
701 else // recoil target proton
702 {
703 G4double eTkin = eRecoil - rM;
704 G4double eTh = 0.01*MeV; // 10.*MeV;
705
706 if( eTkin > eTh )
707 {
710 }
712 }
713 return;
714}
715
716
717////////////////////////////////////////////////////////////
718//
719// Excited barion decay to meson and barion,
720// mass distributions and charge exchange are free parameters
721
723{
724 G4bool finB = false;
725 G4int pdgB(0), i(0), qM(0), qB(0); // pdgM(0),
726 G4double mM(0.), mB(0.), eM(0.), eB(0.), pM(0.), pB(0.);
727 G4double mm1(0.), mm22(0.), M1(0.), M2(0.), mX(0.);
728
729 mX = lvX.m();
730
733
734 // G4double deltaM = 1.*MeV; // 30.*MeV; // 10.*MeV; // 100.*MeV; // 20.*MeV; //
735 G4double deltaMr[4] = { 0.*MeV, 0.*MeV, 100.*MeV, 0.*MeV};
736
737 G4ThreeVector dir(0.,0.,0.);
738 G4ThreeVector bst(0.,0.,0.);
739 G4LorentzVector lvM(0.,0.,0.,0.);
740 G4LorentzVector lvB(0.,0.,0.,0.);
741
742 for( i = 0; i < fClustNumber; ++i) // check resonance
743 {
744 if( mX >= fBarMass[i] )
745 {
746 pdgB = fBarPDG[i];
747 // mB = G4ParticleTable::GetParticleTable()->FindParticle(pdgB)->GetPDGMass();
748 break;
749 }
750 }
751 if( i == fClustNumber || i == fClustNumber-1 ) // low mass, p || n
752 {
753 if ( qX == 2 || qX == 0) { pdgB = 2212; qB = 1;} // p for 2, 0
754
755 else if( qX == 1 || qX == -1) { pdgB = 2112; qB = 0;} // n for 1, -1
756
757 return FinalBarion( lvX, qB, pdgB);
758 }
759 else if( mX < fBarMass[i] + deltaMr[i] || mX < mN + mPi )
760 {
761 finB = true; // final barion -> out
762
763 if ( qX == 1 && pdgB != 2212) pdgB = pdgB - 10;
764 else if( qX == 0 && pdgB != 2212) pdgB = pdgB - 110;
765 else if( qX == 0 && pdgB == 2212) pdgB = pdgB - 100;
766
767 if( finB ) return FinalBarion( lvX, qX, pdgB ); // out
768 }
769 // no barion resonance, try 1->2 decay in COM frame
770
771 // try meson mass
772
773 mm1 = mPi + 1.*MeV; // pi+
774 mm22 = mX - mN; // mX-n
775
776 if( mm22 <= mm1 ) // out with p or n
777 {
778 if( qX == 2 || qX == 0) { pdgB = 2212; qB = 1;} // p
779 else if( qX == 1 || qX == -1) { pdgB = 2112; qB = 0;} // n
780
781 return FinalBarion(lvX, qB, pdgB);
782 }
783 else // try decay -> meson(cluster) + barion(cluster)
784 {
785 // G4double sigmaM = 50.*MeV; // 100.*MeV; // 200.*MeV; // 400.*MeV; // 800.*MeV; //
786 G4double rand = G4UniformRand();
787
788 // mM = mm1*mm22/( mm1 + rand*(mm22 - mm1) );
789 // mM = mm1*mm22/sqrt( mm1*mm1 + rand*(mm22*mm22 - mm1*mm1) );
790 // mM = -sigmaM*log( (1.- rand)*exp(-mm22/sigmaM) + rand*exp(-mm1/sigmaM) );
791 mM = mm1 + rand*(mm22-mm1);
792
793
794 for( i = 0; i < fClustNumber; ++i)
795 {
796 if( mM >= fMesMass[i] )
797 {
798 // pdgM = fMesPDG[i];
799 // mM = G4ParticleTable::GetParticleTable()->FindParticle(pdgM)->GetPDGMass();
800 break;
801 }
802 }
803 M1 = G4ParticleTable::GetParticleTable()->FindParticle(2112)->GetPDGMass()+2.*MeV; // n
804 M2 = mX - mM;
805
806 if( M2 <= M1 ) //
807 {
808 if ( qX == 2 || qX == 0) { pdgB = 2212; qB = 1;} // p
809 else if( qX == 1 || qX == -1) { pdgB = 2112; qB = 0;} // n
810
811 return FinalBarion(lvX, qB, pdgB);
812 }
813 mB = M1 + G4UniformRand()*(M2-M1);
814 // mB = -sigmaM*log( (1.- rand)*exp(-M2/sigmaM) + rand*exp(-M1/sigmaM) );
815
816 bst = lvX.boostVector();
817
818 // dir = G4RandomDirection(); // ???
819 // dir = G4ThreeVector(0.,0.,1.);
820 dir = bst.orthogonal().unit(); // ??
821 // G4double cost = exp(-G4UniformRand());
822 // G4double sint = sqrt((1.-cost)*(1.+cost));
823 // G4double phi = twopi*G4UniformRand();
824 // dir = G4ThreeVector(sint*cos(phi), sint*sin(phi), cost);
825
826 eM = 0.5*(mX*mX + mM*mM - mB*mB)/mX;
827 pM = sqrt(eM*eM - mM*mM);
828 lvM = G4LorentzVector( pM*dir, eM);
829 lvM.boost(bst);
830
831 eB = 0.5*(mX*mX + mB*mB - mM*mM)/mX;
832 pB = sqrt(eB*eB - mB*mB);
833 lvB = G4LorentzVector(-pB*dir, eB);
834 lvB.boost(bst);
835
836 // G4cout<<mM<<"/"<<mB<<", ";
837
838 // charge exchange
839
840 if ( qX == 2 ) { qM = 1; qB = 1;}
841 else if( qX == 1 ) { qM = 0; qB = 1;}
842 else if( qX == 0 ) { qM = 0; qB = 0;}
843 else if( qX == -1 ) { qM = -1; qB = 0;}
844
845 // if ( qM == 0 ) pdgM = pdgM - 100;
846 // else if( qM == -1 ) pdgM = -pdgM;
847
848 MesonDecay( lvM, qM); // pdgM ); //
849
850 // else
851 ClusterDecay( lvB, qB ); // continue
852 }
853}
854
855
856////////////////////////////////////////////////////////////
857//
858// Excited barion decay to meson and barion,
859// mass distributions and charge exchange are free parameters
860
862{
863 G4bool finB = false;
864 G4int pdgM(0), pdgB(0), i(0), qM(0), qB(0);
865 G4double mM(0.), mB(0.), eM(0.), eB(0.), pM(0.), pB(0.);
866 G4double mm1(0.), mm22(0.), M1(0.), M2(0.), mX(0.), Tkin(0.);
867
868 mX = lvX.m();
869 Tkin = lvX.e() - mX;
870
871 // if( mX < 1120*MeV && mX > 1020*MeV ) // phi(1020)->K+K-
872 if( mX < 1080*MeV && mX > 990*MeV && Tkin < 600*MeV ) // phi(1020)->K+K-
873 {
874 return FinalMeson( lvX, qB, 333);
875 }
877
878 G4double deltaMr[4] = { 0.*MeV, 0.*MeV, 100.*MeV, 0.*MeV};
879
880 G4ThreeVector dir(0.,0.,0.);
881 G4ThreeVector bst(0.,0.,0.);
882 G4LorentzVector lvM(0.,0.,0.,0.);
883 G4LorentzVector lvB(0.,0.,0.,0.);
884
885 for( i = 0; i < fClustNumber; ++i) // check resonance
886 {
887 if( mX >= fMesMass[i] )
888 {
889 pdgB = fMesPDG[i];
890 // mB = G4ParticleTable::GetParticleTable()->FindParticle(pdgB)->GetPDGMass();
891 break;
892 }
893 }
894 if( i == fClustNumber ) // || i == fClustNumber-1 ) // low mass, p || n
895 {
896 if ( qX == 1) { pdgB = 211; qB = 1;} // pi+
897 else if( qX == 0 ) { pdgB = 111; qB = 0;} // pi0
898 else if( qX == -1) { pdgB = -211; qB = -1;} // pi-
899
900 return FinalMeson( lvX, qB, pdgB);
901 }
902 else if( mX < fMesMass[i] + deltaMr[i] ) // || mX < mPi + mPi ) //
903 {
904 finB = true; // final barion -> out
905 pdgB = fMesPDG[i];
906
907 // if ( qX == 1 && pdgB != 2212) pdgB = pdgB - 10;
908
909 if( qX == 0 ) pdgB = pdgB - 100;
910 else if( qX == -1 ) pdgB = -pdgB;
911
912 if( finB ) return FinalMeson( lvX, qX, pdgB ); // out
913 }
914 // no resonance, try 1->2 decay in COM frame
915
916 // try meson
917
918 mm1 = mPi + 1.*MeV; // pi+
919 mm22 = mX - mPi - 1.*MeV; // mX-n
920
921 if( mm22 <= mm1 ) // out
922 {
923 if ( qX == 1) { pdgB = 211; qB = 1;} // pi+
924 else if( qX == 0 ) { pdgB = 111; qB = 0;} // pi0
925 else if( qX == -1) { pdgB = -211; qB = -1;} // pi-
926
927 return FinalMeson(lvX, qB, pdgB);
928 }
929 else // try decay -> pion + meson(cluster)
930 {
931 // G4double sigmaM = 50.*MeV; // 100.*MeV; // 200.*MeV; // 400.*MeV; // 800.*MeV; //
932 G4double rand = G4UniformRand();
933
934 if ( qX == 1 ) { qM = 1; qB = 0;}
935 else if( qX == 0 ) { qM = -1; qB = 1;} // { qM = 0; qB = 0;} //
936 else if( qX == -1 ) { qM = -1; qB = 0;}
937 /*
938 mM = mPi;
939 if(qM == 0) mM = G4ParticleTable::GetParticleTable()->FindParticle(111)->GetPDGMass(); //pi0
940 pdgM = fMesPDG[fClustNumber-1];
941 */
942 // mm1*mm22/( mm1 + rand*(mm22 - mm1) );
943 // mM = mm1*mm22/sqrt( mm1*mm1 + rand*(mm22*mm22 - mm1*mm1) );
944 // mM = -sigmaM*log( (1.- rand)*exp(-mm22/sigmaM) + rand*exp(-mm1/sigmaM) );
945 mM = mm1 + rand*(mm22-mm1);
946 // mM = mm1 + 0.9*(mm22-mm1);
947
948
949 for( i = 0; i < fClustNumber; ++i)
950 {
951 if( mM >= fMesMass[i] )
952 {
953 pdgM = fMesPDG[i];
954 // mM = G4ParticleTable::GetParticleTable()->FindParticle(pdgM)->GetPDGMass();
955 break;
956 }
957 }
958 if( i == fClustNumber || i == fClustNumber-1 ) // low mass, p || n
959 {
960 if ( qX == 1) { pdgB = 211; qB = 1;} // pi+
961 else if( qX == 0 ) { pdgB = 111; qB = 0;} // pi0
962 else if( qX == -1) { pdgB = -211; qB = -1;} // pi-
963
964 return FinalMeson( lvX, qB, pdgB);
965 }
966 else if( mX < fMesMass[i] + deltaMr[i] ) // || mX < mPi + mPi ) //
967 {
968 finB = true; // final barion -> out
969 pdgB = fMesPDG[i];
970
971 // if ( qX == 1 && pdgB != 2212) pdgB = pdgB - 10;
972
973 if( qX == 0 ) pdgB = pdgB - 100;
974 else if( qX == -1 ) pdgB = -pdgB;
975
976 if( finB ) return FinalMeson( lvX, qX, pdgB ); // out
977 }
978
980 M2 = mX - mM;
981
982 if( M2 <= M1 ) //
983 {
984 if ( qX == 1) { pdgB = 211; qB = 1;} // pi+
985 else if( qX == 0 ) { pdgB = 111; qB = 0;} // pi0
986 else if( qX == -1) { pdgB = -211; qB = -1;} // pi-
987
988 return FinalMeson(lvX, qB, pdgB);
989 }
990 mB = M1 + G4UniformRand()*(M2-M1);
991 // mB = -sigmaM*log( (1.- rand)*exp(-M2/sigmaM) + rand*exp(-M1/sigmaM) );
992 // mB = M1 + 0.9*(M2-M1);
993
994 bst = lvX.boostVector();
995
996 // dir = G4RandomDirection();
997 dir = bst.orthogonal().unit();
998
999 eM = 0.5*(mX*mX + mM*mM - mB*mB)/mX;
1000 pM = sqrt(eM*eM - mM*mM);
1001 lvM = G4LorentzVector( pM*dir, eM);
1002 lvM.boost(bst);
1003
1004 eB = 0.5*(mX*mX + mB*mB - mM*mM)/mX;
1005 pB = sqrt(eB*eB - mB*mB);
1006 lvB = G4LorentzVector(-pB*dir, eB);
1007 lvB.boost(bst);
1008
1009 // G4cout<<mM<<"/"<<mB<<", ";
1010
1011 // charge exchange
1012
1013 // if ( qX == 2 ) { qM = 1; qB = 1;}
1014
1015 if ( qM == 0 ) pdgM = pdgM - 100;
1016 else if( qM == -1 ) pdgM = -pdgM;
1017
1018 MesonDecay( lvM, qM ); //
1019
1020 MesonDecay( lvB, qB ); // continue
1021 }
1022}
1023
1024///////////////////////////////////////////////////////////////////////
1025//
1026// return final momentum x in the reaction lvX + mI -> mF + mP with momenta p-x, x
1027
1029{
1030 G4double result(0.), delta(0.);
1031 // G4double mI2 = mI*mI;
1032 G4double mF2 = mF*mF;
1033 G4double mP2 = mP*mP;
1034 G4double eX = lvX.e();
1035 // G4double mX = lvX.m();
1036 G4double pX = lvX.vect().mag();
1037 G4double pX2 = pX*pX;
1038 G4double sI = eX + mI;
1039 G4double sI2 = sI*sI;
1040 G4double B = sI2 - mF2 -pX2 + mP2;
1041 G4double B2 = B*B;
1042 G4double a = 4.*(sI2-pX2);
1043 G4double b = -4.*B*pX;
1044 G4double c = 4.*sI2*mP2 - B2;
1045 G4double delta2 = b*b -4.*a*c;
1046
1047 if( delta2 >= 0. ) delta = sqrt(delta2);
1048
1049 result = 0.5*(-b-delta)/a;
1050 // result = 0.5*(-b+delta)/a;
1051
1052 return result;
1053}
1054
1055/////////////////////////////////////////////////////////////////
1056//
1057//
1058
1060{
1061 G4int Z = targetNucleus.GetZ_asInt();
1062 G4int A = targetNucleus.GetA_asInt();
1063
1064 G4double kF(250.*MeV);
1065 G4double kp = 365.*MeV;
1066 G4double kn = 231.*MeV;
1067 G4double t1 = 0.479;
1068 G4double t2 = 0.526;
1069 G4double ZpA = G4double(Z)/G4double(A);
1070 G4double NpA = 1. - ZpA;
1071
1072 if ( Z == 1 && A == 1 ) { kF = 0.; } // hydrogen ???
1073 else if ( Z == 1 && A == 2 ) { kF = 87.*MeV; }
1074 else if ( Z == 2 && A == 3 ) { kF = 134.*MeV; }
1075 else if ( Z == 6 && A == 12 ) { kF = 221.*MeV; }
1076 else if ( Z == 14 && A == 28 ) { kF = 239.*MeV; }
1077 else if ( Z == 26 && A == 56 ) { kF = 257.*MeV; }
1078 else if ( Z == 82 && A == 208 ) { kF = 265.*MeV; }
1079 else
1080 {
1081 kF = kp*ZpA*( 1 - pow( G4double(A), -t1 ) ) + kn*NpA*( 1 - pow( G4double(A), -t2 ) );
1082 }
1083 return kF;
1084}
1085
1086/////////////////////////////////////////////////////////////////
1087//
1088// sample nucleon momentum of Fermi motion for 1p1h and 2p2h modes
1089
1091{
1092 G4int A = targetNucleus.GetA_asInt();
1093 G4double kF = FermiMomentum( targetNucleus);
1094 G4double mom(0.), kCut = 0.5*GeV; // kCut = 1.*GeV; // kCut = 2.*GeV; // kCut = 4.*GeV; //
1095 // G4double cof = 2./GeV;
1096 // G4double ksi = kF*kF*cof*cof/pi/pi;
1097 G4double th = 1.; // 1. - 6.*ksi; //
1098
1099 if( G4UniformRand() < th || A < 3 ) // 1p1h
1100 {
1101 mom = kF*pow( G4UniformRand(), 1./3.);
1102 }
1103 else // 2p2h
1104 {
1105 mom = kF*kCut;
1106 mom /= kCut - G4UniformRand()*(kCut - kF);
1107 f2p2h = true;
1108 }
1109 return mom;
1110}
1111
1112
1113//////////////////////////////////////////////////////
1114//
1115// Excitation energy according Bodek
1116
1118{
1119 G4double eX(10.*MeV), a1(0.), a2(0.), e1(0.), e2(0.), aa = G4double(A);
1120 G4int i(0);
1121 const G4int maxBin = 12;
1122
1123 G4double refA[maxBin] = { 2., 6., 12., 16., 27., 28., 40., 50., 56., 58., 197., 208. };
1124
1125 G4double pEx[maxBin] = { 0., 12.2, 10.1, 10.9, 21.6, 12.4, 17.8, 17., 19., 16.8, 19.5, 14.7 };
1126
1127 G4double nEx[maxBin] = { 0., 12.2, 10., 10.2, 21.6, 12.4, 21.8, 17., 19., 16.8, 19.5, 16.9 };
1128
1129 G4DataVector dE(12,0.);
1130
1131 if(fP) for( i = 0; i < maxBin; ++i ) dE[i] = pEx[i];
1132 else dE[i] = nEx[i];
1133
1134 for( i = 0; i < maxBin; ++i )
1135 {
1136 if( aa <= refA[i] ) break;
1137 }
1138 if( i >= maxBin ) eX = dE[maxBin-1];
1139 else if( i <= 0 ) eX = dE[0];
1140 else
1141 {
1142 a1 = refA[i-1];
1143 a2 = refA[i];
1144 e1 = dE[i-1];
1145 e2 = dE[i];
1146 if (a1 == a2 || e1 == e2 ) eX = dE[i];
1147 else eX = e1 + (e2-e1)*(aa-a1)/(a2-a1);
1148 }
1149 return eX;
1150}
1151
1152
1153
1154///////////////////////////////////////////////////////
1155//
1156// Two G-function sampling for the nucleon momentum
1157
1159{
1160 f2p2h = false;
1161 G4double /* distr(0.), tail(0.), */ shift(1.), xx(1.), mom(0.), th(0.1);
1162 G4double kF = FermiMomentum( nucl);
1163 G4double momMax = 2.*kF; // 1.*GeV; // 1.*GeV; //
1164 G4double aa = 5.5;
1165 G4double ll = 6.0; // 6.5; //
1166
1167 G4int A = nucl.GetA_asInt();
1168
1169 if( A <= 12) th = 0.1;
1170 else
1171 {
1172 // th = 0.1/(1.+log(G4double(A)/12.));
1173 th = 1.2/( G4double(A) + 1.35*log(G4double(A)/12.) );
1174 }
1175 shift = 0.99; // 0.95; //
1176 xx = mom/shift/kF;
1177
1178 G4double rr = G4UniformRand();
1179
1180 if( rr > th )
1181 {
1182 aa = 5.5;
1183
1184 if( A <= 12 ) ll = 6.0;
1185 else
1186 {
1187 ll = 6.0 + 1.35*log(G4double(A)/12.);
1188 }
1189 xx = RandGamma::shoot(aa,ll);
1190 shift = 0.99;
1191 mom = xx*shift*kF;
1192 }
1193 else
1194 {
1195 f2p2h = true;
1196 aa = 6.5;
1197 ll = 6.5;
1198 xx = RandGamma::shoot(aa,ll);
1199 shift = 2.5;
1200 mom = xx*shift*kF;
1201 }
1202 if( mom > momMax ) mom = G4UniformRand()*momMax;
1203 if( mom > 2.*kF ) f2p2h = true;
1204
1205 // mom = 0.;
1206
1207 return mom;
1208}
1209
1210
1211///////////////////////////////////// experimental arrays and get functions ////////////////////////////////////////
1212//
1213// Return index of nu/anu energy array corresponding to the neutrino energy
1214
1216{
1217 G4int i, eIndex = 0;
1218
1219 for( i = 0; i < fIndex; i++)
1220 {
1221 if( energy <= fNuMuEnergy[i]*GeV )
1222 {
1223 eIndex = i;
1224 break;
1225 }
1226 }
1227 if( i >= fIndex ) eIndex = fIndex;
1228 // G4cout<<"eIndex = "<<eIndex<<G4endl;
1229 return eIndex;
1230}
1231
1232/////////////////////////////////////////////////////
1233//
1234// nu_mu QE/Tot ratio for index-1, index linear over energy
1235
1237{
1238 G4double ratio(0.);
1239 // GetMinNuMuEnergy()
1240 if( index <= 0 || energy < fNuMuEnergy[0] ) ratio = 0.;
1241 else if (index >= fIndex) ratio = fNuMuQeTotRat[fIndex-1]*fOnePionEnergy[fIndex-1]*GeV/energy;
1242 else
1243 {
1244 G4double x1 = fNuMuEnergy[index-1]*GeV;
1245 G4double x2 = fNuMuEnergy[index]*GeV;
1246 G4double y1 = fNuMuQeTotRat[index-1];
1247 G4double y2 = fNuMuQeTotRat[index];
1248
1249 if(x1 >= x2) return fNuMuQeTotRat[index];
1250 else
1251 {
1252 G4double angle = (y2-y1)/(x2-x1);
1253 ratio = y1 + (energy-x1)*angle;
1254 }
1255 }
1256 return ratio;
1257}
1258
1259////////////////////////////////////////////////////////
1260
1262{
1263 0.112103, 0.117359, 0.123119, 0.129443, 0.136404,
1264 0.144084, 0.152576, 0.161991, 0.172458, 0.184126,
1265 0.197171, 0.211801, 0.228261, 0.24684, 0.267887,
1266 0.291816, 0.319125, 0.350417, 0.386422, 0.428032,
1267 0.47634, 0.532692, 0.598756, 0.676612, 0.768868,
1268 0.878812, 1.01062, 1.16963, 1.36271, 1.59876,
1269 1.88943, 2.25002, 2.70086, 3.26916, 3.99166,
1270 4.91843, 6.11836, 7.6872, 9.75942, 12.5259,
1271 16.2605, 21.3615, 28.4141, 38.2903, 52.3062,
1272 72.4763, 101.93, 145.6, 211.39, 312.172
1273};
1274
1275////////////////////////////////////////////////////////
1276
1278{
1279 // 1., 1., 1., 1., 1., 1., 1., 1., 1., 1.,
1280 // 1., 1., 1., 1., 1., 1., 1., 1., 1., 1.,
1281 // 1., 1., 1., 0.982311,
1282 0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98,
1283 0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98,
1284 0.97, 0.96, 0.95, 0.93,
1285 0.917794, 0.850239, 0.780412, 0.709339, 0.638134, 0.568165,
1286 0.500236, 0.435528, 0.375015, 0.319157, 0.268463, 0.2232, 0.183284,
1287 0.148627, 0.119008, 0.0940699, 0.0733255, 0.0563819, 0.0427312, 0.0319274,
1288 0.0235026, 0.0170486, 0.0122149, 0.00857825, 0.00594018, 0.00405037
1289};
1290
1291/////////////////////////////////////////////////////
1292//
1293// Return index of one pion array corresponding to the neutrino energy
1294
1296{
1297 G4int i, eIndex = 0;
1298
1299 for( i = 0; i < fOnePionIndex; i++)
1300 {
1301 if( energy <= fOnePionEnergy[i]*GeV )
1302 {
1303 eIndex = i;
1304 break;
1305 }
1306 }
1307 if( i >= fOnePionIndex ) eIndex = fOnePionIndex;
1308 // G4cout<<"eIndex = "<<eIndex<<G4endl;
1309 return eIndex;
1310}
1311
1312/////////////////////////////////////////////////////
1313//
1314// nu_mu 1pi/Tot ratio for index-1, index linear over energy
1315
1317{
1318 G4double ratio(0.);
1319
1320 if( index <= 0 || energy < fOnePionEnergy[0] ) ratio = 0.;
1321 else if ( index >= fOnePionIndex ) ratio = fOnePionProb[fOnePionIndex-1]*fOnePionEnergy[fOnePionIndex-1]*GeV/energy;
1322 else
1323 {
1324 G4double x1 = fOnePionEnergy[index-1]*GeV;
1325 G4double x2 = fOnePionEnergy[index]*GeV;
1326 G4double y1 = fOnePionProb[index-1];
1327 G4double y2 = fOnePionProb[index];
1328
1329 if( x1 >= x2) return fOnePionProb[index];
1330 else
1331 {
1332 G4double angle = (y2-y1)/(x2-x1);
1333 ratio = y1 + (energy-x1)*angle;
1334 }
1335 }
1336 return ratio;
1337}
1338
1339////////////////////////////////////////////////////////////////////////////////////////////////////
1340
1342{
1343
1344 0.275314, 0.293652, 0.31729, 0.33409, 0.351746, 0.365629, 0.380041, 0.400165, 0.437941, 0.479237,
1345 0.504391, 0.537803, 0.588487, 0.627532, 0.686839, 0.791905, 0.878332, 0.987405, 1.08162, 1.16971,
1346 1.2982, 1.40393, 1.49854, 1.64168, 1.7524, 1.87058, 2.02273, 2.15894, 2.3654, 2.55792, 2.73017,
1347 3.03005, 3.40733, 3.88128, 4.53725, 5.16786, 5.73439, 6.53106, 7.43879, 8.36214, 9.39965, 10.296,
1348 11.5735, 13.1801, 15.2052, 17.5414, 19.7178, 22.7462, 25.9026, 29.4955, 33.5867, 39.2516, 46.4716,
1349 53.6065, 63.4668, 73.2147, 85.5593, 99.9854
1350};
1351
1352
1353////////////////////////////////////////////////////////////////////////////////////////////////////
1354
1356{
1357 0.0019357, 0.0189361, 0.0378722, 0.0502758, 0.0662559, 0.0754581, 0.0865008, 0.0987275, 0.124112,
1358 0.153787, 0.18308, 0.213996, 0.245358, 0.274425, 0.301536, 0.326612, 0.338208, 0.337806, 0.335948,
1359 0.328092, 0.313557, 0.304965, 0.292169, 0.28481, 0.269474, 0.254138, 0.247499, 0.236249, 0.221654,
1360 0.205492, 0.198781, 0.182216, 0.162251, 0.142878, 0.128631, 0.116001, 0.108435, 0.0974843, 0.082092,
1361 0.0755204, 0.0703121, 0.0607066, 0.0554278, 0.0480401, 0.0427023, 0.0377123, 0.0323248, 0.0298584,
1362 0.0244296, 0.0218526, 0.019121, 0.016477, 0.0137309, 0.0137963, 0.0110371, 0.00834028, 0.00686127, 0.00538226
1363};
1364
1365
1366//
1367//
1368///////////////////////////
double B(double temperature)
double A(double temperature)
G4double G4Log(G4double x)
Definition: G4Log.hh:226
CLHEP::HepLorentzVector G4LorentzVector
std::vector< G4ReactionProduct * > G4ReactionProductVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector orthogonal() const
double mag() const
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
static double shoot()
void SetNumberOfHoles(G4int valueTot, G4int valueP=0)
Definition: G4Fragment.hh:367
void SetExcEnergyAndMomentum(G4double eexc, const G4LorentzVector &)
Definition: G4Fragment.hh:285
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetLocalEnergyDeposit(G4double aE)
const G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
G4KineticTrackVector * Decay()
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:98
static const G4double fBarMass[4]
void CoherentPion(G4LorentzVector &lvP, G4int pdgP, G4Nucleus &targetNucleus)
static G4double fNuMuQarrayKR[50][51][51]
virtual void ModelDescription(std::ostream &) const
G4GeneratorPrecompoundInterface * fPrecoInterface
G4double FermiMomentum(G4Nucleus &targetNucleus)
static G4double fNuMuXarrayKR[50][51]
G4double NucleonMomentum(G4Nucleus &targetNucleus)
G4int GetOnePionIndex(G4double energy)
static const G4double fResMass[6]
static const G4int fMesPDG[4]
static const G4double fNuMuEnergy[50]
G4ExcitationHandler * fDeExcitation
G4double SampleXkr(G4double energy)
G4ParticleDefinition * theMuonMinus
G4double SampleQkr(G4double energy, G4double xx)
G4double GetNuMuQeTotRat(G4int index, G4double energy)
G4int GetEnergyIndex(G4double energy)
G4double GetXkr(G4int iEnergy, G4double prob)
void MesonDecay(G4LorentzVector &lvX, G4int qX)
G4double GgSampleNM(G4Nucleus &nucl)
G4double GetNuMuOnePionProb(G4int index, G4double energy)
void FinalMeson(G4LorentzVector &lvM, G4int qM, G4int pdgM)
G4double GetQkr(G4int iE, G4int jX, G4double prob)
static const G4double fMesMass[4]
G4PreCompoundModel * fPreCompound
static const G4int fBarPDG[4]
static G4double fNuMuXdistrKR[50][50]
static G4double fNuMuQdistrKR[50][51][50]
virtual G4bool IsApplicable(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
static const G4double fNuMuQeTotRat[50]
G4NeutrinoNucleusModel(const G4String &name="neutrino-nucleus")
G4double FinalMomentum(G4double mI, G4double mF, G4double mP, G4LorentzVector lvX)
G4ParticleDefinition * theMuonPlus
void ClusterDecay(G4LorentzVector &lvX, G4int qX)
void FinalBarion(G4LorentzVector &lvB, G4int qB, G4int pdgB)
static const G4double fOnePionProb[58]
void RecoilDeexcitation(G4Fragment &fragment)
static const G4double fNuMuEnergyLogVector[50]
static const G4double fOnePionEnergy[58]
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4double AtomicMass(const G4double A, const G4double Z) const
Definition: G4Nucleus.cc:254
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment) final
static G4Proton * Proton()
Definition: G4Proton.cc:92
void SetDeExcitation(G4VPreCompoundModel *ptr)
Definition: DoubConv.h:17