Geant4  9.6.p02
 All Classes 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$
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 // 15-10-10 Fixed 4-momentum balance if deexcitation is active (L.Pandola)
110 // 30-05-12 Call new ApplySecondaryBiasing so 2ries may be unique (D. Sawkey)
111 // 30-05-12 Fix bug in forced biasing: now called on first step (D. Sawkey)
112 //
113 // Class Description:
114 //
115 // It is the unified energy loss process it calculates the continuous
116 // energy loss for charged particles using a set of Energy Loss
117 // models valid for different energy regions. There are a possibility
118 // to create and access to dE/dx and range tables, or to calculate
119 // that information on fly.
120 // -------------------------------------------------------------------
121 //
122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
123 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
124 
125 #include "G4VEnergyLossProcess.hh"
126 #include "G4PhysicalConstants.hh"
127 #include "G4SystemOfUnits.hh"
128 #include "G4ProcessManager.hh"
129 #include "G4LossTableManager.hh"
130 #include "G4LossTableBuilder.hh"
131 #include "G4Step.hh"
132 #include "G4ParticleDefinition.hh"
133 #include "G4VEmModel.hh"
134 #include "G4VEmFluctuationModel.hh"
135 #include "G4DataVector.hh"
136 #include "G4PhysicsLogVector.hh"
137 #include "G4VParticleChange.hh"
138 #include "G4Gamma.hh"
139 #include "G4Electron.hh"
140 #include "G4Positron.hh"
141 #include "G4Proton.hh"
142 #include "G4ProcessManager.hh"
143 #include "G4UnitsTable.hh"
144 #include "G4GenericIon.hh"
145 #include "G4ProductionCutsTable.hh"
146 #include "G4Region.hh"
147 #include "G4RegionStore.hh"
148 #include "G4PhysicsTableHelper.hh"
149 #include "G4SafetyHelper.hh"
151 #include "G4EmConfigurator.hh"
152 #include "G4VAtomDeexcitation.hh"
153 #include "G4EmBiasingManager.hh"
154 
155 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
156 
158  G4ProcessType type):
159  G4VContinuousDiscreteProcess(name, type),
160  secondaryParticle(0),
161  nSCoffRegions(0),
162  idxSCoffRegions(0),
163  nProcesses(0),
164  theDEDXTable(0),
165  theDEDXSubTable(0),
166  theDEDXunRestrictedTable(0),
167  theIonisationTable(0),
168  theIonisationSubTable(0),
169  theRangeTableForLoss(0),
170  theCSDARangeTable(0),
171  theSecondaryRangeTable(0),
172  theInverseRangeTable(0),
173  theLambdaTable(0),
174  theSubLambdaTable(0),
175  theDensityFactor(0),
176  theDensityIdx(0),
177  baseParticle(0),
178  minSubRange(0.1),
179  lossFluctuationFlag(true),
180  rndmStepFlag(false),
181  tablesAreBuilt(false),
182  integral(true),
183  isIon(false),
184  isIonisation(true),
185  useSubCutoff(false),
186  useDeexcitation(false),
187  particle(0),
188  currentCouple(0),
189  nWarnings(0),
190  mfpKinEnergy(0.0)
191 {
192  SetVerboseLevel(1);
193 
194  // low energy limit
195  lowestKinEnergy = 1.*eV;
196 
197  // Size of tables assuming spline
198  minKinEnergy = 0.1*keV;
199  maxKinEnergy = 10.0*TeV;
200  nBins = 77;
201  maxKinEnergyCSDA = 1.0*GeV;
202  nBinsCSDA = 35;
203 
204  // default linear loss limit for spline
205  linLossLimit = 0.01;
206 
207  // default dRoverRange and finalRange
208  SetStepFunction(0.2, 1.0*mm);
209 
210  // default lambda factor
211  lambdaFactor = 0.8;
212 
213  // cross section biasing
214  biasFactor = 1.0;
215 
216  // particle types
217  theElectron = G4Electron::Electron();
218  thePositron = G4Positron::Positron();
219  theGamma = G4Gamma::Gamma();
220  theGenericIon = 0;
221 
222  // run time objects
225  modelManager = new G4EmModelManager();
227  ->GetSafetyHelper();
228  aGPILSelection = CandidateForSelection;
229 
230  // initialise model
231  (G4LossTableManager::Instance())->Register(this);
232  fluctModel = 0;
233  atomDeexcitation = 0;
234 
235  biasManager = 0;
236  biasFlag = false;
237  weightFlag = false;
238 
239  scTracks.reserve(5);
240  secParticles.reserve(5);
241 }
242 
243 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
244 
246 {
247  if(1 < verboseLevel) {
248  G4cout << "G4VEnergyLossProcess destruct " << GetProcessName()
249  << " " << this << " " << baseParticle << G4endl;
250  }
251  Clean();
252 
253  if ( !baseParticle ) {
254  if(theDEDXTable) {
255  if(theIonisationTable == theDEDXTable) { theIonisationTable = 0; }
256  theDEDXTable->clearAndDestroy();
257  delete theDEDXTable;
258  if(theDEDXSubTable) {
259  if(theIonisationSubTable == theDEDXSubTable)
260  theIonisationSubTable = 0;
261  theDEDXSubTable->clearAndDestroy();
262  delete theDEDXSubTable;
263  }
264  }
265  if(theIonisationTable) {
266  theIonisationTable->clearAndDestroy();
267  delete theIonisationTable;
268  }
269  if(theIonisationSubTable) {
270  theIonisationSubTable->clearAndDestroy();
271  delete theIonisationSubTable;
272  }
273  if(theDEDXunRestrictedTable && theCSDARangeTable) {
274  theDEDXunRestrictedTable->clearAndDestroy();
275  delete theDEDXunRestrictedTable;
276  }
277  if(theCSDARangeTable) {
278  theCSDARangeTable->clearAndDestroy();
279  delete theCSDARangeTable;
280  }
281  if(theRangeTableForLoss) {
282  theRangeTableForLoss->clearAndDestroy();
283  delete theRangeTableForLoss;
284  }
285  if(theInverseRangeTable) {
286  theInverseRangeTable->clearAndDestroy();
287  delete theInverseRangeTable;
288  }
289  if(theLambdaTable) {
290  theLambdaTable->clearAndDestroy();
291  delete theLambdaTable;
292  }
293  if(theSubLambdaTable) {
294  theSubLambdaTable->clearAndDestroy();
295  delete theSubLambdaTable;
296  }
297  }
298 
299  delete modelManager;
300  delete biasManager;
301  (G4LossTableManager::Instance())->DeRegister(this);
302 }
303 
304 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
305 
306 void G4VEnergyLossProcess::Clean()
307 {
308  if(1 < verboseLevel) {
309  G4cout << "G4VEnergyLossProcess::Clear() for " << GetProcessName()
310  << G4endl;
311  }
312  delete [] idxSCoffRegions;
313 
314  tablesAreBuilt = false;
315 
316  scProcesses.clear();
317  nProcesses = 0;
318 }
319 
320 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
321 
323  const G4Material*,
324  G4double cut)
325 {
326  return cut;
327 }
328 
329 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
330 
332  G4VEmFluctuationModel* fluc,
333  const G4Region* region)
334 {
335  modelManager->AddEmModel(order, p, fluc, region);
336  if(p) { p->SetParticleChange(pParticleChange, fluc); }
337 }
338 
339 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
340 
342  G4double emin, G4double emax)
343 {
344  modelManager->UpdateEmModel(nam, emin, emax);
345 }
346 
347 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
348 
350 {
351  G4int n = emModels.size();
352  if(index >= n) { for(G4int i=n; i<=index; ++i) {emModels.push_back(0);} }
353  emModels[index] = p;
354 }
355 
356 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
357 
359 {
360  G4VEmModel* p = 0;
361  if(index >= 0 && index < G4int(emModels.size())) { p = emModels[index]; }
362  return p;
363 }
364 
365 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
366 
368 {
369  return modelManager->GetModel(idx, ver);
370 }
371 
372 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
373 
375 {
376  return modelManager->NumberOfModels();
377 }
378 
379 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
380 
381 void
383 {
384  if(1 < verboseLevel) {
385  G4cout << "G4VEnergyLossProcess::PreparePhysicsTable for "
386  << GetProcessName() << " for " << part.GetParticleName()
387  << " " << this << G4endl;
388  }
389 
390  currentCouple = 0;
391  preStepLambda = 0.0;
392  mfpKinEnergy = DBL_MAX;
393  fRange = DBL_MAX;
394  preStepKinEnergy = 0.0;
395  chargeSqRatio = 1.0;
396  massRatio = 1.0;
397  reduceFactor = 1.0;
398  fFactor = 1.0;
399 
401 
402  // Are particle defined?
403  if( !particle ) { particle = &part; }
404 
405  if(part.GetParticleType() == "nucleus") {
406 
407  G4String pname = part.GetParticleName();
408  if(pname != "deuteron" && pname != "triton" &&
409  pname != "alpha+" && pname != "helium" &&
410  pname != "hydrogen") {
411 
412  theGenericIon = G4GenericIon::GenericIon();
413  isIon = true;
414  // process is shared between all ions inheriting G4GenericIon
415  // for all excluding He3 and alpha
416  if(pname != "He3" && pname != "alpha") { particle = theGenericIon; }
417  }
418  }
419 
420  if( particle != &part ) {
421  if(isIon) {
422  lManager->RegisterIon(&part, this);
423  } else {
424  lManager->RegisterExtraParticle(&part, this);
425  }
426  if(1 < verboseLevel) {
427  G4cout << "### G4VEnergyLossProcess::PreparePhysicsTable() interrupted for "
428  << part.GetParticleName() << " isIon= " << isIon
429  << " particle " << particle << " GenericIon " << theGenericIon
430  << G4endl;
431  }
432  return;
433  }
434 
435  Clean();
436  lManager->PreparePhysicsTable(&part, this);
437  G4LossTableBuilder* bld = lManager->GetTableBuilder();
438 
439  // Base particle and set of models can be defined here
440  InitialiseEnergyLossProcess(particle, baseParticle);
441 
442  const G4ProductionCutsTable* theCoupleTable=
444  size_t n = theCoupleTable->GetTableSize();
445 
446  theDEDXAtMaxEnergy.resize(n, 0.0);
447  theRangeAtMaxEnergy.resize(n, 0.0);
448  theEnergyOfCrossSectionMax.resize(n, 0.0);
449  theCrossSectionMax.resize(n, DBL_MAX);
450 
451  // Tables preparation
452  if (!baseParticle) {
453 
454  theDEDXTable = G4PhysicsTableHelper::PreparePhysicsTable(theDEDXTable);
455  bld->InitialiseBaseMaterials(theDEDXTable);
456 
457  if (lManager->BuildCSDARange()) {
458  theDEDXunRestrictedTable =
459  G4PhysicsTableHelper::PreparePhysicsTable(theDEDXunRestrictedTable);
460  theCSDARangeTable =
462  }
463 
464  theLambdaTable = G4PhysicsTableHelper::PreparePhysicsTable(theLambdaTable);
465  bld->InitialiseBaseMaterials(theLambdaTable);
466 
467  if(isIonisation) {
468  theRangeTableForLoss =
469  G4PhysicsTableHelper::PreparePhysicsTable(theRangeTableForLoss);
470  theInverseRangeTable =
471  G4PhysicsTableHelper::PreparePhysicsTable(theInverseRangeTable);
472  }
473 
474  if (nSCoffRegions) {
475  theDEDXSubTable =
477  theSubLambdaTable =
479  }
480  }
481 
482  theDensityFactor = bld->GetDensityFactors();
483  theDensityIdx = bld->GetCoupleIndexes();
484 
485  // forced biasing
486  if(biasManager) {
487  biasManager->Initialise(part,GetProcessName(),verboseLevel);
488  biasFlag = false;
489  }
490 
491  G4double initialCharge = particle->GetPDGCharge();
492  G4double initialMass = particle->GetPDGMass();
493 
494  if (baseParticle) {
495  massRatio = (baseParticle->GetPDGMass())/initialMass;
496  G4double q = initialCharge/baseParticle->GetPDGCharge();
497  chargeSqRatio = q*q;
498  if(chargeSqRatio > 0.0) { reduceFactor = 1.0/(chargeSqRatio*massRatio); }
499  }
500 
501  // initialisation of models
502  G4int nmod = modelManager->NumberOfModels();
503  for(G4int i=0; i<nmod; ++i) {
504  G4VEmModel* mod = modelManager->GetModel(i);
505  if(mod->HighEnergyLimit() > maxKinEnergy) {
506  mod->SetHighEnergyLimit(maxKinEnergy);
507  }
508  }
509 
510  theCuts = modelManager->Initialise(particle, secondaryParticle,
511  minSubRange, verboseLevel);
512 
513  // Sub Cutoff
514  if (nSCoffRegions>0) {
515  theSubCuts = modelManager->SubCutoff();
516 
517  if(nSCoffRegions>0) { idxSCoffRegions = new G4bool[n]; }
518  for (size_t j=0; j<n; ++j) {
519 
520  const G4MaterialCutsCouple* couple =
521  theCoupleTable->GetMaterialCutsCouple(j);
522  const G4ProductionCuts* pcuts = couple->GetProductionCuts();
523 
524  if(nSCoffRegions>0) {
525  G4bool reg = false;
526  for(G4int i=0; i<nSCoffRegions; ++i) {
527  if( pcuts == scoffRegions[i]->GetProductionCuts()) { reg = true; }
528  }
529  idxSCoffRegions[j] = reg;
530  }
531  }
532  }
533 
534  if(1 < verboseLevel) {
535  G4cout << "G4VEnergyLossProcess::Initialise() is done "
536  << " for local " << particle->GetParticleName()
537  << " isIon= " << isIon;
538  if(baseParticle) { G4cout << "; base: " << baseParticle->GetParticleName(); }
539  G4cout << " chargeSqRatio= " << chargeSqRatio
540  << " massRatio= " << massRatio
541  << " reduceFactor= " << reduceFactor << G4endl;
542  if (nSCoffRegions) {
543  G4cout << " SubCutoff Regime is ON for regions: " << G4endl;
544  for (G4int i=0; i<nSCoffRegions; ++i) {
545  const G4Region* r = scoffRegions[i];
546  G4cout << " " << r->GetName() << G4endl;
547  }
548  }
549  }
550 }
551 
552 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
553 
555 {
556  if(1 < verboseLevel) {
557  G4cout << "### G4VEnergyLossProcess::BuildPhysicsTable() for "
558  << GetProcessName()
559  << " and particle " << part.GetParticleName()
560  << "; local: " << particle->GetParticleName();
561  if(baseParticle) { G4cout << "; base: " << baseParticle->GetParticleName(); }
562  G4cout << " TablesAreBuilt= " << tablesAreBuilt
563  << " isIon= " << isIon << " " << this << G4endl;
564  }
565 
566  if(&part == particle) {
567  if(!tablesAreBuilt) {
569  }
570  if(!baseParticle) {
571  // needs to be done only once
572  safetyHelper->InitialiseHelper();
573  }
574  }
575 
576  // explicitly defined printout by particle name
577  G4String num = part.GetParticleName();
578  if(1 < verboseLevel ||
579  (0 < verboseLevel && (num == "e-" ||
580  num == "e+" || num == "mu+" ||
581  num == "mu-" || num == "proton"||
582  num == "pi+" || num == "pi-" ||
583  num == "kaon+" || num == "kaon-" ||
584  num == "alpha" || num == "anti_proton" ||
585  num == "GenericIon")))
586  {
587  particle = &part;
589  }
590 
591  // Added tracking cut to avoid tracking artifacts
592  // identify deexcitation flag
593  if(isIonisation) {
594  fParticleChange.SetLowEnergyLimit(lowestKinEnergy);
595  atomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
596  if(atomDeexcitation) {
597  if(atomDeexcitation->IsPIXEActive()) { useDeexcitation = true; }
598  }
599  }
600 
601  if(1 < verboseLevel) {
602  G4cout << "### G4VEnergyLossProcess::BuildPhysicsTable() done for "
603  << GetProcessName()
604  << " and particle " << part.GetParticleName();
605  if(isIonisation) { G4cout << " isIonisation flag = 1"; }
606  G4cout << G4endl;
607  }
608 }
609 
610 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
611 
613 {
614  if(1 < verboseLevel) {
615  G4cout << "G4VEnergyLossProcess::BuildDEDXTable() of type " << tType
616  << " for " << GetProcessName()
617  << " and particle " << particle->GetParticleName()
618  << G4endl;
619  }
620  G4PhysicsTable* table = 0;
621  G4double emax = maxKinEnergy;
622  G4int bin = nBins;
623 
624  if(fTotal == tType) {
625  emax = maxKinEnergyCSDA;
626  bin = nBinsCSDA;
627  table = theDEDXunRestrictedTable;
628  } else if(fRestricted == tType) {
629  table = theDEDXTable;
630  if(theIonisationTable)
631  table = G4PhysicsTableHelper::PreparePhysicsTable(theIonisationTable);
632  } else if(fSubRestricted == tType) {
633  table = theDEDXSubTable;
634  if(theIonisationSubTable)
635  table = G4PhysicsTableHelper::PreparePhysicsTable(theIonisationSubTable);
636  } else {
637  G4cout << "G4VEnergyLossProcess::BuildDEDXTable WARNING: wrong type "
638  << tType << G4endl;
639  }
640 
641  // Access to materials
642  const G4ProductionCutsTable* theCoupleTable=
644  size_t numOfCouples = theCoupleTable->GetTableSize();
645 
646  if(1 < verboseLevel) {
647  G4cout << numOfCouples << " materials"
648  << " minKinEnergy= " << minKinEnergy
649  << " maxKinEnergy= " << emax
650  << " nbin= " << bin
651  << " EmTableType= " << tType
652  << " table= " << table << " " << this
653  << G4endl;
654  }
655  if(!table) return table;
656 
657  G4LossTableBuilder* bld = (G4LossTableManager::Instance())->GetTableBuilder();
658  G4bool splineFlag = (G4LossTableManager::Instance())->SplineFlag();
659  G4PhysicsLogVector* aVector = 0;
660  G4PhysicsLogVector* bVector = 0;
661 
662  for(size_t i=0; i<numOfCouples; ++i) {
663 
664  if(1 < verboseLevel) {
665  G4cout << "G4VEnergyLossProcess::BuildDEDXVector flagTable= "
666  << table->GetFlag(i) << " Flag= " << bld->GetFlag(i) << G4endl;
667  }
668  if(bld->GetFlag(i)) {
669 
670  // create physics vector and fill it
671  const G4MaterialCutsCouple* couple =
672  theCoupleTable->GetMaterialCutsCouple(i);
673  delete (*table)[i];
674  if(!bVector) {
675  aVector = new G4PhysicsLogVector(minKinEnergy, emax, bin);
676  bVector = aVector;
677  } else {
678  aVector = new G4PhysicsLogVector(*bVector);
679  }
680  aVector->SetSpline(splineFlag);
681 
682  modelManager->FillDEDXVector(aVector, couple, tType);
683  if(splineFlag) { aVector->FillSecondDerivatives(); }
684 
685  // Insert vector for this material into the table
686  G4PhysicsTableHelper::SetPhysicsVector(table, i, aVector);
687  }
688  }
689 
690  if(1 < verboseLevel) {
691  G4cout << "G4VEnergyLossProcess::BuildDEDXTable(): table is built for "
692  << particle->GetParticleName()
693  << " and process " << GetProcessName()
694  << G4endl;
695  // if(2 < verboseLevel) G4cout << (*table) << G4endl;
696  }
697 
698  return table;
699 }
700 
701 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
702 
704 {
705  G4PhysicsTable* table = 0;
706 
707  if(fRestricted == tType) {
708  table = theLambdaTable;
709  } else if(fSubRestricted == tType) {
710  table = theSubLambdaTable;
711  } else {
712  G4cout << "G4VEnergyLossProcess::BuildLambdaTable WARNING: wrong type "
713  << tType << G4endl;
714  }
715 
716  if(1 < verboseLevel) {
717  G4cout << "G4VEnergyLossProcess::BuildLambdaTable() of type "
718  << tType << " for process "
719  << GetProcessName() << " and particle "
720  << particle->GetParticleName()
721  << " EmTableType= " << tType
722  << " table= " << table
723  << G4endl;
724  }
725  if(!table) {return table;}
726 
727  // Access to materials
728  const G4ProductionCutsTable* theCoupleTable=
730  size_t numOfCouples = theCoupleTable->GetTableSize();
731 
732  G4LossTableBuilder* bld = (G4LossTableManager::Instance())->GetTableBuilder();
733  G4bool splineFlag = (G4LossTableManager::Instance())->SplineFlag();
734  G4PhysicsLogVector* aVector = 0;
735  G4double scale = std::log(maxKinEnergy/minKinEnergy);
736 
737  for(size_t i=0; i<numOfCouples; ++i) {
738 
739  if (bld->GetFlag(i)) {
740 
741  // create physics vector and fill it
742  const G4MaterialCutsCouple* couple =
743  theCoupleTable->GetMaterialCutsCouple(i);
744  delete (*table)[i];
745  G4double emin = MinPrimaryEnergy(particle,couple->GetMaterial(),(*theCuts)[i]);
746  if(0.0 >= emin) { emin = eV; }
747  else if(maxKinEnergy <= emin) { emin = 0.5*maxKinEnergy; }
748  G4int bin = G4int(nBins*std::log(maxKinEnergy/emin)/scale + 0.5);
749  if(bin < 3) { bin = 3; }
750  aVector = new G4PhysicsLogVector(emin, maxKinEnergy, bin);
751  aVector->SetSpline(splineFlag);
752 
753  modelManager->FillLambdaVector(aVector, couple, true, tType);
754  if(splineFlag) { aVector->FillSecondDerivatives(); }
755 
756  // Insert vector for this material into the table
757  G4PhysicsTableHelper::SetPhysicsVector(table, i, aVector);
758  }
759  }
760 
761  if(1 < verboseLevel) {
762  G4cout << "Lambda table is built for "
763  << particle->GetParticleName()
764  << G4endl;
765  }
766 
767  return table;
768 }
769 
770 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
771 
773 {
774  if(0 < verboseLevel) {
775  G4cout << G4endl << GetProcessName() << ": for "
776  << particle->GetParticleName()
777  << " SubType= " << GetProcessSubType()
778  << G4endl
779  << " dE/dx and range tables from "
780  << G4BestUnit(minKinEnergy,"Energy")
781  << " to " << G4BestUnit(maxKinEnergy,"Energy")
782  << " in " << nBins << " bins" << G4endl
783  << " Lambda tables from threshold to "
784  << G4BestUnit(maxKinEnergy,"Energy")
785  << " in " << nBins << " bins, spline: "
786  << (G4LossTableManager::Instance())->SplineFlag()
787  << G4endl;
788  if(theRangeTableForLoss && isIonisation) {
789  G4cout << " finalRange(mm)= " << finalRange/mm
790  << ", dRoverRange= " << dRoverRange
791  << ", integral: " << integral
792  << ", fluct: " << lossFluctuationFlag
793  << ", linLossLimit= " << linLossLimit
794  << G4endl;
795  }
796  PrintInfo();
797  modelManager->DumpModelList(verboseLevel);
798  if(theCSDARangeTable && isIonisation) {
799  G4cout << " CSDA range table up"
800  << " to " << G4BestUnit(maxKinEnergyCSDA,"Energy")
801  << " in " << nBinsCSDA << " bins" << G4endl;
802  }
803  if(nSCoffRegions>0 && isIonisation) {
804  G4cout << " Subcutoff sampling in " << nSCoffRegions
805  << " regions" << G4endl;
806  }
807  if(2 < verboseLevel) {
808  G4cout << " DEDXTable address= " << theDEDXTable << G4endl;
809  if(theDEDXTable && isIonisation) G4cout << (*theDEDXTable) << G4endl;
810  G4cout << "non restricted DEDXTable address= "
811  << theDEDXunRestrictedTable << G4endl;
812  if(theDEDXunRestrictedTable && isIonisation) {
813  G4cout << (*theDEDXunRestrictedTable) << G4endl;
814  }
815  if(theDEDXSubTable && isIonisation) {
816  G4cout << (*theDEDXSubTable) << G4endl;
817  }
818  G4cout << " CSDARangeTable address= " << theCSDARangeTable
819  << G4endl;
820  if(theCSDARangeTable && isIonisation) {
821  G4cout << (*theCSDARangeTable) << G4endl;
822  }
823  G4cout << " RangeTableForLoss address= " << theRangeTableForLoss
824  << G4endl;
825  if(theRangeTableForLoss && isIonisation) {
826  G4cout << (*theRangeTableForLoss) << G4endl;
827  }
828  G4cout << " InverseRangeTable address= " << theInverseRangeTable
829  << G4endl;
830  if(theInverseRangeTable && isIonisation) {
831  G4cout << (*theInverseRangeTable) << G4endl;
832  }
833  G4cout << " LambdaTable address= " << theLambdaTable << G4endl;
834  if(theLambdaTable && isIonisation) {
835  G4cout << (*theLambdaTable) << G4endl;
836  }
837  G4cout << " SubLambdaTable address= " << theSubLambdaTable << G4endl;
838  if(theSubLambdaTable && isIonisation) {
839  G4cout << (*theSubLambdaTable) << G4endl;
840  }
841  }
842  }
843 }
844 
845 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
846 
848 {
849  G4RegionStore* regionStore = G4RegionStore::GetInstance();
850  const G4Region* reg = r;
851  if (!reg) {reg = regionStore->GetRegion("DefaultRegionForTheWorld", false);}
852 
853  // the region is in the list
854  if (nSCoffRegions) {
855  for (G4int i=0; i<nSCoffRegions; ++i) {
856  if (reg == scoffRegions[i]) {
857  return;
858  }
859  }
860  }
861 
862  // new region
863  if(val) {
864  useSubCutoff = true;
865  scoffRegions.push_back(reg);
866  ++nSCoffRegions;
867  } else {
868  useSubCutoff = false;
869  }
870 }
871 
872 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
873 
875 {
876  // reset parameters for the new track
878  mfpKinEnergy = DBL_MAX;
879 
880  // reset ion
881  if(isIon) {
882  chargeSqRatio = 0.5;
883 
884  G4double newmass = track->GetDefinition()->GetPDGMass();
885  if(baseParticle) {
886  massRatio = baseParticle->GetPDGMass()/newmass;
887  } else {
888  massRatio = proton_mass_c2/newmass;
889  }
890  }
891  // forced biasing only for primary particles
892  if(biasManager) {
893  if(0 == track->GetParentID()) {
894  // primary particle
895  biasFlag = true;
896  biasManager->ResetForcedInteraction();
897  }
898  }
899 }
900 
901 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
902 
905  G4GPILSelection* selection)
906 {
907  G4double x = DBL_MAX;
908  *selection = aGPILSelection;
909  if(isIonisation) {
910  fRange = GetScaledRangeForScaledEnergy(preStepScaledEnergy)*reduceFactor;
911 
912  x = fRange;
913  G4double y = x*dRoverRange;
914  G4double finR = finalRange;
915  if(rndmStepFlag) {
916  finR = std::min(finR,currentCouple->GetProductionCuts()->GetProductionCut(1));
917  }
918  if(x > finR) { x = y + finR*(1.0 - dRoverRange)*(2.0 - finR/fRange); }
919  /*
920  if(particle->GetPDGMass() > 0.9*GeV)
921  G4cout<<GetProcessName()<<": e= "<<preStepKinEnergy
922  <<" range= "<<fRange << " idx= " << basedCoupleIndex
923  << " y= " << y << " finR= " << finR
924  << " limit= " << x <<G4endl;
925  */
926  }
927  //G4cout<<GetProcessName()<<": e= "<<preStepKinEnergy <<" stepLimit= "<<x<<G4endl;
928  return x;
929 }
930 
931 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
932 
934  const G4Track& track,
935  G4double previousStepSize,
937 {
938  // condition is set to "Not Forced"
939  *condition = NotForced;
940  G4double x = DBL_MAX;
941 
942  // initialisation of material, mass, charge, model at the beginning of the step
943  /*
944  if(!theDensityFactor || !theDensityIdx) {
945  G4cout << "G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength 1: "
946  << theDensityFactor << " " << theDensityIdx
947  << G4endl;
948  G4cout << track.GetDefinition()->GetParticleName()
949  << " e(MeV)= " << track.GetKineticEnergy()
950  << " mat " << track.GetMaterialCutsCouple()->GetMaterial()->GetName()
951  << G4endl;
952  }
953  */
954  DefineMaterial(track.GetMaterialCutsCouple());
955  preStepKinEnergy = track.GetKineticEnergy();
956  preStepScaledEnergy = preStepKinEnergy*massRatio;
957  SelectModel(preStepScaledEnergy);
958 
959  if(!currentModel->IsActive(preStepScaledEnergy)) { return x; }
960 
961  // change effective charge of an ion on fly
962  if(isIon) {
963  G4double q2 = currentModel->ChargeSquareRatio(track);
964  if(q2 != chargeSqRatio) {
965  chargeSqRatio = q2;
966  fFactor = q2*biasFactor*(*theDensityFactor)[currentCoupleIndex];
967  reduceFactor = 1.0/(fFactor*massRatio);
968  }
969  }
970  // if(particle->GetPDGMass() > 0.9*GeV)
971  //G4cout << "q2= " << chargeSqRatio << " massRatio= " << massRatio << G4endl;
972  // initialisation for sampling of the interaction length
973  //if(previousStepSize <= 0.0) { theNumberOfInteractionLengthLeft = -1.0; }
974  //if(theNumberOfInteractionLengthLeft < 0.0) { mfpKinEnergy = DBL_MAX; }
975 
976  // forced biasing only for primary particles
977  if(biasManager) {
978  if(0 == track.GetParentID()) {
979  if(biasFlag && biasManager->ForcedInteractionRegion(currentCoupleIndex)) {
980  return biasManager->GetStepLimit(currentCoupleIndex, previousStepSize);
981  }
982  }
983  }
984 
985  // compute mean free path
986  if(preStepScaledEnergy < mfpKinEnergy) {
987  if (integral) { ComputeLambdaForScaledEnergy(preStepScaledEnergy); }
988  else { preStepLambda = GetLambdaForScaledEnergy(preStepScaledEnergy); }
989 
990  // zero cross section
991  if(preStepLambda <= 0.0) {
994  }
995  }
996 
997  // non-zero cross section
998  if(preStepLambda > 0.0) {
1000 
1001  // beggining of tracking (or just after DoIt of this process)
1003 
1004  } else if(currentInteractionLength < DBL_MAX) {
1005 
1006  // subtract NumberOfInteractionLengthLeft using previous step
1008  // SubtractNumberOfInteractionLengthLeft(previousStepSize);
1011  //theNumberOfInteractionLengthLeft = perMillion;
1012  }
1013  }
1014 
1015  // new mean free path and step limit
1016  currentInteractionLength = 1.0/preStepLambda;
1018 
1019 #ifdef G4VERBOSE
1020  if (verboseLevel>2){
1021  // if(particle->GetPDGMass() > 0.9*GeV){
1022  G4cout << "G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength ";
1023  G4cout << "[ " << GetProcessName() << "]" << G4endl;
1024  G4cout << " for " << track.GetDefinition()->GetParticleName()
1025  << " in Material " << currentMaterial->GetName()
1026  << " Ekin(MeV)= " << preStepKinEnergy/MeV
1027  <<G4endl;
1028  G4cout << "MeanFreePath = " << currentInteractionLength/cm << "[cm]"
1029  << "InteractionLength= " << x/cm <<"[cm] " <<G4endl;
1030  }
1031 #endif
1032  }
1033  return x;
1034 }
1035 
1036 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1037 
1039  const G4Step& step)
1040 {
1042  // The process has range table - calculate energy loss
1043  if(!isIonisation || !currentModel->IsActive(preStepScaledEnergy)) {
1044  return &fParticleChange;
1045  }
1046 
1047  // Get the actual (true) Step length
1048  G4double length = step.GetStepLength();
1049  if(length <= 0.0) { return &fParticleChange; }
1050  G4double eloss = 0.0;
1051 
1052  /*
1053  if(1 < verboseLevel) {
1054  const G4ParticleDefinition* d = track.GetParticleDefinition();
1055  G4cout << "AlongStepDoIt for "
1056  << GetProcessName() << " and particle "
1057  << d->GetParticleName()
1058  << " eScaled(MeV)= " << preStepScaledEnergy/MeV
1059  << " range(mm)= " << fRange/mm
1060  << " s(mm)= " << length/mm
1061  << " rf= " << reduceFactor
1062  << " q^2= " << chargeSqRatio
1063  << " md= " << d->GetPDGMass()
1064  << " status= " << track.GetTrackStatus()
1065  << G4endl;
1066  }
1067  */
1068 
1069  const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
1070 
1071  // define new weight for primary and secondaries
1073  if(weightFlag) {
1074  weight /= biasFactor;
1076  }
1077 
1078  // stopping
1079  if (length >= fRange) {
1080  eloss = preStepKinEnergy;
1081  if (useDeexcitation) {
1082  atomDeexcitation->AlongStepDeexcitation(scTracks, step,
1083  eloss, currentCoupleIndex);
1084  if(scTracks.size() > 0) { FillSecondariesAlongStep(eloss, weight); }
1085  if(eloss < 0.0) { eloss = 0.0; }
1086  }
1089  return &fParticleChange;
1090  }
1091 
1092  // Short step
1093  eloss = GetDEDXForScaledEnergy(preStepScaledEnergy)*length;
1094 
1095  // Long step
1096  if(eloss > preStepKinEnergy*linLossLimit) {
1097 
1098  //G4double x = GetScaledRangeForScaledEnergy(preStepScaledEnergy)
1099  // - length/reduceFactor;
1100  G4double x = (fRange - length)/reduceFactor;
1101  eloss = preStepKinEnergy - ScaledKinEnergyForLoss(x)/massRatio;
1102 
1103  /*
1104  if(-1 < verboseLevel)
1105  G4cout << "Long STEP: rPre(mm)= "
1106  << GetScaledRangeForScaledEnergy(preStepScaledEnergy)/mm
1107  << " rPost(mm)= " << x/mm
1108  << " ePre(MeV)= " << preStepScaledEnergy/MeV
1109  << " eloss(MeV)= " << eloss/MeV
1110  << " eloss0(MeV)= "
1111  << GetDEDXForScaledEnergy(preStepScaledEnergy)*length/MeV
1112  << " lim(MeV)= " << preStepKinEnergy*linLossLimit/MeV
1113  << G4endl;
1114  */
1115  }
1116 
1117  /*
1118  G4double eloss0 = eloss;
1119  if(-1 < verboseLevel ) {
1120  G4cout << "Before fluct: eloss(MeV)= " << eloss/MeV
1121  << " e-eloss= " << preStepKinEnergy-eloss
1122  << " step(mm)= " << length/mm
1123  << " range(mm)= " << fRange/mm
1124  << " fluct= " << lossFluctuationFlag
1125  << G4endl;
1126  }
1127  */
1128 
1129  G4double cut = (*theCuts)[currentCoupleIndex];
1130  G4double esec = 0.0;
1131 
1132  // SubCutOff
1133  if(useSubCutoff) {
1134  if(idxSCoffRegions[currentCoupleIndex]) {
1135 
1136  G4bool yes = false;
1137  G4StepPoint* prePoint = step.GetPreStepPoint();
1138 
1139  // Check boundary
1140  if(prePoint->GetStepStatus() == fGeomBoundary) { yes = true; }
1141 
1142  // Check PrePoint
1143  else {
1144  G4double preSafety = prePoint->GetSafety();
1145  G4double rcut = currentCouple->GetProductionCuts()->GetProductionCut(1);
1146 
1147  // recompute presafety
1148  if(preSafety < rcut) {
1149  preSafety = safetyHelper->ComputeSafety(prePoint->GetPosition());
1150  }
1151 
1152  if(preSafety < rcut) { yes = true; }
1153 
1154  // Check PostPoint
1155  else {
1156  G4double postSafety = preSafety - length;
1157  if(postSafety < rcut) {
1158  postSafety =
1159  safetyHelper->ComputeSafety(step.GetPostStepPoint()->GetPosition());
1160  if(postSafety < rcut) { yes = true; }
1161  }
1162  }
1163  }
1164 
1165  // Decided to start subcut sampling
1166  if(yes) {
1167 
1168  cut = (*theSubCuts)[currentCoupleIndex];
1169  eloss -= GetSubDEDXForScaledEnergy(preStepScaledEnergy)*length;
1170  esec = SampleSubCutSecondaries(scTracks, step,
1171  currentModel,currentCoupleIndex);
1172  // add bremsstrahlung sampling
1173  /*
1174  if(nProcesses > 0) {
1175  for(G4int i=0; i<nProcesses; ++i) {
1176  (scProcesses[i])->SampleSubCutSecondaries(
1177  scTracks, step, (scProcesses[i])->
1178  SelectModelForMaterial(preStepKinEnergy, currentCoupleIndex),
1179  currentCoupleIndex);
1180  }
1181  }
1182  */
1183  }
1184  }
1185  }
1186 
1187  // Corrections, which cannot be tabulated
1188  if(isIon) {
1189  G4double eadd = 0.0;
1190  G4double eloss_before = eloss;
1191  currentModel->CorrectionsAlongStep(currentCouple, dynParticle,
1192  eloss, eadd, length);
1193  if(eloss < 0.0) { eloss = 0.5*eloss_before; }
1194  }
1195 
1196  // Sample fluctuations
1197  if (lossFluctuationFlag) {
1198  G4VEmFluctuationModel* fluc = currentModel->GetModelOfFluctuations();
1199  if(fluc &&
1200  (eloss + esec + lowestKinEnergy) < preStepKinEnergy) {
1201 
1202  G4double tmax =
1203  std::min(currentModel->MaxSecondaryKinEnergy(dynParticle),cut);
1204  G4double emean = eloss;
1205  eloss = fluc->SampleFluctuations(currentMaterial,dynParticle,
1206  tmax,length,emean);
1207  /*
1208  if(-1 < verboseLevel)
1209  G4cout << "After fluct: eloss(MeV)= " << eloss/MeV
1210  << " fluc= " << (eloss-eloss0)/MeV
1211  << " ChargeSqRatio= " << chargeSqRatio
1212  << " massRatio= " << massRatio
1213  << " tmax= " << tmax
1214  << G4endl;
1215  */
1216  }
1217  }
1218 
1219  // deexcitation
1220  if (useDeexcitation) {
1221  G4double esecfluo = preStepKinEnergy - esec;
1222  G4double de = esecfluo;
1223  //G4double eloss0 = eloss;
1224  /*
1225  G4cout << "### 1: E(keV)= " << preStepKinEnergy/keV
1226  << " Efluomax(keV)= " << de/keV
1227  << " Eloss(keV)= " << eloss/keV << G4endl;
1228  */
1229  atomDeexcitation->AlongStepDeexcitation(scTracks, step,
1230  de, currentCoupleIndex);
1231 
1232  // sum of de-excitation energies
1233  esecfluo -= de;
1234 
1235  // subtracted from energy loss
1236  if(eloss >= esecfluo) {
1237  esec += esecfluo;
1238  eloss -= esecfluo;
1239  } else {
1240  esec += esecfluo;
1241  eloss = 0.0;
1242  }
1243  /*
1244  if(esecfluo > 0.0) {
1245  G4cout << "### 2: E(keV)= " << preStepKinEnergy/keV
1246  << " Esec(keV)= " << esec/keV
1247  << " Esecf(kV)= " << esecfluo/keV
1248  << " Eloss0(kV)= " << eloss0/keV
1249  << " Eloss(keV)= " << eloss/keV
1250  << G4endl;
1251  }
1252  */
1253  }
1254  if(scTracks.size() > 0) { FillSecondariesAlongStep(eloss, weight); }
1255 
1256  // Energy balanse
1257  G4double finalT = preStepKinEnergy - eloss - esec;
1258  if (finalT <= lowestKinEnergy) {
1259  eloss += finalT;
1260  finalT = 0.0;
1261  } else if(isIon) {
1263  currentModel->GetParticleCharge(track.GetParticleDefinition(),
1264  currentMaterial,finalT));
1265  }
1266 
1267  if(eloss < 0.0) { eloss = 0.0; }
1270 
1271  if(1 < verboseLevel) {
1272  G4double del = finalT + eloss + esec - preStepKinEnergy;
1273  G4cout << "Final value eloss(MeV)= " << eloss/MeV
1274  << " preStepKinEnergy= " << preStepKinEnergy
1275  << " postStepKinEnergy= " << finalT
1276  << " de(keV)= " << del/keV
1277  << " lossFlag= " << lossFluctuationFlag
1278  << " status= " << track.GetTrackStatus()
1279  << G4endl;
1280  }
1281 
1282  return &fParticleChange;
1283 }
1284 
1285 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1286 
1287 void
1288 G4VEnergyLossProcess::FillSecondariesAlongStep(G4double&, G4double& weight)
1289 {
1290  // weight may be changed by biasing manager
1291  if(biasManager) {
1292  if(biasManager->SecondaryBiasingRegion(currentCoupleIndex)) {
1293  weight *=
1294  biasManager->ApplySecondaryBiasing(scTracks, currentCoupleIndex);
1295  }
1296  }
1297 
1298  // fill secondaries
1299  G4int n = scTracks.size();
1301 
1302  for(G4int i=0; i<n; ++i) {
1303  G4Track* t = scTracks[i];
1304  if(t) {
1305  //esec += t->GetKineticEnergy();
1306  t->SetWeight(weight);
1308  //G4cout << "Secondary(along step) has weight " << t->GetWeight()
1309  //<< ", kenergy " << t->GetKineticEnergy()/MeV << " MeV" <<G4endl;
1310  }
1311  }
1312  scTracks.clear();
1313 }
1314 
1315 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1316 
1317 G4double
1318 G4VEnergyLossProcess::SampleSubCutSecondaries(std::vector<G4Track*>& tracks,
1319  const G4Step& step,
1320  G4VEmModel* model,
1321  G4int idx)
1322 {
1323  // Fast check weather subcutoff can work
1324  G4double esec = 0.0;
1325  G4double subcut = (*theSubCuts)[idx];
1326  G4double cut = (*theCuts)[idx];
1327  if(cut <= subcut) { return esec; }
1328 
1329  const G4Track* track = step.GetTrack();
1330  const G4DynamicParticle* dp = track->GetDynamicParticle();
1331  G4double e = dp->GetKineticEnergy()*massRatio;
1332  G4double cross = (*theDensityFactor)[idx]*chargeSqRatio
1333  *(((*theSubLambdaTable)[(*theDensityIdx)[idx]])->Value(e));
1334  G4double length = step.GetStepLength();
1335 
1336  // negligible probability to get any interaction
1337  if(length*cross < perMillion) { return esec; }
1338  /*
1339  if(-1 < verboseLevel)
1340  G4cout << "<<< Subcutoff for " << GetProcessName()
1341  << " cross(1/mm)= " << cross*mm << ">>>"
1342  << " e(MeV)= " << preStepScaledEnergy
1343  << " matIdx= " << currentCoupleIndex
1344  << G4endl;
1345  */
1346 
1347  // Sample subcutoff secondaries
1348  G4StepPoint* preStepPoint = step.GetPreStepPoint();
1349  G4StepPoint* postStepPoint = step.GetPostStepPoint();
1350  G4ThreeVector prepoint = preStepPoint->GetPosition();
1351  G4ThreeVector dr = postStepPoint->GetPosition() - prepoint;
1352  G4double pretime = preStepPoint->GetGlobalTime();
1353  G4double dt = postStepPoint->GetGlobalTime() - pretime;
1354  //G4double dt = length/preStepPoint->GetVelocity();
1355  G4double fragment = 0.0;
1356 
1357  do {
1358  G4double del = -std::log(G4UniformRand())/cross;
1359  fragment += del/length;
1360  if (fragment > 1.0) break;
1361 
1362  // sample secondaries
1363  secParticles.clear();
1364  model->SampleSecondaries(&secParticles,track->GetMaterialCutsCouple(),
1365  dp,subcut,cut);
1366 
1367  // position of subcutoff particles
1368  G4ThreeVector r = prepoint + fragment*dr;
1369  std::vector<G4DynamicParticle*>::iterator it;
1370  for(it=secParticles.begin(); it!=secParticles.end(); ++it) {
1371 
1372  G4bool addSec = true;
1373  /*
1374  // do not track very low-energy delta-electrons
1375  if(theSecondaryRangeTable && (*it)->GetParticleDefinition() == theElectron) {
1376  G4double ekin = (*it)->GetKineticEnergy();
1377  G4double rg = ((*theSecondaryRangeTable)[idx]->Value(ekin));
1378  // if(rg < currentMinSafety) {
1379  if(rg < safetyHelper->ComputeSafety(r)) {
1380  extraEdep += ekin;
1381  delete (*it);
1382  addSec = false;
1383  }
1384  }
1385  */
1386  if(addSec) {
1387  G4Track* t = new G4Track((*it), pretime + fragment*dt, r);
1389  tracks.push_back(t);
1390  esec += t->GetKineticEnergy();
1391  if (t->GetParticleDefinition() == thePositron) {
1392  esec += 2.0*electron_mass_c2;
1393  }
1394 
1395  /*
1396  if(-1 < verboseLevel)
1397  G4cout << "New track " << t->GetParticleDefinition()->GetParticleName()
1398  << " e(keV)= " << t->GetKineticEnergy()/keV
1399  << " fragment= " << fragment
1400  << G4endl;
1401  */
1402  }
1403  }
1404  } while (fragment <= 1.0);
1405  return esec;
1406 }
1407 
1408 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1409 
1411  const G4Step& step)
1412 {
1413  // In all cases clear number of interaction lengths
1415  mfpKinEnergy = DBL_MAX;
1416 
1418  G4double finalT = track.GetKineticEnergy();
1419  if(finalT <= lowestKinEnergy) { return &fParticleChange; }
1420 
1421  G4double postStepScaledEnergy = finalT*massRatio;
1422 
1423  if(!currentModel->IsActive(postStepScaledEnergy)) { return &fParticleChange; }
1424  /*
1425  if(-1 < verboseLevel) {
1426  G4cout << GetProcessName()
1427  << "::PostStepDoIt: E(MeV)= " << finalT/MeV
1428  << G4endl;
1429  }
1430  */
1431 
1432  // forced process - should happen only once per track
1433  if(biasFlag) {
1434  if(biasManager->ForcedInteractionRegion(currentCoupleIndex)) {
1435  biasFlag = false;
1436  }
1437  }
1438 
1439  // Integral approach
1440  if (integral) {
1441  G4double lx = GetLambdaForScaledEnergy(postStepScaledEnergy);
1442  /*
1443  if(preStepLambda<lx && 1 < verboseLevel && nWarnings<200) {
1444  G4cout << "WARNING: for " << particle->GetParticleName()
1445  << " and " << GetProcessName()
1446  << " E(MeV)= " << finalT/MeV
1447  << " preLambda= " << preStepLambda
1448  << " < " << lx << " (postLambda) "
1449  << G4endl;
1450  ++nWarnings;
1451  }
1452  */
1453  if(lx <= 0.0) {
1454  return &fParticleChange;
1455  } else if(preStepLambda*G4UniformRand() > lx) {
1456  return &fParticleChange;
1457  }
1458  }
1459 
1460  SelectModel(postStepScaledEnergy);
1461 
1462  // define new weight for primary and secondaries
1464  if(weightFlag) {
1465  weight /= biasFactor;
1467  }
1468 
1469  const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
1470  G4double tcut = (*theCuts)[currentCoupleIndex];
1471 
1472  // sample secondaries
1473  secParticles.clear();
1474  //G4cout << "Energy of primary: " << dynParticle->GetKineticEnergy()/MeV<<G4endl;
1475  currentModel->SampleSecondaries(&secParticles, currentCouple, dynParticle, tcut);
1476 
1477  // bremsstrahlung splitting or Russian roulette
1478  if(biasManager) {
1479  if(biasManager->SecondaryBiasingRegion(currentCoupleIndex)) {
1480  G4double eloss = 0.0;
1481  weight *= biasManager->ApplySecondaryBiasing(secParticles,
1482  track, currentModel,
1483  &fParticleChange, eloss,
1484  currentCoupleIndex, tcut,
1485  step.GetPostStepPoint()->GetSafety());
1486  if(eloss > 0.0) {
1489  }
1490  }
1491  }
1492 
1493  // save secondaries
1494  G4int num = secParticles.size();
1495  if(num > 0) {
1496 
1498 
1499  for (G4int i=0; i<num; ++i) {
1500  if(secParticles[i]) {
1501  G4Track* t = new G4Track(secParticles[i], track.GetGlobalTime(),
1502  track.GetPosition());
1504  t->SetWeight(weight);
1505  //G4cout << "Secondary(post step) has weight " << t->GetWeight()
1506  //<< ", kenergy " << t->GetKineticEnergy()/MeV << " MeV" <<G4endl;
1508  }
1509  }
1510  }
1511 
1514  if(particle->GetProcessManager()->GetAtRestProcessVector()->size() > 0)
1517  }
1518 
1519  /*
1520  if(-1 < verboseLevel) {
1521  G4cout << "::PostStepDoIt: Sample secondary; Efin= "
1522  << fParticleChange.GetProposedKineticEnergy()/MeV
1523  << " MeV; model= (" << currentModel->LowEnergyLimit()
1524  << ", " << currentModel->HighEnergyLimit() << ")"
1525  << " preStepLambda= " << preStepLambda
1526  << " dir= " << track.GetMomentumDirection()
1527  << " status= " << track.GetTrackStatus()
1528  << G4endl;
1529  }
1530  */
1531  return &fParticleChange;
1532 }
1533 
1534 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1535 
1537  const G4ParticleDefinition* part, const G4String& directory,
1538  G4bool ascii)
1539 {
1540  G4bool res = true;
1541  if ( baseParticle || part != particle ) return res;
1542 
1543  if(!StoreTable(part,theDEDXTable,ascii,directory,"DEDX"))
1544  {res = false;}
1545 
1546  if(!StoreTable(part,theDEDXunRestrictedTable,ascii,directory,"DEDXnr"))
1547  {res = false;}
1548 
1549  if(!StoreTable(part,theDEDXSubTable,ascii,directory,"SubDEDX"))
1550  {res = false;}
1551 
1552  if(!StoreTable(part,theIonisationTable,ascii,directory,"Ionisation"))
1553  {res = false;}
1554 
1555  if(!StoreTable(part,theIonisationSubTable,ascii,directory,"SubIonisation"))
1556  {res = false;}
1557 
1558  if(isIonisation &&
1559  !StoreTable(part,theCSDARangeTable,ascii,directory,"CSDARange"))
1560  {res = false;}
1561 
1562  if(isIonisation &&
1563  !StoreTable(part,theRangeTableForLoss,ascii,directory,"Range"))
1564  {res = false;}
1565 
1566  if(isIonisation &&
1567  !StoreTable(part,theInverseRangeTable,ascii,directory,"InverseRange"))
1568  {res = false;}
1569 
1570  if(!StoreTable(part,theLambdaTable,ascii,directory,"Lambda"))
1571  {res = false;}
1572 
1573  if(!StoreTable(part,theSubLambdaTable,ascii,directory,"SubLambda"))
1574  {res = false;}
1575 
1576  if ( res ) {
1577  if(0 < verboseLevel) {
1578  G4cout << "Physics tables are stored for " << particle->GetParticleName()
1579  << " and process " << GetProcessName()
1580  << " in the directory <" << directory
1581  << "> " << G4endl;
1582  }
1583  } else {
1584  G4cout << "Fail to store Physics Tables for "
1585  << particle->GetParticleName()
1586  << " and process " << GetProcessName()
1587  << " in the directory <" << directory
1588  << "> " << G4endl;
1589  }
1590  return res;
1591 }
1592 
1593 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1594 
1595 G4bool
1597  const G4String& directory,
1598  G4bool ascii)
1599 {
1600  G4bool res = true;
1601  const G4String particleName = part->GetParticleName();
1602 
1603  if(1 < verboseLevel) {
1604  G4cout << "G4VEnergyLossProcess::RetrievePhysicsTable() for "
1605  << particleName << " and process " << GetProcessName()
1606  << "; tables_are_built= " << tablesAreBuilt
1607  << G4endl;
1608  }
1609  if(particle == part) {
1610 
1611  if ( !baseParticle ) {
1612 
1613  G4bool fpi = true;
1614  if(!RetrieveTable(part,theDEDXTable,ascii,directory,"DEDX",fpi))
1615  {fpi = false;}
1616 
1617  // ionisation table keeps individual dEdx and not sum of sub-processes
1618  if(!RetrieveTable(part,theDEDXTable,ascii,directory,"Ionisation",false))
1619  {fpi = false;}
1620 
1621  if(!RetrieveTable(part,theRangeTableForLoss,ascii,directory,"Range",fpi))
1622  {res = false;}
1623 
1624  if(!RetrieveTable(part,theDEDXunRestrictedTable,ascii,directory,"DEDXnr",false))
1625  {res = false;}
1626 
1627  if(!RetrieveTable(part,theCSDARangeTable,ascii,directory,"CSDARange",false))
1628  {res = false;}
1629 
1630  if(!RetrieveTable(part,theInverseRangeTable,ascii,directory,"InverseRange",fpi))
1631  {res = false;}
1632 
1633  if(!RetrieveTable(part,theLambdaTable,ascii,directory,"Lambda",true))
1634  {res = false;}
1635 
1636  G4bool yes = false;
1637  if(nSCoffRegions > 0) {yes = true;}
1638 
1639  if(!RetrieveTable(part,theDEDXSubTable,ascii,directory,"SubDEDX",yes))
1640  {res = false;}
1641 
1642  if(!RetrieveTable(part,theSubLambdaTable,ascii,directory,"SubLambda",yes))
1643  {res = false;}
1644 
1645  if(!fpi) yes = false;
1646  if(!RetrieveTable(part,theIonisationSubTable,ascii,directory,"SubIonisation",yes))
1647  {res = false;}
1648  }
1649  }
1650 
1651  return res;
1652 }
1653 
1654 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1655 
1656 G4bool G4VEnergyLossProcess::StoreTable(const G4ParticleDefinition* part,
1657  G4PhysicsTable* aTable, G4bool ascii,
1658  const G4String& directory,
1659  const G4String& tname)
1660 {
1661  G4bool res = true;
1662  if ( aTable ) {
1663  const G4String name = GetPhysicsTableFileName(part,directory,tname,ascii);
1664  if( !aTable->StorePhysicsTable(name,ascii)) res = false;
1665  }
1666  return res;
1667 }
1668 
1669 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1670 
1671 G4bool
1672 G4VEnergyLossProcess::RetrieveTable(const G4ParticleDefinition* part,
1673  G4PhysicsTable* aTable,
1674  G4bool ascii,
1675  const G4String& directory,
1676  const G4String& tname,
1677  G4bool mandatory)
1678 {
1679  G4bool isRetrieved = false;
1680  G4String filename = GetPhysicsTableFileName(part,directory,tname,ascii);
1681  if(aTable) {
1682  if(aTable->ExistPhysicsTable(filename)) {
1683  if(G4PhysicsTableHelper::RetrievePhysicsTable(aTable,filename,ascii)) {
1684  isRetrieved = true;
1685  if((G4LossTableManager::Instance())->SplineFlag()) {
1686  size_t n = aTable->length();
1687  for(size_t i=0; i<n; ++i) {
1688  if((*aTable)[i]) { (*aTable)[i]->SetSpline(true); }
1689  }
1690  }
1691  if (0 < verboseLevel) {
1692  G4cout << tname << " table for " << part->GetParticleName()
1693  << " is Retrieved from <" << filename << ">"
1694  << G4endl;
1695  }
1696  }
1697  }
1698  }
1699  if(mandatory && !isRetrieved) {
1700  if(0 < verboseLevel) {
1701  G4cout << tname << " table for " << part->GetParticleName()
1702  << " from file <"
1703  << filename << "> is not Retrieved"
1704  << G4endl;
1705  }
1706  return false;
1707  }
1708  return true;
1709 }
1710 
1711 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1712 
1714  const G4MaterialCutsCouple *couple,
1715  const G4DynamicParticle* dp,
1716  G4double length)
1717 {
1718  DefineMaterial(couple);
1719  G4double ekin = dp->GetKineticEnergy();
1720  SelectModel(ekin*massRatio);
1721  G4double tmax = currentModel->MaxSecondaryKinEnergy(dp);
1722  tmax = std::min(tmax,(*theCuts)[currentCoupleIndex]);
1723  G4double d = 0.0;
1724  G4VEmFluctuationModel* fm = currentModel->GetModelOfFluctuations();
1725  if(fm) { d = fm->Dispersion(currentMaterial,dp,tmax,length); }
1726  return d;
1727 }
1728 
1729 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1730 
1732  G4double kineticEnergy, const G4MaterialCutsCouple* couple)
1733 {
1734  // Cross section per volume is calculated
1735  DefineMaterial(couple);
1736  G4double cross = 0.0;
1737  if(theLambdaTable) {
1738  cross = (*theDensityFactor)[currentCoupleIndex]*
1739  ((*theLambdaTable)[basedCoupleIndex])->Value(kineticEnergy);
1740  } else {
1741  SelectModel(kineticEnergy);
1742  cross = currentModel->CrossSectionPerVolume(currentMaterial,
1743  particle, kineticEnergy,
1744  (*theCuts)[currentCoupleIndex]);
1745  }
1746  if(cross < 0.0) { cross = 0.0; }
1747  return cross;
1748 }
1749 
1750 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1751 
1753 {
1754  DefineMaterial(track.GetMaterialCutsCouple());
1755  preStepLambda = GetLambdaForScaledEnergy(track.GetKineticEnergy()*massRatio);
1756  G4double x = DBL_MAX;
1757  if(0.0 < preStepLambda) { x = 1.0/preStepLambda; }
1758  return x;
1759 }
1760 
1761 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1762 
1764  G4double x, G4double y,
1765  G4double& z)
1766 {
1767  G4GPILSelection sel;
1768  return AlongStepGetPhysicalInteractionLength(track, x, y, z, &sel);
1769 }
1770 
1771 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1772 
1774  const G4Track& track,
1775  G4double,
1777 
1778 {
1779  *condition = NotForced;
1780  return MeanFreePath(track);
1781 }
1782 
1783 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1784 
1786  const G4Track&,
1788 {
1789  return DBL_MAX;
1790 }
1791 
1792 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1793 
1796  G4double /*cut*/)
1797 {
1798  G4PhysicsVector* v = new G4PhysicsLogVector(minKinEnergy, maxKinEnergy, nBins);
1799  v->SetSpline((G4LossTableManager::Instance())->SplineFlag());
1800  return v;
1801 }
1802 
1803 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1804 
1807 {
1808  G4bool add = true;
1809  if(p->GetProcessName() != "eBrem") { add = false; }
1810  if(add && nProcesses > 0) {
1811  for(G4int i=0; i<nProcesses; ++i) {
1812  if(p == scProcesses[i]) {
1813  add = false;
1814  break;
1815  }
1816  }
1817  }
1818  if(add) {
1819  scProcesses.push_back(p);
1820  ++nProcesses;
1821  if (1 < verboseLevel) {
1822  G4cout << "### The process " << p->GetProcessName()
1823  << " is added to the list of collaborative processes of "
1824  << GetProcessName() << G4endl;
1825  }
1826  }
1827 }
1828 
1829 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1830 
1832 {
1833  if(fTotal == tType && theDEDXunRestrictedTable != p && !baseParticle) {
1834  if(theDEDXunRestrictedTable) {
1835  theDEDXunRestrictedTable->clearAndDestroy();
1836  delete theDEDXunRestrictedTable;
1837  }
1838  theDEDXunRestrictedTable = p;
1839  if(p) {
1840  size_t n = p->length();
1841  G4PhysicsVector* pv = (*p)[0];
1842  G4double emax = maxKinEnergyCSDA;
1843 
1844  for (size_t i=0; i<n; ++i) {
1845  G4double dedx = 0.0;
1846  pv = (*p)[i];
1847  if(pv) { dedx = pv->Value(emax); }
1848  else {
1849  pv = (*p)[(*theDensityIdx)[i]];
1850  if(pv) { dedx = pv->Value(emax)*(*theDensityFactor)[i]; }
1851  }
1852  theDEDXAtMaxEnergy[i] = dedx;
1853  //G4cout << "i= " << i << " emax(MeV)= " << emax/MeV<< " dedx= "
1854  //<< dedx << G4endl;
1855  }
1856  }
1857 
1858  } else if(fRestricted == tType && theDEDXTable != p) {
1859  //G4cout << "G4VEnergyLossProcess::SetDEDXTable " << particle->GetParticleName()
1860  // << " old table " << theDEDXTable << " new table " << p
1861  // << " ion " << theIonisationTable << " bp " << baseParticle << G4endl;
1862  if(theDEDXTable && !baseParticle) {
1863  if(theDEDXTable == theIonisationTable) { theIonisationTable = 0; }
1864  theDEDXTable->clearAndDestroy();
1865  delete theDEDXTable;
1866  }
1867  theDEDXTable = p;
1868  } else if(fSubRestricted == tType && theDEDXSubTable != p) {
1869  if(theDEDXSubTable && !baseParticle) {
1870  if(theDEDXSubTable == theIonisationSubTable) { theIonisationSubTable = 0; }
1871  theDEDXSubTable->clearAndDestroy();
1872  delete theDEDXSubTable;
1873  }
1874  theDEDXSubTable = p;
1875  } else if(fIsIonisation == tType && theIonisationTable != p) {
1876  if(theIonisationTable && theIonisationTable != theDEDXTable && !baseParticle) {
1877  theIonisationTable->clearAndDestroy();
1878  delete theIonisationTable;
1879  }
1880  theIonisationTable = p;
1881  } else if(fIsSubIonisation == tType && theIonisationSubTable != p) {
1882  if(theIonisationSubTable && theIonisationSubTable != theDEDXSubTable && !baseParticle) {
1883  theIonisationSubTable->clearAndDestroy();
1884  delete theIonisationSubTable;
1885  }
1886  theIonisationSubTable = p;
1887  }
1888 }
1889 
1890 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1891 
1893 {
1894  if(theCSDARangeTable != p) { theCSDARangeTable = p; }
1895 
1896  if(p) {
1897  size_t n = p->length();
1898  G4PhysicsVector* pv;
1899  G4double emax = maxKinEnergyCSDA;
1900 
1901  for (size_t i=0; i<n; ++i) {
1902  pv = (*p)[i];
1903  G4double rmax = 0.0;
1904  if(pv) { rmax = pv->Value(emax); }
1905  else {
1906  pv = (*p)[(*theDensityIdx)[i]];
1907  if(pv) { rmax = pv->Value(emax)/(*theDensityFactor)[i]; }
1908  }
1909  theRangeAtMaxEnergy[i] = rmax;
1910  //G4cout << "i= " << i << " Emax(MeV)= " << emax/MeV << " Rmax= "
1911  //<< rmax<< G4endl;
1912  }
1913  }
1914 }
1915 
1916 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1917 
1919 {
1920  if(theRangeTableForLoss != p) {
1921  theRangeTableForLoss = p;
1922  if(1 < verboseLevel) {
1923  G4cout << "### Set Range table " << p
1924  << " for " << particle->GetParticleName()
1925  << " and process " << GetProcessName() << G4endl;
1926  }
1927  }
1928 }
1929 
1930 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1931 
1933 {
1934  if(theSecondaryRangeTable != p) {
1935  theSecondaryRangeTable = p;
1936  if(1 < verboseLevel) {
1937  G4cout << "### Set SecondaryRange table " << p
1938  << " for " << particle->GetParticleName()
1939  << " and process " << GetProcessName() << G4endl;
1940  }
1941  }
1942 }
1943 
1944 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1945 
1947 {
1948  if(theInverseRangeTable != p) {
1949  theInverseRangeTable = p;
1950  if(1 < verboseLevel) {
1951  G4cout << "### Set InverseRange table " << p
1952  << " for " << particle->GetParticleName()
1953  << " and process " << GetProcessName() << G4endl;
1954  }
1955  }
1956 }
1957 
1958 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1959 
1961 {
1962  if(1 < verboseLevel) {
1963  G4cout << "### Set Lambda table " << p
1964  << " for " << particle->GetParticleName()
1965  << " and process " << GetProcessName() << G4endl;
1966  }
1967  if(theLambdaTable != p) { theLambdaTable = p; }
1968  tablesAreBuilt = true;
1969 
1970  if(theLambdaTable) {
1971  size_t n = theLambdaTable->length();
1972  G4PhysicsVector* pv = (*theLambdaTable)[0];
1973  G4double e, ss, smax, emax;
1974 
1975  size_t i;
1976 
1977  // first loop on existing vectors
1978  for (i=0; i<n; ++i) {
1979  pv = (*theLambdaTable)[i];
1980  if(pv) {
1981  size_t nb = pv->GetVectorLength();
1982  emax = DBL_MAX;
1983  smax = 0.0;
1984  if(nb > 0) {
1985  for (size_t j=0; j<nb; ++j) {
1986  e = pv->Energy(j);
1987  ss = (*pv)(j);
1988  if(ss > smax) {
1989  smax = ss;
1990  emax = e;
1991  }
1992  }
1993  }
1994  theEnergyOfCrossSectionMax[i] = emax;
1995  theCrossSectionMax[i] = smax;
1996  if(1 < verboseLevel) {
1997  G4cout << "For " << particle->GetParticleName()
1998  << " Max CS at i= " << i << " emax(MeV)= " << emax/MeV
1999  << " lambda= " << smax << G4endl;
2000  }
2001  }
2002  }
2003  // second loop using base materials
2004  for (i=0; i<n; ++i) {
2005  pv = (*theLambdaTable)[i];
2006  if(!pv){
2007  G4int j = (*theDensityIdx)[i];
2008  theEnergyOfCrossSectionMax[i] = theEnergyOfCrossSectionMax[j];
2009  theCrossSectionMax[i] = (*theDensityFactor)[i]*theCrossSectionMax[j];
2010  }
2011  }
2012  }
2013 }
2014 
2015 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2016 
2018 {
2019  if(theSubLambdaTable != p) {
2020  theSubLambdaTable = p;
2021  if(1 < verboseLevel) {
2022  G4cout << "### Set SebLambda table " << p
2023  << " for " << particle->GetParticleName()
2024  << " and process " << GetProcessName() << G4endl;
2025  }
2026  }
2027 }
2028 
2029 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2030 
2032 {
2033  const G4Element* elm = 0;
2034  if(currentModel) { elm = currentModel->GetCurrentElement(); }
2035  return elm;
2036 }
2037 
2038 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2039 
2041  G4bool flag)
2042 {
2043  if(f > 0.0) {
2044  biasFactor = f;
2045  weightFlag = flag;
2046  if(1 < verboseLevel) {
2047  G4cout << "### SetCrossSectionBiasingFactor: for "
2048  << " process " << GetProcessName()
2049  << " biasFactor= " << f << " weightFlag= " << flag
2050  << G4endl;
2051  }
2052  }
2053 }
2054 
2055 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2056 
2057 void
2059  const G4String& region,
2060  G4bool flag)
2061 {
2062  if(!biasManager) { biasManager = new G4EmBiasingManager(); }
2063  if(1 < verboseLevel) {
2064  G4cout << "### ActivateForcedInteraction: for "
2065  << " process " << GetProcessName()
2066  << " length(mm)= " << length/mm
2067  << " in G4Region <" << region
2068  << "> weightFlag= " << flag
2069  << G4endl;
2070  }
2071  weightFlag = flag;
2072  biasManager->ActivateForcedInteraction(length, region);
2073 }
2074 
2075 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2076 
2077 void
2079  G4double factor,
2080  G4double energyLimit)
2081 {
2082  if (0.0 <= factor) {
2083 
2084  // Range cut can be applied only for e-
2085  if(0.0 == factor && secondaryParticle != G4Electron::Electron())
2086  { return; }
2087 
2088  if(!biasManager) { biasManager = new G4EmBiasingManager(); }
2089  biasManager->ActivateSecondaryBiasing(region, factor, energyLimit);
2090  if(1 < verboseLevel) {
2091  G4cout << "### ActivateSecondaryBiasing: for "
2092  << " process " << GetProcessName()
2093  << " factor= " << factor
2094  << " in G4Region <" << region
2095  << "> energyLimit(MeV)= " << energyLimit/MeV
2096  << G4endl;
2097  }
2098  }
2099 }
2100 
2101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2102