Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GammaGeneralProcess Class Reference

#include <G4GammaGeneralProcess.hh>

+ Inheritance diagram for G4GammaGeneralProcess:

Public Member Functions

 G4GammaGeneralProcess ()
 
virtual ~G4GammaGeneralProcess ()
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 
void AddEmProcess (G4VEmProcess *)
 
void AddMMProcess (G4GammaConversionToMuons *)
 
void AddHadProcess (G4HadronicProcess *)
 
void ProcessDescription (std::ostream &outFile) const override
 
void PreparePhysicsTable (const G4ParticleDefinition &) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void StartTracking (G4Track *) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
 
G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
 
const G4StringGetProcessName () const
 
G4int GetProcessSubType () const
 
G4VEmProcessGetEmProcess (const G4String &name) override
 
 G4GammaGeneralProcess (G4GammaGeneralProcess &)=delete
 
G4GammaGeneralProcessoperator= (const G4GammaGeneralProcess &right)=delete
 
- Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEmProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p) override=0
 
virtual void PrintInfo ()
 
virtual void ProcessDescription (std::ostream &outFile) const override
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
virtual void StartTracking (G4Track *) override
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
 
virtual G4VEmProcessGetEmProcess (const G4String &name)
 
G4double CrossSectionPerVolume (G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
 
G4double MeanFreePath (const G4Track &track)
 
G4double GetLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple)
 
G4double GetLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy)
 
void SetLambdaBinning (G4int nbins)
 
void SetMinKinEnergy (G4double e)
 
void SetMinKinEnergyPrim (G4double e)
 
void SetMaxKinEnergy (G4double e)
 
G4PhysicsTableLambdaTable () const
 
G4PhysicsTableLambdaTablePrim () const
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
G4VEmModelSelectModelForMaterial (G4double kinEnergy, size_t idxRegion) const
 
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=nullptr)
 
void SetEmModel (G4VEmModel *, G4int index=0)
 
G4VEmModelEmModel (size_t index=0) const
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
G4int GetNumberOfModels () const
 
G4int GetNumberOfRegionModels (size_t couple_index) const
 
G4VEmModelGetRegionModel (G4int idx, size_t couple_index) const
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
const G4VEmModelGetCurrentModel () const
 
const G4ElementGetCurrentElement () const
 
void SetCrossSectionBiasingFactor (G4double f, G4bool flag=true)
 
G4double CrossSectionBiasingFactor () const
 
void ActivateForcedInteraction (G4double length=0.0, const G4String &r="", G4bool flag=true)
 
void ActivateSecondaryBiasing (const G4String &region, G4double factor, G4double energyLimit)
 
void SetEmMasterProcess (const G4VEmProcess *)
 
void SetIntegral (G4bool val)
 
void SetBuildTableFlag (G4bool val)
 
void CurrentSetup (const G4MaterialCutsCouple *, G4double energy)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

void InitialiseProcess (const G4ParticleDefinition *) override
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4double ComputeGeneralLambda (size_t idxe, size_t idxt)
 
G4double GetProbability (size_t idxt)
 
void SelectedProcess (const G4Step &track, const G4VProcess *ptr)
 
G4VParticleChangeSampleEmSecondaries (const G4Track &, const G4Step &, G4VEmProcess *)
 
G4VParticleChangeSampleHadSecondaries (const G4Track &, const G4Step &, G4HadronicProcess *)
 
- Protected Member Functions inherited from G4VEmProcess
virtual void StreamProcessInfo (std::ostream &) const
 
virtual void InitialiseProcess (const G4ParticleDefinition *)=0
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
G4VEmModelSelectModel (G4double kinEnergy, size_t index)
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
 
void DefineMaterial (const G4MaterialCutsCouple *couple)
 
G4int LambdaBinning () const
 
G4double MinKinEnergy () const
 
G4double MaxKinEnergy () const
 
G4double PolarAngleLimit () const
 
G4bool IsIntegral () const
 
G4double RecalculateLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple)
 
G4ParticleChangeForGammaGetParticleChange ()
 
void SetParticle (const G4ParticleDefinition *p)
 
void SetSecondaryParticle (const G4ParticleDefinition *p)
 
size_t CurrentMaterialCutsCoupleIndex () const
 
const G4MaterialCutsCoupleMaterialCutsCouple () const
 
G4bool ApplyCuts () const
 
G4double GetGammaEnergyCut ()
 
G4double GetElectronEnergyCut ()
 
void SetStartFromNullFlag (G4bool val)
 
void SetSplineFlag (G4bool val)
 
const G4ElementGetTargetElement () const
 
const G4IsotopeGetTargetIsotope () const
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4HadronicProcesstheGammaNuclear
 
const G4VProcessselectedProc
 
- Protected Attributes inherited from G4VEmProcess
G4LossTableManagerlManager
 
G4EmParameterstheParameters
 
G4EmBiasingManagerbiasManager
 
const G4ParticleDefinitiontheGamma
 
