Geant4 Cross Reference |
1 // 1 2 // ******************************************* 3 // * License and Disclaimer 4 // * 5 // * The Geant4 software is copyright of th 6 // * the Geant4 Collaboration. It is provided 7 // * conditions of the Geant4 Software License 8 // * LICENSE and available at http://cern.ch/ 9 // * include a list of copyright holders. 10 // * 11 // * Neither the authors of this software syst 12 // * institutes,nor the agencies providing fin 13 // * work make any representation or warran 14 // * regarding this software system or assum 15 // * use. Please see the license in the file 16 // * for the full disclaimer and the limitatio 17 // * 18 // * This code implementation is the result 19 // * technical work of the GEANT4 collaboratio 20 // * By using, copying, modifying or distri 21 // * any work based on the software) you ag 22 // * use in resulting scientific publicati 23 // * acceptance of all terms of the Geant4 Sof 24 // ******************************************* 25 // 26 // $Id: G4NuTauNucleusNcModel.cc 91806 2015-08 27 // 28 // Geant4 Header : G4NuTauNucleusNcModel 29 // 30 // Author : V.Grichine 12.2.19 31 // 32 33 #include "G4NuTauNucleusNcModel.hh" 34 #include "G4NeutrinoNucleusModel.hh" 35 36 // #include "G4NuMuResQX.hh" 37 38 #include "G4SystemOfUnits.hh" 39 #include "G4ParticleTable.hh" 40 #include "G4ParticleDefinition.hh" 41 #include "G4IonTable.hh" 42 #include "Randomize.hh" 43 #include "G4RandomDirection.hh" 44 45 // #include "G4Integrator.hh" 46 #include "G4DataVector.hh" 47 #include "G4PhysicsTable.hh" 48 #include "G4KineticTrack.hh" 49 #include "G4DecayKineticTracks.hh" 50 #include "G4KineticTrackVector.hh" 51 #include "G4Fragment.hh" 52 #include "G4ReactionProductVector.hh" 53 54 55 #include "G4NeutrinoTau.hh" 56 #include "G4AntiNeutrinoTau.hh" 57 #include "G4Nucleus.hh" 58 #include "G4LorentzVector.hh" 59 60 using namespace std; 61 using namespace CLHEP; 62 63 #ifdef G4MULTITHREADED 64 G4Mutex G4NuTauNucleusNcModel::numuNucleus 65 #endif 66 67 68 G4NuTauNucleusNcModel::G4NuTauNucleusNcModel(c 69 : G4NeutrinoNucleusModel(name) 70 { 71 SetMinEnergy( 0.0*GeV ); 72 SetMaxEnergy( 100.*TeV ); 73 SetMinEnergy(1.e-6*eV); 74 75 theNuTau = G4NeutrinoTau::NeutrinoTau(); 76 theANuTau = G4AntiNeutrinoTau::AntiNeutrino 77 78 fMnumu = 0.; 79 fData = fMaster = false; 80 InitialiseModel(); 81 82 } 83 84 85 G4NuTauNucleusNcModel::~G4NuTauNucleusNcModel( 86 {} 87 88 89 void G4NuTauNucleusNcModel::ModelDescription(s 90 { 91 92 outFile << "G4NuTauNucleusNcModel is a tau 93 << "model which uses the standard 94 << "transfer parameterization. Th 95 96 } 97 98 ////////////////////////////////////////////// 99 // 100 // Read data from G4PARTICLEXSDATA (locally PA 101 102 void G4NuTauNucleusNcModel::InitialiseModel() 103 { 104 G4String pName = "nu_mu"; 105 106 G4int nSize(0), i(0), j(0), k(0); 107 108 if(!fData) 109 { 110 #ifdef G4MULTITHREADED 111 G4MUTEXLOCK(&numuNucleusModel); 112 if(!fData) 113 { 114 #endif 115 fMaster = true; 116 #ifdef G4MULTITHREADED 117 } 118 G4MUTEXUNLOCK(&numuNucleusModel); 119 #endif 120 } 121 122 if(fMaster) 123 { 124 const char* path = G4FindDataDir("G4PARTIC 125 std::ostringstream ost1, ost2, ost3, ost4; 126 ost1 << path << "/" << "neutrino" << "/" < 127 128 std::ifstream filein1( ost1.str().c_str() 129 130 // filein.open("$PARTICLEXSDATA/"); 131 132 filein1>>nSize; 133 134 for( k = 0; k < fNbin; ++k ) 135 { 136 for( i = 0; i <= fNbin; ++i ) 137 { 138 filein1 >> fNuMuXarrayKR[k][i]; 139 // G4cout<< fNuMuXarrayKR[k][i] << " 140 } 141 } 142 // G4cout<<G4endl<<G4endl; 143 144 ost2 << path << "/" << "neutrino" << "/" < 145 std::ifstream filein2( ost2.str().c_str() 146 147 filein2>>nSize; 148 149 for( k = 0; k < fNbin; ++k ) 150 { 151 for( i = 0; i < fNbin; ++i ) 152 { 153 filein2 >> fNuMuXdistrKR[k][i]; 154 // G4cout<< fNuMuXdistrKR[k][i] << " 155 } 156 } 157 // G4cout<<G4endl<<G4endl; 158 159 ost3 << path << "/" << "neutrino" << "/" < 160 std::ifstream filein3( ost3.str().c_str() 161 162 filein3>>nSize; 163 164 for( k = 0; k < fNbin; ++k ) 165 { 166 for( i = 0; i <= fNbin; ++i ) 167 { 168 for( j = 0; j <= fNbin; ++j ) 169 { 170 filein3 >> fNuMuQarrayKR[k][i][j]; 171 // G4cout<< fNuMuQarrayKR[k][i][j] < 172 } 173 } 174 } 175 // G4cout<<G4endl<<G4endl; 176 177 ost4 << path << "/" << "neutrino" << "/" < 178 std::ifstream filein4( ost4.str().c_str() 179 180 filein4>>nSize; 181 182 for( k = 0; k < fNbin; ++k ) 183 { 184 for( i = 0; i <= fNbin; ++i ) 185 { 186 for( j = 0; j < fNbin; ++j ) 187 { 188 filein4 >> fNuMuQdistrKR[k][i][j]; 189 // G4cout<< fNuMuQdistrKR[k][i][j] < 190 } 191 } 192 } 193 fData = true; 194 } 195 } 196 197 ////////////////////////////////////////////// 198 199 G4bool G4NuTauNucleusNcModel::IsApplicable(con 200 G4Nucleus & ) 201 { 202 G4bool result = false; 203 G4String pName = aPart.GetDefinition()->GetP 204 G4double energy = aPart.GetTotalEnergy(); 205 206 if( pName == "nu_tau" // || pName == "anti_ 207 && 208 energy > fMinNuEnergy 209 { 210 result = true; 211 } 212 213 return result; 214 } 215 216 /////////////////////////////////////////// Cl 217 // 218 // 219 220 G4HadFinalState* G4NuTauNucleusNcModel::ApplyY 221 const G4HadProjectile& aTrack, G4Nucleus& 222 { 223 theParticleChange.Clear(); 224 fProton = f2p2h = fBreak = false; 225 const G4HadProjectile* aParticle = &aTrack; 226 G4double energy = aParticle->GetTotalEnergy( 227 228 G4String pName = aParticle->GetDefinition() 229 230 if( energy < fMinNuEnergy ) 231 { 232 theParticleChange.SetEnergyChange(energy); 233 theParticleChange.SetMomentumChange(aTrack 234 return &theParticleChange; 235 } 236 SampleLVkr( aTrack, targetNucleus); 237 238 if( fBreak == true || fEmu < fMnumu ) // ~5* 239 { 240 // G4cout<<"ni, "; 241 theParticleChange.SetEnergyChange(energy); 242 theParticleChange.SetMomentumChange(aTrack 243 return &theParticleChange; 244 } 245 246 // LVs of initial state 247 248 G4LorentzVector lvp1 = aParticle->Get4Moment 249 G4LorentzVector lvt1( 0., 0., 0., fM1 ); 250 G4double mPip = G4ParticleTable::GetParticle 251 252 // 1-pi by fQtransfer && nu-energy 253 G4LorentzVector lvpip1( 0., 0., 0., mPip ); 254 G4LorentzVector lvsum, lv2, lvX; 255 G4ThreeVector eP; 256 G4double cost(1.), sint(0.), phi(0.), muMom( 257 G4DynamicParticle* aLept = nullptr; // lepto 258 259 G4int Z = targetNucleus.GetZ_asInt(); 260 G4int A = targetNucleus.GetA_asInt(); 261 G4double mTarg = targetNucleus.AtomicMass(A 262 G4int pdgP(0), qB(0); 263 // G4double mSum = G4ParticleTable::GetParti 264 265 G4int iPi = GetOnePionIndex(energy); 266 G4double p1pi = GetNuMuOnePionProb( iPi, ene 267 268 if( p1pi > G4UniformRand() && fCosTheta > 0. 269 { 270 // lvsum = lvp1 + lvpip1; 271 lvsum = lvp1 + lvt1; 272 // cost = fCosThetaPi; 273 cost = fCosTheta; 274 sint = std::sqrt( (1.0 - cost)*(1.0 + cost 275 phi = G4UniformRand()*CLHEP::twopi; 276 eP = G4ThreeVector( sint*std::cos(phi), 277 278 // muMom = sqrt(fEmuPi*fEmuPi-fMnumu*fMnum 279 muMom = sqrt(fEmu*fEmu-fMnumu*fMnumu); 280 281 eP *= muMom; 282 283 // lv2 = G4LorentzVector( eP, fEmuPi ); 284 lv2 = G4LorentzVector( eP, fEmu ); 285 lv2 = fLVl; 286 287 lvX = lvsum - lv2; 288 lvX = fLVh; 289 massX2 = lvX.m2(); 290 G4double massX = lvX.m(); 291 G4double massR = fLVt.m(); 292 293 // if ( massX2 <= 0. ) // vmg: very rarely 294 if ( massX2 <= fM1*fM1 ) // 9-3-20 vmg: ve 295 if ( lvX.e() <= fM1 ) // 9-3-20 vmg: ver 296 { 297 theParticleChange.SetEnergyChange(energy 298 theParticleChange.SetMomentumChange(aTra 299 return &theParticleChange; 300 } 301 fW2 = massX2; 302 303 if( pName == "nu_tau" ) aLept = n 304 else if( pName == "anti_nu_tau") aLept = n 305 else 306 { 307 theParticleChange.SetEnergyChange(energy 308 theParticleChange.SetMomentumChange(aTra 309 return &theParticleChange; 310 } 311 312 pdgP = 111; 313 314 G4double eCut; // = fMpi + 0.5*(fMpi*fMpi 315 316 if( A > 1 ) 317 { 318 eCut = (fMpi + mTarg)*(fMpi + mTarg) - ( 319 eCut /= 2.*massR; 320 eCut += massX; 321 } 322 else eCut = fM1 + fMpi; 323 324 if ( lvX.e() > eCut ) // && sqrt( GetW2() 325 { 326 CoherentPion( lvX, pdgP, targetNucleus); 327 } 328 else 329 { 330 theParticleChange.SetEnergyChange(energy 331 theParticleChange.SetMomentumChange(aTra 332 return &theParticleChange; 333 } 334 theParticleChange.AddSecondary( aLept, fSe 335 336 return &theParticleChange; 337 } 338 else // lepton part in lab 339 { 340 lvsum = lvp1 + lvt1; 341 cost = fCosTheta; 342 sint = std::sqrt( (1.0 - cost)*(1.0 + cost 343 phi = G4UniformRand()*CLHEP::twopi; 344 eP = G4ThreeVector( sint*std::cos(phi), 345 346 muMom = sqrt(fEmu*fEmu-fMnumu*fMnumu); 347 348 eP *= muMom; 349 350 lv2 = G4LorentzVector( eP, fEmu ); 351 352 lvX = lvsum - lv2; 353 354 massX2 = lvX.m2(); 355 356 if ( massX2 <= 0. ) // vmg: very rarely ~ 357 { 358 theParticleChange.SetEnergyChange(energy 359 theParticleChange.SetMomentumChange(aTra 360 return &theParticleChange; 361 } 362 fW2 = massX2; 363 364 aLept = new G4DynamicParticle( theNuTau, l 365 366 theParticleChange.AddSecondary( aLept, fSe 367 } 368 369 // hadron part 370 371 fRecoil = nullptr; 372 fCascade = false; 373 fString = false; 374 375 if( A == 1 ) 376 { 377 qB = 1; 378 379 // if( G4UniformRand() > 0.1 ) // > 0.999 380 { 381 ClusterDecay( lvX, qB ); 382 } 383 return &theParticleChange; 384 } 385 G4Nucleus recoil; 386 G4double ratio = G4double(Z)/G4double(A); 387 388 if( ratio > G4UniformRand() ) // proton is e 389 { 390 fProton = true; 391 recoil = G4Nucleus(A-1,Z-1); 392 fRecoil = &recoil; 393 fMt = G4ParticleTable::GetParticleTable()- 394 + G4ParticleTable::GetParticleTable( 395 } 396 else // excited neutron 397 { 398 fProton = false; 399 recoil = G4Nucleus(A-1,Z); 400 fRecoil = &recoil; 401 fMt = G4ParticleTable::GetParticleTable()- 402 + G4ParticleTable::GetParticleTable( 403 } 404 // G4int index = GetEnergyIndex(energy 405 G4int nepdg = aParticle->GetDefinition()->Ge 406 407 G4double qeTotRat; // = GetNuMuQeTotRat(ind 408 qeTotRat = CalculateQEratioA( Z, A, energy, 409 410 G4ThreeVector dX = (lvX.vect()).unit(); 411 G4double eX = lvX.e(); // excited nucleon 412 G4double mX = sqrt(massX2); 413 414 if( qeTotRat > G4UniformRand() || mX <= fMt 415 { 416 fString = false; 417 418 G4double rM; 419 if( fProton ) 420 { 421 fPDGencoding = 2212; 422 fMr = proton_mass_c2; 423 recoil = G4Nucleus(A-1,Z-1); 424 fRecoil = &recoil; 425 rM = recoil.AtomicMass(A-1,Z-1); 426 } 427 else 428 { 429 fPDGencoding = 2112; 430 fMr = G4ParticleTable::GetParticleTabl 431 FindParticle(fPDGencoding)->GetPDGMass(); // 432 recoil = G4Nucleus(A-1,Z); 433 fRecoil = &recoil; 434 rM = recoil.AtomicMass(A-1,Z); 435 } 436 G4double eTh = fMr+0.5*(fMr*fMr-mX*mX)/rM; 437 438 if(eX <= eTh) // vmg, very rarely out of k 439 { 440 theParticleChange.SetEnergyChange(energy 441 theParticleChange.SetMomentumChange(aTra 442 return &theParticleChange; 443 } 444 FinalBarion( lvX, 0, fPDGencoding ); // p( 445 } 446 else // if ( eX < 9500000.*GeV ) // < 25.*Ge 447 { 448 if ( fProton && pName == "nu_tau" ) 449 else if( !fProton && pName == "nu_tau" ) 450 451 ClusterDecay( lvX, qB ); 452 } 453 return &theParticleChange; 454 } 455 456 457 ////////////////////////////////////////////// 458 ////////////////////////////////////////////// 459 ////////////////////////////////////////////// 460 461 ////////////////////////////////////////////// 462 // 463 // sample x, then Q2 464 465 void G4NuTauNucleusNcModel::SampleLVkr(const G 466 { 467 fBreak = false; 468 G4int A = targetNucleus.GetA_asInt(), iTer(0 469 G4int Z = targetNucleus.GetZ_asInt(); 470 G4double e3(0.), pMu2(0.), pX2(0.), nMom(0.) 471 G4double cost(1.), sint(0.), phi(0.), muMom( 472 G4ThreeVector eP, bst; 473 const G4HadProjectile* aParticle = &aTrack; 474 G4LorentzVector lvp1 = aParticle->Get4Moment 475 nMom = NucleonMomentum( targetNucleus ); 476 477 if( A == 1 || nMom == 0. ) // hydrogen, no F 478 { 479 fNuEnergy = aParticle->GetTotalEnergy(); 480 iTer = 0; 481 482 do 483 { 484 fXsample = SampleXkr(fNuEnergy); 485 fQtransfer = SampleQkr(fNuEnergy, fXsamp 486 fQ2 = fQtransfer*fQtransfer; 487 488 if( fXsample > 0. ) 489 { 490 fW2 = fM1*fM1 - fQ2 + fQ2/fXsample; // 491 fEmu = fNuEnergy - fQ2/2./fM1/fXsample 492 } 493 else 494 { 495 fW2 = fM1*fM1; 496 fEmu = fNuEnergy; 497 } 498 e3 = fNuEnergy + fM1 - fEmu; 499 500 // if( e3 < sqrt(fW2) ) G4cout<<"energy 501 502 pMu2 = fEmu*fEmu - fMnumu*fMnumu; 503 pX2 = e3*e3 - fW2; 504 505 fCosTheta = fNuEnergy*fNuEnergy + pMu2 506 fCosTheta /= 2.*fNuEnergy*sqrt(pMu2); 507 iTer++; 508 } 509 while( ( abs(fCosTheta) > 1. || fEmu < fMn 510 511 if( iTer >= iTerMax ) { fBreak = true; ret 512 513 if( abs(fCosTheta) > 1.) // vmg: due to bi 514 { 515 G4cout<<"H2: fCosTheta = "<<fCosTheta<<" 516 // fCosTheta = -1. + 2.*G4UniformRand(); 517 if(fCosTheta < -1.) fCosTheta = -1.; 518 if(fCosTheta > 1.) fCosTheta = 1.; 519 } 520 // LVs 521 522 G4LorentzVector lvt1 = G4LorentzVector( 0 523 G4LorentzVector lvsum = lvp1 + lvt1; 524 525 cost = fCosTheta; 526 sint = std::sqrt( (1.0 - cost)*(1.0 + cost 527 phi = G4UniformRand()*CLHEP::twopi; 528 eP = G4ThreeVector( sint*std::cos(phi), 529 muMom = sqrt(fEmu*fEmu-fMnumu*fMnumu); 530 eP *= muMom; 531 fLVl = G4LorentzVector( eP, fEmu ); 532 533 fLVh = lvsum - fLVl; 534 fLVt = G4LorentzVector( 0., 0., 0., 0. ); 535 } 536 else // Fermi motion, Q2 in nucleon rest fra 537 { 538 G4ThreeVector nMomDir = nMom*G4RandomDirec 539 540 if( !f2p2h ) // 1p1h 541 { 542 G4Nucleus recoil(A-1,Z); 543 rM = sqrt( recoil.AtomicMass(A-1,Z)*reco 544 hM = tM - rM; 545 546 fLVt = G4LorentzVector( nMomDir, sqrt( r 547 fLVh = G4LorentzVector(-nMomDir, sqrt( h 548 } 549 else // 2p2h 550 { 551 G4Nucleus recoil(A-2,Z-1); 552 rM = recoil.AtomicMass(A-2,Z-1)+sqrt(nMo 553 hM = tM - rM; 554 555 fLVt = G4LorentzVector( nMomDir, sqrt( r 556 fLVh = G4LorentzVector(-nMomDir, sqrt( h 557 } 558 // G4cout<<hM<<", "; 559 // bst = fLVh.boostVector(); // 9-3-20 560 561 // lvp1.boost(-bst); // 9-3-20 -> nucleon 562 563 fNuEnergy = lvp1.e(); 564 iTer = 0; 565 566 do 567 { 568 fXsample = SampleXkr(fNuEnergy); 569 fQtransfer = SampleQkr(fNuEnergy, fXsamp 570 fQ2 = fQtransfer*fQtransfer; 571 572 if( fXsample > 0. ) 573 { 574 fW2 = fM1*fM1 - fQ2 + fQ2/fXsample; // 575 fEmu = fNuEnergy - fQ2/2./fM1/fXsample 576 } 577 else 578 { 579 fW2 = fM1*fM1; 580 fEmu = fNuEnergy; 581 } 582 583 // if(fEmu < 0.) G4cout<<"fEmu = "<<fEmu 584 585 e3 = fNuEnergy + fM1 - fEmu; 586 587 // if( e3 < sqrt(fW2) ) G4cout<<"energy 588 589 pMu2 = fEmu*fEmu - fMnumu*fMnumu; 590 pX2 = e3*e3 - fW2; 591 592 fCosTheta = fNuEnergy*fNuEnergy + pMu2 593 fCosTheta /= 2.*fNuEnergy*sqrt(pMu2); 594 iTer++; 595 } 596 while( ( abs(fCosTheta) > 1. || fEmu < fMn 597 598 if( iTer >= iTerMax ) { fBreak = true; ret 599 600 if( abs(fCosTheta) > 1.) // vmg: due to bi 601 { 602 G4cout<<"FM: fCosTheta = "<<fCosTheta<<" 603 // fCosTheta = -1. + 2.*G4UniformRand(); 604 if(fCosTheta < -1.) fCosTheta = -1.; 605 if(fCosTheta > 1.) fCosTheta = 1.; 606 } 607 // LVs 608 G4LorentzVector lvt1 = G4LorentzVector( 0 609 G4LorentzVector lvsum = lvp1 + lvt1; 610 611 cost = fCosTheta; 612 sint = std::sqrt( (1.0 - cost)*(1.0 + cost 613 phi = G4UniformRand()*CLHEP::twopi; 614 eP = G4ThreeVector( sint*std::cos(phi), 615 muMom = sqrt(fEmu*fEmu-fMnumu*fMnumu); 616 eP *= muMom; 617 fLVl = G4LorentzVector( eP, fEmu ); 618 fLVh = lvsum - fLVl; 619 // back to lab system 620 // fLVl.boost(bst); // 9-3-20 621 // fLVh.boost(bst); // 9-3-20 622 } 623 //G4cout<<iTer<<", "<<fBreak<<"; "; 624 } 625 626 // 627 // 628 /////////////////////////// 629