Geant4-11
Public Member Functions | Private Attributes
G4VAtomDeexcitation Class Referenceabstract

#include <G4VAtomDeexcitation.hh>

Inheritance diagram for G4VAtomDeexcitation:
G4UAtomicDeexcitation

Public Member Functions

void AlongStepDeexcitation (std::vector< G4Track * > &tracks, const G4Step &step, G4double &eLoss, G4int coupleIndex)
 
G4bool CheckAugerActiveRegion (G4int coupleIndex)
 
G4bool CheckDeexcitationActiveRegion (G4int coupleIndex)
 
virtual G4double ComputeShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr)=0
 
 G4VAtomDeexcitation (const G4String &modname="Deexcitation")
 
 G4VAtomDeexcitation (G4VAtomDeexcitation &)=delete
 
virtual void GenerateParticles (std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4double gammaCut, G4double eCut)=0
 
void GenerateParticles (std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
 
virtual const G4AtomicShellGetAtomicShell (G4int Z, G4AtomicShellEnumerator shell)=0
 
const std::vector< G4bool > & GetListOfActiveAtoms () const
 
const G4StringGetName () const
 
virtual G4double GetShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr)=0
 
G4int GetVerboseLevel () const
 
void InitialiseAtomicDeexcitation ()
 
virtual void InitialiseForExtraAtom (G4int Z)=0
 
virtual void InitialiseForNewRun ()=0
 
G4bool IsAugerActive () const
 
G4bool IsAugerCascadeActive () const
 
G4bool IsFluoActive () const
 
G4bool IsPIXEActive () const
 
G4VAtomDeexcitationoperator= (const G4VAtomDeexcitation &right)=delete
 
void SetAuger (G4bool)
 
void SetAugerCascade (G4bool)
 
