Geant4 11.4.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 (const G4String &pname="GammaGeneralProc")
 ~G4GammaGeneralProcess () override
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 G4VProcessGetCreatorProcess () const override
const G4VProcessGetSelectedProcess () const
const G4StringGetSubProcessName () const
G4int GetSubProcessSubType () const
G4VEmProcessGetEmProcess (const G4String &name) override
G4HadronicProcessGetGammaNuclear () const
 G4GammaGeneralProcess (G4GammaGeneralProcess &)=delete
G4GammaGeneralProcessoperator= (const G4GammaGeneralProcess &right)=delete
Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 ~G4VEmProcess () override
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
G4double GetLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy)
G4double GetCrossSection (const G4double kinEnergy, const G4MaterialCutsCouple *couple) override
G4double ComputeCrossSectionPerAtom (G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
G4double MeanFreePath (const G4Track &track)
void SetLambdaBinning (G4int nbins)
void SetMinKinEnergy (G4double e)
void SetMinKinEnergyPrim (G4double e)
void SetMaxKinEnergy (G4double e)
G4PhysicsTableLambdaTable () const
G4PhysicsTableLambdaTablePrim () const
void SetLambdaTable (G4PhysicsTable *)
void SetLambdaTablePrim (G4PhysicsTable *)
std::vector< G4double > * EnergyOfCrossSectionMax () const
void SetEnergyOfCrossSectionMax (std::vector< G4double > *)
G4CrossSectionType CrossSectionType () const
void SetCrossSectionType (G4CrossSectionType val)
const G4ParticleDefinitionParticle () const
const G4ParticleDefinitionSecondaryParticle () const
G4VEmModelSelectModelForMaterial (G4double kinEnergy, std::size_t idxCouple) const
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel (G4VEmModel *, G4int index=0)
G4int NumberOfModels () const
G4VEmModelEmModel (std::size_t index=0) const
const G4VEmModelGetCurrentModel () const
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) 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 SetBuildTableFlag (G4bool val)
void CurrentSetup (const G4MaterialCutsCouple *, G4double energy)
G4bool UseBaseMaterial () const
void BuildLambdaTable ()
void StreamInfo (std::ostream &outFile, const G4ParticleDefinition &, G4bool rst=false) const
 G4VEmProcess (G4VEmProcess &)=delete
G4VEmProcessoperator= (const G4VEmProcess &right)=delete
Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 G4VDiscreteProcess (G4VDiscreteProcess &)
virtual ~G4VDiscreteProcess ()
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
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 &)
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
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
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)
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 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
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 (std::size_t idxe, std::size_t idxt)
G4double GetProbability (std::size_t idxt)
void SelectedProcess (const G4Step &step, G4VProcess *ptr)
void SelectEmProcess (const G4Step &, G4VEmProcess *)
void SelectHadProcess (const G4Track &, const G4Step &, G4HadronicProcess *)
G4double TotalCrossSectionPerVolume ()
Protected Member Functions inherited from G4VEmProcess
virtual void StreamProcessInfo (std::ostream &) const
G4VEmModelSelectModel (G4double kinEnergy, std::size_t)
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
G4ParticleChangeForGammaGetParticleChange ()
void SetParticle (const G4ParticleDefinition *p)
void SetSecondaryParticle (const G4ParticleDefinition *p)
std::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
G4int DensityIndex (G4int idx) const
G4double DensityFactor (G4int idx) const
Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
void ClearNumberOfInteractionLengthLeft ()

Protected Attributes

