Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ANuMuNucleusCcModel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id: G4ANuMuNucleusCcModel.cc 91806 2015-08-06 12:20:45Z gcosmo $
27//
28// Geant4 Header : G4ANuMuNucleusCcModel
29//
30// Author : V.Grichine 12.2.19
31//
32
33#include <iostream>
34#include <fstream>
35#include <sstream>
36
38// #include "G4NuMuNuclCcDistrKR.hh"
39
40// #include "G4NuMuResQX.hh"
41
42#include "G4SystemOfUnits.hh"
43#include "G4ParticleTable.hh"
45#include "G4IonTable.hh"
46#include "Randomize.hh"
47#include "G4RandomDirection.hh"
48// #include "G4Threading.hh"
49
50// #include "G4Integrator.hh"
51#include "G4DataVector.hh"
52#include "G4PhysicsTable.hh"
53/*
54#include "G4CascadeInterface.hh"
55// #include "G4BinaryCascade.hh"
56#include "G4TheoFSGenerator.hh"
57#include "G4LundStringFragmentation.hh"
58#include "G4ExcitedStringDecay.hh"
59#include "G4FTFModel.hh"
60// #include "G4BinaryCascade.hh"
61#include "G4HadFinalState.hh"
62#include "G4HadSecondary.hh"
63#include "G4HadronicInteractionRegistry.hh"
64// #include "G4INCLXXInterface.hh"
65#include "G4QGSModel.hh"
66#include "G4QGSMFragmentation.hh"
67#include "G4QGSParticipants.hh"
68*/
69#include "G4KineticTrack.hh"
72#include "G4Fragment.hh"
73#include "G4NucleiProperties.hh"
75
77#include "G4PreCompoundModel.hh"
79
80
81// #include "G4MuonMinus.hh"
82#include "G4MuonPlus.hh"
83#include "G4Nucleus.hh"
84#include "G4LorentzVector.hh"
85
86using namespace std;
87using namespace CLHEP;
88
89#ifdef G4MULTITHREADED
90 G4Mutex G4ANuMuNucleusCcModel::numuNucleusModel = G4MUTEX_INITIALIZER;
91#endif
92
93
96{
97 fData = fMaster = false;
99}
100
101
103{}
104
105
106void G4ANuMuNucleusCcModel::ModelDescription(std::ostream& outFile) const
107{
108
109 outFile << "G4ANuMuNucleusCcModel is a neutrino-nucleus (charge current) scattering\n"
110 << "model which uses the standard model \n"
111 << "transfer parameterization. The model is fully relativistic\n";
112
113}
114
115/////////////////////////////////////////////////////////
116//
117// Read data from G4PARTICLEXSDATA (locally PARTICLEXSDATA)
118
120{
121 G4String pName = "anti_nu_mu";
122
123 G4int nSize(0), i(0), j(0), k(0);
124
125 if(!fData)
126 {
127#ifdef G4MULTITHREADED
128 G4MUTEXLOCK(&numuNucleusModel);
129 if(!fData)
130 {
131#endif
132 fMaster = true;
133#ifdef G4MULTITHREADED
134 }
135 G4MUTEXUNLOCK(&numuNucleusModel);
136#endif
137 }
138
139 if(fMaster)
140 {
141 const char* path = G4FindDataDir("G4PARTICLEXSDATA");
142 std::ostringstream ost1, ost2, ost3, ost4;
143 ost1 << path << "/" << "neutrino" << "/" << pName << "/xarraycckr";
144
145 std::ifstream filein1( ost1.str().c_str() );
146
147 // filein.open("$PARTICLEXSDATA/");
148
149 filein1>>nSize;
150
151 for( k = 0; k < fNbin; ++k )
152 {
153 for( i = 0; i <= fNbin; ++i )
154 {
155 filein1 >> fNuMuXarrayKR[k][i];
156 // G4cout<< fNuMuXarrayKR[k][i] << " ";
157 }
158 }
159 // G4cout<<G4endl<<G4endl;
160
161 ost2 << path << "/" << "neutrino" << "/" << pName << "/xdistrcckr";
162 std::ifstream filein2( ost2.str().c_str() );
163
164 filein2>>nSize;
165
166 for( k = 0; k < fNbin; ++k )
167 {
168 for( i = 0; i < fNbin; ++i )
169 {
170 filein2 >> fNuMuXdistrKR[k][i];
171 // G4cout<< fNuMuXdistrKR[k][i] << " ";
172 }
173 }
174 // G4cout<<G4endl<<G4endl;
175
176 ost3 << path << "/" << "neutrino" << "/" << pName << "/q2arraycckr";
177 std::ifstream filein3( ost3.str().c_str() );
178
179 filein3>>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 filein3 >> fNuMuQarrayKR[k][i][j];
188 // G4cout<< fNuMuQarrayKR[k][i][j] << " ";
189 }
190 }
191 }
192 // G4cout<<G4endl<<G4endl;
193
194 ost4 << path << "/" << "neutrino" << "/" << pName << "/q2distrcckr";
195 std::ifstream filein4( ost4.str().c_str() );
196
197 filein4>>nSize;
198
199 for( k = 0; k < fNbin; ++k )
200 {
201 for( i = 0; i <= fNbin; ++i )
202 {
203 for( j = 0; j < fNbin; ++j )
204 {
205 filein4 >> fNuMuQdistrKR[k][i][j];
206 // G4cout<< fNuMuQdistrKR[k][i][j] << " ";
207 }
208 }
209 }
210 fData = true;
211 }
212}
213
214/////////////////////////////////////////////////////////
215
217 G4Nucleus & )
218{
219 G4bool result = false;
220 G4String pName = aPart.GetDefinition()->GetParticleName();
221 G4double energy = aPart.GetTotalEnergy();
222
223 if( pName == "anti_nu_mu"
224 &&
225 energy > fMinNuEnergy )
226 {
227 result = true;
228 }
229
230 return result;
231}
232
233/////////////////////////////////////////// ClusterDecay ////////////////////////////////////////////////////////////
234//
235//
236
238 const G4HadProjectile& aTrack, G4Nucleus& targetNucleus)
239{
241 fProton = f2p2h = fBreak = false;
242 fCascade = fString = false;
243 fLVh = fLVl = fLVt = fLVcpi = G4LorentzVector(0.,0.,0.,0.);
244
245 const G4HadProjectile* aParticle = &aTrack;
246 G4double energy = aParticle->GetTotalEnergy();
247
248 G4String pName = aParticle->GetDefinition()->GetParticleName();
249
250 if( energy < fMinNuEnergy )
251 {
254 return &theParticleChange;
255 }
256
257 SampleLVkr( aTrack, targetNucleus);
258
259 if( fBreak == true || fEmu < fMu ) // ~5*10^-6
260 {
261 // G4cout<<"ni, ";
264 return &theParticleChange;
265 }
266
267 // LVs of initial state
268
269 G4LorentzVector lvp1 = aParticle->Get4Momentum();
270 G4LorentzVector lvt1( 0., 0., 0., fM1 );
272
273 // 1-pi by fQtransfer && nu-energy
274 G4LorentzVector lvpip1( 0., 0., 0., mPip );
275 G4LorentzVector lvsum, lv2, lvX;
276 G4ThreeVector eP;
277 G4double cost(1.), sint(0.), phi(0.), muMom(0.), massX2(0.), massX(0.), massR(0.), eCut(0.);
278 G4DynamicParticle* aLept = nullptr; // lepton lv
279
280 G4int Z = targetNucleus.GetZ_asInt();
281 G4int A = targetNucleus.GetA_asInt();
282 G4double mTarg = targetNucleus.AtomicMass(A,Z);
283 G4int pdgP(0), qB(0);
284 // G4double mSum = G4ParticleTable::GetParticleTable()->FindParticle(2212)->GetPDGMass() + mPip;
285
286 G4int iPi = GetOnePionIndex(energy);
287 G4double p1pi = GetNuMuOnePionProb( iPi, energy);
288
289 if( p1pi > G4UniformRand() && fCosTheta > 0.9 ) // && fQtransfer < 0.95*GeV ) // mu- & coherent pion + nucleus
290 {
291 // lvsum = lvp1 + lvpip1;
292 lvsum = lvp1 + lvt1;
293 // cost = fCosThetaPi;
294 cost = fCosTheta;
295 sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
296 phi = G4UniformRand()*CLHEP::twopi;
297 eP = G4ThreeVector( sint*std::cos(phi), sint*std::sin(phi), cost );
298
299 // muMom = sqrt(fEmuPi*fEmuPi-fMu*fMu);
300 muMom = sqrt(fEmu*fEmu-fMu*fMu);
301
302 eP *= muMom;
303
304 // lv2 = G4LorentzVector( eP, fEmuPi );
305 // lv2 = G4LorentzVector( eP, fEmu );
306 lv2 = fLVl;
307
308 // lvX = lvsum - lv2;
309 lvX = fLVh;
310 massX2 = lvX.m2();
311 massX = lvX.m();
312 massR = fLVt.m();
313
314 if ( massX2 <= 0. ) // vmg: very rarely ~ (1-4)e-6 due to big Q2/x, to be improved
315 {
316 fCascade = true;
319 return &theParticleChange;
320 }
321 fW2 = massX2;
322
323 if( pName == "anti_nu_mu") aLept = new G4DynamicParticle( theMuonPlus, lv2 );
324 else
325 {
328 return &theParticleChange;
329 }
330 if( pName == "anti_nu_mu" ) pdgP = -211;
331 // else pdgP = -211;
332 // eCut = fMpi + 0.5*(fMpi*fMpi-massX2)/mTarg; // massX -> fMpi
333
334 if( A > 1 )
335 {
336 eCut = (fMpi + mTarg)*(fMpi + mTarg) - (massX + massR)*(massX + massR);
337 eCut /= 2.*massR;
338 eCut += massX;
339 }
340 else eCut = fM1 + fMpi;
341
342 if ( lvX.e() > eCut ) // && sqrt( GetW2() ) < 1.4*GeV ) //
343 {
344 CoherentPion( lvX, pdgP, targetNucleus);
345 }
346 else
347 {
348 fCascade = true;
351 return &theParticleChange;
352 }
354
355 return &theParticleChange;
356 }
357 else // lepton part in lab
358 {
359 lvsum = lvp1 + lvt1;
360 cost = fCosTheta;
361 sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
362 phi = G4UniformRand()*CLHEP::twopi;
363 eP = G4ThreeVector( sint*std::cos(phi), sint*std::sin(phi), cost );
364
365 muMom = sqrt(fEmu*fEmu-fMu*fMu);
366
367 eP *= muMom;
368
369 lv2 = G4LorentzVector( eP, fEmu );
370 lv2 = fLVl;
371 lvX = lvsum - lv2;
372 lvX = fLVh;
373 massX2 = lvX.m2();
374
375 if ( massX2 <= 0. ) // vmg: very rarely ~ (1-4)e-6 due to big Q2/x, to be improved
376 {
377 fCascade = true;
380 return &theParticleChange;
381 }
382 fW2 = massX2;
383
384 if( pName == "anti_nu_mu") aLept = new G4DynamicParticle( theMuonPlus, lv2 );
385 else
386 {
389 return &theParticleChange;
390 }
392 }
393
394 // hadron part
395
396 fRecoil = nullptr;
397
398 if( A == 1 )
399 {
400 if( pName == "anti_nu_mu" ) qB = 2;
401 // else qB = 0;
402
403 // if( G4UniformRand() > 0.1 ) // > 0.9999 ) // > 0.0001 ) //
404 {
405 ClusterDecay( lvX, qB );
406 }
407 return &theParticleChange;
408 }
409 /*
410 // else
411 {
412 if( pName == "nu_mu" ) pdgP = 211;
413 else pdgP = -211;
414
415
416 if ( fQtransfer < 0.95*GeV ) // < 0.35*GeV ) //
417 {
418 if( lvX.m() > mSum ) CoherentPion( lvX, pdgP, targetNucleus);
419 }
420 }
421 return &theParticleChange;
422 }
423 */
424 G4Nucleus recoil;
425 G4double rM(0.), ratio = G4double(Z)/G4double(A);
426
427 if( ratio > G4UniformRand() ) // proton is excited
428 {
429 fProton = true;
430 recoil = G4Nucleus(A-1,Z-1);
431 fRecoil = &recoil;
432 rM = recoil.AtomicMass(A-1,Z-1);
433
434 if( pName == "anti_nu_mu" ) // (0) state -> p + pi-
435 {
438 }
439 else // (0) state -> p + pi-, n + pi0
440 {
441 // fMt = G4ParticleTable::GetParticleTable()->FindParticle(2212)->GetPDGMass()
442 // + G4ParticleTable::GetParticleTable()->FindParticle(-211)->GetPDGMass();
443 }
444 }
445 else // excited neutron
446 {
447 fProton = false;
448 recoil = G4Nucleus(A-1,Z);
449 fRecoil = &recoil;
450 rM = recoil.AtomicMass(A-1,Z);
451
452 if( pName == "anti_nu_mu" ) // (+) state -> n + pi+
453 {
456 }
457 else // (-) state -> n + pi-, // n + pi0
458 {
459 // fMt = G4ParticleTable::GetParticleTable()->FindParticle(2112)->GetPDGMass()
460 // + G4ParticleTable::GetParticleTable()->FindParticle(-211)->GetPDGMass();
461 }
462 }
463 // G4int index = GetEnergyIndex(energy);
464 G4int nepdg = aParticle->GetDefinition()->GetPDGEncoding();
465
466 G4double qeTotRat; // = GetNuMuQeTotRat(index, energy);
467 qeTotRat = CalculateQEratioA( Z, A, energy, nepdg);
468
469 G4ThreeVector dX = (lvX.vect()).unit();
470 G4double eX = lvX.e(); // excited nucleon
471 G4double mX = sqrt(massX2);
472 // G4double pX = sqrt( eX*eX - mX*mX );
473 // G4double sumE = eX + rM;
474
475 if( qeTotRat > G4UniformRand() || mX <= fMt ) // || eX <= 1232.*MeV) // QE
476 {
477 fString = false;
478
479 if( fProton )
480 {
481 fPDGencoding = 2212;
482 fMr = proton_mass_c2;
483 recoil = G4Nucleus(A-1,Z-1);
484 fRecoil = &recoil;
485 rM = recoil.AtomicMass(A-1,Z-1);
486 }
487 else // if( pName == "anti_nu_mu" )
488 {
489 fPDGencoding = 2112;
491 FindParticle(fPDGencoding)->GetPDGMass(); // 939.5654133*MeV;
492 recoil = G4Nucleus(A-1,Z);
493 fRecoil = &recoil;
494 rM = recoil.AtomicMass(A-1,Z);
495 }
496 // sumE = eX + rM;
497 G4double eTh = fMr + 0.5*(fMr*fMr - mX*mX)/rM;
498
499 if( eX <= eTh ) // vmg, very rarely out of kinematics
500 {
501 fString = true;
504 return &theParticleChange;
505 }
506 // FinalBarion( fLVh, 0, fPDGencoding ); // p(n)+deexcited recoil
507 FinalBarion( lvX, 0, fPDGencoding ); // p(n)+deexcited recoil
508 }
509 else // if ( eX < 9500000.*GeV ) // < 25.*GeV) // < 95.*GeV ) // < 2.5*GeV ) //cluster decay
510 {
511 if ( fProton && pName == "anti_nu_mu" ) qB = 0;
512 else if( !fProton && pName == "anti_nu_mu" ) qB = -1;
513
514 ClusterDecay( lvX, qB );
515 }
516 return &theParticleChange;
517}
518
519
520/////////////////////////////////////////////////////////////////////
521////////////////////////////////////////////////////////////////////
522///////////////////////////////////////////////////////////////////
523
524/////////////////////////////////////////////////
525//
526// sample x, then Q2
527
529{
530 fBreak = false;
531 G4int A = targetNucleus.GetA_asInt(), iTer(0), iTerMax(100);
532 G4int Z = targetNucleus.GetZ_asInt();
533 G4double e3(0.), pMu2(0.), pX2(0.), nMom(0.), rM(0.), hM(0.), tM = targetNucleus.AtomicMass(A,Z);
534 G4double Ex(0.), ei(0.), nm2(0.);
535 G4double cost(1.), sint(0.), phi(0.), muMom(0.);
536 G4ThreeVector eP, bst;
537 const G4HadProjectile* aParticle = &aTrack;
538 G4LorentzVector lvp1 = aParticle->Get4Momentum();
539
540 if( A == 1 ) // hydrogen, no Fermi motion ???
541 {
542 fNuEnergy = aParticle->GetTotalEnergy();
543 iTer = 0;
544
545 do
546 {
550
551 if( fXsample > 0. )
552 {
553 fW2 = fM1*fM1 - fQ2 + fQ2/fXsample; // sample excited hadron mass
555 }
556 else
557 {
558 fW2 = fM1*fM1;
559 fEmu = fNuEnergy;
560 }
561 e3 = fNuEnergy + fM1 - fEmu;
562
563 if( e3 < sqrt(fW2) ) G4cout<<"energyX = "<<e3/GeV<<", fW = "<<sqrt(fW2)/GeV<<G4endl;
564
565 pMu2 = fEmu*fEmu - fMu*fMu;
566
567 if(pMu2 < 0.) { fBreak = true; return; }
568
569 pX2 = e3*e3 - fW2;
570
571 fCosTheta = fNuEnergy*fNuEnergy + pMu2 - pX2;
572 fCosTheta /= 2.*fNuEnergy*sqrt(pMu2);
573 iTer++;
574 }
575 while( ( abs(fCosTheta) > 1. || fEmu < fMu ) && iTer < iTerMax );
576
577 if( iTer >= iTerMax ) { fBreak = true; return; }
578
579 if( abs(fCosTheta) > 1.) // vmg: due to big Q2/x values. To be improved ...
580 {
581 G4cout<<"H2: fCosTheta = "<<fCosTheta<<", fEmu = "<<fEmu<<G4endl;
582 // fCosTheta = -1. + 2.*G4UniformRand();
583 if(fCosTheta < -1.) fCosTheta = -1.;
584 if(fCosTheta > 1.) fCosTheta = 1.;
585 }
586 // LVs
587
588 G4LorentzVector lvt1 = G4LorentzVector( 0., 0., 0., fM1 );
589 G4LorentzVector lvsum = lvp1 + lvt1;
590
591 cost = fCosTheta;
592 sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
593 phi = G4UniformRand()*CLHEP::twopi;
594 eP = G4ThreeVector( sint*std::cos(phi), sint*std::sin(phi), cost );
595 muMom = sqrt(fEmu*fEmu-fMu*fMu);
596 eP *= muMom;
597 fLVl = G4LorentzVector( eP, fEmu );
598
599 fLVh = lvsum - fLVl;
600 fLVt = G4LorentzVector( 0., 0., 0., 0. ); // no recoil
601 }
602 else // Fermi motion, Q2 in nucleon rest frame
603 {
604 G4Nucleus recoil1( A-1, Z );
605 rM = recoil1.AtomicMass(A-1,Z);
606 do
607 {
608 // nMom = NucleonMomentumBR( targetNucleus ); // BR
609 nMom = GgSampleNM( targetNucleus ); // Gg
610 Ex = GetEx(A-1, fProton);
611 ei = tM - sqrt( (rM + Ex)*(rM + Ex) + nMom*nMom );
612 // ei = 0.5*( tM - s2M - 2*eX );
613
614 nm2 = ei*ei - nMom*nMom;
615 iTer++;
616 }
617 while( nm2 < 0. && iTer < iTerMax );
618
619 if( iTer >= iTerMax ) { fBreak = true; return; }
620
621 G4ThreeVector nMomDir = nMom*G4RandomDirection();
622
623 if( !f2p2h || A < 3 ) // 1p1h
624 {
625 // hM = tM - rM;
626
627 fLVt = G4LorentzVector( -nMomDir, sqrt( (rM + Ex)*(rM + Ex) + nMom*nMom ) ); // rM ); //
628 fLVh = G4LorentzVector( nMomDir, ei ); // hM); //
629 }
630 else // 2p2h
631 {
632 G4Nucleus recoil(A-2,Z-1);
633 rM = recoil.AtomicMass(A-2,Z-1)+sqrt(nMom*nMom+fM1*fM1);
634 hM = tM - rM;
635
636 fLVt = G4LorentzVector( nMomDir, sqrt( rM*rM+nMom*nMom ) );
637 fLVh = G4LorentzVector(-nMomDir, sqrt( hM*hM+nMom*nMom ) );
638 }
639 // G4cout<<hM<<", ";
640 // bst = fLVh.boostVector();
641
642 // lvp1.boost(-bst); // -> nucleon rest system, where Q2 transfer is ???
643
644 fNuEnergy = lvp1.e();
645 // G4double mN = fLVh.m(); // better mN = fM1 !? vmg
646 iTer = 0;
647
648 do // no FM!?, 5.4.20 vmg
649 {
653
654 // G4double mR = mN + fM1*(A-1.)*std::exp(-2.0*fQtransfer/mN); // recoil mass in+el
655
656 if( fXsample > 0. )
657 {
658 fW2 = fM1*fM1 - fQ2 + fQ2/fXsample; // sample excited hadron mass
659
660 // fW2 = mN*mN - fQ2 + fQ2/fXsample; // sample excited hadron mass
661 // fEmu = fNuEnergy - fQ2/2./mR/fXsample; // fM1->mN
662
663 fEmu = fNuEnergy - fQ2/2./fM1/fXsample; // fM1->mN
664 }
665 else
666 {
667 // fW2 = mN*mN;
668
669 fW2 = fM1*fM1;
670 fEmu = fNuEnergy;
671 }
672 // if(fEmu < 0.) G4cout<<"fEmu = "<<fEmu<<" hM = "<<hM<<G4endl;
673 // e3 = fNuEnergy + mR - fEmu;
674
675 e3 = fNuEnergy + fM1 - fEmu;
676
677 // if( e3 < sqrt(fW2) ) G4cout<<"energyX = "<<e3/GeV<<", fW = "<<sqrt(fW2)/GeV<<G4endl;
678
679 pMu2 = fEmu*fEmu - fMu*fMu;
680 pX2 = e3*e3 - fW2;
681
682 if(pMu2 < 0.) { fBreak = true; return; }
683
684 fCosTheta = fNuEnergy*fNuEnergy + pMu2 - pX2;
685 fCosTheta /= 2.*fNuEnergy*sqrt(pMu2);
686 iTer++;
687 }
688 while( ( abs(fCosTheta) > 1. || fEmu < fMu ) && iTer < iTerMax );
689
690 if( iTer >= iTerMax ) { fBreak = true; return; }
691
692 if( abs(fCosTheta) > 1.) // vmg: due to big Q2/x values. To be improved ...
693 {
694 G4cout<<"FM: fCosTheta = "<<fCosTheta<<", fEmu = "<<fEmu<<G4endl;
695 // fCosTheta = -1. + 2.*G4UniformRand();
696 if( fCosTheta < -1.) fCosTheta = -1.;
697 if( fCosTheta > 1.) fCosTheta = 1.;
698 }
699 // LVs
700 // G4LorentzVector lvt1 = G4LorentzVector( 0., 0., 0., mN ); // fM1 );
701
702 G4LorentzVector lvt1 = G4LorentzVector( 0., 0., 0., fM1 ); // fM1 );
703 G4LorentzVector lvsum = lvp1 + lvt1;
704
705 cost = fCosTheta;
706 sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
707 phi = G4UniformRand()*CLHEP::twopi;
708 eP = G4ThreeVector( sint*std::cos(phi), sint*std::sin(phi), cost );
709 muMom = sqrt(fEmu*fEmu-fMu*fMu);
710 eP *= muMom;
711 fLVl = G4LorentzVector( eP, fEmu );
712 fLVh = lvsum - fLVl;
713
714 // if( fLVh.e() < mN || fLVh.m2() < 0.) { fBreak = true; return; }
715
716 if( fLVh.e() < fM1 || fLVh.m2() < 0.) { fBreak = true; return; }
717
718 // back to lab system
719
720 // fLVl.boost(bst);
721 // fLVh.boost(bst);
722 }
723 //G4cout<<iTer<<", "<<fBreak<<"; ";
724}
725
726//
727//
728///////////////////////////
const char * G4FindDataDir(const char *)
CLHEP::HepLorentzVector G4LorentzVector
G4ThreeVector G4RandomDirection()
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
std::mutex G4Mutex
Definition: G4Threading.hh:81
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
Hep3Vector unit() const
Hep3Vector vect() const
virtual void ModelDescription(std::ostream &) const
virtual G4bool IsApplicable(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
void SampleLVkr(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
G4ANuMuNucleusCcModel(const G4String &name="ANuMuNuclCcModel")
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)
static G4double fNuMuQarrayKR[50][51][51]
static G4double fNuMuXarrayKR[50][51]
G4int GetOnePionIndex(G4double energy)
G4double SampleXkr(G4double energy)
G4double SampleQkr(G4double energy, G4double xx)
G4double GgSampleNM(G4Nucleus &nucl)
G4double GetNuMuOnePionProb(G4int index, G4double energy)
static G4double fNuMuXdistrKR[50][50]
static G4double fNuMuQdistrKR[50][51][50]
G4double CalculateQEratioA(G4int Z, G4int A, G4double energy, G4int nepdg)
G4ParticleDefinition * theMuonPlus
void ClusterDecay(G4LorentzVector &lvX, G4int qX)
void FinalBarion(G4LorentzVector &lvB, G4int qB, G4int pdgB)
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:357
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
Definition: DoubConv.h:17