Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutrinoNucleusModel Class Referenceabstract

#include <G4NeutrinoNucleusModel.hh>

+ Inheritance diagram for G4NeutrinoNucleusModel:

Public Member Functions

 G4NeutrinoNucleusModel (const G4String &name="neutrino-nucleus")
 
virtual ~G4NeutrinoNucleusModel ()
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double SampleXkr (G4double energy)
 
G4double GetXkr (G4int iEnergy, G4double prob)
 
G4double SampleQkr (G4double energy, G4double xx)
 
G4double GetQkr (G4int iE, G4int jX, G4double prob)
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
 
void ClusterDecay (G4LorentzVector &lvX, G4int qX)
 
void MesonDecay (G4LorentzVector &lvX, G4int qX)
 
void FinalBarion (G4LorentzVector &lvB, G4int qB, G4int pdgB)
 
void RecoilDeexcitation (G4Fragment &fragment)
 
void FinalMeson (G4LorentzVector &lvM, G4int qM, G4int pdgM)
 
void CoherentPion (G4LorentzVector &lvP, G4int pdgP, G4Nucleus &targetNucleus)
 
void SetCutEnergy (G4double ec)
 
G4double GetCutEnergy ()
 
G4double GetNuEnergy ()
 
G4double GetQtransfer ()
 
G4double GetQ2 ()
 
G4double GetXsample ()
 
G4int GetPDGencoding ()
 
G4bool GetCascade ()
 
G4bool GetString ()
 
G4double GetCosTheta ()
 
G4double GetEmu ()
 
G4double GetEx ()
 
G4double GetMuMass ()
 
G4double GetW2 ()
 
G4double GetM1 ()
 
G4double GetMr ()
 
G4double GetTr ()
 
G4double GetDp ()
 
G4bool GetfBreak ()
 
G4bool GetfCascade ()
 
G4bool GetfString ()
 
G4LorentzVector GetLVl ()
 
G4LorentzVector GetLVh ()
 
G4LorentzVector GetLVt ()
 
G4LorentzVector GetLVcpi ()
 
G4double GetMinNuMuEnergy ()
 
G4double ThresholdEnergy (G4double mI, G4double mF, G4double mP)
 
G4double GetQEratioA ()
 
void SetQEratioA (G4double qea)
 
G4double FinalMomentum (G4double mI, G4double mF, G4double mP, G4LorentzVector lvX)
 
G4double FermiMomentum (G4Nucleus &targetNucleus)
 
G4double NucleonMomentum (G4Nucleus &targetNucleus)
 
G4double GetEx (G4int A, G4bool fP)
 
G4double GgSampleNM (G4Nucleus &nucl)
 
G4int GetEnergyIndex (G4double energy)
 
G4double GetNuMuQeTotRat (G4int index, G4double energy)
 
G4int GetOnePionIndex (G4double energy)
 
G4double GetNuMuOnePionProb (G4int index, G4double energy)
 
G4double CalculateQEratioA (G4int Z, G4int A, G4double energy, G4int nepdg)
 
virtual void ModelDescription (std::ostream &) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void InitialiseModel ()
 
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
 
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 

Protected Attributes

G4ParticleDefinitiontheMuonMinus
 
G4ParticleDefinitiontheMuonPlus
 
G4double fSin2tW
 
G4double fCutEnergy
 
G4int fNbin
 
G4int fIndex
 
G4int fEindex
 
G4int fXindex
 
G4int fQindex
 
G4int fOnePionIndex
 
G4int fPDGencoding
 
G4bool fCascade
 
G4bool fString
 
G4bool fProton
 
G4bool f2p2h
 
G4bool fBreak
 
G4double fNuEnergy
 
G4double fQ2
 
G4double fQtransfer
 
G4double fXsample
 
G4double fM1
 
G4double fM2
 
G4double fMt
 
G4double fMu
 
G4double fW2
 
G4double fMpi
 
G4double fW2pi
 
G4double fMinNuEnergy
 
G4double fDp
 
G4double fTr
 
G4double fEmu
 
G4double fEmuPi
 
G4double fEx
 
G4double fMr
 
G4double fCosTheta
 
G4double fCosThetaPi
 
G4double fQEratioA
 
G4LorentzVector fLVh
 
G4LorentzVector fLVl
 
G4LorentzVector fLVt
 
G4LorentzVector fLVcpi
 
G4GeneratorPrecompoundInterfacefPrecoInterface
 
G4PreCompoundModelfPreCompound
 
G4ExcitationHandlerfDeExcitation
 
G4NucleusfRecoil
 
G4int fSecID
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Static Protected Attributes

static const G4int fResNumber = 6
 
static const G4double fResMass [6]
 
static const G4int fClustNumber = 4
 
static const G4double fMesMass [4] = {1260., 980., 770., 139.57}
 
static const G4int fMesPDG [4] = {20213, 9000211, 213, 211}
 
static const G4double fBarMass [4] = {1700., 1600., 1232., 939.57}
 
static const G4int fBarPDG [4] = {12224, 32224, 2224, 2212}
 
static const G4double fNuMuResQ [50][50]
 
static const G4double fNuMuEnergy [50]
 
static const G4double fNuMuQeTotRat [50]
 
static const G4double fOnePionEnergy [58]
 
static const G4double fOnePionProb [58]
 
static const G4double fNuMuEnergyLogVector [50]
 
static G4double fNuMuXarrayKR [50][51] = {{1.0}}
 
static G4double fNuMuXdistrKR [50][50] = {{1.0}}
 
static G4double fNuMuQarrayKR [50][51][51] = {{{1.0}}}
 
static G4double fNuMuQdistrKR [50][51][50] = {{{1.0}}}
 
static const G4double fQEnergy [50]
 
static const G4double fANeMuQEratio [50]
 
static const G4double fNeMuQEratio [50]
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 

Detailed Description

Definition at line 62 of file G4NeutrinoNucleusModel.hh.

Constructor & Destructor Documentation

◆ G4NeutrinoNucleusModel()

G4NeutrinoNucleusModel::G4NeutrinoNucleusModel ( const G4String & name = "neutrino-nucleus")

Definition at line 109 of file G4NeutrinoNucleusModel.cc.

110 : G4HadronicInteraction(name), fSecID(-1)
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
143 fQEratioA = 0.5; // mean value around 1 GeV neutrino beams
144
147
148 // PDG2016: sin^2 theta Weinberg
149
150 fSin2tW = 0.23129; // 0.2312;
151
152 fCutEnergy = 0.; // default value
153
154
155 /*
156 // G4VPreCompoundModel* ptr;
157 // reuse existing pre-compound model as in binary cascade
158
159 fPrecoInterface = new G4GeneratorPrecompoundInterface ;
160
161 if( !fPreCompound )
162 {
163 G4HadronicInteraction* p =
164 G4HadronicInteractionRegistry::Instance()->FindModel("PRECO");
165 G4VPreCompoundModel* fPreCompound = static_cast<G4VPreCompoundModel*>(p);
166
167 if(!fPreCompound)
168 {
169 fPreCompound = new G4PreCompoundModel();
170 }
171 fPrecoInterface->SetDeExcitation(fPreCompound);
172 }
173 fDeExcitation = GetDeExcitation()->GetExcitationHandler();
174 */
175
180
181 fPDGencoding = 0; // unphysical as default
182 fRecoil = nullptr;
183
184 // Creator model ID
186}
CLHEP::HepLorentzVector G4LorentzVector
void SetMinEnergy(G4double anEnergy)
G4HadronicInteraction(const G4String &modelName="HadronicModel")
const G4String & GetModelName() const
void SetMaxEnergy(const G4double anEnergy)
static G4MuonMinus * MuonMinus()
static G4MuonPlus * MuonPlus()
Definition G4MuonPlus.cc:98
G4GeneratorPrecompoundInterface * fPrecoInterface
G4ExcitationHandler * fDeExcitation
G4ParticleDefinition * theMuonMinus
G4PreCompoundModel * fPreCompound
G4ParticleDefinition * theMuonPlus
static G4int GetModelID(const G4int modelIndex)
void SetDeExcitation(G4VPreCompoundModel *ptr)

◆ ~G4NeutrinoNucleusModel()

G4NeutrinoNucleusModel::~G4NeutrinoNucleusModel ( )
virtual

Definition at line 189 of file G4NeutrinoNucleusModel.cc.

190{
192}

Member Function Documentation

◆ ApplyYourself()

◆ CalculateQEratioA()

G4double G4NeutrinoNucleusModel::CalculateQEratioA ( G4int Z,
G4int A,
G4double energy,
G4int nepdg )

Definition at line 1365 of file G4NeutrinoNucleusModel.cc.

1366{
1367 energy /= GeV;
1368 G4double qerata(0.5), rr(0.), x1(0.), x2(0.), y1(0.), y2(0.), aa(0.);
1369 G4int i(0), N(0);
1370
1371 if( A > Z ) N = A-Z;
1372
1373 for( i = 0; i < 50; i++)
1374 {
1375 if( fQEnergy[i] >= energy ) break;
1376 }
1377 if(i <= 0) return 1.;
1378 else if (i >= 49) return 0.;
1379 else
1380 {
1381 x1 = fQEnergy[i-1];
1382 x2 = fQEnergy[i];
1383
1384 if( nepdg == 12 || nepdg == 14 )
1385 {
1386 if( x1 >= x2) return fNeMuQEratio[i];
1387
1388 y1 = fNeMuQEratio[i-1];
1389 y2 = fNeMuQEratio[i];
1390 }
1391 else
1392 {
1393 if( x1 >= x2) return fANeMuQEratio[i];
1394
1395 y1 = fANeMuQEratio[i-1];
1396 y2 = fANeMuQEratio[i];
1397 }
1398 aa = (y2-y1)/(x2-x1);
1399 rr = y1 + (energy-x1)*aa;
1400
1401 if( nepdg == 12 || nepdg == 14 ) qerata = N*rr/( N*rr + A*( 1 - rr ) );
1402 else qerata = Z*rr/( Z*rr + A*( 1 - rr ) );
1403 }
1404 fQEratioA = qerata;
1405
1406 return qerata;
1407}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
static const G4double fQEnergy[50]
static const G4double fNeMuQEratio[50]
static const G4double fANeMuQEratio[50]
#define N
Definition crc32.c:57
G4double energy(const ThreeVector &p, const G4double m)

Referenced by G4ANuElNucleusCcModel::ApplyYourself(), G4ANuElNucleusNcModel::ApplyYourself(), G4ANuMuNucleusCcModel::ApplyYourself(), G4ANuMuNucleusNcModel::ApplyYourself(), G4ANuTauNucleusCcModel::ApplyYourself(), G4ANuTauNucleusNcModel::ApplyYourself(), G4NuElNucleusCcModel::ApplyYourself(), G4NuElNucleusNcModel::ApplyYourself(), G4NuMuNucleusCcModel::ApplyYourself(), G4NuMuNucleusNcModel::ApplyYourself(), G4NuTauNucleusCcModel::ApplyYourself(), and G4NuTauNucleusNcModel::ApplyYourself().

◆ ClusterDecay()

void G4NeutrinoNucleusModel::ClusterDecay ( G4LorentzVector & lvX,
G4int qX )

Definition at line 720 of file G4NeutrinoNucleusModel.cc.

721{
722 G4bool finB = false;
723 G4int pdgB(0), i(0), qM(0), qB(0); // pdgM(0),
724 G4double mM(0.), mB(0.), eM(0.), eB(0.), pM(0.), pB(0.);
725 G4double mm1(0.), mm22(0.), M1(0.), M2(0.), mX(0.);
726
727 mX = lvX.m();
728
731
732 // G4double deltaM = 1.*MeV; // 30.*MeV; // 10.*MeV; // 100.*MeV; // 20.*MeV; //
733 G4double deltaMr[4] = { 0.*MeV, 0.*MeV, 100.*MeV, 0.*MeV};
734
735 G4ThreeVector dir(0.,0.,0.);
736 G4ThreeVector bst(0.,0.,0.);
737 G4LorentzVector lvM(0.,0.,0.,0.);
738 G4LorentzVector lvB(0.,0.,0.,0.);
739
740 for( i = 0; i < fClustNumber; ++i) // check resonance
741 {
742 if( mX >= fBarMass[i] )
743 {
744 pdgB = fBarPDG[i];
745 // mB = G4ParticleTable::GetParticleTable()->FindParticle(pdgB)->GetPDGMass();
746 break;
747 }
748 }
749 if( i == fClustNumber || i == fClustNumber-1 ) // low mass, p || n
750 {
751 if ( qX == 2 || qX == 0) { pdgB = 2212; qB = 1;} // p for 2, 0
752
753 else if( qX == 1 || qX == -1) { pdgB = 2112; qB = 0;} // n for 1, -1
754
755 return FinalBarion( lvX, qB, pdgB);
756 }
757 else if( mX < fBarMass[i] + deltaMr[i] || mX < mN + mPi )
758 {
759 finB = true; // final barion -> out
760
761 if ( qX == 1 && pdgB != 2212) pdgB = pdgB - 10;
762 else if( qX == 0 && pdgB != 2212) pdgB = pdgB - 110;
763 else if( qX == 0 && pdgB == 2212) pdgB = pdgB - 100;
764
765 if( finB ) return FinalBarion( lvX, qX, pdgB ); // out
766 }
767 // no barion resonance, try 1->2 decay in COM frame
768
769 // try meson mass
770
771 mm1 = mPi + 1.*MeV; // pi+
772 mm22 = mX - mN; // mX-n
773
774 if( mm22 <= mm1 ) // out with p or n
775 {
776 if( qX == 2 || qX == 0) { pdgB = 2212; qB = 1;} // p
777 else if( qX == 1 || qX == -1) { pdgB = 2112; qB = 0;} // n
778
779 return FinalBarion(lvX, qB, pdgB);
780 }
781 else // try decay -> meson(cluster) + barion(cluster)
782 {
783 // G4double sigmaM = 50.*MeV; // 100.*MeV; // 200.*MeV; // 400.*MeV; // 800.*MeV; //
784 G4double rand = G4UniformRand();
785
786 // mM = mm1*mm22/( mm1 + rand*(mm22 - mm1) );
787 // mM = mm1*mm22/sqrt( mm1*mm1 + rand*(mm22*mm22 - mm1*mm1) );
788 // mM = -sigmaM*log( (1.- rand)*exp(-mm22/sigmaM) + rand*exp(-mm1/sigmaM) );
789 mM = mm1 + rand*(mm22-mm1);
790
791
792 for( i = 0; i < fClustNumber; ++i)
793 {
794 if( mM >= fMesMass[i] )
795 {
796 // pdgM = fMesPDG[i];
797 // mM = G4ParticleTable::GetParticleTable()->FindParticle(pdgM)->GetPDGMass();
798 break;
799 }
800 }
801 M1 = G4ParticleTable::GetParticleTable()->FindParticle(2112)->GetPDGMass()+2.*MeV; // n
802 M2 = mX - mM;
803
804 if( M2 <= M1 ) //
805 {
806 if ( qX == 2 || qX == 0) { pdgB = 2212; qB = 1;} // p
807 else if( qX == 1 || qX == -1) { pdgB = 2112; qB = 0;} // n
808
809 return FinalBarion(lvX, qB, pdgB);
810 }
811 mB = M1 + G4UniformRand()*(M2-M1);
812 // mB = -sigmaM*log( (1.- rand)*exp(-M2/sigmaM) + rand*exp(-M1/sigmaM) );
813
814 bst = lvX.boostVector();
815
816 // dir = G4RandomDirection(); // ???
817 // dir = G4ThreeVector(0.,0.,1.);
818 dir = bst.orthogonal().unit(); // ??
819 // G4double cost = exp(-G4UniformRand());
820 // G4double sint = sqrt((1.-cost)*(1.+cost));
821 // G4double phi = twopi*G4UniformRand();
822 // dir = G4ThreeVector(sint*cos(phi), sint*sin(phi), cost);
823
824 eM = 0.5*(mX*mX + mM*mM - mB*mB)/mX;
825 pM = sqrt(eM*eM - mM*mM);
826 lvM = G4LorentzVector( pM*dir, eM);
827 lvM.boost(bst);
828
829 eB = 0.5*(mX*mX + mB*mB - mM*mM)/mX;
830 pB = sqrt(eB*eB - mB*mB);
831 lvB = G4LorentzVector(-pB*dir, eB);
832 lvB.boost(bst);
833
834 // G4cout<<mM<<"/"<<mB<<", ";
835
836 // charge exchange
837
838 if ( qX == 2 ) { qM = 1; qB = 1;}
839 else if( qX == 1 ) { qM = 0; qB = 1;}
840 else if( qX == 0 ) { qM = 0; qB = 0;}
841 else if( qX == -1 ) { qM = -1; qB = 0;}
842
843 // if ( qM == 0 ) pdgM = pdgM - 100;
844 // else if( qM == -1 ) pdgM = -pdgM;
845
846 MesonDecay( lvM, qM); // pdgM ); //
847
848 // else
849 ClusterDecay( lvB, qB ); // continue
850 }
851}
bool G4bool
Definition G4Types.hh:86
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
Hep3Vector orthogonal() const
Hep3Vector boostVector() const
static const G4double fBarMass[4]
void MesonDecay(G4LorentzVector &lvX, G4int qX)
static const G4double fMesMass[4]
static const G4int fBarPDG[4]
void ClusterDecay(G4LorentzVector &lvX, G4int qX)
void FinalBarion(G4LorentzVector &lvB, G4int qB, G4int pdgB)
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()

