Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EmTableUtil Class Reference

#include <G4EmTableUtil.hh>

Static Public Member Functions

static const G4DataVectorPrepareEmProcess (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 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 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 void BuildLambdaTable (G4VEnergyLossProcess *proc, const G4ParticleDefinition *part, G4EmModelManager *modelManager, G4LossTableBuilder *bld, G4PhysicsTable *theLambdaTable, const G4DataVector *theCuts, const G4double minKinEnergy, const G4double maxKinEnergy, const G4double scale, const G4int verbose, const G4bool splineFlag)
 
static const G4ParticleDefinitionCheckIon (G4VEnergyLossProcess *proc, const G4ParticleDefinition *part, const G4ParticleDefinition *particle, const G4int verboseLevel, G4bool &isIon)
 
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 void BuildLocalElossProcess (G4VEnergyLossProcess *proc, const G4VEnergyLossProcess *masterProc, const G4ParticleDefinition *part, const G4int nModels)
 
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 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 G4bool StoreMscTable (G4VMultipleScattering *proc, const G4ParticleDefinition *part, const G4String &directory, const G4int nModels, const G4int verb, const G4bool ascii)
 
static G4bool StoreTable (G4VProcess *, const G4ParticleDefinition *, G4PhysicsTable *, const G4String &dir, const G4String &tname, G4int verb, G4bool ascii)
 
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)
 

Detailed Description

Definition at line 48 of file G4EmTableUtil.hh.

Member Function Documentation

◆ BuildDEDXTable()

void G4EmTableUtil::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

Definition at line 442 of file G4EmTableUtil.cc.

453{
454 // Access to materials
455 const G4ProductionCutsTable* theCoupleTable=
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
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}
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void FillDEDXVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4EmTableType t=fRestricted)
static G4bool GetFlag(std::size_t idx)
const G4String & GetParticleName() const
static void SetPhysicsVector(G4PhysicsTable *physTable, std::size_t idx, G4PhysicsVector *vec)
G4bool GetFlag(std::size_t i) const
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
const G4String & GetProcessName() const

Referenced by G4VEnergyLossProcess::BuildDEDXTable().

◆ BuildEmProcess()

void G4EmTableUtil::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

Definition at line 110 of file G4EmTableUtil.cc.

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());
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 {
157 }
158 if(nullptr == v) { proc->SetCrossSectionType(fEmIncreasing); }
159 }
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}
@ fEmOnePeak
@ fEmIncreasing
bool G4bool
Definition G4Types.hh:86
static std::vector< G4double > * FindCrossSectionMax(G4PhysicsTable *)
void SetUseBaseMaterials(G4bool val)
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
G4CrossSectionType CrossSectionType() const
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
void SetLambdaTablePrim(G4PhysicsTable *)
void SetEnergyOfCrossSectionMax(std::vector< G4double > *)
G4PhysicsTable * LambdaTable() const
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()

Referenced by G4VEmProcess::BuildPhysicsTable().

◆ BuildLambdaTable() [1/2]

void G4EmTableUtil::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

Definition at line 188 of file G4EmTableUtil.cc.

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=
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,
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}
@ fIsCrossSectionPrim
G4double G4Log(G4double x)
Definition G4Log.hh:227
double G4double
Definition G4Types.hh:83
void FillLambdaVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4bool startFromNull=true, G4EmTableType t=fRestricted)
const G4Material * GetMaterial() const
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *)
int G4lrint(double ad)
Definition templates.hh:134

Referenced by G4VEmProcess::BuildLambdaTable(), and G4VEnergyLossProcess::BuildLambdaTable().

◆ BuildLambdaTable() [2/2]

void G4EmTableUtil::BuildLambdaTable ( G4VEnergyLossProcess * proc,
const G4ParticleDefinition * part,
G4EmModelManager * modelManager,
G4LossTableBuilder * bld,
G4PhysicsTable * theLambdaTable,
const G4DataVector * theCuts,
const G4double minKinEnergy,
const G4double maxKinEnergy,
const G4double scale,
const G4int verbose,
const G4bool splineFlag )
static

Definition at line 279 of file G4EmTableUtil.cc.

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=
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}
@ fRestricted
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *, G4double cut)

◆ BuildLocalElossProcess()

void G4EmTableUtil::BuildLocalElossProcess ( G4VEnergyLossProcess * proc,
const G4VEnergyLossProcess * masterProc,
const G4ParticleDefinition * part,
const G4int nModels )
static

Definition at line 411 of file G4EmTableUtil.cc.

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());
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}
@ fTotal
@ fIsIonisation
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 * > *)
std::vector< G4TwoPeaksXS * > * TwoPeaksXS() const
void SetCrossSectionType(G4CrossSectionType val)
void SetInverseRangeTable(G4PhysicsTable *p)
G4CrossSectionType CrossSectionType() const
void SetEnergyOfCrossSectionMax(std::vector< G4double > *)
void SetDEDXTable(G4PhysicsTable *p, G4EmTableType tType)
void SetLambdaTable(G4PhysicsTable *p)
G4PhysicsTable * IonisationTable() const
G4PhysicsTable * LambdaTable() const
void SetCSDARangeTable(G4PhysicsTable *pRange)
G4PhysicsTable * DEDXunRestrictedTable() const
G4PhysicsTable * DEDXTable() const

