Geant4  10.00.p02
G4VMultipleScattering.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4VMultipleScattering.cc 81864 2014-06-06 11:30:54Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4VMultipleScattering
34 //
35 // Author: Vladimir Ivanchenko on base of Laszlo Urban code
36 //
37 // Creation date: 25.03.2003
38 //
39 // Modifications:
40 //
41 // 13.04.03 Change printout (V.Ivanchenko)
42 // 04-06-03 Fix compilation warnings (V.Ivanchenko)
43 // 16-07-03 Use G4VMscModel interface (V.Ivanchenko)
44 // 03-11-03 Fix initialisation problem in RetrievePhysicsTable (V.Ivanchenko)
45 // 04-11-03 Update PrintInfoDefinition (V.Ivanchenko)
46 // 01-03-04 SampleCosineTheta signature changed
47 // 22-04-04 SampleCosineTheta signature changed back to original
48 // 27-08-04 Add InitialiseForRun method (V.Ivanchneko)
49 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
50 // 11-03-05 Shift verbose level by 1 (V.Ivantchenko)
51 // 15-04-05 optimize internal interface (V.Ivanchenko)
52 // 15-04-05 remove boundary flag (V.Ivanchenko)
53 // 27-10-05 introduce virtual function MscStepLimitation() (V.Ivanchenko)
54 // 12-04-07 Add verbosity at destruction (V.Ivanchenko)
55 // 27-10-07 Virtual functions moved to source (V.Ivanchenko)
56 // 11-03-08 Set skin value does not effect step limit type (V.Ivanchenko)
57 // 24-06-09 Removed hidden bin in G4PhysicsVector (V.Ivanchenko)
58 // 04-06-13 Adoptation to MT mode (V.Ivanchenko)
59 //
60 // Class Description:
61 //
62 // It is the generic process of multiple scattering it includes common
63 // part of calculations for all charged particles
64 
65 // -------------------------------------------------------------------
66 //
67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
68 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
69 
70 #include "G4VMultipleScattering.hh"
71 #include "G4PhysicalConstants.hh"
72 #include "G4SystemOfUnits.hh"
73 #include "G4LossTableManager.hh"
74 #include "G4MaterialCutsCouple.hh"
75 #include "G4Step.hh"
76 #include "G4ParticleDefinition.hh"
77 #include "G4VEmFluctuationModel.hh"
78 #include "G4UnitsTable.hh"
79 #include "G4ProductionCutsTable.hh"
80 #include "G4Electron.hh"
81 #include "G4GenericIon.hh"
83 #include "G4SafetyHelper.hh"
84 #include "G4ParticleTable.hh"
85 #include "G4ProcessVector.hh"
86 #include "G4ProcessManager.hh"
87 
88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89 
93  numberOfModels(0),
94  firstParticle(0),
95  currParticle(0),
96  stepLimit(fUseSafety),
97  skin(1.0),
98  facrange(0.04),
99  facgeom(2.5),
100  latDisplasment(true),
101  isIon(false)
102 {
103  SetVerboseLevel(1);
105  if("ionmsc" == name) { firstParticle = G4GenericIon::GenericIon(); }
106 
107  geomMin = 0.05*CLHEP::nm;
108  lowestKinEnergy = 10*eV;
109 
110  // default limit on polar angle
111  polarAngleLimit = 0.0;
112 
114  fIonisation = 0;
115 
117  safetyHelper = 0;
118  fPositionChanged = false;
119  isActive = false;
120 
121  currentModel = 0;
124  emManager->Register(this);
125 
126  warn = 0;
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
130 
132 {
133  /*
134  if(1 < verboseLevel) {
135  G4cout << "G4VMultipleScattering destruct " << GetProcessName()
136  << G4endl;
137  }
138  */
139  delete modelManager;
140  emManager->DeRegister(this);
141 }
142 
143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
144 
146  const G4Region* region)
147 {
149  modelManager->AddEmModel(order, p, fm, region);
150  if(p) { p->SetParticleChange(pParticleChange); }
151 }
152 
153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
154 
156 {
157  G4int n = mscModels.size();
158  if(index >= n) { for(G4int i=n; i<=index; ++i) { mscModels.push_back(0); } }
159  mscModels[index] = p;
160 }
161 
162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
163 
165 {
166  G4VMscModel* p = 0;
167  if(index >= 0 && index < G4int(mscModels.size())) { p = mscModels[index]; }
168  return p;
169 }
170 
171 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
172 
173 G4VEmModel*
175 {
176  return modelManager->GetModel(idx, ver);
177 }
178 
179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180 
181 void
183 {
184  G4bool master = true;
185  const G4VMultipleScattering* masterProc =
186  static_cast<const G4VMultipleScattering*>(GetMasterProcess());
187  if(masterProc && masterProc != this) { master = false; }
188 
189  if(!firstParticle) { firstParticle = &part; }
190  if(part.GetParticleType() == "nucleus") {
193  SetRangeFactor(0.2);
194  G4String pname = part.GetParticleName();
195  if(pname != "deuteron" && pname != "triton" &&
196  pname != "alpha+" && pname != "helium" &&
197  pname != "alpha" && pname != "He3" &&
198  pname != "hydrogen") {
199 
200  const G4ParticleDefinition* theGenericIon =
202 
203  if(theGenericIon && firstParticle != theGenericIon) {
204  G4ProcessManager* pm = theGenericIon->GetProcessManager();
206  size_t n = v->size();
207  for(size_t j=0; j<n; ++j) {
208  if((*v)[j] == this) {
209  firstParticle = theGenericIon;
210  isIon = true;
211  break;
212  }
213  }
214  }
215  }
216  }
217 
218  emManager->PreparePhysicsTable(&part, this, master);
219  currParticle = 0;
220 
221  if(1 < verboseLevel) {
222  G4cout << "### G4VMultipleScattering::PrepearPhysicsTable() for "
223  << GetProcessName()
224  << " and particle " << part.GetParticleName()
225  << " local particle " << firstParticle->GetParticleName()
226  << " isIon= " << isIon
227  << G4endl;
228  }
229 
230  if(firstParticle == &part) {
231 
233 
234  // initialisation of models
236  for(G4int i=0; i<numberOfModels; ++i) {
237  G4VMscModel* msc = static_cast<G4VMscModel*>(modelManager->GetModel(i));
238  msc->SetIonisation(0, firstParticle);
239  msc->SetMasterThread(master);
240  if(0 == i) { currentModel = msc; }
241  if(isIon) {
243  msc->SetLateralDisplasmentFlag(false);
244  msc->SetRangeFactor(0.2);
245  } else {
248  msc->SetSkin(Skin());
249  msc->SetRangeFactor(RangeFactor());
250  msc->SetGeomFactor(GeomFactor());
251  }
253  G4double emax =
255  msc->SetHighEnergyLimit(emax);
256  }
257 
259  10.0, verboseLevel);
260 
261  if(!safetyHelper) {
263  ->GetSafetyHelper();
265  }
266  }
267 }
268 
269 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
270 
272 {
273  G4String num = part.GetParticleName();
274  if(1 < verboseLevel) {
275  G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
276  << GetProcessName()
277  << " and particle " << num
278  << " IsMaster= " << G4LossTableManager::Instance()->IsMaster()
279  << G4endl;
280  }
281  G4bool master = true;
282  const G4VMultipleScattering* masterProcess =
283  static_cast<const G4VMultipleScattering*>(GetMasterProcess());
284  if(masterProcess && masterProcess != this) { master = false; }
285 
286  if(firstParticle == &part) {
287  /*
288  G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
289  << GetProcessName()
290  << " and particle " << num
291  << " IsMaster= " << G4LossTableManager::Instance()->IsMaster()
292  << " " << this
293  << G4endl;
294  */
296 
297  if(!master) {
298  // initialisation of models
299  G4bool printing = true;
301  /*
302  G4cout << "### G4VMultipleScattering::SlaveBuildPhysicsTable() for "
303  << GetProcessName()
304  << " and particle " << num
305  << " Nmod= " << numberOfModels << " " << this
306  << G4endl;
307  */
308  for(G4int i=0; i<numberOfModels; ++i) {
309  G4VMscModel* msc =
310  static_cast<G4VMscModel*>(GetModelByIndex(i, printing));
311  G4VMscModel* msc0=
312  static_cast<G4VMscModel*>(masterProcess->GetModelByIndex(i,printing));
313  msc->SetCrossSectionTable(msc0->GetCrossSectionTable(), false);
314  }
315  }
316 
317  }
318 
319  // explicitly defined printout by particle name
320  if(1 < verboseLevel ||
321  (0 < verboseLevel && (num == "e-" ||
322  num == "e+" || num == "mu+" ||
323  num == "mu-" || num == "proton"||
324  num == "pi+" || num == "pi-" ||
325  num == "kaon+" || num == "kaon-" ||
326  num == "alpha" || num == "anti_proton" ||
327  num == "GenericIon")))
328  {
329  G4cout << G4endl << GetProcessName()
330  << ": for " << num
331  << " SubType= " << GetProcessSubType()
332  << G4endl;
333  PrintInfo();
335  }
336 
337  if(1 < verboseLevel) {
338  G4cout << "### G4VMultipleScattering::BuildPhysicsTable() done for "
339  << GetProcessName()
340  << " and particle " << num
341  << G4endl;
342  }
343 }
344 
345 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
346 
348 {
349  if (0 < verboseLevel) {
350  G4cout << G4endl << GetProcessName()
351  << ": for " << firstParticle->GetParticleName()
352  << " SubType= " << GetProcessSubType()
353  << G4endl;
354  PrintInfo();
356  }
357 }
358 
359 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
360 
362 {
363  G4VEnergyLossProcess* eloss = 0;
364  if(track->GetParticleDefinition() != currParticle) {
367  eloss = fIonisation;
368  }
369  /*
370  G4cout << "G4VMultipleScattering::StartTracking Nmod= " << numberOfModels
371  << " " << currParticle->GetParticleName()
372  << " E(MeV)= " << track->GetKineticEnergy()
373  << " Ion= " << eloss << " " << fIonisation << " IsMaster= "
374  << G4LossTableManager::Instance()->IsMaster()
375  << G4endl;
376  */
377  // one model
378  if(1 == numberOfModels) {
379  currentModel->StartTracking(track);
381 
382  // many models
383  } else {
384  for(G4int i=0; i<numberOfModels; ++i) {
385  G4VMscModel* msc = static_cast<G4VMscModel*>(GetModelByIndex(i,true));
386  msc->StartTracking(track);
387  if(eloss) { msc->SetIonisation(fIonisation, currParticle); }
388  }
389  }
390 }
391 
392 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
393 
395  const G4Track& track,
396  G4double,
397  G4double currentMinimalStep,
398  G4double&,
399  G4GPILSelection* selection)
400 {
401  // get Step limit proposed by the process
402  *selection = NotCandidateForSelection;
403  physStepLimit = gPathLength = tPathLength = currentMinimalStep;
404 
405  G4double ekin = track.GetKineticEnergy();
406  /*
407  G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
408  << " " << currParticle->GetParticleName()
409  << " currMod " << currentModel
410  << G4endl;
411  */
412  // isIon flag is used only to select a model
413  if(isIon) {
414  ekin *= proton_mass_c2/track.GetParticleDefinition()->GetPDGMass();
415  }
416 
417  // select new model
418  if(1 < numberOfModels) {
419  currentModel = static_cast<G4VMscModel*>(
420  SelectModel(ekin,track.GetMaterialCutsCouple()->GetIndex()));
421  }
422  // step limit
423  if(currentModel->IsActive(ekin) && gPathLength >= geomMin
424  && ekin >= lowestKinEnergy) {
425  isActive = true;
426  tPathLength =
428  if (tPathLength < physStepLimit) {
429  *selection = CandidateForSelection;
430  }
431  } else { isActive = false; }
432 
433 
434  //if(currParticle->GetPDGMass() > GeV)
435  /*
436  G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
437  << " " << currParticle->GetParticleName()
438  << " gPathLength= " << gPathLength
439  << " tPathLength= " << tPathLength
440  << " currentMinimalStep= " << currentMinimalStep
441  << " isActive " << isActive << G4endl;
442  */
443  return gPathLength;
444 }
445 
446 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
447 
448 G4double
451 {
452  *condition = Forced;
453  //*condition = NotForced;
454  return DBL_MAX;
455 }
456 
457 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
458 
461 {
466  fPositionChanged = false;
467 
468  G4double geomLength = step.GetStepLength();
469 
470  // very small step - no msc
471  if(!isActive) {
472  tPathLength = geomLength;
473 
474  // sample msc
475  } else {
476  G4double range =
478  track.GetMaterialCutsCouple());
479 
481 
482  // protection against wrong t->g->t conversion
483  /*
484  if(currParticle->GetPDGMass() > GeV)
485  G4cout << "G4VMsc::AlongStepDoIt: GeomLength= "
486  << geomLength
487  << " tPathLength= " << tPathLength
488  << " physStepLimit= " << physStepLimit
489  << " dr= " << range - trueLength
490  << " ekin= " << track.GetKineticEnergy() << G4endl;
491  */
492  if (tPathLength > physStepLimit) {
494  }
495 
496  // do not sample scattering at the last or at a small step
497  if(tPathLength + geomMin < range && tPathLength > geomMin) {
498 
499  G4double preSafety = step.GetPreStepPoint()->GetSafety();
500  G4double maxDisp = (tPathLength + geomLength)*0.5;
501  G4double postSafety= preSafety - maxDisp;
502  G4bool safetyRecomputed = false;
503  if(postSafety < maxDisp) {
504  safetyRecomputed = true;
505  postSafety = safetyHelper->ComputeSafety(fNewPosition,maxDisp);
506  }
508  step.GetPostStepPoint()->GetMomentumDirection(), postSafety);
509 
510  G4double r2 = displacement.mag2();
511 
512  //G4cout << "R= " << sqrt(r2) << " postSafety= " << postSafety
513  // << G4endl;
514 
515  // make correction for displacement
516  if(r2 > 0.0) {
517 
518  fPositionChanged = true;
519  G4double fac = 1.0;
520 
521  // displaced point is definitely within the volume
522  if(r2 > postSafety*postSafety) {
523  G4double dispR = std::sqrt(r2);
524  if(!safetyRecomputed) {
525  postSafety = safetyHelper->ComputeSafety(fNewPosition, dispR);
526  }
527 
528  if(dispR > postSafety) {
529  fac = 0.99*postSafety/dispR;
530  }
531  }
532  // compute new endpoint of the Step
533  fNewPosition += fac*displacement;
534  //safetyHelper->ReLocateWithinVolume(fNewPosition);
535  }
536  }
537  }
539  //fParticleChange.ProposePosition(fNewPosition);
540  return &fParticleChange;
541 }
542 
543 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
544 
547 {
548  fParticleChange.Initialize(track);
549 
550  if(fPositionChanged) {
553  }
554 
555  return &fParticleChange;
556 }
557 
558 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
559 
561  const G4Track& track,
562  G4double previousStepSize,
563  G4double currentMinimalStep,
564  G4double& currentSafety)
565 {
567  G4double x = AlongStepGetPhysicalInteractionLength(track,previousStepSize,
568  currentMinimalStep,
569  currentSafety,
570  &selection);
571  return x;
572 }
573 
574 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
575 
577  const G4Track& track,
578  G4double previousStepSize,
579  G4double currentMinimalStep,
580  G4double& currentSafety)
581 {
582  return GetContinuousStepLimit(track,previousStepSize,currentMinimalStep,
583  currentSafety);
584 }
585 
586 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
587 
590 {
591  *condition = Forced;
592  return DBL_MAX;
593 }
594 
595 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
596 
597 G4bool
599  const G4String& directory,
600  G4bool ascii)
601 {
602  G4bool yes = true;
603  if(part != firstParticle) { return yes; }
605  static const G4String ss[4] = {"1","2","3","4"};
606  for(G4int i=0; i<nmod; ++i) {
607  G4VEmModel* msc = modelManager->GetModel(i);
608  yes = true;
609  G4PhysicsTable* table = msc->GetCrossSectionTable();
610  if (table) {
611  G4int j = std::min(i,3);
612  G4String name =
613  GetPhysicsTableFileName(part,directory,"LambdaMod"+ss[j],ascii);
614  yes = table->StorePhysicsTable(name,ascii);
615 
616  if ( yes ) {
617  if ( verboseLevel>0 ) {
618  G4cout << "Physics table are stored for "
619  << part->GetParticleName()
620  << " and process " << GetProcessName()
621  << " with a name <" << name << "> " << G4endl;
622  }
623  } else {
624  G4cout << "Fail to store Physics Table for "
625  << part->GetParticleName()
626  << " and process " << GetProcessName()
627  << " in the directory <" << directory
628  << "> " << G4endl;
629  }
630  }
631  }
632  return yes;
633 }
634 
635 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
636 
637 G4bool
639  const G4String&,
640  G4bool)
641 {
642  return true;
643 }
644 
645 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
646 
648 {
649  for(G4int i=0; i<numberOfModels; ++i) {
650  G4VMscModel* msc = static_cast<G4VMscModel*>(GetModelByIndex(i, true));
651  msc->SetIonisation(p, firstParticle);
652  }
653 }
654 
655 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
656 
G4double condition(const G4ErrorSymMatrix &m)
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:538
G4VEmModel * SelectModel(G4double kinEnergy, size_t idx)
virtual void Initialize(const G4Track &)
G4SafetyHelper * GetSafetyHelper() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
void ReLocateWithinVolume(const G4ThreeVector &pGlobalPoint)
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:271
void InitialiseHelper()
static G4LossTableManager * Instance()
G4EmModelManager * modelManager
G4int verboseLevel
Definition: G4VProcess.hh:368
static const G4double fac
CLHEP::Hep3Vector G4ThreeVector
G4double GetStepLength() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
void SetLateralDisplasmentFlag(G4bool val)
G4VEmModel * GetModel(G4int, G4bool ver=false)
void DeRegister(G4VEnergyLossProcess *p)
G4String name
Definition: TRTMaterials.hh:40
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double GetMeanFreePath(const G4Track &track, G4double, G4ForceCondition *condition)
void SetLateralDisplasmentFlag(G4bool val)
Definition: G4VMscModel.hh:196
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *condition)
void SetEmModel(G4VMscModel *, G4int index=1)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *, const G4Region *)
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:784
void SetStepLimitType(G4MscStepLimitType)
Definition: G4VMscModel.hh:224
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
G4VEnergyLossProcess * fIonisation
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:186
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
G4VMultipleScattering(const G4String &name="msc", G4ProcessType type=fElectromagnetic)
G4StepPoint * GetPreStepPoint() const
const G4ParticleDefinition * currParticle
G4MscStepLimitType StepLimitType() const
G4ParticleChangeForMSC fParticleChange
virtual G4double ComputeTrueStepLength(G4double geomPathLength)
Definition: G4VMscModel.cc:159
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)
Definition: G4VMscModel.cc:145
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=0)
Definition: G4VEmModel.cc:387
const G4ThreeVector & GetPosition() const
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:288
bool G4bool
Definition: G4Types.hh:79
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
void SetRangeFactor(G4double)
Definition: G4VMscModel.hh:210
void ProposeTrueStepLength(G4double truePathLength)
const G4ParticleDefinition * GetParticleDefinition() const
G4bool IsMaster() const
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
const G4String & GetParticleType() const
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:395
void Register(G4VEnergyLossProcess *p)
Definition: G4Step.hh:76
const G4int n
G4LossTableManager * emManager
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4bool IsActive(G4double kinEnergy)
Definition: G4VEmModel.hh:711
void DumpModelList(G4int verb)
G4bool LateralDisplasmentFlag() const
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4double minSubRange, G4int verb)
G4double ComputeSafety(const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)
static G4TransportationManager * GetTransportationManager()
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:669
G4int size() const
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
static const double eV
Definition: G4SIunits.hh:194
#define fm
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
Definition: G4VMscModel.hh:335
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:93
G4int NumberOfModels() const
void ProposePosition(const G4ThreeVector &finalPosition)
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
G4double ContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
const G4ParticleDefinition * firstParticle
G4StepPoint * GetPostStepPoint() const
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double GetSafety() const
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
void SetGeomFactor(G4double)
Definition: G4VMscModel.hh:217
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
virtual void PrintInfo()=0
double G4double
Definition: G4Types.hh:76
void SetSkin(G4double)
Definition: G4VMscModel.hh:203
void BuildPhysicsTable(const G4ParticleDefinition &)
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4ForceCondition
void SetRangeFactor(G4double val)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection)
#define DBL_MAX
Definition: templates.hh:83
G4double MaxKinEnergy() const
G4VMscModel * EmModel(G4int index=1) const
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
std::vector< G4VMscModel * > mscModels
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:718
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
G4GPILSelection
void SetStepLimitType(G4MscStepLimitType val)
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
void PreparePhysicsTable(const G4ParticleDefinition &)
G4ProcessType
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)
Definition: G4VMscModel.cc:138
void SetIonisation(G4VEnergyLossProcess *)