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 // neutron_hp -- source file 27 // J.P. Wellisch, Nov-1996 28 // A prototype of the low energy neutron trans 29 // 30 // 070523 Try to limit sum of secondary photon 31 // in the of nDiscrete = 1 an nPartial 32 // T. Koi 33 // 070606 Add Partial case by T. Koi 34 // 070618 fix memory leaking by T. Koi 35 // 080801 fix memory leaking by T. Koi 36 // 080801 Correcting data disorder which happe 37 // and InitAnglurar methods was called 38 // 090514 Fix bug in IC electron emission case 39 // Contribution from Chao Zhang (Chao.Z 40 // But it looks like never cause real e 41 // 101111 Change warning message for "repFlag 42 // 43 // there is a lot of unused (and undebugged) c 44 // P. Arce, June-2014 Conversion neutron_hp to 45 // 46 #include "G4ParticleHPPhotonDist.hh" 47 48 #include "G4Electron.hh" 49 #include "G4ParticleHPLegendreStore.hh" 50 #include "G4PhysicalConstants.hh" 51 #include "G4Poisson.hh" 52 #include "G4SystemOfUnits.hh" 53 54 #include <numeric> 55 56 G4bool G4ParticleHPPhotonDist::InitMean(std::i 57 { 58 G4bool result = true; 59 if (aDataFile >> repFlag) { 60 aDataFile >> targetMass; 61 if (repFlag == 1) { 62 // multiplicities 63 aDataFile >> nDiscrete; 64 const std::size_t msize = nDiscrete > 0 65 disType = new G4int[msize]; 66 energy = new G4double[msize]; 67 // actualMult = new G4int[msize]; 68 theYield = new G4ParticleHPVector[msize] 69 for (std::size_t i = 0; i < msize; ++i) 70 aDataFile >> disType[i] >> energy[i]; 71 energy[i] *= eV; 72 theYield[i].Init(aDataFile, eV); 73 } 74 } 75 else if (repFlag == 2) { 76 aDataFile >> theInternalConversionFlag; 77 aDataFile >> theBaseEnergy; 78 theBaseEnergy *= eV; 79 aDataFile >> theInternalConversionFlag; 80 aDataFile >> nGammaEnergies; 81 const std::size_t esize = nGammaEnergies 82 theLevelEnergies = new G4double[esize]; 83 theTransitionProbabilities = new G4doubl 84 if (theInternalConversionFlag == 2) 85 thePhotonTransitionFraction = new G4do 86 for (std::size_t ii = 0; ii < esize; ++i 87 if (theInternalConversionFlag == 1) { 88 aDataFile >> theLevelEnergies[ii] >> 89 theLevelEnergies[ii] *= eV; 90 } 91 else if (theInternalConversionFlag == 92 aDataFile >> theLevelEnergies[ii] >> 93 >> thePhotonTransitionFraction[ii] 94 theLevelEnergies[ii] *= eV; 95 } 96 else { 97 throw G4HadronicException(__FILE__, 98 "G4Particl 99 } 100 } 101 } 102 else { 103 G4cout << "Data representation in G4Part 104 throw G4HadronicException( 105 __FILE__, __LINE__, "G4ParticleHPPhoto 106 } 107 } 108 else { 109 result = false; 110 } 111 return result; 112 } 113 114 void G4ParticleHPPhotonDist::InitAngular(std:: 115 { 116 G4int i, ii; 117 // angular distributions 118 aDataFile >> isoFlag; 119 if (isoFlag != 1) { 120 if (repFlag == 2) 121 G4cout << "G4ParticleHPPhotonDist: repFl 122 "G4ND3.x, then please report t 123 << G4endl; 124 aDataFile >> tabulationType >> nDiscrete2 125 if (theGammas != nullptr && nDiscrete2 != 126 G4cout << "080731c G4ParticleHPPhotonDis 127 "wrong in your NDL files. Plea 128 "messages after the update, th 129 << G4endl; 130 131 // The order of cross section (InitPartial 132 // (InitAngular here) data are different, 133 // consistent data order. 134 std::vector<G4double> vct_gammas_par; 135 std::vector<G4double> vct_shells_par; 136 std::vector<G4int> vct_primary_par; 137 std::vector<G4int> vct_distype_par; 138 std::vector<G4ParticleHPVector*> vct_pXS_p 139 if (theGammas != nullptr && theShells != n 140 // copy the cross section data 141 for (i = 0; i < nDiscrete; ++i) { 142 vct_gammas_par.push_back(theGammas[i]) 143 vct_shells_par.push_back(theShells[i]) 144 vct_primary_par.push_back(isPrimary[i] 145 vct_distype_par.push_back(disType[i]); 146 auto hpv = new G4ParticleHPVector; 147 *hpv = thePartialXsec[i]; 148 vct_pXS_par.push_back(hpv); 149 } 150 } 151 const std::size_t psize = nDiscrete2 > 0 ? 152 if (theGammas == nullptr) theGammas = new 153 if (theShells == nullptr) theShells = new 154 155 for (i = 0; i < nIso; ++i) // isotropic p 156 { 157 aDataFile >> theGammas[i] >> theShells[i 158 theGammas[i] *= eV; 159 theShells[i] *= eV; 160 } 161 const std::size_t tsize = nDiscrete2 - nIs 162 nNeu = new G4int[tsize]; 163 if (tabulationType == 1) theLegendre = new 164 if (tabulationType == 2) theAngular = new 165 for (i = nIso; i < nDiscrete2; ++i) { 166 if (tabulationType == 1) { 167 aDataFile >> theGammas[i] >> theShells 168 theGammas[i] *= eV; 169 theShells[i] *= eV; 170 const std::size_t lsize = nNeu[i - nIs 171 theLegendre[i - nIso] = new G4Particle 172 theLegendreManager.Init(aDataFile); 173 for (ii = 0; ii < nNeu[i - nIso]; ++ii 174 theLegendre[i - nIso][ii].Init(aData 175 } 176 } 177 else if (tabulationType == 2) { 178 aDataFile >> theGammas[i] >> theShells 179 theGammas[i] *= eV; 180 theShells[i] *= eV; 181 const std::size_t asize = nNeu[i - nIs 182 theAngular[i - nIso] = new G4ParticleH 183 for (ii = 0; ii < nNeu[i - nIso]; ++ii 184 theAngular[i - nIso][ii].Init(aDataF 185 } 186 } 187 else { 188 G4cout << "tabulation type: tabulation 189 throw G4HadronicException( 190 __FILE__, __LINE__, "cannot deal wit 191 } 192 } 193 194 if (!vct_gammas_par.empty()) { 195 // Reordering cross section data to corr 196 for (i = 0; i < nDiscrete; ++i) { 197 for (G4int j = 0; j < nDiscrete; ++j) 198 // Checking gamma and shell to ident 199 if (theGammas[i] == vct_gammas_par[j 200 isPrimary[i] = vct_primary_par[j]; 201 disType[i] = vct_distype_par[j]; 202 thePartialXsec[i] = (*(vct_pXS_par 203 } 204 } 205 } 206 // Garbage collection 207 for (auto it = vct_pXS_par.cbegin(); it 208 delete *it; 209 } 210 } 211 } 212 } 213 214 void G4ParticleHPPhotonDist::InitEnergies(std: 215 { 216 G4int i, energyDistributionsNeeded = 0; 217 for (i = 0; i < nDiscrete; ++i) { 218 if (disType[i] == 1) energyDistributionsNe 219 } 220 if (energyDistributionsNeeded == 0) return; 221 aDataFile >> nPartials; 222 const std::size_t dsize = nPartials > 0 ? nP 223 distribution = new G4int[dsize]; 224 probs = new G4ParticleHPVector[dsize]; 225 partials = new G4ParticleHPPartial*[dsize]; 226 G4int nen; 227 G4int dummy; 228 for (i = 0; i < nPartials; ++i) { 229 aDataFile >> dummy; 230 probs[i].Init(aDataFile, eV); 231 aDataFile >> nen; 232 partials[i] = new G4ParticleHPPartial(nen) 233 partials[i]->InitInterpolation(aDataFile); 234 partials[i]->Init(aDataFile); 235 } 236 } 237 238 void G4ParticleHPPhotonDist::InitPartials(std: 239 { 240 if (theXsec != nullptr) theReactionXsec = th 241 242 aDataFile >> nDiscrete >> targetMass; 243 if (nDiscrete != 1) { 244 theTotalXsec.Init(aDataFile, eV); 245 } 246 const std::size_t dsize = nDiscrete > 0 ? nD 247 theGammas = new G4double[dsize]; 248 theShells = new G4double[dsize]; 249 isPrimary = new G4int[dsize]; 250 disType = new G4int[dsize]; 251 thePartialXsec = new G4ParticleHPVector[dsiz 252 for (std::size_t i = 0; i < dsize; ++i) { 253 aDataFile >> theGammas[i] >> theShells[i] 254 theGammas[i] *= eV; 255 theShells[i] *= eV; 256 thePartialXsec[i].Init(aDataFile, eV); 257 } 258 } 259 260 G4ReactionProductVector* G4ParticleHPPhotonDis 261 { 262 // the partial cross-section case is not all 263 if (actualMult.Get() == nullptr) { 264 actualMult.Get() = new std::vector<G4int>( 265 } 266 G4int i, ii, iii; 267 G4int nSecondaries = 0; 268 auto thePhotons = new G4ReactionProductVecto 269 270 if (repFlag == 1) { 271 G4double current = 0; 272 for (i = 0; i < nDiscrete; ++i) { 273 current = theYield[i].GetY(anEnergy); 274 actualMult.Get()->at(i) = (G4int)G4Poiss 275 if (nDiscrete == 1 && current < 1.0001) 276 actualMult.Get()->at(i) = static_cast< 277 if (current < 1) { 278 actualMult.Get()->at(i) = 0; 279 if (G4UniformRand() < current) actua 280 } 281 } 282 nSecondaries += actualMult.Get()->at(i); 283 } 284 for (i = 0; i < nSecondaries; ++i) { 285 auto theOne = new G4ReactionProduct; 286 theOne->SetDefinition(G4Gamma::Gamma()); 287 thePhotons->push_back(theOne); 288 } 289 290 G4int count = 0; 291 292 if (nDiscrete == 1 && nPartials == 1) { 293 if (actualMult.Get()->at(0) > 0) { 294 if (disType[0] == 1) { 295 // continuum 296 G4ParticleHPVector* temp; 297 temp = partials[0]->GetY(anEnergy); 298 G4double maximumE = temp->GetX(temp- 299 300 std::vector<G4double> photons_e_best 301 G4double best = DBL_MAX; 302 G4int maxTry = 1000; 303 for (G4int j = 0; j < maxTry; ++j) { 304 std::vector<G4double> photons_e(ac 305 for (auto it = photons_e.begin(); 306 *it = temp->Sample(); 307 } 308 309 if (std::accumulate(photons_e.cbeg 310 if (std::accumulate(photons_e.cb 311 photons_e_best = std::move(pho 312 continue; 313 } 314 G4int iphot = 0; 315 for (auto it = photons_e.cbegin(); 316 thePhotons->operator[](iphot)->S 317 *it); // Replace index count, 318 // with iphot, which is 319 // bug report 2167 320 ++iphot; 321 } 322 323 break; 324 } 325 delete temp; 326 } 327 else { 328 // discrete 329 thePhotons->operator[](count)->SetKi 330 } 331 ++count; 332 if (count > nSecondaries) 333 throw G4HadronicException(__FILE__, 334 "G4Particl 335 } 336 } 337 else { // nDiscrete != 1 or nPartials != 338 for (i = 0; i < nDiscrete; ++i) { 339 for (ii = 0; ii < actualMult.Get()->at 340 if (disType[i] == 1) { 341 // continuum 342 G4double sum = 0, run = 0; 343 for (iii = 0; iii < nPartials; ++i 344 sum += probs[iii].GetY(anEnergy) 345 G4double random = G4UniformRand(); 346 G4int theP = 0; 347 for (iii = 0; iii < nPartials; ++i 348 run += probs[iii].GetY(anEnergy) 349 theP = iii; 350 if (random < run / sum) break; 351 } 352 353 if (theP == nPartials) theP = nPar 354 sum = 0; 355 G4ParticleHPVector* temp; 356 temp = partials[theP]->GetY(anEner 357 G4double eGamm = temp->Sample(); 358 thePhotons->operator[](count)->Set 359 delete temp; 360 } 361 else { 362 // discrete 363 thePhotons->operator[](count)->Set 364 } 365 ++count; 366 if (count > nSecondaries) 367 throw G4HadronicException(__FILE__ 368 "G4Parti 369 } 370 } 371 } 372 373 // now do the angular distributions... 374 if (isoFlag == 1) { 375 for (i = 0; i < nSecondaries; ++i) { 376 G4double costheta = 2. * G4UniformRand 377 G4double theta = std::acos(costheta); 378 G4double phi = twopi * G4UniformRand() 379 G4double sinth = std::sin(theta); 380 G4double en = thePhotons->operator[](i 381 G4ThreeVector temp(en * sinth * std::c 382 en * std::cos(theta 383 thePhotons->operator[](i)->SetMomentum 384 } 385 } 386 else { 387 for (i = 0; i < nSecondaries; ++i) { 388 G4double currentEnergy = thePhotons->o 389 for (ii = 0; ii < nDiscrete2; ++ii) { 390 if (std::abs(currentEnergy - theGamm 391 } 392 if (ii == nDiscrete2) 393 --ii; // fix for what seems an (fil 394 // data. @@ 395 if (ii < nIso) { 396 // isotropic distribution 397 // 398 // Fix Bugzilla report #1745 399 // G4double theta = pi*G4UniformRand 400 G4double costheta = 2. * G4UniformRa 401 G4double theta = std::acos(costheta) 402 G4double phi = twopi * G4UniformRand 403 G4double sinth = std::sin(theta); 404 G4double en = thePhotons->operator[] 405 // DHW G4ThreeVector tempVector(en*s 406 // en*std::cos(theta) ); 407 G4ThreeVector tempVector(en * sinth 408 en * costhe 409 thePhotons->operator[](i)->SetMoment 410 } 411 else if (tabulationType == 1) { 412 // legendre polynomials 413 G4int it(0); 414 for (iii = 0; iii < nNeu[ii - nIso]; 415 { 416 it = iii; 417 if (theLegendre[ii - nIso][iii].Ge 418 } 419 G4ParticleHPLegendreStore aStore(2); 420 aStore.SetCoeff(1, &(theLegendre[ii 421 if (it > 0) { 422 aStore.SetCoeff(0, &(theLegendre[i 423 } 424 else { 425 aStore.SetCoeff(0, &(theLegendre[i 426 } 427 G4double cosTh = aStore.SampleMax(an 428 G4double theta = std::acos(cosTh); 429 G4double phi = twopi * G4UniformRand 430 G4double sinth = std::sin(theta); 431 G4double en = thePhotons->operator[] 432 G4ThreeVector tempVector(en * sinth 433 en * std::c 434 thePhotons->operator[](i)->SetMoment 435 } 436 else { 437 // tabulation of probabilities. 438 G4int it(0); 439 for (iii = 0; iii < nNeu[ii - nIso]; 440 { 441 it = iii; 442 if (theAngular[ii - nIso][iii].Get 443 } 444 G4double costh = theAngular[ii - nIs 445 G4double theta = std::acos(costh); 446 G4double phi = twopi * G4UniformRand 447 G4double sinth = std::sin(theta); 448 G4double en = thePhotons->operator[] 449 G4ThreeVector tmpVector(en * sinth * 450 en * costh); 451 thePhotons->operator[](i)->SetMoment 452 } 453 } 454 } 455 } 456 else if (repFlag == 2) { 457 auto running = new G4double[nGammaEnergies 458 running[0] = theTransitionProbabilities[0] 459 for (i = 1; i < nGammaEnergies; ++i) { 460 running[i] = running[i - 1] + theTransit 461 } 462 G4double random = G4UniformRand(); 463 G4int it = 0; 464 for (i = 0; i < nGammaEnergies; ++i) { 465 it = i; 466 if (random < running[i] / running[nGamma 467 } 468 delete[] running; 469 G4double totalEnergy = theBaseEnergy - the 470 auto theOne = new G4ReactionProduct; 471 theOne->SetDefinition(G4Gamma::Gamma()); 472 random = G4UniformRand(); 473 if (theInternalConversionFlag == 2 && rand 474 theOne->SetDefinition(G4Electron::Electr 475 // Bug reported Chao Zhang (Chao.Zhang@u 476 // 2009 But never enter at least with G4 477 totalEnergy += 478 G4Electron::Electron()->GetPDGMass(); 479 } 480 theOne->SetTotalEnergy(totalEnergy); 481 if (isoFlag == 1) { 482 G4double costheta = 2. * G4UniformRand() 483 G4double theta = std::acos(costheta); 484 G4double phi = twopi * G4UniformRand(); 485 G4double sinth = std::sin(theta); 486 // Bug reported Chao Zhang (Chao.Zhang@u 487 // 2009 G4double en = theOne->GetTotalEn 488 G4double en = theOne->GetTotalMomentum() 489 // But never cause real effect at least 490 G4ThreeVector temp(en * sinth * std::cos 491 en * std::cos(theta)) 492 theOne->SetMomentum(temp); 493 } 494 else { 495 G4double currentEnergy = theOne->GetTota 496 for (ii = 0; ii < nDiscrete2; ++ii) { 497 if (std::abs(currentEnergy - theGammas 498 } 499 if (ii == nDiscrete2) 500 --ii; // fix for what seems an (file1 501 // data. @@ 502 if (ii < nIso) { 503 // Bug reported Chao Zhang (Chao.Zhang 504 // 2009 505 // isotropic distribution 506 // G4double theta = pi*G4UniformRand() 507 G4double theta = std::acos(2. * G4Unif 508 // But this is alos never cause real e 509 // isoFlag != 1 510 G4double phi = twopi * G4UniformRand() 511 G4double sinth = std::sin(theta); 512 // Bug reported Chao Zhang (Chao.Zhang 513 // 2009 G4double en = theOne->GetTotal 514 G4double en = theOne->GetTotalMomentum 515 // But never cause real effect at leas 516 G4ThreeVector tempVector(en * sinth * 517 en * std::cos 518 theOne->SetMomentum(tempVector); 519 } 520 else if (tabulationType == 1) { 521 // legendre polynomials 522 G4int itt(0); 523 for (iii = 0; iii < nNeu[ii - nIso]; + 524 { 525 itt = iii; 526 if (theLegendre[ii - nIso][iii].GetE 527 } 528 G4ParticleHPLegendreStore aStore(2); 529 aStore.SetCoeff(1, &(theLegendre[ii - 530 // aStore.SetCoeff(0, &(theLegendre[ii 531 // TKDB 110512 532 if (itt > 0) { 533 aStore.SetCoeff(0, &(theLegendre[ii 534 } 535 else { 536 aStore.SetCoeff(0, &(theLegendre[ii 537 } 538 G4double cosTh = aStore.SampleMax(anEn 539 G4double theta = std::acos(cosTh); 540 G4double phi = twopi * G4UniformRand() 541 G4double sinth = std::sin(theta); 542 // Bug reported Chao Zhang (Chao.Zhang 543 // 2009 G4double en = theOne->GetTotal 544 G4double en = theOne->GetTotalMomentum 545 // But never cause real effect at leas 546 G4ThreeVector tempVector(en * sinth * 547 en * std::cos 548 theOne->SetMomentum(tempVector); 549 } 550 else { 551 // tabulation of probabilities. 552 G4int itt(0); 553 for (iii = 0; iii < nNeu[ii - nIso]; + 554 { 555 itt = iii; 556 if (theAngular[ii - nIso][iii].GetEn 557 } 558 G4double costh = theAngular[ii - nIso] 559 G4double theta = std::acos(costh); 560 G4double phi = twopi * G4UniformRand() 561 G4double sinth = std::sin(theta); 562 // Bug reported Chao Zhang (Chao.Zhang 563 // 2009 G4double en = theOne->GetTotal 564 G4double en = theOne->GetTotalMomentum 565 // But never cause real effect at leas 566 G4ThreeVector tmpVector(en * sinth * s 567 theOne->SetMomentum(tmpVector); 568 } 569 } 570 thePhotons->push_back(theOne); 571 } 572 else if (repFlag == 0) { 573 if (thePartialXsec == nullptr) { 574 return thePhotons; 575 } 576 577 // Partial Case 578 579 auto theOne = new G4ReactionProduct; 580 theOne->SetDefinition(G4Gamma::Gamma()); 581 thePhotons->push_back(theOne); 582 583 // Energy 584 585 G4double sum = 0.0; 586 std::vector<G4double> dif(nDiscrete, 0.0); 587 for (G4int j = 0; j < nDiscrete; ++j) { 588 G4double x = thePartialXsec[j].GetXsec(a 589 if (x > 0) { 590 sum += x; 591 } 592 dif[j] = sum; 593 } 594 595 G4double rand = G4UniformRand(); 596 597 G4int iphoton = 0; 598 for (G4int j = 0; j < nDiscrete; ++j) { 599 G4double y = rand * sum; 600 if (dif[j] > y) { 601 iphoton = j; 602 break; 603 } 604 } 605 606 // Statistically suppress the photon accor 607 // Fix proposed by Artem Zontikov, Bug rep 608 if (theReactionXsec != nullptr) { 609 if (thePartialXsec[iphoton].GetXsec(anEn 610 < G4UniformRand()) 611 { 612 delete thePhotons; 613 thePhotons = nullptr; 614 return thePhotons; 615 } 616 } 617 618 // Angle 619 G4double cosTheta = 0.0; // mu 620 621 if (isoFlag == 1) { 622 // Isotropic Case 623 624 cosTheta = 2. * G4UniformRand() - 1; 625 } 626 else { 627 if (iphoton < nIso) { 628 // still Isotropic 629 630 cosTheta = 2. * G4UniformRand() - 1; 631 } 632 else { 633 if (tabulationType == 1) { 634 // Legendre polynomials 635 636 G4int iangle = 0; 637 for (G4int j = 0; j < nNeu[iphoton - 638 iangle = j; 639 if (theLegendre[iphoton - nIso][j] 640 } 641 642 G4ParticleHPLegendreStore aStore(2); 643 aStore.SetCoeff(1, &(theLegendre[iph 644 aStore.SetCoeff(0, &(theLegendre[iph 645 646 cosTheta = aStore.SampleMax(anEnergy 647 } 648 else if (tabulationType == 2) { 649 // tabulation of probabilities. 650 651 G4int iangle = 0; 652 for (G4int j = 0; j < nNeu[iphoton - 653 iangle = j; 654 if (theAngular[iphoton - nIso][j]. 655 } 656 cosTheta = theAngular[iphoton - nIso 657 // no interpolation yet @@ 658 } 659 } 660 } 661 662 // Set 663 G4double phi = twopi * G4UniformRand(); 664 G4double theta = std::acos(cosTheta); 665 G4double sinTheta = std::sin(theta); 666 667 G4double photonE = theGammas[iphoton]; 668 G4ThreeVector direction(sinTheta * std::co 669 G4ThreeVector photonP = photonE * directio 670 thePhotons->operator[](0)->SetMomentum(pho 671 } 672 else { 673 delete thePhotons; 674 thePhotons = nullptr; // no gamma data av 675 } 676 return thePhotons; 677 } 678