G4HadronicProcesstheGammaNuclear = nullptr
G4VProcessselectedProc = nullptr
G4double preStepLogE = 1.0
G4double factor = 1.0
Protected Attributes inherited from G4VEmProcess
const G4MaterialCutsCouplecurrentCouple = nullptr
const G4MaterialcurrentMaterial = nullptr
G4EmBiasingManagerbiasManager = nullptr
std::vector< G4double > * theEnergyOfCrossSectionMax = nullptr
G4double mfpKinEnergy = DBL_MAX
G4double preStepKinEnergy = 0.0
G4double preStepLambda = 0.0
G4int mainSecondaries = 1
G4int secID = _EM
G4int fluoID = _Fluorescence
G4int augerID = _AugerElectron
G4int biasID = _EM
G4int tripletID = _TripletElectron
std::size_t currentCoupleIndex = 0
std::size_t basedCoupleIndex = 0
std::size_t coupleIdxLambda = 0
std::size_t idxLambda = 0
G4bool isTheMaster = false
G4bool baseMat = false
std::vector< G4DynamicParticle * > secParticles
G4ParticleChangeForGamma fParticleChange
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 64 of file G4GammaGeneralProcess.hh.

Constructor & Destructor Documentation

◆ G4GammaGeneralProcess() [1/2]

G4GammaGeneralProcess::G4GammaGeneralProcess ( const G4String & pname = "GammaGeneralProc")
explicit

Definition at line 86 of file G4GammaGeneralProcess.cc.

86 :
88 minPEEnergy(150*CLHEP::keV),
89 minEEEnergy(2*CLHEP::electron_mass_c2),
90 minMMEnergy(100*CLHEP::MeV)
91{
95 if (nullptr == theHandler) {
96 theHandler = new G4EmDataHandler(nTables);
97 }
98}
@ fGammaGeneralProcess
@ fElectromagnetic
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
void SetParticle(const G4ParticleDefinition *p)
void SetVerboseLevel(G4int value)
void SetProcessSubType(G4int)

Referenced by G4GammaGeneralProcess(), and operator=().

◆ ~G4GammaGeneralProcess()

G4GammaGeneralProcess::~G4GammaGeneralProcess ( )
override

Definition at line 102 of file G4GammaGeneralProcess.cc.

103{}

◆ G4GammaGeneralProcess() [2/2]

G4GammaGeneralProcess::G4GammaGeneralProcess ( G4GammaGeneralProcess & )
delete

Member Function Documentation

◆ AddEmProcess()

void G4GammaGeneralProcess::AddEmProcess ( G4VEmProcess * ptr)

Definition at line 114 of file G4GammaGeneralProcess.cc.

115{
116 if(nullptr == ptr) { return; }
117 G4int stype = ptr->GetProcessSubType();
118 if(stype == fRayleigh) { theRayleigh = ptr; }
119 else if(stype == fPhotoElectricEffect) { thePhotoElectric = ptr; }
120 else if(stype == fComptonScattering) { theCompton = ptr; }
121 else if(stype == fGammaConversion) { theConversionEE = ptr; }
122}
@ fGammaConversion
@ fRayleigh
@ fComptonScattering
@ fPhotoElectricEffect
int G4int
Definition G4Types.hh:85
G4int GetProcessSubType() const

◆ AddHadProcess()

void G4GammaGeneralProcess::AddHadProcess ( G4HadronicProcess * ptr)

Definition at line 133 of file G4GammaGeneralProcess.cc.

134{
135 theGammaNuclear = ptr;
136}
G4HadronicProcess * theGammaNuclear

◆ AddMMProcess()

void G4GammaGeneralProcess::AddMMProcess ( G4GammaConversionToMuons * ptr)

Definition at line 126 of file G4GammaGeneralProcess.cc.

127{
128 theConversionMM = ptr;
129}

◆ BuildPhysicsTable()