Referenced by G4ANuElNucleusCcModel::ApplyYourself(), G4ANuElNucleusNcModel::ApplyYourself(), G4ANuMuNucleusCcModel::ApplyYourself(), G4ANuMuNucleusNcModel::ApplyYourself(), G4ANuTauNucleusCcModel::ApplyYourself(), G4ANuTauNucleusNcModel::ApplyYourself(), G4NuElNucleusCcModel::ApplyYourself(), G4NuElNucleusNcModel::ApplyYourself(), G4NuMuNucleusCcModel::ApplyYourself(), G4NuMuNucleusNcModel::ApplyYourself(), G4NuTauNucleusCcModel::ApplyYourself(), G4NuTauNucleusNcModel::ApplyYourself(), and ClusterDecay().

◆ CoherentPion()

void G4NeutrinoNucleusModel::CoherentPion ( G4LorentzVector & lvP,
G4int pdgP,
G4Nucleus & targetNucleus )

Definition at line 595 of file G4NeutrinoNucleusModel.cc.

596{
597 G4int A(0), Z(0), pdg = pdgP;
598 fLVcpi = G4LorentzVector(0.,0.,0.,0.);
599
600 G4double rM(0.), mN(938.), det(0.), det2(0.);
601 G4double mI(0.);
602 mN = G4ParticleTable::GetParticleTable()->FindParticle(2212)->GetPDGMass(); // *0.85; // *0.9; //
603
604 // mN = 1.*139.57 + G4UniformRand()*(938. - 1.*139.57);
605
606 G4ThreeVector vN = lvP.boostVector(), bst(0.,0.,0.);
607 // G4double gN = lvP.e()/lvP.m();
608 // G4LorentzVector lvNu(vN*gN*mN, mN*gN);
609 G4LorentzVector lvNu(0.,0.,0., mN); // lvNu(bst, mN);
610 lvP.boost(-vN); // 9-3-20
611 lvP = lvP - lvNu; // 9-3-20 already 1pi
612 lvP.boost(vN); // 9-3-20
613 lvNu.boost(vN); // 9-3-20
614
615 // G4cout<<vN-lvP.boostVector()<<", ";
616
617 Z = targetNucleus.GetZ_asInt();
618 A = targetNucleus.GetA_asInt();
619 rM = targetNucleus.AtomicMass(A,Z); //->AtomicMass(); //
620
621 // G4cout<<rM<<", ";
622 // G4cout<<A<<", ";
623
624 if( A == 1 )
625 {
626 bst = vN; // lvNu.boostVector(); // 9-3-20
627 // mI = 0.; // 9-3-20
628 rM = mN;
629 }
630 else
631 {
632 G4Nucleus targ(A-1,Z);
633 mI = targ.AtomicMass(A-1,Z);
634 G4LorentzVector lvTar(0.,0.,0.,mI);
635 lvNu = lvNu + lvTar;
636 bst = lvNu.boostVector();
637 // bst = fLVt.boostVector(); // to recoil rest frame
638 // G4cout<<fLVt<<" "<<bst<<G4endl;
639 }
640 lvP.boost(-bst); // 9-3-20
642 G4double eX = lvP.e();
643 G4double mX = lvP.m();
644 // G4cout<<mX-fMr<<", ";
645 G4ThreeVector dX = (lvP.vect()).unit();
646 // G4cout<<dX<<", ";
647 G4double pX = sqrt(eX*eX-mX*mX);
648 // G4cout<<pX<<", ";
649 G4double sumE = eX + rM;
650 G4double B = sumE*sumE + rM*rM - fMr*fMr - pX*pX;
651 G4double a = 4.*(sumE*sumE - pX*pX);
652 G4double b = -4.*B*pX;
653 G4double c = 4.*sumE*sumE*rM*rM - B*B;
654 det2 = b*b-4.*a*c;
655 if(det2 > 0.) det = sqrt(det2);
656 G4double dP = 0.5*(-b - det )/a;
657
658 // dP = FinalMomentum( mI, rM, fMr, lvP);
659 dP = FinalMomentum( rM, rM, fMr, lvP); // 9-3-20
660
661 // G4cout<<dP<<", ";
662 pX -= dP;
663 if( pX < 0. ) pX = 0.;
664
665 eX = sqrt( dP*dP + fMr*fMr );
666 G4LorentzVector lvN( dP*dX, eX );
667
668 if( A >= 1 ) lvN.boost(bst); // 9-3-20 back to lab
669
670 fLVcpi = lvN;
671
673 G4DynamicParticle* dp2 = new G4DynamicParticle( pd2, lvN);
674 theParticleChange.AddSecondary( dp2, fSecID ); // coherent pion
675
676 // recoil nucleus
677
678 G4double eRecoil = sqrt( rM*rM + pX*pX );
679 G4ThreeVector vRecoil(pX*dX);
680 G4LorentzVector lvTarg1( vRecoil, eRecoil);
681 lvTarg1.boost(bst);
682
683 const G4LorentzVector lvTarg = lvTarg1;
684
685 if( A > 1 ) // recoil target nucleus*
686 {
688 G4double exE = fLVt.m() - grM;
689
690 if( exE < 5.*MeV ) exE = 5.*MeV + G4UniformRand()*10.*MeV; // vmg???
691
692 const G4LorentzVector in4v( G4ThreeVector( 0., 0., 0.), grM );
693 G4Fragment fragment( A, Z, in4v); // lvTarg );
694 fragment.SetNumberOfHoles(1);
695 fragment.SetExcEnergyAndMomentum( exE, lvTarg );
696
697 RecoilDeexcitation(fragment);
698 }
699 else // recoil target proton
700 {
701 G4double eTkin = eRecoil - rM;
702 G4double eTh = 0.01*MeV; // 10.*MeV;
703
704 if( eTkin > eTh )
705 {
708 }
710 }
711 return;
712}
G4double B(G4double temperature)
CLHEP::Hep3Vector G4ThreeVector
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetLocalEnergyDeposit(G4double aE)
G4double FinalMomentum(G4double mI, G4double mF, G4double mP, G4LorentzVector lvX)
void RecoilDeexcitation(G4Fragment &fragment)
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105
G4double AtomicMass(const G4double A, const G4double Z, const G4int numberOfLambdas=0) const
Definition G4Nucleus.cc:369
static G4Proton * Proton()
Definition G4Proton.cc:90

Referenced by G4ANuElNucleusCcModel::ApplyYourself(), G4ANuElNucleusNcModel::ApplyYourself(), G4ANuMuNucleusCcModel::ApplyYourself(), G4ANuMuNucleusNcModel::ApplyYourself(), G4ANuTauNucleusCcModel::ApplyYourself(), G4ANuTauNucleusNcModel::ApplyYourself(), G4NuElNucleusCcModel::ApplyYourself(), G4NuElNucleusNcModel::ApplyYourself(), G4NuMuNucleusCcModel::ApplyYourself(), G4NuMuNucleusNcModel::ApplyYourself(), G4NuTauNucleusCcModel::ApplyYourself(), and G4NuTauNucleusNcModel::ApplyYourself().

◆ FermiMomentum()

G4double G4NeutrinoNucleusModel::FermiMomentum ( G4Nucleus & targetNucleus)

Definition at line 1057 of file G4NeutrinoNucleusModel.cc.

1058{
1059 G4int Z = targetNucleus.GetZ_asInt();
1060 G4int A = targetNucleus.GetA_asInt();
1061
1062 G4double kF(250.*MeV);
1063 G4double kp = 365.*MeV;
1064 G4double kn = 231.*MeV;
1065 G4double t1 = 0.479;
1066 G4double t2 = 0.526;
1067 G4double ZpA = G4double(Z)/G4double(A);
1068 G4double NpA = 1. - ZpA;
1069
1070 if ( Z == 1 && A == 1 ) { kF = 0.; } // hydrogen ???
1071 else if ( Z == 1 && A == 2 ) { kF = 87.*MeV; }
1072 else if ( Z == 2 && A == 3 ) { kF = 134.*MeV; }
1073 else if ( Z == 6 && A == 12 ) { kF = 221.*MeV; }
1074 else if ( Z == 14 && A == 28 ) { kF = 239.*MeV; }
1075 else if ( Z == 26 && A == 56 ) { kF = 257.*MeV; }
1076 else if ( Z == 82 && A == 208 ) { kF = 265.*MeV; }
1077 else
1078 {
1079 kF = kp*ZpA*( 1 - pow( G4double(A), -t1 ) ) + kn*NpA*( 1 - pow( G4double(A), -t2 ) );
1080 }
1081 return kF;
1082}

Referenced by GgSampleNM(), and NucleonMomentum().

◆ FinalBarion()

void G4NeutrinoNucleusModel::FinalBarion ( G4LorentzVector & lvB,
G4int qB,
G4int pdgB )

Definition at line 449 of file G4NeutrinoNucleusModel.cc.

450{
451 G4int A(0), Z(0), pdg = pdgB;
452 // G4bool FiNucleon(false);
453
454 // if ( qB == 1 ) pdg = pdgB - 10;
455 // else if ( qB == 0 ) pdg = pdgB - 110;
456 // else if ( qB == -1 ) pdg = pdgB - 1110;
457
458 if( pdg == 2212 || pdg == 2112)
459 {
461 // FiNucleon = true;
462 }
463 else fMr = lvB.m();
464
466 lvB.boost(-bst); // in fLVt rest system
467
468 G4double eX = lvB.e();
469 G4double det(0.), det2(0.), rM(0.), mX = lvB.m();
470 G4ThreeVector dX = (lvB.vect()).unit();
471 G4double pX = sqrt(eX*eX-mX*mX);
472
473 if( fRecoil )
474 {
475 Z = fRecoil->GetZ_asInt();
476 A = fRecoil->GetA_asInt();
477 rM = fRecoil->AtomicMass(A,Z); //->AtomicMass(); //
478 rM = fLVt.m();
479 }
480 else // A=0 nu+p
481 {
482 A = 0;
483 Z = 1;
484 rM = electron_mass_c2;
485 }
486 // G4cout<<A<<", ";
487
488 G4double sumE = eX + rM;
489 G4double B = sumE*sumE + rM*rM - fMr*fMr - pX*pX;
490 G4double a = 4.*(sumE*sumE - pX*pX);
491 G4double b = -4.*B*pX;
492 G4double c = 4.*sumE*sumE*rM*rM - B*B;
493 det2 = b*b-4.*a*c;
494 if( det2 <= 0. ) det = 0.;
495 else det = sqrt(det2);
496 G4double dP = 0.5*(-b - det )/a;
497
498 fDp = dP;
499
500 pX -= dP;
501
502 if(pX < 0.) pX = 0.;
503
504 // if( A == 0 ) G4cout<<pX/MeV<<", ";
505
506 eX = sqrt( pX*pX + fMr*fMr );
507 G4LorentzVector lvN( pX*dX, eX );
508 lvN.boost(bst); // back to lab
509
510 if( pdg == 2212 || pdg == 2112) // nucleons mX >= fMr, dP >= 0
511 {
513 G4DynamicParticle* dp2 = new G4DynamicParticle( pd2, lvN);
515
516 }
517 else // delta resonances
518 {
520 G4KineticTrack ddkt( rePart, 0., G4ThreeVector( 0., 0., 0.), lvN);
521 G4KineticTrackVector* ddktv = ddkt.Decay();
522
524
525 for( unsigned int i = 0; i < ddktv->size(); i++ ) // add products to partchange
526 {
527 G4DynamicParticle * aNew =
528 new G4DynamicParticle( ddktv->operator[](i)->GetDefinition(),
529 ddktv->operator[](i)->Get4Momentum() );
530
531 // G4cout<<" "<<i<<", "<<aNew->GetDefinition()->GetParticleName()<<", "<<aNew->Get4Momentum()<<G4endl;
532
534 delete ddktv->operator[](i);
535 }
536 delete ddktv;
537 }
538 // recoil nucleus
539
540 G4double eRecoil = sqrt( rM*rM + dP*dP );
541 fTr = eRecoil - rM;
542 G4ThreeVector vRecoil(dP*dX);
543 // dP += G4UniformRand()*10.*MeV;
544 G4LorentzVector rec4v(vRecoil, 0.);
545 rec4v.boost(bst); // back to lab
546 fLVt += rec4v;
547 const G4LorentzVector lvTarg = fLVt; // (vRecoil, eRecoil);
548
549
550 if( fRecoil ) // proton*?
551 {
553 G4double exE = fLVt.m() - grM;
554 if( exE < 5.*MeV ) exE = 5.*MeV + G4UniformRand()*10.*MeV;
555
556 const G4LorentzVector in4v( G4ThreeVector( 0., 0., 0.), grM );
557 G4Fragment fragment( A, Z, in4v); // lvTarg );
558 fragment.SetNumberOfHoles(1);
559 fragment.SetExcEnergyAndMomentum( exE, lvTarg );
560
561 RecoilDeexcitation(fragment);
562 }
563 else // momentum?
564 {
566 }
567}
ParticleList decay(Cluster *const c)
Carries out a cluster decay.

Referenced by G4ANuElNucleusCcModel::ApplyYourself(), G4ANuElNucleusNcModel::ApplyYourself(), G4ANuMuNucleusCcModel::ApplyYourself(), G4ANuMuNucleusNcModel::ApplyYourself(), G4ANuTauNucleusCcModel::ApplyYourself(), G4ANuTauNucleusNcModel::ApplyYourself(), G4NuElNucleusCcModel::ApplyYourself(), G4NuElNucleusNcModel::ApplyYourself(), G4NuMuNucleusCcModel::ApplyYourself(), G4NuMuNucleusNcModel::ApplyYourself(), G4NuTauNucleusCcModel::ApplyYourself(), G4NuTauNucleusNcModel::ApplyYourself(), and ClusterDecay().

◆ FinalMeson()

void G4NeutrinoNucleusModel::FinalMeson ( G4LorentzVector & lvM,
G4int qM,
G4int pdgM )

Definition at line 409 of file G4NeutrinoNucleusModel.cc.