void SetDeexcitationActiveRegion (const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
 
void SetFluo (G4bool)
 
void SetPIXE (G4bool)
 
void SetVerboseLevel (G4int)
 
virtual ~G4VAtomDeexcitation ()
 

Private Attributes

std::vector< G4boolactiveAugerMedia
 
std::vector< G4boolactiveDeexcitationMedia
 
std::vector< G4boolactivePIXEMedia
 
std::vector< G4StringactiveRegions
 
std::vector< G4boolactiveZ
 
std::vector< G4boolAugerRegions
 
std::vector< G4booldeRegions
 
G4bool flagAuger = false
 
G4bool flagPIXE = false
 
const G4ParticleDefinitiongamma
 
G4bool ignoreCuts = false
 
G4bool isActive = false
 
G4bool isActiveLocked = false
 
G4bool isAugerLocked = false
 
G4bool isPIXELocked = false
 
G4String name
 
G4int nCouples = 0
 
std::vector< G4boolPIXERegions
 
const G4ProductionCutsTabletheCoupleTable = nullptr
 
std::vector< G4DynamicParticle * > vdyn
 
G4int verbose = 1
 

Detailed Description

Definition at line 64 of file G4VAtomDeexcitation.hh.

Constructor & Destructor Documentation

◆ G4VAtomDeexcitation() [1/2]

G4VAtomDeexcitation::G4VAtomDeexcitation ( const G4String modname = "Deexcitation")
explicit

Definition at line 72 of file G4VAtomDeexcitation.cc.

73 : name(modname)
74{
75 vdyn.reserve(5);
76 theCoupleTable = nullptr;
78}
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
std::vector< G4DynamicParticle * > vdyn
const G4ParticleDefinition * gamma
const G4ProductionCutsTable * theCoupleTable

References G4Gamma::Gamma(), gamma, theCoupleTable, and vdyn.

◆ ~G4VAtomDeexcitation()

G4VAtomDeexcitation::~G4VAtomDeexcitation ( )
virtual

Definition at line 82 of file G4VAtomDeexcitation.cc.

83{}

◆ G4VAtomDeexcitation() [2/2]

G4VAtomDeexcitation::G4VAtomDeexcitation ( G4VAtomDeexcitation )
delete

Member Function Documentation

◆ AlongStepDeexcitation()

void G4VAtomDeexcitation::AlongStepDeexcitation ( std::vector< G4Track * > &  tracks,
const G4Step step,
G4double eLoss,
G4int  coupleIndex 
)

Definition at line 261 of file G4VAtomDeexcitation.cc.

265{
266 G4double truelength = step.GetStepLength();
267 if(!flagPIXE && !activePIXEMedia[coupleIndex]) { return; }
268 if(eLossMax <= 0.0 || truelength <= 0.0) { return; }
269
270 // step parameters
271 const G4StepPoint* preStep = step.GetPreStepPoint();
272 const G4ThreeVector prePos = preStep->GetPosition();
273 const G4ThreeVector delta = step.GetPostStepPoint()->GetPosition() - prePos;
274 const G4double preTime = preStep->GetGlobalTime();
275 const G4double dt = step.GetPostStepPoint()->GetGlobalTime() - preTime;
276
277 // particle parameters
278 const G4Track* track = step.GetTrack();
279 const G4ParticleDefinition* part = track->GetDefinition();
280 G4double ekin = preStep->GetKineticEnergy();
281
282 // media parameters
283 G4double gCut = (*theCoupleTable->GetEnergyCutsVector(0))[coupleIndex];
284 if(ignoreCuts) { gCut = 0.0; }
285 G4double eCut = DBL_MAX;
286 if(CheckAugerActiveRegion(coupleIndex)) {
287 eCut = (*theCoupleTable->GetEnergyCutsVector(1))[coupleIndex];
288 if(ignoreCuts) { eCut = 0.0; }
289 }
290
291 //G4cout<<"!Sample PIXE gCut(MeV)= "<<gCut<<" eCut(MeV)= "<<eCut
292 // <<" Ekin(MeV)= " << ekin/MeV << G4endl;
293
294 const G4Material* material = preStep->GetMaterial();
295 const G4ElementVector* theElementVector = material->GetElementVector();
296 const G4double* theAtomNumDensityVector =
297 material->GetVecNbOfAtomsPerVolume();
298 const G4int nelm = material->GetNumberOfElements();
299
300 // loop over deexcitations
301 for(G4int i=0; i<nelm; ++i) {
302 G4int Z = (*theElementVector)[i]->GetZasInt();
303 if(activeZ[Z] && Z < 93) {
304 G4int nshells =
305 std::min(9,(*theElementVector)[i]->GetNbOfAtomicShells());
306 G4double rho = truelength*theAtomNumDensityVector[i];
307 //G4cout<<" Z "<< Z <<" is active x(mm)= " << truelength/mm << G4endl;
308 for(G4int ii=0; ii<nshells; ++ii) {
310 const G4AtomicShell* shell = GetAtomicShell(Z, as);
311 const G4double bindingEnergy = shell->BindingEnergy();
312
313 if(gCut > bindingEnergy) { break; }
314
315 if(eLossMax > bindingEnergy) {
316 G4double sig = rho*
318
319 // mfp is mean free path in units of step size
320 if(sig > 0.0) {
321 G4double mfp = 1.0/sig;
322 G4double stot = 0.0;
323 //G4cout << " Shell " << ii << " mfp(mm)= " << mfp/mm << G4endl;
324 // sample ionisation points
325 do {
326 stot -= mfp*G4Log(G4UniformRand());
327 if( stot > 1.0 || eLossMax < bindingEnergy) { break; }
328 // sample deexcitation
329 vdyn.clear();
330 GenerateParticles(&vdyn, shell, Z, gCut, eCut);
331 G4int nsec = vdyn.size();
332 if(nsec > 0) {
333 G4ThreeVector r = prePos + stot*delta;
334 G4double time = preTime + stot*dt;
335 for(G4int j=0; j<nsec; ++j) {
336 G4DynamicParticle* dp = vdyn[j];
337 G4double e = dp->GetKineticEnergy();
338
339 // save new secondary if there is enough energy
340 if(eLossMax >= e) {
341 eLossMax -= e;
342 G4Track* t = new G4Track(dp, time, r);
343
344 // defined secondary type
345 if(dp->GetDefinition() == gamma) {
347 } else {
349 }
350 tracks.push_back(t);
351 } else {
352 delete dp;
353 }
354 }
355 }
356 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
357 } while (stot < 1.0);
358 }
359 }
360 }
361 }
362 }
363 return;
364}
G4AtomicShellEnumerator
std::vector< const G4Element * > G4ElementVector
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4UniformRand()
Definition: Randomize.hh:52
G4double BindingEnergy() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const std::vector< G4double > * GetEnergyCutsVector(std::size_t pcIdx) const
G4double GetGlobalTime() const
G4Material * GetMaterial() const
const G4ThreeVector & GetPosition() const
G4double GetKineticEnergy() const
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
G4ParticleDefinition * GetDefinition() const
void SetCreatorModelID(const G4int id)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4bool CheckAugerActiveRegion(G4int coupleIndex)
std::vector< G4bool > activePIXEMedia
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
std::vector< G4bool > activeZ
virtual G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr)=0
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double bindingEnergy(G4int A, G4int Z)
string material
Definition: eplot.py:19
#define DBL_MAX
Definition: templates.hh:62

