Geant4-11
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4MuonicAtomDecay Class Reference

#include <G4MuonicAtomDecay.hh>

Inheritance diagram for G4MuonicAtomDecay:
G4VRestDiscreteProcess G4VProcess

Public Member Functions

virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &theTrack, const G4Step &theStep)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void DumpInfo () const
 
virtual void EndTracking ()
 
 G4MuonicAtomDecay (G4HadronicInteraction *hiptr=nullptr, const G4String &processName="MuonicAtomDecay")
 
G4double GetCurrentInteractionLength () const
 
const G4VProcessGetMasterProcess () const
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)
 
G4double GetNumberOfInteractionLengthLeft () const
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
G4double GetPILfactor () const
 
virtual const G4ProcessManagerGetProcessManager ()
 
const G4StringGetProcessName () const
 
G4int GetProcessSubType () const
 
G4ProcessType GetProcessType () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4int GetVerboseLevel () const
 
G4bool isAlongStepDoItIsEnabled () const
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
G4bool operator!= (const G4VProcess &right) const
 
G4bool operator== (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &theTrack, const G4Step &theStep)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 
virtual void ProcessDescription (std::ostream &outFile) const
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
void SetPILfactor (G4double value)
 
virtual void SetProcessManager (const G4ProcessManager *)
 
void SetProcessSubType (G4int)
 
void SetProcessType (G4ProcessType)
 
void SetVerboseLevel (G4int value)
 
virtual void StartTracking (G4Track *)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual ~G4MuonicAtomDecay ()
 

Static Public Member Functions

static const G4StringGetProcessTypeName (G4ProcessType)
 

Protected Member Functions

void ClearNumberOfInteractionLengthLeft ()
 
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 

Protected Attributes

G4ParticleChange aParticleChange
 
const G4ProcessManageraProcessManager = nullptr
 
G4double currentInteractionLength = -1.0
 
G4bool enableAlongStepDoIt = true
 
G4bool enableAtRestDoIt = true
 
G4bool enablePostStepDoIt = true
 
G4VParticleChangepParticleChange = nullptr
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4String thePhysicsTableFileName
 
G4double thePILfactor = 1.0
 
G4String theProcessName
 
G4int theProcessSubType = -1
 
G4ProcessType theProcessType = fNotDefined
 

Private Member Functions

G4VParticleChangeDecayIt (const G4Track &theTrack, const G4Step &theStep)
 
void DumpState (const G4Track &, const G4String &, G4ExceptionDescription &)
 
void FillResult (G4HadFinalState *aR, const G4Track &aT)
 
 G4MuonicAtomDecay (const G4MuonicAtomDecay &right)
 
G4MuonicAtomDecayoperator= (const G4MuonicAtomDecay &right)
 

Private Attributes

G4HadronicInteractioncmptr
 
const G4double fMuMass
 
G4ProcessTablefProcessTable = nullptr
 
G4VProcessmasterProcessShadow = nullptr
 
G4ParticleChange theTotalResult
 
G4int verboseLevel
 

Detailed Description

Definition at line 50 of file G4MuonicAtomDecay.hh.

Constructor & Destructor Documentation

◆ G4MuonicAtomDecay() [1/2]

G4MuonicAtomDecay::G4MuonicAtomDecay ( G4HadronicInteraction hiptr = nullptr,
const G4String processName = "MuonicAtomDecay" 
)
explicit

Definition at line 67 of file G4MuonicAtomDecay.cc.

70 fMuMass(G4MuonMinus::MuonMinus()->GetPDGMass()),
71 cmptr(hiptr),
73{
74 // This is not a hadronic process; assume it is a kind of decay
75 enableAtRestDoIt = true;
76 enablePostStepDoIt = true; // it is a streach; fixme
77 theProcessSubType = 221; // (see G4DecayProcessType.hh) fixme
78 if (!cmptr) {
79 // cmptr = new G4CascadeInterface(); // Bertini - Pointer owned by InteractionRegistry
80 cmptr = new G4MuMinusCapturePrecompound(); // Precompound
81 }
82}
@ fDecay
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:99
const G4double fMuMass
G4HadronicInteraction * cmptr
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:359
G4bool enablePostStepDoIt
Definition: G4VProcess.hh:361
G4int theProcessSubType
Definition: G4VProcess.hh:349
const char * name(G4int ptype)

References cmptr, G4VProcess::enableAtRestDoIt, G4VProcess::enablePostStepDoIt, and G4VProcess::theProcessSubType.

◆ ~G4MuonicAtomDecay()

G4MuonicAtomDecay::~G4MuonicAtomDecay ( )
virtual

Definition at line 86 of file G4MuonicAtomDecay.cc.

87{}

◆ G4MuonicAtomDecay() [2/2]

G4MuonicAtomDecay::G4MuonicAtomDecay ( const G4MuonicAtomDecay right)
private

Member Function Documentation

◆ AlongStepDoIt()

virtual G4VParticleChange * G4VRestDiscreteProcess::AlongStepDoIt ( const G4Track ,
const G4Step  
)
inlinevirtualinherited

Implements G4VProcess.

Definition at line 87 of file G4VRestDiscreteProcess.hh.

90 { return 0; }

◆ AlongStepGetPhysicalInteractionLength()

virtual G4double G4VRestDiscreteProcess::AlongStepGetPhysicalInteractionLength ( const G4Track ,
G4double  ,
G4double  ,
G4double ,
G4GPILSelection  
)
inlinevirtualinherited

Implements G4VProcess.

Definition at line 78 of file G4VRestDiscreteProcess.hh.

84 { return -1.0; }

◆ AlongStepGPIL()

G4double G4VProcess::AlongStepGPIL ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double proposedSafety,
G4GPILSelection selection 
)
inlineinherited

Definition at line 461 of file G4VProcess.hh.

466{
467 return AlongStepGetPhysicalInteractionLength(track, previousStepSize,
468 currentMinimumStep, proposedSafety, selection);
469}
virtual G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0

References G4VProcess::AlongStepGetPhysicalInteractionLength().

Referenced by G4SteppingManager::DefinePhysicalStepLength(), and G4ITStepProcessor::DoDefinePhysicalStepLength().

◆ AtRestDoIt()

virtual G4VParticleChange * G4MuonicAtomDecay::AtRestDoIt ( const G4Track theTrack,
const G4Step theStep 
)
inlinevirtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 86 of file G4MuonicAtomDecay.hh.

88 {return DecayIt(theTrack, theStep);}
G4VParticleChange * DecayIt(const G4Track &theTrack, const G4Step &theStep)

References DecayIt().

◆ AtRestGetPhysicalInteractionLength()

G4double G4MuonicAtomDecay::AtRestGetPhysicalInteractionLength ( const G4Track aTrack,
G4ForceCondition condition 
)
virtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 113 of file G4MuonicAtomDecay.cc.

115{
117 // check if this is the beginning of tracking
120 }
122}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
virtual G4double GetMeanLifeTime(const G4Track &aTrack, G4ForceCondition *condition)
virtual void ResetNumberOfInteractionLengthLeft()
Definition: G4VProcess.cc:80
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:331

References condition(), GetMeanLifeTime(), NotForced, G4VProcess::ResetNumberOfInteractionLengthLeft(), and G4VProcess::theNumberOfInteractionLengthLeft.

◆ AtRestGPIL()

G4double G4VProcess::AtRestGPIL ( const G4Track track,
G4ForceCondition condition 
)
inlineinherited

Definition at line 472 of file G4VProcess.hh.

474{
476}
G4double thePILfactor
Definition: G4VProcess.hh:352
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)=0

References G4VProcess::AtRestGetPhysicalInteractionLength(), condition(), and G4VProcess::thePILfactor.

Referenced by G4ITStepProcessor::GetAtRestIL(), and G4SteppingManager::InvokeAtRestDoItProcs().

◆ BuildPhysicsTable()

virtual void G4VProcess::BuildPhysicsTable ( const G4ParticleDefinition )
inlinevirtualinherited

◆ BuildWorkerPhysicsTable()

void G4VProcess::BuildWorkerPhysicsTable ( const G4ParticleDefinition part)
virtualinherited

Reimplemented in G4BiasingProcessInterface.

Definition at line 200 of file G4VProcess.cc.

201{
202 BuildPhysicsTable(part);
203}
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:187

References G4VProcess::BuildPhysicsTable().

Referenced by G4BiasingProcessInterface::BuildWorkerPhysicsTable().

◆ ClearNumberOfInteractionLengthLeft()

void G4VProcess::ClearNumberOfInteractionLengthLeft ( )
inlineprotectedinherited

◆ DecayIt()

G4VParticleChange * G4MuonicAtomDecay::DecayIt ( const G4Track theTrack,
const G4Step theStep 
)
private

Definition at line 151 of file G4MuonicAtomDecay.cc.