410{
411 G4int pdg = pdgM;
412 // if ( qM == 0 ) pdg = pdgM - 100;
413 // else if ( qM == -1 ) pdg = -pdgM;
414
415 if( pdg == 211 || pdg == -211 || pdg == 111) // pions
416 {
418 G4DynamicParticle* dp2 = new G4DynamicParticle( pd2, lvM);
420 }
421 else // meson resonances
422 {
424 FindParticle(pdg);
425 G4KineticTrack ddkt( rePart, 0., G4ThreeVector(0.,0.,0.), lvM);
426 G4KineticTrackVector* ddktv = ddkt.Decay();
427
429
430 for( unsigned int i = 0; i < ddktv->size(); i++ ) // add products to partchange
431 {
432 G4DynamicParticle * aNew =
433 new G4DynamicParticle( ddktv->operator[](i)->GetDefinition(),
434 ddktv->operator[](i)->Get4Momentum());
435
436 // G4cout<<" "<<i<<", "<<aNew->GetDefinition()->GetParticleName()<<", "<<aNew->Get4Momentum()<<G4endl;
437
439 delete ddktv->operator[](i);
440 }
441 delete ddktv;
442 }
443}

Referenced by MesonDecay().

◆ FinalMomentum()

G4double G4NeutrinoNucleusModel::FinalMomentum ( G4double mI,
G4double mF,
G4double mP,
G4LorentzVector lvX )

Definition at line 1026 of file G4NeutrinoNucleusModel.cc.

1027{
1028 G4double result(0.), delta(0.);
1029 // G4double mI2 = mI*mI;
1030 G4double mF2 = mF*mF;
1031 G4double mP2 = mP*mP;
1032 G4double eX = lvX.e();
1033 // G4double mX = lvX.m();
1034 G4double pX = lvX.vect().mag();
1035 G4double pX2 = pX*pX;
1036 G4double sI = eX + mI;
1037 G4double sI2 = sI*sI;
1038 G4double B = sI2 - mF2 -pX2 + mP2;
1039 G4double B2 = B*B;
1040 G4double a = 4.*(sI2-pX2);
1041 G4double b = -4.*B*pX;
1042 G4double c = 4.*sI2*mP2 - B2;
1043 G4double delta2 = b*b -4.*a*c;
1044
1045 if( delta2 >= 0. ) delta = sqrt(delta2);
1046
1047 result = 0.5*(-b-delta)/a;
1048 // result = 0.5*(-b+delta)/a;
1049
1050 return result;
1051}
double mag() const

Referenced by CoherentPion().

◆ GetCascade()

G4bool G4NeutrinoNucleusModel::GetCascade ( )
inline

Definition at line 107 of file G4NeutrinoNucleusModel.hh.

107{return fCascade;};

◆ GetCosTheta()

G4double G4NeutrinoNucleusModel::GetCosTheta ( )
inline

Definition at line 110 of file G4NeutrinoNucleusModel.hh.

110{return fCosTheta;};

◆ GetCutEnergy()

G4double G4NeutrinoNucleusModel::GetCutEnergy ( )
inline

Definition at line 99 of file G4NeutrinoNucleusModel.hh.

99{return fCutEnergy;};

◆ GetDp()

G4double G4NeutrinoNucleusModel::GetDp ( )
inline

Definition at line 118 of file G4NeutrinoNucleusModel.hh.

118{return fDp;};

◆ GetEmu()

G4double G4NeutrinoNucleusModel::GetEmu ( )
inline

Definition at line 111 of file G4NeutrinoNucleusModel.hh.

111{return fEmu;};

◆ GetEnergyIndex()

G4int G4NeutrinoNucleusModel::GetEnergyIndex ( G4double energy)

Definition at line 1213 of file G4NeutrinoNucleusModel.cc.

1214{
1215 G4int i, eIndex = 0;
1216
1217 for( i = 0; i < fIndex; i++)
1218 {
1219 if( energy <= fNuMuEnergy[i]*GeV )
1220 {
1221 eIndex = i;
1222 break;
1223 }
1224 }
1225 if( i >= fIndex ) eIndex = fIndex;
1226 // G4cout<<"eIndex = "<<eIndex<<G4endl;
1227 return eIndex;
1228}
static const G4double fNuMuEnergy[50]

◆ GetEx() [1/2]

◆ GetEx() [2/2]

G4double G4NeutrinoNucleusModel::GetEx ( G4int A,
G4bool fP )

Definition at line 1115 of file G4NeutrinoNucleusModel.cc.

1116{
1117 G4double eX(10.*MeV), a1(0.), a2(0.), e1(0.), e2(0.), aa = G4double(A);
1118 G4int i(0);
1119 const G4int maxBin = 12;
1120
1121 G4double refA[maxBin] = { 2., 6., 12., 16., 27., 28., 40., 50., 56., 58., 197., 208. };
1122
1123 G4double pEx[maxBin] = { 0., 12.2, 10.1, 10.9, 21.6, 12.4, 17.8, 17., 19., 16.8, 19.5, 14.7 };
1124
1125 G4double nEx[maxBin] = { 0., 12.2, 10., 10.2, 21.6, 12.4, 21.8, 17., 19., 16.8, 19.5, 16.9 };
1126
1127 G4DataVector dE(12,0.);
1128
1129 if(fP) for( i = 0; i < maxBin; ++i ) dE[i] = pEx[i];
1130 else dE[i] = nEx[i];
1131
1132 for( i = 0; i < maxBin; ++i )
1133 {
1134 if( aa <= refA[i] ) break;
1135 }
1136 if( i >= maxBin ) eX = dE[maxBin-1];
1137 else if( i <= 0 ) eX = dE[0];
1138 else
1139 {
1140 a1 = refA[i-1];
1141 a2 = refA[i];
1142 e1 = dE[i-1];
1143 e2 = dE[i];
1144 if (a1 == a2 || e1 == e2 ) eX = dE[i];
1145 else eX = e1 + (e2-e1)*(aa-a1)/(a2-a1);
1146 }
1147 return eX;
1148}

◆ GetfBreak()

G4bool G4NeutrinoNucleusModel::GetfBreak ( )
inline

Definition at line 120 of file G4NeutrinoNucleusModel.hh.

120{return fBreak;};

◆ GetfCascade()

G4bool G4NeutrinoNucleusModel::GetfCascade ( )
inline

Definition at line 121 of file G4NeutrinoNucleusModel.hh.

121{return fCascade;};

◆ GetfString()

G4bool G4NeutrinoNucleusModel::GetfString ( )
inline

Definition at line 122 of file G4NeutrinoNucleusModel.hh.

122{return fString;};

◆ GetLVcpi()

G4LorentzVector G4NeutrinoNucleusModel::GetLVcpi ( )
inline

Definition at line 127 of file G4NeutrinoNucleusModel.hh.

127{return fLVcpi;};

◆ GetLVh()

G4LorentzVector G4NeutrinoNucleusModel::GetLVh ( )
inline

Definition at line 125 of file G4NeutrinoNucleusModel.hh.

125{return fLVh;};

◆ GetLVl()

G4LorentzVector G4NeutrinoNucleusModel::GetLVl ( )
inline

Definition at line 124 of file G4NeutrinoNucleusModel.hh.

124{return fLVl;};

◆ GetLVt()

G4LorentzVector G4NeutrinoNucleusModel::GetLVt ( )
inline

Definition at line 126 of file G4NeutrinoNucleusModel.hh.

126{return fLVt;};

◆ GetM1()

G4double G4NeutrinoNucleusModel::GetM1 ( )
inline

Definition at line 115 of file G4NeutrinoNucleusModel.hh.

115{return fM1;};

◆ GetMinNuMuEnergy()

G4double G4NeutrinoNucleusModel::GetMinNuMuEnergy ( )
inline

Definition at line 129 of file G4NeutrinoNucleusModel.hh.

129{ return fMu + 0.5*fMu*fMu/fM1 + 4.*CLHEP::MeV; }; // kinematics + accuracy for sqrts

Referenced by G4NeutrinoNucleusModel().

◆ GetMr()

G4double G4NeutrinoNucleusModel::GetMr ( )
inline

Definition at line 116 of file G4NeutrinoNucleusModel.hh.

116{return fMr;};

◆ GetMuMass()

G4double G4NeutrinoNucleusModel::GetMuMass ( )
inline

Definition at line 113 of file G4NeutrinoNucleusModel.hh.

113{return fMu;};

◆ GetNuEnergy()

G4double G4NeutrinoNucleusModel::GetNuEnergy ( )
inline

Definition at line 101 of file G4NeutrinoNucleusModel.hh.

101{return fNuEnergy;};

◆ GetNuMuOnePionProb()

G4double G4NeutrinoNucleusModel::GetNuMuOnePionProb ( G4int index,
G4double energy )

Definition at line 1314 of file G4NeutrinoNucleusModel.cc.

1315{
1316 G4double ratio(0.);
1317
1318 if( index <= 0 || energy < fOnePionEnergy[0] ) ratio = 0.;
1319 else if ( index >= fOnePionIndex ) ratio = fOnePionProb[fOnePionIndex-1]*fOnePionEnergy[fOnePionIndex-1]*GeV/energy;
1320 else
1321 {
1322 G4double x1 = fOnePionEnergy[index-1]*GeV;
1323 G4double x2 = fOnePionEnergy[index]*GeV;
1324 G4double y1 = fOnePionProb[index-1];
1325 G4double y2 = fOnePionProb[index];
1326
1327 if( x1 >= x2) return fOnePionProb[index];
1328 else
1329 {
1330 G4double angle = (y2-y1)/(x2-x1);
1331 ratio = y1 + (energy-x1)*angle;
1332 }
1333 }
1334 return ratio;
1335}
static const G4double fOnePionProb[58]
static const G4double fOnePionEnergy[58]

Referenced by G4ANuElNucleusCcModel::ApplyYourself(), G4ANuElNucleusNcModel::ApplyYourself(), G4ANuMuNucleusCcModel::ApplyYourself(), G4ANuMuNucleusNcModel::ApplyYourself(), G4ANuTauNucleusCcModel::ApplyYourself(), G4ANuTauNucleusNcModel::ApplyYourself(), G4NuElNucleusCcModel::ApplyYourself(), G4NuElNucleusNcModel::ApplyYourself(), G4NuMuNucleusCcModel::ApplyYourself(), G4NuMuNucleusNcModel::ApplyYourself(), G4NuTauNucleusCcModel::ApplyYourself(), and G4NuTauNucleusNcModel::ApplyYourself().

◆ GetNuMuQeTotRat()

G4double G4NeutrinoNucleusModel::GetNuMuQeTotRat ( G4int index,
G4double energy )

Definition at line 1234 of file G4NeutrinoNucleusModel.cc.

1235{
1236 G4double ratio(0.);
1237 // GetMinNuMuEnergy()
1238 if( index <= 0 || energy < fNuMuEnergy[0] ) ratio = 0.;
1239 else if (index >= fIndex) ratio = fNuMuQeTotRat[fIndex-1]*fOnePionEnergy[fIndex-1]*GeV/energy;
1240 else
1241 {
1242 G4double x1 = fNuMuEnergy[index-1]*GeV;
1243 G4double x2 = fNuMuEnergy[index]*GeV;
1244 G4double y1 = fNuMuQeTotRat[index-1];
1245 G4double y2 = fNuMuQeTotRat[index];
1246
1247 if(x1 >= x2) return fNuMuQeTotRat[index];
1248 else
1249 {
1250 G4double angle = (y2-y1)/(x2-x1);
1251 ratio = y1 + (energy-x1)*angle;
1252 }
1253 }
1254 return ratio;
1255}
static const G4double fNuMuQeTotRat[50]

◆ GetOnePionIndex()

G4int G4NeutrinoNucleusModel::GetOnePionIndex ( G4double energy)

◆ GetPDGencoding()

G4int G4NeutrinoNucleusModel::GetPDGencoding ( )
inline

Definition at line 106 of file G4NeutrinoNucleusModel.hh.

106{return fPDGencoding;};

◆ GetQ2()

G4double G4NeutrinoNucleusModel::GetQ2 ( )
inline

Definition at line 103 of file G4NeutrinoNucleusModel.hh.

103{return fQ2;};

◆ GetQEratioA()

G4double G4NeutrinoNucleusModel::GetQEratioA ( )
inline

Definition at line 136 of file G4NeutrinoNucleusModel.hh.

136{ return fQEratioA; };

◆ GetQkr()

G4double G4NeutrinoNucleusModel::GetQkr ( G4int iE,
G4int jX,
G4double prob )

Definition at line 366 of file G4NeutrinoNucleusModel.cc.

367{
368 G4int i(0), nBin=50;
369 G4double qq(0.);
370
371 for( i = 0; i < nBin; ++i )
372 {
373 if( prob <= fNuMuQdistrKR[iE][jX][i] )
374 break;
375 }
376 if(i <= 0 ) // Q-edge
377 {
378 fQindex = 0;
379 qq = fNuMuQarrayKR[iE][jX][0];
380 }
381 if ( i >= nBin )
382 {
383 fQindex = nBin;
384 qq = fNuMuQarrayKR[iE][jX][nBin];
385 }
386 else
387 {
388 fQindex = i;
389 G4double q1 = fNuMuQarrayKR[iE][jX][i];
390 G4double q2 = fNuMuQarrayKR[iE][jX][i+1];
391
392 G4double p1 = 0.;
393
394 if( i > 0 ) p1 = fNuMuQdistrKR[iE][jX][i-1];
395
396 G4double p2 = fNuMuQdistrKR[iE][jX][i];
397
398 if( p2 <= p1 ) qq = q1 + G4UniformRand()*(q2-q1);
399 else qq = q1 + (prob-p1)*(q2-q1)/(p2-p1);
400 }
401 return qq;
402}
static G4double fNuMuQarrayKR[50][51][51]
static G4double fNuMuQdistrKR[50][51][50]

Referenced by SampleQkr().

◆ GetQtransfer()

G4double G4NeutrinoNucleusModel::GetQtransfer ( )
inline

Definition at line 102 of file G4NeutrinoNucleusModel.hh.

102{return fQtransfer;};

◆ GetString()

G4bool G4NeutrinoNucleusModel::GetString ( )
inline

Definition at line 108 of file G4NeutrinoNucleusModel.hh.

108{return fString;};

◆ GetTr()

G4double G4NeutrinoNucleusModel::GetTr ( )
inline

Definition at line 117 of file G4NeutrinoNucleusModel.hh.

117{return fTr;};

◆ GetW2()

G4double G4NeutrinoNucleusModel::GetW2 ( )
inline

Definition at line 114 of file G4NeutrinoNucleusModel.hh.

114{return fW2;};

◆ GetXkr()

G4double G4NeutrinoNucleusModel::GetXkr ( G4int iEnergy,
G4double prob )

Definition at line 264 of file G4NeutrinoNucleusModel.cc.

265{
266 G4int i(0), nBin=50;
267 G4double xx(0.);
268
269 for( i = 0; i < nBin; ++i )
270 {
271 if( prob <= fNuMuXdistrKR[iEnergy][i] )
272 break;
273 }
274 if(i <= 0 ) // X-edge
275 {
276 fXindex = 0;
277 xx = fNuMuXarrayKR[iEnergy][0];
278 }
279 if ( i >= nBin )
280 {
281 fXindex = nBin;
282 xx = fNuMuXarrayKR[iEnergy][nBin];
283 }
284 else
285 {
286 fXindex = i;
287 G4double x1 = fNuMuXarrayKR[iEnergy][i];
288 G4double x2 = fNuMuXarrayKR[iEnergy][i+1];
289
290 G4double p1 = 0.;
291
292 if( i > 0 ) p1 = fNuMuXdistrKR[iEnergy][i-1];
293
294 G4double p2 = fNuMuXdistrKR[iEnergy][i];
295
296 if( p2 <= p1 ) xx = x1 + G4UniformRand()*(x2-x1);
297 else xx = x1 + (prob-p1)*(x2-x1)/(p2-p1);
298 }
299 return xx;
300}
static G4double fNuMuXarrayKR[50][51]
static G4double fNuMuXdistrKR[50][50]

Referenced by SampleXkr().

◆ GetXsample()

G4double G4NeutrinoNucleusModel::GetXsample ( )
inline

Definition at line 104 of file G4NeutrinoNucleusModel.hh.

104{return fXsample;};

◆ GgSampleNM()

G4double G4NeutrinoNucleusModel::GgSampleNM ( G4Nucleus & nucl)

