Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4LossTableManager.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: G4LossTableManager.cc 79268 2014-02-20 16:46:31Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4LossTableManager
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 03.01.2002
38 //
39 // Modifications:
40 //
41 // 20-01-03 Migrade to cut per region (V.Ivanchenko)
42 // 15-02-03 Lambda table can be scaled (V.Ivanchenko)
43 // 17-02-03 Fix problem of store/restore tables (V.Ivanchenko)
44 // 10-03-03 Add Ion registration (V.Ivanchenko)
45 // 25-03-03 Add deregistration (V.Ivanchenko)
46 // 02-04-03 Change messenger (V.Ivanchenko)
47 // 26-04-03 Fix retrieve tables (V.Ivanchenko)
48 // 13-05-03 Add calculation of precise range (V.Ivanchenko)
49 // 23-07-03 Add exchange with G4EnergyLossTables (V.Ivanchenko)
50 // 05-10-03 Add G4VEmProcesses registration and Verbose command (V.Ivanchenko)
51 // 17-10-03 Add SetParameters method (V.Ivanchenko)
52 // 23-10-03 Add control on inactive processes (V.Ivanchenko)
53 // 04-11-03 Add checks in RetrievePhysicsTable (V.Ivanchenko)
54 // 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
55 // 14-01-04 Activate precise range calculation (V.Ivanchenko)
56 // 10-03-04 Fix a problem of Precise Range table (V.Ivanchenko)
57 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivanchenko)
58 // 13-01-04 Fix problem which takes place for inactivate eIoni (V.Ivanchenko)
59 // 25-01-04 Fix initialisation problem for ions (V.Ivanchenko)
60 // 11-03-05 Shift verbose level by 1 (V.Ivantchenko)
61 // 10-01-06 PreciseRange -> CSDARange (V.Ivantchenko)
62 // 20-01-06 Introduce G4EmTableType to remove repeating code (VI)
63 // 23-03-06 Set flag isIonisation (VI)
64 // 10-05-06 Add methods SetMscStepLimitation, FacRange and MscFlag (VI)
65 // 22-05-06 Add methods Set/Get bremsTh (VI)
66 // 05-06-06 Do not clear loss_table map between runs (VI)
67 // 16-01-07 Create new energy loss table for e+,e-,mu+,mu- and
68 // left ionisation table for further usage (VI)
69 // 12-02-07 Add SetSkin, SetLinearLossLimit (V.Ivanchenko)
70 // 18-06-07 Move definition of msc parameters to G4EmProcessOptions (V.Ivanchenko)
71 // 21-02-08 Added G4EmSaturation (V.Ivanchenko)
72 // 12-04-10 Added PreparePhysicsTables and BuildPhysicsTables entries (V.Ivanchenko)
73 // 04-06-13 (V.Ivanchenko) Adaptation for MT mode; new method LocalPhysicsTables;
74 // ions expect G4GenericIon are not included in the map of energy loss
75 // processes for performnc reasons
76 //
77 // Class Description:
78 //
79 // -------------------------------------------------------------------
80 //
81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83 
84 #include "G4LossTableManager.hh"
85 #include "G4SystemOfUnits.hh"
86 #include "G4EnergyLossMessenger.hh"
87 #include "G4PhysicsTable.hh"
88 #include "G4ParticleDefinition.hh"
89 #include "G4MaterialCutsCouple.hh"
90 #include "G4ProcessManager.hh"
91 #include "G4Electron.hh"
92 #include "G4Proton.hh"
93 #include "G4VMultipleScattering.hh"
94 #include "G4VEmProcess.hh"
95 #include "G4ProductionCutsTable.hh"
96 #include "G4PhysicsTableHelper.hh"
97 #include "G4EmCorrections.hh"
98 #include "G4EmSaturation.hh"
99 #include "G4EmConfigurator.hh"
100 #include "G4ElectronIonPair.hh"
101 #include "G4EmTableType.hh"
102 #include "G4LossTableBuilder.hh"
103 #include "G4VAtomDeexcitation.hh"
104 #include "G4Region.hh"
105 #include "G4PhysicalConstants.hh"
106 #include "G4Threading.hh"
107 
108 //G4ThreadLocal G4LossTableManager* G4LossTableManager::theInstance = 0;
109 
110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
111 
113 {
114  /*
115  if(!theInstance) {
116  theInstance = new G4LossTableManager;
117  }
118  return theInstance;
119  */
121  return instance.Instance();
122 }
123 
124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
125 
127 {
128  //G4cout << "### G4LossTableManager::~G4LossTableManager()" << G4endl;
129  for (G4int i=0; i<n_loss; ++i) {
130  if( loss_vector[i] ) { delete loss_vector[i]; }
131  }
132  size_t msc = msc_vector.size();
133  for (size_t j=0; j<msc; ++j) {
134  if( msc_vector[j] ) { delete msc_vector[j]; }
135  }
136  size_t emp = emp_vector.size();
137  for (size_t k=0; k<emp; ++k) {
138  if( emp_vector[k] ) { delete emp_vector[k]; }
139  }
140  size_t mod = mod_vector.size();
141  size_t fmod = fmod_vector.size();
142  for (size_t a=0; a<mod; ++a) {
143  if( mod_vector[a] ) {
144  for (size_t b=0; b<fmod; ++b) {
145  if((G4VEmModel*)(fmod_vector[b]) == mod_vector[a]) {
146  fmod_vector[b] = 0;
147  }
148  }
149  delete mod_vector[a];
150  }
151  }
152  for (size_t b=0; b<fmod; ++b) {
153  if( fmod_vector[b] ) { delete fmod_vector[b]; }
154  }
155  Clear();
156  delete theMessenger;
157  delete tableBuilder;
158  delete emCorrections;
159  delete emSaturation;
160  delete emConfigurator;
161  delete emElectronIonPair;
162  delete atomDeexcitation;
163 }
164 
165 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
166 
167 G4LossTableManager::G4LossTableManager()
168 {
169  n_loss = 0;
170  run = -1;
171  startInitialisation = false;
172  all_tables_are_built = false;
173  currentLoss = 0;
174  currentParticle = 0;
175  firstParticle = 0;
176  lossFluctuationFlag = true;
177  subCutoffFlag = false;
178  rndmStepFlag = false;
179  minSubRange = 0.0;
180  maxRangeVariation = 1.0;
181  maxFinalStep = 0.0;
182  minKinEnergy = 0.1*keV;
183  maxKinEnergy = 10.0*TeV;
184  nbinsLambda = 77;
185  nbinsPerDecade = 7;
186  maxKinEnergyForMuons = 10.*TeV;
187  integral = true;
188  integralActive = false;
189  buildCSDARange = false;
190  minEnergyActive = false;
191  maxEnergyActive = false;
192  maxEnergyForMuonsActive = false;
193  stepFunctionActive = false;
194  flagLPM = true;
195  splineFlag = true;
196  isMaster = true;
197  bremsTh = DBL_MAX;
198  factorForAngleLimit = 1.0;
199  verbose = 1;
200  theMessenger = new G4EnergyLossMessenger();
201  theElectron = G4Electron::Electron();
202  theGenericIon= 0;
203  tableBuilder = new G4LossTableBuilder();
204  emCorrections= new G4EmCorrections();
205  emSaturation = new G4EmSaturation();
206  emConfigurator = new G4EmConfigurator(verbose);
207  emElectronIonPair = new G4ElectronIonPair();
208  tableBuilder->SetSplineFlag(splineFlag);
209  atomDeexcitation = 0;
211  verbose = 0;
212  isMaster = false;
213  }
214 }
215 
216 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
217 
219 {
220  all_tables_are_built = false;
221  currentLoss = 0;
222  currentParticle = 0;
223  if(n_loss)
224  {
225  dedx_vector.clear();
226  range_vector.clear();
227  inv_range_vector.clear();
228  loss_map.clear();
229  loss_vector.clear();
230  part_vector.clear();
231  base_part_vector.clear();
232  tables_are_built.clear();
233  isActive.clear();
234  n_loss = 0;
235  }
236 }
237 
238 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
239 
241 {
242  if(!p) { return; }
243  for (G4int i=0; i<n_loss; ++i) {
244  if(loss_vector[i] == p) { return; }
245  }
246  if(verbose > 1) {
247  G4cout << "G4LossTableManager::Register G4VEnergyLossProcess : "
248  << p->GetProcessName() << " idx= " << n_loss << G4endl;
249  }
250  ++n_loss;
251  loss_vector.push_back(p);
252  part_vector.push_back(0);
253  base_part_vector.push_back(0);
254  dedx_vector.push_back(0);
255  range_vector.push_back(0);
256  inv_range_vector.push_back(0);
257  tables_are_built.push_back(false);
258  isActive.push_back(true);
259  all_tables_are_built = false;
260  if(!lossFluctuationFlag) { p->SetLossFluctuations(false); }
261  if(subCutoffFlag) { p->ActivateSubCutoff(true); }
262  if(rndmStepFlag) { p->SetRandomStep(true); }
263  if(stepFunctionActive) { p->SetStepFunction(maxRangeVariation,
264  maxFinalStep); }
265  if(integralActive) { p->SetIntegral(integral); }
266  if(minEnergyActive) { p->SetMinKinEnergy(minKinEnergy); }
267  if(maxEnergyActive) { p->SetMaxKinEnergy(maxKinEnergy); }
268 }
269 
270 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
271 
273 {
274  if(!p) { return; }
275  for (G4int i=0; i<n_loss; ++i) {
276  if(loss_vector[i] == p) { loss_vector[i] = 0; }
277  }
278 }
279 
280 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
281 
283 {
284  if(!p) { return; }
285  G4int n = msc_vector.size();
286  for (G4int i=0; i<n; ++i) {
287  if(msc_vector[i] == p) { return; }
288  }
289  if(verbose > 1) {
290  G4cout << "G4LossTableManager::Register G4VMultipleScattering : "
291  << p->GetProcessName() << " idx= " << msc_vector.size() << G4endl;
292  }
293  msc_vector.push_back(p);
294 }
295 
296 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
297 
299 {
300  if(!p) { return; }
301  size_t msc = msc_vector.size();
302  for (size_t i=0; i<msc; ++i) {
303  if(msc_vector[i] == p) { msc_vector[i] = 0; }
304  }
305 }
306 
307 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
308 
310 {
311  if(!p) { return; }
312  G4int n = emp_vector.size();
313  for (G4int i=0; i<n; ++i) {
314  if(emp_vector[i] == p) { return; }
315  }
316  if(verbose > 1) {
317  G4cout << "G4LossTableManager::Register G4VEmProcess : "
318  << p->GetProcessName() << " idx= " << emp_vector.size() << G4endl;
319  }
320  emp_vector.push_back(p);
321 }
322 
323 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
324 
326 {
327  if(!p) { return; }
328  size_t emp = emp_vector.size();
329  for (size_t i=0; i<emp; ++i) {
330  if(emp_vector[i] == p) { emp_vector[i] = 0; }
331  }
332 }
333 
334 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
335 
337 {
338  mod_vector.push_back(p);
339  if(verbose > 1) {
340  G4cout << "G4LossTableManager::Register G4VEmModel : "
341  << p->GetName() << G4endl;
342  }
343 }
344 
345 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
346 
348 {
349  size_t n = mod_vector.size();
350  for (size_t i=0; i<n; ++i) {
351  if(mod_vector[i] == p) { mod_vector[i] = 0; }
352  }
353 }
354 
355 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
356 
358 {
359  fmod_vector.push_back(p);
360  if(verbose > 1) {
361  G4cout << "G4LossTableManager::Register G4VEmFluctuationModel : "
362  << p->GetName() << G4endl;
363  }
364 }
365 
366 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
367 
369 {
370  size_t n = fmod_vector.size();
371  for (size_t i=0; i<n; ++i) {
372  if(fmod_vector[i] == p) { fmod_vector[i] = 0; }
373  }
374 }
375 
376 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
377 
379  const G4ParticleDefinition* part,
381 {
382  if(!p || !part) { return; }
383  for (G4int i=0; i<n_loss; ++i) {
384  if(loss_vector[i] == p) { return; }
385  }
386  if(verbose > 1) {
387  G4cout << "G4LossTableManager::RegisterExtraParticle "
388  << part->GetParticleName() << " G4VEnergyLossProcess : "
389  << p->GetProcessName() << " idx= " << n_loss << G4endl;
390  }
391  ++n_loss;
392  loss_vector.push_back(p);
393  part_vector.push_back(part);
394  base_part_vector.push_back(p->BaseParticle());
395  dedx_vector.push_back(0);
396  range_vector.push_back(0);
397  inv_range_vector.push_back(0);
398  tables_are_built.push_back(false);
399  all_tables_are_built = false;
400 }
401 
402 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
403 
404 void
407  G4bool theMaster)
408 {
409  if (1 < verbose) {
410  G4cout << "G4LossTableManager::PreparePhysicsTable for "
411  << particle->GetParticleName()
412  << " and " << p->GetProcessName() << " run= " << run
413  << " loss_vector " << loss_vector.size() << G4endl;
414  }
415 
416  isMaster = theMaster;
417 
418  if(!startInitialisation) {
419  tableBuilder->SetInitialisationFlag(false);
420  if (1 < verbose) {
421  G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
422  << G4endl;
423  }
424  }
425 
426  // start initialisation for the first run
427  if( -1 == run ) {
428  emConfigurator->PrepareModels(particle, p);
429 
430  // initialise particles for given process
431  for (G4int j=0; j<n_loss; ++j) {
432  if (p == loss_vector[j] && !part_vector[j]) {
433  part_vector[j] = particle;
434  if(particle->GetParticleName() == "GenericIon") {
435  theGenericIon = particle;
436  }
437  }
438  }
439  }
440  startInitialisation = true;
441 }
442 
443 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
444 
445 void
447  G4VEmProcess* p, G4bool theMaster)
448 {
449  if (1 < verbose) {
450  G4cout << "G4LossTableManager::PreparePhysicsTable for "
451  << particle->GetParticleName()
452  << " and " << p->GetProcessName() << G4endl;
453  }
454  isMaster = theMaster;
455 
456  if(!startInitialisation) {
457  tableBuilder->SetInitialisationFlag(false);
458  if (1 < verbose) {
459  G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
460  << G4endl;
461  }
462  }
463 
464  // start initialisation for the first run
465  if( -1 == run ) {
466  emConfigurator->PrepareModels(particle, p);
467  }
468  startInitialisation = true;
469 }
470 
471 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
472 
473 void
476  G4bool theMaster)
477 {
478  if (1 < verbose) {
479  G4cout << "G4LossTableManager::PreparePhysicsTable for "
480  << particle->GetParticleName()
481  << " and " << p->GetProcessName() << G4endl;
482  }
483 
484  isMaster = theMaster;
485 
486  if(!startInitialisation) {
487  tableBuilder->SetInitialisationFlag(false);
488  if (1 < verbose) {
489  G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
490  << G4endl;
491  }
492  }
493 
494  // start initialisation for the first run
495  if( -1 == run ) {
496  emConfigurator->PrepareModels(particle, p);
497  }
498  startInitialisation = true;
499 }
500 
501 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
502 
503 void
505 {
506  if(-1 == run && startInitialisation) {
507  emConfigurator->Clear();
508  }
509 }
510 
511 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
512 
514  const G4ParticleDefinition* aParticle,
516 {
517  if(1 < verbose) {
518  G4cout << "### G4LossTableManager::SlavePhysicsTable() for "
519  << aParticle->GetParticleName()
520  << " and process " << p->GetProcessName()
521  << G4endl;
522  }
523 
524  if(-1 == run && startInitialisation) {
525  emConfigurator->Clear();
526  firstParticle = aParticle;
527  }
528 
529  if(startInitialisation) {
530  ++run;
531  SetVerbose(verbose);
532  if(1 < verbose) {
533  G4cout << "===== G4LossTableManager::SlavePhysicsTable() for run "
534  << run << " =====" << G4endl;
535  }
536  if(atomDeexcitation) {
537  atomDeexcitation->InitialiseAtomicDeexcitation();
538  }
539  currentParticle = 0;
540  startInitialisation = false;
541  for (G4int i=0; i<n_loss; ++i) {
542  if(loss_vector[i]) {
543  tables_are_built[i] = false;
544  } else {
545  tables_are_built[i] = true;
546  part_vector[i] = 0;
547  }
548  }
549  }
550 
551  all_tables_are_built= true;
552  for (G4int i=0; i<n_loss; ++i) {
553  if(p == loss_vector[i]) {
554  tables_are_built[i] = true;
555  isActive[i] = true;
556  /*
557  const G4ProcessManager* pm = p->GetProcessManager();
558  isActive[i] = false;
559  if(pm) { isActive[i] = pm->GetProcessActivation(p); }
560  */
561  part_vector[i] = p->Particle();
562  base_part_vector[i] = p->BaseParticle();
563  dedx_vector[i] = p->DEDXTable();
564  range_vector[i] = p->RangeTableForLoss();
565  inv_range_vector[i] = p->InverseRangeTable();
566  if(0 == run && p->IsIonisationProcess()) {
567  loss_map[part_vector[i]] = loss_vector[i];
568  }
569 
570  if(1 < verbose) {
571  G4cout << i <<". "<< p->GetProcessName();
572  if(part_vector[i]) {
573  G4cout << " for " << part_vector[i]->GetParticleName();
574  }
575  G4cout << " active= " << isActive[i]
576  << " table= " << tables_are_built[i]
577  << " isIonisation= " << p->IsIonisationProcess()
578  << G4endl;
579  }
580  break;
581  } else if(!tables_are_built[i]) {
582  all_tables_are_built = false;
583  }
584  }
585 
586  // Set run time parameters
587  SetParameters(aParticle, p);
588 
589  if(1 < verbose) {
590  G4cout << "### G4LossTableManager::SlavePhysicsTable end"
591  << G4endl;
592  }
593  if(all_tables_are_built) {
594  if(1 < verbose) {
595  G4cout << "%%%%% All dEdx and Range tables for worker are ready for run "
596  << run << " %%%%%" << G4endl;
597  }
598  }
599 }
600 
601 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
602 
604  const G4ParticleDefinition* aParticle,
606 {
607  if(1 < verbose) {
608  G4cout << "### G4LossTableManager::BuildPhysicsTable() for "
609  << aParticle->GetParticleName()
610  << " and process " << p->GetProcessName() << G4endl;
611  }
612  // clear configurator
613  if(-1 == run && startInitialisation) {
614  emConfigurator->Clear();
615  firstParticle = aParticle;
616  }
617  if(startInitialisation) {
618  ++run;
619  if(1 < verbose) {
620  G4cout << "===== G4LossTableManager::BuildPhysicsTable() for run "
621  << run << " =====" << G4endl;
622  }
623  if(atomDeexcitation) {
624  atomDeexcitation->InitialiseAtomicDeexcitation();
625  }
626  currentParticle = 0;
627  all_tables_are_built= true;
628  }
629 
630  // initialisation before any table is built
631  if ( startInitialisation && aParticle == firstParticle ) {
632 
633  startInitialisation = false;
634  if(1 < verbose) {
635  G4cout << "### G4LossTableManager start initilisation for first particle "
636  << firstParticle->GetParticleName()
637  << G4endl;
638  }
639  for (G4int i=0; i<n_loss; ++i) {
640  G4VEnergyLossProcess* el = loss_vector[i];
641 
642  if(el) {
643  isActive[i] = true;
644  /*
645  const G4ProcessManager* pm = el->GetProcessManager();
646  isActive[i] = false;
647  if(pm) { isActive[i] = pm->GetProcessActivation(el); }
648  */
649  base_part_vector[i] = el->BaseParticle();
650  tables_are_built[i] = false;
651  all_tables_are_built= false;
652  if(!isActive[i]) {
653  el->SetIonisation(false);
654  tables_are_built[i] = true;
655  }
656 
657  if(1 < verbose) {
658  G4cout << i <<". "<< el->GetProcessName();
659  if(el->Particle()) {
660  G4cout << " for " << el->Particle()->GetParticleName();
661  }
662  G4cout << " active= " << isActive[i]
663  << " table= " << tables_are_built[i]
664  << " isIonisation= " << el->IsIonisationProcess();
665  if(base_part_vector[i]) {
666  G4cout << " base particle "
667  << base_part_vector[i]->GetParticleName();
668  }
669  G4cout << G4endl;
670  }
671  } else {
672  tables_are_built[i] = true;
673  part_vector[i] = 0;
674  isActive[i] = false;
675  }
676  }
677  }
678 
679  // Set run time parameters
680  SetParameters(aParticle, p);
681 
682  if (all_tables_are_built) { return; }
683 
684  // Build tables for given particle
685  all_tables_are_built = true;
686 
687  for(G4int i=0; i<n_loss; ++i) {
688  if(p == loss_vector[i] && !tables_are_built[i] && !base_part_vector[i]) {
689  const G4ParticleDefinition* curr_part = part_vector[i];
690  if(1 < verbose) {
691  G4cout << "### Build Table for " << p->GetProcessName()
692  << " and " << curr_part->GetParticleName() << G4endl;
693  }
694  G4VEnergyLossProcess* curr_proc = BuildTables(curr_part);
695  if(curr_proc) { CopyTables(curr_part, curr_proc); }
696  }
697  if ( !tables_are_built[i] ) { all_tables_are_built = false; }
698  }
699  if(0 == run && p->IsIonisationProcess()) { loss_map[aParticle] = p; }
700 
701  if(1 < verbose) {
702  G4cout << "### G4LossTableManager::BuildPhysicsTable end: "
703  << "all_tables_are_built= " << all_tables_are_built
704  << G4endl;
705  }
706  if(all_tables_are_built) {
707  if(1 < verbose) {
708  G4cout << "%%%%% All dEdx and Range tables are built for master run= "
709  << run << " %%%%%" << G4endl;
710  }
711  }
712 }
713 
714 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
715 
716 void G4LossTableManager::CopyTables(const G4ParticleDefinition* part,
717  G4VEnergyLossProcess* base_proc)
718 {
719  for (G4int j=0; j<n_loss; ++j) {
720 
721  G4VEnergyLossProcess* proc = loss_vector[j];
722 
723  if (!tables_are_built[j] && part == base_part_vector[j]) {
724  tables_are_built[j] = true;
725  proc->SetDEDXTable(base_proc->DEDXTable(),fRestricted);
726  proc->SetDEDXTable(base_proc->DEDXTableForSubsec(),fSubRestricted);
727  proc->SetDEDXTable(base_proc->DEDXunRestrictedTable(),fTotal);
728  proc->SetCSDARangeTable(base_proc->CSDARangeTable());
729  proc->SetRangeTableForLoss(base_proc->RangeTableForLoss());
730  proc->SetInverseRangeTable(base_proc->InverseRangeTable());
731  proc->SetLambdaTable(base_proc->LambdaTable());
732  proc->SetSubLambdaTable(base_proc->SubLambdaTable());
733  proc->SetIonisation(base_proc->IsIonisationProcess());
734  if(proc->IsIonisationProcess()) {
735  range_vector[j] = base_proc->RangeTableForLoss();
736  inv_range_vector[j] = base_proc->InverseRangeTable();
737  loss_map[part_vector[j]] = proc;
738  }
739  if (1 < verbose) {
740  G4cout << "For " << proc->GetProcessName()
741  << " for " << part_vector[j]->GetParticleName()
742  << " base_part= " << part->GetParticleName()
743  << " tables are assigned "
744  << G4endl;
745  }
746  }
747 
748  if (theElectron == part && theElectron == proc->SecondaryParticle() ) {
749  proc->SetSecondaryRangeTable(base_proc->RangeTableForLoss());
750  }
751  }
752 }
753 
754 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
755 
756 G4VEnergyLossProcess* G4LossTableManager::BuildTables(
757  const G4ParticleDefinition* aParticle)
758 {
759  if(1 < verbose) {
760  G4cout << "G4LossTableManager::BuildTables() for "
761  << aParticle->GetParticleName() << G4endl;
762  }
763 
764  std::vector<G4PhysicsTable*> t_list;
765  std::vector<G4VEnergyLossProcess*> loss_list;
766  loss_list.clear();
767  G4VEnergyLossProcess* em = 0;
768  G4VEnergyLossProcess* p = 0;
769  G4int iem = 0;
770  G4PhysicsTable* dedx = 0;
771  G4int i;
772 
773  G4ProcessVector* pvec =
774  aParticle->GetProcessManager()->GetProcessList();
775  G4int nvec = pvec->size();
776 
777  for (i=0; i<n_loss; ++i) {
778  p = loss_vector[i];
779  if (p) {
780  G4bool yes = (aParticle == part_vector[i]);
781 
782  // possible case of process sharing between particle/anti-particle
783  if(!yes) {
784  G4VProcess* ptr = static_cast<G4VProcess*>(p);
785  for(G4int j=0; j<nvec; ++j) {
786  //G4cout << "j= " << j << " " << (*pvec)[j] << " " << ptr << G4endl;
787  if(ptr == (*pvec)[j]) {
788  yes = true;
789  break;
790  }
791  }
792  }
793  // process belong to this particle
794  if(yes && isActive[i]) {
795  if (p->IsIonisationProcess() || !em) {
796  em = p;
797  iem= i;
798  }
799  // tables may be shared between particle/anti-particle
800  if (!tables_are_built[i]) {
801  dedx = p->BuildDEDXTable(fRestricted);
802  //G4cout << "Build DEDX table for " << p->GetProcessName()
803  // << " idx= " << i << dedx << " " << dedx->length() << G4endl;
804  p->SetDEDXTable(dedx,fRestricted);
805  tables_are_built[i] = true;
806  } else {
807  dedx = p->DEDXTable();
808  }
809  t_list.push_back(dedx);
810  loss_list.push_back(p);
811  }
812  }
813  }
814 
815  G4int n_dedx = t_list.size();
816  if (0 == n_dedx || !em) {
817  G4cout << "G4LossTableManager WARNING: no DEDX processes for "
818  << aParticle->GetParticleName() << G4endl;
819  return 0;
820  }
821  G4int nSubRegions = em->NumberOfSubCutoffRegions();
822 
823  if (1 < verbose) {
824  G4cout << "G4LossTableManager::BuildTables() start to build range tables"
825  << " and the sum of " << n_dedx << " processes"
826  << " iem= " << iem << " em= " << em->GetProcessName()
827  << " buildCSDARange= " << buildCSDARange
828  << " nSubRegions= " << nSubRegions
829  << G4endl;
830  }
831 
832  dedx = em->DEDXTable();
833  em->SetIonisation(true);
834 
835  if (1 < n_dedx) {
836  em->SetDEDXTable(dedx, fIsIonisation);
837  dedx = 0;
839  tableBuilder->BuildDEDXTable(dedx, t_list);
840  em->SetDEDXTable(dedx, fRestricted);
841  }
842  /*
843  if(2==run && "e-" == aParticle->GetParticleName()) {
844  G4cout << "G4LossTableManager::BuildTables for e- " << dedx << G4endl;
845  G4cout << (*dedx) << G4endl;
846  G4cout << "%%%%% Instance ID= " << (*dedx)[0]->GetInstanceID() << G4endl;
847  G4cout << "%%%%% LastValue= " << (*dedx)[0]->GetLastValue() << G4endl;
848  G4cout << "%%%%% 1.2 " << (*(dedx))[0]->Value(1.2) << G4endl;
849  }
850  */
851  dedx_vector[iem] = dedx;
852 
853  G4PhysicsTable* range = em->RangeTableForLoss();
854  if(!range) range = G4PhysicsTableHelper::PreparePhysicsTable(range);
855  range_vector[iem] = range;
856 
857  G4PhysicsTable* invrange = em->InverseRangeTable();
858  if(!invrange) invrange = G4PhysicsTableHelper::PreparePhysicsTable(invrange);
859  inv_range_vector[iem] = invrange;
860 
861  tableBuilder->BuildRangeTable(dedx, range, true);
862  tableBuilder->BuildInverseRangeTable(range, invrange, true);
863 
864  // if(1<verbose) G4cout << *dedx << G4endl;
865 
866  em->SetRangeTableForLoss(range);
867  em->SetInverseRangeTable(invrange);
868 
869  // if(1<verbose) G4cout << *range << G4endl;
870 
871  std::vector<G4PhysicsTable*> listSub;
872  std::vector<G4PhysicsTable*> listCSDA;
873 
874  for (i=0; i<n_dedx; ++i) {
875  p = loss_list[i];
876  if(p != em) { p->SetIonisation(false); }
878  if (0 < nSubRegions) {
879  dedx = p->BuildDEDXTable(fSubRestricted);
880  p->SetDEDXTable(dedx,fSubRestricted);
881  listSub.push_back(dedx);
883  if(p != em) { em->AddCollaborativeProcess(p); }
884  }
885  if(buildCSDARange) {
886  dedx = p->BuildDEDXTable(fTotal);
887  p->SetDEDXTable(dedx,fTotal);
888  listCSDA.push_back(dedx);
889  }
890  }
891 
892  if (0 < nSubRegions) {
893  G4PhysicsTable* dedxSub = em->IonisationTableForSubsec();
894  if (1 < listSub.size()) {
895  em->SetDEDXTable(dedxSub, fIsSubIonisation);
896  dedxSub = 0;
898  tableBuilder->BuildDEDXTable(dedxSub, listSub);
899  em->SetDEDXTable(dedxSub, fSubRestricted);
900  }
901  }
902  if(buildCSDARange) {
903  G4PhysicsTable* dedxCSDA = em->DEDXunRestrictedTable();
904  if (1 < n_dedx) {
905  dedxCSDA = 0;
906  dedxCSDA = G4PhysicsTableHelper::PreparePhysicsTable(dedxCSDA);
907  tableBuilder->BuildDEDXTable(dedxCSDA, listCSDA);
908  em->SetDEDXTable(dedxCSDA,fTotal);
909  }
910  G4PhysicsTable* rCSDA = em->CSDARangeTable();
911  if(!rCSDA) { rCSDA = G4PhysicsTableHelper::PreparePhysicsTable(rCSDA); }
912  tableBuilder->BuildRangeTable(dedxCSDA, rCSDA, true);
913  em->SetCSDARangeTable(rCSDA);
914  }
915 
916  if (1 < verbose) {
917  G4cout << "G4LossTableManager::BuildTables: Tables are built for "
918  << aParticle->GetParticleName()
919  << "; ionisation process: " << em->GetProcessName()
920  << G4endl;
921  }
922  return em;
923 }
924 
925 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
926 
928 {
929  return theMessenger;
930 }
931 
932 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
933 
934 void G4LossTableManager::ParticleHaveNoLoss(
935  const G4ParticleDefinition* aParticle)
936 {
938  ed << "Energy loss process not found for " << aParticle->GetParticleName()
939  << " !";
940  G4Exception("G4LossTableManager::ParticleHaveNoLoss", "em0001",
941  FatalException, ed);
942 
943 }
944 
945 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
946 
948 {
949  return buildCSDARange;
950 }
951 
952 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
953 
955 {
956  lossFluctuationFlag = val;
957  for(G4int i=0; i<n_loss; ++i) {
958  if(loss_vector[i]) { loss_vector[i]->SetLossFluctuations(val); }
959  }
960 }
961 
962 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
963 
965 {
966  subCutoffFlag = val;
967  for(G4int i=0; i<n_loss; ++i) {
968  if(loss_vector[i]) { loss_vector[i]->ActivateSubCutoff(val, r); }
969  }
970 }
971 
972 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
973 
975 {
976  integral = val;
977  integralActive = true;
978  for(G4int i=0; i<n_loss; ++i) {
979  if(loss_vector[i]) { loss_vector[i]->SetIntegral(val); }
980  }
981  size_t emp = emp_vector.size();
982  for (size_t k=0; k<emp; ++k) {
983  if(emp_vector[k]) { emp_vector[k]->SetIntegral(val); }
984  }
985 }
986 
987 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
988 
990 {
991  minSubRange = val;
992  for(G4int i=0; i<n_loss; ++i) {
993  if(loss_vector[i]) { loss_vector[i]->SetMinSubRange(val); }
994  }
995 }
996 
997 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
998 
1000 {
1001  rndmStepFlag = val;
1002  for(G4int i=0; i<n_loss; ++i) {
1003  if(loss_vector[i]) { loss_vector[i]->SetRandomStep(val); }
1004  }
1005 }
1006 
1007 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1008 
1010 {
1011  minEnergyActive = true;
1012  minKinEnergy = val;
1013  for(G4int i=0; i<n_loss; ++i) {
1014  if(loss_vector[i]) { loss_vector[i]->SetMinKinEnergy(val); }
1015  }
1016  size_t emp = emp_vector.size();
1017  for (size_t k=0; k<emp; ++k) {
1018  if(emp_vector[k]) { emp_vector[k]->SetMinKinEnergy(val); }
1019  }
1020 }
1021 
1022 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1023 
1025 {
1026  maxEnergyActive = true;
1027  maxKinEnergy = val;
1028  for(G4int i=0; i<n_loss; ++i) {
1029  if(loss_vector[i]) { loss_vector[i]->SetMaxKinEnergy(val); }
1030  }
1031  size_t emp = emp_vector.size();
1032  for (size_t k=0; k<emp; ++k) {
1033  if(emp_vector[k]) { emp_vector[k]->SetMaxKinEnergy(val); }
1034  }
1035 }
1036 
1037 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1038 
1040 {
1041  for(G4int i=0; i<n_loss; ++i) {
1042  if(loss_vector[i]) { loss_vector[i]->SetMaxKinEnergyForCSDARange(val); }
1043  }
1044 }
1045 
1046 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1047 
1049 {
1050  maxEnergyForMuonsActive = true;
1051  maxKinEnergyForMuons = val;
1052 }
1053 
1054 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1055 
1057 {
1058  for(G4int i=0; i<n_loss; ++i) {
1059  if(loss_vector[i]) { loss_vector[i]->SetDEDXBinning(val); }
1060  }
1061 }
1062 
1063 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1064 
1066 {
1067  for(G4int i=0; i<n_loss; ++i) {
1068  if(loss_vector[i]) { loss_vector[i]->SetDEDXBinningForCSDARange(val); }
1069  }
1070 }
1071 
1072 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1073 
1075 {
1076  G4int n = val/G4int(std::log10(maxKinEnergy/minKinEnergy) + 0.5);
1077  if(n < 5) {
1078  G4cout << "G4LossTableManager::SetLambdaBinning WARNING "
1079  << "too small number of bins " << val << " ignored"
1080  << G4endl;
1081  return;
1082  }
1083  nbinsLambda = val;
1084  nbinsPerDecade = n;
1085  size_t emp = emp_vector.size();
1086  for (size_t k=0; k<emp; ++k) {
1087  if(emp_vector[k]) { emp_vector[k]->SetLambdaBinning(val); }
1088  }
1089 }
1090 
1091 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1092 
1094 {
1095  return nbinsPerDecade;
1096 }
1097 
1098 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1099 
1101 {
1102  verbose = val;
1103  for(G4int i=0; i<n_loss; ++i) {
1104  if(loss_vector[i]) { loss_vector[i]->SetVerboseLevel(val); }
1105  }
1106  size_t msc = msc_vector.size();
1107  for (size_t j=0; j<msc; ++j) {
1108  if(msc_vector[j]) { msc_vector[j]->SetVerboseLevel(val); }
1109  }
1110  size_t emp = emp_vector.size();
1111  for (size_t k=0; k<emp; ++k) {
1112  if(emp_vector[k]) { emp_vector[k]->SetVerboseLevel(val); }
1113  }
1114  emConfigurator->SetVerbose(val);
1115  //tableBuilder->SetVerbose(val);
1116  //emCorrections->SetVerbose(val);
1117  emSaturation->SetVerbose(val);
1118  emElectronIonPair->SetVerbose(val);
1119  if(atomDeexcitation) { atomDeexcitation->SetVerboseLevel(val); }
1120 }
1121 
1122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1123 
1125 {
1126  stepFunctionActive = true;
1127  maxRangeVariation = v1;
1128  maxFinalStep = v2;
1129  for(G4int i=0; i<n_loss; ++i) {
1130  if(loss_vector[i]) { loss_vector[i]->SetStepFunction(v1, v2); }
1131  }
1132 }
1133 
1134 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1135 
1137 {
1138  for(G4int i=0; i<n_loss; ++i) {
1139  if(loss_vector[i]) { loss_vector[i]->SetLinearLossLimit(val); }
1140  }
1141 }
1142 
1143 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1144 
1146 {
1147  buildCSDARange = val;
1148 }
1149 
1150 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1151 
1152 void
1153 G4LossTableManager::SetParameters(const G4ParticleDefinition* aParticle,
1155 {
1156  if(stepFunctionActive) {
1157  p->SetStepFunction(maxRangeVariation, maxFinalStep);
1158  }
1159  if(integralActive) { p->SetIntegral(integral); }
1160  if(minEnergyActive) { p->SetMinKinEnergy(minKinEnergy); }
1161  if(maxEnergyActive) { p->SetMaxKinEnergy(maxKinEnergy); }
1162  // p->SetVerboseLevel(verbose);
1163  if(maxEnergyForMuonsActive) {
1164  G4double dm = std::abs(aParticle->GetPDGMass() - 105.7*MeV);
1165  if(dm < 5.*MeV) { p->SetMaxKinEnergy(maxKinEnergyForMuons); }
1166  }
1167 }
1168 
1169 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1170 
1171 const std::vector<G4VEnergyLossProcess*>&
1173 {
1174  return loss_vector;
1175 }
1176 
1177 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1178 
1179 const std::vector<G4VEmProcess*>& G4LossTableManager::GetEmProcessVector()
1180 {
1181  return emp_vector;
1182 }
1183 
1184 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1185 
1186 const std::vector<G4VMultipleScattering*>&
1188 {
1189  return msc_vector;
1190 }
1191 
1192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1193 
1195 {
1196  flagLPM = val;
1197 }
1198 
1199 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1200 
1202 {
1203  return flagLPM;
1204 }
1205 
1206 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1207 
1209 {
1210  splineFlag = val;
1211  tableBuilder->SetSplineFlag(splineFlag);
1212 }
1213 
1214 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1215 
1217 {
1218  return splineFlag;
1219 }
1220 
1221 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1222 
1224 {
1225  return isMaster;
1226 }
1227 
1228 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1229 
1231 {
1232  bremsTh = val;
1233 }
1234 
1235 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1236 
1238 {
1239  return bremsTh;
1240 }
1241 
1242 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1243 
1245 {
1246  if(val > 0.0) { factorForAngleLimit = val; }
1247 }
1248 
1249 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1250 
1252 {
1253  return factorForAngleLimit;
1254 }
1255 
1256 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1257 
1259 {
1260  return minKinEnergy;
1261 }
1262 
1263 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1264 
1266 {
1267  return maxKinEnergy;
1268 }
1269 
1270 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1271 
1273 {
1274  return emCorrections;
1275 }
1276 
1277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1278 
1280 {
1281  return emSaturation;
1282 }
1283 
1284 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1285 
1287 {
1288  return emConfigurator;
1289 }
1290 
1291 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1292 
1294 {
1295  return emElectronIonPair;
1296 }
1297 
1298 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1299 
1301 {
1302  return atomDeexcitation;
1303 }
1304 
1305 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1306 
1308 {
1309  return tableBuilder;
1310 }
1311 
1312 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1313 
1315 {
1316  atomDeexcitation = p;
1317 }
1318 
1319 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1320 
1323 {
1324  if(aParticle != currentParticle) {
1325  currentParticle = aParticle;
1326  std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
1327  if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
1328  currentLoss = (*pos).second;
1329  } else {
1330  currentLoss = 0;
1331  if(aParticle->GetParticleType() == "nucleus") {
1332  if ((pos = loss_map.find(theGenericIon)) != loss_map.end()) {
1333  currentLoss = (*pos).second;
1334  }
1335  }
1336  }
1337  }
1338  return currentLoss;
1339 }
1340 
1341 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1342 
1344  G4double kineticEnergy,
1345  const G4MaterialCutsCouple *couple)
1346 {
1347  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1348  G4double x = 0.0;
1349  if(currentLoss) { x = currentLoss->GetDEDX(kineticEnergy, couple); }
1350  return x;
1351 }
1352 
1353 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1354 
1356  G4double kineticEnergy,
1357  const G4MaterialCutsCouple *couple)
1358 {
1359  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1360  G4double x = 0.0;
1361  if(currentLoss) { x = currentLoss->GetDEDXForSubsec(kineticEnergy, couple); }
1362  return x;
1363 }
1364 
1365 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1366 
1368  G4double kineticEnergy,
1369  const G4MaterialCutsCouple *couple)
1370 {
1371  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1372  G4double x = DBL_MAX;
1373  if(currentLoss) { x = currentLoss->GetCSDARange(kineticEnergy, couple); }
1374  return x;
1375 }
1376 
1377 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1378 
1380  const G4ParticleDefinition *aParticle,
1381  G4double kineticEnergy,
1382  const G4MaterialCutsCouple *couple)
1383 {
1384  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1385  G4double x = DBL_MAX;
1386  if(currentLoss) { x = currentLoss->GetRangeForLoss(kineticEnergy, couple); }
1387  return x;
1388 }
1389 
1390 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1391 
1393  G4double kineticEnergy,
1394  const G4MaterialCutsCouple *couple)
1395 {
1396  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1397  G4double x = DBL_MAX;
1398  if(currentLoss) { x = currentLoss->GetRange(kineticEnergy, couple); }
1399  return x;
1400 }
1401 
1402 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1403 
1405  G4double range,
1406  const G4MaterialCutsCouple *couple)
1407 {
1408  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1409  G4double x = 0;
1410  if(currentLoss) { x = currentLoss->GetKineticEnergy(range, couple); }
1411  return x;
1412 }
1413 
1414 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1415 
1417  const G4MaterialCutsCouple *couple,
1418  const G4DynamicParticle* dp,
1419  G4double& length)
1420 {
1421  const G4ParticleDefinition* aParticle = dp->GetParticleDefinition();
1422  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1423  G4double x = 0.0;
1424  if(currentLoss) { currentLoss->GetDEDXDispersion(couple, dp, length); }
1425  return x;
1426 }
1427 
1428 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
G4EmConfigurator * EmConfigurator()
void SetRandomStep(G4bool val)
void BuildRangeTable(const G4PhysicsTable *dedxTable, G4PhysicsTable *rangeTable, G4bool isIonisation=false)
void SetDEDXTable(G4PhysicsTable *p, G4EmTableType tType)
void SetIntegral(G4bool val)
G4bool SplineFlag() const
static G4LossTableManager * Instance()
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
void SetIonisation(G4bool val)
void SetBremsstrahlungTh(G4double val)
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4PhysicsTable * SubLambdaTable() const
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double &length)
void DeRegister(G4VEnergyLossProcess *p)
G4PhysicsTable * RangeTableForLoss() const
void SetLambdaBinning(G4int val)
const char * p
Definition: xmltok.h:285
G4double GetSubDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4bool BuildCSDARange() const
void push_back(G4PhysicsVector *)
G4double GetRangeForLoss(G4double &kineticEnergy, const G4MaterialCutsCouple *)
void AddCollaborativeProcess(G4VEnergyLossProcess *)
void SetStepFunction(G4double v1, G4double v2)
G4PhysicsTable * CSDARangeTable() const
G4double GetCSDARange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4EnergyLossMessenger * GetMessenger()
G4double FactorForAngleLimit() const
void SetFactorForAngleLimit(G4double val)
G4PhysicsTable * IonisationTableForSubsec() const
G4double GetDEDXForSubsec(G4double &kineticEnergy, const G4MaterialCutsCouple *)
G4ProcessManager * GetProcessManager() const
const std::vector< G4VEmProcess * > & GetEmProcessVector()
int G4int
Definition: G4Types.hh:78
G4PhysicsTable * BuildLambdaTable(G4EmTableType tType=fRestricted)
void SetVerbose(G4int)
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
const G4String & GetParticleName() const
G4double GetCSDARange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4LossTableBuilder * GetTableBuilder()
void SetBuildCSDARange(G4bool val)
void SetInitialisationFlag(G4bool flag)
G4PhysicsTable * LambdaTable() const
void SetInverseRangeTable(G4PhysicsTable *p)
const G4ParticleDefinition * SecondaryParticle() const
G4GLOB_DLL std::ostream G4cout
G4PhysicsTable * DEDXTable() const
G4double MinKinEnergy() const
void SetLPMFlag(G4bool val)
bool G4bool
Definition: G4Types.hh:79
G4EmCorrections * EmCorrections()
void SetLossFluctuations(G4bool val)
G4ElectronIonPair * ElectronIonPair()
G4EmSaturation * EmSaturation()
void SetMaxEnergyForCSDARange(G4double val)
G4bool IsMaster() const
G4int NumberOfSubCutoffRegions() const
const G4String & GetParticleType() const
void SetMaxKinEnergy(G4double e)
void Register(G4VEnergyLossProcess *p)
const G4ParticleDefinition * BaseParticle() const
const G4ParticleDefinition * GetParticleDefinition() const
const G4int n
void SetDEDXBinningForCSDARange(G4int val)
void BuildInverseRangeTable(const G4PhysicsTable *rangeTable, G4PhysicsTable *invRangeTable, G4bool isIonisation=false)
G4PhysicsTable * DEDXTableForSubsec() const
G4double GetEnergy(const G4ParticleDefinition *aParticle, G4double range, const G4MaterialCutsCouple *couple)
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
void SetSubCutoff(G4bool val, const G4Region *r=0)
G4double BremsstrahlungTh() const
G4double GetRangeFromRestricteDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
void RegisterExtraParticle(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4bool IsWorkerThread()
Definition: G4Threading.cc:104
void BuildDEDXTable(G4PhysicsTable *dedxTable, const std::vector< G4PhysicsTable * > &)
void SetVerbose(G4int value)
void SetMinSubRange(G4double val)
G4int size() const
void SetMinEnergy(G4double val)
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
void SetLambdaTable(G4PhysicsTable *p)
void SetMaxEnergy(G4double val)
void SetStepFunction(G4double v1, G4double v2)
G4PhysicsTable * InverseRangeTable() const
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
G4double GetPDGMass() const
const G4ParticleDefinition * Particle() const
void SetLinearLossLimit(G4double val)
G4int GetNumberOfBinsPerDecade() const
void PrepareModels(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
void SetCSDARangeTable(G4PhysicsTable *pRange)
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
void SetVerbose(G4int val)
void LocalPhysicsTables(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
void SetDEDXBinning(G4int val)
void SetMaxEnergyForMuons(G4double val)
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
void SetIntegral(G4bool val)
void SetLossFluctuations(G4bool val)
const G4String & GetName() const
Definition: G4VEmModel.hh:753
void SetSecondaryRangeTable(G4PhysicsTable *p)
G4VAtomDeexcitation * AtomDeexcitation()
G4PhysicsTable * DEDXunRestrictedTable() const
double G4double
Definition: G4Types.hh:76
void SetSubLambdaTable(G4PhysicsTable *p)
void ActivateSubCutoff(G4bool val, const G4Region *region=0)
void SetRangeTableForLoss(G4PhysicsTable *p)
#define DBL_MAX
Definition: templates.hh:83
G4double MaxKinEnergy() const
G4double GetRange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
void SetAtomDeexcitation(G4VAtomDeexcitation *)
G4PhysicsTable * BuildDEDXTable(G4EmTableType tType=fRestricted)
void SetSplineFlag(G4bool val)
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()
G4ProcessVector * GetProcessList() const
G4double GetRange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
void SetRandomStep(G4bool val)
void SetSplineFlag(G4bool flag)
G4bool IsIonisationProcess() const
void SetMinKinEnergy(G4double e)