Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ANuElNucleusNcModel Class Reference

#include <G4ANuElNucleusNcModel.hh>

+ Inheritance diagram for G4ANuElNucleusNcModel:

Public Member Functions

 G4ANuElNucleusNcModel (const G4String &name="NuMuNuclNcModel")
 
virtual ~G4ANuElNucleusNcModel ()
 
virtual void InitialiseModel ()
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
void SampleLVkr (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double GetMinNuElEnergy ()
 
G4double ThresholdEnergy (G4double mI, G4double mF, G4double mP)
 
virtual void ModelDescription (std::ostream &) const
 
- Public Member Functions inherited from G4NeutrinoNucleusModel
 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 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)
 
virtual void ModelDescription (std::ostream &) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
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 ModelDescription (std::ostream &outFile) const
 
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
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4NeutrinoNucleusModel
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
 
G4LorentzVector fLVh
 
G4LorentzVector fLVl
 
G4LorentzVector fLVt
 
G4LorentzVector fLVcpi
 
G4GeneratorPrecompoundInterfacefPrecoInterface
 
G4PreCompoundModelfPreCompound
 
G4ExcitationHandlerfDeExcitation
 
G4NucleusfRecoil
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 
- Static Protected Attributes inherited from G4NeutrinoNucleusModel
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}}}
 

Detailed Description

Definition at line 60 of file G4ANuElNucleusNcModel.hh.

Constructor & Destructor Documentation

◆ G4ANuElNucleusNcModel()

G4ANuElNucleusNcModel::G4ANuElNucleusNcModel ( const G4String name = "NuMuNuclNcModel")

Definition at line 68 of file G4ANuElNucleusNcModel.cc.

70{
71 SetMinEnergy( 0.0*GeV );
72 SetMaxEnergy( 100.*TeV );
73 SetMinEnergy(1.e-6*eV);
74
75 theNuE = G4NeutrinoE::NeutrinoE();
76
77 fMnumu = 0.;
78 fData = fMaster = false;
80
81}
void SetMinEnergy(G4double anEnergy)
void SetMaxEnergy(const G4double anEnergy)
static G4NeutrinoE * NeutrinoE()
Definition: G4NeutrinoE.cc:84

◆ ~G4ANuElNucleusNcModel()

G4ANuElNucleusNcModel::~G4ANuElNucleusNcModel ( )
virtual

Definition at line 84 of file G4ANuElNucleusNcModel.cc.

85{}

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4ANuElNucleusNcModel::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
virtual

Implements G4NeutrinoNucleusModel.

Definition at line 222 of file G4ANuElNucleusNcModel.cc.

