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

#include <G4VMultipleScattering.hh>

+ Inheritance diagram for G4VMultipleScattering:

Public Member Functions

 G4VMultipleScattering (const G4String &name="msc", G4ProcessType type=fElectromagnetic)
 
virtual ~G4VMultipleScattering ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p) override=0
 
virtual void PrintInfo ()
 
virtual void ProcessDescription (std::ostream &outFile) const override
 
void PreparePhysicsTable (const G4ParticleDefinition &) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
 
G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
 
void StartTracking (G4Track *) override
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
G4double ContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
 
G4VEmModelSelectModel (G4double kinEnergy, size_t idx)
 
void AddEmModel (G4int order, G4VEmModel *, const G4Region *region=nullptr)
 
void SetEmModel (G4VMscModel *, size_t index=0)
 
G4VMscModelEmModel (size_t index=0) const
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
void SetIonisation (G4VEnergyLossProcess *)
 
G4bool LateralDisplasmentFlag () const
 
void SetLateralDisplasmentFlag (G4bool val)
 
G4double Skin () const
 
void SetSkin (G4double val)
 
G4double RangeFactor () const
 
void SetRangeFactor (G4double val)
 
G4double GeomFactor () const
 
G4double PolarAngleLimit () const
 
G4MscStepLimitType StepLimitType () const
 
void SetStepLimitType (G4MscStepLimitType val)
 
G4double LowestKinEnergy () const
 
void SetLowestKinEnergy (G4double val)
 
const G4ParticleDefinitionFirstParticle () const
 
- Public Member Functions inherited from G4VContinuousDiscreteProcess
 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
 
virtual ~G4VContinuousDiscreteProcess ()
 
G4VContinuousDiscreteProcessoperator= (const G4VContinuousDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

virtual void InitialiseProcess (const G4ParticleDefinition *)=0
 
virtual void StreamProcessInfo (std::ostream &) const
 
G4double GetMeanFreePath (const G4Track &track, G4double, G4ForceCondition *condition) override
 
G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety) override
 
G4int NumberOfModels () const
 
- Protected Member Functions inherited from G4VContinuousDiscreteProcess
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
virtual G4double GetContinuousStepLimit (const G4Track &aTrack, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
 
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4ParticleChangeForMSC fParticleChange
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Additional Inherited Members

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

Detailed Description

Definition at line 91 of file G4VMultipleScattering.hh.

Constructor & Destructor Documentation

◆ G4VMultipleScattering()

G4VMultipleScattering::G4VMultipleScattering ( const G4String name = "msc",
G4ProcessType  type = fElectromagnetic 
)

Definition at line 90 of file G4VMultipleScattering.cc.

92 numberOfModels(0),
93 firstParticle(nullptr),
94 currParticle(nullptr),
95 stepLimit(fUseSafety),
96 facrange(0.04),
97 latDisplacement(true),
98 isIon(false),
99 fNewPosition(0.,0.,0.),
100 fNewDirection(0.,0.,1.)
101{
102 theParameters = G4EmParameters::Instance();
105 if("ionmsc" == name) { firstParticle = G4GenericIon::GenericIon(); }
106
107 lowestKinEnergy = 10*CLHEP::eV;
108
109 physStepLimit = gPathLength = tPathLength = 0.0;
110 fIonisation = nullptr;
111
112 geomMin = 0.05*CLHEP::nm;
113 minDisplacement2 = geomMin*geomMin;
114
116 safetyHelper = nullptr;
117 fPositionChanged = false;
118 isActive = false;
119
120 currentModel = nullptr;
121 modelManager = new G4EmModelManager();
122 emManager = G4LossTableManager::Instance();
123 mscModels.reserve(2);
124 emManager->Register(this);
125}
@ fMultipleScattering
@ fUseSafety
@ fElectromagnetic
static G4EmParameters * Instance()
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
static G4LossTableManager * Instance()
void Register(G4VEnergyLossProcess *p)
G4ParticleChangeForMSC fParticleChange
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321

◆ ~G4VMultipleScattering()

G4VMultipleScattering::~G4VMultipleScattering ( )
virtual

Definition at line 129 of file G4VMultipleScattering.cc.

130{
131 /*
132 if(1 < verboseLevel) {
133 G4cout << "G4VMultipleScattering destruct " << GetProcessName()
134 << G4endl;
135 }
136 */
137 delete modelManager;
138 emManager->DeRegister(this);
139}
void DeRegister(G4VEnergyLossProcess *p)

Member Function Documentation

◆ AddEmModel()

void G4VMultipleScattering::AddEmModel ( G4int  order,
G4VEmModel p,
const G4Region region = nullptr 
)

Definition at line 143 of file G4VMultipleScattering.cc.

145{
146 if(!p) { return; }
147 G4VEmFluctuationModel* fm = nullptr;
148 modelManager->AddEmModel(order, p, fm, region);
150}
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *, const G4Region *)
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:456