Definition at line 1156 of file G4NeutrinoNucleusModel.cc.

1157{
1158 f2p2h = false;
1159 G4double /* distr(0.), tail(0.), */ shift(1.), xx(1.), mom(0.), th(0.1);
1160 G4double kF = FermiMomentum( nucl);
1161 G4double momMax = 2.*kF; // 1.*GeV; // 1.*GeV; //
1162 G4double aa = 5.5;
1163 G4double ll = 6.0; // 6.5; //
1164
1165 G4int A = nucl.GetA_asInt();
1166
1167 if( A <= 12) th = 0.1;
1168 else
1169 {
1170 // th = 0.1/(1.+log(G4double(A)/12.));
1171 th = 1.2/( G4double(A) + 1.35*log(G4double(A)/12.) );
1172 }
1173 shift = 0.99; // 0.95; //
1174 xx = mom/shift/kF;
1175
1176 G4double rr = G4UniformRand();
1177
1178 if( rr > th )
1179 {
1180 aa = 5.5;
1181
1182 if( A <= 12 ) ll = 6.0;
1183 else
1184 {
1185 ll = 6.0 + 1.35*log(G4double(A)/12.);
1186 }
1187 xx = RandGamma::shoot(aa,ll);
1188 shift = 0.99;
1189 mom = xx*shift*kF;
1190 }
1191 else
1192 {
1193 f2p2h = true;
1194 aa = 6.5;
1195 ll = 6.5;
1196 xx = RandGamma::shoot(aa,ll);
1197 shift = 2.5;
1198 mom = xx*shift*kF;
1199 }
1200 if( mom > momMax ) mom = G4UniformRand()*momMax;
1201 if( mom > 2.*kF ) f2p2h = true;
1202
1203 // mom = 0.;
1204
1205 return mom;
1206}
static double shoot()
G4double FermiMomentum(G4Nucleus &targetNucleus)

Referenced by G4ANuElNucleusCcModel::SampleLVkr(), G4ANuMuNucleusCcModel::SampleLVkr(), G4ANuTauNucleusCcModel::SampleLVkr(), G4NuElNucleusCcModel::SampleLVkr(), G4NuMuNucleusCcModel::SampleLVkr(), and G4NuTauNucleusCcModel::SampleLVkr().

◆ IsApplicable()

G4bool G4NeutrinoNucleusModel::IsApplicable ( const G4HadProjectile & aTrack,
G4Nucleus & targetNucleus )
virtual

Reimplemented from G4HadronicInteraction.

Reimplemented in G4ANuElNucleusCcModel, G4ANuElNucleusNcModel, G4ANuMuNucleusCcModel, G4ANuMuNucleusNcModel, G4ANuTauNucleusCcModel, G4ANuTauNucleusNcModel, G4NuElNucleusCcModel, G4NuElNucleusNcModel, G4NuMuNucleusCcModel, G4NuMuNucleusNcModel, G4NuTauNucleusCcModel, and G4NuTauNucleusNcModel.

Definition at line 206 of file G4NeutrinoNucleusModel.cc.

208{
209 G4bool result = false;
210 G4String pName = aPart.GetDefinition()->GetParticleName();
211 G4double energy = aPart.GetTotalEnergy();
212
213 if( pName == "nu_mu" // || pName == "anti_nu_mu" )
214 &&
215 energy > fMinNuEnergy )
216 {
217 result = true;
218 }
219
220 return result;
221}

◆ MesonDecay()

void G4NeutrinoNucleusModel::MesonDecay ( G4LorentzVector & lvX,
G4int qX )

Definition at line 859 of file G4NeutrinoNucleusModel.cc.

860{
861 G4bool finB = false;
862 G4int pdgM(0), pdgB(0), i(0), qM(0), qB(0);
863 G4double mM(0.), mB(0.), eM(0.), eB(0.), pM(0.), pB(0.);
864 G4double mm1(0.), mm22(0.), M1(0.), M2(0.), mX(0.), Tkin(0.);
865
866 mX = lvX.m();
867 Tkin = lvX.e() - mX;
868
869 // if( mX < 1120*MeV && mX > 1020*MeV ) // phi(1020)->K+K-
870 if( mX < 1080*MeV && mX > 990*MeV && Tkin < 600*MeV ) // phi(1020)->K+K-
871 {
872 return FinalMeson( lvX, qB, 333);
873 }
875
876 G4double deltaMr[4] = { 0.*MeV, 0.*MeV, 100.*MeV, 0.*MeV};
877
878 G4ThreeVector dir(0.,0.,0.);
879 G4ThreeVector bst(0.,0.,0.);
880 G4LorentzVector lvM(0.,0.,0.,0.);
881 G4LorentzVector lvB(0.,0.,0.,0.);
882
883 for( i = 0; i < fClustNumber; ++i) // check resonance
884 {
885 if( mX >= fMesMass[i] )
886 {
887 pdgB = fMesPDG[i];
888 // mB = G4ParticleTable::GetParticleTable()->FindParticle(pdgB)->GetPDGMass();
889 break;
890 }
891 }
892 if( i == fClustNumber ) // || i == fClustNumber-1 ) // low mass, p || n
893 {
894 if ( qX == 1) { pdgB = 211; qB = 1;} // pi+
895 else if( qX == 0 ) { pdgB = 111; qB = 0;} // pi0
896 else if( qX == -1) { pdgB = -211; qB = -1;} // pi-
897
898 return FinalMeson( lvX, qB, pdgB);
899 }
900 else if( mX < fMesMass[i] + deltaMr[i] ) // || mX < mPi + mPi ) //
901 {
902 finB = true; // final barion -> out
903 pdgB = fMesPDG[i];
904
905 // if ( qX == 1 && pdgB != 2212) pdgB = pdgB - 10;
906
907 if( qX == 0 ) pdgB = pdgB - 100;
908 else if( qX == -1 ) pdgB = -pdgB;
909
910 if( finB ) return FinalMeson( lvX, qX, pdgB ); // out
911 }
912 // no resonance, try 1->2 decay in COM frame
913
914 // try meson
915
916 mm1 = mPi + 1.*MeV; // pi+
917 mm22 = mX - mPi - 1.*MeV; // mX-n
918
919 if( mm22 <= mm1 ) // out
920 {
921 if ( qX == 1) { pdgB = 211; qB = 1;} // pi+
922 else if( qX == 0 ) { pdgB = 111; qB = 0;} // pi0
923 else if( qX == -1) { pdgB = -211; qB = -1;} // pi-
924
925 return FinalMeson(lvX, qB, pdgB);
926 }
927 else // try decay -> pion + meson(cluster)
928 {
929 // G4double sigmaM = 50.*MeV; // 100.*MeV; // 200.*MeV; // 400.*MeV; // 800.*MeV; //
930 G4double rand = G4UniformRand();
931
932 if ( qX == 1 ) { qM = 1; qB = 0;}
933 else if( qX == 0 ) { qM = -1; qB = 1;} // { qM = 0; qB = 0;} //
934 else if( qX == -1 ) { qM = -1; qB = 0;}
935 /*
936 mM = mPi;
937 if(qM == 0) mM = G4ParticleTable::GetParticleTable()->FindParticle(111)->GetPDGMass(); //pi0
938 pdgM = fMesPDG[fClustNumber-1];
939 */
940 // mm1*mm22/( mm1 + rand*(mm22 - mm1) );
941 // mM = mm1*mm22/sqrt( mm1*mm1 + rand*(mm22*mm22 - mm1*mm1) );
942 // mM = -sigmaM*log( (1.- rand)*exp(-mm22/sigmaM) + rand*exp(-mm1/sigmaM) );
943 mM = mm1 + rand*(mm22-mm1);
944 // mM = mm1 + 0.9*(mm22-mm1);
945
946
947 for( i = 0; i < fClustNumber; ++i)
948 {
949 if( mM >= fMesMass[i] )
950 {
951 pdgM = fMesPDG[i];
952 // mM = G4ParticleTable::GetParticleTable()->FindParticle(pdgM)->GetPDGMass();
953 break;
954 }
955 }
956 if( i == fClustNumber || i == fClustNumber-1 ) // low mass, p || n
957 {
958 if ( qX == 1) { pdgB = 211; qB = 1;} // pi+
959 else if( qX == 0 ) { pdgB = 111; qB = 0;} // pi0
960 else if( qX == -1) { pdgB = -211; qB = -1;} // pi-
961
962 return FinalMeson( lvX, qB, pdgB);
963 }
964 else if( mX < fMesMass[i] + deltaMr[i] ) // || mX < mPi + mPi ) //
965 {
966 finB = true; // final barion -> out
967 pdgB = fMesPDG[i];
968
969 // if ( qX == 1 && pdgB != 2212) pdgB = pdgB - 10;
970
971 if( qX == 0 ) pdgB = pdgB - 100;
972 else if( qX == -1 ) pdgB = -pdgB;
973
974 if( finB ) return FinalMeson( lvX, qX, pdgB ); // out
975 }
976
978 M2 = mX - mM;
979
980 if( M2 <= M1 ) //
981 {
982 if ( qX == 1) { pdgB = 211; qB = 1;} // pi+
983 else if( qX == 0 ) { pdgB = 111; qB = 0;} // pi0
984 else if( qX == -1) { pdgB = -211; qB = -1;} // pi-
985
986 return FinalMeson(lvX, qB, pdgB);
987 }
988 mB = M1 + G4UniformRand()*(M2-M1);
989 // mB = -sigmaM*log( (1.- rand)*exp(-M2/sigmaM) + rand*exp(-M1/sigmaM) );
990 // mB = M1 + 0.9*(M2-M1);
991
992 bst = lvX.boostVector();
993
994 // dir = G4RandomDirection();
995 dir = bst.orthogonal().unit();
996
997 eM = 0.5*(mX*mX + mM*mM - mB*mB)/mX;
998 pM = sqrt(eM*eM - mM*mM);
999 lvM = G4LorentzVector( pM*dir, eM);
1000 lvM.boost(bst);
1001
1002 eB = 0.5*(mX*mX + mB*mB - mM*mM)/mX;
1003 pB = sqrt(eB*eB - mB*mB);
1004 lvB = G4LorentzVector(-pB*dir, eB);
1005 lvB.boost(bst);
1006
1007 // G4cout<<mM<<"/"<<mB<<", ";
1008
1009 // charge exchange
1010
1011 // if ( qX == 2 ) { qM = 1; qB = 1;}
1012
1013 if ( qM == 0 ) pdgM = pdgM - 100;
1014 else if( qM == -1 ) pdgM = -pdgM;
1015
1016 MesonDecay( lvM, qM ); //
1017
1018 MesonDecay( lvB, qB ); // continue
1019 }
1020}
static const G4int fMesPDG[4]
void FinalMeson(G4LorentzVector &lvM, G4int qM, G4int pdgM)

Referenced by ClusterDecay(), and MesonDecay().

◆ ModelDescription()

void G4NeutrinoNucleusModel::ModelDescription ( std::ostream & outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Reimplemented in G4ANuElNucleusCcModel, G4ANuElNucleusNcModel, G4ANuMuNucleusCcModel, G4ANuMuNucleusNcModel, G4ANuTauNucleusCcModel, G4ANuTauNucleusNcModel, G4NuElNucleusCcModel, G4NuElNucleusNcModel, G4NuMuNucleusCcModel, G4NuMuNucleusNcModel, G4NuTauNucleusCcModel, and G4NuTauNucleusNcModel.

Definition at line 195 of file G4NeutrinoNucleusModel.cc.

196{
197
198 outFile << "G4NeutrinoNucleusModel is a neutrino-nucleus general\n"
199 << "model which uses the standard model \n"
200 << "transfer parameterization. The model is fully relativistic\n";
201
202}

◆ NucleonMomentum()

G4double G4NeutrinoNucleusModel::NucleonMomentum ( G4Nucleus & targetNucleus)

Definition at line 1088 of file G4NeutrinoNucleusModel.cc.

1089{
1090 G4int A = targetNucleus.GetA_asInt();
1091 G4double kF = FermiMomentum( targetNucleus);
1092 G4double mom(0.), kCut = 0.5*GeV; // kCut = 1.*GeV; // kCut = 2.*GeV; // kCut = 4.*GeV; //
1093 // G4double cof = 2./GeV;
1094 // G4double ksi = kF*kF*cof*cof/pi/pi;
1095 G4double th = 1.; // 1. - 6.*ksi; //
1096
1097 if( G4UniformRand() < th || A < 3 ) // 1p1h
1098 {
1099 mom = kF*pow( G4UniformRand(), 1./3.);
1100 }
1101 else // 2p2h
1102 {
1103 mom = kF*kCut;
1104 mom /= kCut - G4UniformRand()*(kCut - kF);
1105 f2p2h = true;
1106 }
1107 return mom;
1108}

Referenced by G4ANuElNucleusNcModel::SampleLVkr(), G4ANuMuNucleusNcModel::SampleLVkr(), G4ANuTauNucleusNcModel::SampleLVkr(), G4NuElNucleusNcModel::SampleLVkr(), G4NuMuNucleusNcModel::SampleLVkr(), and G4NuTauNucleusNcModel::SampleLVkr().

◆ RecoilDeexcitation()

void G4NeutrinoNucleusModel::RecoilDeexcitation ( G4Fragment & fragment)

Definition at line 574 of file G4NeutrinoNucleusModel.cc.

575{
576 G4ReactionProductVector* products = fPreCompound->DeExcite(fragment);
577
578 if( products != nullptr )
579 {
580 for( auto & prod : *products ) // prod is the pointer to final hadronic particle
581 {
582 theParticleChange.AddSecondary(new G4DynamicParticle( prod->GetDefinition(),
583 prod->GetTotalEnergy(),
584 prod->GetMomentum() ), fSecID );
585 delete prod;
586 }
587 delete products;
588 }
589}
std::vector< G4ReactionProduct * > G4ReactionProductVector
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment) final

Referenced by CoherentPion(), and FinalBarion().

◆ SampleQkr()

G4double G4NeutrinoNucleusModel::SampleQkr ( G4double energy,
G4double xx )

Definition at line 306 of file G4NeutrinoNucleusModel.cc.

307{
308 G4int nBin(50), iE=fEindex, jX=fXindex;
309 G4double qq(0.), qq1(0.), qq2(0.);
310 G4double prob = G4UniformRand();
311
312 // first E
313
314 if( iE <= 0 )
315 {
316 qq1 = GetQkr( 0, jX, prob);
317 }
318 else if ( iE >= nBin-1)
319 {
320 qq1 = GetQkr( nBin-1, jX, prob);
321 }
322 else
323 {
324 G4double q1 = GetQkr(iE-1,jX, prob);
325 G4double q2 = GetQkr(iE,jX, prob);
326
329 G4double e = G4Log(energy);
330
331 if( e2 <= e1) qq1 = q1 + G4UniformRand()*(q2-q1);
332 else qq1 = q1 + (e-e1)*(q2-q1)/(e2-e1); // lin in energy log-scale
333 }
334
335 // then X
336
337 if( jX <= 0 )
338 {
339 qq2 = GetQkr( iE, 0, prob);
340 }
341 else if ( jX >= nBin)
342 {
343 qq2 = GetQkr( iE, nBin, prob);
344 }
345 else
346 {
347 G4double q1 = GetQkr(iE,jX-1, prob);
348 G4double q2 = GetQkr(iE,jX, prob);
349
350 G4double e1 = G4Log(fNuMuXarrayKR[iE][jX-1]);
351 G4double e2 = G4Log(fNuMuXarrayKR[iE][jX]);
352 G4double e = G4Log(xx);
353
354 if( e2 <= e1) qq2 = q1 + G4UniformRand()*(q2-q1);
355 else qq2 = q1 + (e-e1)*(q2-q1)/(e2-e1); // lin in energy log-scale
356 }
357 qq = 0.5*(qq1+qq2);
358
359 return qq;
360}
G4double G4Log(G4double x)
Definition G4Log.hh:227
G4double GetQkr(G4int iE, G4int jX, G4double prob)
static const G4double fNuMuEnergyLogVector[50]

