Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4VMultipleScattering
33//
34// Author: Vladimir Ivanchenko on base of Laszlo Urban code
35//
36// Creation date: 25.03.2003
37//
38// Modifications:
39//
40// 13.04.03 Change printout (V.Ivanchenko)
41// 04-06-03 Fix compilation warnings (V.Ivanchenko)
42// 16-07-03 Use G4VMscModel interface (V.Ivanchenko)
43// 03-11-03 Fix initialisation problem in RetrievePhysicsTable (V.Ivanchenko)
44// 04-11-03 Update PrintInfoDefinition (V.Ivanchenko)
45// 01-03-04 SampleCosineTheta signature changed
46// 22-04-04 SampleCosineTheta signature changed back to original
47// 27-08-04 Add InitialiseForRun method (V.Ivanchneko)
48// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
49// 11-03-05 Shift verbose level by 1 (V.Ivantchenko)
50// 15-04-05 optimize internal interface (V.Ivanchenko)
51// 15-04-05 remove boundary flag (V.Ivanchenko)
52// 27-10-05 introduce virtual function MscStepLimitation() (V.Ivanchenko)
53// 12-04-07 Add verbosity at destruction (V.Ivanchenko)
54// 27-10-07 Virtual functions moved to source (V.Ivanchenko)
55// 11-03-08 Set skin value does not effect step limit type (V.Ivanchenko)
56// 24-06-09 Removed hidden bin in G4PhysicsVector (V.Ivanchenko)
57// 04-06-13 Adoptation to MT mode (V.Ivanchenko)
58//
59// Class Description:
60//
61// It is the generic process of multiple scattering it includes common
62// part of calculations for all charged particles
63
64// -------------------------------------------------------------------
65//
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
67//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
68
71#include "G4SystemOfUnits.hh"
72#include "G4LossTableManager.hh"
74#include "G4Step.hh"
77#include "G4UnitsTable.hh"
79#include "G4Electron.hh"
80#include "G4GenericIon.hh"
82#include "G4SafetyHelper.hh"
83#include "G4ParticleTable.hh"
84#include "G4ProcessVector.hh"
85#include "G4ProcessManager.hh"
86#include <iostream>
87
88//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89
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}
126
127//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
128
130{
131 /*
132 if(1 < verboseLevel) {
133 G4cout << "G4VMultipleScattering destruct " << GetProcessName()
134 << G4endl;
135 }
136 */
137 delete modelManager;
138 emManager->DeRegister(this);
139}
140
141//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
142
144 const G4Region* region)
145{
146 if(!p) { return; }
147 G4VEmFluctuationModel* fm = nullptr;
148 modelManager->AddEmModel(order, p, fm, region);
150}
151
152//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
153
155{
156 for(auto & msc : mscModels) { if(msc == p) { return; } }
157 mscModels.push_back(p);
158}
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
161
163{
164 return (index < mscModels.size()) ? mscModels[index] : nullptr;
165}
166
167//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
168
169void
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}
273
274//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
275
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}
341
342//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
343
344void G4VMultipleScattering::StreamInfo(std::ostream& outFile,
345 const G4ParticleDefinition& part, G4bool rst) const
346{
347 G4String indent = (rst ? " " : "");
348 outFile << G4endl << indent << GetProcessName() << ": ";
349 if (!rst) outFile << " for " << part.GetParticleName();
350 outFile << " SubType= " << GetProcessSubType() << G4endl;
351 //StreamProcessInfo(outFile);
352 modelManager->DumpModelList(outFile, verboseLevel);
353}
354
355//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
356
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}
388
389//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
390
392 const G4Track& track,
393 G4double,
394 G4double currentMinimalStep,
395 G4double&,
396 G4GPILSelection* selection)
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}
443
444//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
445
449{
451 return DBL_MAX;
452}
453
454//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
455
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}
541
542//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
543
546{
548 return &fParticleChange;
549}
550
551//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
552
554 const G4Track& track,
555 G4double previousStepSize,
556 G4double currentMinimalStep,
557 G4double& currentSafety)
558{
560 G4double x = AlongStepGetPhysicalInteractionLength(track,previousStepSize,
561 currentMinimalStep,
562 currentSafety,
563 &selection);
564 return x;
565}
566
567//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
568
570 const G4Track& track,
571 G4double previousStepSize,
572 G4double currentMinimalStep,
573 G4double& currentSafety)
574{
575 return GetContinuousStepLimit(track,previousStepSize,currentMinimalStep,
576 currentSafety);
577}
578
579//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
580
583{
584 *condition = Forced;
585 return DBL_MAX;
586}
587
588//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
589
590G4bool
592 const G4String& directory,
593 G4bool ascii)
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}
631
632//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
633
634G4bool
636 const G4String&,
637 G4bool)
638{
639 return true;
640}
641
642//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
643
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}
651
652//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
653
654void G4VMultipleScattering::ProcessDescription(std::ostream& outFile) const
655{
656 if(firstParticle) {
657 StreamInfo(outFile, *firstParticle, true);
658 }
659}
660
661//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
662
@ fMultipleScattering
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
@ Forced
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
@ fMinimal
@ fUseSafety
G4ProcessType
@ fElectromagnetic
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double mag2() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *, const G4Region *)
void DumpModelList(std::ostream &out, G4int verb)
G4int NumberOfModels() const
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4double minSubRange, G4int verb)
G4VEmModel * GetModel(G4int idx, G4bool ver=false)
static G4EmParameters * Instance()
G4double MscMuHadRangeFactor() const
G4double MscThetaLimit() 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
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
static G4LossTableManager * Instance()
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
void DeRegister(G4VEnergyLossProcess *p)
void Register(G4VEnergyLossProcess *p)
virtual void Initialize(const G4Track &)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
void ProposePosition(const G4ThreeVector &finalPosition)
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleType() const
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
G4double ComputeSafety(const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)
void ReLocateWithinVolume(const G4ThreeVector &pGlobalPoint)
void InitialiseHelper()
const G4ThreeVector & GetPosition() const
const G4ThreeVector & GetMomentumDirection() const
Definition: G4Step.hh:62
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:464
void SetPolarAngleLimit(G4double)
Definition: G4VEmModel.hh:792
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:729
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:456
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
G4bool IsActive(G4double kinEnergy) const
Definition: G4VEmModel.hh:785
G4PhysicsTable * GetCrossSectionTable()
Definition: G4VEmModel.hh:860
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:288
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:220
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)=0
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
Definition: G4VMscModel.hh:395
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 SetIonisation(G4VEnergyLossProcess *)
G4double GetMeanFreePath(const G4Track &track, G4double, G4ForceCondition *condition) override
G4ParticleChangeForMSC fParticleChange
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=nullptr)
void StartTracking(G4Track *) override
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
void PreparePhysicsTable(const G4ParticleDefinition &) override
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4double ContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety)
G4VMultipleScattering(const G4String &name="msc", G4ProcessType type=fElectromagnetic)
void SetEmModel(G4VMscModel *, size_t index=0)
virtual void ProcessDescription(std::ostream &outFile) const override
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
virtual void InitialiseProcess(const G4ParticleDefinition *)=0
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection) override
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *condition) override
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
G4VEmModel * SelectModel(G4double kinEnergy, size_t idx)
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety) override
G4VMscModel * EmModel(size_t index=0) const
void ProposeTrueStepLength(G4double truePathLength)
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:518
G4int verboseLevel
Definition: G4VProcess.hh:356
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
G4int GetProcessSubType() const
Definition: G4VProcess.hh:400
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:182
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
#define DBL_MAX
Definition: templates.hh:62