const G4ParticleDefinitiontheElectron
 
G4ParticleChangeForGamma fParticleChange
 
std::vector< G4DynamicParticle * > secParticles
 
const G4MaterialCutsCouplecurrentCouple
 
const G4MaterialcurrentMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t currentCoupleIndex
 
size_t basedCoupleIndex
 
G4int mainSecondaries
 
G4int secID
 
G4int fluoID
 
G4int augerID
 
G4int biasID
 
G4bool isTheMaster
 
G4double mfpKinEnergy
 
G4double preStepKinEnergy
 
G4double preStepLogKinEnergy
 
G4double preStepLambda
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Detailed Description

Definition at line 66 of file G4GammaGeneralProcess.hh.

Constructor & Destructor Documentation

◆ G4GammaGeneralProcess() [1/2]

G4GammaGeneralProcess::G4GammaGeneralProcess ( )
explicit

Definition at line 85 of file G4GammaGeneralProcess.cc.

85 :
86 G4VEmProcess("GammaGeneralProc", fElectromagnetic),
87 minPEEnergy(50*CLHEP::keV),
88 minEEEnergy(2*CLHEP::electron_mass_c2),
89 minMMEnergy(100*CLHEP::MeV),
90 peLambda(0.0),
91 nLowE(40),
92 nHighE(50),
93 splineFlag(false)
94{
95 thePhotoElectric = theCompton = theConversionEE = theRayleigh = nullptr;
96 theGammaNuclear = nullptr;
97 theConversionMM = nullptr;
98 selectedProc = nullptr;
99
100 idxEnergy = 0;
101 preStepLogE = factor = 1.0;
102
106}
@ fGammaGeneralProcess
@ fElectromagnetic
const G4VProcess * selectedProc
G4HadronicProcess * theGammaNuclear
const G4ParticleDefinition * theGamma
void SetParticle(const G4ParticleDefinition *p)
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406

◆ ~G4GammaGeneralProcess()

G4GammaGeneralProcess::~G4GammaGeneralProcess ( )
virtual

Definition at line 110 of file G4GammaGeneralProcess.cc.

111{
112 //std::cout << "G4GammaGeneralProcess::~G4GammaGeneralProcess " << isTheMaster
113 // << " " << theHandler << G4endl;
114 if(isTheMaster) {
115 delete theHandler;
116 theHandler = nullptr;
117 }
118}
G4bool isTheMaster

◆ G4GammaGeneralProcess() [2/2]

G4GammaGeneralProcess::G4GammaGeneralProcess ( G4GammaGeneralProcess )
delete

Member Function Documentation

◆ AddEmProcess()

void G4GammaGeneralProcess::AddEmProcess ( G4VEmProcess ptr)

Definition at line 129 of file G4GammaGeneralProcess.cc.

130{
131 if(!ptr) { return; }
132 G4int stype = ptr->GetProcessSubType();
133 if(stype == fRayleigh) { theRayleigh = ptr; }
134 else if(stype == fPhotoElectricEffect) { thePhotoElectric = ptr; }
135 else if(stype == fComptonScattering) { theCompton = ptr; }
136 else if(stype == fGammaConversion) { theConversionEE = ptr; }
137}
@ fGammaConversion
@ fRayleigh
@ fComptonScattering
@ fPhotoElectricEffect
int G4int
Definition: G4Types.hh:85
G4int GetProcessSubType() const
Definition: G4VProcess.hh:400

◆ AddHadProcess()

void G4GammaGeneralProcess::AddHadProcess ( G4HadronicProcess ptr)

Definition at line 148 of file G4GammaGeneralProcess.cc.

149{
150 theGammaNuclear = ptr;
151}

◆ AddMMProcess()

void G4GammaGeneralProcess::AddMMProcess ( G4GammaConversionToMuons ptr)

Definition at line 141 of file G4GammaGeneralProcess.cc.

142{
143 theConversionMM = ptr;
144}

◆ BuildPhysicsTable()

void G4GammaGeneralProcess::BuildPhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4VEmProcess.

Definition at line 254 of file G4GammaGeneralProcess.cc.