void G4GammaGeneralProcess::BuildPhysicsTable ( const G4ParticleDefinition & part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 245 of file G4GammaGeneralProcess.cc.

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

◆ ComputeGeneralLambda()

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

Definition at line 195 of file G4GammaGeneralProcess.hh.

196{
197 idxEnergy = idxe;
198 return factor*theHandler->GetVector(idxt, basedCoupleIndex)
199 ->LogVectorValue(preStepKinEnergy, preStepLogE);
200}
std::size_t basedCoupleIndex
G4double preStepKinEnergy

Referenced by TotalCrossSectionPerVolume().

◆ GetCreatorProcess()

const G4VProcess * G4GammaGeneralProcess::GetCreatorProcess ( ) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 777 of file G4GammaGeneralProcess.cc.

778{
779 return selectedProc;
780}

◆ GetEmProcess()

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

Reimplemented from G4VEmProcess.

Definition at line 760 of file G4GammaGeneralProcess.cc.

761{
762 G4VEmProcess* proc = nullptr;
763 if(name == thePhotoElectric->GetProcessName()) {
764 proc = thePhotoElectric;
765 } else if(name == theCompton->GetProcessName()) {
766 proc = theCompton;
767 } else if(name == theConversionEE->GetProcessName()) {
768 proc = theConversionEE;
769 } else if(theRayleigh != nullptr && name == theRayleigh->GetProcessName()) {
770 proc = theRayleigh;
771 }
772 return proc;
773}

◆ GetGammaNuclear()

G4HadronicProcess * G4GammaGeneralProcess::GetGammaNuclear ( ) const
inline

Definition at line 238 of file G4GammaGeneralProcess.hh.

239{
240 return theGammaNuclear;
241}

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 722 of file G4GammaGeneralProcess.cc.

725{
727 return MeanFreePath(track);
728}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
G4double MeanFreePath(const G4Track &track)

◆ GetProbability()

G4double G4GammaGeneralProcess::GetProbability ( std::size_t idxt)
inlineprotected

Definition at line 204 of file G4GammaGeneralProcess.hh.

205{
206 return theHandler->GetVector(idxt, basedCoupleIndex)
207 ->LogVectorValue(preStepKinEnergy, preStepLogE);
208}

Referenced by PostStepDoIt().

◆ GetSelectedProcess()

const G4VProcess * G4GammaGeneralProcess::GetSelectedProcess ( ) const
inline

Definition at line 230 of file G4GammaGeneralProcess.hh.

231{
232 return selectedProc;
233}

◆ GetSubProcessName()

const G4String & G4GammaGeneralProcess::GetSubProcessName ( ) const

Definition at line 744 of file G4GammaGeneralProcess.cc.

745{
746 return (selectedProc) ? selectedProc->GetProcessName()
748}

◆ GetSubProcessSubType()

G4int G4GammaGeneralProcess::GetSubProcessSubType ( ) const

Definition at line 752 of file G4GammaGeneralProcess.cc.

753{
754 return (selectedProc) ? selectedProc->GetProcessSubType()
756}

◆ InitialiseProcess()

void G4GammaGeneralProcess::InitialiseProcess ( const G4ParticleDefinition * )
overrideprotectedvirtual

Implements G4VEmProcess.

Definition at line 186 of file G4GammaGeneralProcess.cc.

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

Referenced by PreparePhysicsTable().

◆ IsApplicable()

G4bool G4GammaGeneralProcess::IsApplicable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 107 of file G4GammaGeneralProcess.cc.

108{
109 return true;
110}

◆ operator=()

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

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 568 of file G4GammaGeneralProcess.cc.

570{
571 // In all cases clear number of interaction lengths
573 selectedProc = nullptr;
575 /*
576 G4cout << "PostStep: preStepLambda= " << preStepLambda
577 << " PE= " << peLambda << " q= " << q << " idxE= " << idxEnergy
578 << G4endl;
579 */
580 switch (idxEnergy) {
581 case 0:
582 if(preStepLambda*q <= peLambda) {
583 SelectEmProcess(step, thePhotoElectric);
584 } else {
585 if(theT[1] && preStepLambda*q < (preStepLambda - peLambda)*GetProbability(1) + peLambda) {
586 SelectEmProcess(step, theRayleigh);
587 } else {
588 SelectEmProcess(step, theCompton);
589 }
590 }
591 break;
592
593 case 1:
594 if(q <= GetProbability(3)) {
595 SelectEmProcess(step, thePhotoElectric);
596 } else if(q <= GetProbability(4)) {
597 SelectEmProcess(step, theCompton);
598 } else if(theRayleigh) {
599 SelectEmProcess(step, theRayleigh);
600 } else {
601 SelectEmProcess(step, thePhotoElectric);
602 }
603 break;
604
605 case 2:
606 if(q <= GetProbability(7)) {
607 SelectEmProcess(step, theConversionEE);
608 } else if(q <= GetProbability(8)) {
609 SelectEmProcess(step, theCompton);
610 } else if(q <= GetProbability(9)) {
611 SelectEmProcess(step, thePhotoElectric);
612 } else if(theGammaNuclear) {
613 SelectHadProcess(track, step, theGammaNuclear);
614 } else {
615 SelectEmProcess(step, theConversionEE);
616 }
617 break;
618
619 case 3:
620 if(q + GetProbability(11) <= 1.0) {
621 SelectEmProcess(step, theConversionEE);
622 } else if(q + GetProbability(12) <= 1.0) {
623 SelectEmProcess(step, theCompton);
624 } else if(q + GetProbability(13) <= 1.0) {
625 SelectEmProcess(step, thePhotoElectric);
626 } else if(theGammaNuclear && q + GetProbability(14) <= 1.0) {
627 SelectHadProcess(track, step, theGammaNuclear);
628 } else if(theConversionMM) {
629 SelectedProcess(step, theConversionMM);
630 } else {
631 SelectEmProcess(step, theConversionEE);
632 }
633 break;
634 }
635 // sample secondaries
636 if(selectedProc != nullptr) {
637 return selectedProc->PostStepDoIt(track, step);
638 }
639 // no interaction - exception case
640 fParticleChange.InitializeForPostStep(track);
641 return &fParticleChange;
642}
#define G4UniformRand()
Definition Randomize.hh:52
void SelectedProcess(const G4Step &step, G4VProcess *ptr)
void SelectHadProcess(const G4Track &, const G4Step &, G4HadronicProcess *)
void SelectEmProcess(const G4Step &, G4VEmProcess *)
G4double GetProbability(std::size_t idxt)
G4double preStepLambda
G4ParticleChangeForGamma fParticleChange
G4double theNumberOfInteractionLengthLeft

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 460 of file G4GammaGeneralProcess.cc.

464{
466 G4double x = DBL_MAX;
467
469 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
470
471 // compute mean free path
472 G4bool recompute = false;
473 if(couple != currentCouple) {
474 currentCouple = couple;
476 currentMaterial = couple->GetMaterial();
477 factor = 1.0;
478 if(baseMat) {
481 }
482 recompute = true;
483 }
484 if(energy != preStepKinEnergy) {
487 recompute = true;
488 }
489 if(recompute) {
491
492 // zero cross section
493 if(preStepLambda <= 0.0) {
496 }
497 }
498
499 // non-zero cross section
500 if(preStepLambda > 0.0) {
501
503
504 // beggining of tracking (or just after DoIt of this process)
507
508 } else if(currentInteractionLength < DBL_MAX) {
509
511 previousStepSize/currentInteractionLength;
514 }
515
516 // new mean free path and step limit for the next step
519 }
520 /*
521 G4cout << "PostStepGetPhysicalInteractionLength: e= " << energy
522 << " idxe= " << idxEnergy << " xs= " << preStepLambda
523 << " x= " << x << G4endl;
524 */
525 return x;
526}
bool G4bool
Definition G4Types.hh:86
G4double GetLogKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
std::size_t currentCoupleIndex
G4double DensityFactor(G4int idx) const
const G4MaterialCutsCouple * currentCouple
const G4Material * currentMaterial
G4double currentInteractionLength
G4double theInitialNumberOfInteractionLength
G4double energy(const ThreeVector &p, const G4double m)

