Geant4 9.6.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)=0
 
virtual void PrintInfo ()=0
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
void PrintInfoDefinition ()
 
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
 
void StartTracking (G4Track *)
 
G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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=0)
 
void SetModel (G4VMscModel *, G4int index=1)
 
G4VMscModelModel (G4int index=1)
 
void SetEmModel (G4VMscModel *, G4int index=1)
 
G4VMscModelEmModel (G4int index=1)
 
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
 
void SetGeomFactor (G4double val)
 
G4double PolarAngleLimit () const
 
void SetPolarAngleLimit (G4double val)
 
G4MscStepLimitType StepLimitType () const
 
void SetStepLimitType (G4MscStepLimitType val)
 
const G4ParticleDefinitionFirstParticle () const
 
- Public Member Functions inherited from G4VContinuousDiscreteProcess
 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
 
virtual ~G4VContinuousDiscreteProcess ()
 
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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Protected Member Functions

virtual void InitialiseProcess (const G4ParticleDefinition *)=0
 
G4double GetMeanFreePath (const G4Track &track, G4double, G4ForceCondition *condition)
 
G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
 
- 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 previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4GPILSelection valueGPILSelectionMSC
 
G4ParticleChangeForMSC fParticleChange
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

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 86 of file G4VMultipleScattering.cc.

86 :
88 numberOfModels(0),
89 firstParticle(0),
90 currParticle(0),
91 stepLimit(fUseSafety),
92 skin(1.0),
93 facrange(0.04),
94 facgeom(2.5),
95 latDisplasment(true),
96 isIon(false)
97{
100 if("ionmsc" == name) { firstParticle = G4GenericIon::GenericIon(); }
101
102 geomMin = 1.e-6*CLHEP::mm;
103 lowestKinEnergy = 1*eV;
104
105 // default limit on polar angle
106 polarAngleLimit = 0.0;
107
108 physStepLimit = gPathLength = tPathLength = 0.0;
109 fIonisation = 0;
110
112 safetyHelper = 0;
113 fPositionChanged = false;
114 isActive = false;
115
116 modelManager = new G4EmModelManager();
117 emManager = G4LossTableManager::Instance();
118 emManager->Register(this);
119
120 warn = 0;
121}
@ fMultipleScattering
@ fUseSafety
@ fElectromagnetic
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
static G4LossTableManager * Instance()
void Register(G4VEnergyLossProcess *p)
G4ParticleChangeForMSC fParticleChange
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:408
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283

◆ ~G4VMultipleScattering()

G4VMultipleScattering::~G4VMultipleScattering ( )
virtual

Definition at line 125 of file G4VMultipleScattering.cc.

126{
127 if(1 < verboseLevel) {
128 G4cout << "G4VMultipleScattering destruct " << GetProcessName()
129 << G4endl;
130 }
131 delete modelManager;
132 emManager->DeRegister(this);
133}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void DeRegister(G4VEnergyLossProcess *p)
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379

Member Function Documentation

◆ AddEmModel()

◆ AlongStepDoIt()

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 403 of file G4VMultipleScattering.cc.