255{
256 if(1 < verboseLevel) {
257 G4cout << "### G4VEmProcess::BuildPhysicsTable() for "
258 << GetProcessName()
259 << " and particle " << part.GetParticleName()
260 << G4endl;
261 }
262 if(thePhotoElectric != nullptr) {
263 if(!isTheMaster) {
264 thePhotoElectric->SetEmMasterProcess(theHandler->GetMasterProcess(0));
265 }
266 thePhotoElectric->BuildPhysicsTable(part);
267 }
268 if(theCompton != nullptr) {
269 if(!isTheMaster) {
270 theCompton->SetEmMasterProcess(theHandler->GetMasterProcess(1));
271 }
272 theCompton->BuildPhysicsTable(part);
273 }
274 if(theConversionEE != nullptr) {
275 if(!isTheMaster) {
276 theConversionEE->SetEmMasterProcess(theHandler->GetMasterProcess(2));
277 }
278 theConversionEE->BuildPhysicsTable(part);
279 }
280 if(theRayleigh != nullptr) {
281 if(!isTheMaster) {
282 theRayleigh->SetEmMasterProcess(theHandler->GetMasterProcess(3));
283 }
284 theRayleigh->BuildPhysicsTable(part);
285 }
286 if(theGammaNuclear != nullptr) { theGammaNuclear->BuildPhysicsTable(part); }
287 if(theConversionMM != nullptr) { theConversionMM->BuildPhysicsTable(part); }
288
289 if(isTheMaster) {
290 const G4ProductionCutsTable* theCoupleTable=
292 size_t numOfCouples = theCoupleTable->GetTableSize();
293
295 const std::vector<G4PhysicsTable*>& tables = theHandler->GetTables();
296
299 G4DynamicParticle* dynParticle =
301
302 G4double sigComp(0.), sigPE(0.), sigConv(0.), sigR(0.),
303 sigN(0.), sigM(0.), val(0.);
304
305 for(size_t i=0; i<numOfCouples; ++i) {
306
307 if (bld->GetFlag(i)) {
308
309 G4int idx = (*theDensityIdx)[i];
310 const G4MaterialCutsCouple* couple =
311 theCoupleTable->GetMaterialCutsCouple(i);
312 const G4Material* material = couple->GetMaterial();
313
314 // energy interval 0
315 size_t nn = (*(tables[0]))[idx]->GetVectorLength();
316 if(1 < verboseLevel) {
317 G4cout << "======= Zone 0 ======= N= " << nn
318 << " for " << material->GetName() << G4endl;
319 }
320 for(size_t j=0; j<nn; ++j) {
321 G4double e = (*(tables[0]))[idx]->Energy(j);
322 G4double loge = G4Log(e);
323 sigComp = (theCompton) ? theCompton->GetLambda(e, couple, loge) : 0.0;
324 sigR = (theRayleigh) ? theRayleigh->GetLambda(e, couple, loge) : 0.0;
325 G4double sum = sigComp + sigR;
326 if(1 < verboseLevel) {
327 G4cout << j << ". E= " << e << " xs= " << sum
328 << " compt= " << sigComp << " Rayl= " << sigR << G4endl;
329 }
330 (*(tables[0]))[idx]->PutValue(j, sum);
331 if(theT[1]) {
332 val = (sum > 0.0) ? sigComp/sum : 0.0;
333 (*(tables[1]))[idx]->PutValue(j, val);
334 }
335 }
336
337 // energy interval 1
338 nn = (*(tables[2]))[idx]->GetVectorLength();
339 if(1 < verboseLevel) {
340 G4cout << "======= Zone 1 ======= N= " << nn << G4endl;
341 }
342 for(size_t j=0; j<nn; ++j) {
343 G4double e = (*(tables[2]))[idx]->Energy(j);
344 G4double loge = G4Log(e);
345 sigComp = (theCompton) ? theCompton->GetLambda(e, couple, loge) : 0.0;
346 sigR = (theRayleigh) ? theRayleigh->GetLambda(e, couple, loge) : 0.0;
347 sigPE = (thePhotoElectric)
348 ? thePhotoElectric->GetLambda(e, couple, loge) : 0.0;
349 sigN = 0.0;
350 if(gn) {
351 dynParticle->SetKineticEnergy(e);
352 sigN = gn->ComputeCrossSection(dynParticle, material);
353 }
354 G4double sum = sigComp + sigR + sigPE + sigN;
355 if(1 < verboseLevel) {
356 G4cout << j << ". E= " << e << " xs= " << sum
357 << " compt= " << sigComp << " conv= " << sigConv
358 << " PE= " << sigPE << " Rayl= " << sigR
359 << " GN= " << sigN << G4endl;
360 }
361 (*(tables[2]))[idx]->PutValue(j, sum);
362 val = (sum > 0.0) ? sigPE/sum : 0.0;
363 (*(tables[3]))[idx]->PutValue(j, val);
364 if(theT[4]) {
365 val = (sum > 0.0) ? (sigComp + sigPE)/sum : 0.0;
366 (*(tables[4]))[idx]->PutValue(j, val);
367 }
368 if(theT[5]) {
369 val = (sum > 0.0) ? (sigComp + sigPE + sigR)/sum : 0.0;
370 (*(tables[5]))[idx]->PutValue(j, val);
371 }
372 }
373
374 // energy interval 2
375 nn = (*(tables[6]))[idx]->GetVectorLength();
376 if(1 < verboseLevel) {
377 G4cout << "======= Zone 2 ======= N= " << nn << G4endl;
378 }
379 for(size_t j=0; j<nn; ++j) {
380 G4double e = (*(tables[6]))[idx]->Energy(j);
381 G4double loge = G4Log(e);
382 sigComp = (theCompton) ? theCompton->GetLambda(e, couple, loge) : 0.0;
383 sigConv = (theConversionEE)
384 ? theConversionEE->GetLambda(e, couple, loge) : 0.0;
385 sigPE = (thePhotoElectric)
386 ? thePhotoElectric->GetLambda(e, couple, loge) : 0.0;
387 sigN = 0.0;
388 if(gn) {
389 dynParticle->SetKineticEnergy(e);
390 sigN = gn->ComputeCrossSection(dynParticle, material);
391 }
392 G4double sum = sigComp + sigConv + sigPE + sigN;
393 if(1 < verboseLevel) {
394 G4cout << j << ". E= " << e << " xs= " << sum
395 << " compt= " << sigComp << " conv= " << sigConv
396 << " PE= " << sigPE
397 << " GN= " << sigN << G4endl;
398 }
399 (*(tables[6]))[idx]->PutValue(j, sum);
400 val = (sum > 0.0) ? sigConv/sum : 0.0;
401 (*(tables[7]))[idx]->PutValue(j, val);
402 val = (sum > 0.0) ? (sigConv + sigComp)/sum : 0.0;
403 (*(tables[8]))[idx]->PutValue(j, val);
404 if(theT[9]) {
405 val = (sum > 0.0) ? (sigConv + sigComp + sigPE)/sum : 0.0;
406 (*(tables[9]))[idx]->PutValue(j, val);
407 }
408 }
409
410 // energy interval 3
411 nn = (*(tables[10]))[idx]->GetVectorLength();
412 if(1 < verboseLevel) {
413 G4cout << "======= Zone 3 ======= N= " << nn
414 << " for " << material->GetName() << G4endl;
415 }
416 for(size_t j=0; j<nn; ++j) {
417 G4double e = (*(tables[10]))[idx]->Energy(j);
418 G4double loge = G4Log(e);
419 sigComp = (theCompton) ? theCompton->GetLambda(e, couple, loge) : 0.0;
420 sigConv = (theConversionEE)
421 ? theConversionEE->GetLambda(e, couple, loge) : 0.0;
422 sigPE = (thePhotoElectric)
423 ? thePhotoElectric->GetLambda(e, couple, loge) : 0.0;
424 sigN = 0.0;
425 if(gn) {
426 dynParticle->SetKineticEnergy(e);
427 sigN = gn->ComputeCrossSection(dynParticle, material);
428 }
429 sigM = 0.0;
430 if(theConversionMM) {
431 val = theConversionMM->ComputeMeanFreePath(e, material);
432 sigM = (val < DBL_MAX) ? 1./val : 0.0;
433 }
434 G4double sum = sigComp + sigConv + sigPE + sigN + sigM;
435 if(1 < verboseLevel) {
436 G4cout << j << ". E= " << e << " xs= " << sum
437 << " compt= " << sigComp << " conv= " << sigConv
438 << " PE= " << sigPE
439 << " GN= " << sigN << G4endl;
440 }
441 (*(tables[10]))[idx]->PutValue(j, sum);
442 val = (sum > 0.0) ? 1.0 - sigConv/sum : 1.0;
443 (*(tables[11]))[idx]->PutValue(j, val);
444 val = (sum > 0.0) ? 1.0 - (sigConv + sigComp)/sum : 1.0;
445 (*(tables[12]))[idx]->PutValue(j, val);
446 if(theT[13]) {
447 val = (sum > 0.0) ? 1.0 - (sigConv + sigComp + sigPE)/sum : 1.0;
448 (*(tables[13]))[idx]->PutValue(j, val);
449 }
450 if(theT[14]) {
451 val = (sum > 0.0)
452 ? 1.0 - (sigConv + sigComp + sigPE + sigN)/sum : 1.0;
453 (*(tables[14]))[idx]->PutValue(j, val);
454 }
455 }
456 for(size_t k=0; k<nTables; ++k) {
457 if(theT[k] && splineFlag) {
458 (*(tables[k]))[idx]->FillSecondDerivatives();
459 }
460 }
461 }
462 }
463 delete dynParticle;
464 }
465
466 if(1 < verboseLevel) {
467 G4cout << "### G4VEmProcess::BuildPhysicsTable() done for "
468 << GetProcessName()
469 << " and particle " << part.GetParticleName()
470 << G4endl;
471 }
472}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
void SetKineticEnergy(G4double aEnergy)
const G4VEmProcess * GetMasterProcess(size_t idx) const
const std::vector< G4PhysicsTable * > & GetTables() const
G4double ComputeMeanFreePath(G4double GammaEnergy, const G4Material *aMaterial)
void BuildPhysicsTable(const G4ParticleDefinition &) override
const G4String & GetProcessName() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4CrossSectionDataStore * GetCrossSectionDataStore()
G4bool GetFlag(size_t idx)
G4LossTableBuilder * GetTableBuilder()
const G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:175
const G4String & GetParticleName() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4LossTableManager * lManager
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
void SetEmMasterProcess(const G4VEmProcess *)
G4double GetLambda(G4double kinEnergy, const G4MaterialCutsCouple *couple)
G4int verboseLevel
Definition: G4VProcess.hh:356
#define DBL_MAX
Definition: templates.hh:62