224{
226 fProton = f2p2h = fBreak = false;
227 const G4HadProjectile* aParticle = &aTrack;
228 G4double energy = aParticle->GetTotalEnergy();
229
230 G4String pName = aParticle->GetDefinition()->GetParticleName();
231
232 if( energy < fMinNuEnergy )
233 {
236 return &theParticleChange;
237 }
238 SampleLVkr( aTrack, targetNucleus);
239
240 if( fBreak == true || fEmu < fMnumu ) // ~5*10^-6
241 {
242 // G4cout<<"ni, ";
245 return &theParticleChange;
246 }
247
248 // LVs of initial state
249
250 G4LorentzVector lvp1 = aParticle->Get4Momentum();
251 G4LorentzVector lvt1( 0., 0., 0., fM1 );
253
254 // 1-pi by fQtransfer && nu-energy
255 G4LorentzVector lvpip1( 0., 0., 0., mPip );
256 G4LorentzVector lvsum, lv2, lvX;
257 G4ThreeVector eP;
258 G4double cost(1.), sint(0.), phi(0.), muMom(0.), massX2(0.);
259 G4DynamicParticle* aLept = nullptr; // lepton lv
260
261 G4int Z = targetNucleus.GetZ_asInt();
262 G4int A = targetNucleus.GetA_asInt();
263 G4double mTarg = targetNucleus.AtomicMass(A,Z);
264 G4int pdgP(0), qB(0);
265 // G4double mSum = G4ParticleTable::GetParticleTable()->FindParticle(2212)->GetPDGMass() + mPip;
266
267 G4int iPi = GetOnePionIndex(energy);
268 G4double p1pi = GetNuMuOnePionProb( iPi, energy);
269
270 if( p1pi > G4UniformRand() && fCosTheta > 0.9 ) // && fQtransfer < 0.95*GeV ) // mu- & coherent pion + nucleus
271 {
272 // lvsum = lvp1 + lvpip1;
273 lvsum = lvp1 + lvt1;
274 // cost = fCosThetaPi;
275 cost = fCosTheta;
276 sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
277 phi = G4UniformRand()*CLHEP::twopi;
278 eP = G4ThreeVector( sint*std::cos(phi), sint*std::sin(phi), cost );
279
280 // muMom = sqrt(fEmuPi*fEmuPi-fMnumu*fMnumu);
281 muMom = sqrt(fEmu*fEmu-fMnumu*fMnumu);
282
283 eP *= muMom;
284
285 // lv2 = G4LorentzVector( eP, fEmuPi );
286 lv2 = G4LorentzVector( eP, fEmu );
287 lv2 = fLVl;
288
289 lvX = lvsum - lv2;
290 lvX = fLVh;
291 massX2 = lvX.m2();
292 G4double massX = lvX.m();
293 G4double massR = fLVt.m();
294
295 // if ( massX2 <= 0. ) // vmg: very rarely ~ (1-4)e-6 due to big Q2/x, to be improved
296 if ( massX2 <= fM1*fM1 ) // 9-3-20 vmg: very rarely ~ (1-4)e-6 due to big Q2/x, to be improved
297 if ( lvX.e() <= fM1 ) // 9-3-20 vmg: very rarely ~ (1-4)e-6 due to big Q2/x, to be improved
298 {
301 return &theParticleChange;
302 }
303 fW2 = massX2;
304
305 if( pName == "anti_nu_e" ) aLept = new G4DynamicParticle( theNuE, lv2 );
306 // else if( pName == "anti_nu_mu") aLept = new G4DynamicParticle( theANuMu, lv2 );
307 else
308 {
311 return &theParticleChange;
312 }
313
314 pdgP = 111;
315
316 G4double eCut; // = fMpi + 0.5*(fMpi*fMpi - massX2)/mTarg; // massX -> fMpi
317
318 if( A > 1 )
319 {
320 eCut = (fMpi + mTarg)*(fMpi + mTarg) - (massX + massR)*(massX + massR);
321 eCut /= 2.*massR;
322 eCut += massX;
323 }
324 else eCut = fM1 + fMpi;
325
326 if ( lvX.e() > eCut ) // && sqrt( GetW2() ) < 1.4*GeV ) //
327 {
328 CoherentPion( lvX, pdgP, targetNucleus);
329 }
330 else
331 {
334 return &theParticleChange;
335 }
337
338 return &theParticleChange;
339 }
340 else // lepton part in lab
341 {
342 lvsum = lvp1 + lvt1;
343 cost = fCosTheta;
344 sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
345 phi = G4UniformRand()*CLHEP::twopi;
346 eP = G4ThreeVector( sint*std::cos(phi), sint*std::sin(phi), cost );
347
348 muMom = sqrt(fEmu*fEmu-fMnumu*fMnumu);
349
350 eP *= muMom;
351
352 lv2 = G4LorentzVector( eP, fEmu );
353
354 lvX = lvsum - lv2;
355
356 massX2 = lvX.m2();
357
358 if ( massX2 <= 0. ) // vmg: very rarely ~ (1-4)e-6 due to big Q2/x, to be improved
359 {
362 return &theParticleChange;
363 }
364 fW2 = massX2;
365
366 aLept = new G4DynamicParticle( theNuE, lv2 );
367
369 }
370
371 // hadron part
372
373 fRecoil = nullptr;
374 fCascade = false;
375 fString = false;
376
377 if( A == 1 )
378 {
379 qB = 1;
380
381 // if( G4UniformRand() > 0.1 ) // > 0.9999 ) // > 0.0001 ) //
382 {
383 ClusterDecay( lvX, qB );
384 }
385 return &theParticleChange;
386 }
387 G4Nucleus recoil;
388 G4double rM(0.), ratio = G4double(Z)/G4double(A);
389
390 if( ratio > G4UniformRand() ) // proton is excited
391 {
392 fProton = true;
393 recoil = G4Nucleus(A-1,Z-1);
394 fRecoil = &recoil;
395 rM = recoil.AtomicMass(A-1,Z-1);
396
399 }
400 else // excited neutron
401 {
402 fProton = false;
403 recoil = G4Nucleus(A-1,Z);
404 fRecoil = &recoil;
405 rM = recoil.AtomicMass(A-1,Z);
406
409 }
410 G4int index = GetEnergyIndex(energy);
411 G4double qeTotRat = GetNuMuQeTotRat(index, energy);
412
413 G4ThreeVector dX = (lvX.vect()).unit();
414 G4double eX = lvX.e(); // excited nucleon
415 G4double mX = sqrt(massX2);
416
417 if( qeTotRat > G4UniformRand() || mX <= fMt ) // || eX <= 1232.*MeV) // QE
418 {
419 fString = false;
420
421 if( fProton )
422 {
423 fPDGencoding = 2212;
424 fMr = proton_mass_c2;
425 recoil = G4Nucleus(A-1,Z-1);
426 fRecoil = &recoil;
427 rM = recoil.AtomicMass(A-1,Z-1);
428 }
429 else
430 {
431 fPDGencoding = 2112;
433 FindParticle(fPDGencoding)->GetPDGMass(); // 939.5654133*MeV;
434 recoil = G4Nucleus(A-1,Z);
435 fRecoil = &recoil;
436 rM = recoil.AtomicMass(A-1,Z);
437 }
438 G4double eTh = fMr+0.5*(fMr*fMr-mX*mX)/rM;
439
440 if(eX <= eTh) // vmg, very rarely out of kinematics
441 {
444 return &theParticleChange;
445 }
446 FinalBarion( lvX, 0, fPDGencoding ); // p(n)+deexcited recoil
447 }
448 else // if ( eX < 9500000.*GeV ) // < 25.*GeV) // < 95.*GeV ) // < 2.5*GeV ) //cluster decay
449 {
450 if ( fProton && pName == "anti_nu_e" ) qB = 1;
451 else if( !fProton && pName == "anti_nu_e" ) qB = 0;
452
453 ClusterDecay( lvX, qB );
454 }
455 return &theParticleChange;
456}
double A(double temperature)
CLHEP::HepLorentzVector G4LorentzVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector vect() const
void SampleLVkr(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
void SetEnergyChange(G4double anEnergy)
void SetMomentumChange(const G4ThreeVector &aV)
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const
G4double GetTotalEnergy() const
void CoherentPion(G4LorentzVector &lvP, G4int pdgP, G4Nucleus &targetNucleus)
G4int GetOnePionIndex(G4double energy)
G4double GetNuMuQeTotRat(G4int index, G4double energy)
G4int GetEnergyIndex(G4double energy)
G4double GetNuMuOnePionProb(G4int index, G4double energy)
void ClusterDecay(G4LorentzVector &lvX, G4int qX)
void FinalBarion(G4LorentzVector &lvB, G4int qB, G4int pdgB)
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
G4double AtomicMass(const G4double A, const G4double Z) const
Definition: G4Nucleus.cc:254
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
G4double energy(const ThreeVector &p, const G4double m)

◆ GetMinNuElEnergy()

G4double G4ANuElNucleusNcModel::GetMinNuElEnergy ( )
inline

Definition at line 81 of file G4ANuElNucleusNcModel.hh.

81{ return fMnumu + 0.5*fMnumu*fMnumu/fM1 + 0.05*CLHEP::keV; }; // kinematics + accuracy for sqrts

Referenced by IsApplicable().

◆ InitialiseModel()

void G4ANuElNucleusNcModel::InitialiseModel ( )
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 101 of file G4ANuElNucleusNcModel.cc.

102{
103 G4String pName = "anti_nu_e";
104
105 G4int nSize(0), i(0), j(0), k(0);
106
107 if(!fData)
108 {
109#ifdef G4MULTITHREADED
110 G4MUTEXLOCK(&numuNucleusModel);
111 if(!fData)
112 {
113#endif
114 fMaster = true;
115#ifdef G4MULTITHREADED
116 }
117 G4MUTEXUNLOCK(&numuNucleusModel);
118#endif
119 }
120
121 if(fMaster)
122 {
123 char* path = getenv("G4PARTICLEXSDATA");
124 std::ostringstream ost1, ost2, ost3, ost4;
125 ost1 << path << "/" << "neutrino" << "/" << pName << "/xarraynckr";
126
127 std::ifstream filein1( ost1.str().c_str() );
128
129 // filein.open("$PARTICLEXSDATA/");
130
131 filein1>>nSize;
132
133 for( k = 0; k < fNbin; ++k )
134 {
135 for( i = 0; i <= fNbin; ++i )
136 {
137 filein1 >> fNuMuXarrayKR[k][i];
138 // G4cout<< fNuMuXarrayKR[k][i] << " ";
139 }
140 }
141 // G4cout<<G4endl<<G4endl;
142
143 ost2 << path << "/" << "neutrino" << "/" << pName << "/xdistrnckr";
144 std::ifstream filein2( ost2.str().c_str() );
145
146 filein2>>nSize;
147
148 for( k = 0; k < fNbin; ++k )
149 {
150 for( i = 0; i < fNbin; ++i )
151 {
152 filein2 >> fNuMuXdistrKR[k][i];
153 // G4cout<< fNuMuXdistrKR[k][i] << " ";
154 }
155 }
156 // G4cout<<G4endl<<G4endl;
157
158 ost3 << path << "/" << "neutrino" << "/" << pName << "/q2arraynckr";
159 std::ifstream filein3( ost3.str().c_str() );
160
161 filein3>>nSize;
162
163 for( k = 0; k < fNbin; ++k )
164 {
165 for( i = 0; i <= fNbin; ++i )
166 {
167 for( j = 0; j <= fNbin; ++j )
168 {
169 filein3 >> fNuMuQarrayKR[k][i][j];
170 // G4cout<< fNuMuQarrayKR[k][i][j] << " ";
171 }
172 }
173 }
174 // G4cout<<G4endl<<G4endl;
175
176 ost4 << path << "/" << "neutrino" << "/" << pName << "/q2distrnckr";
177 std::ifstream filein4( ost4.str().c_str() );
178
179 filein4>>nSize;
180
181 for( k = 0; k < fNbin; ++k )
182 {
183 for( i = 0; i <= fNbin; ++i )
184 {
185 for( j = 0; j < fNbin; ++j )
186 {
187 filein4 >> fNuMuQdistrKR[k][i][j];
188 // G4cout<< fNuMuQdistrKR[k][i][j] << " ";
189 }
190 }
191 }
192 fData = true;
193 }
194}
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
static G4double fNuMuQarrayKR[50][51][51]
static G4double fNuMuXarrayKR[50][51]
static G4double fNuMuXdistrKR[50][50]
static G4double fNuMuQdistrKR[50][51][50]

Referenced by G4ANuElNucleusNcModel().

◆ IsApplicable()

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

Reimplemented from G4NeutrinoNucleusModel.

Definition at line 198 of file G4ANuElNucleusNcModel.cc.

200{
201 G4bool result = false;
202 G4String pName = aPart.GetDefinition()->GetParticleName();
203 G4double energy = aPart.GetTotalEnergy();
205
206 if( pName == "anti_nu_e"
207 &&
208 energy > fMinNuEnergy )
209 {
210 result = true;
211 }
212 G4int Z = targetNucleus.GetZ_asInt();
213 Z *= 1;
214
215 return result;
216}
bool G4bool
Definition: G4Types.hh:86

◆ ModelDescription()

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

Reimplemented from G4NeutrinoNucleusModel.

Definition at line 88 of file G4ANuElNucleusNcModel.cc.

89{
90
91 outFile << "G4ANuElNucleusNcModel is a neutrino-nucleus (neutral current) scattering\n"
92 << "model which uses the standard model \n"
93 << "transfer parameterization. The model is fully relativistic\n";
94
95}

◆ SampleLVkr()

void G4ANuElNucleusNcModel::SampleLVkr ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)