Referenced by G4AdjointhMultipleScattering::InitialiseProcess(), G4eAdjointMultipleScattering::InitialiseProcess(), G4MuMultipleScattering::InitialiseProcess(), G4eMultipleScattering::InitialiseProcess(), G4hMultipleScattering::InitialiseProcess(), and G4EmConfigurator::PrepareModels().

◆ AlongStepDoIt()

G4VParticleChange * G4VMultipleScattering::AlongStepDoIt ( const G4Track track,
const G4Step step 
)
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 457 of file G4VMultipleScattering.cc.

458{
461 fNewPosition = step.GetPostStepPoint()->GetPosition();
462 fParticleChange.ProposePosition(fNewPosition);
463 fPositionChanged = false;
464
465 G4double geomLength = step.GetStepLength();
466
467 // very small step - no msc
468 if(!isActive) {
469 tPathLength = geomLength;
470
471 // sample msc
472 } else {
473 G4double range =
474 currentModel->GetRange(currParticle,track.GetKineticEnergy(),
475 track.GetMaterialCutsCouple());
476
477 tPathLength = currentModel->ComputeTrueStepLength(geomLength);
478
479 /*
480 if(currParticle->GetPDGMass() > 0.9*GeV)
481 G4cout << "G4VMsc::AlongStepDoIt: GeomLength= "
482 << geomLength
483 << " tPathLength= " << tPathLength
484 << " physStepLimit= " << physStepLimit
485 << " dr= " << range - tPathLength
486 << " ekin= " << track.GetKineticEnergy() << G4endl;
487 */
488 // protection against wrong t->g->t conversion
489 tPathLength = std::min(tPathLength, physStepLimit);
490
491 // do not sample scattering at the last or at a small step
492 if(tPathLength < range && tPathLength > geomMin) {
493
494 static const G4double minSafety = 1.20*CLHEP::nm;
495 static const G4double sFact = 0.99;
496
497 G4ThreeVector displacement = currentModel->SampleScattering(
498 step.GetPostStepPoint()->GetMomentumDirection(),minSafety);
499
500 G4double r2 = displacement.mag2();
501 //G4cout << " R= " << sqrt(r2) << " Rmin= " << sqrt(minDisplacement2)
502 // << " flag= " << fDispBeyondSafety << G4endl;
503 if(r2 > minDisplacement2) {
504
505 fPositionChanged = true;
506 G4double dispR = std::sqrt(r2);
507 G4double postSafety =
508 sFact*safetyHelper->ComputeSafety(fNewPosition, dispR);
509 //G4cout<<" R= "<< dispR<<" postSafety= "<<postSafety<<G4endl;
510
511 // far away from geometry boundary
512 if(postSafety > 0.0 && dispR <= postSafety) {
513 fNewPosition += displacement;
514
515 //near the boundary
516 } else {
517 // displaced point is definitely within the volume
518 //G4cout<<" R= "<<dispR<<" postSafety= "<<postSafety<<G4endl;
519 if(dispR < postSafety) {
520 fNewPosition += displacement;
521
522 // reduced displacement
523 } else if(postSafety > geomMin) {
524 fNewPosition += displacement*(postSafety/dispR);
525
526 // very small postSafety
527 } else {
528 fPositionChanged = false;
529 }
530 }
531 if(fPositionChanged) {
532 safetyHelper->ReLocateWithinVolume(fNewPosition);
533 fParticleChange.ProposePosition(fNewPosition);
534 }
535 }
536 }
537 }
539 return &fParticleChange;
540}
double G4double
Definition: G4Types.hh:83
double mag2() const
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void ProposePosition(const G4ThreeVector &finalPosition)
G4double ComputeSafety(const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)
void ReLocateWithinVolume(const G4ThreeVector &pGlobalPoint)
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual G4double ComputeTrueStepLength(G4double geomPathLength)=0
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:330
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)=0
void ProposeTrueStepLength(G4double truePathLength)