◆ ComputeGeneralLambda()

G4double G4GammaGeneralProcess::ComputeGeneralLambda ( size_t  idxe,
size_t  idxt 
)
inlineprotected

Definition at line 189 of file G4GammaGeneralProcess.hh.

190{
191 idxEnergy = idxe;
192 return factor*theHandler->GetVector(idxt, basedCoupleIndex)
193 ->LogVectorValue(preStepKinEnergy, preStepLogE);
194}
const G4PhysicsVector * GetVector(size_t itable, size_t ivec) const
G4double LogVectorValue(const G4double theEnergy, const G4double theLogEnergy) const
size_t basedCoupleIndex
G4double preStepKinEnergy

◆ GetEmProcess()

G4VEmProcess * G4GammaGeneralProcess::GetEmProcess ( const G4String name)
overridevirtual

Reimplemented from G4VEmProcess.

Definition at line 797 of file G4GammaGeneralProcess.cc.

798{
799 G4VEmProcess* proc = nullptr;
800 if(thePhotoElectric != nullptr && name == thePhotoElectric->GetProcessName()) {
801 proc = thePhotoElectric;
802 } else if(theCompton != nullptr && name == theCompton->GetProcessName()) {
803 proc = theCompton;
804 } else if(theConversionEE != nullptr && name == theConversionEE->GetProcessName()) {
805 proc = theConversionEE;
806 } else if(theRayleigh != nullptr && name == theRayleigh->GetProcessName()) {
807 proc = theRayleigh;
808 }
809 return proc;
810}
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