404{
406 fNewPosition = step.GetPostStepPoint()->GetPosition();
407 fParticleChange.ProposePosition(fNewPosition);
408 fPositionChanged = false;
409
410 G4double geomLength = step.GetStepLength();
411
412 // very small step - no msc
413 if(!isActive) {
414 tPathLength = geomLength;
415
416 // sample msc
417 } else {
418 G4double range =
419 currentModel->GetRange(currParticle,track.GetKineticEnergy(),
420 track.GetMaterialCutsCouple());
421
422 G4double trueLength = currentModel->ComputeTrueStepLength(geomLength);
423
424
425 // protection against wrong t->g->t conversion
426 // if(trueLength > tPathLength)
427 /*
428 if(currParticle->GetPDGMass() > GeV)
429 G4cout << "G4VMsc::AlongStepDoIt: GeomLength= "
430 << geomLength
431 << " trueLenght= " << trueLength
432 << " tPathLength= " << tPathLength
433 << " dr= " << range - trueLength
434 << " ekin= " << track.GetKineticEnergy() << G4endl;
435 */
436 if (trueLength <= physStepLimit) {
437 tPathLength = trueLength;
438 } else {
439 tPathLength = physStepLimit - 0.5*geomMin;
440 }
441
442 // do not sample scattering at the last or at a small step
443 if(tPathLength + geomMin < range && tPathLength > geomMin) {
444
445 G4double preSafety = step.GetPreStepPoint()->GetSafety();
446 G4double postSafety= preSafety - geomLength;
447 G4bool safetyRecomputed = false;
448 if( postSafety < geomMin ) {
449 safetyRecomputed = true;
450 postSafety = currentModel->ComputeSafety(fNewPosition,0.0);
451 }
452 G4ThreeVector displacement =
453 currentModel->SampleScattering(track.GetDynamicParticle(),postSafety);
454
455 G4double r2 = displacement.mag2();
456
457 //G4cout << "R= " << sqrt(r2) << " postSafety= " << postSafety << G4endl;
458
459 // make correction for displacement
460 if(r2 > 0.0) {
461
462 fPositionChanged = true;
463 G4double fac = 1.0;
464
465 // displaced point is definitely within the volume
466 if(r2 > postSafety*postSafety) {
467 if(!safetyRecomputed) {
468 postSafety = currentModel->ComputeSafety(fNewPosition, 0.0);
469 }
470 // add a factor which ensure numerical stability
471 if(r2 > postSafety*postSafety) { fac = 0.99*postSafety/std::sqrt(r2); }
472 }
473 // compute new endpoint of the Step
474 fNewPosition += fac*displacement;
475 //safetyHelper->ReLocateWithinVolume(fNewPosition);
476 }
477 }
478 }
480 //fParticleChange.ProposePosition(fNewPosition);
481 return &fParticleChange;
482}
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
double mag2() const
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void ProposePosition(const G4ThreeVector &finalPosition)
G4double GetSafety() const
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4DynamicParticle * GetDynamicParticle() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual G4ThreeVector & SampleScattering(const G4DynamicParticle *, G4double safety)
Definition: G4VMscModel.cc:131
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:288
G4double ComputeSafety(const G4ThreeVector &position, G4double limit)
Definition: G4VMscModel.hh:238
virtual G4double ComputeTrueStepLength(G4double geomPathLength)
Definition: G4VMscModel.cc:152
void ProposeTrueStepLength(G4double truePathLength)

◆ AlongStepGetPhysicalInteractionLength()

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 346 of file G4VMultipleScattering.cc.

352{
353 // get Step limit proposed by the process
354 *selection = NotCandidateForSelection;
355 physStepLimit = gPathLength = tPathLength = currentMinimalStep;
356
357 G4double ekin = track.GetKineticEnergy();
358 // isIon flag is used only to select a model
359 if(isIon) {
360 ekin *= proton_mass_c2/track.GetParticleDefinition()->GetPDGMass();
361 }
362
363 // select new model
364 if(1 < numberOfModels) {
365 currentModel = static_cast<G4VMscModel*>(
367 }
368
369 // step limit
370 if(currentModel->IsActive(ekin) && gPathLength >= geomMin
371 && ekin >= lowestKinEnergy) {
372 isActive = true;
373 tPathLength = currentModel->ComputeTruePathLengthLimit(track, gPathLength);
374 if (tPathLength < physStepLimit) {
375 *selection = CandidateForSelection;
376 }
377 } else { isActive = false; }
378 /*
379 if(currParticle->GetPDGMass() > GeV)
380 G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
381 << " gPathLength= " << gPathLength
382 << " tPathLength= " << tPathLength
383 << " currentMinimalStep= " << currentMinimalStep
384 << " isActive " << isActive << G4endl;
385 */
386 return gPathLength;
387}
@ CandidateForSelection
@ NotCandidateForSelection
const G4ParticleDefinition * GetParticleDefinition() const
G4bool IsActive(G4double kinEnergy)
Definition: G4VEmModel.hh:613
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)
Definition: G4VMscModel.cc:138
G4VEmModel * SelectModel(G4double kinEnergy, size_t idx)

Referenced by GetContinuousStepLimit().

◆ BuildPhysicsTable()

void G4VMultipleScattering::BuildPhysicsTable ( const G4ParticleDefinition part)
virtual

Reimplemented from G4VProcess.

Definition at line 267 of file G4VMultipleScattering.cc.

268{
269 G4String num = part.GetParticleName();
270 if(1 < verboseLevel) {
271 G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
272 << GetProcessName()
273 << " and particle " << num
274 << G4endl;
275 }
276
277 emManager->BuildPhysicsTable(firstParticle);
278
279 // explicitly defined printout by particle name
280 if(1 < verboseLevel ||
281 (0 < verboseLevel && (num == "e-" ||
282 num == "e+" || num == "mu+" ||
283 num == "mu-" || num == "proton"||
284 num == "pi+" || num == "pi-" ||
285 num == "kaon+" || num == "kaon-" ||
286 num == "alpha" || num == "anti_proton" ||
287 num == "GenericIon")))
288 {
290 << ": for " << num
291 << " SubType= " << GetProcessSubType()
292 << G4endl;
293 PrintInfo();
294 modelManager->DumpModelList(verboseLevel);
295 }
296
297 if(1 < verboseLevel) {
298 G4cout << "### G4VMultipleScattering::BuildPhysicsTable() done for "
299 << GetProcessName()
300 << " and particle " << num
301 << G4endl;
302 }
303}
void DumpModelList(G4int verb)
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
const G4String & GetParticleName() const
virtual void PrintInfo()=0
G4int GetProcessSubType() const
Definition: G4VProcess.hh:397

