Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4VEnergyLossProcess.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: G4VEnergyLossProcess.cc 79268 2014-02-20 16:46:31Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4VEnergyLossProcess
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 03.01.2002
38 //
39 // Modifications:
40 //
41 // 13-11-02 Minor fix - use normalised direction (V.Ivanchenko)
42 // 04-12-02 Minor change in PostStepDoIt (V.Ivanchenko)
43 // 23-12-02 Change interface in order to move to cut per region (V.Ivanchenko)
44 // 26-12-02 Secondary production moved to derived classes (V.Ivanchenko)
45 // 04-01-03 Fix problem of very small steps for ions (V.Ivanchenko)
46 // 20-01-03 Migrade to cut per region (V.Ivanchenko)
47 // 24-01-03 Temporarily close a control on usage of couples (V.Ivanchenko)
48 // 24-01-03 Make models region aware (V.Ivanchenko)
49 // 05-02-03 Fix compilation warnings (V.Ivanchenko)
50 // 06-02-03 Add control on tmax in PostStepDoIt (V.Ivanchenko)
51 // 13-02-03 SubCutoffProcessors defined for regions (V.Ivanchenko)
52 // 15-02-03 Lambda table can be scaled (V.Ivanchenko)
53 // 17-02-03 Fix problem of store/restore tables (V.Ivanchenko)
54 // 18-02-03 Add control on CutCouple usage (V.Ivanchenko)
55 // 26-02-03 Simplify control on GenericIons (V.Ivanchenko)
56 // 06-03-03 Control on GenericIons using SubType+ update verbose (V.Ivanchenko)
57 // 10-03-03 Add Ion registration (V.Ivanchenko)
58 // 22-03-03 Add Initialisation of cash (V.Ivanchenko)
59 // 26-03-03 Remove finalRange modification (V.Ivanchenko)
60 // 09-04-03 Fix problem of negative range limit for non integral (V.Ivanchenko)
61 // 26-04-03 Fix retrieve tables (V.Ivanchenko)
62 // 06-05-03 Set defalt finalRange = 1 mm (V.Ivanchenko)
63 // 12-05-03 Update range calculations + lowKinEnergy (V.Ivanchenko)
64 // 13-05-03 Add calculation of precise range (V.Ivanchenko)
65 // 23-05-03 Remove tracking cuts (V.Ivanchenko)
66 // 03-06-03 Fix initialisation problem for STD ionisation (V.Ivanchenko)
67 // 21-07-03 Add UpdateEmModel method (V.Ivanchenko)
68 // 03-11-03 Fix initialisation problem in RetrievePhysicsTable (V.Ivanchenko)
69 // 04-11-03 Add checks in RetrievePhysicsTable (V.Ivanchenko)
70 // 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
71 // 21-01-04 Migrade to G4ParticleChangeForLoss (V.Ivanchenko)
72 // 27-02-04 Fix problem of loss in low presure gases, cleanup precise range
73 // calculation, use functions ForLoss in AlongStepDoIt (V.Ivanchenko)
74 // 10-03-04 Fix a problem of Precise Range table (V.Ivanchenko)
75 // 19-03-04 Fix a problem energy below lowestKinEnergy (V.Ivanchenko)
76 // 31-03-04 Fix a problem of retrieve tables (V.Ivanchenko)
77 // 21-07-04 Check weather AtRest are active or not (V.Ivanchenko)
78 // 03-08-04 Add pointer of DEDX table to all processes (V.Ivanchenko)
79 // 06-08-04 Clear up names of member functions (V.Ivanchenko)
80 // 06-08-04 Clear up names of member functions (V.Ivanchenko)
81 // 27-08-04 Add NeedBuildTables method (V.Ivanchneko)
82 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
83 // 11-03-05 Shift verbose level by 1 (V.Ivantchenko)
84 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
85 // 11-04-05 Use MaxSecondaryEnergy from a model (V.Ivanchenko)
86 // 25-07-05 Add extra protection PostStep for non-integral mode (V.Ivanchenko)
87 // 12-08-05 Integral=false; SetStepFunction(0.2, 0.1*mm) (mma)
88 // 18-08-05 Return back both AlongStep and PostStep from 7.0 (V.Ivanchenko)
89 // 02-09-05 Default StepFunction 0.2 1 mm + integral (V.Ivanchenko)
90 // 04-09-05 default lambdaFactor 0.8 (V.Ivanchenko)
91 // 05-10-05 protection against 0 energy loss added (L.Urban)
92 // 17-10-05 protection above has been removed (L.Urban)
93 // 06-01-06 reset currentCouple when StepFunction is changed (V.Ivanchenko)
94 // 10-01-06 PreciseRange -> CSDARange (V.Ivantchenko)
95 // 18-01-06 Clean up subcutoff including recalculation of presafety (VI)
96 // 20-01-06 Introduce G4EmTableType and reducing number of methods (VI)
97 // 22-03-06 Add control on warning printout AlongStep (VI)
98 // 23-03-06 Use isIonisation flag (V.Ivanchenko)
99 // 07-06-06 Do not reflect AlongStep in subcutoff regime (V.Ivanchenko)
100 // 14-01-07 add SetEmModel(index) and SetFluctModel() (mma)
101 // 16-01-07 add IonisationTable and IonisationSubTable (V.Ivanchenko)
102 // 16-02-07 set linLossLimit=1.e-6 (V.Ivanchenko)
103 // 13-03-07 use SafetyHelper instead of navigator (V.Ivanchenko)
104 // 10-04-07 use unique SafetyHelper (V.Ivanchenko)
105 // 12-04-07 Add verbosity at destruction (V.Ivanchenko)
106 // 25-04-07 move initialisation of safety helper to BuildPhysicsTable (VI)
107 // 27-10-07 Virtual functions moved to source (V.Ivanchenko)
108 // 24-06-09 Removed hidden bin in G4PhysicsVector (V.Ivanchenko)
109 // 01-25-09 (Xin Dong) Phase II change for Geant4 multi-threading:
110 // New methods SlavePreparePhysicsTable, SlaveBuildPhysicsTable
111 // Worker threads share physics tables with the master thread for
112 // this kind of process. This member function is used by worker
113 // threads to achieve the partial effect of the master thread when
114 // it builds physcis tables.
115 // 15-10-10 Fixed 4-momentum balance if deexcitation is active (L.Pandola)
116 // 30-05-12 Call new ApplySecondaryBiasing so 2ries may be unique (D. Sawkey)
117 // 30-05-12 Fix bug in forced biasing: now called on first step (D. Sawkey)
118 // 04-06-13 Adoptation to MT mode, adding internal cache to GetRangeForLoss,
119 // more accurate initialisation for ions (V.Ivanchenko)
120 //
121 // Class Description:
122 //
123 // It is the unified energy loss process it calculates the continuous
124 // energy loss for charged particles using a set of Energy Loss
125 // models valid for different energy regions. There are a possibility
126 // to create and access to dE/dx and range tables, or to calculate
127 // that information on fly.
128 // -------------------------------------------------------------------
129 //
130 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
131 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
132 
133 #include "G4VEnergyLossProcess.hh"
134 #include "G4PhysicalConstants.hh"
135 #include "G4SystemOfUnits.hh"
136 #include "G4ProcessManager.hh"
137 #include "G4LossTableManager.hh"
138 #include "G4LossTableBuilder.hh"
139 #include "G4Step.hh"
140 #include "G4ParticleDefinition.hh"
141 #include "G4ParticleTable.hh"
142 #include "G4VEmModel.hh"
143 #include "G4VEmFluctuationModel.hh"
144 #include "G4DataVector.hh"
145 #include "G4PhysicsLogVector.hh"
146 #include "G4VParticleChange.hh"
147 #include "G4Gamma.hh"
148 #include "G4Electron.hh"
149 #include "G4Positron.hh"
150 #include "G4ProcessManager.hh"
151 #include "G4UnitsTable.hh"
152 #include "G4ProductionCutsTable.hh"
153 #include "G4Region.hh"
154 #include "G4RegionStore.hh"
155 #include "G4PhysicsTableHelper.hh"
156 #include "G4SafetyHelper.hh"
158 #include "G4EmConfigurator.hh"
159 #include "G4VAtomDeexcitation.hh"
160 #include "G4EmBiasingManager.hh"
161 #include "G4Log.hh"
162 
163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
164 
166  G4ProcessType type):
167  G4VContinuousDiscreteProcess(name, type),
168  secondaryParticle(0),
169  nSCoffRegions(0),
170  idxSCoffRegions(0),
171  nProcesses(0),
172  theDEDXTable(0),
173  theDEDXSubTable(0),
174  theDEDXunRestrictedTable(0),
175  theIonisationTable(0),
176  theIonisationSubTable(0),
177  theRangeTableForLoss(0),
178  theCSDARangeTable(0),
179  theSecondaryRangeTable(0),
180  theInverseRangeTable(0),
181  theLambdaTable(0),
182  theSubLambdaTable(0),
183  theDensityFactor(0),
184  theDensityIdx(0),
185  baseParticle(0),
186  minSubRange(0.1),
187  lossFluctuationFlag(true),
188  rndmStepFlag(false),
189  tablesAreBuilt(false),
190  integral(true),
191  isIon(false),
192  isIonisation(true),
193  useSubCutoff(false),
194  useDeexcitation(false),
195  particle(0),
196  currentCouple(0),
197  nWarnings(0),
198  mfpKinEnergy(0.0)
199 {
200  SetVerboseLevel(1);
201 
202  // low energy limit
203  lowestKinEnergy = 1.*eV;
204  preStepKinEnergy = 0.0;
205  preStepRangeEnergy = 0.0;
206  computedRange = DBL_MAX;
207 
208  // Size of tables assuming spline
209  minKinEnergy = 0.1*keV;
210  maxKinEnergy = 10.0*TeV;
211  nBins = 77;
212  maxKinEnergyCSDA = 1.0*GeV;
213  nBinsCSDA = 35;
214 
215  // default linear loss limit for spline
216  linLossLimit = 0.01;
217 
218  // default dRoverRange and finalRange
219  SetStepFunction(0.2, 1.0*mm);
220 
221  // default lambda factor
222  lambdaFactor = 0.8;
223 
224  // cross section biasing
225  biasFactor = 1.0;
226 
227  // particle types
228  theElectron = G4Electron::Electron();
229  thePositron = G4Positron::Positron();
230  theGamma = G4Gamma::Gamma();
231  theGenericIon = 0;
232 
233  // run time objects
236  modelManager = new G4EmModelManager();
238  ->GetSafetyHelper();
239  aGPILSelection = CandidateForSelection;
240 
241  // initialise model
242  lManager = G4LossTableManager::Instance();
243  lManager->Register(this);
244  fluctModel = 0;
245  atomDeexcitation = 0;
246 
247  biasManager = 0;
248  biasFlag = false;
249  weightFlag = false;
250  isMaster = true;
251  lastIdx = 0;
252 
253  idxDEDX = idxDEDXSub = idxDEDXunRestricted = idxIonisation =
254  idxIonisationSub = idxRange = idxCSDA = idxSecRange =
255  idxInverseRange = idxLambda = idxSubLambda = 0;
256 
257  scTracks.reserve(5);
258  secParticles.reserve(5);
259 
260  secID = biasID = subsecID = -1;
261 }
262 
263 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
264 
266 {
267  /*
268  G4cout << "** G4VEnergyLossProcess::~G4VEnergyLossProcess() for "
269  << GetProcessName()
270  << " isMaster: " << isMaster << G4endl;
271  */
272  Clean();
273 
274  if ( !baseParticle && isMaster ) {
275  //G4cout << " isIonisation " << isIonisation << " "
276  // << theDEDXTable << G4endl;
277 
278  if(theDEDXTable) {
279  if(theIonisationTable == theDEDXTable) { theIonisationTable = 0; }
280  delete theDEDXTable;
281  if(theDEDXSubTable) {
282  if(theIonisationSubTable == theDEDXSubTable)
283  { theIonisationSubTable = 0; }
284  delete theDEDXSubTable;
285  }
286  }
287  delete theIonisationTable;
288  delete theIonisationSubTable;
289  if(theDEDXunRestrictedTable && isIonisation) {
290  delete theDEDXunRestrictedTable;
291  }
292  if(theCSDARangeTable && isIonisation) {
293  delete theCSDARangeTable;
294  }
295  if(theRangeTableForLoss && isIonisation) {
296  delete theRangeTableForLoss;
297  }
298  if(theInverseRangeTable && isIonisation) {
299  delete theInverseRangeTable;
300  }
301  delete theLambdaTable;
302  delete theSubLambdaTable;
303  }
304 
305  delete modelManager;
306  delete biasManager;
307  lManager->DeRegister(this);
308  //G4cout << "** all removed" << G4endl;
309 }
310 
311 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
312 
313 void G4VEnergyLossProcess::Clean()
314 {
315  /*
316  if(1 < verboseLevel) {
317  G4cout << "G4VEnergyLossProcess::Clear() for " << GetProcessName()
318  << G4endl;
319  }
320  */
321  delete [] idxSCoffRegions;
322 
323  tablesAreBuilt = false;
324 
325  scProcesses.clear();
326  nProcesses = 0;
327 
328  idxDEDX = idxDEDXSub = idxDEDXunRestricted = idxIonisation =
329  idxIonisationSub = idxRange = idxCSDA = idxSecRange =
330  idxInverseRange = idxLambda = idxSubLambda = 0;
331 }
332 
333 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
334 
336  const G4Material*,
337  G4double cut)
338 {
339  return cut;
340 }
341 
342 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
343 
345  G4VEmFluctuationModel* fluc,
346  const G4Region* region)
347 {
348  modelManager->AddEmModel(order, p, fluc, region);
349  if(p) { p->SetParticleChange(pParticleChange, fluc); }
350 }
351 
352 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
353 
355  G4double emin, G4double emax)
356 {
357  modelManager->UpdateEmModel(nam, emin, emax);
358 }
359 
360 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
361 
363 {
364  G4int n = emModels.size();
365  if(index >= n) { for(G4int i=n; i<=index; ++i) {emModels.push_back(0);} }
366  emModels[index] = p;
367 }
368 
369 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
370 
372 {
373  G4VEmModel* p = 0;
374  if(index >= 0 && index < G4int(emModels.size())) { p = emModels[index]; }
375  return p;
376 }
377 
378 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
379 
381 {
382  return modelManager->GetModel(idx, ver);
383 }
384 
385 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
386 
388 {
389  return modelManager->NumberOfModels();
390 }
391 
392 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
393 
394 void
396 {
397  if(1 < verboseLevel) {
398  G4cout << "G4VEnergyLossProcess::PreparePhysicsTable for "
399  << GetProcessName() << " for " << part.GetParticleName()
400  << " " << this << G4endl;
401  }
402 
403  if(GetMasterProcess() != this) { isMaster = false; }
404 
405  currentCouple = 0;
406  preStepLambda = 0.0;
407  mfpKinEnergy = DBL_MAX;
408  fRange = DBL_MAX;
409  preStepKinEnergy = 0.0;
410  preStepRangeEnergy = 0.0;
411  chargeSqRatio = 1.0;
412  massRatio = 1.0;
413  reduceFactor = 1.0;
414  fFactor = 1.0;
415  lastIdx = 0;
416 
417  // Are particle defined?
418  if( !particle ) { particle = &part; }
419 
420  if(part.GetParticleType() == "nucleus") {
421 
422  G4String pname = part.GetParticleName();
423  if(pname != "deuteron" && pname != "triton" &&
424  pname != "alpha+" && pname != "helium" &&
425  pname != "hydrogen") {
426 
427  if(!theGenericIon) {
428  theGenericIon =
430  }
431  isIon = true;
432  if(theGenericIon && particle != theGenericIon) {
433  G4ProcessManager* pm = theGenericIon->GetProcessManager();
435  size_t n = v->size();
436  for(size_t j=0; j<n; ++j) {
437  if((*v)[j] == this) {
438  particle = theGenericIon;
439  break;
440  }
441  }
442  }
443  }
444  }
445 
446  if( particle != &part ) {
447  if(!isIon) {
448  lManager->RegisterExtraParticle(&part, this);
449  }
450  if(1 < verboseLevel) {
451  G4cout << "### G4VEnergyLossProcess::PreparePhysicsTable()"
452  << " interrupted for "
453  << part.GetParticleName() << " isIon= " << isIon
454  << " particle " << particle << " GenericIon " << theGenericIon
455  << G4endl;
456  }
457  return;
458  }
459 
460  Clean();
461  lManager->PreparePhysicsTable(&part, this, isMaster);
462  G4LossTableBuilder* bld = lManager->GetTableBuilder();
463 
464  // Base particle and set of models can be defined here
465  InitialiseEnergyLossProcess(particle, baseParticle);
466 
467  const G4ProductionCutsTable* theCoupleTable=
469  size_t n = theCoupleTable->GetTableSize();
470 
471  theDEDXAtMaxEnergy.resize(n, 0.0);
472  theRangeAtMaxEnergy.resize(n, 0.0);
473  theEnergyOfCrossSectionMax.resize(n, 0.0);
474  theCrossSectionMax.resize(n, DBL_MAX);
475 
476  // Tables preparation
477  if (isMaster && !baseParticle) {
478 
479  if(theDEDXTable && isIonisation) {
480  if(theIonisationTable && theDEDXTable != theIonisationTable) {
481  theDEDXTable->clearAndDestroy();
482  theDEDXTable = theIonisationTable;
483  }
484  if(theDEDXSubTable && theIonisationSubTable &&
485  theDEDXSubTable != theIonisationSubTable) {
486  theDEDXSubTable->clearAndDestroy();
487  theDEDXSubTable = theIonisationSubTable;
488  }
489  }
490 
491  theDEDXTable = G4PhysicsTableHelper::PreparePhysicsTable(theDEDXTable);
492  bld->InitialiseBaseMaterials(theDEDXTable);
493 
494  if(theDEDXSubTable) {
495  theDEDXSubTable =
497  }
498 
499  if (lManager->BuildCSDARange()) {
500  theDEDXunRestrictedTable =
501  G4PhysicsTableHelper::PreparePhysicsTable(theDEDXunRestrictedTable);
502  theCSDARangeTable =
504  }
505 
506  theLambdaTable = G4PhysicsTableHelper::PreparePhysicsTable(theLambdaTable);
507 
508  if(isIonisation) {
509  theRangeTableForLoss =
510  G4PhysicsTableHelper::PreparePhysicsTable(theRangeTableForLoss);
511  theInverseRangeTable =
512  G4PhysicsTableHelper::PreparePhysicsTable(theInverseRangeTable);
513  }
514 
515  if (nSCoffRegions) {
516  theDEDXSubTable =
518  theSubLambdaTable =
520  }
521  }
522 
523  // forced biasing
524  if(biasManager) {
525  biasManager->Initialise(part,GetProcessName(),verboseLevel);
526  biasFlag = false;
527  }
528 
529  G4double initialCharge = particle->GetPDGCharge();
530  G4double initialMass = particle->GetPDGMass();
531 
532  if (baseParticle) {
533  massRatio = (baseParticle->GetPDGMass())/initialMass;
534  G4double q = initialCharge/baseParticle->GetPDGCharge();
535  chargeSqRatio = q*q;
536  if(chargeSqRatio > 0.0) { reduceFactor = 1.0/(chargeSqRatio*massRatio); }
537  }
538 
539  // defined ID of secondary particles
540  if(isMaster) {
541  G4String nam1 = GetProcessName();
542  G4String nam4 = nam1 + "_split";
543  G4String nam5 = nam1 + "_subcut";
544  secID = G4PhysicsModelCatalog::Register(nam1);
545  biasID = G4PhysicsModelCatalog::Register(nam4);
546  subsecID= G4PhysicsModelCatalog::Register(nam5);
547  }
548 
549  // initialisation of models
550  G4int nmod = modelManager->NumberOfModels();
551  for(G4int i=0; i<nmod; ++i) {
552  G4VEmModel* mod = modelManager->GetModel(i);
553  mod->SetMasterThread(isMaster);
554  if(mod->HighEnergyLimit() > maxKinEnergy) {
555  mod->SetHighEnergyLimit(maxKinEnergy);
556  }
557  }
558 
559  theCuts = modelManager->Initialise(particle, secondaryParticle,
560  minSubRange, verboseLevel);
561 
562  // Sub Cutoff
563  if (nSCoffRegions>0) {
564  theSubCuts = modelManager->SubCutoff();
565 
566  if(nSCoffRegions>0) { idxSCoffRegions = new G4bool[n]; }
567  for (size_t j=0; j<n; ++j) {
568 
569  const G4MaterialCutsCouple* couple =
570  theCoupleTable->GetMaterialCutsCouple(j);
571  const G4ProductionCuts* pcuts = couple->GetProductionCuts();
572 
573  if(nSCoffRegions>0) {
574  G4bool reg = false;
575  for(G4int i=0; i<nSCoffRegions; ++i) {
576  if( pcuts == scoffRegions[i]->GetProductionCuts()) { reg = true; }
577  }
578  idxSCoffRegions[j] = reg;
579  }
580  }
581  }
582 
583  if(1 < verboseLevel) {
584  G4cout << "G4VEnergyLossProcess::PrepearPhysicsTable() is done "
585  << " for local " << particle->GetParticleName()
586  << " isIon= " << isIon;
587  if(baseParticle) { G4cout << "; base: " << baseParticle->GetParticleName(); }
588  G4cout << " chargeSqRatio= " << chargeSqRatio
589  << " massRatio= " << massRatio
590  << " reduceFactor= " << reduceFactor << G4endl;
591  if (nSCoffRegions) {
592  G4cout << " SubCutoff Regime is ON for regions: " << G4endl;
593  for (G4int i=0; i<nSCoffRegions; ++i) {
594  const G4Region* r = scoffRegions[i];
595  G4cout << " " << r->GetName() << G4endl;
596  }
597  }
598  }
599 }
600 
601 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
602 
604 {
605  G4bool verb = false;
606  if(1 < verboseLevel || verb) {
607 
608  //if(1 < verboseLevel) {
609  G4cout << "### G4VEnergyLossProcess::BuildPhysicsTable() for "
610  << GetProcessName()
611  << " and particle " << part.GetParticleName()
612  << "; local: " << particle->GetParticleName();
613  if(baseParticle) {
614  G4cout << "; base: " << baseParticle->GetParticleName();
615  }
616  G4cout << " TablesAreBuilt= " << tablesAreBuilt
617  << " isIon= " << isIon << " " << this << G4endl;
618  }
619 
620  G4bool master = true;
621  const G4VEnergyLossProcess* masterProcess =
622  static_cast<const G4VEnergyLossProcess*>(GetMasterProcess());
623  if(masterProcess != this) { master = false; }
624 
625  if(&part == particle) {
626 
627  G4LossTableBuilder* bld = lManager->GetTableBuilder();
628  if(master) {
629  theDensityFactor = bld->GetDensityFactors();
630  theDensityIdx = bld->GetCoupleIndexes();
631  lManager->BuildPhysicsTable(particle, this);
632 
633  } else {
634 
635  // define density factors for worker thread
636  bld->InitialiseBaseMaterials(masterProcess->DEDXTable());
637  theDensityFactor = bld->GetDensityFactors();
638  theDensityIdx = bld->GetCoupleIndexes();
639 
640  // copy table pointers from master thread
641  SetDEDXTable(masterProcess->DEDXTable(),fRestricted);
643  SetDEDXTable(masterProcess->DEDXunRestrictedTable(),fTotal);
644  SetDEDXTable(masterProcess->IonisationTable(),fIsIonisation);
646  SetRangeTableForLoss(masterProcess->RangeTableForLoss());
647  SetCSDARangeTable(masterProcess->CSDARangeTable());
649  SetInverseRangeTable(masterProcess->InverseRangeTable());
650  SetLambdaTable(masterProcess->LambdaTable());
651  SetSubLambdaTable(masterProcess->SubLambdaTable());
652  isIonisation = masterProcess->IsIonisationProcess();
653 
654  tablesAreBuilt = true;
655  // local initialisation of models
656  G4bool printing = true;
657  G4int numberOfModels = modelManager->NumberOfModels();
658  for(G4int i=0; i<numberOfModels; ++i) {
659  G4VEmModel* mod = GetModelByIndex(i, printing);
660  G4VEmModel* mod0= masterProcess->GetModelByIndex(i,printing);
661  mod->InitialiseLocal(particle, mod0);
662  }
663 
664  lManager->LocalPhysicsTables(particle, this);
665  }
666 
667  // needs to be done only once
668  safetyHelper->InitialiseHelper();
669  }
670  // explicitly defined printout by particle name
671  G4String num = part.GetParticleName();
672  if(1 < verboseLevel ||
673  (0 < verboseLevel && (num == "e-" ||
674  num == "e+" || num == "mu+" ||
675  num == "mu-" || num == "proton"||
676  num == "pi+" || num == "pi-" ||
677  num == "kaon+" || num == "kaon-" ||
678  num == "alpha" || num == "anti_proton" ||
679  num == "GenericIon")))
680  {
681  PrintInfoDefinition(part);
682  }
683 
684  // Added tracking cut to avoid tracking artifacts
685  // identify deexcitation flag
686  if(isIonisation) {
687  fParticleChange.SetLowEnergyLimit(lowestKinEnergy);
688  atomDeexcitation = lManager->AtomDeexcitation();
689  if(atomDeexcitation) {
690  if(atomDeexcitation->IsPIXEActive()) { useDeexcitation = true; }
691  }
692  }
693 
694  //if(1 < verboseLevel || verb) {
695  if(1 < verboseLevel) {
696  G4cout << "### G4VEnergyLossProcess::BuildPhysicsTable() done for "
697  << GetProcessName()
698  << " and particle " << part.GetParticleName();
699  if(isIonisation) { G4cout << " isIonisation flag = 1"; }
700  G4cout << G4endl;
701  }
702 }
703 
704 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
705 
707 {
708  G4bool verb = false;
709  if(1 < verboseLevel || verb) {
710  G4cout << "G4VEnergyLossProcess::BuildDEDXTable() of type " << tType
711  << " for " << GetProcessName()
712  << " and particle " << particle->GetParticleName()
713  << G4endl;
714  }
715  G4PhysicsTable* table = 0;
716  G4double emax = maxKinEnergy;
717  G4int bin = nBins;
718 
719  if(fTotal == tType) {
720  emax = maxKinEnergyCSDA;
721  bin = nBinsCSDA;
722  table = theDEDXunRestrictedTable;
723  } else if(fRestricted == tType) {
724  table = theDEDXTable;
725  } else if(fSubRestricted == tType) {
726  table = theDEDXSubTable;
727  } else {
728  G4cout << "G4VEnergyLossProcess::BuildDEDXTable WARNING: wrong type "
729  << tType << G4endl;
730  }
731 
732  // Access to materials
733  const G4ProductionCutsTable* theCoupleTable=
735  size_t numOfCouples = theCoupleTable->GetTableSize();
736 
737  if(1 < verboseLevel || verb) {
738  G4cout << numOfCouples << " materials"
739  << " minKinEnergy= " << minKinEnergy
740  << " maxKinEnergy= " << emax
741  << " nbin= " << bin
742  << " EmTableType= " << tType
743  << " table= " << table << " " << this
744  << G4endl;
745  }
746  if(!table) { return table; }
747 
748  G4LossTableBuilder* bld = lManager->GetTableBuilder();
749  G4bool splineFlag = lManager->SplineFlag();
750  G4PhysicsLogVector* aVector = 0;
751  G4PhysicsLogVector* bVector = 0;
752 
753  for(size_t i=0; i<numOfCouples; ++i) {
754 
755  if(1 < verboseLevel || verb) {
756  G4cout << "G4VEnergyLossProcess::BuildDEDXVector Idx= " << i
757  << " flagTable= " << table->GetFlag(i)
758  << " Flag= " << bld->GetFlag(i) << G4endl;
759  }
760  if(bld->GetFlag(i)) {
761 
762  // create physics vector and fill it
763  const G4MaterialCutsCouple* couple =
764  theCoupleTable->GetMaterialCutsCouple(i);
765  delete (*table)[i];
766  if(!bVector) {
767  aVector = new G4PhysicsLogVector(minKinEnergy, emax, bin);
768  bVector = aVector;
769  } else {
770  aVector = new G4PhysicsLogVector(*bVector);
771  }
772  aVector->SetSpline(splineFlag);
773 
774  modelManager->FillDEDXVector(aVector, couple, tType);
775  if(splineFlag) { aVector->FillSecondDerivatives(); }
776 
777  // Insert vector for this material into the table
778  G4PhysicsTableHelper::SetPhysicsVector(table, i, aVector);
779  }
780  }
781 
782  if(1 < verboseLevel || verb) {
783  G4cout << "G4VEnergyLossProcess::BuildDEDXTable(): table is built for "
784  << particle->GetParticleName()
785  << " and process " << GetProcessName()
786  << G4endl;
787  //if(2 < verboseLevel) G4cout << (*table) << G4endl;
788  }
789 
790  return table;
791 }
792 
793 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
794 
796 {
797  G4PhysicsTable* table = 0;
798 
799  if(fRestricted == tType) {
800  table = theLambdaTable;
801  } else if(fSubRestricted == tType) {
802  table = theSubLambdaTable;
803  } else {
804  G4cout << "G4VEnergyLossProcess::BuildLambdaTable WARNING: wrong type "
805  << tType << G4endl;
806  }
807 
808  if(1 < verboseLevel) {
809  G4cout << "G4VEnergyLossProcess::BuildLambdaTable() of type "
810  << tType << " for process "
811  << GetProcessName() << " and particle "
812  << particle->GetParticleName()
813  << " EmTableType= " << tType
814  << " table= " << table
815  << G4endl;
816  }
817  if(!table) {return table;}
818 
819  // Access to materials
820  const G4ProductionCutsTable* theCoupleTable=
822  size_t numOfCouples = theCoupleTable->GetTableSize();
823 
824  G4LossTableBuilder* bld = lManager->GetTableBuilder();
825  theDensityFactor = bld->GetDensityFactors();
826  theDensityIdx = bld->GetCoupleIndexes();
827 
828  G4bool splineFlag = lManager->SplineFlag();
829  G4PhysicsLogVector* aVector = 0;
830  G4double scale = G4Log(maxKinEnergy/minKinEnergy);
831 
832  for(size_t i=0; i<numOfCouples; ++i) {
833 
834  if (bld->GetFlag(i)) {
835 
836  // create physics vector and fill it
837  const G4MaterialCutsCouple* couple =
838  theCoupleTable->GetMaterialCutsCouple(i);
839  delete (*table)[i];
840 
841  G4bool startNull = true;
842  G4double emin =
843  MinPrimaryEnergy(particle,couple->GetMaterial(),(*theCuts)[i]);
844  if(minKinEnergy > emin) {
845  emin = minKinEnergy;
846  startNull = false;
847  }
848 
849  G4double emax = maxKinEnergy;
850  if(emax <= emin) { emax = 2*emin; }
851  G4int bin = G4lrint(nBins*G4Log(emax/emin)/scale);
852  if(bin < 3) { bin = 3; }
853  aVector = new G4PhysicsLogVector(emin, emax, bin);
854  aVector->SetSpline(splineFlag);
855 
856  modelManager->FillLambdaVector(aVector, couple, startNull, tType);
857  if(splineFlag) { aVector->FillSecondDerivatives(); }
858 
859  // Insert vector for this material into the table
860  G4PhysicsTableHelper::SetPhysicsVector(table, i, aVector);
861  }
862  }
863 
864  if(1 < verboseLevel) {
865  G4cout << "Lambda table is built for "
866  << particle->GetParticleName()
867  << G4endl;
868  }
869 
870  return table;
871 }
872 
873 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
874 
875 void
877 {
878  if(0 < verboseLevel) {
879  G4cout << std::setprecision(6);
880  G4cout << G4endl << GetProcessName() << ": for "
881  << part.GetParticleName()
882  << " SubType= " << GetProcessSubType()
883  << G4endl;
884  G4cout << " dE/dx and range tables from "
885  << G4BestUnit(minKinEnergy,"Energy")
886  << " to " << G4BestUnit(maxKinEnergy,"Energy")
887  << " in " << nBins << " bins" << G4endl
888  << " Lambda tables from threshold to "
889  << G4BestUnit(maxKinEnergy,"Energy")
890  << " in " << nBins << " bins, spline: "
891  << lManager->SplineFlag()
892  << G4endl;
893  if(theRangeTableForLoss && isIonisation) {
894  G4cout << " finalRange(mm)= " << finalRange/mm
895  << ", dRoverRange= " << dRoverRange
896  << ", integral: " << integral
897  << ", fluct: " << lossFluctuationFlag
898  << ", linLossLimit= " << linLossLimit
899  << G4endl;
900  }
901  PrintInfo();
902  modelManager->DumpModelList(verboseLevel);
903  if(theCSDARangeTable && isIonisation) {
904  G4cout << " CSDA range table up"
905  << " to " << G4BestUnit(maxKinEnergyCSDA,"Energy")
906  << " in " << nBinsCSDA << " bins" << G4endl;
907  }
908  if(nSCoffRegions>0 && isIonisation) {
909  G4cout << " Subcutoff sampling in " << nSCoffRegions
910  << " regions" << G4endl;
911  }
912  if(2 < verboseLevel) {
913  G4cout << " DEDXTable address= " << theDEDXTable << G4endl;
914  if(theDEDXTable && isIonisation) G4cout << (*theDEDXTable) << G4endl;
915  G4cout << "non restricted DEDXTable address= "
916  << theDEDXunRestrictedTable << G4endl;
917  if(theDEDXunRestrictedTable && isIonisation) {
918  G4cout << (*theDEDXunRestrictedTable) << G4endl;
919  }
920  if(theDEDXSubTable && isIonisation) {
921  G4cout << (*theDEDXSubTable) << G4endl;
922  }
923  G4cout << " CSDARangeTable address= " << theCSDARangeTable
924  << G4endl;
925  if(theCSDARangeTable && isIonisation) {
926  G4cout << (*theCSDARangeTable) << G4endl;
927  }
928  G4cout << " RangeTableForLoss address= " << theRangeTableForLoss
929  << G4endl;
930  if(theRangeTableForLoss && isIonisation) {
931  G4cout << (*theRangeTableForLoss) << G4endl;
932  }
933  G4cout << " InverseRangeTable address= " << theInverseRangeTable
934  << G4endl;
935  if(theInverseRangeTable && isIonisation) {
936  G4cout << (*theInverseRangeTable) << G4endl;
937  }
938  G4cout << " LambdaTable address= " << theLambdaTable << G4endl;
939  if(theLambdaTable && isIonisation) {
940  G4cout << (*theLambdaTable) << G4endl;
941  }
942  G4cout << " SubLambdaTable address= " << theSubLambdaTable
943  << G4endl;
944  if(theSubLambdaTable && isIonisation) {
945  G4cout << (*theSubLambdaTable) << G4endl;
946  }
947  }
948  }
949 }
950 
951 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
952 
954 {
955  G4RegionStore* regionStore = G4RegionStore::GetInstance();
956  const G4Region* reg = r;
957  if (!reg) {
958  reg = regionStore->GetRegion("DefaultRegionForTheWorld", false);
959  }
960 
961  // the region is in the list
962  if (nSCoffRegions) {
963  for (G4int i=0; i<nSCoffRegions; ++i) {
964  if (reg == scoffRegions[i]) {
965  return;
966  }
967  }
968  }
969 
970  // new region
971  if(val) {
972  useSubCutoff = true;
973  scoffRegions.push_back(reg);
974  ++nSCoffRegions;
975  } else {
976  useSubCutoff = false;
977  }
978 }
979 
980 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
981 
983 {
984  /*
985  G4cout << track->GetDefinition()->GetParticleName()
986  << " e(MeV)= " << track->GetKineticEnergy()
987  << " baseParticle " << baseParticle << " proc " << this;
988  if(particle) G4cout << " " << particle->GetParticleName();
989  G4cout << " isIon= " << isIon << " dedx " << theDEDXTable <<G4endl;
990  */
991  // reset parameters for the new track
993  mfpKinEnergy = DBL_MAX;
994  preStepRangeEnergy = 0.0;
995 
996  // reset ion
997  if(isIon) {
998  chargeSqRatio = 0.5;
999 
1000  G4double newmass = track->GetDefinition()->GetPDGMass();
1001  if(baseParticle) {
1002  massRatio = baseParticle->GetPDGMass()/newmass;
1003  } else if(theGenericIon) {
1004  massRatio = proton_mass_c2/newmass;
1005  } else {
1006  massRatio = 1.0;
1007  }
1008  }
1009  // forced biasing only for primary particles
1010  if(biasManager) {
1011  if(0 == track->GetParentID()) {
1012  // primary particle
1013  biasFlag = true;
1014  biasManager->ResetForcedInteraction();
1015  }
1016  }
1017 }
1018 
1019 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1020 
1023  G4GPILSelection* selection)
1024 {
1025  G4double x = DBL_MAX;
1026  *selection = aGPILSelection;
1027  if(isIonisation) {
1028  fRange = GetScaledRangeForScaledEnergy(preStepScaledEnergy)*reduceFactor;
1029  x = fRange;
1030  G4double finR = finalRange;
1031  if(rndmStepFlag) {
1032  finR = std::min(finR,
1033  currentCouple->GetProductionCuts()->GetProductionCut(1));
1034  }
1035  if(fRange > finR) {
1036  x = fRange*dRoverRange + finR*(1.0 - dRoverRange)*(2.0 - finR/fRange);
1037  }
1038  /*
1039  if(particle->GetPDGMass() > 0.9*GeV)
1040  G4cout<<GetProcessName()<<": e= "<<preStepKinEnergy
1041  <<" range= "<<fRange << " idx= " << basedCoupleIndex
1042  << " finR= " << finR
1043  << " limit= " << x <<G4endl;
1044  G4cout << "massRatio= " << massRatio << " Q^2= " << chargeSqRatio
1045  << " finR= " << finR << " dRoverRange= " << dRoverRange
1046  << " finalRange= " << finalRange << G4endl;
1047  */
1048  }
1049  //G4cout<<GetProcessName()<<": e= "<<preStepKinEnergy
1050  //<<" stepLimit= "<<x<<G4endl;
1051  return x;
1052 }
1053 
1054 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1055 
1057  const G4Track& track,
1058  G4double previousStepSize,
1060 {
1061  // condition is set to "Not Forced"
1062  *condition = NotForced;
1063  G4double x = DBL_MAX;
1064 
1065  // initialisation of material, mass, charge, model
1066  // at the beginning of the step
1067  DefineMaterial(track.GetMaterialCutsCouple());
1068  preStepKinEnergy = track.GetKineticEnergy();
1069  preStepScaledEnergy = preStepKinEnergy*massRatio;
1070  SelectModel(preStepScaledEnergy);
1071 
1072  if(!currentModel->IsActive(preStepScaledEnergy)) {
1074  return x;
1075  }
1076 
1077  // change effective charge of an ion on fly
1078  if(isIon) {
1079  G4double q2 = currentModel->ChargeSquareRatio(track);
1080  if(q2 != chargeSqRatio && q2 > 0.0) {
1081  chargeSqRatio = q2;
1082  fFactor = q2*biasFactor*(*theDensityFactor)[currentCoupleIndex];
1083  reduceFactor = 1.0/(fFactor*massRatio);
1084  }
1085  }
1086  // if(particle->GetPDGMass() > 0.9*GeV)
1087  //G4cout << "q2= "<<chargeSqRatio << " massRatio= " << massRatio << G4endl;
1088  // initialisation for sampling of the interaction length
1089  //if(previousStepSize <= 0.0) { theNumberOfInteractionLengthLeft = -1.0; }
1090  //if(theNumberOfInteractionLengthLeft < 0.0) { mfpKinEnergy = DBL_MAX; }
1091 
1092  // forced biasing only for primary particles
1093  if(biasManager) {
1094  if(0 == track.GetParentID()) {
1095  if(biasFlag &&
1096  biasManager->ForcedInteractionRegion(currentCoupleIndex)) {
1097  return biasManager->GetStepLimit(currentCoupleIndex, previousStepSize);
1098  }
1099  }
1100  }
1101 
1102  // compute mean free path
1103  if(preStepScaledEnergy < mfpKinEnergy) {
1104  if (integral) { ComputeLambdaForScaledEnergy(preStepScaledEnergy); }
1105  else { preStepLambda = GetLambdaForScaledEnergy(preStepScaledEnergy); }
1106 
1107  // zero cross section
1108  if(preStepLambda <= 0.0) {
1111  }
1112  }
1113 
1114  // non-zero cross section
1115  if(preStepLambda > 0.0) {
1117 
1118  // beggining of tracking (or just after DoIt of this process)
1119  // ResetNumberOfInteractionLengthLeft();
1122 
1123  } else if(currentInteractionLength < DBL_MAX) {
1124 
1125  // subtract NumberOfInteractionLengthLeft using previous step
1127  previousStepSize/currentInteractionLength;
1128  // SubtractNumberOfInteractionLengthLeft(previousStepSize);
1131  //theNumberOfInteractionLengthLeft = perMillion;
1132  }
1133  }
1134 
1135  // new mean free path and step limit
1136  currentInteractionLength = 1.0/preStepLambda;
1138  }
1139 #ifdef G4VERBOSE
1140  if (verboseLevel>2){
1141  // if(particle->GetPDGMass() > 0.9*GeV){
1142  G4cout << "G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength ";
1143  G4cout << "[ " << GetProcessName() << "]" << G4endl;
1144  G4cout << " for " << track.GetDefinition()->GetParticleName()
1145  << " in Material " << currentMaterial->GetName()
1146  << " Ekin(MeV)= " << preStepKinEnergy/MeV
1147  << " " << track.GetMaterial()->GetName()
1148  <<G4endl;
1149  G4cout << "MeanFreePath = " << currentInteractionLength/cm << "[cm]"
1150  << "InteractionLength= " << x/cm <<"[cm] " <<G4endl;
1151  }
1152 #endif
1153  return x;
1154 }
1155 
1156 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1157 
1159  const G4Step& step)
1160 {
1162  // The process has range table - calculate energy loss
1163  if(!isIonisation || !currentModel->IsActive(preStepScaledEnergy)) {
1164  return &fParticleChange;
1165  }
1166 
1167  // Get the actual (true) Step length
1168  G4double length = step.GetStepLength();
1169  if(length <= 0.0) { return &fParticleChange; }
1170  G4double eloss = 0.0;
1171 
1172  /*
1173  if(-1 < verboseLevel) {
1174  const G4ParticleDefinition* d = track.GetParticleDefinition();
1175  G4cout << "AlongStepDoIt for "
1176  << GetProcessName() << " and particle "
1177  << d->GetParticleName()
1178  << " eScaled(MeV)= " << preStepScaledEnergy/MeV
1179  << " range(mm)= " << fRange/mm
1180  << " s(mm)= " << length/mm
1181  << " rf= " << reduceFactor
1182  << " q^2= " << chargeSqRatio
1183  << " md= " << d->GetPDGMass()
1184  << " status= " << track.GetTrackStatus()
1185  << " " << track.GetMaterial()->GetName()
1186  << G4endl;
1187  }
1188  */
1189 
1190  const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
1191 
1192  // define new weight for primary and secondaries
1194  if(weightFlag) {
1195  weight /= biasFactor;
1197  }
1198 
1199  // stopping
1200  if (length >= fRange) {
1201  eloss = preStepKinEnergy;
1202  if (useDeexcitation) {
1203  atomDeexcitation->AlongStepDeexcitation(scTracks, step,
1204  eloss, currentCoupleIndex);
1205  if(scTracks.size() > 0) { FillSecondariesAlongStep(eloss, weight); }
1206  if(eloss < 0.0) { eloss = 0.0; }
1207  }
1210  return &fParticleChange;
1211  }
1212  //G4cout << theDEDXTable << " idx= " << basedCoupleIndex
1213  // << " " << GetProcessName() << " "<< currentMaterial->GetName()<<G4endl;
1214  //if(particle->GetParticleName() == "e-")G4cout << (*theDEDXTable) <<G4endl;
1215  // Short step
1216  eloss = GetDEDXForScaledEnergy(preStepScaledEnergy)*length;
1217 
1218  //G4cout << "eloss= " << eloss << G4endl;
1219 
1220  // Long step
1221  if(eloss > preStepKinEnergy*linLossLimit) {
1222 
1223  G4double x = (fRange - length)/reduceFactor;
1224  //G4cout << "x= " << x << " " << theInverseRangeTable << G4endl;
1225  eloss = preStepKinEnergy - ScaledKinEnergyForLoss(x)/massRatio;
1226 
1227  /*
1228  if(-1 < verboseLevel)
1229  G4cout << "Long STEP: rPre(mm)= "
1230  << GetScaledRangeForScaledEnergy(preStepScaledEnergy)/mm
1231  << " rPost(mm)= " << x/mm
1232  << " ePre(MeV)= " << preStepScaledEnergy/MeV
1233  << " eloss(MeV)= " << eloss/MeV
1234  << " eloss0(MeV)= "
1235  << GetDEDXForScaledEnergy(preStepScaledEnergy)*length/MeV
1236  << " lim(MeV)= " << preStepKinEnergy*linLossLimit/MeV
1237  << G4endl;
1238  */
1239  }
1240 
1241  /*
1242  G4double eloss0 = eloss;
1243  if(-1 < verboseLevel ) {
1244  G4cout << "Before fluct: eloss(MeV)= " << eloss/MeV
1245  << " e-eloss= " << preStepKinEnergy-eloss
1246  << " step(mm)= " << length/mm
1247  << " range(mm)= " << fRange/mm
1248  << " fluct= " << lossFluctuationFlag
1249  << G4endl;
1250  }
1251  */
1252 
1253  G4double cut = (*theCuts)[currentCoupleIndex];
1254  G4double esec = 0.0;
1255 
1256  // SubCutOff
1257  if(useSubCutoff) {
1258  if(idxSCoffRegions[currentCoupleIndex]) {
1259 
1260  G4bool yes = false;
1261  G4StepPoint* prePoint = step.GetPreStepPoint();
1262 
1263  // Check boundary
1264  if(prePoint->GetStepStatus() == fGeomBoundary) { yes = true; }
1265 
1266  // Check PrePoint
1267  else {
1268  G4double preSafety = prePoint->GetSafety();
1269  G4double rcut =
1270  currentCouple->GetProductionCuts()->GetProductionCut(1);
1271 
1272  // recompute presafety
1273  if(preSafety < rcut) {
1274  preSafety = safetyHelper->ComputeSafety(prePoint->GetPosition());
1275  }
1276 
1277  if(preSafety < rcut) { yes = true; }
1278 
1279  // Check PostPoint
1280  else {
1281  G4double postSafety = preSafety - length;
1282  if(postSafety < rcut) {
1283  postSafety = safetyHelper->ComputeSafety(
1284  step.GetPostStepPoint()->GetPosition());
1285  if(postSafety < rcut) { yes = true; }
1286  }
1287  }
1288  }
1289 
1290  // Decided to start subcut sampling
1291  if(yes) {
1292 
1293  cut = (*theSubCuts)[currentCoupleIndex];
1294  eloss -= GetSubDEDXForScaledEnergy(preStepScaledEnergy)*length;
1295  esec = SampleSubCutSecondaries(scTracks, step,
1296  currentModel,currentCoupleIndex);
1297  // add bremsstrahlung sampling
1298  /*
1299  if(nProcesses > 0) {
1300  for(G4int i=0; i<nProcesses; ++i) {
1301  (scProcesses[i])->SampleSubCutSecondaries(
1302  scTracks, step, (scProcesses[i])->
1303  SelectModelForMaterial(preStepKinEnergy, currentCoupleIndex),
1304  currentCoupleIndex);
1305  }
1306  }
1307  */
1308  }
1309  }
1310  }
1311 
1312  // Corrections, which cannot be tabulated
1313  if(isIon) {
1314  G4double eadd = 0.0;
1315  G4double eloss_before = eloss;
1316  currentModel->CorrectionsAlongStep(currentCouple, dynParticle,
1317  eloss, eadd, length);
1318  if(eloss < 0.0) { eloss = 0.5*eloss_before; }
1319  }
1320 
1321  // Sample fluctuations
1322  if (lossFluctuationFlag) {
1323  G4VEmFluctuationModel* fluc = currentModel->GetModelOfFluctuations();
1324  if(fluc &&
1325  (eloss + esec + lowestKinEnergy) < preStepKinEnergy) {
1326 
1327  G4double tmax =
1328  std::min(currentModel->MaxSecondaryKinEnergy(dynParticle),cut);
1329  eloss = fluc->SampleFluctuations(currentCouple,dynParticle,
1330  tmax,length,eloss);
1331  /*
1332  if(-1 < verboseLevel)
1333  G4cout << "After fluct: eloss(MeV)= " << eloss/MeV
1334  << " fluc= " << (eloss-eloss0)/MeV
1335  << " ChargeSqRatio= " << chargeSqRatio
1336  << " massRatio= " << massRatio
1337  << " tmax= " << tmax
1338  << G4endl;
1339  */
1340  }
1341  }
1342 
1343  // deexcitation
1344  if (useDeexcitation) {
1345  G4double esecfluo = preStepKinEnergy - esec;
1346  G4double de = esecfluo;
1347  //G4double eloss0 = eloss;
1348  /*
1349  G4cout << "### 1: E(keV)= " << preStepKinEnergy/keV
1350  << " Efluomax(keV)= " << de/keV
1351  << " Eloss(keV)= " << eloss/keV << G4endl;
1352  */
1353  atomDeexcitation->AlongStepDeexcitation(scTracks, step,
1354  de, currentCoupleIndex);
1355 
1356  // sum of de-excitation energies
1357  esecfluo -= de;
1358 
1359  // subtracted from energy loss
1360  if(eloss >= esecfluo) {
1361  esec += esecfluo;
1362  eloss -= esecfluo;
1363  } else {
1364  esec += esecfluo;
1365  eloss = 0.0;
1366  }
1367  /*
1368  if(esecfluo > 0.0) {
1369  G4cout << "### 2: E(keV)= " << preStepKinEnergy/keV
1370  << " Esec(keV)= " << esec/keV
1371  << " Esecf(kV)= " << esecfluo/keV
1372  << " Eloss0(kV)= " << eloss0/keV
1373  << " Eloss(keV)= " << eloss/keV
1374  << G4endl;
1375  }
1376  */
1377  }
1378  if(scTracks.size() > 0) { FillSecondariesAlongStep(eloss, weight); }
1379 
1380  // Energy balanse
1381  G4double finalT = preStepKinEnergy - eloss - esec;
1382  if (finalT <= lowestKinEnergy) {
1383  eloss += finalT;
1384  finalT = 0.0;
1385  } else if(isIon) {
1387  currentModel->GetParticleCharge(track.GetParticleDefinition(),
1388  currentMaterial,finalT));
1389  }
1390 
1391  if(eloss < 0.0) { eloss = 0.0; }
1394  /*
1395  if(-1 < verboseLevel) {
1396  G4double del = finalT + eloss + esec - preStepKinEnergy;
1397  G4cout << "Final value eloss(MeV)= " << eloss/MeV
1398  << " preStepKinEnergy= " << preStepKinEnergy
1399  << " postStepKinEnergy= " << finalT
1400  << " de(keV)= " << del/keV
1401  << " lossFlag= " << lossFluctuationFlag
1402  << " status= " << track.GetTrackStatus()
1403  << G4endl;
1404  }
1405  */
1406  return &fParticleChange;
1407 }
1408 
1409 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1410 
1411 void
1412 G4VEnergyLossProcess::FillSecondariesAlongStep(G4double&, G4double& weight)
1413 {
1414  G4int n0 = scTracks.size();
1415 
1416  // weight may be changed by biasing manager
1417  if(biasManager) {
1418  if(biasManager->SecondaryBiasingRegion(currentCoupleIndex)) {
1419  weight *=
1420  biasManager->ApplySecondaryBiasing(scTracks, currentCoupleIndex);
1421  }
1422  }
1423 
1424  // fill secondaries
1425  G4int n = scTracks.size();
1427 
1428  for(G4int i=0; i<n; ++i) {
1429  G4Track* t = scTracks[i];
1430  if(t) {
1431  t->SetWeight(weight);
1433  if(i < n0) { t->SetCreatorModelIndex(secID); }
1434  else { t->SetCreatorModelIndex(biasID); }
1435  //G4cout << "Secondary(along step) has weight " << t->GetWeight()
1436  //<< ", kenergy " << t->GetKineticEnergy()/MeV << " MeV" <<G4endl;
1437  }
1438  }
1439  scTracks.clear();
1440 }
1441 
1442 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1443 
1444 G4double
1445 G4VEnergyLossProcess::SampleSubCutSecondaries(std::vector<G4Track*>& tracks,
1446  const G4Step& step,
1447  G4VEmModel* model,
1448  G4int idx)
1449 {
1450  // Fast check weather subcutoff can work
1451  G4double esec = 0.0;
1452  G4double subcut = (*theSubCuts)[idx];
1453  G4double cut = (*theCuts)[idx];
1454  if(cut <= subcut) { return esec; }
1455 
1456  const G4Track* track = step.GetTrack();
1457  const G4DynamicParticle* dp = track->GetDynamicParticle();
1458  G4double e = dp->GetKineticEnergy()*massRatio;
1459  G4double cross = (*theDensityFactor)[idx]*chargeSqRatio
1460  *(((*theSubLambdaTable)[(*theDensityIdx)[idx]])->Value(e, idxSubLambda));
1461  G4double length = step.GetStepLength();
1462 
1463  // negligible probability to get any interaction
1464  if(length*cross < perMillion) { return esec; }
1465  /*
1466  if(-1 < verboseLevel)
1467  G4cout << "<<< Subcutoff for " << GetProcessName()
1468  << " cross(1/mm)= " << cross*mm << ">>>"
1469  << " e(MeV)= " << preStepScaledEnergy
1470  << " matIdx= " << currentCoupleIndex
1471  << G4endl;
1472  */
1473 
1474  // Sample subcutoff secondaries
1475  G4StepPoint* preStepPoint = step.GetPreStepPoint();
1476  G4StepPoint* postStepPoint = step.GetPostStepPoint();
1477  G4ThreeVector prepoint = preStepPoint->GetPosition();
1478  G4ThreeVector dr = postStepPoint->GetPosition() - prepoint;
1479  G4double pretime = preStepPoint->GetGlobalTime();
1480  G4double dt = postStepPoint->GetGlobalTime() - pretime;
1481  //G4double dt = length/preStepPoint->GetVelocity();
1482  G4double fragment = 0.0;
1483 
1484  do {
1485  G4double del = -std::log(G4UniformRand())/cross;
1486  fragment += del/length;
1487  if (fragment > 1.0) break;
1488 
1489  // sample secondaries
1490  secParticles.clear();
1491  model->SampleSecondaries(&secParticles,track->GetMaterialCutsCouple(),
1492  dp,subcut,cut);
1493 
1494  // position of subcutoff particles
1495  G4ThreeVector r = prepoint + fragment*dr;
1496  std::vector<G4DynamicParticle*>::iterator it;
1497  for(it=secParticles.begin(); it!=secParticles.end(); ++it) {
1498 
1499  G4Track* t = new G4Track((*it), pretime + fragment*dt, r);
1501  t->SetCreatorModelIndex(subsecID);
1502  tracks.push_back(t);
1503  esec += t->GetKineticEnergy();
1504  if (t->GetParticleDefinition() == thePositron) {
1505  esec += 2.0*electron_mass_c2;
1506  }
1507 
1508  /*
1509  if(-1 < verboseLevel)
1510  G4cout << "New track "
1511  << t->GetParticleDefinition()->GetParticleName()
1512  << " e(keV)= " << t->GetKineticEnergy()/keV
1513  << " fragment= " << fragment
1514  << G4endl;
1515  */
1516  }
1517  } while (fragment <= 1.0);
1518  return esec;
1519 }
1520 
1521 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1522 
1524  const G4Step& step)
1525 {
1526  // In all cases clear number of interaction lengths
1528  mfpKinEnergy = DBL_MAX;
1529 
1531  G4double finalT = track.GetKineticEnergy();
1532  if(finalT <= lowestKinEnergy) { return &fParticleChange; }
1533 
1534  G4double postStepScaledEnergy = finalT*massRatio;
1535 
1536  if(!currentModel->IsActive(postStepScaledEnergy)) {
1537  return &fParticleChange;
1538  }
1539  /*
1540  if(-1 < verboseLevel) {
1541  G4cout << GetProcessName()
1542  << "::PostStepDoIt: E(MeV)= " << finalT/MeV
1543  << G4endl;
1544  }
1545  */
1546 
1547  // forced process - should happen only once per track
1548  if(biasFlag) {
1549  if(biasManager->ForcedInteractionRegion(currentCoupleIndex)) {
1550  biasFlag = false;
1551  }
1552  }
1553 
1554  // Integral approach
1555  if (integral) {
1556  G4double lx = GetLambdaForScaledEnergy(postStepScaledEnergy);
1557  /*
1558  if(preStepLambda<lx && 1 < verboseLevel && nWarnings<200) {
1559  G4cout << "WARNING: for " << particle->GetParticleName()
1560  << " and " << GetProcessName()
1561  << " E(MeV)= " << finalT/MeV
1562  << " preLambda= " << preStepLambda
1563  << " < " << lx << " (postLambda) "
1564  << G4endl;
1565  ++nWarnings;
1566  }
1567  */
1568  if(lx <= 0.0) {
1569  return &fParticleChange;
1570  } else if(preStepLambda*G4UniformRand() > lx) {
1571  return &fParticleChange;
1572  }
1573  }
1574 
1575  SelectModel(postStepScaledEnergy);
1576 
1577  // define new weight for primary and secondaries
1579  if(weightFlag) {
1580  weight /= biasFactor;
1582  }
1583 
1584  const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
1585  G4double tcut = (*theCuts)[currentCoupleIndex];
1586 
1587  // sample secondaries
1588  secParticles.clear();
1589  //G4cout<< "Eprimary: "<<dynParticle->GetKineticEnergy()/MeV<<G4endl;
1590  currentModel->SampleSecondaries(&secParticles, currentCouple,
1591  dynParticle, tcut);
1592 
1593  G4int num0 = secParticles.size();
1594 
1595  // bremsstrahlung splitting or Russian roulette
1596  if(biasManager) {
1597  if(biasManager->SecondaryBiasingRegion(currentCoupleIndex)) {
1598  G4double eloss = 0.0;
1599  weight *= biasManager->ApplySecondaryBiasing(
1600  secParticles,
1601  track, currentModel,
1602  &fParticleChange, eloss,
1603  currentCoupleIndex, tcut,
1604  step.GetPostStepPoint()->GetSafety());
1605  if(eloss > 0.0) {
1608  }
1609  }
1610  }
1611 
1612  // save secondaries
1613  G4int num = secParticles.size();
1614  if(num > 0) {
1615 
1617  G4double time = track.GetGlobalTime();
1618 
1619  for (G4int i=0; i<num; ++i) {
1620  if(secParticles[i]) {
1621  G4Track* t = new G4Track(secParticles[i], time, track.GetPosition());
1623  t->SetWeight(weight);
1624  if(i < num0) { t->SetCreatorModelIndex(secID); }
1625  else { t->SetCreatorModelIndex(biasID); }
1626 
1627  //G4cout << "Secondary(post step) has weight " << t->GetWeight()
1628  // << ", kenergy " << t->GetKineticEnergy()/MeV << " MeV"
1629  // << " time= " << time/ns << " ns " << G4endl;
1631  }
1632  }
1633  }
1634 
1637  if(particle->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
1640  }
1641 
1642  /*
1643  if(-1 < verboseLevel) {
1644  G4cout << "::PostStepDoIt: Sample secondary; Efin= "
1645  << fParticleChange.GetProposedKineticEnergy()/MeV
1646  << " MeV; model= (" << currentModel->LowEnergyLimit()
1647  << ", " << currentModel->HighEnergyLimit() << ")"
1648  << " preStepLambda= " << preStepLambda
1649  << " dir= " << track.GetMomentumDirection()
1650  << " status= " << track.GetTrackStatus()
1651  << G4endl;
1652  }
1653  */
1654  return &fParticleChange;
1655 }
1656 
1657 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1658 
1660  const G4ParticleDefinition* part, const G4String& directory,
1661  G4bool ascii)
1662 {
1663  G4bool res = true;
1664  if ( baseParticle || part != particle ) return res;
1665 
1666  if(!StoreTable(part,theDEDXTable,ascii,directory,"DEDX"))
1667  {res = false;}
1668 
1669  if(!StoreTable(part,theDEDXunRestrictedTable,ascii,directory,"DEDXnr"))
1670  {res = false;}
1671 
1672  if(!StoreTable(part,theDEDXSubTable,ascii,directory,"SubDEDX"))
1673  {res = false;}
1674 
1675  if(!StoreTable(part,theIonisationTable,ascii,directory,"Ionisation"))
1676  {res = false;}
1677 
1678  if(!StoreTable(part,theIonisationSubTable,ascii,directory,"SubIonisation"))
1679  {res = false;}
1680 
1681  if(isIonisation &&
1682  !StoreTable(part,theCSDARangeTable,ascii,directory,"CSDARange"))
1683  {res = false;}
1684 
1685  if(isIonisation &&
1686  !StoreTable(part,theRangeTableForLoss,ascii,directory,"Range"))
1687  {res = false;}
1688 
1689  if(isIonisation &&
1690  !StoreTable(part,theInverseRangeTable,ascii,directory,"InverseRange"))
1691  {res = false;}
1692 
1693  if(!StoreTable(part,theLambdaTable,ascii,directory,"Lambda"))
1694  {res = false;}
1695 
1696  if(!StoreTable(part,theSubLambdaTable,ascii,directory,"SubLambda"))
1697  {res = false;}
1698 
1699  if ( res ) {
1700  if(0 < verboseLevel) {
1701  G4cout << "Physics tables are stored for "
1702  << particle->GetParticleName()
1703  << " and process " << GetProcessName()
1704  << " in the directory <" << directory
1705  << "> " << G4endl;
1706  }
1707  } else {
1708  G4cout << "Fail to store Physics Tables for "
1709  << particle->GetParticleName()
1710  << " and process " << GetProcessName()
1711  << " in the directory <" << directory
1712  << "> " << G4endl;
1713  }
1714  return res;
1715 }
1716 
1717 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1718 
1719 G4bool
1721  const G4String& directory,
1722  G4bool ascii)
1723 {
1724  G4bool res = true;
1725  const G4String particleName = part->GetParticleName();
1726 
1727  if(1 < verboseLevel) {
1728  G4cout << "G4VEnergyLossProcess::RetrievePhysicsTable() for "
1729  << particleName << " and process " << GetProcessName()
1730  << "; tables_are_built= " << tablesAreBuilt
1731  << G4endl;
1732  }
1733  if(particle == part) {
1734 
1735  if ( !baseParticle ) {
1736 
1737  G4bool fpi = true;
1738  if(!RetrieveTable(part,theDEDXTable,ascii,directory,"DEDX",fpi))
1739  {fpi = false;}
1740 
1741  // ionisation table keeps individual dEdx and not sum of sub-processes
1742  if(!RetrieveTable(part,theDEDXTable,ascii,directory,"Ionisation",false))
1743  {fpi = false;}
1744 
1745  if(!RetrieveTable(part,theRangeTableForLoss,ascii,directory,"Range",fpi))
1746  {res = false;}
1747 
1748  if(!RetrieveTable(part,theDEDXunRestrictedTable,ascii,directory,
1749  "DEDXnr",false))
1750  {res = false;}
1751 
1752  if(!RetrieveTable(part,theCSDARangeTable,ascii,directory,
1753  "CSDARange",false))
1754  {res = false;}
1755 
1756  if(!RetrieveTable(part,theInverseRangeTable,ascii,directory,
1757  "InverseRange",fpi))
1758  {res = false;}
1759 
1760  if(!RetrieveTable(part,theLambdaTable,ascii,directory,"Lambda",true))
1761  {res = false;}
1762 
1763  G4bool yes = false;
1764  if(nSCoffRegions > 0) {yes = true;}
1765 
1766  if(!RetrieveTable(part,theDEDXSubTable,ascii,directory,"SubDEDX",yes))
1767  {res = false;}
1768 
1769  if(!RetrieveTable(part,theSubLambdaTable,ascii,directory,
1770  "SubLambda",yes))
1771  {res = false;}
1772 
1773  if(!fpi) yes = false;
1774  if(!RetrieveTable(part,theIonisationSubTable,ascii,directory,
1775  "SubIonisation",yes))
1776  {res = false;}
1777  }
1778  }
1779 
1780  return res;
1781 }
1782 
1783 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1784 
1785 G4bool G4VEnergyLossProcess::StoreTable(const G4ParticleDefinition* part,
1786  G4PhysicsTable* aTable, G4bool ascii,
1787  const G4String& directory,
1788  const G4String& tname)
1789 {
1790  G4bool res = true;
1791  if ( aTable ) {
1792  const G4String name = GetPhysicsTableFileName(part,directory,tname,ascii);
1793  if( !aTable->StorePhysicsTable(name,ascii)) res = false;
1794  }
1795  return res;
1796 }
1797 
1798 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1799 
1800 G4bool
1801 G4VEnergyLossProcess::RetrieveTable(const G4ParticleDefinition* part,
1802  G4PhysicsTable* aTable,
1803  G4bool ascii,
1804  const G4String& directory,
1805  const G4String& tname,
1806  G4bool mandatory)
1807 {
1808  G4bool isRetrieved = false;
1809  G4String filename = GetPhysicsTableFileName(part,directory,tname,ascii);
1810  if(aTable) {
1811  if(aTable->ExistPhysicsTable(filename)) {
1812  if(G4PhysicsTableHelper::RetrievePhysicsTable(aTable,filename,ascii)) {
1813  isRetrieved = true;
1814  if(lManager->SplineFlag()) {
1815  size_t n = aTable->length();
1816  for(size_t i=0; i<n; ++i) {
1817  if((*aTable)[i]) { (*aTable)[i]->SetSpline(true); }
1818  }
1819  }
1820  if (0 < verboseLevel) {
1821  G4cout << tname << " table for " << part->GetParticleName()
1822  << " is Retrieved from <" << filename << ">"
1823  << G4endl;
1824  }
1825  }
1826  }
1827  }
1828  if(mandatory && !isRetrieved) {
1829  if(0 < verboseLevel) {
1830  G4cout << tname << " table for " << part->GetParticleName()
1831  << " from file <"
1832  << filename << "> is not Retrieved"
1833  << G4endl;
1834  }
1835  return false;
1836  }
1837  return true;
1838 }
1839 
1840 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1841 
1843  const G4MaterialCutsCouple *couple,
1844  const G4DynamicParticle* dp,
1845  G4double length)
1846 {
1847  DefineMaterial(couple);
1848  G4double ekin = dp->GetKineticEnergy();
1849  SelectModel(ekin*massRatio);
1850  G4double tmax = currentModel->MaxSecondaryKinEnergy(dp);
1851  tmax = std::min(tmax,(*theCuts)[currentCoupleIndex]);
1852  G4double d = 0.0;
1853  G4VEmFluctuationModel* fm = currentModel->GetModelOfFluctuations();
1854  if(fm) { d = fm->Dispersion(currentMaterial,dp,tmax,length); }
1855  return d;
1856 }
1857 
1858 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1859 
1861  G4double kineticEnergy, const G4MaterialCutsCouple* couple)
1862 {
1863  // Cross section per volume is calculated
1864  DefineMaterial(couple);
1865  G4double cross = 0.0;
1866  if(theLambdaTable) {
1867  cross = GetLambdaForScaledEnergy(kineticEnergy*massRatio);
1868  } else {
1869  SelectModel(kineticEnergy*massRatio);
1870  cross = biasFactor*(*theDensityFactor)[currentCoupleIndex]
1871  *(currentModel->CrossSectionPerVolume(currentMaterial,
1872  particle, kineticEnergy,
1873  (*theCuts)[currentCoupleIndex]));
1874  }
1875  if(cross < 0.0) { cross = 0.0; }
1876  return cross;
1877 }
1878 
1879 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1880 
1882 {
1883  DefineMaterial(track.GetMaterialCutsCouple());
1884  preStepLambda = GetLambdaForScaledEnergy(track.GetKineticEnergy()*massRatio);
1885  G4double x = DBL_MAX;
1886  if(0.0 < preStepLambda) { x = 1.0/preStepLambda; }
1887  return x;
1888 }
1889 
1890 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1891 
1893  G4double x, G4double y,
1894  G4double& z)
1895 {
1896  G4GPILSelection sel;
1897  return AlongStepGetPhysicalInteractionLength(track, x, y, z, &sel);
1898 }
1899 
1900 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1901 
1903  const G4Track& track,
1904  G4double,
1906 
1907 {
1908  *condition = NotForced;
1909  return MeanFreePath(track);
1910 }
1911 
1912 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1913 
1915  const G4Track&,
1917 {
1918  return DBL_MAX;
1919 }
1920 
1921 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1922 
1925  G4double)
1926 {
1927  G4PhysicsVector* v =
1928  new G4PhysicsLogVector(minKinEnergy, maxKinEnergy, nBins);
1929  v->SetSpline(lManager->SplineFlag());
1930  return v;
1931 }
1932 
1933 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1934 
1937 {
1938  G4bool add = true;
1939  if(p->GetProcessName() != "eBrem") { add = false; }
1940  if(add && nProcesses > 0) {
1941  for(G4int i=0; i<nProcesses; ++i) {
1942  if(p == scProcesses[i]) {
1943  add = false;
1944  break;
1945  }
1946  }
1947  }
1948  if(add) {
1949  scProcesses.push_back(p);
1950  ++nProcesses;
1951  if (1 < verboseLevel) {
1952  G4cout << "### The process " << p->GetProcessName()
1953  << " is added to the list of collaborative processes of "
1954  << GetProcessName() << G4endl;
1955  }
1956  }
1957 }
1958 
1959 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1960 
1961 void
1963 {
1964  G4bool verb = false;
1965  if(fTotal == tType) {
1966  theDEDXunRestrictedTable = p;
1967  if(p) {
1968  size_t n = p->length();
1969  G4PhysicsVector* pv = (*p)[0];
1970  G4double emax = maxKinEnergyCSDA;
1971 
1972  for (size_t i=0; i<n; ++i) {
1973  G4double dedx = 0.0;
1974  pv = (*p)[i];
1975  if(pv) {
1976  dedx = pv->Value(emax, idxDEDXunRestricted);
1977  } else {
1978  pv = (*p)[(*theDensityIdx)[i]];
1979  if(pv) {
1980  dedx =
1981  pv->Value(emax, idxDEDXunRestricted)*(*theDensityFactor)[i];
1982  }
1983  }
1984  theDEDXAtMaxEnergy[i] = dedx;
1985  //G4cout << "i= " << i << " emax(MeV)= " << emax/MeV<< " dedx= "
1986  //<< dedx << G4endl;
1987  }
1988  }
1989 
1990  } else if(fRestricted == tType) {
1991  if(verb) {
1992  G4cout<< "G4VEnergyLossProcess::SetDEDXTable "
1993  << particle->GetParticleName()
1994  << " oldTable " << theDEDXTable << " newTable " << p
1995  << " ion " << theIonisationTable
1996  << " IsMaster " << isMaster
1997  << " " << GetProcessName() << G4endl;
1998  G4cout << (*p) << G4endl;
1999  }
2000  theDEDXTable = p;
2001  } else if(fSubRestricted == tType) {
2002  theDEDXSubTable = p;
2003  } else if(fIsIonisation == tType) {
2004  if(verb) {
2005  G4cout<< "G4VEnergyLossProcess::SetIonisationTable "
2006  << particle->GetParticleName()
2007  << " oldTable " << theDEDXTable << " newTable " << p
2008  << " ion " << theIonisationTable
2009  << " IsMaster " << isMaster
2010  << " " << GetProcessName() << G4endl;
2011  }
2012  theIonisationTable = p;
2013  } else if(fIsSubIonisation == tType) {
2014  theIonisationSubTable = p;
2015  }
2016 }
2017 
2018 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2019 
2021 {
2022  theCSDARangeTable = p;
2023 
2024  if(p) {
2025  size_t n = p->length();
2026  G4PhysicsVector* pv;
2027  G4double emax = maxKinEnergyCSDA;
2028 
2029  for (size_t i=0; i<n; ++i) {
2030  pv = (*p)[i];
2031  G4double rmax = 0.0;
2032  if(pv) { rmax = pv->Value(emax, idxCSDA); }
2033  else {
2034  pv = (*p)[(*theDensityIdx)[i]];
2035  if(pv) { rmax = pv->Value(emax, idxCSDA)/(*theDensityFactor)[i]; }
2036  }
2037  theRangeAtMaxEnergy[i] = rmax;
2038  //G4cout << "i= " << i << " Emax(MeV)= " << emax/MeV << " Rmax= "
2039  //<< rmax<< G4endl;
2040  }
2041  }
2042 }
2043 
2044 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2045 
2047 {
2048  theRangeTableForLoss = p;
2049  if(1 < verboseLevel) {
2050  G4cout << "### Set Range table " << p
2051  << " for " << particle->GetParticleName()
2052  << " and process " << GetProcessName() << G4endl;
2053  }
2054 }
2055 
2056 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2057 
2059 {
2060  theSecondaryRangeTable = p;
2061  if(1 < verboseLevel) {
2062  G4cout << "### Set SecondaryRange table " << p
2063  << " for " << particle->GetParticleName()
2064  << " and process " << GetProcessName() << G4endl;
2065  }
2066 }
2067 
2068 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2069 
2071 {
2072  theInverseRangeTable = p;
2073  if(1 < verboseLevel) {
2074  G4cout << "### Set InverseRange table " << p
2075  << " for " << particle->GetParticleName()
2076  << " and process " << GetProcessName() << G4endl;
2077  }
2078 }
2079 
2080 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2081 
2083 {
2084  if(1 < verboseLevel) {
2085  G4cout << "### Set Lambda table " << p
2086  << " for " << particle->GetParticleName()
2087  << " and process " << GetProcessName() << G4endl;
2088  }
2089  theLambdaTable = p;
2090  tablesAreBuilt = true;
2091 
2092  G4LossTableBuilder* bld = lManager->GetTableBuilder();
2093  theDensityFactor = bld->GetDensityFactors();
2094  theDensityIdx = bld->GetCoupleIndexes();
2095 
2096  if(theLambdaTable) {
2097  size_t n = theLambdaTable->length();
2098  G4PhysicsVector* pv = (*theLambdaTable)[0];
2099  G4double e, ss, smax, emax;
2100 
2101  size_t i;
2102 
2103  // first loop on existing vectors
2104  for (i=0; i<n; ++i) {
2105  pv = (*theLambdaTable)[i];
2106  if(pv) {
2107  size_t nb = pv->GetVectorLength();
2108  emax = DBL_MAX;
2109  smax = 0.0;
2110  if(nb > 0) {
2111  for (size_t j=0; j<nb; ++j) {
2112  e = pv->Energy(j);
2113  ss = (*pv)(j);
2114  if(ss > smax) {
2115  smax = ss;
2116  emax = e;
2117  }
2118  }
2119  }
2120  theEnergyOfCrossSectionMax[i] = emax;
2121  theCrossSectionMax[i] = smax;
2122  if(1 < verboseLevel) {
2123  G4cout << "For " << particle->GetParticleName()
2124  << " Max CS at i= " << i << " emax(MeV)= " << emax/MeV
2125  << " lambda= " << smax << G4endl;
2126  }
2127  }
2128  }
2129  // second loop using base materials
2130  for (i=0; i<n; ++i) {
2131  pv = (*theLambdaTable)[i];
2132  if(!pv){
2133  G4int j = (*theDensityIdx)[i];
2134  theEnergyOfCrossSectionMax[i] = theEnergyOfCrossSectionMax[j];
2135  theCrossSectionMax[i] = (*theDensityFactor)[i]*theCrossSectionMax[j];
2136  }
2137  }
2138  }
2139 }
2140 
2141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2142 
2144 {
2145  theSubLambdaTable = p;
2146  if(1 < verboseLevel) {
2147  G4cout << "### Set SebLambda table " << p
2148  << " for " << particle->GetParticleName()
2149  << " and process " << GetProcessName() << G4endl;
2150  }
2151 }
2152 
2153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2154 
2156 {
2157  const G4Element* elm = 0;
2158  if(currentModel) { elm = currentModel->GetCurrentElement(); }
2159  return elm;
2160 }
2161 
2162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2163 
2165  G4bool flag)
2166 {
2167  if(f > 0.0) {
2168  biasFactor = f;
2169  weightFlag = flag;
2170  if(1 < verboseLevel) {
2171  G4cout << "### SetCrossSectionBiasingFactor: for "
2172  << " process " << GetProcessName()
2173  << " biasFactor= " << f << " weightFlag= " << flag
2174  << G4endl;
2175  }
2176  }
2177 }
2178 
2179 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2180 
2181 void
2183  const G4String& region,
2184  G4bool flag)
2185 {
2186  if(!biasManager) { biasManager = new G4EmBiasingManager(); }
2187  if(1 < verboseLevel) {
2188  G4cout << "### ActivateForcedInteraction: for "
2189  << " process " << GetProcessName()
2190  << " length(mm)= " << length/mm
2191  << " in G4Region <" << region
2192  << "> weightFlag= " << flag
2193  << G4endl;
2194  }
2195  weightFlag = flag;
2196  biasManager->ActivateForcedInteraction(length, region);
2197 }
2198 
2199 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2200 
2201 void
2203  G4double factor,
2204  G4double energyLimit)
2205 {
2206  if (0.0 <= factor) {
2207 
2208  // Range cut can be applied only for e-
2209  if(0.0 == factor && secondaryParticle != G4Electron::Electron())
2210  { return; }
2211 
2212  if(!biasManager) { biasManager = new G4EmBiasingManager(); }
2213  biasManager->ActivateSecondaryBiasing(region, factor, energyLimit);
2214  if(1 < verboseLevel) {
2215  G4cout << "### ActivateSecondaryBiasing: for "
2216  << " process " << GetProcessName()
2217  << " factor= " << factor
2218  << " in G4Region <" << region
2219  << "> energyLimit(MeV)= " << energyLimit/MeV
2220  << G4endl;
2221  }
2222  }
2223 }
2224 
2225 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2226 
G4VEnergyLossProcess(const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
G4double condition(const G4ErrorSymMatrix &m)
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:448
G4ParticleDefinition * GetDefinition() const
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:538
void SetDEDXTable(G4PhysicsTable *p, G4EmTableType tType)
G4int GetParentID() const
G4bool SplineFlag() const
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:245
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
void InitializeForPostStep(const G4Track &)
const G4String & GetName() const
G4SafetyHelper * GetSafetyHelper() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
const std::vector< G4double > * GetDensityFactors()
static void SetPhysicsVector(G4PhysicsTable *physTable, size_t idx, G4PhysicsVector *vec)
void InitialiseHelper()
static G4LossTableManager * Instance()
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="")
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmax, G4double length, G4double meanLoss)=0
void AlongStepDeexcitation(std::vector< G4Track * > &tracks, const G4Step &step, G4double &eLoss, G4int coupleIndex)
G4double ApplySecondaryBiasing(std::vector< G4DynamicParticle * > &, const G4Track &track, G4VEmModel *currentModel, G4ParticleChangeForGamma *pParticleChange, G4double &eloss, G4int coupleIdx, G4double tcut, G4double safety=0.0)
G4int verboseLevel
Definition: G4VProcess.hh:368
tuple bin
Definition: plottest35.py:22
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
Definition: G4VEmModel.cc:339
void UpdateEmModel(const G4String &, G4double, G4double)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
G4double GetKineticEnergy() const
void PrintInfoDefinition(const G4ParticleDefinition &part)
G4double GetStepLength() const
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
void BuildPhysicsTable(const G4ParticleDefinition &)
const G4DynamicParticle * GetDynamicParticle() const
G4bool SecondaryBiasingRegion(G4int coupleIdx)
G4bool IsPIXEActive() const
G4PhysicsTable * SubLambdaTable() const
G4double z
Definition: TRTMaterials.hh:39
void SetLowEnergyLimit(G4double elimit)
G4VEmModel * GetModel(G4int, G4bool ver=false)
G4double GetProductionCut(G4int index) const
void DeRegister(G4VEnergyLossProcess *p)
G4PhysicsTable * RangeTableForLoss() const
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
const char * p
Definition: xmltok.h:285
G4StepStatus GetStepStatus() const
G4bool BuildCSDARange() const
const G4String & GetName() const
Definition: G4Material.hh:176
G4PhysicsTable * SecondaryRangeTable() const
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void AddCollaborativeProcess(G4VEnergyLossProcess *)
const G4ThreeVector & GetPosition() const
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:206
void ActivateForcedInteraction(G4double length=0.0, const G4String &region="", G4bool flag=true)
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
void SetStepFunction(G4double v1, G4double v2)
G4bool ForcedInteractionRegion(G4int coupleIdx)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4PhysicsTable * CSDARangeTable() const
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
G4ParticleChangeForLoss fParticleChange
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
void SelectModel(G4double kinEnergy)
G4bool GetFlag(size_t idx) const
const XML_Char * name
void SetTouchableHandle(const G4TouchableHandle &apValue)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *, const G4Region *)
G4VEmModel * GetModelByIndex(G4int idx=0, G4bool ver=false) const
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4PhysicsTable * IonisationTableForSubsec() const
G4double GetParentWeight() const
#define G4BestUnit(a, b)
#define G4_USE_G4BESTUNIT_FOR_VERBOSE 1
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:571
const G4DataVector * SubCutoff() const
size_t GetVectorLength() const
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *, G4double cut)
void InitializeForAlongStep(const G4Track &)
G4double GetProposedKineticEnergy() const
G4PhysicsTable * IonisationTable() const
void FillLambdaVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4bool startFromNull=true, G4EmTableType t=fRestricted)
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
void FillSecondDerivatives()
G4PhysicsTable * BuildLambdaTable(G4EmTableType tType=fRestricted)
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:186
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
const G4int smax
const G4String & GetParticleName() const
G4LossTableBuilder * GetTableBuilder()
void SetWeight(G4double aValue)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetCreatorModelIndex(G4int idx)
Definition: G4Track.hh:245
G4bool ExistPhysicsTable(const G4String &fileName) const
void SetSpline(G4bool)
void Initialise(const G4ParticleDefinition &part, const G4String &procName, G4int verbose)
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, G4double tmax, G4double length)=0
static G4RegionStore * GetInstance()
G4StepPoint * GetPreStepPoint() const
G4PhysicsTable * LambdaTable() const
void ProposeWeight(G4double finalWeight)
void SetSecondaryWeightByProcess(G4bool)
void SetInverseRangeTable(G4PhysicsTable *p)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4double GetKineticEnergy() const
G4double SampleSubCutSecondaries(std::vector< G4Track * > &, const G4Step &, G4VEmModel *model, G4int matIdx)
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
void FillDEDXVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4EmTableType t=fRestricted)
G4PhysicsTable * DEDXTable() const
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=0)
Definition: G4VEmModel.cc:387
const G4ThreeVector & GetPosition() const
const G4Element * GetCurrentElement() const
bool G4bool
Definition: G4Types.hh:79
G4double GetStepLimit(G4int coupleIdx, G4double previousStep)
G4double currentInteractionLength
Definition: G4VProcess.hh:297
void SetProposedKineticEnergy(G4double proposedKinEnergy)
const G4ParticleDefinition * GetParticleDefinition() const
const XML_Char XML_Content * model
float proton_mass_c2
Definition: hepunit.py:275
const G4String & GetParticleType() const
void Register(G4VEnergyLossProcess *p)
float electron_mass_c2
Definition: hepunit.py:274
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &)
Definition: G4Step.hh:76
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
const G4int n
G4double GetGlobalTime() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4PhysicsTable * DEDXTableForSubsec() const
G4double Energy(size_t index) const
G4double Value(G4double theEnergy, size_t &lastidx) const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
void AddSecondary(G4Track *aSecondary)
virtual G4double ChargeSquareRatio(const G4Track &)
Definition: G4VEmModel.cc:314
const G4TouchableHandle & GetTouchableHandle() const
string pname
Definition: eplot.py:33
G4Material * GetMaterial() const
void RegisterExtraParticle(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
G4bool IsActive(G4double kinEnergy)
Definition: G4VEmModel.hh:711
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()
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:669
void PreparePhysicsTable(const G4ParticleDefinition &)
size_t length() const
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4int size() const
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:331
static G4ProductionCutsTable * GetProductionCutsTable()
void SetLambdaTable(G4PhysicsTable *p)
#define fm
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4int NumberOfModels() const
virtual void PrintInfo()=0
G4PhysicsTable * InverseRangeTable() const
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
G4double GetPDGMass() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ParticleTable * GetParticleTable()
int G4lrint(double ad)
Definition: templates.hh:163
void InitialiseBaseMaterials(G4PhysicsTable *table)
void SetNumberOfSecondaries(G4int totSecondaries)
G4double GetLocalEnergyDeposit() const
G4double MeanFreePath(const G4Track &track)
G4StepPoint * GetPostStepPoint() const
void UpdateEmModel(const G4String &, G4double, G4double)
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
void SetEmModel(G4VEmModel *, G4int index=1)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double ContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
G4double GetSafety() const
void SetProposedCharge(G4double theCharge)
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
void SetCSDARangeTable(G4PhysicsTable *pRange)
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
void LocalPhysicsTables(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
G4EmTableType
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
G4double GetGlobalTime() const
void SetSecondaryRangeTable(G4PhysicsTable *p)
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
G4TrackStatus GetTrackStatus() const
G4VAtomDeexcitation * AtomDeexcitation()
G4PhysicsTable * DEDXunRestrictedTable() const
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:300
G4bool GetFlag(size_t i) const
static G4int Register(G4String &)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetSubLambdaTable(G4PhysicsTable *p)
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *, G4double cut)
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4ForceCondition
G4Track * GetTrack() const
G4double GetPDGCharge() const
void ActivateSubCutoff(G4bool val, const G4Region *region=0)
void SetRangeTableForLoss(G4PhysicsTable *p)
static G4bool RetrievePhysicsTable(G4PhysicsTable *physTable, const G4String &fileName, G4bool ascii)
G4ProductionCuts * GetProductionCuts() const
#define DBL_MAX
Definition: templates.hh:83
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
G4PhysicsTable * BuildDEDXTable(G4EmTableType tType=fRestricted)
G4VEmModel * EmModel(G4int index=1) const
void clearAndDestroy()
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:437
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *)=0
G4GPILSelection
const G4Material * GetMaterial() const
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple)
const std::vector< G4int > * GetCoupleIndexes()
G4ProcessType
float perMillion
Definition: hepunit.py:241
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:440
G4bool IsIonisationProcess() const