Geant4 10.7.0
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 char* path = getenv("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 & targetNucleus)
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 G4int Z = targetNucleus.GetZ_asInt();
230 Z *= 1;
231
232 return result;
233}
234
235/////////////////////////////////////////// ClusterDecay ////////////////////////////////////////////////////////////
236//
237//
238
240 const G4HadProjectile& aTrack, G4Nucleus& targetNucleus)
241{
243 fProton = f2p2h = fBreak = false;
244 fCascade = fString = false;
245 fLVh = fLVl = fLVt = fLVcpi = G4LorentzVector(0.,0.,0.,0.);
246
247 const G4HadProjectile* aParticle = &aTrack;
248 G4double energy = aParticle->GetTotalEnergy();
249
250 G4String pName = aParticle->GetDefinition()->GetParticleName();
251
252 if( energy < fMinNuEnergy )
253 {
256 return &theParticleChange;
257 }
258
259 SampleLVkr( aTrack, targetNucleus);
260
261 if( fBreak == true || fEmu < fMu ) // ~5*10^-6
262 {
263 // G4cout<<"ni, ";
266 return &theParticleChange;
267 }
268
269 // LVs of initial state
270
271 G4LorentzVector lvp1 = aParticle->Get4Momentum();
272 G4LorentzVector lvt1( 0., 0., 0., fM1 );
274
275 // 1-pi by fQtransfer && nu-energy
276 G4LorentzVector lvpip1( 0., 0., 0., mPip );
277 G4LorentzVector lvsum, lv2, lvX;
278 G4ThreeVector eP;
279 G4double cost(1.), sint(0.), phi(0.), muMom(0.), massX2(0.), massX(0.), massR(0.), eCut(0.);
280 G4DynamicParticle* aLept = nullptr; // lepton lv
281
282 G4int Z = targetNucleus.GetZ_asInt();
283 G4int A = targetNucleus.GetA_asInt();
284 G4double mTarg = targetNucleus.AtomicMass(A,Z);
285 G4int pdgP(0), qB(0);
286 // G4double mSum = G4ParticleTable::GetParticleTable()->FindParticle(2212)->GetPDGMass() + mPip;
287
288 G4int iPi = GetOnePionIndex(energy);
289 G4double p1pi = GetNuMuOnePionProb( iPi, energy);
290
291 if( p1pi > G4UniformRand() && fCosTheta > 0.9 ) // && fQtransfer < 0.95*GeV ) // mu- & coherent pion + nucleus
292 {
293 // lvsum = lvp1 + lvpip1;
294 lvsum = lvp1 + lvt1;
295 // cost = fCosThetaPi;
296 cost = fCosTheta;
297 sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
298 phi = G4UniformRand()*CLHEP::twopi;
299 eP = G4ThreeVector( sint*std::cos(phi), sint*std::sin(phi), cost );
300
301 // muMom = sqrt(fEmuPi*fEmuPi-fMu*fMu);
302 muMom = sqrt(fEmu*fEmu-fMu*fMu);
303
304 eP *= muMom;
305
306 // lv2 = G4LorentzVector( eP, fEmuPi );
307 // lv2 = G4LorentzVector( eP, fEmu );
308 lv2 = fLVl;
309
310 // lvX = lvsum - lv2;
311 lvX = fLVh;
312 massX2 = lvX.m2();
313 massX = lvX.m();
314 massR = fLVt.m();
315
316 if ( massX2 <= 0. ) // vmg: very rarely ~ (1-4)e-6 due to big Q2/x, to be improved
317 {
318 fCascade = true;
321 return &theParticleChange;
322 }
323 fW2 = massX2;
324
325 if( pName == "anti_nu_mu") aLept = new G4DynamicParticle( theMuonPlus, lv2 );
326 else
327 {
330 return &theParticleChange;
331 }
332 if( pName == "anti_nu_mu" ) pdgP = -211;
333 // else pdgP = -211;
334 // eCut = fMpi + 0.5*(fMpi*fMpi-massX2)/mTarg; // massX -> fMpi
335
336 if( A > 1 )
337 {
338 eCut = (fMpi + mTarg)*(fMpi + mTarg) - (massX + massR)*(massX + massR);
339 eCut /= 2.*massR;
340 eCut += massX;
341 }
342 else eCut = fM1 + fMpi;
343
344 if ( lvX.e() > eCut ) // && sqrt( GetW2() ) < 1.4*GeV ) //
345 {
346 CoherentPion( lvX, pdgP, targetNucleus);
347 }
348 else
349 {
350 fCascade = true;
353 return &theParticleChange;
354 }
356
357 return &theParticleChange;
358 }
359 else // lepton part in lab
360 {
361 lvsum = lvp1 + lvt1;
362 cost = fCosTheta;
363 sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
364 phi = G4UniformRand()*CLHEP::twopi;
365 eP = G4ThreeVector( sint*std::cos(phi), sint*std::sin(phi), cost );
366
367 muMom = sqrt(fEmu*fEmu-fMu*fMu);
368
369 eP *= muMom;
370
371 lv2 = G4LorentzVector( eP, fEmu );
372 lv2 = fLVl;
373 lvX = lvsum - lv2;
374 lvX = fLVh;
375 massX2 = lvX.m2();
376
377 if ( massX2 <= 0. ) // vmg: very rarely ~ (1-4)e-6 due to big Q2/x, to be improved
378 {
379 fCascade = true;
382 return &theParticleChange;
383 }
384 fW2 = massX2;
385
386 if( pName == "anti_nu_mu") aLept = new G4DynamicParticle( theMuonPlus, lv2 );
387 else
388 {
391 return &theParticleChange;
392 }
394 }
395
396 // hadron part
397
398 fRecoil = nullptr;
399
400 if( A == 1 )
401 {
402 if( pName == "anti_nu_mu" ) qB = 2;
403 // else qB = 0;
404
405 // if( G4UniformRand() > 0.1 ) // > 0.9999 ) // > 0.0001 ) //
406 {
407 ClusterDecay( lvX, qB );
408 }
409 return &theParticleChange;
410 }
411 /*
412 // else
413 {
414 if( pName == "nu_mu" ) pdgP = 211;
415 else pdgP = -211;
416
417
418 if ( fQtransfer < 0.95*GeV ) // < 0.35*GeV ) //
419 {
420 if( lvX.m() > mSum ) CoherentPion( lvX, pdgP, targetNucleus);
421 }
422 }
423 return &theParticleChange;
424 }
425 */
426 G4Nucleus recoil;
427 G4double rM(0.), ratio = G4double(Z)/G4double(A);
428
429 if( ratio > G4UniformRand() ) // proton is excited
430 {
431 fProton = true;
432 recoil = G4Nucleus(A-1,Z-1);
433 fRecoil = &recoil;
434 rM = recoil.AtomicMass(A-1,Z-1);
435
436 if( pName == "anti_nu_mu" ) // (0) state -> p + pi-
437 {
440 }
441 else // (0) state -> p + pi-, n + pi0
442 {
443 // fMt = G4ParticleTable::GetParticleTable()->FindParticle(2212)->GetPDGMass()
444 // + G4ParticleTable::GetParticleTable()->FindParticle(-211)->GetPDGMass();
445 }
446 }
447 else // excited neutron
448 {
449 fProton = false;
450 recoil = G4Nucleus(A-1,Z);
451 fRecoil = &recoil;
452 rM = recoil.AtomicMass(A-1,Z);
453
454 if( pName == "anti_nu_mu" ) // (+) state -> n + pi+
455 {
458 }
459 else // (-) state -> n + pi-, // n + pi0
460 {
461 // fMt = G4ParticleTable::GetParticleTable()->FindParticle(2112)->GetPDGMass()
462 // + G4ParticleTable::GetParticleTable()->FindParticle(-211)->GetPDGMass();
463 }
464 }
465 G4int index = GetEnergyIndex(energy);
466 G4double qeTotRat = GetNuMuQeTotRat(index, energy);
467
468 G4ThreeVector dX = (lvX.vect()).unit();
469 G4double eX = lvX.e(); // excited nucleon
470 G4double mX = sqrt(massX2);
471 // G4double pX = sqrt( eX*eX - mX*mX );
472 // G4double sumE = eX + rM;
473
474 if( qeTotRat > G4UniformRand() || mX <= fMt ) // || eX <= 1232.*MeV) // QE
475 {
476 fString = false;
477
478 if( fProton )
479 {
480 fPDGencoding = 2212;
481 fMr = proton_mass_c2;
482 recoil = G4Nucleus(A-1,Z-1);
483 fRecoil = &recoil;
484 rM = recoil.AtomicMass(A-1,Z-1);
485 }
486 else // if( pName == "anti_nu_mu" )
487 {
488 fPDGencoding = 2112;
490 FindParticle(fPDGencoding)->GetPDGMass(); // 939.5654133*MeV;
491 recoil = G4Nucleus(A-1,Z);
492 fRecoil = &recoil;
493 rM = recoil.AtomicMass(A-1,Z);
494 }
495 // sumE = eX + rM;
496 G4double eTh = fMr + 0.5*(fMr*fMr - mX*mX)/rM;
497
498 if( eX <= eTh ) // vmg, very rarely out of kinematics
499 {
500 fString = true;
503 return &theParticleChange;
504 }
505 // FinalBarion( fLVh, 0, fPDGencoding ); // p(n)+deexcited recoil
506 FinalBarion( lvX, 0, fPDGencoding ); // p(n)+deexcited recoil
507 }
508 else // if ( eX < 9500000.*GeV ) // < 25.*GeV) // < 95.*GeV ) // < 2.5*GeV ) //cluster decay
509 {
510 if ( fProton && pName == "anti_nu_mu" ) qB = 0;
511 else if( !fProton && pName == "anti_nu_mu" ) qB = -1;
512
513 ClusterDecay( lvX, qB );
514 }
515 return &theParticleChange;
516}
517
518
519/////////////////////////////////////////////////////////////////////
520////////////////////////////////////////////////////////////////////
521///////////////////////////////////////////////////////////////////
522
523/////////////////////////////////////////////////
524//
525// sample x, then Q2
526
528{
529 fBreak = false;
530 G4int A = targetNucleus.GetA_asInt(), iTer(0), iTerMax(100);
531 G4int Z = targetNucleus.GetZ_asInt();
532 G4double e3(0.), pMu2(0.), pX2(0.), nMom(0.), rM(0.), hM(0.), tM = targetNucleus.AtomicMass(A,Z);
533 G4double Ex(0.), ei(0.), nm2(0.);
534 G4double cost(1.), sint(0.), phi(0.), muMom(0.);
535 G4ThreeVector eP, bst;
536 const G4HadProjectile* aParticle = &aTrack;
537 G4LorentzVector lvp1 = aParticle->Get4Momentum();
538
539 if( A == 1 ) // hydrogen, no Fermi motion ???
540 {
541 fNuEnergy = aParticle->GetTotalEnergy();
542 iTer = 0;
543
544 do
545 {
549
550 if( fXsample > 0. )
551 {
552 fW2 = fM1*fM1 - fQ2 + fQ2/fXsample; // sample excited hadron mass
554 }
555 else
556 {
557 fW2 = fM1*fM1;
558 fEmu = fNuEnergy;
559 }
560 e3 = fNuEnergy + fM1 - fEmu;
561
562 if( e3 < sqrt(fW2) ) G4cout<<"energyX = "<<e3/GeV<<", fW = "<<sqrt(fW2)/GeV<<G4endl;
563
564 pMu2 = fEmu*fEmu - fMu*fMu;
565
566 if(pMu2 < 0.) { fBreak = true; return; }
567
568 pX2 = e3*e3 - fW2;
569
570 fCosTheta = fNuEnergy*fNuEnergy + pMu2 - pX2;
571 fCosTheta /= 2.*fNuEnergy*sqrt(pMu2);
572 iTer++;
573 }
574 while( ( abs(fCosTheta) > 1. || fEmu < fMu ) && iTer < iTerMax );
575
576 if( iTer >= iTerMax ) { fBreak = true; return; }
577
578 if( abs(fCosTheta) > 1.) // vmg: due to big Q2/x values. To be improved ...
579 {
580 G4cout<<"H2: fCosTheta = "<<fCosTheta<<", fEmu = "<<fEmu<<G4endl;
581 // fCosTheta = -1. + 2.*G4UniformRand();
582 if(fCosTheta < -1.) fCosTheta = -1.;
583 if(fCosTheta > 1.) fCosTheta = 1.;
584 }
585 // LVs
586
587 G4LorentzVector lvt1 = G4LorentzVector( 0., 0., 0., fM1 );
588 G4LorentzVector lvsum = lvp1 + lvt1;
589
590 cost = fCosTheta;
591 sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
592 phi = G4UniformRand()*CLHEP::twopi;
593 eP = G4ThreeVector( sint*std::cos(phi), sint*std::sin(phi), cost );
594 muMom = sqrt(fEmu*fEmu-fMu*fMu);
595 eP *= muMom;
596 fLVl = G4LorentzVector( eP, fEmu );
597
598 fLVh = lvsum - fLVl;
599 fLVt = G4LorentzVector( 0., 0., 0., 0. ); // no recoil
600 }
601 else // Fermi motion, Q2 in nucleon rest frame
602 {
603 G4Nucleus recoil1( A-1, Z );
604 rM = recoil1.AtomicMass(A-1,Z);
605 do
606 {
607 // nMom = NucleonMomentumBR( targetNucleus ); // BR
608 nMom = GgSampleNM( targetNucleus ); // Gg
609 Ex = GetEx(A-1, fProton);
610 ei = tM - sqrt( (rM + Ex)*(rM + Ex) + nMom*nMom );
611 // ei = 0.5*( tM - s2M - 2*eX );
612
613 nm2 = ei*ei - nMom*nMom;
614 iTer++;
615 }
616 while( nm2 < 0. && iTer < iTerMax );
617
618 if( iTer >= iTerMax ) { fBreak = true; return; }
619
620 G4ThreeVector nMomDir = nMom*G4RandomDirection();
621
622 if( !f2p2h || A < 3 ) // 1p1h
623 {
624 // hM = tM - rM;
625
626 fLVt = G4LorentzVector( -nMomDir, sqrt( (rM + Ex)*(rM + Ex) + nMom*nMom ) ); // rM ); //
627 fLVh = G4LorentzVector( nMomDir, ei ); // hM); //
628 }
629 else // 2p2h
630 {
631 G4Nucleus recoil(A-2,Z-1);
632 rM = recoil.AtomicMass(A-2,Z-1)+sqrt(nMom*nMom+fM1*fM1);
633 hM = tM - rM;
634
635 fLVt = G4LorentzVector( nMomDir, sqrt( rM*rM+nMom*nMom ) );
636 fLVh = G4LorentzVector(-nMomDir, sqrt( hM*hM+nMom*nMom ) );
637 }
638 // G4cout<<hM<<", ";
639 // bst = fLVh.boostVector();
640
641 // lvp1.boost(-bst); // -> nucleon rest system, where Q2 transfer is ???
642
643 fNuEnergy = lvp1.e();
644 // G4double mN = fLVh.m(); // better mN = fM1 !? vmg
645 iTer = 0;
646
647 do // no FM!?, 5.4.20 vmg
648 {
652
653 // G4double mR = mN + fM1*(A-1.)*std::exp(-2.0*fQtransfer/mN); // recoil mass in+el
654
655 if( fXsample > 0. )
656 {
657 fW2 = fM1*fM1 - fQ2 + fQ2/fXsample; // sample excited hadron mass
658
659 // fW2 = mN*mN - fQ2 + fQ2/fXsample; // sample excited hadron mass
660 // fEmu = fNuEnergy - fQ2/2./mR/fXsample; // fM1->mN
661
662 fEmu = fNuEnergy - fQ2/2./fM1/fXsample; // fM1->mN
663 }
664 else
665 {
666 // fW2 = mN*mN;
667
668 fW2 = fM1*fM1;
669 fEmu = fNuEnergy;
670 }
671 // if(fEmu < 0.) G4cout<<"fEmu = "<<fEmu<<" hM = "<<hM<<G4endl;
672 // e3 = fNuEnergy + mR - fEmu;
673
674 e3 = fNuEnergy + fM1 - fEmu;
675
676 // if( e3 < sqrt(fW2) ) G4cout<<"energyX = "<<e3/GeV<<", fW = "<<sqrt(fW2)/GeV<<G4endl;
677
678 pMu2 = fEmu*fEmu - fMu*fMu;
679 pX2 = e3*e3 - fW2;
680
681 if(pMu2 < 0.) { fBreak = true; return; }
682
683 fCosTheta = fNuEnergy*fNuEnergy + pMu2 - pX2;
684 fCosTheta /= 2.*fNuEnergy*sqrt(pMu2);
685 iTer++;
686 }
687 while( ( abs(fCosTheta) > 1. || fEmu < fMu ) && iTer < iTerMax );
688
689 if( iTer >= iTerMax ) { fBreak = true; return; }
690
691 if( abs(fCosTheta) > 1.) // vmg: due to big Q2/x values. To be improved ...
692 {
693 G4cout<<"FM: fCosTheta = "<<fCosTheta<<", fEmu = "<<fEmu<<G4endl;
694 // fCosTheta = -1. + 2.*G4UniformRand();
695 if( fCosTheta < -1.) fCosTheta = -1.;
696 if( fCosTheta > 1.) fCosTheta = 1.;
697 }
698 // LVs
699 // G4LorentzVector lvt1 = G4LorentzVector( 0., 0., 0., mN ); // fM1 );
700
701 G4LorentzVector lvt1 = G4LorentzVector( 0., 0., 0., fM1 ); // fM1 );
702 G4LorentzVector lvsum = lvp1 + lvt1;
703
704 cost = fCosTheta;
705 sint = std::sqrt( (1.0 - cost)*(1.0 + cost) );
706 phi = G4UniformRand()*CLHEP::twopi;
707 eP = G4ThreeVector( sint*std::cos(phi), sint*std::sin(phi), cost );
708 muMom = sqrt(fEmu*fEmu-fMu*fMu);
709 eP *= muMom;
710 fLVl = G4LorentzVector( eP, fEmu );
711 fLVh = lvsum - fLVl;
712
713 // if( fLVh.e() < mN || fLVh.m2() < 0.) { fBreak = true; return; }
714
715 if( fLVh.e() < fM1 || fLVh.m2() < 0.) { fBreak = true; return; }
716
717 // back to lab system
718
719 // fLVl.boost(bst);
720 // fLVh.boost(bst);
721 }
722 //G4cout<<iTer<<", "<<fBreak<<"; ";
723}
724
725//
726//
727///////////////////////////
double A(double temperature)
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
#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="ANuMuNucleCcModel")
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 GetNuMuQeTotRat(G4int index, G4double energy)
G4int GetEnergyIndex(G4double energy)
G4double GgSampleNM(G4Nucleus &nucl)
G4double GetNuMuOnePionProb(G4int index, G4double energy)
static G4double fNuMuXdistrKR[50][50]
static G4double fNuMuQdistrKR[50][51][50]
G4ParticleDefinition * theMuonPlus
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()
Definition: DoubConv.h:17