153{
154
155 // mainly based on G4HadronStoppingProcess & G4Decay
156 // if primary is not Alive then do nothing
157 theTotalResult.Clear(); // G4ParticleChange*
159 theTotalResult.ProposeWeight(aTrack.GetWeight());
160 if(aTrack.GetTrackStatus() != fAlive &&
161 aTrack.GetTrackStatus() != fStopButAlive) {
162 return &theTotalResult;
163 }
164
165 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
166 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
167 G4MuonicAtom const* muatom = static_cast<G4MuonicAtom const*>(aParticleDef);
168 G4Ions const* baseion = muatom->GetBaseIon();
169 G4int Z = baseion->GetAtomicNumber();
170 G4double Zd = Z;
171 G4double KEnergy = G4MuonicAtomHelper::GetKShellEnergy(Zd); // fixme check
172 G4HadProjectile thePro(aTrack); // G4HadProjectile, here the muonic atom
173 thePro.SetBoundEnergy(KEnergy);
174
175 G4ForceCondition* condition = nullptr; // it is unused in the following call anyway
176 G4double meanlife = GetMeanLifeTime(aTrack, condition);
177
178 G4HadFinalState* result = nullptr; // result before converting to G4VParticleChange*
179 // G4int nSecondaries = 0;
180 // save track time and start from zero time
181 // G4double time0 = aTrack.GetGlobalTime(); FillResult does it
182 // see G4Decay DecayIt
183 // see time0 down below
184 thePro.SetGlobalTime(0.0);
185
186 // do we need G4double fRemainderLifeTime; ???
187
188 G4double maDTime = theNumberOfInteractionLengthLeft*meanlife; //fixme check
189#ifdef G4VERBOSE
190 if (GetVerboseLevel()>1) {
191 G4cout << "G4MuonicAtomDecay::DecayIt time set to: "<< maDTime/ns << "[ns]" << G4endl;
192 }
193#endif
194
195 // decide on DIO or Capture
196
197 G4double lambdac = 1./muatom->GetDIOLifeTime();
198 G4double lambdad = 1./muatom->GetNCLifeTime();
199 G4double lambda = lambdac + lambdad;
200
201 if ( G4UniformRand()*lambda < lambdac) {
202 // if ( false ) { // force NC for testing
203
204 // DIO
205 // result = dmptr->ApplyYourself(thePro, *nucleus); // not quite the reaction;
206 // using G4PhaseSpaceDecayChannel
207
208#ifdef G4VERBOSE
209 if (GetVerboseLevel()>0) {
210 G4cout << "G4MuonicAtomDecay::DecayIt: selected DIO mode" << G4endl;
211 }
212#endif
213
214 // decay table; we use it only for the DIO which is more of a decay
215 // code mostly copied from G4Decay
216
217 G4DecayProducts* products = nullptr;
218 G4DecayTable *decaytable = aParticleDef->GetDecayTable();
219 G4VDecayChannel* decaychannel = nullptr;
220 G4double massParent = aParticle->GetMass();
221 decaychannel = decaytable->SelectADecayChannel(massParent);
222 if (decaychannel == nullptr) {
223 // decay channel not found
225 ed << "Can not determine decay channel for "
226 << aParticleDef->GetParticleName() << G4endl
227 << " mass of dynamic particle: " << massParent/GeV << " (GEV)" << G4endl
228 << " dacay table has " << decaytable->entries() << " entries" << G4endl;
229 G4double checkedmass=massParent;
230 if (massParent < 0.) {
231 checkedmass=aParticleDef->GetPDGMass();
232 ed << "Using PDG mass ("<<checkedmass/GeV << "(GeV)) in IsOKWithParentMass" << G4endl;
233 }
234 for (G4int ic =0;ic <decaytable->entries();++ic) {
235 G4VDecayChannel * dc= decaytable->GetDecayChannel(ic);
236 ed << ic << ": BR " << dc->GetBR() << ", IsOK? "
237 << dc->IsOKWithParentMass(checkedmass)
238 << ", --> ";
239 G4int ndaughters=dc->GetNumberOfDaughters();
240 for (G4int id=0;id<ndaughters;++id) {
241 if (id>0) ed << " + "; // seperator, except for first
242 ed << dc->GetDaughterName(id);
243 }
244 ed << G4endl;
245 }
246 G4Exception("G4MuonicAtomDecay::DecayIt", "DECAY003", FatalException,ed);
247 return &theTotalResult;
248 } else {
249 // execute DecayIt()
250#ifdef G4VERBOSE
251 G4int temp = decaychannel->GetVerboseLevel();
252 if (GetVerboseLevel()>1) {
253 G4cout << "G4MuonicAtomDecay::DecayIt : selected decay channel addr:"
254 << decaychannel <<G4endl;
255 decaychannel->SetVerboseLevel(GetVerboseLevel());
256 }
257#endif
258 products = decaychannel->DecayIt(aParticle->GetMass());
259 if(products == nullptr) {
261 ed << "No products are generated for "
262 << aParticleDef->GetParticleName();
263 G4Exception("G4MuonicAtomDecay::DecayIt","DECAY003",FatalException,ed);
264 return &theTotalResult;
265 }
266#ifdef G4VERBOSE
267 if (GetVerboseLevel()>1) {
268 decaychannel->SetVerboseLevel(temp);
269 }
270 if (GetVerboseLevel()>2) {
271 if (! products->IsChecked() ) products->DumpInfo();
272 }
273#endif
274 }
275
276 // get parent particle information ...................................
277 G4double ParentEnergy = aParticle->GetTotalEnergy();
278 G4double ParentMass = aParticle->GetMass();
279 if (ParentEnergy < ParentMass) {
280 if (GetVerboseLevel()>0) {
281 G4cout << "G4MuonicAtomDecay::DecayIt : Total Energy is less than its mass" << G4endl;
282 G4cout << " Particle: " << aParticle->GetDefinition()->GetParticleName();
283 G4cout << " Energy:" << ParentEnergy/MeV << "[MeV]";
284 G4cout << " Mass:" << ParentMass/MeV << "[MeV]";
285 G4cout << G4endl;
286 }
287 G4Exception( "G4MuonicAtomDecay::DecayIt ",
288 "DECAY102",JustWarning,
289 "Total Energy is less than its mass");
290 ParentEnergy = ParentMass;
291 }
292
293 G4ThreeVector ParentDirection(aParticle->GetMomentumDirection());
294
295 //boost all decay products to laboratory frame
296 G4double energyDeposit = 0.0;
297 G4double finalGlobalTime = aTrack.GetGlobalTime();
298 G4double finalLocalTime = aTrack.GetLocalTime();
299 if (aTrack.GetTrackStatus() == fStopButAlive ){
300 // AtRest case
301 finalGlobalTime += maDTime;
302 finalLocalTime += maDTime;
303 energyDeposit += aParticle->GetKineticEnergy();
304 } else {
305 // PostStep case
306 products->Boost( ParentEnergy, ParentDirection);
307 }
308
309 //add products in theTotalResult
310 G4int numberOfSecondaries = products->entries();
311 theTotalResult.SetNumberOfSecondaries(numberOfSecondaries);
312
313#ifdef G4VERBOSE
314 if (GetVerboseLevel()>1) {
315 G4cout << "G4MuonicAtomDecay::DecayIt : Decay vertex :";
316 G4cout << " Time: " << finalGlobalTime/ns << "[ns]";
317 G4cout << " X:" << (aTrack.GetPosition()).x() /cm << "[cm]";
318 G4cout << " Y:" << (aTrack.GetPosition()).y() /cm << "[cm]";
319 G4cout << " Z:" << (aTrack.GetPosition()).z() /cm << "[cm]";
320 G4cout << G4endl;
321 G4cout << "G4MuonicAtomDecay::DecayIt : decay products in Lab. Frame" << G4endl;
322 products->DumpInfo();
323 }
324#endif
325 G4int index;
326 G4ThreeVector currentPosition;
327 const G4TouchableHandle thand = aTrack.GetTouchableHandle();
328 for (index=0; index < numberOfSecondaries; index++)
329 {
330 // get current position of the track
331 currentPosition = aTrack.GetPosition();
332 // create a new track object
333 G4Track* secondary = new G4Track( products->PopProducts(),
334 finalGlobalTime ,
335 currentPosition );
336 // switch on good for tracking flag
337 secondary->SetGoodForTrackingFlag();
338 secondary->SetTouchableHandle(thand);
339 // add the secondary track in the List
340 theTotalResult.AddSecondary(secondary);
341 }
342 delete products;
343
344 // Kill the parent particle
347 theTotalResult.ProposeLocalTime( finalLocalTime );
348
349 // Clear NumberOfInteractionLengthLeft
351
352 return &theTotalResult ;
353
354 } else { //either or
355
356 // nuclearCapture
357
358 // model
359 // need to be able to choose between preco or bertini; no good way to do it?
360 // hardcoded in the constructor for now
361
362#ifdef G4VERBOSE
363 if (GetVerboseLevel()>0) {
364 G4cout << "G4MuonicAtomDecay::DecayIt: selected NC mode" << G4endl;
365 }
366#endif
367
368 G4int A = baseion->GetAtomicMass();
369 // G4Nucleus* nucleus = GetTargetNucleusPointer(); // from G4HadronicProcess
370 G4Nucleus nucleus;
371 nucleus.SetParameters(A, Z);
372
373 // we define a local projectile here which will be the orbiting muon
374 // we shall assume it is at rest; fixme
375
376 // G4HadProjectile, here the muon
378 G4ThreeVector(0.,0.,0.)));
379 theMuPro.SetBoundEnergy(KEnergy);
380 theMuPro.SetGlobalTime(0.0);
381
382 G4int reentryCount = 0; // this may be in the model already; check fixme <---
383 do {
384 // sample final state
385 // nuclear interaction should keep G4HadFinalState object
386 // model should define time of each secondary particle
387 try {
388 result = cmptr->ApplyYourself(theMuPro, nucleus); // muon and muonic atom nucleus
389 ++reentryCount;
390 }
391 catch(G4HadronicException & aR) {
393 ed << "Call for " << cmptr->GetModelName() << G4endl;
394 ed << " Z= "
395 << nucleus.GetZ_asInt()
396 << " A= " << nucleus.GetA_asInt() << G4endl;
397 DumpState(aTrack,"ApplyYourself",ed);
398 ed << " ApplyYourself failed" << G4endl;
399 G4Exception("G4MuonicAtomDecay::DecayIt", "HAD_MAD_101",
400 FatalException, ed);
401 }
402
403 // Check the result for catastrophic energy non-conservation
404 // result = CheckResult(theMuPro, nucleus, result);
405
406 if(reentryCount>100) {
408 ed << "Call for " << cmptr->GetModelName() << G4endl;
409 ed << " Z= "
410 << nucleus.GetZ_asInt()
411 << " A= " << nucleus.GetA_asInt() << G4endl;
412 DumpState(aTrack,"ApplyYourself",ed);
413 ed << " ApplyYourself does not completed after 100 attempts" << G4endl;
414 G4Exception("G4MuonicAtomDecay::DecayIt", "HAD_MAD_102",
415 FatalException, ed);
416 }
417 // Loop checking, 06-Aug-2015, Vladimir Ivanchenko
418 } while(result == nullptr);
419
420 // add delay time of capture (inter + intra)
421 G4int nsec = result->GetNumberOfSecondaries();
422 for(G4int i=0; i<nsec; ++i) {
423 G4HadSecondary* sec = result->GetSecondary(i);
424 G4double ctime = sec->GetTime();
425 sec->SetTime(maDTime + ctime); // we add time0 in the next stage
426#ifdef G4VERBOSE
427 if (GetVerboseLevel()>1) {
428 G4cout << "G4MuonicAtomDecay::DecayIt time set to: "
429 << (maDTime + ctime)/ns << "[ns]" << G4endl;
430 }
431#endif
432 }
433
434 FillResult(result,aTrack);
435
437 return &theTotalResult;
438 }
439}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
static constexpr double cm
Definition: G4SIunits.hh:99
CLHEP::Hep3Vector G4ThreeVector
@ fAlive
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
void DumpInfo() const
G4int entries() const
G4DynamicParticle * PopProducts()
G4bool IsChecked() const
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
G4VDecayChannel * GetDecayChannel(G4int index) const
G4VDecayChannel * SelectADecayChannel(G4double parentMass=-1.)
Definition: G4DecayTable.cc:82
G4int entries() const
G4double GetMass() const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
std::size_t GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
void SetTime(G4double aT)
G4double GetTime() const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
const G4String & GetModelName() const
Definition: G4Ions.hh:52
G4int GetVerboseLevel() const
void DumpState(const G4Track &, const G4String &, G4ExceptionDescription &)
G4ParticleChange theTotalResult
void FillResult(G4HadFinalState *aR, const G4Track &aT)
static G4double GetKShellEnergy(G4double A)
G4double GetNCLifeTime() const
G4Ions const * GetBaseIon() const
Definition: G4MuonicAtom.hh:96
G4double GetDIOLifeTime() const
G4int GetA_asInt() const
Definition: G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:105
void SetParameters(const G4double A, const G4double Z, const G4int numberOfLambdas=0)
Definition: G4Nucleus.cc:307
void AddSecondary(G4Track *aSecondary)
void ProposeLocalTime(G4double t)
virtual void Initialize(const G4Track &)
G4int GetAtomicNumber() const
G4int GetAtomicMass() const
G4DecayTable * GetDecayTable() const
const G4String & GetParticleName() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
void SetGoodForTrackingFlag(G4bool value=true)
G4double GetBR() const
void SetVerboseLevel(G4int value)
G4int GetVerboseLevel() const
G4int GetNumberOfDaughters() const
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
virtual G4bool IsOKWithParentMass(G4double parentMass)
const G4String & GetDaughterName(G4int anIndex) const
void ProposeTrackStatus(G4TrackStatus status)
void ProposeWeight(G4double finalWeight)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:424
#define ns
Definition: xmlparse.cc:614

