Geant4  10.01.p01
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 88981 2015-03-17 10:14:15Z 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 
90 static const G4double minSafety = 1.20*CLHEP::nm;
91 static const G4double geomMin = 0.05*CLHEP::nm;
93 
97  numberOfModels(0),
98  firstParticle(0),
99  currParticle(0),
100  stepLimit(fUseSafety),
101  facrange(0.04),
102  latDisplacement(true),
103  isIon(false),
104  fDispBeyondSafety(false)
105 {
107  SetVerboseLevel(1);
109  if("ionmsc" == name) { firstParticle = G4GenericIon::GenericIon(); }
110 
112 
114  fIonisation = 0;
115 
117  safetyHelper = 0;
118  fPositionChanged = false;
119  isActive = false;
120  actStepLimit = false;
121  actFacRange = false;
122  actLatDisp = false;
123 
124  currentModel = 0;
127  emManager->Register(this);
128 }
129 
130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
131 
133 {
134  /*
135  if(1 < verboseLevel) {
136  G4cout << "G4VMultipleScattering destruct " << GetProcessName()
137  << G4endl;
138  }
139  */
140  delete modelManager;
141  emManager->DeRegister(this);
142 }
143 
144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
145 
147  const G4Region* region)
148 {
149  G4VEmFluctuationModel* fm = 0;
150  modelManager->AddEmModel(order, p, fm, region);
151  if(p) { p->SetParticleChange(pParticleChange); }
152 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
155 
157 {
158  G4int n = mscModels.size();
159  if(index >= n) { for(G4int i=n; i<=index; ++i) { mscModels.push_back(0); } }
160  mscModels[index] = p;
161 }
162 
163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
164 
166 {
167  G4VMscModel* p = 0;
168  if(index >= 0 && index < G4int(mscModels.size())) { p = mscModels[index]; }
169  return p;
170 }
171 
172 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
173 
174 G4VEmModel*
176 {
177  return modelManager->GetModel(idx, ver);
178 }
179 
180 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
181 
182 void
184 {
185  G4bool master = true;
186  const G4VMultipleScattering* masterProc =
187  static_cast<const G4VMultipleScattering*>(GetMasterProcess());
188  if(masterProc && masterProc != this) { master = false; }
189 
190  if(!firstParticle) { firstParticle = &part; }
191  if(part.GetParticleType() == "nucleus") {
194  SetRangeFactor(0.2);
195  G4String pname = part.GetParticleName();
196  if(pname != "deuteron" && pname != "triton" &&
197  pname != "alpha+" && pname != "helium" &&
198  pname != "alpha" && pname != "He3" &&
199  pname != "hydrogen") {
200 
201  const G4ParticleDefinition* theGenericIon =
203 
204  if(theGenericIon && firstParticle != theGenericIon) {
205  G4ProcessManager* pm = theGenericIon->GetProcessManager();
207  size_t n = v->size();
208  for(size_t j=0; j<n; ++j) {
209  if((*v)[j] == this) {
210  firstParticle = theGenericIon;
211  isIon = true;
212  break;
213  }
214  }
215  }
216  }
217  }
218 
219  emManager->PreparePhysicsTable(&part, this, master);
220  currParticle = 0;
221 
222  if(1 < verboseLevel) {
223  G4cout << "### G4VMultipleScattering::PrepearPhysicsTable() for "
224  << GetProcessName()
225  << " and particle " << part.GetParticleName()
226  << " local particle " << firstParticle->GetParticleName()
227  << " isIon= " << isIon
228  << G4endl;
229  }
230 
231  if(firstParticle == &part) {
232 
233  // initialise process
235  if(part.GetPDGMass() > MeV) {
236  if(!actStepLimit) { stepLimit = fMinimal; }
237  if(!actFacRange) { facrange = 0.2; }
238  if(!actLatDisp) {
240  }
241  } else {
244  if(!actLatDisp) {
246  }
247  }
248  if(latDisplacement) {
250  }
251  if(master) { SetVerboseLevel(theParameters->Verbose()); }
253 
254  // initialisation of models
256  for(G4int i=0; i<numberOfModels; ++i) {
257  G4VMscModel* msc = static_cast<G4VMscModel*>(modelManager->GetModel(i));
258  msc->SetIonisation(0, firstParticle);
259  msc->SetMasterThread(master);
260  if(0 == i) { currentModel = msc; }
263  msc->SetSkin(theParameters->MscSkin());
264  msc->SetRangeFactor(facrange);
267  G4double emax =
269  msc->SetHighEnergyLimit(emax);
270  }
271 
273  10.0, verboseLevel);
274 
275  if(!safetyHelper) {
277  ->GetSafetyHelper();
279  }
280  }
281 }
282 
283 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
284 
286 {
287  G4String num = part.GetParticleName();
288  if(1 < verboseLevel) {
289  G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
290  << GetProcessName()
291  << " and particle " << num
292  << " IsMaster= " << G4LossTableManager::Instance()->IsMaster()
293  << G4endl;
294  }
295  G4bool master = true;
296  const G4VMultipleScattering* masterProcess =
297  static_cast<const G4VMultipleScattering*>(GetMasterProcess());
298  if(masterProcess && masterProcess != this) { master = false; }
299 
300  if(firstParticle == &part) {
301  /*
302  G4cout << "### G4VMultipleScattering::BuildPhysicsTable() for "
303  << GetProcessName()
304  << " and particle " << num
305  << " IsMaster= " << G4LossTableManager::Instance()->IsMaster()
306  << " " << this
307  << G4endl;
308  */
310 
311  if(!master) {
312  // initialisation of models
313  G4bool printing = true;
315  /*
316  G4cout << "### G4VMultipleScattering::SlaveBuildPhysicsTable() for "
317  << GetProcessName()
318  << " and particle " << num
319  << " Nmod= " << numberOfModels << " " << this
320  << G4endl;
321  */
322  for(G4int i=0; i<numberOfModels; ++i) {
323  G4VMscModel* msc =
324  static_cast<G4VMscModel*>(GetModelByIndex(i, printing));
325  G4VMscModel* msc0=
326  static_cast<G4VMscModel*>(masterProcess->GetModelByIndex(i,printing));
327  msc->SetCrossSectionTable(msc0->GetCrossSectionTable(), false);
328  msc->InitialiseLocal(firstParticle, msc0);
329  }
330  }
331  }
332 
333  // explicitly defined printout by particle name
334  if(1 < verboseLevel ||
335  (0 < verboseLevel && (num == "e-" ||
336  num == "e+" || num == "mu+" ||
337  num == "mu-" || num == "proton"||
338  num == "pi+" || num == "pi-" ||
339  num == "kaon+" || num == "kaon-" ||
340  num == "alpha" || num == "anti_proton" ||
341  num == "GenericIon")))
342  {
343  G4cout << G4endl << GetProcessName()
344  << ": for " << num
345  << " SubType= " << GetProcessSubType()
346  << G4endl;
347  PrintInfo();
349  }
350 
351  if(1 < verboseLevel) {
352  G4cout << "### G4VMultipleScattering::BuildPhysicsTable() done for "
353  << GetProcessName()
354  << " and particle " << num
355  << G4endl;
356  }
357 }
358 
359 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
360 
362 {
363  if (0 < verboseLevel) {
364  G4cout << G4endl << GetProcessName()
365  << ": for " << firstParticle->GetParticleName()
366  << " SubType= " << GetProcessSubType()
367  << G4endl;
368  PrintInfo();
370  }
371 }
372 
373 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
374 
376 {
377  G4VEnergyLossProcess* eloss = 0;
378  if(track->GetParticleDefinition() != currParticle) {
381  eloss = fIonisation;
382  }
383  /*
384  G4cout << "G4VMultipleScattering::StartTracking Nmod= " << numberOfModels
385  << " " << currParticle->GetParticleName()
386  << " E(MeV)= " << track->GetKineticEnergy()
387  << " Ion= " << eloss << " " << fIonisation << " IsMaster= "
388  << G4LossTableManager::Instance()->IsMaster()
389  << G4endl;
390  */
391  // one model
392  if(1 == numberOfModels) {
393  currentModel->StartTracking(track);
395 
396  // many models
397  } else {
398  for(G4int i=0; i<numberOfModels; ++i) {
399  G4VMscModel* msc = static_cast<G4VMscModel*>(GetModelByIndex(i,true));
400  msc->StartTracking(track);
401  if(eloss) { msc->SetIonisation(fIonisation, currParticle); }
402  }
403  }
404 }
405 
406 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
407 
409  const G4Track& track,
410  G4double,
411  G4double currentMinimalStep,
412  G4double&,
413  G4GPILSelection* selection)
414 {
415  // get Step limit proposed by the process
416  *selection = NotCandidateForSelection;
417  physStepLimit = gPathLength = tPathLength = currentMinimalStep;
418 
419  G4double ekin = track.GetKineticEnergy();
420  /*
421  G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
422  << " " << currParticle->GetParticleName()
423  << " currMod " << currentModel
424  << G4endl;
425  */
426  // isIon flag is used only to select a model
427  if(isIon) {
428  ekin *= proton_mass_c2/track.GetParticleDefinition()->GetPDGMass();
429  }
430 
431  // select new model
432  if(1 < numberOfModels) {
433  currentModel = static_cast<G4VMscModel*>(
434  SelectModel(ekin,track.GetMaterialCutsCouple()->GetIndex()));
435  }
436  // msc is active is model is active, energy above the limit,
437  // and step size is above the limit;
438  // if it is active msc may limit the step
440  && ekin >= lowestKinEnergy) {
441  isActive = true;
442  tPathLength =
444  if (tPathLength < physStepLimit) {
445  *selection = CandidateForSelection;
446  }
447  } else { isActive = false; }
448 
449  //if(currParticle->GetPDGMass() > GeV)
450  /*
451  G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
452  << " " << currParticle->GetParticleName()
453  << " gPathLength= " << gPathLength
454  << " tPathLength= " << tPathLength
455  << " currentMinimalStep= " << currentMinimalStep
456  << " isActive " << isActive << G4endl;
457  */
458  return gPathLength;
459 }
460 
461 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
462 
463 G4double
466 {
467  *condition = Forced;
468  //*condition = NotForced;
469  return DBL_MAX;
470 }
471 
472 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
473 
476 {
481  fPositionChanged = false;
482 
483  G4double geomLength = step.GetStepLength();
484 
485  // very small step - no msc
486  if(!isActive) {
487  tPathLength = geomLength;
488 
489  // sample msc
490  } else {
491  G4double range =
493  track.GetMaterialCutsCouple());
494 
496 
497  // protection against wrong t->g->t conversion
498  /*
499  if(currParticle->GetPDGMass() > 0.9*GeV)
500  G4cout << "G4VMsc::AlongStepDoIt: GeomLength= "
501  << geomLength
502  << " tPathLength= " << tPathLength
503  << " physStepLimit= " << physStepLimit
504  << " dr= " << range - tPathLength
505  << " ekin= " << track.GetKineticEnergy() << G4endl;
506  */
508 
509  // do not sample scattering at the last or at a small step
510  if(tPathLength < range && tPathLength > geomMin) {
511 
514 
515  G4double r2 = displacement.mag2();
516  //G4cout << " R= " << sqrt(r2) << " Rmin= " << sqrt(minDisplacement2)
517  // << " flag= " << fDispBeyondSafety << G4endl;
518  if(r2 > minDisplacement2) {
519 
520  fPositionChanged = true;
521  const G4double sFact = 0.99;
522  G4double postSafety =
523  sFact*(step.GetPreStepPoint()->GetSafety() - geomLength);
524  //G4cout<<" R= "<<sqrt(r2)<<" postSafety= "<<postSafety<<G4endl;
525 
526  // far away from geometry boundary
527  if(postSafety > 0.0 && r2 <= postSafety*postSafety) {
528  fNewPosition += displacement;
529 
530  } else {
531  G4double dispR = std::sqrt(r2);
532  postSafety =
533  sFact*safetyHelper->ComputeSafety(fNewPosition, dispR);
534 
535  // displaced point is definitely within the volume
536  //G4cout<<" R= "<<dispR<<" postSafety= "<<postSafety<<G4endl;
537  if(dispR < postSafety) {
538  fNewPosition += displacement;
539 
540  // optional extra mechanism is applied only if a particle
541  // is stopped by the boundary
542  } else if(fDispBeyondSafety && 0.0 == postSafety) {
543  fNewPosition += displacement;
544  G4double maxshift =
545  std::min(2.0*dispR, physStepLimit-tPathLength);
546  G4double dist = 0.0;
547  G4double safety = postSafety + dispR;
549  /*
550  G4cout << "##MSC before Recheck maxshift= " << maxshift
551  << " postsafety= " << postSafety
552  << " Ekin= " << track.GetKineticEnergy()
553  << " " << track.GetDefinition()->GetParticleName()
554  << G4endl;
555  */
556  // check if it is possible to shift to the boundary
559  maxshift,
560  &dist,
561  &safety))
562  {
563  /*
564  G4cout << "##MSC after Recheck dist= " << dist
565  << " postsafety= " << postSafety
566  << " t= " << tPathLength
567  << " g= " << geomLength
568  << " p= " << physStepLimit
569  << G4endl;
570  */
571  G4double tnew = tPathLength*(1.0 + dist/geomLength);
572  if(tnew >= 0.0 && tnew < physStepLimit) {
573  tPathLength = tnew;
574  fNewPosition += dist*fNewDirection;
575  } else {
576  fNewPosition += displacement*(postSafety/dispR - 1.0);
577  }
578  }
579  else
580  // shift on boundary is not possible
581  {
582  fNewPosition += displacement*(postSafety/dispR - 1.0);
583  }
584  // reduced displacement
585  } else if(postSafety > geomMin) {
586  fNewPosition += displacement*(postSafety/dispR);
587 
588  // very small postSafety
589  } else {
590  fPositionChanged = false;
591  }
592  }
593  //safetyHelper->ReLocateWithinVolume(fNewPosition);
594  }
595  }
596  }
598  //fParticleChange.ProposePosition(fNewPosition);
599  return &fParticleChange;
600 }
601 
602 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
603 
606 {
608 
609  if(fPositionChanged) {
612  }
613 
614  return &fParticleChange;
615 }
616 
617 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
618 
620  const G4Track& track,
621  G4double previousStepSize,
622  G4double currentMinimalStep,
623  G4double& currentSafety)
624 {
626  G4double x = AlongStepGetPhysicalInteractionLength(track,previousStepSize,
627  currentMinimalStep,
628  currentSafety,
629  &selection);
630  return x;
631 }
632 
633 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
634 
636  const G4Track& track,
637  G4double previousStepSize,
638  G4double currentMinimalStep,
639  G4double& currentSafety)
640 {
641  return GetContinuousStepLimit(track,previousStepSize,currentMinimalStep,
642  currentSafety);
643 }
644 
645 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
646 
649 {
650  *condition = Forced;
651  return DBL_MAX;
652 }
653 
654 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
655 
656 G4bool
658  const G4String& directory,
659  G4bool ascii)
660 {
661  G4bool yes = true;
662  if(part != firstParticle) { return yes; }
663  const G4VMultipleScattering* masterProcess =
664  static_cast<const G4VMultipleScattering*>(GetMasterProcess());
665  if(masterProcess && masterProcess != this) { return yes; }
666 
668  static const G4String ss[4] = {"1","2","3","4"};
669  for(G4int i=0; i<nmod; ++i) {
670  G4VEmModel* msc = modelManager->GetModel(i);
671  yes = true;
672  G4PhysicsTable* table = msc->GetCrossSectionTable();
673  if (table) {
674  G4int j = std::min(i,3);
675  G4String name =
676  GetPhysicsTableFileName(part,directory,"LambdaMod"+ss[j],ascii);
677  yes = table->StorePhysicsTable(name,ascii);
678 
679  if ( yes ) {
680  if ( verboseLevel>0 ) {
681  G4cout << "Physics table are stored for "
682  << part->GetParticleName()
683  << " and process " << GetProcessName()
684  << " with a name <" << name << "> " << G4endl;
685  }
686  } else {
687  G4cout << "Fail to store Physics Table for "
688  << part->GetParticleName()
689  << " and process " << GetProcessName()
690  << " in the directory <" << directory
691  << "> " << G4endl;
692  }
693  }
694  }
695  return yes;
696 }
697 
698 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
699 
700 G4bool
702  const G4String&,
703  G4bool)
704 {
705  return true;
706 }
707 
708 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
709 
711 {
712  for(G4int i=0; i<numberOfModels; ++i) {
713  G4VMscModel* msc = static_cast<G4VMscModel*>(GetModelByIndex(i, true));
714  msc->SetIonisation(p, firstParticle);
715  }
716 }
717 
718 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
719 
G4double condition(const G4ErrorSymMatrix &m)
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
const G4ThreeVector * GetMomentumDirection() const
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:538
static const G4double geomMin
static const double MeV
Definition: G4SIunits.hh:193
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:284
void InitialiseHelper()
G4int WorkerVerbose() const
G4double MaxKinEnergy() const
static G4LossTableManager * Instance()
G4EmModelManager * modelManager
G4int verboseLevel
Definition: G4VProcess.hh:368
CLHEP::Hep3Vector G4ThreeVector
G4double GetStepLength() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:615
void SetLateralDisplasmentFlag(G4bool val)
G4VEmModel * GetModel(G4int, G4bool ver=false)
void DeRegister(G4VEnergyLossProcess *p)
G4double MscGeomFactor() const
G4String name
Definition: TRTMaterials.hh:40
G4double MscThetaLimit() const
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:219
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4double GetMeanFreePath(const G4Track &track, G4double, G4ForceCondition *condition)
void SetLateralDisplasmentFlag(G4bool val)
Definition: G4VMscModel.hh:197
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:807
void SetStepLimitType(G4MscStepLimitType)
Definition: G4VMscModel.hh:225
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
static const G4double minDisplacement2
const G4String & GetParticleName() const
G4bool LatDisplacementBeyondSafety() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:706
G4VMultipleScattering(const G4String &name="msc", G4ProcessType type=fElectromagnetic)
G4StepPoint * GetPreStepPoint() const
const G4ParticleDefinition * currParticle
G4ParticleChangeForMSC fParticleChange
virtual G4double ComputeTrueStepLength(G4double geomPathLength)
Definition: G4VMscModel.cc:159
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4int Verbose() const
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)
Definition: G4VMscModel.cc:145
G4bool MuHadLateralDisplacement() const
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=0)
Definition: G4VEmModel.cc:402
const G4ThreeVector & GetPosition() const
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.hh:289
bool G4bool
Definition: G4Types.hh:79
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
void SetRangeFactor(G4double)
Definition: G4VMscModel.hh:211
void ProposeTrueStepLength(G4double truePathLength)
G4bool LateralDisplacement() const
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:410
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:734
void DumpModelList(G4int verb)
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:692
static const G4double minSafety
G4int size() const
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
static const double eV
Definition: G4SIunits.hh:194
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
Definition: G4VMscModel.hh:336
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)
static G4EmParameters * Instance()
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double MscRangeFactor() const
G4double GetSafety() const
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
G4MscStepLimitType stepLimit
void SetGeomFactor(G4double)
Definition: G4VMscModel.hh:218
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
G4MscStepLimitType MscStepLimitType() const
virtual void PrintInfo()=0
double G4double
Definition: G4Types.hh:76
void SetSkin(G4double)
Definition: G4VMscModel.hh:204
void BuildPhysicsTable(const G4ParticleDefinition &)
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4ForceCondition
G4bool RecheckDistanceToCurrentBoundary(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double *prDistance, G4double *prNewSafety=0) const
void SetRangeFactor(G4double val)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection)
#define DBL_MAX
Definition: templates.hh:83
G4VMscModel * EmModel(G4int index=1) const
G4double MscSkin() const
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
std::vector< G4VMscModel * > mscModels
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:741
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 *)