◆ GetMeanFreePath()

G4double G4GammaGeneralProcess::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
overrideprotectedvirtual

Reimplemented from G4VEmProcess.

Definition at line 758 of file G4GammaGeneralProcess.cc.

761{
763 return MeanFreePath(track);
764}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
G4double MeanFreePath(const G4Track &track)

◆ GetProbability()

G4double G4GammaGeneralProcess::GetProbability ( size_t  idxt)
inlineprotected

Definition at line 198 of file G4GammaGeneralProcess.hh.

199{
200 return (theT[idxt]) ? theHandler->GetVector(idxt, basedCoupleIndex)
201 ->LogVectorValue(preStepKinEnergy, preStepLogE) : 1.0;
202}

Referenced by PostStepDoIt().

◆ GetProcessName()

const G4String & G4GammaGeneralProcess::GetProcessName ( ) const

Definition at line 780 of file G4GammaGeneralProcess.cc.

781{
782 //G4cout << "GetProcessName(): " << selectedProc << G4endl;
785}

Referenced by BuildPhysicsTable(), PreparePhysicsTable(), and RetrievePhysicsTable().

◆ GetProcessSubType()

G4int G4GammaGeneralProcess::GetProcessSubType ( ) const

Definition at line 789 of file G4GammaGeneralProcess.cc.

◆ InitialiseProcess()

void G4GammaGeneralProcess::InitialiseProcess ( const G4ParticleDefinition )
overrideprotectedvirtual

Implements G4VEmProcess.

Definition at line 186 of file G4GammaGeneralProcess.cc.

187{
188 if(isTheMaster) {
189
190 // tables are created and its size is defined only once
191 if(!theHandler) {
192 theHandler = new G4EmDataHandler(nTables);
193 if(theRayleigh) { theT[1] = theT[4] = true; }
194 if(theGammaNuclear) { theT[9] = theT[13] = true; }
195 if(theConversionMM) { theT[14] = true; }
196
197 theHandler->SetMasterProcess(thePhotoElectric);
198 theHandler->SetMasterProcess(theCompton);
199 theHandler->SetMasterProcess(theConversionEE);
200 theHandler->SetMasterProcess(theRayleigh);
201 }
202 auto bld = lManager->GetTableBuilder();
203
204 const G4ProductionCutsTable* theCoupleTable=
206 size_t numOfCouples = theCoupleTable->GetTableSize();
207
211 size_t nbin1 = std::max(5, nd*G4lrint(std::log10(minPEEnergy/mine)));
212 size_t nbin2 = std::max(5, nd*G4lrint(std::log10(maxe/minMMEnergy)));
213
214 G4PhysicsVector* vec = nullptr;
215 G4PhysicsLogVector aVector(mine,minPEEnergy,nbin1);
216 G4PhysicsLogVector bVector(minPEEnergy,minEEEnergy,nLowE);
217 G4PhysicsLogVector cVector(minEEEnergy,minMMEnergy,nHighE);
218 G4PhysicsLogVector dVector(minMMEnergy,maxe,nbin2);
219 if(splineFlag) {
220 aVector.SetSpline(splineFlag);
221 bVector.SetSpline(splineFlag);
222 cVector.SetSpline(splineFlag);
223 dVector.SetSpline(splineFlag);
224 }
225
226 for(size_t i=0; i<nTables; ++i) {
227 //G4cout << "## PreparePhysTable " << i << "." << G4endl;
228 if(theT[i]) {
229 G4PhysicsTable* table = theHandler->MakeTable(i);
230 //G4cout << " make table " << table << G4endl;
231 for(size_t j=0; j<numOfCouples; ++j) {
232 vec = (*table)[j];
233 if (bld->GetFlag(j) && !vec) {
234 //G4cout << " i= " << i << " j= " << j << " make new vector" << G4endl;
235 if(i<=1) {
236 vec = new G4PhysicsVector(aVector);
237 } else if(i<=5) {
238 vec = new G4PhysicsVector(bVector);
239 } else if(i<=9) {
240 vec = new G4PhysicsVector(cVector);
241 } else {
242 vec = new G4PhysicsVector(dVector);
243 }
245 }
246 }
247 }
248 }
249 }
250}
void SetMasterProcess(const G4VEmProcess *)
G4PhysicsTable * MakeTable(size_t idx)
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4double MaxKinEnergy() const
static void SetPhysicsVector(G4PhysicsTable *physTable, std::size_t idx, G4PhysicsVector *vec)
G4EmParameters * theParameters
int G4lrint(double ad)
Definition: templates.hh:134

