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