References _ePIXE, _GammaPIXE, activePIXEMedia, activeZ, G4AtomicShell::BindingEnergy(), G4InuclSpecialFunctions::bindingEnergy(), CheckAugerActiveRegion(), DBL_MAX, flagPIXE, G4Log(), G4UniformRand, gamma, GenerateParticles(), GetAtomicShell(), G4DynamicParticle::GetDefinition(), G4Track::GetDefinition(), G4ProductionCutsTable::GetEnergyCutsVector(), G4StepPoint::GetGlobalTime(), G4DynamicParticle::GetKineticEnergy(), G4StepPoint::GetKineticEnergy(), G4StepPoint::GetMaterial(), G4StepPoint::GetPosition(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), GetShellIonisationCrossSectionPerAtom(), G4Step::GetStepLength(), G4Step::GetTrack(), ignoreCuts, eplot::material, G4INCL::Math::min(), G4Track::SetCreatorModelID(), theCoupleTable, vdyn, and Z.

Referenced by G4VEnergyLossProcess::AlongStepDoIt().

◆ CheckAugerActiveRegion()

G4bool G4VAtomDeexcitation::CheckAugerActiveRegion ( G4int  coupleIndex)
inline

Definition at line 266 of file G4VAtomDeexcitation.hh.

267{
268 return (idx < nCouples) ? activeAugerMedia[idx] : false;
269}
std::vector< G4bool > activeAugerMedia

References activeAugerMedia, and nCouples.

Referenced by AlongStepDeexcitation(), and GenerateParticles().

◆ CheckDeexcitationActiveRegion()

G4bool G4VAtomDeexcitation::CheckDeexcitationActiveRegion ( G4int  coupleIndex)
inline

◆ ComputeShellIonisationCrossSectionPerAtom()

virtual G4double G4VAtomDeexcitation::ComputeShellIonisationCrossSectionPerAtom ( const G4ParticleDefinition ,
G4int  Z,
G4AtomicShellEnumerator  shell,
G4double  kinE,
const G4Material mat = nullptr 
)
pure virtual

◆ GenerateParticles() [1/2]

virtual void G4VAtomDeexcitation::GenerateParticles ( std::vector< G4DynamicParticle * > *  secVect,
const G4AtomicShell ,
G4int  Z,
G4double  gammaCut,
G4double  eCut 
)
pure virtual

Implemented in G4UAtomicDeexcitation.

◆ GenerateParticles() [2/2]

void G4VAtomDeexcitation::GenerateParticles ( std::vector< G4DynamicParticle * > *  secVect,
const G4AtomicShell as,
G4int  Z,
G4int  coupleIndex 
)

Definition at line 235 of file G4VAtomDeexcitation.cc.