◆ PreparePhysicsTable()

void G4GammaGeneralProcess::PreparePhysicsTable ( const G4ParticleDefinition & part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 140 of file G4GammaGeneralProcess.cc.

141{
142 SetParticle(&part);
143 preStepLambda = 0.0;
144 idxEnergy = 0;
145 currentCouple = nullptr;
146
147 G4LossTableManager* man = G4LossTableManager::Instance();
148
149 // initialise base material for the current run
150 G4LossTableBuilder* bld = man->GetTableBuilder();
153
154 if(1 < verboseLevel) {
155 G4cout << "G4GammaGeneralProcess::PreparePhysicsTable() for "
156 << GetProcessName()
157 << " and particle " << part.GetParticleName()
158 << " isMaster: " << isTheMaster << G4endl;
159 }
160
161 // 3 sub-processes must be always defined
162 if(thePhotoElectric == nullptr || theCompton == nullptr ||
163 theConversionEE == nullptr) {
165 ed << "### G4GeneralGammaProcess is initialized incorrectly"
166 << "\n Photoelectric: " << thePhotoElectric
167 << "\n Compton: " << theCompton
168 << "\n Conversion: " << theConversionEE;
169 G4Exception("G4GeneralGammaProcess","em0004",
170 FatalException, ed,"");
171 return;
172 }
173
174 thePhotoElectric->PreparePhysicsTable(part);
175 theCompton->PreparePhysicsTable(part);
176 theConversionEE->PreparePhysicsTable(part);
177 if (nullptr != theRayleigh) { theRayleigh->PreparePhysicsTable(part); }
178 if (nullptr != theGammaNuclear) { theGammaNuclear->PreparePhysicsTable(part); }
179 if (nullptr != theConversionMM) { theConversionMM->PreparePhysicsTable(part); }
180
181 InitialiseProcess(&part);
182}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
void InitialiseProcess(const G4ParticleDefinition *) override
static G4bool GetBaseMaterialFlag()
void InitialiseBaseMaterials(const G4PhysicsTable *table=nullptr)

◆ ProcessDescription()

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

Reimplemented from G4VProcess.

Definition at line 732 of file G4GammaGeneralProcess.cc.

733{
734 thePhotoElectric->ProcessDescription(out);
735 theCompton->ProcessDescription(out);
736 theConversionEE->ProcessDescription(out);
737 if(theRayleigh) { theRayleigh->ProcessDescription(out); }
738 if(theGammaNuclear) { theGammaNuclear->ProcessDescription(out); }
739 if(theConversionMM) { theConversionMM->ProcessDescription(out); }
740}

◆ RetrievePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 686 of file G4GammaGeneralProcess.cc.

689{
690 if (!isTheMaster) { return true; }
691 if(1 < verboseLevel) {
692 G4cout << "G4GammaGeneralProcess::RetrievePhysicsTable() for "
693 << part->GetParticleName() << " and process "
694 << GetProcessName() << G4endl;
695 }
696 G4bool yes = true;
697 if(!thePhotoElectric->RetrievePhysicsTable(part, directory, ascii))
698 { yes = false; }
699 if(!theCompton->RetrievePhysicsTable(part, directory, ascii))
700 { yes = false; }
701 if(!theConversionEE->RetrievePhysicsTable(part, directory, ascii))
702 { yes = false; }
703 if(theRayleigh != nullptr &&
704 !theRayleigh->RetrievePhysicsTable(part, directory, ascii))
705 { yes = false; }
706
707 for(std::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 G4bool spline = (i <= 1 || i >= 10);
713 if(!theHandler->RetrievePhysicsTable(i, part, fnam, ascii, spline))
714 { yes = false; }
715 }
716 }
717 return yes;
718}
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)