◆ AlongStepGetPhysicalInteractionLength()

G4double G4VMultipleScattering::AlongStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimalStep,
G4double currentSafety,
G4GPILSelection selection 
)
overridevirtual

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 391 of file G4VMultipleScattering.cc.

397{
398 // get Step limit proposed by the process
399 *selection = NotCandidateForSelection;
400 physStepLimit = gPathLength = tPathLength = currentMinimalStep;
401
402 G4double ekin = track.GetKineticEnergy();
403 /*
404 G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
405 << " " << currParticle->GetParticleName()
406 << " currMod " << currentModel
407 << G4endl;
408 */
409 // isIon flag is used only to select a model
410 if(isIon) {
411 ekin *= proton_mass_c2/track.GetParticleDefinition()->GetPDGMass();
412 }
413
414 // select new model
415 if(1 < numberOfModels) {
416 currentModel = static_cast<G4VMscModel*>(
418 }
419 // msc is active is model is active, energy above the limit,
420 // and step size is above the limit;
421 // if it is active msc may limit the step
422 if(currentModel->IsActive(ekin) && tPathLength > geomMin
423 && ekin >= lowestKinEnergy) {
424 isActive = true;
425 tPathLength =
426 currentModel->ComputeTruePathLengthLimit(track, gPathLength);
427 if (tPathLength < physStepLimit) {
428 *selection = CandidateForSelection;
429 }
430 } else { isActive = false; }
431
432 //if(currParticle->GetPDGMass() > GeV)
433 /*
434 G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
435 << " " << currParticle->GetParticleName()
436 << " gPathLength= " << gPathLength
437 << " tPathLength= " << tPathLength
438 << " currentMinimalStep= " << currentMinimalStep
439 << " isActive " << isActive << G4endl;
440 */
441 return gPathLength;
442}
@ CandidateForSelection
@ NotCandidateForSelection
const G4ParticleDefinition * GetParticleDefinition() const
G4bool IsActive(G4double kinEnergy) const
Definition: G4VEmModel.hh:785
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)=0
G4VEmModel * SelectModel(G4double kinEnergy, size_t idx)

Referenced by GetContinuousStepLimit().

◆ BuildPhysicsTable()

void G4VMultipleScattering::BuildPhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 276 of file G4VMultipleScattering.cc.