238{
239 G4double gCut = DBL_MAX;
240 if(ignoreCuts) {
241 gCut = 0.0;
242 } else if (nullptr != theCoupleTable) {
243 gCut = (*(theCoupleTable->GetEnergyCutsVector(0)))[idx];
244 }
245 if(gCut < as->BindingEnergy()) {
246 G4double eCut = DBL_MAX;
247 if(CheckAugerActiveRegion(idx)) {
248 if(ignoreCuts) {
249 eCut = 0.0;
250 } else if (nullptr != theCoupleTable) {
251 eCut = (*(theCoupleTable->GetEnergyCutsVector(1)))[idx];
252 }
253 }
254 GenerateParticles(v, as, Z, gCut, eCut);
255 }
256}

References CheckAugerActiveRegion(), DBL_MAX, GenerateParticles(), G4ProductionCutsTable::GetEnergyCutsVector(), ignoreCuts, theCoupleTable, and Z.

Referenced by AlongStepDeexcitation(), G4ECDecay::DecayIt(), G4ITDecay::DecayIt(), GenerateParticles(), G4DNABornIonisationModel1::SampleSecondaries(), G4DNABornIonisationModel2::SampleSecondaries(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNARelativisticIonisationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4MicroElecInelasticModel::SampleSecondaries(), G4MicroElecInelasticModel_new::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeIonisationModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), and G4PEEffectFluoModel::SampleSecondaries().

◆ GetAtomicShell()

virtual const G4AtomicShell * G4VAtomDeexcitation::GetAtomicShell ( G4int  Z,
G4AtomicShellEnumerator  shell 
)
pure virtual

◆ GetListOfActiveAtoms()

const std::vector< G4bool > & G4VAtomDeexcitation::GetListOfActiveAtoms ( ) const
inline

Definition at line 244 of file G4VAtomDeexcitation.hh.

245{
246 return activeZ;
247}

References activeZ.

◆ GetName()

const G4String & G4VAtomDeexcitation::GetName ( ) const
inline

Definition at line 238 of file G4VAtomDeexcitation.hh.

239{
240 return name;
241}

References name.

Referenced by SetDeexcitationActiveRegion().

◆ GetShellIonisationCrossSectionPerAtom()

virtual G4double G4VAtomDeexcitation::GetShellIonisationCrossSectionPerAtom ( const G4ParticleDefinition ,
G4int  Z,
G4AtomicShellEnumerator  shell,
G4double  kinE,
const G4Material mat = nullptr 
)
pure virtual

◆ GetVerboseLevel()

G4int G4VAtomDeexcitation::GetVerboseLevel ( ) const
inline

Definition at line 254 of file G4VAtomDeexcitation.hh.

255{
256 return verbose;
257}

References verbose.

◆ InitialiseAtomicDeexcitation()

void G4VAtomDeexcitation::InitialiseAtomicDeexcitation ( )

Definition at line 87 of file G4VAtomDeexcitation.cc.