Referenced by G4ANuElNucleusCcModel::SampleLVkr(), G4ANuElNucleusNcModel::SampleLVkr(), G4ANuMuNucleusCcModel::SampleLVkr(), G4ANuMuNucleusNcModel::SampleLVkr(), G4ANuTauNucleusCcModel::SampleLVkr(), G4ANuTauNucleusNcModel::SampleLVkr(), G4NuElNucleusCcModel::SampleLVkr(), G4NuElNucleusNcModel::SampleLVkr(), G4NuMuNucleusCcModel::SampleLVkr(), G4NuMuNucleusNcModel::SampleLVkr(), G4NuTauNucleusCcModel::SampleLVkr(), and G4NuTauNucleusNcModel::SampleLVkr().

◆ SampleXkr()

G4double G4NeutrinoNucleusModel::SampleXkr ( G4double energy)

Definition at line 225 of file G4NeutrinoNucleusModel.cc.

226{
227 G4int i(0), nBin(50);
228 G4double xx(0.), prob = G4UniformRand();
229
230 for( i = 0; i < nBin; ++i )
231 {
232 if( energy <= fNuMuEnergyLogVector[i] ) break;
233 }
234 if( i <= 0) // E-edge
235 {
236 fEindex = 0;
237 xx = GetXkr( 0, prob);
238 }
239 else if ( i >= nBin)
240 {
241 fEindex = nBin-1;
242 xx = GetXkr( nBin-1, prob);
243 }
244 else
245 {
246 fEindex = i;
247 G4double x1 = GetXkr(i-1,prob);
248 G4double x2 = GetXkr(i,prob);
249
252 G4double e = G4Log(energy);
253
254 if( e2 <= e1) xx = x1 + G4UniformRand()*(x2-x1);
255 else xx = x1 + (e-e1)*(x2-x1)/(e2-e1); // lin in energy log-scale
256 }
257 return xx;
258}
G4double GetXkr(G4int iEnergy, G4double prob)

Referenced by G4ANuElNucleusCcModel::SampleLVkr(), G4ANuElNucleusNcModel::SampleLVkr(), G4ANuMuNucleusCcModel::SampleLVkr(), G4ANuMuNucleusNcModel::SampleLVkr(), G4ANuTauNucleusCcModel::SampleLVkr(), G4ANuTauNucleusNcModel::SampleLVkr(), G4NuElNucleusCcModel::SampleLVkr(), G4NuElNucleusNcModel::SampleLVkr(), G4NuMuNucleusCcModel::SampleLVkr(), G4NuMuNucleusNcModel::SampleLVkr(), G4NuTauNucleusCcModel::SampleLVkr(), and G4NuTauNucleusNcModel::SampleLVkr().

◆ SetCutEnergy()

void G4NeutrinoNucleusModel::SetCutEnergy ( G4double ec)
inline

Definition at line 98 of file G4NeutrinoNucleusModel.hh.

98{fCutEnergy=ec;};

◆ SetQEratioA()

void G4NeutrinoNucleusModel::SetQEratioA ( G4double qea)
inline

Definition at line 137 of file G4NeutrinoNucleusModel.hh.

137{ fQEratioA = qea; };

◆ ThresholdEnergy()

G4double G4NeutrinoNucleusModel::ThresholdEnergy ( G4double mI,
G4double mF,
G4double mP )
inline

Definition at line 131 of file G4NeutrinoNucleusModel.hh.

132 {
133 G4double w = std::sqrt(fW2);
134 return w + 0.5*( (mP+mF)*(mP+mF)-(w+mI)*(w+mI) )/mI;
135 };

Member Data Documentation

◆ f2p2h

◆ fANeMuQEratio

const G4double G4NeutrinoNucleusModel::fANeMuQEratio
staticprotected
Initial value:
=
{
1, 1, 1, 1, 1, 1, 1, 0.97506, 0.920938, 0.847671, 0.762973, 0.677684, 0.597685,
0.52538, 0.461466, 0.405329, 0.356154, 0.312944, 0.274984, 0.241341, 0.211654, 0.185322,
0.161991, 0.141339, 0.123078, 0.106952, 0.0927909, 0.0803262, 0.0693698, 0.0598207, 0.0514545,
0.044193, 0.0378696, 0.0324138, 0.0276955, 0.0236343, 0.0201497, 0.0171592, 0.014602, 0.0124182,
0.0105536, 0.00896322, 0.00761004, 0.00645821, 0.00547859, 0.00464595, 0.00393928,
0.00333961, 0.00283086, 0.00239927
}

Definition at line 1418 of file G4NeutrinoNucleusModel.hh.