◆ SelectedProcess()

void G4GammaGeneralProcess::SelectedProcess ( const G4Step & step,
G4VProcess * ptr )
inlineprotected

Definition at line 213 of file G4GammaGeneralProcess.hh.

214{
215 selectedProc = ptr;
217}
void SetProcessDefinedStep(const G4VProcess *aValue)
G4StepPoint * GetPostStepPoint() const

Referenced by PostStepDoIt(), SelectEmProcess(), and SelectHadProcess().

◆ SelectEmProcess()

void G4GammaGeneralProcess::SelectEmProcess ( const G4Step & step,
G4VEmProcess * proc )
inlineprotected

Definition at line 222 of file G4GammaGeneralProcess.hh.

223{
225 SelectedProcess(step, proc);
226}
void CurrentSetup(const G4MaterialCutsCouple *, G4double energy)

Referenced by PostStepDoIt().

◆ SelectHadProcess()

void G4GammaGeneralProcess::SelectHadProcess ( const G4Track & track,
const G4Step & step,
G4HadronicProcess * proc )
protected

Definition at line 646 of file G4GammaGeneralProcess.cc.

648{
649 SelectedProcess(step, proc);
652}
G4CrossSectionDataStore * GetCrossSectionDataStore()

Referenced by PostStepDoIt().