◆ ContinuousStepLimit()

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

Definition at line 514 of file G4VMultipleScattering.cc.

519{
520 return GetContinuousStepLimit(track,previousStepSize,currentMinimalStep,
521 currentSafety);
522}
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)

◆ EmModel()

G4VMscModel * G4VMultipleScattering::EmModel ( G4int  index = 1)

Definition at line 186 of file G4VMultipleScattering.cc.

187{
188 G4VMscModel* p = 0;
189 if(index >= 0 && index < G4int(mscModels.size())) { p = mscModels[index]; }
190 return p;
191}
int G4int
Definition: G4Types.hh:66

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

◆ FirstParticle()

const G4ParticleDefinition * G4VMultipleScattering::FirstParticle ( ) const
inline

Definition at line 402 of file G4VMultipleScattering.hh.

403{
404 return firstParticle;
405}

◆ GeomFactor()

G4double G4VMultipleScattering::GeomFactor ( ) const
inline

◆ GetContinuousStepLimit()

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

Implements G4VContinuousDiscreteProcess.

Definition at line 499 of file G4VMultipleScattering.cc.

504{
506 G4double x = AlongStepGetPhysicalInteractionLength(track,previousStepSize,
507 currentMinimalStep,
508 currentSafety, &selection);
509 return x;
510}
G4GPILSelection
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection)

Referenced by ContinuousStepLimit().

◆ GetMeanFreePath()

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

Implements G4VContinuousDiscreteProcess.

Definition at line 526 of file G4VMultipleScattering.cc.

528{
529 *condition = Forced;
530 return DBL_MAX;
531}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced
#define DBL_MAX
Definition: templates.hh:83

◆ GetModelByIndex()

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

Definition at line 196 of file G4VMultipleScattering.cc.

197{
198 return modelManager->GetModel(idx, ver);
199}
G4VEmModel * GetModel(G4int, G4bool ver=false)

◆ InitialiseProcess()

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

◆ IsApplicable()

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

◆ LateralDisplasmentFlag()

◆ Model()

G4VMscModel * G4VMultipleScattering::Model ( G4int  index = 1)

Definition at line 162 of file G4VMultipleScattering.cc.

163{
164 ++warn;
165 if(warn < 10) {
166 G4cout << "### G4VMultipleScattering::Model is obsolete method "
167 << "and will be removed for the next release." << G4endl;
168 G4cout << " Please, use EmModel" << G4endl;
169 }
170 G4VMscModel* p = 0;
171 if(index >= 0 && index < G4int(mscModels.size())) { p = mscModels[index]; }
172 return p;
173}

◆ PolarAngleLimit()

G4double G4VMultipleScattering::PolarAngleLimit ( ) const
inline

Definition at line 371 of file G4VMultipleScattering.hh.

372{
373 return polarAngleLimit;
374}

Referenced by G4MuMultipleScattering::PrintInfo().

◆ PostStepDoIt()

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 487 of file G4VMultipleScattering.cc.

488{
490 if(fPositionChanged) {
491 safetyHelper->ReLocateWithinVolume(fNewPosition);
492 fParticleChange.ProposePosition(fNewPosition);
493 }
494 return &fParticleChange;
495}
virtual void Initialize(const G4Track &)
void ReLocateWithinVolume(const G4ThreeVector &pGlobalPoint)

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4VContinuousDiscreteProcess.

Definition at line 392 of file G4VMultipleScattering.cc.

394{
395 *condition = Forced;
396 //*condition = NotForced;
397 return DBL_MAX;
398}

◆ PreparePhysicsTable()

void G4VMultipleScattering::PreparePhysicsTable ( const G4ParticleDefinition part)
virtual

Reimplemented from G4VProcess.

Definition at line 204 of file G4VMultipleScattering.cc.