Definition at line 467 of file G4ANuElNucleusNcModel.cc.

468{
469 fBreak = false;
470 G4int A = targetNucleus.GetA_asInt(), iTer(0), iTerMax(100);
471 G4int Z = targetNucleus.GetZ_asInt();
472 G4double e3(0.), pMu2(0.), pX2(0.), nMom(0.), rM(0.), hM(0.), tM = targetNucleus.AtomicMass(A,Z);
473 G4double cost(1.), sint(0.), phi(0.), muMom(0.);
474 G4ThreeVector eP, bst;
475 const G4HadProjectile* aParticle = &aTrack;
476 G4LorentzVector lvp1 = aParticle->Get4Momentum();
477 nMom = NucleonMomentum( targetNucleus );
478
479 if( A == 1 || nMom == 0. ) // hydrogen, no Fermi motion ???
480 {
481 fNuEnergy = aParticle->GetTotalEnergy();
482 iTer = 0;
483
484 do
485 {
489
490 if( fXsample > 0. )
491 {
492 fW2 = fM1*fM1 - fQ2 + fQ2/fXsample; // sample excited hadron mass
494 }
495 else
496 {
497 fW2 = fM1*fM1;
498 fEmu = fNuEnergy;
499 }
500 e3 = fNuEnergy + fM1 - fEmu;
501
502 // if( e3 < sqrt(fW2) ) G4cout<<"energyX = "<<e3/GeV<<", fW = "<<sqrt(fW2)/GeV<<G4endl; // vmg ~10^-5 for NC
503
504 pMu2 = fEmu*fEmu - fMnumu*fMnumu;
505 pX2 = e3*e3 - fW2;
506
507 fCosTheta = fNuEnergy*fNuEnergy + pMu2 - pX2;
508 fCosTheta /= 2.*fNuEnergy*sqrt(pMu2);
509 iTer++;
510 }
511 while( ( abs(fCosTheta) > 1. || fEmu < fMnumu ) && iTer < iTerMax );
512
513 if( iTer >= iTerMax ) { fBreak = true; return; }
514
515 if( abs(fCosTheta) > 1.) // vmg: due to big Q2/x values. To be improved ...
516 {
517 G4cout<<"H2: fCosTheta = "<<fCosTheta<<", fEmu = "<<fEmu<<G4endl;
518 // fCosTheta = -1. + 2.*G4UniformRand();
519 if(fCosTheta < -1.) fCosTheta = -1.;
520 if(fCosTheta > 1.) fCosTheta = 1.;
521 }
522 // LVs
523
524 G4LorentzVector lvt1 = G4LorentzVector( 0., 0., 0., fM1 );
525 G4LorentzVector lvsum = lvp1 + lvt1;
526
527 cost = fCosTheta;
528 sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
529 phi = G4UniformRand()*CLHEP::twopi;
530 eP = G4ThreeVector( sint*std::cos(phi), sint*std::sin(phi), cost );
531 muMom = sqrt(fEmu*fEmu-fMnumu*fMnumu);
532 eP *= muMom;
533 fLVl = G4LorentzVector( eP, fEmu );
534
535 fLVh = lvsum - fLVl;
536 fLVt = G4LorentzVector( 0., 0., 0., 0. ); // no recoil
537 }
538 else // Fermi motion, Q2 in nucleon rest frame
539 {
540 G4ThreeVector nMomDir = nMom*G4RandomDirection();
541
542 if( !f2p2h ) // 1p1h
543 {
544 G4Nucleus recoil(A-1,Z);
545 rM = sqrt( recoil.AtomicMass(A-1,Z)*recoil.AtomicMass(A-1,Z) + nMom*nMom );
546 hM = tM - rM;
547
548 fLVt = G4LorentzVector( nMomDir, sqrt( rM*rM+nMom*nMom ) );
549 fLVh = G4LorentzVector(-nMomDir, sqrt( hM*hM+nMom*nMom ) );
550 }
551 else // 2p2h
552 {
553 G4Nucleus recoil(A-2,Z-1);
554 rM = recoil.AtomicMass(A-2,Z-1)+sqrt(nMom*nMom+fM1*fM1);
555 hM = tM - rM;
556
557 fLVt = G4LorentzVector( nMomDir, sqrt( rM*rM+nMom*nMom ) );
558 fLVh = G4LorentzVector(-nMomDir, sqrt( hM*hM+nMom*nMom ) );
559 }
560 // G4cout<<hM<<", ";
561 // bst = fLVh.boostVector(); // 9-3-20
562
563 // lvp1.boost(-bst); // 9-3-20 -> nucleon rest system, where Q2 transfer is ???
564
565 fNuEnergy = lvp1.e();
566 iTer = 0;
567
568 do
569 {
573
574 if( fXsample > 0. )
575 {
576 fW2 = fM1*fM1 - fQ2 + fQ2/fXsample; // sample excited hadron mass
578 }
579 else
580 {
581 fW2 = fM1*fM1;
582 fEmu = fNuEnergy;
583 }
584
585 // if(fEmu < 0.) G4cout<<"fEmu = "<<fEmu<<" hM = "<<hM<<G4endl;
586
587 e3 = fNuEnergy + fM1 - fEmu;
588
589 // if( e3 < sqrt(fW2) ) G4cout<<"energyX = "<<e3/GeV<<", fW = "<<sqrt(fW2)/GeV<<G4endl;
590
591 pMu2 = fEmu*fEmu - fMnumu*fMnumu;
592 pX2 = e3*e3 - fW2;
593
594 fCosTheta = fNuEnergy*fNuEnergy + pMu2 - pX2;
595 fCosTheta /= 2.*fNuEnergy*sqrt(pMu2);
596 iTer++;
597 }
598 while( ( abs(fCosTheta) > 1. || fEmu < fMnumu ) && iTer < iTerMax );
599
600 if( iTer >= iTerMax ) { fBreak = true; return; }
601
602 if( abs(fCosTheta) > 1.) // vmg: due to big Q2/x values. To be improved ...
603 {
604 G4cout<<"FM: fCosTheta = "<<fCosTheta<<", fEmu = "<<fEmu<<G4endl;
605 // fCosTheta = -1. + 2.*G4UniformRand();
606 if(fCosTheta < -1.) fCosTheta = -1.;
607 if(fCosTheta > 1.) fCosTheta = 1.;
608 }
609 // LVs
610 G4LorentzVector lvt1 = G4LorentzVector( 0., 0., 0., fM1 );
611 G4LorentzVector lvsum = lvp1 + lvt1;
612
613 cost = fCosTheta;
614 sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
615 phi = G4UniformRand()*CLHEP::twopi;
616 eP = G4ThreeVector( sint*std::cos(phi), sint*std::sin(phi), cost );
617 muMom = sqrt(fEmu*fEmu-fMnumu*fMnumu);
618 eP *= muMom;
619 fLVl = G4LorentzVector( eP, fEmu );
620 fLVh = lvsum - fLVl;
621 // back to lab system
622 // fLVl.boost(bst); // 9-3-20
623 // fLVh.boost(bst); // 9-3-20
624 }
625 //G4cout<<iTer<<", "<<fBreak<<"; ";
626}
G4ThreeVector G4RandomDirection()
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double NucleonMomentum(G4Nucleus &targetNucleus)
G4double SampleXkr(G4double energy)
G4double SampleQkr(G4double energy, G4double xx)

Referenced by ApplyYourself().

◆ ThresholdEnergy()

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

Definition at line 83 of file G4ANuElNucleusNcModel.hh.

84 {
85 G4double w = std::sqrt(fW2);
86 return w + 0.5*( (mP+mF)*(mP+mF)-(w+mI)*(w+mI) )/mI;
87 };

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