◆ StartTracking()

void G4GammaGeneralProcess::StartTracking ( G4Track * )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 453 of file G4GammaGeneralProcess.cc.

454{
456}

◆ StorePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 656 of file G4GammaGeneralProcess.cc.

659{
660 G4bool yes = true;
661 if(!isTheMaster) { return yes; }
662 if(!thePhotoElectric->StorePhysicsTable(part, directory, ascii))
663 { yes = false; }
664 if(!theCompton->StorePhysicsTable(part, directory, ascii))
665 { yes = false; }
666 if(!theConversionEE->StorePhysicsTable(part, directory, ascii))
667 { yes = false; }
668 if(theRayleigh != nullptr &&
669 !theRayleigh->StorePhysicsTable(part, directory, ascii))
670 { yes = false; }
671
672 for(std::size_t i=0; i<nTables; ++i) {
673 if(theT[i]) {
674 G4String nam = (0==i || 2==i || 6==i || 10==i)
675 ? "LambdaGeneral" + nameT[i] : "ProbGeneral" + nameT[i];
676 G4String fnam = GetPhysicsTableFileName(part,directory,nam,ascii);
677 if(!theHandler->StorePhysicsTable(i, part, fnam, ascii)) { yes = false; }
678 }
679 }
680 return yes;
681}

◆ TotalCrossSectionPerVolume()

G4double G4GammaGeneralProcess::TotalCrossSectionPerVolume ( )
protected

Definition at line 530 of file G4GammaGeneralProcess.cc.

531{
532 G4double cross = 0.0;
533 /*
534 G4cout << "#Total: " << preStepKinEnergy << " " << minPEEnergy << " "
535 << minEEEnergy << " " << minMMEnergy<< G4endl;
536 G4cout << " idxE= " << idxEnergy
537 << " idxC= " << currentCoupleIndex << G4endl;
538 */
539 if(preStepKinEnergy < minPEEnergy) {
540 cross = ComputeGeneralLambda(0, 0);
541 //G4cout << "XS1: " << cross << G4endl;
542 peLambda = thePhotoElectric->GetLambda(preStepKinEnergy, currentCouple, preStepLogE);
543 cross += peLambda;
544 //G4cout << "XS2: " << peLambda << G4endl;
545
546 } else if(preStepKinEnergy < minEEEnergy) {
547 cross = ComputeGeneralLambda(1, 2);
548 //G4cout << "XS3: " << cross << G4endl;
549
550 } else if(preStepKinEnergy < minMMEnergy) {
551 cross = ComputeGeneralLambda(2, 6);
552 //G4cout << "XS4: " << cross << G4endl;
553
554 } else {
555 cross = ComputeGeneralLambda(3, 10);
556 //G4cout << "XS5: " << cross << G4endl;
557 }
558 /*
559 G4cout << "xs= " << cross << " idxE= " << idxEnergy
560 << " idxC= " << currentCoupleIndex
561 << " E= " << preStepKinEnergy << G4endl;
562 */
563 return cross;
564}
G4double ComputeGeneralLambda(std::size_t idxe, std::size_t idxt)

Referenced by PostStepGetPhysicalInteractionLength().

Member Data Documentation

◆ factor

G4double G4GammaGeneralProcess::factor = 1.0
protected

◆ preStepLogE

G4double G4GammaGeneralProcess::preStepLogE = 1.0
protected

◆ selectedProc

G4VProcess* G4GammaGeneralProcess::selectedProc = nullptr
protected

◆ theGammaNuclear

G4HadronicProcess* G4GammaGeneralProcess::theGammaNuclear = nullptr
protected

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