88{
90 theParameters->DefineRegParamForDeex(this);
91
92 // Define list of couples
95
96 // needed for unit tests
97 size_t nn = std::max(nCouples, 1);
98 if(activeDeexcitationMedia.size() != nn) {
99 activeDeexcitationMedia.resize(nn, false);
100 activeAugerMedia.resize(nn, false);
101 activePIXEMedia.resize(nn, false);
102 }
103 if(activeZ.size() != 93) { activeZ.resize(93, false); }
104
105 // initialisation of flags and options
106 // normally there is no locksed flags
107 if(!isActiveLocked) { isActive = theParameters->Fluo(); }
108 if(!isAugerLocked) { flagAuger = theParameters->Auger(); }
109 if(!isPIXELocked) { flagPIXE = theParameters->Pixe(); }
110 ignoreCuts = theParameters->DeexcitationIgnoreCut();
111
112 // Define list of regions
113 size_t nRegions = deRegions.size();
114 // check if deexcitation is active for the given run
115 if(!isActive && 0 == nRegions) { return; }
116
117 // if no active regions add a world
118 if(0 == nRegions) {
120 nRegions = deRegions.size();
121 }
122
123 if(0 < verbose) {
124 G4cout << G4endl;
125 G4cout << "### === Deexcitation model " << name
126 << " is activated for " << nRegions;
127 if(1 == nRegions) { G4cout << " region:" << G4endl; }
128 else { G4cout << " regions:" << G4endl;}
129 }
130
131 // Identify active media
132 const G4RegionStore* regionStore = G4RegionStore::GetInstance();
133 for(size_t j=0; j<nRegions; ++j) {
134 const G4Region* reg = regionStore->GetRegion(activeRegions[j], false);
135 if(nullptr != reg && 0 < nCouples) {
136 const G4ProductionCuts* rpcuts = reg->GetProductionCuts();
137 if(0 < verbose) {
138 G4cout << " " << activeRegions[j]
139 << " " << deRegions[j] << " " << AugerRegions[j]
140 << " " << PIXERegions[j] << G4endl;
141 }
142 for(G4int i=0; i<nCouples; ++i) {
143 const G4MaterialCutsCouple* couple =
145 if (couple->GetProductionCuts() == rpcuts) {
149 }
150 }
151 }
152 }
154 //G4cout << nelm << G4endl;
155 for(G4int k=0; k<nelm; ++k) {
156 G4int Z = (*(G4Element::GetElementTable()))[k]->GetZasInt();
157 if(Z > 5 && Z < 93) {
158 activeZ[Z] = true;
159 //G4cout << "!!! Active de-excitation Z= " << Z << G4endl;
160 }
161 }
162
163 // Initialise derived class
165
166 if(0 < verbose && flagAuger) {
167 G4cout << "### === Auger flag: " << flagAuger
168 << G4endl;
169 }
170 if(0 < verbose) {
171 G4cout << "### === Ignore cuts flag: " << ignoreCuts
172 << G4endl;
173 }
174 if(0 < verbose && flagPIXE) {
175 G4cout << "### === PIXE model for hadrons: "
176 << theParameters->PIXECrossSectionModel()
177 << G4endl;
178 G4cout << "### === PIXE model for e+-: "
179 << theParameters->PIXEElectronCrossSectionModel()
180 << G4endl;
181 }
182}
static const G4double reg
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:397
static size_t GetNumberOfElements()
Definition: G4Element.cc:404
static G4EmParameters * Instance()
const G4String & PIXECrossSectionModel()
const G4String & PIXEElectronCrossSectionModel()
G4bool Fluo() const
G4bool Pixe() const
void DefineRegParamForDeex(G4VAtomDeexcitation *) const
G4bool DeexcitationIgnoreCut() const
G4bool Auger() const
G4ProductionCuts * GetProductionCuts() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
std::vector< G4String > activeRegions
std::vector< G4bool > AugerRegions
std::vector< G4bool > deRegions
std::vector< G4bool > PIXERegions
void SetDeexcitationActiveRegion(const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
virtual void InitialiseForNewRun()=0
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References activeAugerMedia, activeDeexcitationMedia, activePIXEMedia, activeRegions, activeZ, G4EmParameters::Auger(), AugerRegions, G4EmParameters::DeexcitationIgnoreCut(), G4EmParameters::DefineRegParamForDeex(), deRegions, flagAuger, flagPIXE, G4EmParameters::Fluo(), G4cout, G4endl, G4Element::GetElementTable(), G4RegionStore::GetInstance(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4Element::GetNumberOfElements(), G4MaterialCutsCouple::GetProductionCuts(), G4ProductionCutsTable::GetProductionCutsTable(), G4RegionStore::GetRegion(), G4ProductionCutsTable::GetTableSize(), ignoreCuts, InitialiseForNewRun(), G4EmParameters::Instance(), isActive, isActiveLocked, isAugerLocked, isPIXELocked, G4INCL::Math::max(), name, nCouples, G4InuclParticleNames::nn, G4EmParameters::Pixe(), G4EmParameters::PIXECrossSectionModel(), G4EmParameters::PIXEElectronCrossSectionModel(), PIXERegions, reg, SetDeexcitationActiveRegion(), theCoupleTable, verbose, and Z.

Referenced by LBE::ConstructGeneral(), and G4LossTableManager::ResetParameters().

◆ InitialiseForExtraAtom()

virtual void G4VAtomDeexcitation::InitialiseForExtraAtom ( G4int  Z)
pure virtual

Implemented in G4UAtomicDeexcitation.

◆ InitialiseForNewRun()

virtual void G4VAtomDeexcitation::InitialiseForNewRun ( )
pure virtual

Implemented in G4UAtomicDeexcitation.

Referenced by InitialiseAtomicDeexcitation().

◆ IsAugerActive()

G4bool G4VAtomDeexcitation::IsAugerActive ( ) const
inline

Definition at line 213 of file G4VAtomDeexcitation.hh.

214{
215 return flagAuger;
216}

References flagAuger.

Referenced by G4UAtomicDeexcitation::GenerateAuger().

◆ IsAugerCascadeActive()

G4bool G4VAtomDeexcitation::IsAugerCascadeActive ( ) const
inline

◆ IsFluoActive()

G4bool G4VAtomDeexcitation::IsFluoActive ( ) const
inline

◆ IsPIXEActive()

G4bool G4VAtomDeexcitation::IsPIXEActive ( ) const
inline

◆ operator=()

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

◆ SetAuger()

void G4VAtomDeexcitation::SetAuger ( G4bool  val)
inline

Definition at line 208 of file G4VAtomDeexcitation.hh.

209{
210 if(!isAugerLocked) { flagAuger = val; isAugerLocked = true; }
211}

References flagAuger, and isAugerLocked.

Referenced by SetAugerCascade().

◆ SetAugerCascade()

void G4VAtomDeexcitation::SetAugerCascade ( G4bool  val)
inline

Definition at line 218 of file G4VAtomDeexcitation.hh.

219{
220 SetAuger(val);
221}

References SetAuger().

◆ SetDeexcitationActiveRegion()

void G4VAtomDeexcitation::SetDeexcitationActiveRegion ( const G4String rname,
G4bool  valDeexcitation,
G4bool  valAuger,
G4bool  valPIXE 
)

Definition at line 187 of file G4VAtomDeexcitation.cc.

191{
192 // no PIXE in parallel world
193 if(rname == "DefaultRegionForParallelWorld") { return; }
194
195 G4String ss = rname;
196 /*
197 G4cout << "### G4VAtomDeexcitation::SetDeexcitationActiveRegion " << ss
198 << " " << valDeexcitation << " " << valAuger
199 << " " << valPIXE << G4endl;
200 */
201 if(ss == "world" || ss == "World" || ss == "WORLD") {
202 ss = "DefaultRegionForTheWorld";
203 }
204 size_t n = deRegions.size();
205 for(size_t i=0; i<n; ++i) {
206
207 // Region already exist
208 if(ss == activeRegions[i]) {
209 deRegions[i] = valDeexcitation;
210 AugerRegions[i] = valAuger;
211 PIXERegions[i] = valPIXE;
212 return;
213 }
214 }
215 // New region
216 activeRegions.push_back(ss);
217 deRegions.push_back(valDeexcitation);
218 AugerRegions.push_back(valAuger);
219 PIXERegions.push_back(valPIXE);
220
221 // if de-excitation defined for the world volume
222 // it should be active for all G4Regions
223 if(ss == "DefaultRegionForTheWorld") {
225 G4int nn = regions->size();
226 for(G4int i=0; i<nn; ++i) {
227 if(ss == (*regions)[i]->GetName()) { continue; }
228 SetDeexcitationActiveRegion((*regions)[i]->GetName(), valDeexcitation,
229 valAuger, valPIXE);
230
231 }
232 }
233}
const G4String & GetName() const

References activeRegions, AugerRegions, deRegions, G4RegionStore::GetInstance(), GetName(), CLHEP::detail::n, G4InuclParticleNames::nn, PIXERegions, and SetDeexcitationActiveRegion().

Referenced by G4EmLowEParameters::DefineRegParamForDeex(), InitialiseAtomicDeexcitation(), and SetDeexcitationActiveRegion().

◆ SetFluo()

void G4VAtomDeexcitation::SetFluo ( G4bool  val)
inline

Definition at line 198 of file G4VAtomDeexcitation.hh.

199{
200 if(!isActiveLocked) { isActive = val; isActiveLocked = true; }
201}

References isActive, and isActiveLocked.

Referenced by G4EmDNAPhysics_stationary_option2::ConstructProcess().

◆ SetPIXE()

void G4VAtomDeexcitation::SetPIXE ( G4bool  val)
inline

Definition at line 228 of file G4VAtomDeexcitation.hh.

229{
230 if(!isPIXELocked) { flagPIXE = val; isPIXELocked = true; }
231}

References flagPIXE, and isPIXELocked.

◆ SetVerboseLevel()

void G4VAtomDeexcitation::SetVerboseLevel ( G4int  val)
inline

Definition at line 249 of file G4VAtomDeexcitation.hh.

250{
251 verbose = val;
252}

References verbose.

Referenced by G4LossTableManager::ResetParameters().

Field Documentation

◆ activeAugerMedia

std::vector<G4bool> G4VAtomDeexcitation::activeAugerMedia
private

Definition at line 183 of file G4VAtomDeexcitation.hh.

Referenced by CheckAugerActiveRegion(), and InitialiseAtomicDeexcitation().

◆ activeDeexcitationMedia

std::vector<G4bool> G4VAtomDeexcitation::activeDeexcitationMedia
private

◆ activePIXEMedia

std::vector<G4bool> G4VAtomDeexcitation::activePIXEMedia
private

Definition at line 184 of file G4VAtomDeexcitation.hh.

Referenced by AlongStepDeexcitation(), and InitialiseAtomicDeexcitation().

◆ activeRegions

std::vector<G4String> G4VAtomDeexcitation::activeRegions
private

◆ activeZ

std::vector<G4bool> G4VAtomDeexcitation::activeZ
private

◆ AugerRegions

std::vector<G4bool> G4VAtomDeexcitation::AugerRegions
private

◆ deRegions

std::vector<G4bool> G4VAtomDeexcitation::deRegions
private

◆ flagAuger

G4bool G4VAtomDeexcitation::flagAuger = false
private

◆ flagPIXE

G4bool G4VAtomDeexcitation::flagPIXE = false
private

◆ gamma

const G4ParticleDefinition* G4VAtomDeexcitation::gamma
private

Definition at line 166 of file G4VAtomDeexcitation.hh.

Referenced by AlongStepDeexcitation(), and G4VAtomDeexcitation().

◆ ignoreCuts

G4bool G4VAtomDeexcitation::ignoreCuts = false
private

◆ isActive

G4bool G4VAtomDeexcitation::isActive = false
private

Definition at line 172 of file G4VAtomDeexcitation.hh.

Referenced by InitialiseAtomicDeexcitation(), IsFluoActive(), and SetFluo().

◆ isActiveLocked

G4bool G4VAtomDeexcitation::isActiveLocked = false
private

Definition at line 177 of file G4VAtomDeexcitation.hh.

Referenced by InitialiseAtomicDeexcitation(), and SetFluo().

◆ isAugerLocked

G4bool G4VAtomDeexcitation::isAugerLocked = false
private

Definition at line 178 of file G4VAtomDeexcitation.hh.

Referenced by InitialiseAtomicDeexcitation(), and SetAuger().

◆ isPIXELocked

G4bool G4VAtomDeexcitation::isPIXELocked = false
private

Definition at line 179 of file G4VAtomDeexcitation.hh.

Referenced by InitialiseAtomicDeexcitation(), and SetPIXE().

◆ name

G4String G4VAtomDeexcitation::name
private

◆ nCouples

G4int G4VAtomDeexcitation::nCouples = 0
private

◆ PIXERegions

std::vector<G4bool> G4VAtomDeexcitation::PIXERegions
private

◆ theCoupleTable

const G4ProductionCutsTable* G4VAtomDeexcitation::theCoupleTable = nullptr
private

◆ vdyn

std::vector<G4DynamicParticle*> G4VAtomDeexcitation::vdyn
private

Definition at line 188 of file G4VAtomDeexcitation.hh.

Referenced by AlongStepDeexcitation(), and G4VAtomDeexcitation().

◆ verbose

G4int G4VAtomDeexcitation::verbose = 1
private

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