277{
278 G4String num = part.GetParticleName();
279 G4bool master = emManager->IsMaster();
280 if(1 < verboseLevel) {
281 G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
282 << GetProcessName()
283 << " and particle " << num << " isIon: " << isIon
284 << " IsMaster: " << master << G4endl;
285 }
286 const G4VMultipleScattering* masterProcess =
287 static_cast<const G4VMultipleScattering*>(GetMasterProcess());
288
289 if(firstParticle == &part) {
290 /*
291 G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
292 << GetProcessName()
293 << " and particle " << num
294 << " IsMaster= " << G4LossTableManager::Instance()->IsMaster()
295 << " " << this
296 << G4endl;
297 */
298 emManager->BuildPhysicsTable(firstParticle);
299
300 if(!master) {
301 // initialisation of models
302 /*
303 G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
304 << GetProcessName()
305 << " and particle " << num
306 << " Nmod= " << numberOfModels << " " << this
307 << G4endl;
308 */
309 for(G4int i=0; i<numberOfModels; ++i) {
310 G4VMscModel* msc = static_cast<G4VMscModel*>(GetModelByIndex(i));
311 if(!msc) { continue; }
312 G4VMscModel* msc0=
313 static_cast<G4VMscModel*>(masterProcess->GetModelByIndex(i));
314 msc->SetCrossSectionTable(msc0->GetCrossSectionTable(), false);
315 msc->InitialiseLocal(firstParticle, msc0);
316 }
317 }
318 }
319
320 // explicitly defined printout by particle name
321 if(1 < verboseLevel ||
322 (0 < verboseLevel && (num == "e-" ||
323 num == "e+" || num == "mu+" ||
324 num == "mu-" || num == "proton"||
325 num == "pi+" || num == "pi-" ||
326 num == "kaon+" || num == "kaon-" ||
327 num == "alpha" || num == "anti_proton" ||
328 num == "GenericIon" || num == "alpha+" ||
329 num == "alpha++" )))
330 {
331 StreamInfo(G4cout, part);
332 }
333
334 if(1 < verboseLevel) {
335 G4cout << "### G4VMultipleScattering::BuildPhysicsTable() done for "
336 << GetProcessName()
337 << " and particle " << num
338 << G4endl;
339 }
340}
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
const G4String & GetParticleName() const
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:464
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:860
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:220
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:518
G4int verboseLevel
Definition: G4VProcess.hh:356
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382

◆ ContinuousStepLimit()

G4double G4VMultipleScattering::ContinuousStepLimit ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimalStep,
G4double currentSafety 
)

Definition at line 569 of file G4VMultipleScattering.cc.

574{
575 return GetContinuousStepLimit(track,previousStepSize,currentMinimalStep,
576 currentSafety);
577}
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety) override

◆ EmModel()

G4VMscModel * G4VMultipleScattering::EmModel ( size_t  index = 0) const

◆ FirstParticle()

const G4ParticleDefinition * G4VMultipleScattering::FirstParticle ( ) const
inline

Definition at line 406 of file G4VMultipleScattering.hh.

407{
408 return firstParticle;
409}

◆ GeomFactor()

G4double G4VMultipleScattering::GeomFactor ( ) const
inline

◆ GetContinuousStepLimit()

G4double G4VMultipleScattering::GetContinuousStepLimit ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimalStep,
G4double currentSafety 
)
overrideprotectedvirtual

Implements G4VContinuousDiscreteProcess.

Definition at line 553 of file G4VMultipleScattering.cc.

558{
560 G4double x = AlongStepGetPhysicalInteractionLength(track,previousStepSize,
561 currentMinimalStep,
562 currentSafety,
563 &selection);
564 return x;
565}
G4GPILSelection
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection) override

Referenced by ContinuousStepLimit().

◆ GetMeanFreePath()

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

Implements G4VContinuousDiscreteProcess.

Definition at line 581 of file G4VMultipleScattering.cc.

583{
584 *condition = Forced;
585 return DBL_MAX;
586}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced
#define DBL_MAX
Definition: templates.hh:62

◆ GetModelByIndex()

G4VEmModel * G4VMultipleScattering::GetModelByIndex ( G4int  idx = 0,
G4bool  ver = false 
) const
inline

Definition at line 421 of file G4VMultipleScattering.hh.

422{
423 return modelManager->GetModel(idx, ver);
424}
G4VEmModel * GetModel(G4int idx, G4bool ver=false)

Referenced by BuildPhysicsTable(), PreparePhysicsTable(), SetIonisation(), and StartTracking().

◆ InitialiseProcess()

virtual void G4VMultipleScattering::InitialiseProcess ( const G4ParticleDefinition )
protectedpure virtual

◆ IsApplicable()

virtual G4bool G4VMultipleScattering::IsApplicable ( const G4ParticleDefinition p)
overridepure virtual

◆ LateralDisplasmentFlag()

◆ LowestKinEnergy()

G4double G4VMultipleScattering::LowestKinEnergy ( ) const
inline

Definition at line 392 of file G4VMultipleScattering.hh.

393{
394 return lowestKinEnergy;
395}

◆ NumberOfModels()