Referenced by PreparePhysicsTable().

◆ IsApplicable()

G4bool G4GammaGeneralProcess::IsApplicable ( const G4ParticleDefinition )
overridevirtual

Implements G4VEmProcess.

Definition at line 122 of file G4GammaGeneralProcess.cc.

123{
124 return true;
125}

◆ operator=()

G4GammaGeneralProcess & G4GammaGeneralProcess::operator= ( const G4GammaGeneralProcess right)
delete

◆ PostStepDoIt()

G4VParticleChange * G4GammaGeneralProcess::PostStepDoIt ( const G4Track track,
const G4Step step 
)
overridevirtual

Reimplemented from G4VEmProcess.

Definition at line 583 of file G4GammaGeneralProcess.cc.

585{
586 // In all cases clear number of interaction lengths
588 G4double q = G4UniformRand();
590 G4double p;
591 /*
592 G4cout << "PostStep: preStepLambda= " << preStepLambda << " x= " << x
593 << " PE= " << peLambda << " q= " << q << " idxE= " << idxEnergy
594 << G4endl;
595 */
596 switch (idxEnergy) {
597 case 0:
598 if(x <= peLambda) {
599 return SampleEmSecondaries(track, step, thePhotoElectric);
600 } else {
601 p = GetProbability(1);
602 if(x <= peLambda + (preStepLambda - peLambda)*p) {
603 return SampleEmSecondaries(track, step, theCompton);
604 } else if(theRayleigh != nullptr) {
605 return SampleEmSecondaries(track, step, theRayleigh);
606 }
607 }
608 break;
609
610 case 1:
611 p = GetProbability(3);
612 if(q <= p) {
613 return SampleEmSecondaries(track, step, thePhotoElectric);
614 }
615 p = GetProbability(4);
616 if(q <= p) {
617 return SampleEmSecondaries(track, step, theCompton);
618 }
619 p = GetProbability(5);
620 if(q <= p) {
621 if(theRayleigh != nullptr) {
622 return SampleEmSecondaries(track, step, theRayleigh);
623 }
624 } else if(theGammaNuclear != nullptr) {
625 return SampleHadSecondaries(track, step, theGammaNuclear);
626 }
627 break;
628
629 case 2:
630 p = GetProbability(7);
631 if(q <= p) {
632 return SampleEmSecondaries(track, step, theConversionEE);
633 }
634 p = GetProbability(8);
635 if(q <= p) {
636 return SampleEmSecondaries(track, step, theCompton);
637 }
638 p = GetProbability(9);
639 if(q <= p) {
640 return SampleEmSecondaries(track, step, thePhotoElectric);
641 } else if(theGammaNuclear != nullptr) {
642 return SampleHadSecondaries(track, step, theGammaNuclear);
643 }
644 break;
645
646 case 3:
647 p = 1.0 - GetProbability(11);
648 if(q <= p) {
649 return SampleEmSecondaries(track, step, theConversionEE);
650 }
651 p = 1.0 - GetProbability(12);
652 if(q <= p) {
653 return SampleEmSecondaries(track, step, theCompton);
654 }
655 p = 1.0 - GetProbability(13);
656 if(q <= p) {
657 return SampleEmSecondaries(track, step, thePhotoElectric);
658 }
659 p = 1.0 - GetProbability(14);
660 if(q <= p) {
661 if(theGammaNuclear != nullptr) {
662 return SampleHadSecondaries(track, step, theGammaNuclear);
663 }
664 } else if(theConversionMM != nullptr) {
665 SelectedProcess(step, theConversionMM);
666 return theConversionMM->PostStepDoIt(track, step);
667 }
668 break;
669 }
670 // no interaction
672 return &fParticleChange;
673}
#define G4UniformRand()
Definition: Randomize.hh:52
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4double GetProbability(size_t idxt)
G4VParticleChange * SampleHadSecondaries(const G4Track &, const G4Step &, G4HadronicProcess *)
void SelectedProcess(const G4Step &track, const G4VProcess *ptr)
G4VParticleChange * SampleEmSecondaries(const G4Track &, const G4Step &, G4VEmProcess *)
void InitializeForPostStep(const G4Track &)
G4double preStepLambda
G4ParticleChangeForGamma fParticleChange
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:331

