Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4TransportationWithMsc.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// G4TransportationWithMsc
27//
28// Class Description:
29//
30// It is a generic process of transportation with multiple scattering included
31// in the step limitation and propagation.
32//
33// Original author: Jonas Hahnfeld, 2022
34
35// -------------------------------------------------------------------
36//
37//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
38//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
39
41
42#include "G4LossTableBuilder.hh"
43#include "G4LossTableManager.hh"
44#include "G4EmConfigurator.hh"
45#include "G4VMscModel.hh"
46
47#include "G4DynamicParticle.hh"
48#include "G4Step.hh"
49#include "G4StepPoint.hh"
50#include "G4StepStatus.hh"
51#include "G4Track.hh"
52
53#include "G4Electron.hh"
54
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
58
59static constexpr G4double kLowestKinEnergy = 10 * CLHEP::eV;
60static constexpr G4double kGeomMin = 0.05 * CLHEP::nm;
61static constexpr G4double kMinDisplacement2 = kGeomMin * kGeomMin;
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64
66 G4int verbosity)
67 : G4Transportation(verbosity, "TransportationWithMsc")
68 , fType(type)
69{
71
72 fEmManager = G4LossTableManager::Instance();
73 fModelManager = new G4EmModelManager;
74
75 G4ThreeVector zero;
76 fSubStepDynamicParticle =
78 fSubStepTrack = new G4Track(fSubStepDynamicParticle, 0, zero);
79 fSubStep = new G4Step;
80 fSubStepTrack->SetStep(fSubStep);
81}
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
84
86{
87 delete fModelManager;
88 // fSubStepDynamicParticle is owned and also deleted by fSubStepTrack!
89 delete fSubStepTrack;
90 delete fSubStep;
91}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
94
96 const G4Region* region)
97{
99 {
100 G4Exception("G4TransportationWithMsc::AddMscModel", "em0051",
102 "not allowed unless type == MultipleScattering");
103 }
104
105 fModelManager->AddEmModel(order, mscModel, nullptr, region);
107}
108
109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
110
112 const G4ParticleDefinition& part)
113{
114 if(nullptr == fFirstParticle)
115 {
116 fFirstParticle = ∂
117 G4VMultipleScattering* ptr = nullptr;
118 auto emConfigurator = fEmManager->EmConfigurator();
119 emConfigurator->PrepareModels(&part, ptr, this);
120 }
121
122 if(fFirstParticle == &part)
123 {
124 G4bool master = fEmManager->IsMaster();
125 G4LossTableBuilder* bld = fEmManager->GetTableBuilder();
126 G4bool baseMat = bld->GetBaseMaterialFlag();
127 const auto* theParameters = G4EmParameters::Instance();
128
129 if(master)
130 {
131 SetVerboseLevel(theParameters->Verbose());
132 }
133 else
134 {
135 SetVerboseLevel(theParameters->WorkerVerbose());
136 }
137
138 const G4int numberOfModels = fModelManager->NumberOfModels();
139 for(G4int i = 0; i < numberOfModels; ++i)
140 {
141 auto msc = static_cast<G4VMscModel*>(fModelManager->GetModel(i));
142 msc->SetMasterThread(master);
143 msc->SetPolarAngleLimit(theParameters->MscThetaLimit());
144 G4double emax =
145 std::min(msc->HighEnergyLimit(), theParameters->MaxKinEnergy());
146 msc->SetHighEnergyLimit(emax);
147 msc->SetUseBaseMaterials(baseMat);
148 }
149
150 fModelManager->Initialise(fFirstParticle, G4Electron::Electron(),
152 }
153}
154
155//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
156
158 const G4ParticleDefinition& part)
159{
160 if(fFirstParticle == &part)
161 {
162 fEmManager->BuildPhysicsTable(fFirstParticle);
163
164 if(!fEmManager->IsMaster())
165 {
166 const auto masterProcess =
167 static_cast<const G4TransportationWithMsc*>(GetMasterProcess());
168
169 // Initialisation of models.
170 const G4int numberOfModels = fModelManager->NumberOfModels();
171 for(G4int i = 0; i < numberOfModels; ++i)
172 {
173 auto msc = static_cast<G4VMscModel*>(fModelManager->GetModel(i));
174 auto msc0 =
175 static_cast<G4VMscModel*>(masterProcess->fModelManager->GetModel(i));
176 msc->SetCrossSectionTable(msc0->GetCrossSectionTable(), false);
177 msc->InitialiseLocal(fFirstParticle, msc0);
178 }
179 }
180 }
181
183 {
184 G4cout << G4endl;
185 G4cout << GetProcessName() << ": for " << part.GetParticleName();
186 if(fMultipleSteps)
187 {
188 G4cout << " (multipleSteps: 1)";
189 }
190 G4cout << G4endl;
191 fModelManager->DumpModelList(G4cout, verboseLevel);
192 }
193}
194
195//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
196
198{
199 auto* currParticle = track->GetParticleDefinition();
200 auto* ionisation = fEmManager->GetEnergyLossProcess(currParticle);
201
202 fSubStepDynamicParticle->SetDefinition(currParticle);
203
204 const G4int numberOfModels = fModelManager->NumberOfModels();
205 for(G4int i = 0; i < numberOfModels; ++i)
206 {
207 auto msc = static_cast<G4VMscModel*>(fModelManager->GetModel(i));
208 msc->StartTracking(track);
209 msc->SetIonisation(ionisation, currParticle);
210 }
211
212 // Ensure that field propagation state is also cleared / prepared
214}
215
216//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
217
219 const G4Track& track, G4double previousStepSize, G4double currentMinimumStep,
220 G4double& proposedSafety, G4GPILSelection* selection)
221{
222 *selection = NotCandidateForSelection;
223
224 const G4double physStepLimit = currentMinimumStep;
225
226 switch(fType)
227 {
229 // Select the MSC model for the current kinetic energy.
230 G4VMscModel* mscModel = nullptr;
231 const G4double ekin = track.GetKineticEnergy();
232 const auto* couple = track.GetMaterialCutsCouple();
233 const auto* particleDefinition = track.GetParticleDefinition();
234 if(physStepLimit > kGeomMin)
235 {
236 G4double ekinForSelection = ekin;
237 G4double pdgMass = particleDefinition->GetPDGMass();
238 if(pdgMass > CLHEP::GeV)
239 {
240 ekinForSelection *= proton_mass_c2 / pdgMass;
241 }
242
243 if(ekinForSelection >= kLowestKinEnergy)
244 {
245 mscModel = static_cast<G4VMscModel*>(
246 fModelManager->SelectModel(ekinForSelection, couple->GetIndex()));
247 if(mscModel == nullptr)
248 {
249 G4Exception("G4TransportationWithMsc::AlongStepGPIL", "em0052",
250 FatalException, "no MSC model found");
251 }
252 if(!mscModel->IsActive(ekinForSelection))
253 {
254 mscModel = nullptr;
255 }
256 }
257 }
258
259 // Call the MSC model to potentially limit the step and convert to
260 // geometric path length.
261 if(mscModel != nullptr)
262 {
263 mscModel->SetCurrentCouple(couple);
264
265 // Use the provided track for the first step.
266 const G4Track* currentTrackPtr = &track;
267
268 G4double currentSafety = proposedSafety;
269 G4double currentEnergy = ekin;
270
271 G4double stepLimitLeft = physStepLimit;
272 G4double totalGeometryStepLength = 0, totalTruePathLength = 0;
273 G4bool firstStep = true, continueStepping = fMultipleSteps;
274
275 do
276 {
277 G4double gPathLength = stepLimitLeft;
278 G4double tPathLength =
279 mscModel->ComputeTruePathLengthLimit(*currentTrackPtr, gPathLength);
280 G4bool mscLimitsStep = (tPathLength < stepLimitLeft);
281 if(!fMultipleSteps && mscLimitsStep)
282 {
283 // MSC limits the step.
284 *selection = CandidateForSelection;
285 }
286
287 if(!firstStep)
288 {
289 // Move the navigator to where the previous step ended.
292 }
293
294 G4GPILSelection transportSelection;
295 G4double geometryStepLength =
297 *currentTrackPtr, previousStepSize, gPathLength, currentSafety,
298 &transportSelection);
299 if(geometryStepLength < gPathLength)
300 {
301 // Transportation limits the step, ie the track hit a boundary.
302 *selection = CandidateForSelection;
303 continueStepping = false;
304 }
305 if(fTransportEndKineticEnergy != currentEnergy)
306 {
307 // Field propagation changed the energy, it's not possible to
308 // estimate the continuous energy loss and continue stepping.
309 continueStepping = false;
310 }
311
312 if(firstStep)
313 {
314 proposedSafety = currentSafety;
315 }
316 totalGeometryStepLength += geometryStepLength;
317
318 // Sample MSC direction change and displacement.
319 const G4double range =
320 mscModel->GetRange(particleDefinition, currentEnergy, couple);
321
322 tPathLength = mscModel->ComputeTrueStepLength(geometryStepLength);
323
324 // Protect against wrong t->g->t conversion.
325 tPathLength = std::min(tPathLength, stepLimitLeft);
326
327 totalTruePathLength += tPathLength;
328 if(*selection != CandidateForSelection && !mscLimitsStep)
329 {
330 // If neither MSC nor transportation limits the step, we got the
331 // distance we want - make sure we exit the loop.
332 continueStepping = false;
333 }
334 else if(tPathLength >= range)
335 {
336 // The particle will stop, exit the loop.
337 continueStepping = false;
338 }
339 else
340 {
341 stepLimitLeft -= tPathLength;
342 }
343
344 // Do not sample scattering at the last or at a small step.
345 if(tPathLength < range && tPathLength > kGeomMin)
346 {
347 static constexpr G4double minSafety = 1.20 * CLHEP::nm;
348 static constexpr G4double sFact = 0.99;
349
350 // The call to SampleScattering() *may* directly fill in the changed
351 // direction into fParticleChange, so we have to:
352 // 1) Make sure the momentum direction is initialized.
354 // 2) Call SampleScattering(), which *may* change it.
355 const G4ThreeVector displacement =
356 mscModel->SampleScattering(fTransportEndMomentumDir, minSafety);
357 // 3) Get the changed direction and inform G4Transportation.
358 fMomentumChanged = true;
360
361 const G4double r2 = displacement.mag2();
362 if(r2 > kMinDisplacement2)
363 {
364 G4bool positionChanged = true;
365 G4double dispR = std::sqrt(r2);
366 G4double postSafety = sFact * fpSafetyHelper->ComputeSafety(
367 fTransportEndPosition, dispR);
368
369 // Far away from geometry boundary
370 if(postSafety > 0.0 && dispR <= postSafety)
371 {
372 fTransportEndPosition += displacement;
373
374 // Near the boundary
375 }
376 else
377 {
378 // displaced point is definitely within the volume
379 if(dispR < postSafety)
380 {
381 fTransportEndPosition += displacement;
382
383 // reduced displacement
384 }
385 else if(postSafety > kGeomMin)
386 {
387 fTransportEndPosition += displacement * (postSafety / dispR);
388
389 // very small postSafety
390 }
391 else
392 {
393 positionChanged = false;
394 }
395 }
396 if(positionChanged)
397 {
399 }
400 }
401 }
402
403 if(continueStepping)
404 {
405 // Update safety according to the geometry distance.
406 if(currentSafety < fEndPointDistance)
407 {
408 currentSafety = 0;
409 }
410 else
411 {
412 currentSafety -= fEndPointDistance;
413 }
414
415 // Update the kinetic energy according to the continuous loss.
416 currentEnergy = mscModel->GetEnergy(particleDefinition,
417 range - tPathLength, couple);
418
419 // From now on, use the track that we can update below.
420 currentTrackPtr = fSubStepTrack;
421
422 fSubStepDynamicParticle->SetKineticEnergy(currentEnergy);
423 fSubStepDynamicParticle->SetMomentumDirection(
425 fSubStepTrack->SetPosition(fTransportEndPosition);
426
427 G4StepPoint& subPreStepPoint = *fSubStep->GetPreStepPoint();
428 subPreStepPoint.SetMaterialCutsCouple(couple);
429 subPreStepPoint.SetPosition(fTransportEndPosition);
430 subPreStepPoint.SetSafety(currentSafety);
431 subPreStepPoint.SetStepStatus(fAlongStepDoItProc);
432 }
433 firstStep = false;
434 } while(continueStepping);
435
436 // Note: currentEnergy is only updated if continueStepping is true.
437 // In case field propagation changed the energy, this flag is
438 // immediately set to false and currentEnergy is still equal to the
439 // initial kinetic energy stored in ekin.
440 if(currentEnergy != ekin)
441 {
442 // If field propagation didn't change the energy and we potentially
443 // did multiple steps, reset the energy that G4Transportation will
444 // propose to not subtract the energy loss twice.
446 // Also ask for the range again with the initial energy so it is
447 // correctly cached in the G4VEnergyLossProcess.
448 // FIXME: Asking for a range should never change the cached values!
449 (void) mscModel->GetRange(particleDefinition, ekin, couple);
450 }
451
452 fParticleChange.ProposeTrueStepLength(totalTruePathLength);
453
454 return totalGeometryStepLength;
455 }
456 }
457 }
458
459 // If we get here, no scattering has happened.
461 track, previousStepSize, currentMinimumStep, proposedSafety, selection);
462}
463
464//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
@ fAlongStepDoItProc
Definition: G4StepStatus.hh:47
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
void SetMomentumDirection(const G4ThreeVector &aDirection)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(G4double aEnergy)
static G4Electron * Definition()
Definition: G4Electron.cc:48
static G4Electron * Electron()
Definition: G4Electron.cc:93
void PrepareModels(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fm, const G4Region *r)
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4int verb)
G4VEmModel * SelectModel(G4double energy, std::size_t index)
void DumpModelList(std::ostream &out, G4int verb)
G4int NumberOfModels() const
G4VEmModel * GetModel(G4int idx, G4bool ver=false) const
G4bool IsPrintLocked() const
static G4EmParameters * Instance()
static G4LossTableManager * Instance()
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
G4LossTableBuilder * GetTableBuilder()
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
G4EmConfigurator * EmConfigurator()
virtual void LocateGlobalPointWithinVolume(const G4ThreeVector &position)
Definition: G4Navigator.cc:600
const G4ThreeVector * GetMomentumDirection() const
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
G4double ComputeSafety(const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)
void ReLocateWithinVolume(const G4ThreeVector &pGlobalPoint)
void SetSafety(const G4double aValue)
void SetStepStatus(const G4StepStatus aValue)
void SetMaterialCutsCouple(const G4MaterialCutsCouple *)
void SetPosition(const G4ThreeVector &aValue)
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
void SetPosition(const G4ThreeVector &aValue)
void SetStep(const G4Step *aValue)
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void AddMscModel(G4VMscModel *mscModel, G4int order=0, const G4Region *region=nullptr)
void BuildPhysicsTable(const G4ParticleDefinition &part) override
void StartTracking(G4Track *track) override
G4TransportationWithMsc(ScatteringType type, G4int verbosity=0)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection) override
void PreparePhysicsTable(const G4ParticleDefinition &part) override
G4double AlongStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
G4double fTransportEndKineticEnergy
G4ThreeVector fTransportEndPosition
G4ParticleChangeForTransport fParticleChange
void StartTracking(G4Track *aTrack)
G4Navigator * fLinearNavigator
G4ThreeVector fTransportEndMomentumDir
G4SafetyHelper * fpSafetyHelper
void SetCrossSectionTable(G4PhysicsTable *, G4bool isLocal)
Definition: G4VEmModel.cc:398
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:718
void SetCurrentCouple(const G4MaterialCutsCouple *)
Definition: G4VEmModel.hh:468
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:390
G4bool IsActive(G4double kinEnergy) const
Definition: G4VEmModel.hh:774
virtual void StartTracking(G4Track *)
Definition: G4VEmModel.cc:214
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)=0
virtual G4double ComputeTrueStepLength(G4double geomPathLength)=0
G4double GetEnergy(const G4ParticleDefinition *part, G4double range, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:223
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
Definition: G4VMscModel.cc:188
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)=0
void ProposeTrueStepLength(G4double truePathLength)
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:416
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:522
G4int verboseLevel
Definition: G4VProcess.hh:360
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386