G4int G4VMultipleScattering::NumberOfModels ( ) const
inlineprotected

Definition at line 413 of file G4VMultipleScattering.hh.

414{
415 return modelManager->NumberOfModels();
416}
G4int NumberOfModels() const

◆ PolarAngleLimit()

G4double G4VMultipleScattering::PolarAngleLimit ( ) const
inline

Definition at line 369 of file G4VMultipleScattering.hh.

370{
371 return theParameters->MscThetaLimit();
372}
G4double MscThetaLimit() const

Referenced by G4MuMultipleScattering::StreamProcessInfo().

◆ PostStepDoIt()

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 545 of file G4VMultipleScattering.cc.

546{
548 return &fParticleChange;
549}
virtual void Initialize(const G4Track &)

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 447 of file G4VMultipleScattering.cc.

449{
451 return DBL_MAX;
452}
@ NotForced

◆ PreparePhysicsTable()

void G4VMultipleScattering::PreparePhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 170 of file G4VMultipleScattering.cc.

171{
172 if(1 < verboseLevel) {
173 G4cout << "### G4VMultipleScattering::PrepearPhysicsTable() for "
174 << GetProcessName()
175 << " and particle " << part.GetParticleName()
176 << G4endl;
177 }
178 G4bool master = emManager->IsMaster();
179
180 if(!firstParticle) { firstParticle = &part; }
181 if(part.GetParticleType() == "nucleus") {
182 stepLimit = fMinimal;
183 latDisplacement = false;
184 facrange = 0.2;
185 G4String pname = part.GetParticleName();
186 if(pname != "deuteron" && pname != "triton" &&
187 pname != "alpha+" && pname != "helium" &&
188 pname != "alpha" && pname != "He3" &&
189 pname != "hydrogen") {
190
191 const G4ParticleDefinition* theGenericIon =
193 if(&part == theGenericIon) { isIon = true; }
194
195 if(theGenericIon && firstParticle != theGenericIon) {
196 G4ProcessManager* pm = theGenericIon->GetProcessManager();
198 size_t n = v->size();
199 for(size_t j=0; j<n; ++j) {
200 if((*v)[j] == this) {
201 firstParticle = theGenericIon;
202 isIon = true;
203 break;
204 }
205 }
206 }
207 }
208 }
209
210 emManager->PreparePhysicsTable(&part, this, master);
211 currParticle = nullptr;
212
213 if(1 < verboseLevel) {
214 G4cout << "### G4VMultipleScattering::PrepearPhysicsTable() for "
215 << GetProcessName()
216 << " and particle " << part.GetParticleName()
217 << " local particle " << firstParticle->GetParticleName()
218 << " isIon: " << isIon << " isMaster: " << master
219 << G4endl;
220 }
221
222 if(firstParticle == &part) {
223
224 // initialise process
225 InitialiseProcess(firstParticle);
226
227 // heavy particles and not ions
228 if(!isIon) {
229 if(part.GetPDGMass() > MeV) {
230 stepLimit = theParameters->MscMuHadStepLimitType();
231 facrange = theParameters->MscMuHadRangeFactor();
232 latDisplacement = theParameters->MuHadLateralDisplacement();
233 } else {
234 stepLimit = theParameters->MscStepLimitType();
235 facrange = theParameters->MscRangeFactor();
236 latDisplacement = theParameters->LateralDisplacement();
237 }
238 }
239 if(master) { SetVerboseLevel(theParameters->Verbose()); }
240 else { SetVerboseLevel(theParameters->WorkerVerbose()); }
241
242 // initialisation of models
243 numberOfModels = modelManager->NumberOfModels();
244 /*
245 G4cout << "### G4VMultipleScattering::PreparePhysicsTable() for "
246 << GetProcessName()
247 << " and particle " << part.GetParticleName()
248 << " Nmod= " << numberOfModels << " " << this
249 << G4endl;
250 */
251 for(G4int i=0; i<numberOfModels; ++i) {
252 G4VMscModel* msc = static_cast<G4VMscModel*>(GetModelByIndex(i));
253 if(!msc) { continue; }
254 msc->SetIonisation(nullptr, firstParticle);
255 msc->SetMasterThread(master);
256 currentModel = msc;
257 msc->SetPolarAngleLimit(theParameters->MscThetaLimit());
258 G4double emax =
259 std::min(msc->HighEnergyLimit(),theParameters->MaxKinEnergy());
260 msc->SetHighEnergyLimit(emax);
261 }
262
263 modelManager->Initialise(firstParticle, G4Electron::Electron(),
264 10.0, verboseLevel);
265
266 if(!safetyHelper) {
268 ->GetSafetyHelper();
269 safetyHelper->InitialiseHelper();
270 }
271 }
272}
@ fMinimal
static G4Electron * Electron()
Definition: G4Electron.cc:93
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4double minSubRange, G4int verb)
G4double MscMuHadRangeFactor() const
G4MscStepLimitType MscMuHadStepLimitType() const
G4int Verbose() const
G4int WorkerVerbose() const
G4MscStepLimitType MscStepLimitType() const
G4double MaxKinEnergy() const
G4bool LateralDisplacement() const
G4bool MuHadLateralDisplacement() const
G4double MscRangeFactor() const
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleType() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
void InitialiseHelper()
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:792
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:729
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
Definition: G4VMscModel.hh:395
virtual void InitialiseProcess(const G4ParticleDefinition *)=0