References A, G4ParticleChange::AddSecondary(), G4HadronicInteraction::ApplyYourself(), G4DecayProducts::Boost(), G4VParticleChange::Clear(), G4VProcess::ClearNumberOfInteractionLengthLeft(), cm, cmptr, condition(), G4VDecayChannel::DecayIt(), G4DecayProducts::DumpInfo(), DumpState(), G4DecayProducts::entries(), G4DecayTable::entries(), fAlive, FatalException, FillResult(), fStopAndKill, fStopButAlive, G4cout, G4endl, G4Exception(), G4UniformRand, G4Nucleus::GetA_asInt(), G4ParticleDefinition::GetAtomicMass(), G4ParticleDefinition::GetAtomicNumber(), G4MuonicAtom::GetBaseIon(), G4VDecayChannel::GetBR(), G4VDecayChannel::GetDaughterName(), G4DecayTable::GetDecayChannel(), G4ParticleDefinition::GetDecayTable(), G4DynamicParticle::GetDefinition(), G4MuonicAtom::GetDIOLifeTime(), G4Track::GetDynamicParticle(), G4Track::GetGlobalTime(), G4DynamicParticle::GetKineticEnergy(), G4MuonicAtomHelper::GetKShellEnergy(), G4Track::GetLocalTime(), G4DynamicParticle::GetMass(), GetMeanLifeTime(), G4HadronicInteraction::GetModelName(), G4DynamicParticle::GetMomentumDirection(), G4MuonicAtom::GetNCLifeTime(), G4VDecayChannel::GetNumberOfDaughters(), G4HadFinalState::GetNumberOfSecondaries(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGMass(), G4Track::GetPosition(), G4HadFinalState::GetSecondary(), G4HadSecondary::GetTime(), G4DynamicParticle::GetTotalEnergy(), G4Track::GetTouchableHandle(), G4Track::GetTrackStatus(), G4VDecayChannel::GetVerboseLevel(), GetVerboseLevel(), G4Track::GetWeight(), G4Nucleus::GetZ_asInt(), GeV, G4ParticleChange::Initialize(), G4DecayProducts::IsChecked(), G4VDecayChannel::IsOKWithParentMass(), JustWarning, G4InuclParticleNames::lambda, MeV, G4MuonMinus::MuonMinus(), ns, G4DecayProducts::PopProducts(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChange::ProposeLocalTime(), G4VParticleChange::ProposeTrackStatus(), G4VParticleChange::ProposeWeight(), G4DecayTable::SelectADecayChannel(), G4HadProjectile::SetBoundEnergy(), G4HadProjectile::SetGlobalTime(), G4Track::SetGoodForTrackingFlag(), G4VParticleChange::SetNumberOfSecondaries(), G4Nucleus::SetParameters(), G4HadSecondary::SetTime(), G4Track::SetTouchableHandle(), G4VDecayChannel::SetVerboseLevel(), G4VProcess::theNumberOfInteractionLengthLeft, theTotalResult, and Z.

Referenced by AtRestDoIt(), and PostStepDoIt().

◆ DumpInfo()

void G4VProcess::DumpInfo ( ) const
virtualinherited

Reimplemented in G4AdjointAlongStepWeightCorrection, G4AdjointForcedInteractionForGamma, G4AdjointhMultipleScattering, G4ContinuousGainOfEnergy, G4eAdjointMultipleScattering, G4eInverseBremsstrahlung, G4eInverseCompton, G4eInverseIonisation, G4InversePEEffect, G4IonInverseIonisation, G4PolarizedAnnihilation, G4PolarizedBremsstrahlung, G4PolarizedCompton, G4PolarizedGammaConversion, G4PolarizedIonisation, G4PolarizedPhotoElectric, G4Cerenkov, G4ForwardXrayTR, G4GammaXTRadiator, G4GaussXTRadiator, G4RegularXTRadiator, G4Scintillation, G4StrawTubeXTRadiator, G4SynchrotronRadiation, G4TransitionRadiation, G4TransparentRegXTRadiator, G4VTransitionRadiation, G4VXTRenergyLoss, G4XTRGammaRadModel, G4XTRRegularRadModel, and G4XTRTransparentRegRadModel.

Definition at line 167 of file G4VProcess.cc.

168{
169 G4cout << "Process Name " << theProcessName ;
170 G4cout << " : Type[" << GetProcessTypeName(theProcessType) << "]";
171 G4cout << " : SubType[" << theProcessSubType << "]"<< G4endl;
172}
static const G4String & GetProcessTypeName(G4ProcessType)
Definition: G4VProcess.cc:134
G4ProcessType theProcessType
Definition: G4VProcess.hh:346
G4String theProcessName
Definition: G4VProcess.hh:341

References G4cout, G4endl, G4VProcess::GetProcessTypeName(), G4VProcess::theProcessName, G4VProcess::theProcessSubType, and G4VProcess::theProcessType.

Referenced by G4ProcessTable::DumpInfo(), export_G4VProcess(), G4Scintillation::ProcessDescription(), G4Cerenkov::ProcessDescription(), and G4ProcessManagerMessenger::SetNewValue().

◆ DumpState()

void G4MuonicAtomDecay::DumpState ( const G4Track aTrack,
const G4String method,
G4ExceptionDescription ed 
)
private

Definition at line 558 of file G4MuonicAtomDecay.cc.

561{
562 ed << "Unrecoverable error in the method " << method << " of "
563 << GetProcessName() << G4endl;
564 ed << "TrackID= "<< aTrack.GetTrackID() << " ParentID= "
565 << aTrack.GetParentID()
566 << " " << aTrack.GetParticleDefinition()->GetParticleName()
567 << G4endl;
568 ed << "Ekin(GeV)= " << aTrack.GetKineticEnergy()/CLHEP::GeV
569 << "; direction= " << aTrack.GetMomentumDirection() << G4endl;
570 ed << "Position(mm)= " << aTrack.GetPosition()/CLHEP::mm << ";";
571
572 if (aTrack.GetMaterial()) {
573 ed << " material " << aTrack.GetMaterial()->GetName();
574 }
575 ed << G4endl;
576
577 if (aTrack.GetVolume()) {
578 ed << "PhysicalVolume <" << aTrack.GetVolume()->GetName()
579 << ">" << G4endl;
580 }
581}
const G4String & GetName() const
Definition: G4Material.hh:173
G4int GetTrackID() const
const G4ParticleDefinition * GetParticleDefinition() const
G4VPhysicalVolume * GetVolume() const
const G4ThreeVector & GetPosition() const
G4Material * GetMaterial() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4int GetParentID() const
const G4String & GetName() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
static constexpr double mm
Definition: SystemOfUnits.h:96
static constexpr double GeV

References G4endl, G4Track::GetKineticEnergy(), G4Track::GetMaterial(), G4Track::GetMomentumDirection(), G4VPhysicalVolume::GetName(), G4Material::GetName(), G4Track::GetParentID(), G4Track::GetParticleDefinition(), G4ParticleDefinition::GetParticleName(), G4Track::GetPosition(), G4VProcess::GetProcessName(), G4Track::GetTrackID(), G4Track::GetVolume(), CLHEP::GeV, and CLHEP::mm.

Referenced by DecayIt(), and FillResult().

◆ EndTracking()

void G4VProcess::EndTracking ( )
virtualinherited

◆ FillResult()

void G4MuonicAtomDecay::FillResult ( G4HadFinalState aR,
const G4Track aT 
)
private

Definition at line 450 of file G4MuonicAtomDecay.cc.

451{
452 // based on G4HadronicProcess::FillResult
453
455
457 G4ThreeVector it(0., 0., 1.);
458
459 G4double efinal = aR->GetEnergyChange();
460 if(efinal < 0.0) { efinal = 0.0; }
461
462 // check status of primary
463 if(aR->GetStatusChange() == stopAndKill) {
466
467 // check its final energy
468 } else if(0.0 == efinal) {
471 ->GetAtRestProcessVector()->size() > 0)
473 else { theTotalResult.ProposeTrackStatus(fStopAndKill); } // check fixme
474
475 // primary is not killed apply rotation and Lorentz transformation
476 } else {
479 G4double newE = efinal + mass;
480 G4double newP = std::sqrt(efinal*(efinal + 2*mass));
481 G4ThreeVector newPV = newP*aR->GetMomentumChange();
482 G4LorentzVector newP4(newE, newPV);
483 newP4.rotate(rotation, it);
484 newP4 *= aR->GetTrafoToLab();
485 theTotalResult.ProposeMomentumDirection(newP4.vect().unit());
486 newE = newP4.e() - mass;
487#ifdef G4VERBOSE
488 if (GetVerboseLevel()>1 && newE <= 0.0) {
490 DumpState(aT,"Primary has zero energy after interaction",ed);
491 G4Exception("G4MuonicAtomDecay::FillResults", "HAD_MAD_103", JustWarning, ed);
492 }
493#endif
494 if(newE < 0.0) { newE = 0.0; }
496 }
497 //G4cout << "FillResult: Efinal= " << efinal << " status= "
498 // << theTotalResult.GetTrackStatus()
499 // << " fKill= " << fStopAndKill << G4endl;
500
501 // check secondaries: apply rotation and Lorentz transformation
502 G4int nSec = aR->GetNumberOfSecondaries();
504 G4double weight = aT.GetWeight();
505
506 if (nSec > 0) {
507 G4double time0 = aT.GetGlobalTime();
508 for (G4int i = 0; i < nSec; ++i) {
510 theM.rotate(rotation, it);
511 theM *= aR->GetTrafoToLab();
512 aR->GetSecondary(i)->GetParticle()->Set4Momentum(theM);
513
514 // time of interaction starts from zero
515 G4double time = aR->GetSecondary(i)->GetTime();
516 if (time < 0.0) { time = 0.0; }
517
518 // take into account global time
519 time += time0;
520
521 G4Track* track = new G4Track(aR->GetSecondary(i)->GetParticle(),
522 time, aT.GetPosition());
524 G4double newWeight = weight*aR->GetSecondary(i)->GetWeight();
525 // G4cout << "#### ParticleDebug "
526 // <<GetProcessName()<<" "
527 //<<aR->GetSecondary(i)->GetParticle()->GetDefinition()->GetParticleName()<<" "
528 // <<aScaleFactor<<" "
529 // <<XBiasSurvivalProbability()<<" "
530 // <<XBiasSecondaryWeight()<<" "
531 // <<aT.GetWeight()<<" "
532 // <<aR->GetSecondary(i)->GetWeight()<<" "
533 // <<aR->GetSecondary(i)->GetParticle()->Get4Momentum()<<" "
534 // <<G4endl;
535 track->SetWeight(newWeight);
538#ifdef G4VERBOSE
539 if (GetVerboseLevel()>1) {
540 G4double e = track->GetKineticEnergy();
541 if (e <= 0.0) {
543 DumpState(aT,"Secondary has zero energy",ed);
544 ed << "Secondary " << track->GetDefinition()->GetParticleName()
545 << G4endl;
546 G4Exception("G4MuonicAtomDecay::FillResults", "HAD_MAD_103",
547 JustWarning,ed);
548 }
549 }
550#endif
551 }
552 }
553 aR->Clear();
554}
@ stopAndKill
HepLorentzVector & rotate(double, const Hep3Vector &)
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
G4double GetEnergyChange() const
const G4LorentzRotation & GetTrafoToLab() const
G4HadFinalStateStatus GetStatusChange() const
G4double GetLocalEnergyDeposit() const
const G4ThreeVector & GetMomentumChange() const
G4DynamicParticle * GetParticle()
G4double GetWeight() const
G4int GetCreatorModelID() const
void ProposeEnergy(G4double finalEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4ProcessManager * GetProcessManager() const
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
G4double GetWeight() const
void SetWeight(G4double aValue)
G4double GetGlobalTime() const
G4ParticleDefinition * GetDefinition() const
const G4TouchableHandle & GetTouchableHandle() const
void SetCreatorModelID(const G4int id)
static constexpr double twopi
Definition: SystemOfUnits.h:56

References G4ParticleChange::AddSecondary(), G4HadFinalState::Clear(), DumpState(), CLHEP::HepLorentzVector::e(), fAlive, fStopAndKill, fStopButAlive, G4endl, G4Exception(), G4UniformRand, G4DynamicParticle::Get4Momentum(), G4ProcessManager::GetAtRestProcessVector(), G4HadSecondary::GetCreatorModelID(), G4Track::GetDefinition(), G4HadFinalState::GetEnergyChange(), G4Track::GetGlobalTime(), G4Track::GetKineticEnergy(), G4HadFinalState::GetLocalEnergyDeposit(), G4HadFinalState::GetMomentumChange(), G4HadFinalState::GetNumberOfSecondaries(), G4HadSecondary::GetParticle(), G4Track::GetParticleDefinition(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGMass(), G4Track::GetPosition(), G4ParticleDefinition::GetProcessManager(), G4HadFinalState::GetSecondary(), G4HadFinalState::GetStatusChange(), G4HadSecondary::GetTime(), G4Track::GetTouchableHandle(), G4HadFinalState::GetTrafoToLab(), GetVerboseLevel(), G4HadSecondary::GetWeight(), G4Track::GetWeight(), JustWarning, G4ParticleChange::ProposeEnergy(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChange::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), CLHEP::HepLorentzVector::rotate(), G4DynamicParticle::Set4Momentum(), G4Track::SetCreatorModelID(), G4VParticleChange::SetNumberOfSecondaries(), G4Track::SetTouchableHandle(), G4Track::SetWeight(), G4ProcessVector::size(), stopAndKill, theTotalResult, CLHEP::twopi, CLHEP::Hep3Vector::unit(), and CLHEP::HepLorentzVector::vect().

Referenced by DecayIt().

◆ GetCurrentInteractionLength()

G4double G4VProcess::GetCurrentInteractionLength ( ) const
inlineinherited

◆ GetMasterProcess()

const G4VProcess * G4VProcess::GetMasterProcess ( ) const
inlineinherited

◆ GetMeanFreePath()

G4double G4MuonicAtomDecay::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Implements G4VRestDiscreteProcess.

Definition at line 585 of file G4MuonicAtomDecay.cc.

586{
587 // based on G4Decay::GetMeanFreePath check; fixme
588
589 // get particle
590 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
591 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
592 G4double aMass = aParticle->GetMass();
593 G4double aLife = aParticleDef->GetPDGLifeTime();
594
595 // returns the mean free path in GEANT4 internal units
596 G4double pathlength;
597 G4double aCtau = c_light * aLife;
598
599 // check if the particle is stable?
600 if (aParticleDef->GetPDGStable()) {
601 pathlength = DBL_MAX;
602
603 //check if the particle has very short life time ?
604 } else if (aCtau < DBL_MIN) {
605 pathlength = DBL_MIN;
606
607 } else {
608 //calculate the mean free path
609 // by using normalized kinetic energy (= Ekin/mass)
610 G4double rKineticEnergy = aParticle->GetKineticEnergy()/aMass;
611 const G4double HighestValue = 20.0; //
612 if ( rKineticEnergy > HighestValue) {
613 // gamma >> 1
614 pathlength = ( rKineticEnergy + 1.0)* aCtau;
615 } else if ( rKineticEnergy < DBL_MIN ) {
616 // too slow particle
617#ifdef G4VERBOSE
618 if (GetVerboseLevel()>1) {
619 G4cout << "G4MuonicAtomDecay::GetMeanFreePath() !!particle stops!!";
620 G4cout << aParticleDef->GetParticleName() << G4endl;
621 G4cout << "KineticEnergy:" << aParticle->GetKineticEnergy()/GeV <<"[GeV]";
622 }
623#endif
624 pathlength = DBL_MIN;
625 } else {
626 // beta <1
627 pathlength = (aParticle->GetTotalMomentum())/aMass*aCtau ;
628 }
629 }
630 return pathlength;
631}
G4double GetTotalMomentum() const
G4bool GetPDGStable() const
G4double GetPDGLifeTime() const
const G4DynamicParticle * GetDynamicParticle() const
float c_light
Definition: hepunit.py:256
#define DBL_MIN
Definition: templates.hh:54
#define DBL_MAX
Definition: templates.hh:62

References source.hepunit::c_light, DBL_MAX, DBL_MIN, G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMass(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGLifeTime(), G4ParticleDefinition::GetPDGStable(), G4DynamicParticle::GetTotalMomentum(), GetVerboseLevel(), and GeV.

◆ GetMeanLifeTime()

G4double G4MuonicAtomDecay::GetMeanLifeTime ( const G4Track aTrack,
G4ForceCondition condition 
)
virtual

Implements G4VRestDiscreteProcess.

Definition at line 135 of file G4MuonicAtomDecay.cc.

137{
138 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
139 G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
140 G4MuonicAtom* muatom = static_cast<G4MuonicAtom*>(aParticleDef);
141 G4double meanlife = muatom->GetPDGLifeTime();
142#ifdef G4VERBOSE
143 if (GetVerboseLevel()>1) {
144 G4cout << "mean life time: "<< meanlife/ns << "[ns]" << G4endl;
145 }
146#endif
147 return meanlife;
148}

References G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4ParticleDefinition::GetPDGLifeTime(), GetVerboseLevel(), and ns.

Referenced by AtRestGetPhysicalInteractionLength(), and DecayIt().

◆ GetNumberOfInteractionLengthLeft()

G4double G4VProcess::GetNumberOfInteractionLengthLeft ( ) const
inlineinherited

Definition at line 431 of file G4VProcess.hh.

432{
434}

References G4VProcess::theNumberOfInteractionLengthLeft.

◆ GetPhysicsTableFileName()

const G4String & G4VProcess::GetPhysicsTableFileName ( const G4ParticleDefinition particle,
const G4String directory,
const G4String tableName,
G4bool  ascii = false 
)
inherited

◆ GetPILfactor()

G4double G4VProcess::GetPILfactor ( ) const
inlineinherited

Definition at line 455 of file G4VProcess.hh.

456{
457 return thePILfactor;
458}

References G4VProcess::thePILfactor.

Referenced by export_G4VProcess().

◆ GetProcessManager()

const G4ProcessManager * G4VProcess::GetProcessManager ( )
inlinevirtualinherited

Reimplemented in G4BiasingProcessInterface, and G4WrapperProcess.

Definition at line 494 of file G4VProcess.hh.

495{
496 return aProcessManager;
497}
const G4ProcessManager * aProcessManager
Definition: G4VProcess.hh:319

References G4VProcess::aProcessManager.

Referenced by G4BiasingProcessInterface::GetProcessManager(), and G4WrapperProcess::GetProcessManager().

◆ GetProcessName()

const G4String & G4VProcess::GetProcessName ( ) const
inlineinherited

Definition at line 382 of file G4VProcess.hh.

383{
384 return theProcessName;
385}

References G4VProcess::theProcessName.

Referenced by G4VEnergyLossProcess::ActivateForcedInteraction(), G4VEmProcess::ActivateForcedInteraction(), G4ProcessManager::ActivateProcess(), G4VEmProcess::ActivateSecondaryBiasing(), G4VEnergyLossProcess::ActivateSecondaryBiasing(), G4ParallelGeometriesLimiterProcess::AddParallelWorld(), G4IonQMDPhysics::AddProcess(), G4IonINCLXXPhysics::AddProcess(), G4ProcessManager::AddProcess(), G4ProcessPlacer::AddProcessAs(), G4ITSteppingVerbose::AlongStepDoItAllDone(), G4SteppingVerbose::AlongStepDoItAllDone(), G4SteppingVerboseWithUnits::AlongStepDoItAllDone(), G4ITSteppingVerbose::AlongStepDoItOneByOne(), G4SteppingVerbose::AlongStepDoItOneByOne(), G4SteppingVerboseWithUnits::AlongStepDoItOneByOne(), G4VContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4ParallelWorldProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AlongStepGetPhysicalInteractionLength(), G4VRestContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4VContinuousProcess::AlongStepGetPhysicalInteractionLength(), G4BOptnLeadingParticle::ApplyFinalStateBiasing(), G4ITSteppingVerbose::AtRestDoItInvoked(), G4SteppingVerbose::AtRestDoItInvoked(), G4SteppingVerboseWithUnits::AtRestDoItInvoked(), G4ITSteppingVerbose::AtRestDoItOneByOne(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VITRestProcess::AtRestGetPhysicalInteractionLength(), G4VRestProcess::AtRestGetPhysicalInteractionLength(), G4HadronicProcess::BiasCrossSectionByFactor(), G4VXTRenergyLoss::BuildAngleForEnergyBank(), G4VEnergyLossProcess::BuildDEDXTable(), G4VUserPhysicsList::BuildIntegralPhysicsTable(), G4VEmProcess::BuildLambdaTable(), G4VEnergyLossProcess::BuildLambdaTable(), G4DNABrownianTransportation::BuildPhysicsTable(), G4GammaGeneralProcess::BuildPhysicsTable(), G4VEmProcess::BuildPhysicsTable(), G4VEnergyLossProcess::BuildPhysicsTable(), G4VMultipleScattering::BuildPhysicsTable(), G4LossTableManager::BuildPhysicsTable(), G4LossTableManager::BuildTables(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4ProcessManager::CheckOrderingParameters(), G4HadronicProcess::CheckResult(), G4StackChecker::ClassifyNewTrack(), G4BOptrForceCollision::ConfigureForWorker(), G4RunManagerKernel::ConfirmCoupledTransportation(), G4FastSimulationPhysics::ConstructProcess(), G4GenericBiasingPhysics::ConstructProcess(), G4IonElasticPhysics::ConstructProcess(), G4LossTableManager::CopyTables(), G4RichTrajectory::CreateAttValues(), G4RichTrajectoryPoint::CreateAttValues(), G4VPhononProcess::CreateSecondary(), G4EmExtraParameters::DefineRegParamForEM(), G4EmExtraParameters::DefineRegParamForLoss(), G4HadronicProcessStore::DeRegisterExtraProcess(), G4ITSteppingVerbose::DPSLAlongStep(), G4SteppingVerbose::DPSLAlongStep(), G4SteppingVerboseWithUnits::DPSLAlongStep(), G4ITSteppingVerbose::DPSLPostStep(), G4SteppingVerbose::DPSLPostStep(), G4SteppingVerboseWithUnits::DPSLPostStep(), G4HadronicProcessStore::Dump(), G4HadronicProcess::DumpState(), DumpState(), G4ExceptionHandler::DumpTrackInfo(), export_G4VProcess(), G4EmCalculator::FindEmModel(), G4VEmProcess::FindLambdaMax(), G4BiasingProcessInterface::G4BiasingProcessInterface(), G4Cerenkov::G4Cerenkov(), G4ErrorEnergyLoss::G4ErrorEnergyLoss(), G4ErrorTrackLengthTarget::G4ErrorTrackLengthTarget(), G4FastSimulationManagerProcess::G4FastSimulationManagerProcess(), G4ImportanceProcess::G4ImportanceProcess(), G4MaxTimeCuts::G4MaxTimeCuts(), G4MicroElecSurface::G4MicroElecSurface(), G4MinEkineCuts::G4MinEkineCuts(), G4OpAbsorption::G4OpAbsorption(), G4OpBoundaryProcess::G4OpBoundaryProcess(), G4OpMieHG::G4OpMieHG(), G4OpRayleigh::G4OpRayleigh(), G4OpWLS::G4OpWLS(), G4OpWLS2::G4OpWLS2(), G4ParallelWorldProcess::G4ParallelWorldProcess(), G4ParallelWorldScoringProcess::G4ParallelWorldScoringProcess(), G4Scintillation::G4Scintillation(), G4ScoreSplittingProcess::G4ScoreSplittingProcess(), G4SpecialCuts::G4SpecialCuts(), G4StepLimiter::G4StepLimiter(), G4UCNAbsorption::G4UCNAbsorption(), G4UCNBoundaryProcess::G4UCNBoundaryProcess(), G4UCNLoss::G4UCNLoss(), G4UCNMultiScattering::G4UCNMultiScattering(), G4UserSpecialCuts::G4UserSpecialCuts(), G4WeightCutOffProcess::G4WeightCutOffProcess(), G4WeightWindowProcess::G4WeightWindowProcess(), G4HadronicProcess::GetElementCrossSection(), G4VEmProcess::GetEmProcess(), G4GammaGeneralProcess::GetEmProcess(), G4WeightWindowProcess::GetName(), G4ProcessManager::GetProcess(), G4ProcessManager::GetProcessVectorIndex(), G4GammaGeneralProcess::GetSubProcessName(), G4ProcessManager::InActivateProcess(), G4hhIonisation::InitialiseEnergyLossProcess(), G4ProcessTable::Insert(), G4ITStepProcessor::InvokeAlongStepDoItProcs(), G4SteppingManager::InvokeAlongStepDoItProcs(), G4SteppingManager::InvokeAtRestDoItProcs(), G4SteppingManager::InvokePSDIP(), G4LossTableManager::LocalPhysicsTables(), G4ErrorPropagator::MakeOneStep(), G4VEmProcess::PostStepDoIt(), G4ITSteppingVerbose::PostStepDoItAllDone(), G4SteppingVerbose::PostStepDoItAllDone(), G4SteppingVerboseWithUnits::PostStepDoItAllDone(), G4ITSteppingVerbose::PostStepDoItOneByOne(), G4SteppingVerbose::PostStepDoItOneByOne(), G4SteppingVerboseWithUnits::PostStepDoItOneByOne(), G4VITDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4DNASecondOrderReaction::PostStepGetPhysicalInteractionLength(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VITRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4ITSteppingVerbose::PostStepVerbose(), G4EmConfigurator::PrepareModels(), G4HadronStoppingProcess::PreparePhysicsTable(), G4GammaGeneralProcess::PreparePhysicsTable(), G4VEmProcess::PreparePhysicsTable(), G4VEnergyLossProcess::PreparePhysicsTable(), G4VMultipleScattering::PreparePhysicsTable(), G4LossTableManager::PreparePhysicsTable(), G4HadronicProcessStore::Print(), G4HadronicProcessStore::PrintHtml(), G4AnnihiToMuPair::PrintInfoDefinition(), G4GammaConversionToMuons::PrintInfoDefinition(), G4hImpactIonisation::PrintInfoDefinition(), G4ProcessPlacer::PrintProcVec(), G4VEnergyLossProcess::PrintWarning(), G4VEmProcess::PrintWarning(), G4SynchrotronRadiation::ProcessDescription(), G4Decay::ProcessDescription(), G4DecayWithSpin::ProcessDescription(), G4PionDecayMakeSpin::ProcessDescription(), G4UnknownDecay::ProcessDescription(), G4ChannelingOptrChangeCrossSection::ProposeOccurenceBiasingOperation(), G4StackManager::PushOneTrack(), G4HadronicProcessStore::Register(), G4LossTableManager::Register(), G4LossTableManager::RegisterExtraParticle(), G4HadronicProcessStore::RegisterExtraProcess(), G4HadronicProcessStore::RegisterParticle(), G4WrapperProcess::RegisterProcess(), G4PhysicsListHelper::RegisterProcess(), G4ProcessTable::Remove(), G4ParallelGeometriesLimiterProcess::RemoveParallelWorld(), G4ProcessManager::RemoveProcess(), G4ProcessPlacer::RemoveProcess(), G4GammaGeneralProcess::RetrievePhysicsTable(), G4VEmProcess::RetrievePhysicsTable(), G4VEnergyLossProcess::RetrievePhysicsTable(), G4VEmProcess::SetCrossSectionBiasingFactor(), G4VEnergyLossProcess::SetCrossSectionBiasingFactor(), G4VEnergyLossProcess::SetCSDARangeTable(), G4VEnergyLossProcess::SetInverseRangeTable(), G4VEnergyLossProcess::SetLambdaTable(), G4ProcessTableMessenger::SetNewValue(), G4ProcessTable::SetProcessActivation(), G4ProcessManager::SetProcessOrdering(), G4ProcessManager::SetProcessOrderingToFirst(), G4ProcessManager::SetProcessOrderingToLast(), G4ProcessManager::SetProcessOrderingToSecond(), G4VEnergyLossProcess::SetRangeTableForLoss(), G4VEnergyLossProcess::SetSecondaryRangeTable(), G4FastSimulationManagerProcess::SetWorldVolume(), G4ITSteppingVerbose::ShowStep(), G4SteppingVerbose::ShowStep(), G4SteppingVerboseWithUnits::ShowStep(), G4ChannelingOptrChangeCrossSection::StartRun(), G4ITSteppingVerbose::StepInfo(), G4SteppingVerbose::StepInfo(), G4SteppingVerboseWithUnits::StepInfo(), G4ITSteppingVerbose::StepInfoForLeadingTrack(), G4VEmProcess::StorePhysicsTable(), G4VMultipleScattering::StorePhysicsTable(), G4VEnergyLossProcess::StreamInfo(), G4VEmProcess::StreamInfo(), G4VMultipleScattering::StreamInfo(), G4EmCalculator::UpdateParticle(), G4ParallelWorldScoringProcess::Verbose(), G4ScoreSplittingProcess::Verbose(), G4ITSteppingVerbose::VerboseTrack(), G4SteppingVerbose::VerboseTrack(), and G4SteppingVerboseWithUnits::VerboseTrack().

◆ GetProcessSubType()

G4int G4VProcess::GetProcessSubType ( ) const
inlineinherited

◆ GetProcessType()

G4ProcessType G4VProcess::GetProcessType ( ) const
inlineinherited

◆ GetProcessTypeName()

const G4String & G4VProcess::GetProcessTypeName ( G4ProcessType  aType)
staticinherited

Definition at line 134 of file G4VProcess.cc.

135{
136 switch (aType)
137 {
138 case fNotDefined: return typeNotDefined; break;
139 case fTransportation: return typeTransportation; break;
140 case fElectromagnetic: return typeElectromagnetic; break;
141 case fOptical: return typeOptical; break;
142 case fHadronic: return typeHadronic; break;
144 case fDecay: return typeDecay; break;
145 case fGeneral: return typeGeneral; break;
146 case fParameterisation: return typeParameterisation; break;
147 case fUserDefined: return typeUserDefined; break;
148 case fPhonon: return typePhonon; break;
149 default: ;
150 }
151 return noType;
152}
@ fOptical
@ fPhonon
@ fParameterisation
@ fGeneral
@ fElectromagnetic
@ fHadronic
@ fUserDefined
@ fTransportation
@ fPhotolepton_hadron
@ fNotDefined
static const G4String typeNotDefined
Definition: G4VProcess.cc:119
static const G4String typeParameterisation
Definition: G4VProcess.cc:127
static const G4String typePhotolepton_hadron
Definition: G4VProcess.cc:124
static const G4String typeElectromagnetic
Definition: G4VProcess.cc:121
static const G4String noType
Definition: G4VProcess.cc:130
static const G4String typeUserDefined
Definition: G4VProcess.cc:128
static const G4String typeDecay
Definition: G4VProcess.cc:125
static const G4String typeTransportation
Definition: G4VProcess.cc:120
static const G4String typeHadronic
Definition: G4VProcess.cc:123
static const G4String typeOptical
Definition: G4VProcess.cc:122
static const G4String typeGeneral
Definition: G4VProcess.cc:126
static const G4String typePhonon
Definition: G4VProcess.cc:129

References fDecay, fElectromagnetic, fGeneral, fHadronic, fNotDefined, fOptical, fParameterisation, fPhonon, fPhotolepton_hadron, fTransportation, fUserDefined, anonymous_namespace{G4VProcess.cc}::noType, anonymous_namespace{G4VProcess.cc}::typeDecay, anonymous_namespace{G4VProcess.cc}::typeElectromagnetic, anonymous_namespace{G4VProcess.cc}::typeGeneral, anonymous_namespace{G4VProcess.cc}::typeHadronic, anonymous_namespace{G4VProcess.cc}::typeNotDefined, anonymous_namespace{G4VProcess.cc}::typeOptical, anonymous_namespace{G4VProcess.cc}::typeParameterisation, anonymous_namespace{G4VProcess.cc}::typePhonon, anonymous_namespace{G4VProcess.cc}::typePhotolepton_hadron, anonymous_namespace{G4VProcess.cc}::typeTransportation, and anonymous_namespace{G4VProcess.cc}::typeUserDefined.

Referenced by G4RichTrajectory::CreateAttValues(), G4RichTrajectoryPoint::CreateAttValues(), G4ProcessManager::DumpInfo(), G4VProcess::DumpInfo(), G4ProcessTableMessenger::G4ProcessTableMessenger(), G4ProcessTableMessenger::GetProcessType(), G4ProcessTableMessenger::GetProcessTypeName(), and G4ProcessTableMessenger::SetNumberOfProcessType().

◆ GetTotalNumberOfInteractionLengthTraversed()

G4double G4VProcess::GetTotalNumberOfInteractionLengthTraversed ( ) const
inlineinherited

◆ GetVerboseLevel()

G4int G4MuonicAtomDecay::GetVerboseLevel ( ) const
inline

Definition at line 74 of file G4MuonicAtomDecay.hh.

74{return verboseLevel;}

References verboseLevel.

Referenced by DecayIt(), FillResult(), GetMeanFreePath(), and GetMeanLifeTime().

◆ isAlongStepDoItIsEnabled()

G4bool G4VProcess::isAlongStepDoItIsEnabled ( ) const
inlineinherited

Definition at line 506 of file G4VProcess.hh.

507{
508 return enableAlongStepDoIt;
509}
G4bool enableAlongStepDoIt
Definition: G4VProcess.hh:360

References G4VProcess::enableAlongStepDoIt.

Referenced by G4ProcessManager::CheckOrderingParameters().

◆ IsApplicable()

G4bool G4MuonicAtomDecay::IsApplicable ( const G4ParticleDefinition a)
virtual

Reimplemented from G4VProcess.

Definition at line 91 of file G4MuonicAtomDecay.cc.

92{
93 return ( a.GetParticleType() == "MuonicAtom" );
94}
const G4String & GetParticleType() const

References G4ParticleDefinition::GetParticleType().

◆ isAtRestDoItIsEnabled()

G4bool G4VProcess::isAtRestDoItIsEnabled ( ) const
inlineinherited

Definition at line 500 of file G4VProcess.hh.

501{
502 return enableAtRestDoIt;
503}

References G4VProcess::enableAtRestDoIt.

Referenced by G4ProcessManager::CheckOrderingParameters().

◆ isPostStepDoItIsEnabled()

G4bool G4VProcess::isPostStepDoItIsEnabled ( ) const
inlineinherited

Definition at line 512 of file G4VProcess.hh.

513{
514 return enablePostStepDoIt;
515}

References G4VProcess::enablePostStepDoIt.

Referenced by G4ProcessManager::CheckOrderingParameters().

◆ operator!=()

G4bool G4VProcess::operator!= ( const G4VProcess right) const
inherited

Definition at line 161 of file G4VProcess.cc.

162{
163 return (this != &right);
164}

◆ operator=()

G4MuonicAtomDecay & G4MuonicAtomDecay::operator= ( const G4MuonicAtomDecay right)
private

◆ operator==()

G4bool G4VProcess::operator== ( const G4VProcess right) const
inherited

Definition at line 155 of file G4VProcess.cc.

156{
157 return (this == &right);
158}

◆ PostStepDoIt()

virtual G4VParticleChange * G4MuonicAtomDecay::PostStepDoIt ( const G4Track theTrack,
const G4Step theStep 
)
inlinevirtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 90 of file G4MuonicAtomDecay.hh.

92 {return DecayIt(theTrack, theStep);}

References DecayIt().

◆ PostStepGetPhysicalInteractionLength()

G4double G4MuonicAtomDecay::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
virtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 126 of file G4MuonicAtomDecay.cc.

128{
130 return DBL_MAX; // this will need to be changed in future; fixme
131}

References condition(), DBL_MAX, and NotForced.

◆ PostStepGPIL()

G4double G4VProcess::PostStepGPIL ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
inlineinherited

Definition at line 479 of file G4VProcess.hh.

482{
483 return thePILfactor *
484 PostStepGetPhysicalInteractionLength(track, previousStepSize, condition);
485}
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0

References condition(), G4VProcess::PostStepGetPhysicalInteractionLength(), and G4VProcess::thePILfactor.

Referenced by G4SteppingManager::DefinePhysicalStepLength(), and G4ITStepProcessor::DoDefinePhysicalStepLength().

◆ PreparePhysicsTable()

virtual void G4VProcess::PreparePhysicsTable ( const G4ParticleDefinition )
inlinevirtualinherited

◆ PrepareWorkerPhysicsTable()

void G4VProcess::PrepareWorkerPhysicsTable ( const G4ParticleDefinition part)
virtualinherited

Reimplemented in G4BiasingProcessInterface.

Definition at line 206 of file G4VProcess.cc.

207{
209}
virtual void PreparePhysicsTable(const G4ParticleDefinition &)
Definition: G4VProcess.hh:194

References G4VProcess::PreparePhysicsTable().

Referenced by G4BiasingProcessInterface::PrepareWorkerPhysicsTable().

◆ ProcessDescription()

void G4MuonicAtomDecay::ProcessDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VProcess.

Definition at line 443 of file G4MuonicAtomDecay.cc.

444{
445 outFile << "MuonicAtom process where Muon decays in orbit or is captured by the nucleus." <<G4endl;
446}

References G4endl.

◆ ResetNumberOfInteractionLengthLeft()

void G4VProcess::ResetNumberOfInteractionLengthLeft ( )
virtualinherited

◆ RetrievePhysicsTable()

virtual G4bool G4VProcess::RetrievePhysicsTable ( const G4ParticleDefinition ,
const G4String ,
G4bool   
)
inlinevirtualinherited

◆ SetMasterProcess()

void G4VProcess::SetMasterProcess ( G4VProcess masterP)
virtualinherited

◆ SetPILfactor()

void G4VProcess::SetPILfactor ( G4double  value)
inlineinherited

Definition at line 449 of file G4VProcess.hh.

450{
451 if (value>0.) { thePILfactor = value; }
452}

References G4VProcess::thePILfactor.

Referenced by export_G4VProcess().

◆ SetProcessManager()

void G4VProcess::SetProcessManager ( const G4ProcessManager procMan)
inlinevirtualinherited

◆ SetProcessSubType()

void G4VProcess::SetProcessSubType ( G4int  value)
inlineinherited

Definition at line 406 of file G4VProcess.hh.

407{
408 theProcessSubType = value;
409}

References G4VProcess::theProcessSubType.

Referenced by G4DNAElectronHoleRecombination::Create(), G4DNASecondOrderReaction::Create(), G4AnnihiToMuPair::G4AnnihiToMuPair(), G4BiasingProcessInterface::G4BiasingProcessInterface(), G4Cerenkov::G4Cerenkov(), G4ComptonScattering::G4ComptonScattering(), G4CoulombScattering::G4CoulombScattering(), G4CoupledTransportation::G4CoupledTransportation(), G4Decay::G4Decay(), G4DecayWithSpin::G4DecayWithSpin(), G4DNAAttachment::G4DNAAttachment(), G4DNABrownianTransportation::G4DNABrownianTransportation(), G4DNAChargeDecrease::G4DNAChargeDecrease(), G4DNAChargeIncrease::G4DNAChargeIncrease(), G4DNAElastic::G4DNAElastic(), G4DNAElectronSolvation::G4DNAElectronSolvation(), G4DNAExcitation::G4DNAExcitation(), G4DNAIonisation::G4DNAIonisation(), G4DNAMolecularDissociation::G4DNAMolecularDissociation(), G4DNAScavengerProcess::G4DNAScavengerProcess(), G4DNAVibExcitation::G4DNAVibExcitation(), G4eBremsstrahlung::G4eBremsstrahlung(), G4eeToHadrons::G4eeToHadrons(), G4eIonisation::G4eIonisation(), G4ePairProduction::G4ePairProduction(), G4eplusAnnihilation::G4eplusAnnihilation(), G4FastSimulationManagerProcess::G4FastSimulationManagerProcess(), G4GammaConversion::G4GammaConversion(), G4GammaConversionToMuons::G4GammaConversionToMuons(), G4GammaGeneralProcess::G4GammaGeneralProcess(), G4HadronicProcess::G4HadronicProcess(), G4hhIonisation::G4hhIonisation(), G4hIonisation::G4hIonisation(), G4ionIonisation::G4ionIonisation(), G4ITTransportation::G4ITTransportation(), G4JAEAElasticScattering::G4JAEAElasticScattering(), G4MicroElecElastic::G4MicroElecElastic(), G4MicroElecInelastic::G4MicroElecInelastic(), G4MicroElecLOPhononScattering::G4MicroElecLOPhononScattering(), G4MicroElecSurface::G4MicroElecSurface(), G4mplIonisation::G4mplIonisation(), G4MuBremsstrahlung::G4MuBremsstrahlung(), G4MuIonisation::G4MuIonisation(), G4MuonMinusAtomicCapture::G4MuonMinusAtomicCapture(), G4MuPairProduction::G4MuPairProduction(), G4NeutronKiller::G4NeutronKiller(), G4NuclearStopping::G4NuclearStopping(), G4OpAbsorption::G4OpAbsorption(), G4OpBoundaryProcess::G4OpBoundaryProcess(), G4OpMieHG::G4OpMieHG(), G4OpRayleigh::G4OpRayleigh(), G4OpWLS::G4OpWLS(), G4OpWLS2::G4OpWLS2(), G4ParallelWorldProcess::G4ParallelWorldProcess(), G4PhotoElectricEffect::G4PhotoElectricEffect(), G4PionDecayMakeSpin::G4PionDecayMakeSpin(), G4PolarizedCompton::G4PolarizedCompton(), G4PolarizedGammaConversion::G4PolarizedGammaConversion(), G4PolarizedIonisation::G4PolarizedIonisation(), G4PolarizedPhotoElectric::G4PolarizedPhotoElectric(), G4RadioactiveDecay::G4RadioactiveDecay(), G4RayleighScattering::G4RayleighScattering(), G4Scintillation::G4Scintillation(), G4StepLimiter::G4StepLimiter(), G4SynchrotronRadiation::G4SynchrotronRadiation(), G4SynchrotronRadiationInMat::G4SynchrotronRadiationInMat(), G4TransitionRadiation::G4TransitionRadiation(), G4Transportation::G4Transportation(), G4UCNAbsorption::G4UCNAbsorption(), G4UCNBoundaryProcess::G4UCNBoundaryProcess(), G4UCNLoss::G4UCNLoss(), G4UCNMultiScattering::G4UCNMultiScattering(), G4UnknownDecay::G4UnknownDecay(), G4UserSpecialCuts::G4UserSpecialCuts(), G4VMultipleScattering::G4VMultipleScattering(), G4VTransitionRadiation::G4VTransitionRadiation(), G4VXTRenergyLoss::G4VXTRenergyLoss(), and G4Decay::SetExtDecayer().

◆ SetProcessType()

void G4VProcess::SetProcessType ( G4ProcessType  aType)
inlineinherited

Definition at line 394 of file G4VProcess.hh.

395{
396 theProcessType = aType;
397}

References G4VProcess::theProcessType.

Referenced by G4MaxTimeCuts::G4MaxTimeCuts(), and G4MinEkineCuts::G4MinEkineCuts().

◆ SetVerboseLevel()

void G4MuonicAtomDecay::SetVerboseLevel ( G4int  value)
inline

Definition at line 71 of file G4MuonicAtomDecay.hh.

71{verboseLevel = value;}

References verboseLevel.

◆ StartTracking()

void G4VProcess::StartTracking ( G4Track )
virtualinherited

Reimplemented in G4ParallelGeometriesLimiterProcess, G4ImportanceProcess, G4WeightCutOffProcess, G4WeightWindowProcess, G4VITProcess, G4DNASecondOrderReaction, G4WrapperProcess, G4FastSimulationManagerProcess, G4ParallelWorldProcess, G4ParallelWorldScoringProcess, G4ScoreSplittingProcess, G4GammaGeneralProcess, G4Decay, G4AdjointProcessEquivalentToDirectProcess, G4eAdjointMultipleScattering, G4DNAElectronHoleRecombination, G4DNAScavengerProcess, G4VEmProcess, G4VEnergyLossProcess, G4VMultipleScattering, G4ITTransportation, G4DNABrownianTransportation, G4CoupledTransportation, G4Transportation, G4BiasingProcessInterface, and G4VPhononProcess.

Definition at line 87 of file G4VProcess.cc.

88{
92#ifdef G4VERBOSE
93 if (verboseLevel>2)
94 {
95 G4cout << "G4VProcess::StartTracking() - [" << theProcessName << "]"
96 << G4endl;
97 }
98#endif
99}

References G4VProcess::currentInteractionLength, G4cout, G4endl, G4VProcess::theInitialNumberOfInteractionLength, G4VProcess::theNumberOfInteractionLengthLeft, G4VProcess::theProcessName, and G4VProcess::verboseLevel.

Referenced by G4DNASecondOrderReaction::StartTracking(), G4WrapperProcess::StartTracking(), G4AdjointProcessEquivalentToDirectProcess::StartTracking(), G4DNAElectronHoleRecombination::StartTracking(), G4DNAScavengerProcess::StartTracking(), G4ITTransportation::StartTracking(), G4Transportation::StartTracking(), G4BiasingProcessInterface::StartTracking(), and G4VPhononProcess::StartTracking().

◆ StorePhysicsTable()

virtual G4bool G4VProcess::StorePhysicsTable ( const G4ParticleDefinition ,
const G4String ,
G4bool   
)
inlinevirtualinherited

◆ SubtractNumberOfInteractionLengthLeft()

void G4VProcess::SubtractNumberOfInteractionLengthLeft ( G4double  prevStepSize)
inlineprotectedinherited

Definition at line 524 of file G4VProcess.hh.

525{
527 {
530 {
532 }
533 }
534 else
535 {
536#ifdef G4VERBOSE
537 if (verboseLevel>0)
538 {
539 G4cerr << "G4VProcess::SubtractNumberOfInteractionLengthLeft()";
540 G4cerr << " [" << theProcessName << "]" <<G4endl;
541 G4cerr << " currentInteractionLength = "
542 << currentInteractionLength << " [mm]";
543 G4cerr << " previousStepSize = " << prevStepSize << " [mm]";
544 G4cerr << G4endl;
545 }
546#endif
547 G4String msg = "Negative currentInteractionLength for ";
548 msg += theProcessName;
549 G4Exception("G4VProcess::SubtractNumberOfInteractionLengthLeft()",
550 "ProcMan201", EventMustBeAborted, msg);
551 }
552}
@ EventMustBeAborted
G4GLOB_DLL std::ostream G4cerr
static constexpr double perMillion

References G4VProcess::currentInteractionLength, EventMustBeAborted, G4cerr, G4endl, G4Exception(), CLHEP::perMillion, G4VProcess::theNumberOfInteractionLengthLeft, G4VProcess::theProcessName, and G4VProcess::verboseLevel.

Referenced by G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), and G4Decay::PostStepGetPhysicalInteractionLength().

Field Documentation

◆ aParticleChange

G4ParticleChange G4VProcess::aParticleChange
protectedinherited

◆ aProcessManager

const G4ProcessManager* G4VProcess::aProcessManager = nullptr
protectedinherited

Definition at line 319 of file G4VProcess.hh.

Referenced by G4VProcess::GetProcessManager(), and G4VProcess::SetProcessManager().

◆ cmptr

G4HadronicInteraction* G4MuonicAtomDecay::cmptr
private

Definition at line 125 of file G4MuonicAtomDecay.hh.

Referenced by DecayIt(), and G4MuonicAtomDecay().

◆ currentInteractionLength

G4double G4VProcess::currentInteractionLength = -1.0
protectedinherited

◆ enableAlongStepDoIt

G4bool G4VProcess::enableAlongStepDoIt = true
protectedinherited

◆ enableAtRestDoIt

G4bool G4VProcess::enableAtRestDoIt = true
protectedinherited

◆ enablePostStepDoIt

G4bool G4VProcess::enablePostStepDoIt = true
protectedinherited

◆ fMuMass

const G4double G4MuonicAtomDecay::fMuMass
private

Definition at line 123 of file G4MuonicAtomDecay.hh.

◆ fProcessTable

G4ProcessTable* G4VProcess::fProcessTable = nullptr
privateinherited

Definition at line 374 of file G4VProcess.hh.

Referenced by G4VProcess::G4VProcess(), and G4VProcess::~G4VProcess().

◆ masterProcessShadow

G4VProcess* G4VProcess::masterProcessShadow = nullptr
privateinherited

Definition at line 370 of file G4VProcess.hh.

Referenced by G4VProcess::GetMasterProcess(), and G4VProcess::SetMasterProcess().

◆ pParticleChange

G4VParticleChange* G4VProcess::pParticleChange = nullptr
protectedinherited

Definition at line 321 of file G4VProcess.hh.

Referenced by G4VMultipleScattering::AddEmModel(), G4VEmProcess::AddEmModel(), G4VEnergyLossProcess::AddEmModel(), G4ImportanceProcess::AlongStepDoIt(), G4WeightCutOffProcess::AlongStepDoIt(), G4WeightWindowProcess::AlongStepDoIt(), G4VContinuousDiscreteProcess::AlongStepDoIt(), G4VContinuousProcess::AlongStepDoIt(), G4VRestContinuousDiscreteProcess::AlongStepDoIt(), G4VRestContinuousProcess::AlongStepDoIt(), G4ParallelWorldProcess::AlongStepDoIt(), G4ParallelWorldScoringProcess::AlongStepDoIt(), G4VITRestProcess::AtRestDoIt(), G4VRestContinuousDiscreteProcess::AtRestDoIt(), G4VRestContinuousProcess::AtRestDoIt(), G4VRestDiscreteProcess::AtRestDoIt(), G4VRestProcess::AtRestDoIt(), G4ParallelWorldProcess::AtRestDoIt(), G4ParallelWorldScoringProcess::AtRestDoIt(), G4ScoreSplittingProcess::AtRestDoIt(), G4VITRestDiscreteProcess::AtRestDoIt(), G4eplusAnnihilation::AtRestDoIt(), G4DNAElectronHoleRecombination::Create(), G4DNASecondOrderReaction::Create(), G4VEnergyLossProcess::FillSecondariesAlongStep(), G4Decay::G4Decay(), G4DNAMolecularDissociation::G4DNAMolecularDissociation(), G4DNAScavengerProcess::G4DNAScavengerProcess(), G4ImportanceProcess::G4ImportanceProcess(), G4ITTransportation::G4ITTransportation(), G4ParallelWorldProcess::G4ParallelWorldProcess(), G4ParallelWorldScoringProcess::G4ParallelWorldScoringProcess(), G4RadioactiveDecay::G4RadioactiveDecay(), G4ScoreSplittingProcess::G4ScoreSplittingProcess(), G4Transportation::G4Transportation(), G4UnknownDecay::G4UnknownDecay(), G4VEmProcess::G4VEmProcess(), G4VEnergyLossProcess::G4VEnergyLossProcess(), G4VMultipleScattering::G4VMultipleScattering(), G4VProcess::G4VProcess(), G4VXTRenergyLoss::G4VXTRenergyLoss(), G4WeightCutOffProcess::G4WeightCutOffProcess(), G4WeightWindowProcess::G4WeightWindowProcess(), G4VITDiscreteProcess::PostStepDoIt(), G4VContinuousDiscreteProcess::PostStepDoIt(), G4VDiscreteProcess::PostStepDoIt(), G4VRestContinuousDiscreteProcess::PostStepDoIt(), G4VRestDiscreteProcess::PostStepDoIt(), G4ParallelWorldProcess::PostStepDoIt(), G4ParallelWorldScoringProcess::PostStepDoIt(), G4ScoreSplittingProcess::PostStepDoIt(), G4NeutronKiller::PostStepDoIt(), G4VITRestDiscreteProcess::PostStepDoIt(), G4LowECapture::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4Cerenkov::PostStepDoIt(), and G4VTransitionRadiation::PostStepDoIt().

◆ theInitialNumberOfInteractionLength

G4double G4VProcess::theInitialNumberOfInteractionLength = -1.0
protectedinherited

◆ theNumberOfInteractionLengthLeft

G4double G4VProcess::theNumberOfInteractionLengthLeft = -1.0
protectedinherited

Definition at line 331 of file G4VProcess.hh.

Referenced by G4AdjointForcedInteractionForGamma::AlongStepDoIt(), AtRestGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestContinuousProcess::AtRestGetPhysicalInteractionLength(), G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), G4VRestProcess::AtRestGetPhysicalInteractionLength(), G4Decay::AtRestGetPhysicalInteractionLength(), G4VProcess::ClearNumberOfInteractionLengthLeft(), DecayIt(), G4VProcess::EndTracking(), G4VProcess::GetNumberOfInteractionLengthLeft(), G4VProcess::GetTotalNumberOfInteractionLengthTraversed(), G4GammaGeneralProcess::PostStepDoIt(), G4VEmProcess::PostStepDoIt(), G4VEnergyLossProcess::PostStepDoIt(), G4VContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestContinuousDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4VRestDiscreteProcess::PostStepGetPhysicalInteractionLength(), G4GammaGeneralProcess::PostStepGetPhysicalInteractionLength(), G4Decay::PostStepGetPhysicalInteractionLength(), G4AdjointForcedInteractionForGamma::PostStepGetPhysicalInteractionLength(), G4PolarizedAnnihilation::PostStepGetPhysicalInteractionLength(), G4PolarizedCompton::PostStepGetPhysicalInteractionLength(), G4PolarizedIonisation::PostStepGetPhysicalInteractionLength(), G4VEmProcess::PostStepGetPhysicalInteractionLength(), G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength(), G4VProcess::ResetNumberOfInteractionLengthLeft(), G4VProcess::StartTracking(), G4GammaGeneralProcess::StartTracking(), G4VEmProcess::StartTracking(), G4VEnergyLossProcess::StartTracking(), and G4VProcess::SubtractNumberOfInteractionLengthLeft().

◆ thePhysicsTableFileName

G4String G4VProcess::thePhysicsTableFileName
protectedinherited

Definition at line 344 of file G4VProcess.hh.

Referenced by G4VProcess::GetPhysicsTableFileName().

◆ thePILfactor

G4double G4VProcess::thePILfactor = 1.0
protectedinherited

◆ theProcessName

G4String G4VProcess::theProcessName
protectedinherited

◆ theProcessSubType

G4int G4VProcess::theProcessSubType = -1
protectedinherited

◆ theProcessType

G4ProcessType G4VProcess::theProcessType = fNotDefined
protectedinherited

◆ theTotalResult

G4ParticleChange G4MuonicAtomDecay::theTotalResult
private

Definition at line 121 of file G4MuonicAtomDecay.hh.

Referenced by DecayIt(), and FillResult().

◆ verboseLevel

G4int G4MuonicAtomDecay::verboseLevel
private

Definition at line 128 of file G4MuonicAtomDecay.hh.

Referenced by GetVerboseLevel(), and SetVerboseLevel().


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