1479 : public G4HadronicInteraction
1480{
1481public:
1482
1483 G4NeutrinoNucleusModel(const G4String& name = "neutrino-nucleus");
1484
1485 virtual ~G4NeutrinoNucleusModel();
1486
1487 virtual G4bool IsApplicable(const G4HadProjectile & aTrack,
1488 G4Nucleus & targetNucleus);
1489
1490 G4double SampleXkr(G4double energy);
1491 G4double GetXkr(G4int iEnergy, G4double prob);
1492 G4double SampleQkr(G4double energy, G4double xx);
1493 G4double GetQkr(G4int iE, G4int jX, G4double prob);
1494
1495 virtual G4HadFinalState * ApplyYourself(const G4HadProjectile & aTrack,
1496 G4Nucleus & targetNucleus)=0;
1497
1498 //////// fragmentation functions /////////////////////////
1499
1500 void ClusterDecay( G4LorentzVector & lvX, G4int qX);
1501
1502 void MesonDecay( G4LorentzVector & lvX, G4int qX);
1503
1504 void FinalBarion( G4LorentzVector & lvB, G4int qB, G4int pdgB);
1505
1506 void RecoilDeexcitation( G4Fragment& fragment);
1507
1508 void FinalMeson( G4LorentzVector & lvM, G4int qM, G4int pdgM);
1509
1510 void CoherentPion( G4LorentzVector & lvP, G4int pdgP, G4Nucleus & targetNucleus);
1511
1512
1513 // set/get class fields
1514
1515 void SetCutEnergy(G4double ec){fCutEnergy=ec;};
1516 G4double GetCutEnergy(){return fCutEnergy;};
1517
1518 G4double GetNuEnergy(){return fNuEnergy;};
1519 G4double GetQtransfer(){return fQtransfer;};
1520 G4double GetQ2(){return fQ2;};
1521 G4double GetXsample(){return fXsample;};
1522
1524 G4bool GetCascade(){return fCascade;};
1525 G4bool GetString(){return fString;};
1526
1527 G4double GetCosTheta(){return fCosTheta;};
1528 G4double GetEmu(){return fEmu;};
1529 G4double GetEx(){return fEx;};
1530 G4double GetMuMass(){return fMu;};
1531 G4double GetW2(){return fW2;};
1532 G4double GetM1(){return fM1;};
1533 G4double GetMr(){return fMr;};
1534 G4double GetTr(){return fTr;};
1535 G4double GetDp(){return fDp;};
1536
1537 G4bool GetfBreak() {return fBreak;};
1538 G4bool GetfCascade(){return fCascade;};
1539 G4bool GetfString() {return fString;};
1540
1541 G4LorentzVector GetLVl(){return fLVl;};
1542 G4LorentzVector GetLVh(){return fLVh;};
1543 G4LorentzVector GetLVt(){return fLVt;};
1544 G4LorentzVector GetLVcpi(){return fLVcpi;};
1545
1546 G4double GetMinNuMuEnergy(){ return fMu + 0.5*fMu*fMu/fM1 + 4.*CLHEP::MeV; }; // kinematics + accuracy for sqrts
1547
1548 G4double ThresholdEnergy(G4double mI, G4double mF, G4double mP) // for cluster decay
1549 {
1550 G4double w = std::sqrt(fW2);
1551 return w + 0.5*( (mP+mF)*(mP+mF)-(w+mI)*(w+mI) )/mI;
1552 };
1553 G4double GetQEratioA(){ return fQEratioA; };
1554 void SetQEratioA( G4double qea ){ fQEratioA = qea; };
1555
1556
1557 G4double FinalMomentum(G4double mI, G4double mF, G4double mP, G4LorentzVector lvX); // for cluster decay
1558
1559 // nucleon binding
1560
1561 G4double FermiMomentum( G4Nucleus & targetNucleus);
1562 G4double NucleonMomentum( G4Nucleus & targetNucleus);
1563
1564 G4double GetEx( G4int A, G4bool fP );
1566
1568 G4double GetNuMuQeTotRat(G4int index, G4double energy);
1569
1572
1573 G4double CalculateQEratioA( G4int Z, G4int A, G4double energy, G4int nepdg);
1574
1575 virtual void ModelDescription(std::ostream&) const;
1576
1577protected:
1578
1581
1582 G4double fSin2tW; // sin^2theta_Weinberg
1583 G4double fCutEnergy; // minimal recoil electron energy detected
1584
1587
1589
1591
1593
1595
1599
1601
1602 G4int fSecID; // Creator model ID for the secondaries created by this model
1603
1604 static const G4int fResNumber;
1605 static const G4double fResMass[6]; // [fResNumber];
1606
1607 static const G4int fClustNumber;
1608
1609 static const G4double fMesMass[4];
1610 static const G4int fMesPDG[4];
1611
1612 static const G4double fBarMass[4];
1613 static const G4int fBarPDG[4];
1614
1615 static const G4double fNuMuResQ[50][50];
1616
1617
1618 static const G4double fNuMuEnergy[50];
1619 static const G4double fNuMuQeTotRat[50];
1620 static const G4double fOnePionEnergy[58];
1621 static const G4double fOnePionProb[58];
1622
1623 static const G4double fNuMuEnergyLogVector[50];
1624
1625 // KR sample distributions, X at E_nu and Q2 at E_nu and X
1626
1627 static G4double fNuMuXarrayKR[50][51];
1628 static G4double fNuMuXdistrKR[50][50];
1629 static G4double fNuMuQarrayKR[50][51][51];
1630 static G4double fNuMuQdistrKR[50][51][50];
1631
1632 // QEratio(Z,A,Enu)
1633
1634 static const G4double fQEnergy[50];
1635 static const G4double fANeMuQEratio[50];
1636 static const G4double fNeMuQEratio[50];
1637
1638};
1639
1640
1641
1642#endif
void CoherentPion(G4LorentzVector &lvP, G4int pdgP, G4Nucleus &targetNucleus)
virtual void ModelDescription(std::ostream &) const
G4double NucleonMomentum(G4Nucleus &targetNucleus)
G4int GetOnePionIndex(G4double energy)
static const G4double fResMass[6]
G4double SampleXkr(G4double energy)
G4double SampleQkr(G4double energy, G4double xx)
G4double GetNuMuQeTotRat(G4int index, G4double energy)
G4int GetEnergyIndex(G4double energy)
G4double GgSampleNM(G4Nucleus &nucl)
G4double GetNuMuOnePionProb(G4int index, G4double energy)
G4double ThresholdEnergy(G4double mI, G4double mF, G4double mP)
G4double CalculateQEratioA(G4int Z, G4int A, G4double energy, G4int nepdg)
static const G4double fNuMuResQ[50][50]
virtual G4bool IsApplicable(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4NeutrinoNucleusModel(const G4String &name="neutrino-nucleus")
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0

Referenced by CalculateQEratioA().

◆ fBarMass

const G4double G4NeutrinoNucleusModel::fBarMass = {1700., 1600., 1232., 939.57}
staticprotected

Definition at line 92 of file G4NeutrinoNucleusModel.hh.

Referenced by ClusterDecay().

◆ fBarPDG

const G4int G4NeutrinoNucleusModel::fBarPDG = {12224, 32224, 2224, 2212}
staticprotected

Definition at line 93 of file G4NeutrinoNucleusModel.hh.

Referenced by ClusterDecay().

◆ fBreak

◆ fCascade

◆ fClustNumber

const G4int G4NeutrinoNucleusModel::fClustNumber = 4
staticprotected

Definition at line 190 of file G4NeutrinoNucleusModel.hh.

Referenced by ClusterDecay(), and MesonDecay().

◆ fCosTheta

◆ fCosThetaPi

G4double G4NeutrinoNucleusModel::fCosThetaPi
protected

Definition at line 175 of file G4NeutrinoNucleusModel.hh.

Referenced by G4NeutrinoNucleusModel().

◆ fCutEnergy

G4double G4NeutrinoNucleusModel::fCutEnergy
protected

Definition at line 166 of file G4NeutrinoNucleusModel.hh.

Referenced by G4NeutrinoNucleusModel(), GetCutEnergy(), and SetCutEnergy().

◆ fDeExcitation

G4ExcitationHandler* G4NeutrinoNucleusModel::fDeExcitation
protected

Definition at line 181 of file G4NeutrinoNucleusModel.hh.

Referenced by G4NeutrinoNucleusModel().

◆ fDp

G4double G4NeutrinoNucleusModel::fDp
protected

Definition at line 173 of file G4NeutrinoNucleusModel.hh.

Referenced by FinalBarion(), G4NeutrinoNucleusModel(), and GetDp().

◆ fEindex

G4int G4NeutrinoNucleusModel::fEindex
protected

Definition at line 168 of file G4NeutrinoNucleusModel.hh.

Referenced by G4NeutrinoNucleusModel(), SampleQkr(), and SampleXkr().

◆ fEmu

◆ fEmuPi

G4double G4NeutrinoNucleusModel::fEmuPi
protected

Definition at line 175 of file G4NeutrinoNucleusModel.hh.

Referenced by G4NeutrinoNucleusModel().

◆ fEx

G4double G4NeutrinoNucleusModel::fEx
protected

Definition at line 175 of file G4NeutrinoNucleusModel.hh.

Referenced by G4NeutrinoNucleusModel(), and GetEx().

◆ fIndex

G4int G4NeutrinoNucleusModel::fIndex
protected

◆ fLVcpi

◆ fLVh

◆ fLVl

◆ fLVt

◆ fM1

G4double G4NeutrinoNucleusModel::fM1
protected

Definition at line 173 of file G4NeutrinoNucleusModel.hh.

Referenced by G4ANuElNucleusCcModel::ApplyYourself(), G4ANuElNucleusNcModel::ApplyYourself(), G4ANuMuNucleusCcModel::ApplyYourself(), G4ANuMuNucleusNcModel::ApplyYourself(), G4ANuTauNucleusCcModel::ApplyYourself(), G4ANuTauNucleusNcModel::ApplyYourself(), G4NuElNucleusCcModel::ApplyYourself(), G4NuElNucleusNcModel::ApplyYourself(), G4NuMuNucleusCcModel::ApplyYourself(), G4NuMuNucleusNcModel::ApplyYourself(), G4NuTauNucleusCcModel::ApplyYourself(), G4NuTauNucleusNcModel::ApplyYourself(), G4NeutrinoNucleusModel(), GetM1(), G4ANuElNucleusCcModel::GetMinNuElEnergy(), G4ANuElNucleusNcModel::GetMinNuElEnergy(), G4NuElNucleusCcModel::GetMinNuElEnergy(), G4NuElNucleusNcModel::GetMinNuElEnergy(), G4ANuMuNucleusCcModel::GetMinNuMuEnergy(), G4ANuMuNucleusNcModel::GetMinNuMuEnergy(), G4ANuTauNucleusCcModel::GetMinNuMuEnergy(), G4ANuTauNucleusNcModel::GetMinNuMuEnergy(), GetMinNuMuEnergy(), G4NuMuNucleusCcModel::GetMinNuMuEnergy(), G4NuMuNucleusNcModel::GetMinNuMuEnergy(), G4NuTauNucleusCcModel::GetMinNuMuEnergy(), G4NuTauNucleusNcModel::GetMinNuMuEnergy(), G4ANuElNucleusCcModel::SampleLVkr(), G4ANuElNucleusNcModel::SampleLVkr(), G4ANuMuNucleusCcModel::SampleLVkr(), G4ANuMuNucleusNcModel::SampleLVkr(), G4ANuTauNucleusCcModel::SampleLVkr(), G4ANuTauNucleusNcModel::SampleLVkr(), G4NuElNucleusCcModel::SampleLVkr(), G4NuElNucleusNcModel::SampleLVkr(), G4NuMuNucleusCcModel::SampleLVkr(), G4NuMuNucleusNcModel::SampleLVkr(), G4NuTauNucleusCcModel::SampleLVkr(), and G4NuTauNucleusNcModel::SampleLVkr().

◆ fM2

G4double G4NeutrinoNucleusModel::fM2
protected

Definition at line 173 of file G4NeutrinoNucleusModel.hh.

Referenced by G4NeutrinoNucleusModel().

◆ fMesMass

const G4double G4NeutrinoNucleusModel::fMesMass = {1260., 980., 770., 139.57}
staticprotected

Definition at line 86 of file G4NeutrinoNucleusModel.hh.

Referenced by ClusterDecay(), and MesonDecay().

◆ fMesPDG

const G4int G4NeutrinoNucleusModel::fMesPDG = {20213, 9000211, 213, 211}
staticprotected

Definition at line 87 of file G4NeutrinoNucleusModel.hh.

Referenced by MesonDecay().

◆ fMinNuEnergy

◆ fMpi

◆ fMr

◆ fMt

◆ fMu

◆ fNbin

◆ fNeMuQEratio

const G4double G4NeutrinoNucleusModel::fNeMuQEratio
staticprotected
Initial value:
=
{
1, 1, 1, 1, 1, 1, 1, 0.977592, 0.926073, 0.858783, 0.783874, 0.706868, 0.63113, 0.558681,
0.490818, 0.428384, 0.371865, 0.321413, 0.276892, 0.237959, 0.204139, 0.1749, 0.149706, 0.128047,
0.109456, 0.093514, 0.0798548, 0.0681575, 0.0581455, 0.0495804, 0.0422578, 0.036002, 0.0306614,
0.0261061, 0.0222231, 0.0189152, 0.0160987, 0.0137011, 0.0116604, 0.00992366, 0.00844558, 0.00718766,
0.00611714, 0.00520618, 0.00443105, 0.00377158, 0.00321062, 0.0027335, 0.00232774, 0.00198258
}

Definition at line 1428 of file G4NeutrinoNucleusModel.hh.

1489 : public G4HadronicInteraction
1490{
1491public:
1492
1493 G4NeutrinoNucleusModel(const G4String& name = "neutrino-nucleus");
1494
1495 virtual ~G4NeutrinoNucleusModel();
1496
1497 virtual G4bool IsApplicable(const G4HadProjectile & aTrack,
1498 G4Nucleus & targetNucleus);
1499
1500 G4double SampleXkr(G4double energy);
1501 G4double GetXkr(G4int iEnergy, G4double prob);
1502 G4double SampleQkr(G4double energy, G4double xx);
1503 G4double GetQkr(G4int iE, G4int jX, G4double prob);
1504
1505 virtual G4HadFinalState * ApplyYourself(const G4HadProjectile & aTrack,
1506 G4Nucleus & targetNucleus)=0;
1507
1508 //////// fragmentation functions /////////////////////////
1509
1510 void ClusterDecay( G4LorentzVector & lvX, G4int qX);
1511
1512 void MesonDecay( G4LorentzVector & lvX, G4int qX);
1513
1514 void FinalBarion( G4LorentzVector & lvB, G4int qB, G4int pdgB);
1515
1516 void RecoilDeexcitation( G4Fragment& fragment);
1517
1518 void FinalMeson( G4LorentzVector & lvM, G4int qM, G4int pdgM);
1519
1520 void CoherentPion( G4LorentzVector & lvP, G4int pdgP, G4Nucleus & targetNucleus);
1521
1522
1523 // set/get class fields
1524
1525 void SetCutEnergy(G4double ec){fCutEnergy=ec;};
1526 G4double GetCutEnergy(){return fCutEnergy;};
1527
1528 G4double GetNuEnergy(){return fNuEnergy;};
1529 G4double GetQtransfer(){return fQtransfer;};
1530 G4double GetQ2(){return fQ2;};
1531 G4double GetXsample(){return fXsample;};
1532
1534 G4bool GetCascade(){return fCascade;};
1535 G4bool GetString(){return fString;};
1536
1537 G4double GetCosTheta(){return fCosTheta;};
1538 G4double GetEmu(){return fEmu;};
1539 G4double GetEx(){return fEx;};
1540 G4double GetMuMass(){return fMu;};
1541 G4double GetW2(){return fW2;};
1542 G4double GetM1(){return fM1;};
1543 G4double GetMr(){return fMr;};
1544 G4double GetTr(){return fTr;};
1545 G4double GetDp(){return fDp;};
1546
1547 G4bool GetfBreak() {return fBreak;};
1548 G4bool GetfCascade(){return fCascade;};
1549 G4bool GetfString() {return fString;};
1550
1551 G4LorentzVector GetLVl(){return fLVl;};
1552 G4LorentzVector GetLVh(){return fLVh;};
1553 G4LorentzVector GetLVt(){return fLVt;};
1554 G4LorentzVector GetLVcpi(){return fLVcpi;};
1555
1556 G4double GetMinNuMuEnergy(){ return fMu + 0.5*fMu*fMu/fM1 + 4.*CLHEP::MeV; }; // kinematics + accuracy for sqrts
1557
1558 G4double ThresholdEnergy(G4double mI, G4double mF, G4double mP) // for cluster decay
1559 {
1560 G4double w = std::sqrt(fW2);
1561 return w + 0.5*( (mP+mF)*(mP+mF)-(w+mI)*(w+mI) )/mI;
1562 };
1563 G4double GetQEratioA(){ return fQEratioA; };
1564 void SetQEratioA( G4double qea ){ fQEratioA = qea; };
1565
1566
1567 G4double FinalMomentum(G4double mI, G4double mF, G4double mP, G4LorentzVector lvX); // for cluster decay
1568
1569 // nucleon binding
1570
1571 G4double FermiMomentum( G4Nucleus & targetNucleus);
1572 G4double NucleonMomentum( G4Nucleus & targetNucleus);
1573
1574 G4double GetEx( G4int A, G4bool fP );
1576
1578 G4double GetNuMuQeTotRat(G4int index, G4double energy);
1579
1582
1583 G4double CalculateQEratioA( G4int Z, G4int A, G4double energy, G4int nepdg);
1584
1585 virtual void ModelDescription(std::ostream&) const;
1586
1587protected:
1588
1591
1592 G4double fSin2tW; // sin^2theta_Weinberg
1593 G4double fCutEnergy; // minimal recoil electron energy detected
1594
1597
1599
1601
1603
1605
1609
1611
1612 G4int fSecID; // Creator model ID for the secondaries created by this model
1613
1614 static const G4int fResNumber;
1615 static const G4double fResMass[6]; // [fResNumber];
1616
1617 static const G4int fClustNumber;
1618
1619 static const G4double fMesMass[4];
1620 static const G4int fMesPDG[4];
1621
1622 static const G4double fBarMass[4];
1623 static const G4int fBarPDG[4];
1624
1625 static const G4double fNuMuResQ[50][50];
1626
1627
1628 static const G4double fNuMuEnergy[50];
1629 static const G4double fNuMuQeTotRat[50];
1630 static const G4double fOnePionEnergy[58];
1631 static const G4double fOnePionProb[58];
1632
1633 static const G4double fNuMuEnergyLogVector[50];
1634
1635 // KR sample distributions, X at E_nu and Q2 at E_nu and X
1636
1637 static G4double fNuMuXarrayKR[50][51];
1638 static G4double fNuMuXdistrKR[50][50];
1639 static G4double fNuMuQarrayKR[50][51][51];
1640 static G4double fNuMuQdistrKR[50][51][50];
1641
1642 // QEratio(Z,A,Enu)
1643
1644 static const G4double fQEnergy[50];
1645 static const G4double fANeMuQEratio[50];
1646 static const G4double fNeMuQEratio[50];
1647
1648};
1649
1650
1651
1652#endif

Referenced by CalculateQEratioA().

◆ fNuEnergy

◆ fNuMuEnergy

const G4double G4NeutrinoNucleusModel::fNuMuEnergy
staticprotected
Initial value:
=
{
0.112103, 0.117359, 0.123119, 0.129443, 0.136404,
0.144084, 0.152576, 0.161991, 0.172458, 0.184126,
0.197171, 0.211801, 0.228261, 0.24684, 0.267887,
0.291816, 0.319125, 0.350417, 0.386422, 0.428032,
0.47634, 0.532692, 0.598756, 0.676612, 0.768868,
0.878812, 1.01062, 1.16963, 1.36271, 1.59876,
1.88943, 2.25002, 2.70086, 3.26916, 3.99166,
4.91843, 6.11836, 7.6872, 9.75942, 12.5259,
16.2605, 21.3615, 28.4141, 38.2903, 52.3062,
72.4763, 101.93, 145.6, 211.39, 312.172
}

Definition at line 1259 of file G4NeutrinoNucleusModel.hh.

1320 : public G4HadronicInteraction
1321{
1322public:
1323
1324 G4NeutrinoNucleusModel(const G4String& name = "neutrino-nucleus");
1325
1326 virtual ~G4NeutrinoNucleusModel();
1327
1328 virtual G4bool IsApplicable(const G4HadProjectile & aTrack,
1329 G4Nucleus & targetNucleus);
1330
1331 G4double SampleXkr(G4double energy);
1332 G4double GetXkr(G4int iEnergy, G4double prob);
1333 G4double SampleQkr(G4double energy, G4double xx);
1334 G4double GetQkr(G4int iE, G4int jX, G4double prob);
1335
1336 virtual G4HadFinalState * ApplyYourself(const G4HadProjectile & aTrack,
1337 G4Nucleus & targetNucleus)=0;
1338
1339 //////// fragmentation functions /////////////////////////
1340
1341 void ClusterDecay( G4LorentzVector & lvX, G4int qX);
1342
1343 void MesonDecay( G4LorentzVector & lvX, G4int qX);
1344
1345 void FinalBarion( G4LorentzVector & lvB, G4int qB, G4int pdgB);
1346
1347 void RecoilDeexcitation( G4Fragment& fragment);
1348
1349 void FinalMeson( G4LorentzVector & lvM, G4int qM, G4int pdgM);
1350
1351 void CoherentPion( G4LorentzVector & lvP, G4int pdgP, G4Nucleus & targetNucleus);
1352
1353
1354 // set/get class fields
1355
1356 void SetCutEnergy(G4double ec){fCutEnergy=ec;};
1357 G4double GetCutEnergy(){return fCutEnergy;};
1358
1359 G4double GetNuEnergy(){return fNuEnergy;};
1360 G4double GetQtransfer(){return fQtransfer;};
1361 G4double GetQ2(){return fQ2;};
1362 G4double GetXsample(){return fXsample;};
1363
1365 G4bool GetCascade(){return fCascade;};
1366 G4bool GetString(){return fString;};
1367
1368 G4double GetCosTheta(){return fCosTheta;};
1369 G4double GetEmu(){return fEmu;};
1370 G4double GetEx(){return fEx;};
1371 G4double GetMuMass(){return fMu;};
1372 G4double GetW2(){return fW2;};
1373 G4double GetM1(){return fM1;};
1374 G4double GetMr(){return fMr;};
1375 G4double GetTr(){return fTr;};
1376 G4double GetDp(){return fDp;};
1377
1378 G4bool GetfBreak() {return fBreak;};
1379 G4bool GetfCascade(){return fCascade;};
1380 G4bool GetfString() {return fString;};
1381
1382 G4LorentzVector GetLVl(){return fLVl;};
1383 G4LorentzVector GetLVh(){return fLVh;};
1384 G4LorentzVector GetLVt(){return fLVt;};
1385 G4LorentzVector GetLVcpi(){return fLVcpi;};
1386
1387 G4double GetMinNuMuEnergy(){ return fMu + 0.5*fMu*fMu/fM1 + 4.*CLHEP::MeV; }; // kinematics + accuracy for sqrts
1388
1389 G4double ThresholdEnergy(G4double mI, G4double mF, G4double mP) // for cluster decay
1390 {
1391 G4double w = std::sqrt(fW2);
1392 return w + 0.5*( (mP+mF)*(mP+mF)-(w+mI)*(w+mI) )/mI;
1393 };
1394 G4double GetQEratioA(){ return fQEratioA; };
1395 void SetQEratioA( G4double qea ){ fQEratioA = qea; };
1396
1397
1398 G4double FinalMomentum(G4double mI, G4double mF, G4double mP, G4LorentzVector lvX); // for cluster decay
1399
1400 // nucleon binding
1401
1402 G4double FermiMomentum( G4Nucleus & targetNucleus);
1403 G4double NucleonMomentum( G4Nucleus & targetNucleus);
1404
1405 G4double GetEx( G4int A, G4bool fP );
1407
1409 G4double GetNuMuQeTotRat(G4int index, G4double energy);
1410
1413
1414 G4double CalculateQEratioA( G4int Z, G4int A, G4double energy, G4int nepdg);
1415
1416 virtual void ModelDescription(std::ostream&) const;
1417
1418protected:
1419
1422
1423 G4double fSin2tW; // sin^2theta_Weinberg
1424 G4double fCutEnergy; // minimal recoil electron energy detected
1425
1428
1430
1432
1434
1436
1440
1442
1443 G4int fSecID; // Creator model ID for the secondaries created by this model
1444
1445 static const G4int fResNumber;
1446 static const G4double fResMass[6]; // [fResNumber];
1447
1448 static const G4int fClustNumber;
1449
1450 static const G4double fMesMass[4];
1451 static const G4int fMesPDG[4];
1452
1453 static const G4double fBarMass[4];
1454 static const G4int fBarPDG[4];
1455
1456 static const G4double fNuMuResQ[50][50];
1457
1458
1459 static const G4double fNuMuEnergy[50];
1460 static const G4double fNuMuQeTotRat[50];
1461 static const G4double fOnePionEnergy[58];
1462 static const G4double fOnePionProb[58];
1463
1464 static const G4double fNuMuEnergyLogVector[50];
1465
1466 // KR sample distributions, X at E_nu and Q2 at E_nu and X
1467
1468 static G4double fNuMuXarrayKR[50][51];
1469 static G4double fNuMuXdistrKR[50][50];
1470 static G4double fNuMuQarrayKR[50][51][51];
1471 static G4double fNuMuQdistrKR[50][51][50];
1472
1473 // QEratio(Z,A,Enu)
1474
1475 static const G4double fQEnergy[50];
1476 static const G4double fANeMuQEratio[50];
1477 static const G4double fNeMuQEratio[50];
1478
1479};
1480
1481
1482
1483#endif

Referenced by GetEnergyIndex(), and GetNuMuQeTotRat().

◆ fNuMuEnergyLogVector

const G4double G4NeutrinoNucleusModel::fNuMuEnergyLogVector
staticprotected
Initial value:
= {
115.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,
745.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,
4806.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,
30988.8, 35765.7, 41279, 47642.2, 54986.3, 63462.4, 73245.2, 84536, 97567.2, 112607, 129966 }

Definition at line 95 of file G4NeutrinoNucleusModel.hh.

98 {fCutEnergy=ec;};

Referenced by SampleQkr(), and SampleXkr().

◆ fNuMuQarrayKR

◆ fNuMuQdistrKR

◆ fNuMuQeTotRat

const G4double G4NeutrinoNucleusModel::fNuMuQeTotRat
staticprotected
Initial value:
=
{
0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98,
0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98, 0.98,
0.97, 0.96, 0.95, 0.93,
0.917794, 0.850239, 0.780412, 0.709339, 0.638134, 0.568165,
0.500236, 0.435528, 0.375015, 0.319157, 0.268463, 0.2232, 0.183284,
0.148627, 0.119008, 0.0940699, 0.0733255, 0.0563819, 0.0427312, 0.0319274,
0.0235026, 0.0170486, 0.0122149, 0.00857825, 0.00594018, 0.00405037
}

Definition at line 1275 of file G4NeutrinoNucleusModel.hh.

1336 : public G4HadronicInteraction
1337{
1338public:
1339
1340 G4NeutrinoNucleusModel(const G4String& name = "neutrino-nucleus");
1341
1342 virtual ~G4NeutrinoNucleusModel();
1343
1344 virtual G4bool IsApplicable(const G4HadProjectile & aTrack,
1345 G4Nucleus & targetNucleus);
1346
1347 G4double SampleXkr(G4double energy);
1348 G4double GetXkr(G4int iEnergy, G4double prob);
1349 G4double SampleQkr(G4double energy, G4double xx);
1350 G4double GetQkr(G4int iE, G4int jX, G4double prob);
1351
1352 virtual G4HadFinalState * ApplyYourself(const G4HadProjectile & aTrack,
1353 G4Nucleus & targetNucleus)=0;
1354
1355 //////// fragmentation functions /////////////////////////
1356
1357 void ClusterDecay( G4LorentzVector & lvX, G4int qX);
1358
1359 void MesonDecay( G4LorentzVector & lvX, G4int qX);
1360
1361 void FinalBarion( G4LorentzVector & lvB, G4int qB, G4int pdgB);
1362
1363 void RecoilDeexcitation( G4Fragment& fragment);
1364
1365 void FinalMeson( G4LorentzVector & lvM, G4int qM, G4int pdgM);
1366
1367 void CoherentPion( G4LorentzVector & lvP, G4int pdgP, G4Nucleus & targetNucleus);
1368
1369
1370 // set/get class fields
1371
1372 void SetCutEnergy(G4double ec){fCutEnergy=ec;};
1373 G4double GetCutEnergy(){return fCutEnergy;};
1374
1375 G4double GetNuEnergy(){return fNuEnergy;};
1376 G4double GetQtransfer(){return fQtransfer;};
1377 G4double GetQ2(){return fQ2;};
1378 G4double GetXsample(){return fXsample;};
1379
1381 G4bool GetCascade(){return fCascade;};
1382 G4bool GetString(){return fString;};
1383
1384 G4double GetCosTheta(){return fCosTheta;};
1385 G4double GetEmu(){return fEmu;};
1386 G4double GetEx(){return fEx;};
1387 G4double GetMuMass(){return fMu;};
1388 G4double GetW2(){return fW2;};
1389 G4double GetM1(){return fM1;};
1390 G4double GetMr(){return fMr;};
1391 G4double GetTr(){return fTr;};
1392 G4double GetDp(){return fDp;};
1393
1394 G4bool GetfBreak() {return fBreak;};
1395 G4bool GetfCascade(){return fCascade;};
1396 G4bool GetfString() {return fString;};
1397
1398 G4LorentzVector GetLVl(){return fLVl;};
1399 G4LorentzVector GetLVh(){return fLVh;};
1400 G4LorentzVector GetLVt(){return fLVt;};
1401 G4LorentzVector GetLVcpi(){return fLVcpi;};
1402
1403 G4double GetMinNuMuEnergy(){ return fMu + 0.5*fMu*fMu/fM1 + 4.*CLHEP::MeV; }; // kinematics + accuracy for sqrts
1404
1405 G4double ThresholdEnergy(G4double mI, G4double mF, G4double mP) // for cluster decay
1406 {
1407 G4double w = std::sqrt(fW2);
1408 return w + 0.5*( (mP+mF)*(mP+mF)-(w+mI)*(w+mI) )/mI;
1409 };
1410 G4double GetQEratioA(){ return fQEratioA; };
1411 void SetQEratioA( G4double qea ){ fQEratioA = qea; };
1412
1413
1414 G4double FinalMomentum(G4double mI, G4double mF, G4double mP, G4LorentzVector lvX); // for cluster decay
1415
1416 // nucleon binding
1417
1418 G4double FermiMomentum( G4Nucleus & targetNucleus);
1419 G4double NucleonMomentum( G4Nucleus & targetNucleus);
1420
1421 G4double GetEx( G4int A, G4bool fP );
1423
1425 G4double GetNuMuQeTotRat(G4int index, G4double energy);
1426
1429
1430 G4double CalculateQEratioA( G4int Z, G4int A, G4double energy, G4int nepdg);
1431
1432 virtual void ModelDescription(std::ostream&) const;
1433
1434protected:
1435
1438
1439 G4double fSin2tW; // sin^2theta_Weinberg
1440 G4double fCutEnergy; // minimal recoil electron energy detected
1441
1444
1446
1448
1450
1452
1456
1458
1459 G4int fSecID; // Creator model ID for the secondaries created by this model
1460
1461 static const G4int fResNumber;
1462 static const G4double fResMass[6]; // [fResNumber];
1463
1464 static const G4int fClustNumber;
1465
1466 static const G4double fMesMass[4];
1467 static const G4int fMesPDG[4];
1468
1469 static const G4double fBarMass[4];
1470 static const G4int fBarPDG[4];
1471
1472 static const G4double fNuMuResQ[50][50];
1473
1474
1475 static const G4double fNuMuEnergy[50];
1476 static const G4double fNuMuQeTotRat[50];
1477 static const G4double fOnePionEnergy[58];
1478 static const G4double fOnePionProb[58];
1479
1480 static const G4double fNuMuEnergyLogVector[50];
1481
1482 // KR sample distributions, X at E_nu and Q2 at E_nu and X
1483
1484 static G4double fNuMuXarrayKR[50][51];
1485 static G4double fNuMuXdistrKR[50][50];
1486 static G4double fNuMuQarrayKR[50][51][51];
1487 static G4double fNuMuQdistrKR[50][51][50];
1488
1489 // QEratio(Z,A,Enu)
1490
1491 static const G4double fQEnergy[50];
1492 static const G4double fANeMuQEratio[50];
1493 static const G4double fNeMuQEratio[50];
1494
1495};
1496
1497
1498
1499#endif

Referenced by GetNuMuQeTotRat().

◆ fNuMuResQ

const G4double G4NeutrinoNucleusModel::fNuMuResQ[50][50]
staticprotected

Definition at line 198 of file G4NeutrinoNucleusModel.hh.

◆ fNuMuXarrayKR

◆ fNuMuXdistrKR

◆ fOnePionEnergy

const G4double G4NeutrinoNucleusModel::fOnePionEnergy
staticprotected
Initial value:
=
{
0.275314, 0.293652, 0.31729, 0.33409, 0.351746, 0.365629, 0.380041, 0.400165, 0.437941, 0.479237,
0.504391, 0.537803, 0.588487, 0.627532, 0.686839, 0.791905, 0.878332, 0.987405, 1.08162, 1.16971,
1.2982, 1.40393, 1.49854, 1.64168, 1.7524, 1.87058, 2.02273, 2.15894, 2.3654, 2.55792, 2.73017,
3.03005, 3.40733, 3.88128, 4.53725, 5.16786, 5.73439, 6.53106, 7.43879, 8.36214, 9.39965, 10.296,
11.5735, 13.1801, 15.2052, 17.5414, 19.7178, 22.7462, 25.9026, 29.4955, 33.5867, 39.2516, 46.4716,
53.6065, 63.4668, 73.2147, 85.5593, 99.9854
}

Definition at line 1339 of file G4NeutrinoNucleusModel.hh.

1400 : public G4HadronicInteraction
1401{
1402public:
1403
1404 G4NeutrinoNucleusModel(const G4String& name = "neutrino-nucleus");
1405
1406 virtual ~G4NeutrinoNucleusModel();
1407
1408 virtual G4bool IsApplicable(const G4HadProjectile & aTrack,
1409 G4Nucleus & targetNucleus);
1410
1411 G4double SampleXkr(G4double energy);
1412 G4double GetXkr(G4int iEnergy, G4double prob);
1413 G4double SampleQkr(G4double energy, G4double xx);
1414 G4double GetQkr(G4int iE, G4int jX, G4double prob);
1415
1416 virtual G4HadFinalState * ApplyYourself(const G4HadProjectile & aTrack,
1417 G4Nucleus & targetNucleus)=0;
1418
1419 //////// fragmentation functions /////////////////////////
1420
1421 void ClusterDecay( G4LorentzVector & lvX, G4int qX);
1422
1423 void MesonDecay( G4LorentzVector & lvX, G4int qX);
1424
1425 void FinalBarion( G4LorentzVector & lvB, G4int qB, G4int pdgB);
1426
1427 void RecoilDeexcitation( G4Fragment& fragment);
1428
1429 void FinalMeson( G4LorentzVector & lvM, G4int qM, G4int pdgM);
1430
1431 void CoherentPion( G4LorentzVector & lvP, G4int pdgP, G4Nucleus & targetNucleus);
1432
1433
1434 // set/get class fields
1435
1436 void SetCutEnergy(G4double ec){fCutEnergy=ec;};
1437 G4double GetCutEnergy(){return fCutEnergy;};
1438
1439 G4double GetNuEnergy(){return fNuEnergy;};
1440 G4double GetQtransfer(){return fQtransfer;};
1441 G4double GetQ2(){return fQ2;};
1442 G4double GetXsample(){return fXsample;};
1443
1445 G4bool GetCascade(){return fCascade;};
1446 G4bool GetString(){return fString;};
1447
1448 G4double GetCosTheta(){return fCosTheta;};
1449 G4double GetEmu(){return fEmu;};
1450 G4double GetEx(){return fEx;};
1451 G4double GetMuMass(){return fMu;};
1452 G4double GetW2(){return fW2;};
1453 G4double GetM1(){return fM1;};
1454 G4double GetMr(){return fMr;};
1455 G4double GetTr(){return fTr;};
1456 G4double GetDp(){return fDp;};
1457
1458 G4bool GetfBreak() {return fBreak;};
1459 G4bool GetfCascade(){return fCascade;};
1460 G4bool GetfString() {return fString;};
1461
1462 G4LorentzVector GetLVl(){return fLVl;};
1463 G4LorentzVector GetLVh(){return fLVh;};
1464 G4LorentzVector GetLVt(){return fLVt;};
1465 G4LorentzVector GetLVcpi(){return fLVcpi;};
1466
1467 G4double GetMinNuMuEnergy(){ return fMu + 0.5*fMu*fMu/fM1 + 4.*CLHEP::MeV; }; // kinematics + accuracy for sqrts
1468
1469 G4double ThresholdEnergy(G4double mI, G4double mF, G4double mP) // for cluster decay
1470 {
1471 G4double w = std::sqrt(fW2);
1472 return w + 0.5*( (mP+mF)*(mP+mF)-(w+mI)*(w+mI) )/mI;
1473 };
1474 G4double GetQEratioA(){ return fQEratioA; };
1475 void SetQEratioA( G4double qea ){ fQEratioA = qea; };
1476
1477
1478 G4double FinalMomentum(G4double mI, G4double mF, G4double mP, G4LorentzVector lvX); // for cluster decay
1479
1480 // nucleon binding
1481
1482 G4double FermiMomentum( G4Nucleus & targetNucleus);
1483 G4double NucleonMomentum( G4Nucleus & targetNucleus);
1484
1485 G4double GetEx( G4int A, G4bool fP );
1487
1489 G4double GetNuMuQeTotRat(G4int index, G4double energy);
1490
1493
1494 G4double CalculateQEratioA( G4int Z, G4int A, G4double energy, G4int nepdg);
1495
1496 virtual void ModelDescription(std::ostream&) const;
1497
1498protected:
1499
1502
1503 G4double fSin2tW; // sin^2theta_Weinberg
1504 G4double fCutEnergy; // minimal recoil electron energy detected
1505
1508
1510
1512
1514
1516
1520
1522
1523 G4int fSecID; // Creator model ID for the secondaries created by this model
1524
1525 static const G4int fResNumber;
1526 static const G4double fResMass[6]; // [fResNumber];
1527
1528 static const G4int fClustNumber;
1529
1530 static const G4double fMesMass[4];
1531 static const G4int fMesPDG[4];
1532
1533 static const G4double fBarMass[4];
1534 static const G4int fBarPDG[4];
1535
1536 static const G4double fNuMuResQ[50][50];
1537
1538
1539 static const G4double fNuMuEnergy[50];
1540 static const G4double fNuMuQeTotRat[50];
1541 static const G4double fOnePionEnergy[58];
1542 static const G4double fOnePionProb[58];
1543
1544 static const G4double fNuMuEnergyLogVector[50];
1545
1546 // KR sample distributions, X at E_nu and Q2 at E_nu and X
1547
1548 static G4double fNuMuXarrayKR[50][51];
1549 static G4double fNuMuXdistrKR[50][50];
1550 static G4double fNuMuQarrayKR[50][51][51];
1551 static G4double fNuMuQdistrKR[50][51][50];
1552
1553 // QEratio(Z,A,Enu)
1554
1555 static const G4double fQEnergy[50];
1556 static const G4double fANeMuQEratio[50];
1557 static const G4double fNeMuQEratio[50];
1558
1559};
1560
1561
1562
1563#endif

Referenced by GetNuMuOnePionProb(), GetNuMuQeTotRat(), and GetOnePionIndex().

◆ fOnePionIndex

G4int G4NeutrinoNucleusModel::fOnePionIndex
protected

◆ fOnePionProb

const G4double G4NeutrinoNucleusModel::fOnePionProb
staticprotected
Initial value:
=
{
0.0019357, 0.0189361, 0.0378722, 0.0502758, 0.0662559, 0.0754581, 0.0865008, 0.0987275, 0.124112,
0.153787, 0.18308, 0.213996, 0.245358, 0.274425, 0.301536, 0.326612, 0.338208, 0.337806, 0.335948,
0.328092, 0.313557, 0.304965, 0.292169, 0.28481, 0.269474, 0.254138, 0.247499, 0.236249, 0.221654,
0.205492, 0.198781, 0.182216, 0.162251, 0.142878, 0.128631, 0.116001, 0.108435, 0.0974843, 0.082092,
0.0755204, 0.0703121, 0.0607066, 0.0554278, 0.0480401, 0.0427023, 0.0377123, 0.0323248, 0.0298584,
0.0244296, 0.0218526, 0.019121, 0.016477, 0.0137309, 0.0137963, 0.0110371, 0.00834028, 0.00686127, 0.00538226
}

Definition at line 1353 of file G4NeutrinoNucleusModel.hh.

1414 : public G4HadronicInteraction
1415{
1416public:
1417
1418 G4NeutrinoNucleusModel(const G4String& name = "neutrino-nucleus");
1419
1420 virtual ~G4NeutrinoNucleusModel();
1421
1422 virtual G4bool IsApplicable(const G4HadProjectile & aTrack,
1423 G4Nucleus & targetNucleus);
1424
1425 G4double SampleXkr(G4double energy);
1426 G4double GetXkr(G4int iEnergy, G4double prob);
1427 G4double SampleQkr(G4double energy, G4double xx);
1428 G4double GetQkr(G4int iE, G4int jX, G4double prob);
1429
1430 virtual G4HadFinalState * ApplyYourself(const G4HadProjectile & aTrack,
1431 G4Nucleus & targetNucleus)=0;
1432
1433 //////// fragmentation functions /////////////////////////
1434
1435 void ClusterDecay( G4LorentzVector & lvX, G4int qX);
1436
1437 void MesonDecay( G4LorentzVector & lvX, G4int qX);
1438
1439 void FinalBarion( G4LorentzVector & lvB, G4int qB, G4int pdgB);
1440
1441 void RecoilDeexcitation( G4Fragment& fragment);
1442
1443 void FinalMeson( G4LorentzVector & lvM, G4int qM, G4int pdgM);
1444
1445 void CoherentPion( G4LorentzVector & lvP, G4int pdgP, G4Nucleus & targetNucleus);
1446
1447
1448 // set/get class fields
1449
1450 void SetCutEnergy(G4double ec){fCutEnergy=ec;};
1451 G4double GetCutEnergy(){return fCutEnergy;};
1452
1453 G4double GetNuEnergy(){return fNuEnergy;};
1454 G4double GetQtransfer(){return fQtransfer;};
1455 G4double GetQ2(){return fQ2;};
1456 G4double GetXsample(){return fXsample;};
1457
1459 G4bool GetCascade(){return fCascade;};
1460 G4bool GetString(){return fString;};
1461
1462 G4double GetCosTheta(){return fCosTheta;};
1463 G4double GetEmu(){return fEmu;};
1464 G4double GetEx(){return fEx;};
1465 G4double GetMuMass(){return fMu;};
1466 G4double GetW2(){return fW2;};
1467 G4double GetM1(){return fM1;};
1468 G4double GetMr(){return fMr;};
1469 G4double GetTr(){return fTr;};
1470 G4double GetDp(){return fDp;};
1471
1472 G4bool GetfBreak() {return fBreak;};
1473 G4bool GetfCascade(){return fCascade;};
1474 G4bool GetfString() {return fString;};
1475
1476 G4LorentzVector GetLVl(){return fLVl;};
1477 G4LorentzVector GetLVh(){return fLVh;};
1478 G4LorentzVector GetLVt(){return fLVt;};
1479 G4LorentzVector GetLVcpi(){return fLVcpi;};
1480
1481 G4double GetMinNuMuEnergy(){ return fMu + 0.5*fMu*fMu/fM1 + 4.*CLHEP::MeV; }; // kinematics + accuracy for sqrts
1482
1483 G4double ThresholdEnergy(G4double mI, G4double mF, G4double mP) // for cluster decay
1484 {
1485 G4double w = std::sqrt(fW2);
1486 return w + 0.5*( (mP+mF)*(mP+mF)-(w+mI)*(w+mI) )/mI;
1487 };
1488 G4double GetQEratioA(){ return fQEratioA; };
1489 void SetQEratioA( G4double qea ){ fQEratioA = qea; };
1490
1491
1492 G4double FinalMomentum(G4double mI, G4double mF, G4double mP, G4LorentzVector lvX); // for cluster decay
1493
1494 // nucleon binding
1495
1496 G4double FermiMomentum( G4Nucleus & targetNucleus);
1497 G4double NucleonMomentum( G4Nucleus & targetNucleus);
1498
1499 G4double GetEx( G4int A, G4bool fP );
1501
1503 G4double GetNuMuQeTotRat(G4int index, G4double energy);
1504
1507
1508 G4double CalculateQEratioA( G4int Z, G4int A, G4double energy, G4int nepdg);
1509
1510 virtual void ModelDescription(std::ostream&) const;
1511
1512protected:
1513
1516
1517 G4double fSin2tW; // sin^2theta_Weinberg
1518 G4double fCutEnergy; // minimal recoil electron energy detected
1519
1522
1524
1526
1528
1530
1534
1536
1537 G4int fSecID; // Creator model ID for the secondaries created by this model
1538
1539 static const G4int fResNumber;
1540 static const G4double fResMass[6]; // [fResNumber];
1541
1542 static const G4int fClustNumber;
1543
1544 static const G4double fMesMass[4];
1545 static const G4int fMesPDG[4];
1546
1547 static const G4double fBarMass[4];
1548 static const G4int fBarPDG[4];
1549
1550 static const G4double fNuMuResQ[50][50];
1551
1552
1553 static const G4double fNuMuEnergy[50];
1554 static const G4double fNuMuQeTotRat[50];
1555 static const G4double fOnePionEnergy[58];
1556 static const G4double fOnePionProb[58];
1557
1558 static const G4double fNuMuEnergyLogVector[50];
1559
1560 // KR sample distributions, X at E_nu and Q2 at E_nu and X
1561
1562 static G4double fNuMuXarrayKR[50][51];
1563 static G4double fNuMuXdistrKR[50][50];
1564 static G4double fNuMuQarrayKR[50][51][51];
1565 static G4double fNuMuQdistrKR[50][51][50];
1566
1567 // QEratio(Z,A,Enu)
1568
1569 static const G4double fQEnergy[50];
1570 static const G4double fANeMuQEratio[50];
1571 static const G4double fNeMuQEratio[50];
1572
1573};
1574
1575
1576
1577#endif

Referenced by GetNuMuOnePionProb().

◆ fPDGencoding

◆ fPrecoInterface

G4GeneratorPrecompoundInterface* G4NeutrinoNucleusModel::fPrecoInterface
protected

Definition at line 179 of file G4NeutrinoNucleusModel.hh.

Referenced by G4NeutrinoNucleusModel(), and ~G4NeutrinoNucleusModel().

◆ fPreCompound

G4PreCompoundModel* G4NeutrinoNucleusModel::fPreCompound
protected

Definition at line 180 of file G4NeutrinoNucleusModel.hh.

Referenced by G4NeutrinoNucleusModel(), and RecoilDeexcitation().

◆ fProton

◆ fQ2

◆ fQEnergy

const G4double G4NeutrinoNucleusModel::fQEnergy
staticprotected
Initial value:
=
{
0.12, 0.1416, 0.167088, 0.197164, 0.232653, 0.274531, 0.323946, 0.382257, 0.451063, 0.532254,
0.62806, 0.741111, 0.874511, 1.03192, 1.21767, 1.43685, 1.69548, 2.00067, 2.36079, 2.78573,
3.28716, 3.87885, 4.57705, 5.40092, 6.37308, 7.52024, 8.87388, 10.4712, 12.356, 14.5801,
17.2045, 20.3013, 23.9555, 28.2675, 33.3557, 39.3597, 46.4444, 54.8044, 64.6692, 76.3097,
90.0454, 106.254, 125.379, 147.947, 174.578, 206.002, 243.082, 286.837, 338.468, 399.392
}

Definition at line 1409 of file G4NeutrinoNucleusModel.hh.

1470 : public G4HadronicInteraction
1471{
1472public:
1473
1474 G4NeutrinoNucleusModel(const G4String& name = "neutrino-nucleus");
1475
1476 virtual ~G4NeutrinoNucleusModel();
1477
1478 virtual G4bool IsApplicable(const G4HadProjectile & aTrack,
1479 G4Nucleus & targetNucleus);
1480
1481 G4double SampleXkr(G4double energy);
1482 G4double GetXkr(G4int iEnergy, G4double prob);
1483 G4double SampleQkr(G4double energy, G4double xx);
1484 G4double GetQkr(G4int iE, G4int jX, G4double prob);
1485
1486 virtual G4HadFinalState * ApplyYourself(const G4HadProjectile & aTrack,
1487 G4Nucleus & targetNucleus)=0;
1488
1489 //////// fragmentation functions /////////////////////////
1490
1491 void ClusterDecay( G4LorentzVector & lvX, G4int qX);
1492
1493 void MesonDecay( G4LorentzVector & lvX, G4int qX);
1494
1495 void FinalBarion( G4LorentzVector & lvB, G4int qB, G4int pdgB);
1496
1497 void RecoilDeexcitation( G4Fragment& fragment);
1498
1499 void FinalMeson( G4LorentzVector & lvM, G4int qM, G4int pdgM);
1500
1501 void CoherentPion( G4LorentzVector & lvP, G4int pdgP, G4Nucleus & targetNucleus);
1502
1503
1504 // set/get class fields
1505
1506 void SetCutEnergy(G4double ec){fCutEnergy=ec;};
1507 G4double GetCutEnergy(){return fCutEnergy;};
1508
1509 G4double GetNuEnergy(){return fNuEnergy;};
1510 G4double GetQtransfer(){return fQtransfer;};
1511 G4double GetQ2(){return fQ2;};
1512 G4double GetXsample(){return fXsample;};
1513
1515 G4bool GetCascade(){return fCascade;};
1516 G4bool GetString(){return fString;};
1517
1518 G4double GetCosTheta(){return fCosTheta;};
1519 G4double GetEmu(){return fEmu;};
1520 G4double GetEx(){return fEx;};
1521 G4double GetMuMass(){return fMu;};
1522 G4double GetW2(){return fW2;};
1523 G4double GetM1(){return fM1;};
1524 G4double GetMr(){return fMr;};
1525 G4double GetTr(){return fTr;};
1526 G4double GetDp(){return fDp;};
1527
1528 G4bool GetfBreak() {return fBreak;};
1529 G4bool GetfCascade(){return fCascade;};
1530 G4bool GetfString() {return fString;};
1531
1532 G4LorentzVector GetLVl(){return fLVl;};
1533 G4LorentzVector GetLVh(){return fLVh;};
1534 G4LorentzVector GetLVt(){return fLVt;};
1535 G4LorentzVector GetLVcpi(){return fLVcpi;};
1536
1537 G4double GetMinNuMuEnergy(){ return fMu + 0.5*fMu*fMu/fM1 + 4.*CLHEP::MeV; }; // kinematics + accuracy for sqrts
1538
1539 G4double ThresholdEnergy(G4double mI, G4double mF, G4double mP) // for cluster decay
1540 {
1541 G4double w = std::sqrt(fW2);
1542 return w + 0.5*( (mP+mF)*(mP+mF)-(w+mI)*(w+mI) )/mI;
1543 };
1544 G4double GetQEratioA(){ return fQEratioA; };
1545 void SetQEratioA( G4double qea ){ fQEratioA = qea; };
1546
1547
1548 G4double FinalMomentum(G4double mI, G4double mF, G4double mP, G4LorentzVector lvX); // for cluster decay
1549
1550 // nucleon binding
1551
1552 G4double FermiMomentum( G4Nucleus & targetNucleus);
1553 G4double NucleonMomentum( G4Nucleus & targetNucleus);
1554
1555 G4double GetEx( G4int A, G4bool fP );
1557
1559 G4double GetNuMuQeTotRat(G4int index, G4double energy);
1560
1563
1564 G4double CalculateQEratioA( G4int Z, G4int A, G4double energy, G4int nepdg);
1565
1566 virtual void ModelDescription(std::ostream&) const;
1567
1568protected:
1569
1572
1573 G4double fSin2tW; // sin^2theta_Weinberg
1574 G4double fCutEnergy; // minimal recoil electron energy detected
1575
1578
1580
1582
1584
1586
1590
1592
1593 G4int fSecID; // Creator model ID for the secondaries created by this model
1594
1595 static const G4int fResNumber;
1596 static const G4double fResMass[6]; // [fResNumber];
1597
1598 static const G4int fClustNumber;
1599
1600 static const G4double fMesMass[4];
1601 static const G4int fMesPDG[4];
1602
1603 static const G4double fBarMass[4];
1604 static const G4int fBarPDG[4];
1605
1606 static const G4double fNuMuResQ[50][50];
1607
1608
1609 static const G4double fNuMuEnergy[50];
1610 static const G4double fNuMuQeTotRat[50];
1611 static const G4double fOnePionEnergy[58];
1612 static const G4double fOnePionProb[58];
1613
1614 static const G4double fNuMuEnergyLogVector[50];
1615
1616 // KR sample distributions, X at E_nu and Q2 at E_nu and X
1617
1618 static G4double fNuMuXarrayKR[50][51];
1619 static G4double fNuMuXdistrKR[50][50];
1620 static G4double fNuMuQarrayKR[50][51][51];
1621 static G4double fNuMuQdistrKR[50][51][50];
1622
1623 // QEratio(Z,A,Enu)
1624
1625 static const G4double fQEnergy[50];
1626 static const G4double fANeMuQEratio[50];
1627 static const G4double fNeMuQEratio[50];
1628
1629};
1630
1631
1632
1633#endif

Referenced by CalculateQEratioA().

◆ fQEratioA

G4double G4NeutrinoNucleusModel::fQEratioA
protected

◆ fQindex

G4int G4NeutrinoNucleusModel::fQindex
protected

Definition at line 168 of file G4NeutrinoNucleusModel.hh.

Referenced by G4NeutrinoNucleusModel(), and GetQkr().

◆ fQtransfer

◆ fRecoil

◆ fResMass

const G4double G4NeutrinoNucleusModel::fResMass
staticprotected
Initial value:
=
{2190., 1920., 1700., 1600., 1440., 1232. }

Definition at line 81 of file G4NeutrinoNucleusModel.hh.

◆ fResNumber

const G4int G4NeutrinoNucleusModel::fResNumber = 6
staticprotected

Definition at line 187 of file G4NeutrinoNucleusModel.hh.

◆ fSecID

◆ fSin2tW

G4double G4NeutrinoNucleusModel::fSin2tW
protected

Definition at line 165 of file G4NeutrinoNucleusModel.hh.

Referenced by G4NeutrinoNucleusModel().

◆ fString

◆ fTr

G4double G4NeutrinoNucleusModel::fTr
protected

Definition at line 173 of file G4NeutrinoNucleusModel.hh.

Referenced by FinalBarion(), G4NeutrinoNucleusModel(), and GetTr().

◆ fW2

G4double G4NeutrinoNucleusModel::fW2
protected

Definition at line 173 of file G4NeutrinoNucleusModel.hh.

Referenced by G4ANuElNucleusCcModel::ApplyYourself(), G4ANuElNucleusNcModel::ApplyYourself(), G4ANuMuNucleusCcModel::ApplyYourself(), G4ANuMuNucleusNcModel::ApplyYourself(), G4ANuTauNucleusCcModel::ApplyYourself(), G4ANuTauNucleusNcModel::ApplyYourself(), G4NuElNucleusCcModel::ApplyYourself(), G4NuElNucleusNcModel::ApplyYourself(), G4NuMuNucleusCcModel::ApplyYourself(), G4NuMuNucleusNcModel::ApplyYourself(), G4NuTauNucleusCcModel::ApplyYourself(), G4NuTauNucleusNcModel::ApplyYourself(), G4NeutrinoNucleusModel(), GetW2(), G4ANuElNucleusCcModel::SampleLVkr(), G4ANuElNucleusNcModel::SampleLVkr(), G4ANuMuNucleusCcModel::SampleLVkr(), G4ANuMuNucleusNcModel::SampleLVkr(), G4ANuTauNucleusCcModel::SampleLVkr(), G4ANuTauNucleusNcModel::SampleLVkr(), G4NuElNucleusCcModel::SampleLVkr(), G4NuElNucleusNcModel::SampleLVkr(), G4NuMuNucleusCcModel::SampleLVkr(), G4NuMuNucleusNcModel::SampleLVkr(), G4NuTauNucleusCcModel::SampleLVkr(), G4NuTauNucleusNcModel::SampleLVkr(), G4ANuElNucleusCcModel::ThresholdEnergy(), G4ANuElNucleusNcModel::ThresholdEnergy(), G4ANuMuNucleusCcModel::ThresholdEnergy(), G4ANuMuNucleusNcModel::ThresholdEnergy(), G4ANuTauNucleusCcModel::ThresholdEnergy(), G4ANuTauNucleusNcModel::ThresholdEnergy(), ThresholdEnergy(), G4NuElNucleusCcModel::ThresholdEnergy(), G4NuElNucleusNcModel::ThresholdEnergy(), G4NuMuNucleusCcModel::ThresholdEnergy(), G4NuMuNucleusNcModel::ThresholdEnergy(), G4NuTauNucleusCcModel::ThresholdEnergy(), and G4NuTauNucleusNcModel::ThresholdEnergy().

◆ fW2pi

G4double G4NeutrinoNucleusModel::fW2pi
protected

Definition at line 173 of file G4NeutrinoNucleusModel.hh.

Referenced by G4NeutrinoNucleusModel().

◆ fXindex

G4int G4NeutrinoNucleusModel::fXindex
protected

Definition at line 168 of file G4NeutrinoNucleusModel.hh.

Referenced by G4NeutrinoNucleusModel(), GetXkr(), and SampleQkr().

◆ fXsample

◆ theMuonMinus

G4ParticleDefinition* G4NeutrinoNucleusModel::theMuonMinus
protected

◆ theMuonPlus

G4ParticleDefinition* G4NeutrinoNucleusModel::theMuonPlus
protected

The documentation for this class was generated from the following files: