Geant4 Cross Reference |
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 // 27 // Geant4 class G4EmTableUtil 28 // 29 // Author V.Ivanchenko 14.03.2022 30 // 31 32 #include "G4EmTableUtil.hh" 33 #include "G4RegionStore.hh" 34 #include "G4ProductionCutsTable.hh" 35 #include "G4EmParameters.hh" 36 #include "G4EmUtility.hh" 37 #include "G4LossTableManager.hh" 38 #include "G4EmTableType.hh" 39 #include "G4PhysicsModelCatalog.hh" 40 #include "G4LowEnergyEmProcessSubType.hh" 41 #include "G4PhysicsTableHelper.hh" 42 #include "G4PhysicsLogVector.hh" 43 #include "G4ProcessManager.hh" 44 #include "G4UIcommand.hh" 45 #include "G4GenericIon.hh" 46 #include <iostream> 47 48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 49 50 const G4DataVector* 51 G4EmTableUtil::PrepareEmProcess(G4VEmProcess* proc, 52 const G4ParticleDefinition* part, 53 const G4ParticleDefinition* secPart, 54 G4EmModelManager* modelManager, 55 const G4double& maxKinEnergy, 56 G4int& secID, G4int& tripletID, 57 G4int& mainSec, const G4int& verb, 58 const G4bool& master) 59 { 60 G4EmParameters* param = G4EmParameters::Instance(); 61 62 // initialisation of models 63 G4double plimit = param->MscThetaLimit(); 64 G4int nModels = modelManager->NumberOfModels(); 65 for(G4int i=0; i<nModels; ++i) { 66 G4VEmModel* mod = modelManager->GetModel(i); 67 if(nullptr == mod) { continue; } 68 mod->SetPolarAngleLimit(plimit); 69 if(mod->HighEnergyLimit() > maxKinEnergy) { 70 mod->SetHighEnergyLimit(maxKinEnergy); 71 } 72 proc->SetEmModel(mod); 73 } 74 75 // defined ID of secondary particles and verbosity 76 G4int stype = proc->GetProcessSubType(); 77 if(stype == fAnnihilation) { 78 secID = _Annihilation; 79 tripletID = _TripletGamma; 80 } else if(stype == fGammaConversion) { 81 secID = _PairProduction; 82 mainSec = 2; 83 } else if(stype == fPhotoElectricEffect) { 84 secID = _PhotoElectron; 85 } else if(stype == fComptonScattering) { 86 secID = _ComptonElectron; 87 } else if(stype >= fLowEnergyElastic) { 88 secID = fDNAUnknownModel; 89 } 90 if(master) { 91 proc->SetVerboseLevel(param->Verbose()); 92 } else { 93 proc->SetVerboseLevel(param->WorkerVerbose()); 94 } 95 96 // model initialisation 97 const G4DataVector* cuts = modelManager->Initialise(part, secPart, verb); 98 99 if(1 < verb) { 100 G4cout << "### G4EmTableUtil::PreparePhysicsTable() done for " 101 << proc->GetProcessName() 102 << " and particle " << part->GetParticleName() 103 << G4endl; 104 } 105 return cuts; 106 } 107 108 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 109 110 void G4EmTableUtil::BuildEmProcess(G4VEmProcess* proc, 111 const G4VEmProcess* masterProc, 112 const G4ParticleDefinition* firstPart, 113 const G4ParticleDefinition* part, 114 const G4int nModels, const G4int verb, 115 const G4bool master, const G4bool isLocked, 116 const G4bool toBuild, G4bool& baseMat) 117 { 118 G4String num = part->GetParticleName(); 119 if(1 < verb) { 120 G4cout << "### G4EmTableUtil::BuildPhysicsTable() for " 121 << proc->GetProcessName() << " and particle " << num 122 << " buildLambdaTable=" << toBuild << " master= " << master 123 << G4endl; 124 } 125 126 if(firstPart == part) { 127 128 // worker initialisation 129 if(!master) { 130 proc->SetLambdaTable(masterProc->LambdaTable()); 131 proc->SetLambdaTablePrim(masterProc->LambdaTablePrim()); 132 proc->SetCrossSectionType(masterProc->CrossSectionType()); 133 proc->SetEnergyOfCrossSectionMax(masterProc->EnergyOfCrossSectionMax()); 134 135 // local initialisation of models 136 baseMat = masterProc->UseBaseMaterial(); 137 G4bool printing = true; 138 for(G4int i=0; i<nModels; ++i) { 139 G4VEmModel* mod = proc->GetModelByIndex(i, printing); 140 G4VEmModel* mod0= masterProc->GetModelByIndex(i, printing); 141 mod->SetUseBaseMaterials(baseMat); 142 mod->InitialiseLocal(part, mod0); 143 } 144 // master thread 145 } else { 146 if(toBuild) { proc->BuildLambdaTable(); } 147 auto fXSType = proc->CrossSectionType(); 148 auto v = proc->EnergyOfCrossSectionMax(); 149 delete v; 150 v = nullptr; 151 if(fXSType == fEmOnePeak) { 152 auto table = proc->LambdaTable(); 153 if(nullptr == table) { 154 v = G4EmUtility::FindCrossSectionMax(proc, part); 155 } else { 156 v = G4EmUtility::FindCrossSectionMax(table); 157 } 158 if(nullptr == v) { proc->SetCrossSectionType(fEmIncreasing); } 159 } 160 proc->SetEnergyOfCrossSectionMax(v); 161 } 162 } 163 // protection against double printout 164 if(isLocked) { return; } 165 166 // explicitly defined printout by particle name 167 if(1 < verb || (0 < verb && (num == "gamma" || num == "e-" || 168 num == "e+" || num == "mu+" || 169 num == "mu-" || num == "proton"|| 170 num == "pi+" || num == "pi-" || 171 num == "kaon+" || num == "kaon-" || 172 num == "alpha" || num == "anti_proton" || 173 num == "GenericIon" || 174 num == "alpha+" || num == "helium" || 175 num == "hydrogen"))) { 176 proc->StreamInfo(G4cout, *part); 177 } 178 179 if(1 < verb) { 180 G4cout << "### G4EmTableUtil::BuildPhysicsTable() done for " 181 << proc->GetProcessName() << " and particle " << num 182 << " baseMat=" << baseMat << G4endl; 183 } 184 } 185 186 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 187 188 void G4EmTableUtil::BuildLambdaTable(G4VEmProcess* proc, 189 const G4ParticleDefinition* part, 190 G4EmModelManager* modelManager, 191 G4LossTableBuilder* bld, 192 G4PhysicsTable* theLambdaTable, 193 G4PhysicsTable* theLambdaTablePrim, 194 const G4double minKinEnergy, 195 const G4double minKinEnergyPrim, 196 const G4double maxKinEnergy, 197 const G4double scale, 198 const G4int verboseLevel, 199 const G4bool startFromNull, 200 const G4bool splineFlag) 201 { 202 if(1 < verboseLevel) { 203 G4cout << "G4EmTableUtil::BuildLambdaTable() for process " 204 << proc->GetProcessName() << " and particle " 205 << part->GetParticleName() << G4endl; 206 } 207 208 // Access to materials 209 const G4ProductionCutsTable* theCoupleTable= 210 G4ProductionCutsTable::GetProductionCutsTable(); 211 std::size_t numOfCouples = theCoupleTable->GetTableSize(); 212 213 G4PhysicsLogVector* aVector = nullptr; 214 G4PhysicsLogVector* aVectorPrim = nullptr; 215 G4PhysicsLogVector* bVectorPrim = nullptr; 216 217 G4double emax1 = std::min(maxKinEnergy, minKinEnergyPrim); 218 219 for(std::size_t i=0; i<numOfCouples; ++i) { 220 221 if (bld->GetFlag(i)) { 222 // create physics vector and fill it 223 const G4MaterialCutsCouple* couple = 224 theCoupleTable->GetMaterialCutsCouple((G4int)i); 225 226 // build main table 227 if(nullptr != theLambdaTable) { 228 delete (*theLambdaTable)[i]; 229 230 // if start from zero then change the scale 231 G4double emin = minKinEnergy; 232 G4bool startNull = false; 233 if(startFromNull) { 234 G4double e = proc->MinPrimaryEnergy(part, couple->GetMaterial()); 235 if(e >= emin) { 236 emin = e; 237 startNull = true; 238 } 239 } 240 G4double emax = emax1; 241 if(emax <= emin) { emax = 2*emin; } 242 G4int bin = G4lrint(scale*G4Log(emax/emin)); 243 bin = std::max(bin, 5); 244 aVector = new G4PhysicsLogVector(emin, emax, bin, splineFlag); 245 modelManager->FillLambdaVector(aVector, couple, startNull); 246 if(splineFlag) { aVector->FillSecondDerivatives(); } 247 G4PhysicsTableHelper::SetPhysicsVector(theLambdaTable, i, aVector); 248 } 249 // build high energy table 250 if(nullptr != theLambdaTablePrim && minKinEnergyPrim < maxKinEnergy) { 251 delete (*theLambdaTablePrim)[i]; 252 253 // start not from zero and always use spline 254 if(nullptr == bVectorPrim) { 255 G4int bin = G4lrint(scale*G4Log(maxKinEnergy/minKinEnergyPrim)); 256 bin = std::max(bin, 5); 257 aVectorPrim = 258 new G4PhysicsLogVector(minKinEnergyPrim, maxKinEnergy, bin, true); 259 bVectorPrim = aVectorPrim; 260 } else { 261 aVectorPrim = new G4PhysicsLogVector(*bVectorPrim); 262 } 263 modelManager->FillLambdaVector(aVectorPrim, couple, false, 264 fIsCrossSectionPrim); 265 aVectorPrim->FillSecondDerivatives(); 266 G4PhysicsTableHelper::SetPhysicsVector(theLambdaTablePrim, i, 267 aVectorPrim); 268 } 269 } 270 } 271 272 if(1 < verboseLevel) { 273 G4cout << "Lambda table is built for " << part->GetParticleName() << G4endl; 274 } 275 } 276 277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 278 279 void G4EmTableUtil::BuildLambdaTable(G4VEnergyLossProcess* proc, 280 const G4ParticleDefinition* part, 281 G4EmModelManager* modelManager, 282 G4LossTableBuilder* bld, 283 G4PhysicsTable* theLambdaTable, 284 const G4DataVector* theCuts, 285 const G4double minKinEnergy, 286 const G4double maxKinEnergy, 287 const G4double scale, 288 const G4int verboseLevel, 289 const G4bool splineFlag) 290 { 291 if(1 < verboseLevel) { 292 G4cout << "G4EmTableUtil::BuildLambdaTable() for process " 293 << proc->GetProcessName() << " and particle " 294 << part->GetParticleName() << G4endl; 295 } 296 297 const G4ProductionCutsTable* theCoupleTable= 298 G4ProductionCutsTable::GetProductionCutsTable(); 299 std::size_t numOfCouples = theCoupleTable->GetTableSize(); 300 301 G4PhysicsLogVector* aVector = nullptr; 302 for(std::size_t i=0; i<numOfCouples; ++i) { 303 if (bld->GetFlag(i)) { 304 // create physics vector and fill it 305 const G4MaterialCutsCouple* couple = 306 theCoupleTable->GetMaterialCutsCouple((G4int)i); 307 308 delete (*theLambdaTable)[i]; 309 G4bool startNull = true; 310 G4double emin = 311 proc->MinPrimaryEnergy(part, couple->GetMaterial(), (*theCuts)[i]); 312 if(minKinEnergy > emin) { 313 emin = minKinEnergy; 314 startNull = false; 315 } 316 317 G4double emax = maxKinEnergy; 318 if(emax <= emin) { emax = 2*emin; } 319 G4int bin = G4lrint(scale*G4Log(emax/emin)); 320 bin = std::max(bin, 5); 321 aVector = new G4PhysicsLogVector(emin, emax, bin, splineFlag); 322 modelManager->FillLambdaVector(aVector, couple, startNull, fRestricted); 323 if(splineFlag) { aVector->FillSecondDerivatives(); } 324 G4PhysicsTableHelper::SetPhysicsVector(theLambdaTable, i, aVector); 325 } 326 } 327 328 if(1 < verboseLevel) { 329 G4cout << "Lambda table is built for " << part->GetParticleName() << G4endl; 330 } 331 } 332 333 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 334 335 const G4ParticleDefinition* 336 G4EmTableUtil::CheckIon(G4VEnergyLossProcess* proc, 337 const G4ParticleDefinition* part, 338 const G4ParticleDefinition* partLocal, 339 const G4int verb, G4bool& isIon) 340 { 341 if(1 < verb) { 342 G4cout << "G4EmTableUtil::CheckIon for " 343 << proc->GetProcessName() << " for " << part->GetParticleName() 344 << " should be called from G4VEnergyLossProcess::PreparePhysicsTable" 345 << G4endl; 346 } 347 const G4ParticleDefinition* particle = partLocal; 348 349 // Are particle defined? 350 if(nullptr == particle) { particle = part; } 351 if(part->GetParticleType() == "nucleus") { 352 G4String pname = part->GetParticleName(); 353 if(pname != "deuteron" && pname != "triton" && 354 pname != "He3" && pname != "alpha+" && pname != "alpha") { 355 356 const G4ParticleDefinition* theGIon = G4GenericIon::GenericIon(); 357 isIon = true; 358 359 // this is a loop to compare pointers of G4GenericIon processes in order 360 // to confirm that for given particle the G4GenericIon physics is used 361 if(particle != theGIon) { 362 G4ProcessManager* pm = theGIon->GetProcessManager(); 363 G4ProcessVector* v = pm->GetAlongStepProcessVector(); 364 G4int n = (G4int)v->size(); 365 for(G4int j=0; j<n; ++j) { 366 if((*v)[j] == proc) { 367 particle = theGIon; 368 break; 369 } 370 } 371 } 372 } 373 } 374 return particle; 375 } 376 377 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 378 379 void G4EmTableUtil::UpdateModels(G4VEnergyLossProcess* proc, 380 G4EmModelManager* modelManager, 381 const G4double maxKinEnergy, 382 const G4int nModels, 383 G4int& secID, G4int& biasID, 384 G4int& mainSec, const G4bool baseMat, 385 const G4bool, const G4bool useAGen) 386 { 387 // defined ID of secondary particles 388 G4int stype = proc->GetProcessSubType(); 389 if(stype == fBremsstrahlung) { 390 secID = _Bremsstrahlung; 391 biasID = _SplitBremsstrahlung; 392 } else if(stype == fPairProdByCharged) { 393 secID = _PairProduction; 394 mainSec = 2; 395 } 396 397 // initialisation of models 398 for(G4int i=0; i<nModels; ++i) { 399 G4VEmModel* mod = modelManager->GetModel(i); 400 mod->SetAngularGeneratorFlag(useAGen); 401 if(mod->HighEnergyLimit() > maxKinEnergy) { 402 mod->SetHighEnergyLimit(maxKinEnergy); 403 } 404 mod->SetUseBaseMaterials(baseMat); 405 } 406 } 407 408 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 409 410 void 411 G4EmTableUtil::BuildLocalElossProcess(G4VEnergyLossProcess* proc, 412 const G4VEnergyLossProcess* masterProc, 413 const G4ParticleDefinition* part, 414 const G4int nModels) 415 { 416 // copy table pointers from master thread 417 proc->SetDEDXTable(masterProc->DEDXTable(),fRestricted); 418 proc->SetDEDXTable(masterProc->DEDXunRestrictedTable(),fTotal); 419 proc->SetDEDXTable(masterProc->IonisationTable(),fIsIonisation); 420 proc->SetRangeTableForLoss(masterProc->RangeTableForLoss()); 421 proc->SetCSDARangeTable(masterProc->CSDARangeTable()); 422 proc->SetInverseRangeTable(masterProc->InverseRangeTable()); 423 proc->SetLambdaTable(masterProc->LambdaTable()); 424 proc->SetCrossSectionType(masterProc->CrossSectionType()); 425 proc->SetEnergyOfCrossSectionMax(masterProc->EnergyOfCrossSectionMax()); 426 proc->SetTwoPeaksXS(masterProc->TwoPeaksXS()); 427 proc->SetIonisation(masterProc->IsIonisationProcess()); 428 G4bool baseMat = masterProc->UseBaseMaterial(); 429 430 // local initialisation of models 431 G4bool printing = true; 432 for(G4int i=0; i<nModels; ++i) { 433 G4VEmModel* mod = proc->GetModelByIndex(i, printing); 434 G4VEmModel* mod0= masterProc->GetModelByIndex(i, printing); 435 mod->SetUseBaseMaterials(baseMat); 436 mod->InitialiseLocal(part, mod0); 437 } 438 } 439 440 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 441 442 void G4EmTableUtil::BuildDEDXTable(G4VEnergyLossProcess* proc, 443 const G4ParticleDefinition* part, 444 G4EmModelManager* modelManager, 445 G4LossTableBuilder* bld, 446 G4PhysicsTable* table, 447 const G4double emin, 448 const G4double emax, 449 const G4int nbins, 450 const G4int verbose, 451 const G4EmTableType tType, 452 const G4bool spline) 453 { 454 // Access to materials 455 const G4ProductionCutsTable* theCoupleTable= 456 G4ProductionCutsTable::GetProductionCutsTable(); 457 std::size_t numOfCouples = theCoupleTable->GetTableSize(); 458 459 if(1 < verbose) { 460 G4cout << numOfCouples << " couples" << " minKinEnergy(MeV)= " << emin 461 << " maxKinEnergy(MeV)= " << emax << " nbins= " << nbins << G4endl; 462 } 463 G4PhysicsLogVector* aVector = nullptr; 464 G4PhysicsLogVector* bVector = nullptr; 465 466 for(std::size_t i=0; i<numOfCouples; ++i) { 467 468 if(1 < verbose) { 469 G4cout << "G4EmTableUtil::BuildDEDXVector idx= " << i 470 << " flagTable=" << table->GetFlag(i) 471 << " flagBuilder=" << bld->GetFlag(i) << G4endl; 472 } 473 if(bld->GetFlag(i)) { 474 475 // create physics vector and fill it 476 const G4MaterialCutsCouple* couple = 477 theCoupleTable->GetMaterialCutsCouple((G4int)i); 478 delete (*table)[i]; 479 if(nullptr != bVector) { 480 aVector = new G4PhysicsLogVector(*bVector); 481 } else { 482 bVector = new G4PhysicsLogVector(emin, emax, nbins, spline); 483 aVector = bVector; 484 } 485 486 modelManager->FillDEDXVector(aVector, couple, tType); 487 if(spline) { aVector->FillSecondDerivatives(); } 488 489 // Insert vector for this material into the table 490 G4PhysicsTableHelper::SetPhysicsVector(table, i, aVector); 491 } 492 } 493 494 if(1 < verbose) { 495 G4cout << "G4EmTableUtil::BuildDEDXTable(): table is built for " 496 << part->GetParticleName() 497 << " and process " << proc->GetProcessName() 498 << G4endl; 499 if(2 < verbose) G4cout << (*table) << G4endl; 500 } 501 } 502 503 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 504 505 void G4EmTableUtil::PrepareMscProcess(G4VMultipleScattering* proc, 506 const G4ParticleDefinition& part, 507 G4EmModelManager* modelManager, 508 G4MscStepLimitType& stepLimit, 509 G4double& facrange, 510 G4bool& latDisplacement, G4bool& master, 511 G4bool& isIon, G4bool& baseMat) 512 { 513 auto param = G4EmParameters::Instance(); 514 G4int verb = (master) ? param->Verbose() : param->WorkerVerbose(); 515 proc->SetVerboseLevel(verb); 516 517 if(part.GetPDGMass() > CLHEP::GeV || 518 part.GetParticleName() == "GenericIon") { isIon = true; } 519 520 if(1 < verb) { 521 G4cout << "### G4EmTableUtil::PrepearPhysicsTable() for " 522 << proc->GetProcessName() 523 << " and particle " << part.GetParticleName() 524 << " isIon: " << isIon << " isMaster: " << master 525 << G4endl; 526 } 527 528 // initialise process 529 proc->InitialiseProcess(&part); 530 531 // heavy particles 532 if(part.GetPDGMass() > CLHEP::MeV) { 533 stepLimit = param->MscMuHadStepLimitType(); 534 facrange = param->MscMuHadRangeFactor(); 535 latDisplacement = param->MuHadLateralDisplacement(); 536 } else { 537 stepLimit = param->MscStepLimitType(); 538 facrange = param->MscRangeFactor(); 539 latDisplacement = param->LateralDisplacement(); 540 } 541 542 // initialisation of models 543 auto numberOfModels = modelManager->NumberOfModels(); 544 for(G4int i=0; i<numberOfModels; ++i) { 545 G4VMscModel* msc = proc->GetModelByIndex(i); 546 msc->SetIonisation(nullptr, &part); 547 msc->SetPolarAngleLimit(param->MscThetaLimit()); 548 G4double emax = std::min(msc->HighEnergyLimit(),param->MaxKinEnergy()); 549 msc->SetHighEnergyLimit(emax); 550 msc->SetUseBaseMaterials(baseMat); 551 } 552 modelManager->Initialise(&part, nullptr, verb); 553 } 554 555 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo..... 556 557 void G4EmTableUtil::BuildMscProcess(G4VMultipleScattering* proc, 558 const G4VMultipleScattering* masterProc, 559 const G4ParticleDefinition& part, 560 const G4ParticleDefinition* firstPart, 561 G4int nModels, G4bool master) 562 { 563 auto param = G4EmParameters::Instance(); 564 G4int verb = param->Verbose(); 565 566 if (firstPart == &part) { 567 G4LossTableBuilder* bld = G4LossTableManager::Instance()->GetTableBuilder(); 568 G4bool baseMat = bld->GetBaseMaterialFlag(); 569 if (master) { 570 for (G4int i=0; i<nModels; ++i) { 571 G4VMscModel* msc = proc->GetModelByIndex(i); 572 msc->SetUseBaseMaterials(baseMat); 573 // table is always built for low mass particles 574 if (part.GetParticleName() != "GenericIon" && 575 (part.GetPDGMass() < CLHEP::GeV || msc->ForceBuildTableFlag())) { 576 G4double emin = 577 std::max(msc->LowEnergyLimit(), msc->LowEnergyActivationLimit()); 578 G4double emax = 579 std::min(msc->HighEnergyLimit(), msc->HighEnergyActivationLimit()); 580 emin = std::max(emin, param->MinKinEnergy()); 581 emax = std::min(emax, param->MaxKinEnergy()); 582 if (emin < emax) { 583 auto table = bld->BuildTableForModel(msc->GetCrossSectionTable(), 584 msc, &part, emin, emax, true); 585 msc->SetCrossSectionTable(table, true); 586 } 587 } 588 } 589 } else { 590 for (G4int i=0; i<nModels; ++i) { 591 G4VMscModel* msc = proc->GetModelByIndex(i); 592 G4VMscModel* msc0 = masterProc->GetModelByIndex(i); 593 msc->SetUseBaseMaterials(baseMat); 594 msc->SetCrossSectionTable(msc0->GetCrossSectionTable(), false); 595 msc->InitialiseLocal(&part, msc0); 596 } 597 } 598 } 599 if(!param->IsPrintLocked()) { 600 const G4String& num = part.GetParticleName(); 601 602 // explicitly defined printout by particle name 603 if(1 < verb || (0 < verb && (num == "e-" || 604 num == "e+" || num == "mu+" || 605 num == "mu-" || num == "proton"|| 606 num == "pi+" || num == "pi-" || 607 num == "kaon+" || num == "kaon-" || 608 num == "alpha" || num == "anti_proton" || 609 num == "GenericIon" || num == "alpha+" || 610 num == "alpha" ))) { 611 proc->StreamInfo(G4cout, part); 612 } 613 } 614 if(1 < verb) { 615 G4cout << "### G4EmTableUtil::BuildPhysicsTable() done for " 616 << proc->GetProcessName() 617 << " and particle " << part.GetParticleName() << G4endl; 618 } 619 } 620 621 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 622 623 G4bool G4EmTableUtil::StoreMscTable(G4VMultipleScattering* proc, 624 const G4ParticleDefinition* part, 625 const G4String& dir, 626 const G4int nModels, const G4int verb, 627 const G4bool ascii) 628 { 629 G4bool ok = true; 630 for(G4int i=0; i<nModels; ++i) { 631 G4VMscModel* msc = proc->GetModelByIndex(i); 632 G4PhysicsTable* table = msc->GetCrossSectionTable(); 633 if (nullptr != table) { 634 G4String ss = G4UIcommand::ConvertToString(i); 635 G4String name = 636 proc->GetPhysicsTableFileName(part, dir, "LambdaMod"+ss, ascii); 637 G4bool yes = table->StorePhysicsTable(name,ascii); 638 639 if ( yes ) { 640 if ( verb > 0 ) { 641 G4cout << "Physics table are stored for " 642 << part->GetParticleName() 643 << " and process " << proc->GetProcessName() 644 << " with a name <" << name << "> " << G4endl; 645 } 646 } else { 647 G4cout << "Fail to store Physics Table for " 648 << part->GetParticleName() 649 << " and process " << proc->GetProcessName() 650 << " in the directory <" << dir 651 << "> " << G4endl; 652 ok = false; 653 } 654 } 655 } 656 return ok; 657 } 658 659 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 660 661 G4bool G4EmTableUtil::StoreTable(G4VProcess* ptr, 662 const G4ParticleDefinition* part, 663 G4PhysicsTable* aTable, 664 const G4String& dir, 665 const G4String& tname, 666 const G4int verb, const G4bool ascii) 667 { 668 G4bool res = true; 669 if (nullptr != aTable) { 670 const G4String& name = 671 ptr->GetPhysicsTableFileName(part, dir, tname, ascii); 672 if ( aTable->StorePhysicsTable(name, ascii) ) { 673 if (1 < verb) G4cout << "Stored: " << name << G4endl; 674 } else { 675 res = false; 676 G4cout << "G4EmTableUtil::StoreTable fail to store: " << name << G4endl; 677 } 678 } 679 return res; 680 } 681 682 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 683 684 G4bool G4EmTableUtil::RetrieveTable(G4VProcess* ptr, 685 const G4ParticleDefinition* part, 686 G4PhysicsTable* aTable, 687 const G4String& dir, const G4String& tname, 688 const G4int verb, const G4bool ascii, 689 const G4bool spline) 690 { 691 G4bool res = true; 692 if (nullptr == aTable) { return res; } 693 if (1 < verb) { 694 G4cout << tname << " table for " << part->GetParticleName() 695 << " will be retrieved " << G4endl; 696 } 697 const G4String& name = 698 ptr->GetPhysicsTableFileName(part, dir, tname, ascii); 699 if(G4PhysicsTableHelper::RetrievePhysicsTable(aTable, name, ascii, spline)) { 700 if(spline) { 701 for(auto & v : *aTable) { 702 if(nullptr != v) { v->FillSecondDerivatives(); } 703 } 704 } 705 if (0 < verb) { 706 G4cout << tname << " table for " << part->GetParticleName() 707 << " is retrieved from <" << name << ">" 708 << G4endl; 709 } 710 } else { 711 res = false; 712 G4cout << "G4EmTableUtil::RetrieveTable fail to retrieve: " << tname 713 << " from " << name << " for " << part->GetParticleName() << G4endl; 714 } 715 return res; 716 } 717 718 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 719 720