◆ PostStepGetPhysicalInteractionLength()

G4double G4GammaGeneralProcess::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
overridevirtual

Reimplemented from G4VEmProcess.

Definition at line 484 of file G4GammaGeneralProcess.cc.

488{
490 G4double x = DBL_MAX;
491
494 const G4Material* mat = currentCouple->GetMaterial();
495
496 // compute mean free path
497 if(mat != currentMaterial || energy != preStepKinEnergy) {
499 basedCoupleIndex = (*theDensityIdx)[currentCoupleIndex];
500 factor = (*theDensityFactor)[currentCoupleIndex];
501 currentMaterial = mat;
503 preStepLogE = track.GetDynamicParticle()->GetLogKineticEnergy();
504 preStepLambda = TotalCrossSectionPerVolume();
505
506 // zero cross section
507 if(preStepLambda <= 0.0) {
510 }
511 }
512
513 // non-zero cross section
514 if(preStepLambda > 0.0) {
515
517
518 // beggining of tracking (or just after DoIt of this process)
521
522 } else if(currentInteractionLength < DBL_MAX) {
523
525 previousStepSize/currentInteractionLength;
528 }
529
530 // new mean free path and step limit for the next step
533 }
534 /*
535 G4cout << "PostStepGetPhysicalInteractionLength: e= " << energy
536 << " idxe= " << idxEnergy << " xs= " << preStepLambda
537 << " x= " << x << G4endl;
538 */
539 return x;
540}
G4double GetLogKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
const G4MaterialCutsCouple * currentCouple
size_t currentCoupleIndex
const G4Material * currentMaterial
G4double currentInteractionLength
Definition: G4VProcess.hh:335
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:338
G4double energy(const ThreeVector &p, const G4double m)

◆ PreparePhysicsTable()

void G4GammaGeneralProcess::PreparePhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4VEmProcess.

Definition at line 155 of file G4GammaGeneralProcess.cc.

156{
157 SetParticle(&part);
158 currentCouple = nullptr;
159 currentMaterial = nullptr;
160 preStepLambda = 0.0;
161 idxEnergy = 0;
162
166
167 if(1 < verboseLevel) {
168 G4cout << "G4GammaGeneralProcess::PreparePhysicsTable() for "
169 << GetProcessName()
170 << " and particle " << part.GetParticleName()
171 << " isMaster: " << isTheMaster << G4endl;
172 }
173
174 if(thePhotoElectric) { thePhotoElectric->PreparePhysicsTable(part); }
175 if(theCompton) { theCompton->PreparePhysicsTable(part); }
176 if(theConversionEE) { theConversionEE->PreparePhysicsTable(part); }
177 if(theRayleigh) { theRayleigh->PreparePhysicsTable(part); }
179 if(theConversionMM) { theConversionMM->PreparePhysicsTable(part); }
180
181 InitialiseProcess(&part);
182}
G4int Verbose() const
G4int WorkerVerbose() const
void InitialiseProcess(const G4ParticleDefinition *) override
void PreparePhysicsTable(const G4ParticleDefinition &) override
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:194

◆ ProcessDescription()

void G4GammaGeneralProcess::ProcessDescription ( std::ostream &  outFile) const
overridevirtual

Reimplemented from G4VEmProcess.

Definition at line 768 of file G4GammaGeneralProcess.cc.

769{
770 if(thePhotoElectric) { thePhotoElectric->ProcessDescription(out); }
771 if(theCompton) { theCompton->ProcessDescription(out); }
772 if(theConversionEE) { theConversionEE->ProcessDescription(out); }
773 if(theRayleigh) { theRayleigh->ProcessDescription(out); }
775 if(theConversionMM) { theConversionMM->ProcessDescription(out); }
776}
void ProcessDescription(std::ostream &outFile) const override
virtual void ProcessDescription(std::ostream &outFile) const override
virtual void ProcessDescription(std::ostream &outfile) const
Definition: G4VProcess.cc:175

◆ RetrievePhysicsTable()

G4bool G4GammaGeneralProcess::RetrievePhysicsTable ( const G4ParticleDefinition part,
const G4String directory,
G4bool  ascii 
)
overridevirtual

Reimplemented from G4VEmProcess.

Definition at line 721 of file G4GammaGeneralProcess.cc.