205{
206 if(!firstParticle) { firstParticle = &part; }
207 if(part.GetParticleType() == "nucleus") {
210 SetRangeFactor(0.2);
211 if(&part == G4GenericIon::GenericIon()) { firstParticle = &part; }
212 isIon = true;
213 }
214
215 emManager->PreparePhysicsTable(&part, this);
216 currParticle = 0;
217
218 if(1 < verboseLevel) {
219 G4cout << "### G4VMultipleScattering::PrepearPhysicsTable() for "
220 << GetProcessName()
221 << " and particle " << part.GetParticleName()
222 << " local particle " << firstParticle->GetParticleName()
223 << " isIon= " << isIon
224 << G4endl;
225 }
226
227 if(firstParticle == &part) {
228
229 InitialiseProcess(firstParticle);
230
231 // initialisation of models
232 numberOfModels = modelManager->NumberOfModels();
233 for(G4int i=0; i<numberOfModels; ++i) {
234 G4VMscModel* msc = static_cast<G4VMscModel*>(modelManager->GetModel(i));
235 msc->SetIonisation(0, firstParticle);
236 if(0 == i) { currentModel = msc; }
237 if(isIon) {
239 msc->SetLateralDisplasmentFlag(false);
240 msc->SetRangeFactor(0.2);
241 } else {
244 msc->SetSkin(Skin());
247 }
248 msc->SetPolarAngleLimit(polarAngleLimit);
249 G4double emax =
250 std::min(msc->HighEnergyLimit(),emManager->MaxKinEnergy());
251 msc->SetHighEnergyLimit(emax);
252 }
253
254 modelManager->Initialise(firstParticle, G4Electron::Electron(),
255 10.0, verboseLevel);
256
257 if(!safetyHelper) {
259 ->GetSafetyHelper();
260 safetyHelper->InitialiseHelper();
261 }
262 }
263}
@ fMinimal
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4int NumberOfModels() const
const G4DataVector * Initialise(const G4ParticleDefinition *, const G4ParticleDefinition *, G4double, G4int)
G4double MaxKinEnergy() const
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
const G4String & GetParticleType() const
void InitialiseHelper()
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:620
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void SetSkin(G4double)
Definition: G4VMscModel.hh:203
void SetRangeFactor(G4double)
Definition: G4VMscModel.hh:210
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
Definition: G4VMscModel.hh:324
void SetLateralDisplasmentFlag(G4bool val)
Definition: G4VMscModel.hh:196
void SetGeomFactor(G4double)
Definition: G4VMscModel.hh:217
void SetStepLimitType(G4MscStepLimitType)
Definition: G4VMscModel.hh:224
void SetRangeFactor(G4double val)
void SetStepLimitType(G4MscStepLimitType val)
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
void SetLateralDisplasmentFlag(G4bool val)
G4bool LateralDisplasmentFlag() const
G4MscStepLimitType StepLimitType() const

◆ PrintInfo()

virtual void G4VMultipleScattering::PrintInfo ( )
pure virtual

◆ PrintInfoDefinition()

void G4VMultipleScattering::PrintInfoDefinition ( )

Definition at line 307 of file G4VMultipleScattering.cc.

308{
309 if (0 < verboseLevel) {
311 << ": for " << firstParticle->GetParticleName()
312 << " SubType= " << GetProcessSubType()
313 << G4endl;
314 PrintInfo();
315 modelManager->DumpModelList(verboseLevel);
316 }
317}

◆ RangeFactor()

◆ RetrievePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 574 of file G4VMultipleScattering.cc.

577{
578 return true;
579}

◆ SelectModel()

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

Definition at line 307 of file G4VMultipleScattering.hh.

308{
309 return modelManager->SelectModel(kinEnergy, coupleIndex);
310}
G4VEmModel * SelectModel(G4double &energy, size_t &index)

Referenced by AlongStepGetPhysicalInteractionLength().

◆ SetEmModel()

void G4VMultipleScattering::SetEmModel ( G4VMscModel p,
G4int  index = 1 
)

Definition at line 177 of file G4VMultipleScattering.cc.

178{
179 G4int n = mscModels.size();
180 if(index >= n) { for(G4int i=n; i<=index; ++i) { mscModels.push_back(0); } }
181 mscModels[index] = p;
182}

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

◆ SetGeomFactor()

void G4VMultipleScattering::SetGeomFactor ( G4double  val)
inline

Definition at line 364 of file G4VMultipleScattering.hh.

365{
366 if(val > 0.0) facgeom = val;
367}

◆ SetIonisation()

void G4VMultipleScattering::SetIonisation ( G4VEnergyLossProcess p)

Definition at line 583 of file G4VMultipleScattering.cc.

584{
585 for(G4int i=0; i<numberOfModels; ++i) {
586 G4VMscModel* msc = static_cast<G4VMscModel*>(modelManager->GetModel(i));
587 msc->SetIonisation(p, firstParticle);
588 }
589}