◆ PrintInfo()

virtual void G4VMultipleScattering::PrintInfo ( )
inlinevirtual

Reimplemented in G4AdjointhMultipleScattering, and G4eAdjointMultipleScattering.

Definition at line 107 of file G4VMultipleScattering.hh.

107{};

◆ ProcessDescription()

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

Reimplemented from G4VProcess.

Reimplemented in G4MuMultipleScattering, G4eMultipleScattering, and G4hMultipleScattering.

Definition at line 654 of file G4VMultipleScattering.cc.

655{
656 if(firstParticle) {
657 StreamInfo(outFile, *firstParticle, true);
658 }
659}

Referenced by G4MuMultipleScattering::ProcessDescription(), G4eMultipleScattering::ProcessDescription(), and G4hMultipleScattering::ProcessDescription().

◆ RangeFactor()

◆ RetrievePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 635 of file G4VMultipleScattering.cc.

638{
639 return true;
640}

◆ SelectModel()

G4VEmModel * G4VMultipleScattering::SelectModel ( G4double  kinEnergy,
size_t  idx 
)
inline

Definition at line 310 of file G4VMultipleScattering.hh.

311{
312 return modelManager->SelectModel(kinEnergy, coupleIndex);
313}
G4VEmModel * SelectModel(G4double &energy, size_t &index)

Referenced by AlongStepGetPhysicalInteractionLength().

◆ SetEmModel()

◆ SetIonisation()

void G4VMultipleScattering::SetIonisation ( G4VEnergyLossProcess p)

Definition at line 644 of file G4VMultipleScattering.cc.

645{
646 for(G4int i=0; i<numberOfModels; ++i) {
647 G4VMscModel* msc = static_cast<G4VMscModel*>(GetModelByIndex(i, true));
648 if(msc) { msc->SetIonisation(p, firstParticle); }
649 }
650}

◆ SetLateralDisplasmentFlag()

void G4VMultipleScattering::SetLateralDisplasmentFlag ( G4bool  val)
inline

Definition at line 324 of file G4VMultipleScattering.hh.

325{
326 latDisplacement = val;
327}

◆ SetLowestKinEnergy()

void G4VMultipleScattering::SetLowestKinEnergy ( G4double  val)
inline

Definition at line 399 of file G4VMultipleScattering.hh.

400{
401 lowestKinEnergy = val;
402}

◆ SetRangeFactor()

void G4VMultipleScattering::SetRangeFactor ( G4double  val)
inline

Definition at line 352 of file G4VMultipleScattering.hh.

353{
354 if(val > 0.0 && val < 1.0) {
355 facrange = val;
356 theParameters->SetMscRangeFactor(val);
357 }
358}
void SetMscRangeFactor(G4double val)

