Geant4 11.2.2
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// 16-07-03 Use G4VMscModel interface (V.Ivanchenko)
41// 03-11-03 Fix initialisation problem in RetrievePhysicsTable (V.Ivanchenko)
42// 04-11-03 Update PrintInfoDefinition (V.Ivanchenko)
43// 01-03-04 SampleCosineTheta signature changed
44// 22-04-04 SampleCosineTheta signature changed back to original
45// 27-08-04 Add InitialiseForRun method (V.Ivanchneko)
46// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
47// 11-03-05 Shift verbose level by 1 (V.Ivantchenko)
48// 15-04-05 optimize internal interface (V.Ivanchenko)
49// 15-04-05 remove boundary flag (V.Ivanchenko)
50// 27-10-05 introduce virtual function MscStepLimitation() (V.Ivanchenko)
51// 12-04-07 Add verbosity at destruction (V.Ivanchenko)
52// 27-10-07 Virtual functions moved to source (V.Ivanchenko)
53// 11-03-08 Set skin value does not effect step limit type (V.Ivanchenko)
54// 24-06-09 Removed hidden bin in G4PhysicsVector (V.Ivanchenko)
55// 04-06-13 Adoptation to MT mode (V.Ivanchenko)
56//
57
58// -------------------------------------------------------------------
59//
60//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62
65#include "G4SystemOfUnits.hh"
66#include "G4LossTableManager.hh"
68#include "G4Step.hh"
71#include "G4UnitsTable.hh"
73#include "G4Electron.hh"
74#include "G4GenericIon.hh"
76#include "G4SafetyHelper.hh"
77#include "G4ParticleTable.hh"
78#include "G4ProcessVector.hh"
79#include "G4ProcessManager.hh"
80#include "G4LossTableBuilder.hh"
81#include "G4EmTableUtil.hh"
82#include <iostream>
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85
88 fNewPosition(0.,0.,0.),
89 fNewDirection(0.,0.,1.)
90{
91 theParameters = G4EmParameters::Instance();
94
95 lowestKinEnergy = 10*CLHEP::eV;
96
97 geomMin = 0.05*CLHEP::nm;
98 minDisplacement2 = geomMin*geomMin;
99
101
102 modelManager = new G4EmModelManager();
103 emManager = G4LossTableManager::Instance();
104 mscModels.reserve(2);
105 emManager->Register(this);
106}
107
108//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
109
111{
112 delete modelManager;
113 emManager->DeRegister(this);
114}
115
116//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
117
119 const G4Region* region)
120{
121 if(nullptr == ptr) { return; }
122 G4VEmFluctuationModel* fm = nullptr;
123 modelManager->AddEmModel(order, ptr, fm, region);
125}
126
127//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
128
130{
131 if(nullptr == ptr) { return; }
132 if(!mscModels.empty()) {
133 for(auto & msc : mscModels) { if(msc == ptr) { return; } }
134 }
135 mscModels.push_back(ptr);
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
139
140void
142{
143 G4bool master = emManager->IsMaster();
144 if(nullptr == firstParticle) { firstParticle = &part; }
145
146 emManager->PreparePhysicsTable(&part, this);
147 currParticle = nullptr;
148
149 if(firstParticle == &part) {
150 baseMat = emManager->GetTableBuilder()->GetBaseMaterialFlag();
151
152 G4EmTableUtil::PrepareMscProcess(this, part, modelManager,
153 stepLimit, facrange,
154 latDisplacement, master,
155 isIon, baseMat);
156
157 numberOfModels = modelManager->NumberOfModels();
158 currentModel = GetModelByIndex(0);
159
160 if(nullptr == safetyHelper) {
162 ->GetSafetyHelper();
163 safetyHelper->InitialiseHelper();
164 }
165 }
166}
167
168//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
169
171{
172 G4bool master = emManager->IsMaster();
173
174 if(firstParticle == &part) {
175 emManager->BuildPhysicsTable(firstParticle);
176 }
177 const G4VMultipleScattering* ptr = this;
178 if(!master) {
179 ptr = static_cast<const G4VMultipleScattering*>(GetMasterProcess());
180 }
181
182 G4EmTableUtil::BuildMscProcess(this, ptr, part, firstParticle,
183 numberOfModels, master);
184}
185
186//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
187
188void G4VMultipleScattering::StreamInfo(std::ostream& outFile,
189 const G4ParticleDefinition& part, G4bool rst) const
190{
191 G4String indent = (rst ? " " : "");
192 outFile << G4endl << indent << GetProcessName() << ": ";
193 if (!rst) outFile << " for " << part.GetParticleName();
194 outFile << " SubType= " << GetProcessSubType() << G4endl;
195 modelManager->DumpModelList(outFile, verboseLevel);
196}
197
198//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
199
201{
202 G4VEnergyLossProcess* eloss = nullptr;
203 if(track->GetParticleDefinition() != currParticle) {
204 currParticle = track->GetParticleDefinition();
205 fIonisation = emManager->GetEnergyLossProcess(currParticle);
206 eloss = fIonisation;
207 }
208 for(G4int i=0; i<numberOfModels; ++i) {
210 msc->StartTracking(track);
211 if(nullptr != eloss) {
212 msc->SetIonisation(eloss, currParticle);
213 }
214 }
215}
216
217//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
218
220 const G4Track& track,
221 G4double,
222 G4double currentMinimalStep,
223 G4double&,
224 G4GPILSelection* selection)
225{
226 // get Step limit proposed by the process
227 *selection = NotCandidateForSelection;
228 physStepLimit = gPathLength = tPathLength = currentMinimalStep;
229
230 G4double ekin = track.GetKineticEnergy();
231 /*
232 G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
233 << " " << currParticle->GetParticleName()
234 << " currMod " << currentModel
235 << G4endl;
236 */
237 // isIon flag is used only to select a model
238 if(isIon) {
239 ekin *= proton_mass_c2/track.GetParticleDefinition()->GetPDGMass();
240 }
241 const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
242
243 // select new model, static cast is possible in this class
244 if(1 < numberOfModels) {
245 currentModel =
246 static_cast<G4VMscModel*>(SelectModel(ekin,couple->GetIndex()));
247 }
248 currentModel->SetCurrentCouple(couple);
249 // msc is active is model is active, energy above the limit,
250 // and step size is above the limit;
251 // if it is active msc may limit the step
252 if(currentModel->IsActive(ekin) && tPathLength > geomMin
253 && ekin >= lowestKinEnergy) {
254 isActive = true;
255 tPathLength =
256 currentModel->ComputeTruePathLengthLimit(track, gPathLength);
257 if (tPathLength < physStepLimit) {
258 *selection = CandidateForSelection;
259 }
260 } else {
261 isActive = false;
262 gPathLength = DBL_MAX;
263 }
264
265 //if(currParticle->GetPDGMass() > GeV)
266 /*
267 G4cout << "MSC::AlongStepGPIL: Ekin= " << ekin
268 << " " << currParticle->GetParticleName()
269 << " gPathLength= " << gPathLength
270 << " tPathLength= " << tPathLength
271 << " currentMinimalStep= " << currentMinimalStep
272 << " isActive " << isActive << G4endl;
273 */
274 return gPathLength;
275}
276
277//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
278
286
287//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
288
291{
292 fParticleChange.InitialiseMSC(track, step);
293 fNewPosition = fParticleChange.GetProposedPosition();
294 fPositionChanged = false;
295
296 G4double geomLength = step.GetStepLength();
297
298 // very small step - no msc
299 if(!isActive) {
300 tPathLength = geomLength;
301
302 // sample msc
303 } else {
304 G4double range =
305 currentModel->GetRange(currParticle,track.GetKineticEnergy(),
306 track.GetMaterialCutsCouple());
307
308 tPathLength = currentModel->ComputeTrueStepLength(geomLength);
309
310 /*
311 if(currParticle->GetPDGMass() > 0.9*GeV)
312 G4cout << "G4VMsc::AlongStepDoIt: GeomLength= "
313 << geomLength
314 << " tPathLength= " << tPathLength
315 << " physStepLimit= " << physStepLimit
316 << " dr= " << range - tPathLength
317 << " ekin= " << track.GetKineticEnergy() << G4endl;
318 */
319 // protection against wrong t->g->t conversion
320 tPathLength = std::min(tPathLength, physStepLimit);
321
322 // do not sample scattering at the last or at a small step
323 if(tPathLength < range && tPathLength > geomMin) {
324
325 static const G4double minSafety = 1.20*CLHEP::nm;
326 static const G4double sFact = 0.99;
327
328 G4ThreeVector displacement = currentModel->SampleScattering(
329 step.GetPostStepPoint()->GetMomentumDirection(),minSafety);
330
331 G4double r2 = displacement.mag2();
332 //G4cout << " R= " << sqrt(r2) << " Rmin= " << sqrt(minDisplacement2)
333 // << " flag= " << fDispBeyondSafety << G4endl;
334 if(r2 > minDisplacement2) {
335
336 fPositionChanged = true;
337 G4double dispR = std::sqrt(r2);
338 G4double postSafety =
339 sFact*safetyHelper->ComputeSafety(fNewPosition, dispR);
340 //G4cout<<" R= "<< dispR<<" postSafety= "<<postSafety<<G4endl;
341
342 // far away from geometry boundary
343 if(postSafety > 0.0 && dispR <= postSafety) {
344 fNewPosition += displacement;
345
346 //near the boundary
347 } else {
348 // displaced point is definitely within the volume
349 //G4cout<<" R= "<<dispR<<" postSafety= "<<postSafety<<G4endl;
350 if(dispR < postSafety) {
351 fNewPosition += displacement;
352
353 // reduced displacement
354 } else if(postSafety > geomMin) {
355 fNewPosition += displacement*(postSafety/dispR);
356
357 // very small postSafety
358 } else {
359 fPositionChanged = false;
360 }
361 }
362 if(fPositionChanged) {
363 safetyHelper->ReLocateWithinVolume(fNewPosition);
364 fParticleChange.ProposePosition(fNewPosition);
365 }
366 }
367 }
368 }
370 return &fParticleChange;
371}
372
373//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
374
376 const G4Track& track,
377 G4double previousStepSize,
378 G4double currentMinimalStep,
379 G4double& currentSafety)
380{
382 G4double x = AlongStepGetPhysicalInteractionLength(track,previousStepSize,
383 currentMinimalStep,
384 currentSafety,
385 &selection);
386 return x;
387}
388
389//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
390
392 const G4Track& track,
393 G4double previousStepSize,
394 G4double currentMinimalStep,
395 G4double& currentSafety)
396{
397 return GetContinuousStepLimit(track,previousStepSize,currentMinimalStep,
398 currentSafety);
399}
400
401//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
402
409
410//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
411
412G4bool
414 const G4String& directory,
415 G4bool ascii)
416{
417 G4bool yes = true;
418 if(part != firstParticle || !emManager->IsMaster()) { return yes; }
419
420 return G4EmTableUtil::StoreMscTable(this, part, directory,
421 numberOfModels, verboseLevel,
422 ascii);
423}
424
425//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
426
427G4bool
429 const G4String&,
430 G4bool)
431{
432 return true;
433}
434
435//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
436
437void G4VMultipleScattering::ProcessDescription(std::ostream& outFile) const
438{
439 if(nullptr != firstParticle) {
440 StreamInfo(outFile, *firstParticle, true);
441 }
442}
443
444//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
445
@ fMultipleScattering
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
@ NotForced
@ Forced
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
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:67
double mag2() const
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fm, const G4Region *r)
void DumpModelList(std::ostream &out, G4int verb)
G4int NumberOfModels() const
static G4EmParameters * Instance()
static void PrepareMscProcess(G4VMultipleScattering *proc, const G4ParticleDefinition &part, G4EmModelManager *modelManager, G4MscStepLimitType &stepLimit, G4double &facrange, G4bool &latDisplacement, G4bool &master, G4bool &isIon, G4bool &baseMat)
static void BuildMscProcess(G4VMultipleScattering *proc, const G4VMultipleScattering *masterProc, const G4ParticleDefinition &part, const G4ParticleDefinition *firstPart, G4int nModels, G4bool master)
static G4bool StoreMscTable(G4VMultipleScattering *proc, const G4ParticleDefinition *part, const G4String &directory, const G4int nModels, const G4int verb, const G4bool ascii)
static G4LossTableManager * Instance()
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
G4LossTableBuilder * GetTableBuilder()
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
void DeRegister(G4VEnergyLossProcess *p)
void Register(G4VEnergyLossProcess *p)
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
const G4ThreeVector & GetProposedPosition() const
void InitialiseMSC(const G4Track &, const G4Step &step)
void ProposePosition(const G4ThreeVector &finalPosition)
const G4String & GetParticleName() const
G4double ComputeSafety(const G4ThreeVector &pGlobalPoint, G4double maxRadius=DBL_MAX)
void ReLocateWithinVolume(const G4ThreeVector &pGlobalPoint)
const G4ThreeVector & GetMomentumDirection() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
void SetCurrentCouple(const G4MaterialCutsCouple *)
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
G4bool IsActive(G4double kinEnergy) const
virtual void StartTracking(G4Track *)
virtual G4double ComputeTruePathLengthLimit(const G4Track &track, G4double &stepLimit)=0
void SetIonisation(G4VEnergyLossProcess *, const G4ParticleDefinition *part)
virtual G4double ComputeTrueStepLength(G4double geomPathLength)=0
G4double GetRange(const G4ParticleDefinition *part, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
virtual G4ThreeVector & SampleScattering(const G4ThreeVector &, G4double safety)=0
void AddEmModel(G4int order, G4VMscModel *, const G4Region *region=nullptr)
G4double GetMeanFreePath(const G4Track &track, G4double, G4ForceCondition *condition) override
G4ParticleChangeForMSC fParticleChange
G4VMscModel * GetModelByIndex(G4int idx, G4bool ver=false) const
void StartTracking(G4Track *) override
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 ProcessDescription(std::ostream &outFile) const override
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
void SetEmModel(G4VMscModel *, G4int idx=0)
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimalStep, G4double &currentSafety, G4GPILSelection *selection) override
void StreamInfo(std::ostream &outFile, const G4ParticleDefinition &, G4bool rst=false) const
G4double PostStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4ForceCondition *condition) override
G4VParticleChange * AlongStepDoIt(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
void ProposeTrueStepLength(G4double truePathLength)
void SetVerboseLevel(G4int value)
const G4VProcess * GetMasterProcess() const
G4int verboseLevel
void SetProcessSubType(G4int)
G4VParticleChange * pParticleChange
G4int GetProcessSubType() const
const G4String & GetProcessName() const
#define DBL_MAX
Definition templates.hh:62