724{
725 if(1 < verboseLevel) {
726 G4cout << "G4GammaGeneralProcess::RetrievePhysicsTable() for "
727 << part->GetParticleName() << " and process "
728 << GetProcessName() << G4endl;
729 }
730 G4bool yes = true;
731 if(thePhotoElectric != nullptr &&
732 !thePhotoElectric->RetrievePhysicsTable(part, directory, ascii))
733 { yes = false; }
734 if(theCompton != nullptr &&
735 !theCompton->RetrievePhysicsTable(part, directory, ascii))
736 { yes = false; }
737 if(theConversionEE != nullptr &&
738 !theConversionEE->RetrievePhysicsTable(part, directory, ascii))
739 { yes = false; }
740 if(theRayleigh != nullptr &&
741 !theRayleigh->RetrievePhysicsTable(part, directory, ascii))
742 { yes = false; }
743
744 for(size_t i=0; i<nTables; ++i) {
745 if(theT[i]) {
746 G4String nam = (0==i || 2==i || 6==i || 10==i)
747 ? "LambdaGeneral" + nameT[i] : "ProbGeneral" + nameT[i];
748 G4String fnam = GetPhysicsTableFileName(part,directory,nam,ascii);
749 if(!theHandler->RetrievePhysicsTable(i, part, fnam, ascii))
750 { yes = false; }
751 }
752 }
753 return yes;
754}
bool G4bool
Definition: G4Types.hh:86
G4bool RetrievePhysicsTable(size_t idx, const G4ParticleDefinition *part, const G4String &fname, G4bool ascii)
virtual G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:182

◆ SampleEmSecondaries()

G4VParticleChange * G4GammaGeneralProcess::SampleEmSecondaries ( const G4Track track,
const G4Step step,
G4VEmProcess proc 
)
inlineprotected

Definition at line 215 of file G4GammaGeneralProcess.hh.

217{
219 SelectedProcess(step, proc);
220 return proc->PostStepDoIt(track, step);
221}
void CurrentSetup(const G4MaterialCutsCouple *, G4double energy)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override

Referenced by PostStepDoIt().

◆ SampleHadSecondaries()

G4VParticleChange * G4GammaGeneralProcess::SampleHadSecondaries ( const G4Track track,
const G4Step step,
G4HadronicProcess proc 
)
protected

Definition at line 677 of file G4GammaGeneralProcess.cc.

679{
680 SelectedProcess(step, proc);
682 track.GetMaterial());
683 return proc->PostStepDoIt(track, step);
684}
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4Material * GetMaterial() const

Referenced by PostStepDoIt().

◆ SelectedProcess()

void G4GammaGeneralProcess::SelectedProcess ( const G4Step track,
const G4VProcess ptr 
)
inlineprotected

Definition at line 207 of file G4GammaGeneralProcess.hh.

208{
209 selectedProc = ptr;
210 step.GetPostStepPoint()->SetProcessDefinedStep(ptr);
211}

Referenced by PostStepDoIt(), SampleEmSecondaries(), and SampleHadSecondaries().

◆ StartTracking()

void G4GammaGeneralProcess::StartTracking ( G4Track )
overridevirtual

Reimplemented from G4VEmProcess.

Definition at line 476 of file G4GammaGeneralProcess.cc.

477{
479 currentMaterial = nullptr;
480}

◆ StorePhysicsTable()

G4bool G4GammaGeneralProcess::StorePhysicsTable ( const G4ParticleDefinition part,
const G4String directory,
G4bool  ascii = false 
)
overridevirtual

Reimplemented from G4VEmProcess.

Definition at line 688 of file G4GammaGeneralProcess.cc.

691{
692 G4bool yes = true;
693 if(!isTheMaster) { return yes; }
694 if(thePhotoElectric != nullptr &&
695 !thePhotoElectric->StorePhysicsTable(part, directory, ascii))
696 { yes = false; }
697 if(theCompton != nullptr &&
698 !theCompton->StorePhysicsTable(part, directory, ascii))
699 { yes = false; }
700 if(theConversionEE != nullptr &&
701 !theConversionEE->StorePhysicsTable(part, directory, ascii))
702 { yes = false; }
703 if(theRayleigh != nullptr &&
704 !theRayleigh->StorePhysicsTable(part, directory, ascii))
705 { yes = false; }
706
707 for(size_t i=0; i<nTables; ++i) {
708 if(theT[i]) {
709 G4String nam = (0==i || 2==i || 6==i || 10==i)
710 ? "LambdaGeneral" + nameT[i] : "ProbGeneral" + nameT[i];
711 G4String fnam = GetPhysicsTableFileName(part,directory,nam,ascii);
712 if(!theHandler->StorePhysicsTable(i, part, fnam, ascii)) { yes = false; }
713 }
714 }
715 return yes;
716}
G4bool StorePhysicsTable(size_t idx, const G4ParticleDefinition *part, const G4String &fname, G4bool ascii)
virtual G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override

Member Data Documentation

◆ selectedProc

const G4VProcess* G4GammaGeneralProcess::selectedProc
protected

◆ theGammaNuclear


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