Referenced by SetStepLimitType().

◆ SetSkin()

void G4VMultipleScattering::SetSkin ( G4double  val)
inline

Definition at line 338 of file G4VMultipleScattering.hh.

339{
340 theParameters->SetMscSkin(val);
341}
void SetMscSkin(G4double val)

◆ SetStepLimitType()

◆ Skin()

G4double G4VMultipleScattering::Skin ( ) const
inline

◆ StartTracking()

void G4VMultipleScattering::StartTracking ( G4Track track)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 357 of file G4VMultipleScattering.cc.

358{
359 G4VEnergyLossProcess* eloss = nullptr;
360 if(track->GetParticleDefinition() != currParticle) {
361 currParticle = track->GetParticleDefinition();
362 fIonisation = emManager->GetEnergyLossProcess(currParticle);
363 eloss = fIonisation;
364 }
365 /*
366 G4cout << "G4VMultipleScattering::StartTracking Nmod= " << numberOfModels
367 << " " << currParticle->GetParticleName()
368 << " E(MeV)= " << track->GetKineticEnergy()
369 << " Ion= " << eloss << " " << fIonisation << " IsMaster= "
370 << G4LossTableManager::Instance()->IsMaster()
371 << G4endl;
372 */
373 for(G4int i=0; i<numberOfModels; ++i) {
374 /*
375 G4cout << "Next model " << i << " " << msc
376 << " Emin= " << msc->LowEnergyLimit()
377 << " Emax= " << msc->HighEnergyLimit()
378 << " Eact= " << msc->LowEnergyActivationLimit() << G4endl;
379 */
380 G4VEmModel* msc = GetModelByIndex(i);
381 msc->StartTracking(track);
382 if(eloss) {
383 G4VMscModel* mscmod = static_cast<G4VMscModel*>(msc);
384 if(mscmod) { mscmod->SetIonisation(fIonisation, currParticle); }
385 }
386 }
387}
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:288

Referenced by G4eAdjointMultipleScattering::StartTracking().

◆ StepLimitType()

◆ StorePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 591 of file G4VMultipleScattering.cc.

594{
595 G4bool yes = true;
596 if(part != firstParticle) { return yes; }
597 const G4VMultipleScattering* masterProcess =
598 static_cast<const G4VMultipleScattering*>(GetMasterProcess());
599 if(masterProcess && masterProcess != this) { return yes; }
600
601 G4int nmod = modelManager->NumberOfModels();
602 static const G4String ss[4] = {"1","2","3","4"};
603 for(G4int i=0; i<nmod; ++i) {
604 G4VEmModel* msc = modelManager->GetModel(i);
605 yes = true;
606 G4PhysicsTable* table = msc->GetCrossSectionTable();
607 if (table) {
608 G4int j = std::min(i,3);
609 G4String name =
610 GetPhysicsTableFileName(part,directory,"LambdaMod"+ss[j],ascii);
611 yes = table->StorePhysicsTable(name,ascii);
612
613 if ( yes ) {
614 if ( verboseLevel>0 ) {
615 G4cout << "Physics table are stored for "
616 << part->GetParticleName()
617 << " and process " << GetProcessName()
618 << " with a name <" << name << "> " << G4endl;
619 }
620 } else {
621 G4cout << "Fail to store Physics Table for "
622 << part->GetParticleName()
623 << " and process " << GetProcessName()
624 << " in the directory <" << directory
625 << "> " << G4endl;
626 }
627 }
628 }
629 return yes;
630}
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:182
const char * name(G4int ptype)

◆ StreamProcessInfo()

virtual void G4VMultipleScattering::StreamProcessInfo ( std::ostream &  ) const
inlineprotectedvirtual

Reimplemented in G4MuMultipleScattering, G4eMultipleScattering, and G4hMultipleScattering.

Definition at line 115 of file G4VMultipleScattering.hh.

115{};

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForMSC G4VMultipleScattering::fParticleChange
protected

Definition at line 285 of file G4VMultipleScattering.hh.

Referenced by AlongStepDoIt(), G4VMultipleScattering(), and PostStepDoIt().


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