◆ SetLateralDisplasmentFlag()

void G4VMultipleScattering::SetLateralDisplasmentFlag ( G4bool  val)
inline

Definition at line 321 of file G4VMultipleScattering.hh.

322{
323 latDisplasment = val;
324}

Referenced by G4AdjointhMultipleScattering::InitialiseProcess(), and PreparePhysicsTable().

◆ SetModel()

void G4VMultipleScattering::SetModel ( G4VMscModel p,
G4int  index = 1 
)

Definition at line 147 of file G4VMultipleScattering.cc.

148{
149 ++warn;
150 if(warn < 10) {
151 G4cout << "### G4VMultipleScattering::SetModel is obsolete method "
152 << "and will be removed for the next release." << G4endl;
153 G4cout << " Please, use SetEmModel" << G4endl;
154 }
155 G4int n = mscModels.size();
156 if(index >= n) { for(G4int i=n; i<=index; ++i) {mscModels.push_back(0);} }
157 mscModels[index] = p;
158}

◆ SetPolarAngleLimit()

void G4VMultipleScattering::SetPolarAngleLimit ( G4double  val)
inline

Definition at line 378 of file G4VMultipleScattering.hh.

379{
380 if(val < 0.0) { polarAngleLimit = 0.0; }
381 else if(val > CLHEP::pi) { polarAngleLimit = CLHEP::pi; }
382 else { polarAngleLimit = val; }
383}

◆ SetRangeFactor()

void G4VMultipleScattering::SetRangeFactor ( G4double  val)
inline

Definition at line 350 of file G4VMultipleScattering.hh.

351{
352 if(val > 0.0) facrange = val;
353}

Referenced by PreparePhysicsTable().

◆ SetSkin()

void G4VMultipleScattering::SetSkin ( G4double  val)
inline

Definition at line 335 of file G4VMultipleScattering.hh.

336{
337 if(val < 1.0) { skin = 0.0; }
338 else { skin = val; }
339}

◆ SetStepLimitType()

◆ Skin()

◆ StartTracking()

void G4VMultipleScattering::StartTracking ( G4Track track)
virtual

Reimplemented from G4VProcess.

Definition at line 321 of file G4VMultipleScattering.cc.

322{
323 G4VEnergyLossProcess* eloss = 0;
324 if(track->GetParticleDefinition() != currParticle) {
325 currParticle = track->GetParticleDefinition();
326 fIonisation = emManager->GetEnergyLossProcess(currParticle);
327 eloss = fIonisation;
328 }
329 // one model
330 if(1 == numberOfModels) {
331 currentModel->StartTracking(track);
332 if(eloss) { currentModel->SetIonisation(fIonisation, currParticle); }
333
334 // many models
335 } else {
336 for(G4int i=0; i<numberOfModels; ++i) {
337 G4VMscModel* msc = static_cast<G4VMscModel*>(modelManager->GetModel(i));
338 msc->StartTracking(track);
339 if(eloss) { msc->SetIonisation(fIonisation, currParticle); }
340 }
341 }
342}
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:211

◆ StepLimitType()

◆ StorePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 536 of file G4VMultipleScattering.cc.

539{
540 G4bool yes = true;
541 if(part != firstParticle) { return yes; }
542 G4int nmod = modelManager->NumberOfModels();
543 const G4String ss[4] = {"1","2","3","4"};
544 for(G4int i=0; i<nmod; ++i) {
545 G4VEmModel* msc = modelManager->GetModel(i);
546 yes = true;
547 G4PhysicsTable* table = msc->GetCrossSectionTable();
548 if (table) {
549 G4int j = std::min(i,3);
550 G4String name =
551 GetPhysicsTableFileName(part,directory,"LambdaMod"+ss[j],ascii);
552 yes = table->StorePhysicsTable(name,ascii);
553
554 if ( yes ) {
555 if ( verboseLevel>0 ) {
556 G4cout << "Physics table are stored for " << part->GetParticleName()
557 << " and process " << GetProcessName()
558 << " with a name <" << name << "> " << G4endl;
559 }
560 } else {
561 G4cout << "Fail to store Physics Table for " << part->GetParticleName()
562 << " and process " << GetProcessName()
563 << " in the directory <" << directory
564 << "> " << G4endl;
565 }
566 }
567 }
568 return yes;
569}
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:660
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:214

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForMSC G4VMultipleScattering::fParticleChange
protected

Definition at line 283 of file G4VMultipleScattering.hh.

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

◆ valueGPILSelectionMSC

G4GPILSelection G4VMultipleScattering::valueGPILSelectionMSC
protected

Definition at line 282 of file G4VMultipleScattering.hh.


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