Referenced by G4VEnergyLossProcess::BuildPhysicsTable().

◆ BuildMscProcess()

void G4EmTableUtil::BuildMscProcess ( G4VMultipleScattering * proc,
const G4VMultipleScattering * masterProc,
const G4ParticleDefinition & part,
const G4ParticleDefinition * firstPart,
G4int nModels,
G4bool master )
static

Definition at line 557 of file G4EmTableUtil.cc.

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}
static G4EmParameters * Instance()
static G4bool GetBaseMaterialFlag()
G4PhysicsTable * BuildTableForModel(G4PhysicsTable *table, G4VEmModel *model, const G4ParticleDefinition *, G4double emin, G4double emax, G4bool spline)
static G4LossTableManager * Instance()
G4LossTableBuilder * GetTableBuilder()
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
G4double HighEnergyActivationLimit() const
G4PhysicsTable * GetCrossSectionTable()
G4double LowEnergyActivationLimit() const
G4bool ForceBuildTableFlag() const
G4VMscModel * GetModelByIndex(G4int idx, G4bool ver=false) const
void StreamInfo(std::ostream &outFile, const G4ParticleDefinition &, G4bool rst=false) const

Referenced by G4VMultipleScattering::BuildPhysicsTable().

◆ CheckIon()

const G4ParticleDefinition * G4EmTableUtil::CheckIon ( G4VEnergyLossProcess * proc,
const G4ParticleDefinition * part,
const G4ParticleDefinition * particle,
const G4int verboseLevel,
G4bool & isIon )
static

Definition at line 336 of file G4EmTableUtil.cc.

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}
static G4GenericIon * GenericIon()
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleType() const
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const

Referenced by G4VEnergyLossProcess::PreparePhysicsTable().

◆ PrepareEmProcess()

const G4DataVector * G4EmTableUtil::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

Definition at line 51 of file G4EmTableUtil.cc.

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}
@ fDNAUnknownModel
@ fGammaConversion
@ fComptonScattering
@ fAnnihilation
@ fPhotoElectricEffect
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4int verb)
G4int NumberOfModels() const
G4VEmModel * GetModel(G4int idx, G4bool ver=false) const
G4double MscThetaLimit() const
G4int Verbose() const
G4int WorkerVerbose() const
void SetPolarAngleLimit(G4double)
void SetHighEnergyLimit(G4double)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetVerboseLevel(G4int value)
G4int GetProcessSubType() const

Referenced by G4VEmProcess::PreparePhysicsTable().

◆ PrepareMscProcess()

void G4EmTableUtil::PrepareMscProcess ( G4VMultipleScattering * proc,
const G4ParticleDefinition & part,
G4EmModelManager * modelManager,
G4MscStepLimitType & stepLimit,
G4double & facrange,
G4bool & latDisplacement,
G4bool & master,
G4bool & isIon,
G4bool & baseMat )
static

Definition at line 505 of file G4EmTableUtil.cc.

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}
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
virtual void InitialiseProcess(const G4ParticleDefinition *)=0

Referenced by G4VMultipleScattering::PreparePhysicsTable().

◆ RetrieveTable()

G4bool G4EmTableUtil::RetrieveTable ( G4VProcess * ptr,
const G4ParticleDefinition * part,
G4PhysicsTable * aTable,
const G4String & dir,
const G4String & tname,
const G4int verb,
const G4bool ascii,
const G4bool spline )
static

Definition at line 684 of file G4EmTableUtil.cc.

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}
static G4bool RetrievePhysicsTable(G4PhysicsTable *physTable, const G4String &fileName, G4bool ascii, G4bool spline)
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
const char * name(G4int ptype)

Referenced by G4VEmProcess::RetrievePhysicsTable(), and G4VEnergyLossProcess::RetrievePhysicsTable().

◆ StoreMscTable()

G4bool G4EmTableUtil::StoreMscTable ( G4VMultipleScattering * proc,
const G4ParticleDefinition * part,
const G4String & directory,
const G4int nModels,
const G4int verb,
const G4bool ascii )
static

Definition at line 623 of file G4EmTableUtil.cc.

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}
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
static G4String ConvertToString(G4bool boolVal)

Referenced by G4VMultipleScattering::StorePhysicsTable().

◆ StoreTable()

G4bool G4EmTableUtil::StoreTable ( G4VProcess * ptr,
const G4ParticleDefinition * part,
G4PhysicsTable * aTable,
const G4String & dir,
const G4String & tname,
G4int verb,
G4bool ascii )
static

Definition at line 661 of file G4EmTableUtil.cc.

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}

Referenced by G4VEmProcess::StorePhysicsTable(), and G4VEnergyLossProcess::StorePhysicsTable().

◆ UpdateModels()

void G4EmTableUtil::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

Definition at line 379 of file G4EmTableUtil.cc.

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}
@ fBremsstrahlung
@ fPairProdByCharged
void SetAngularGeneratorFlag(G4bool)

Referenced by G4VEnergyLossProcess::PreparePhysicsTable().


The documentation for this class was generated from the following files: