Geant4  10.01
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 85424 2014-10-29 08:23:44Z 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 "G4EmParameters.hh"
87 #include "G4EnergyLossMessenger.hh"
88 #include "G4PhysicsTable.hh"
89 #include "G4ParticleDefinition.hh"
90 #include "G4MaterialCutsCouple.hh"
91 #include "G4ProcessManager.hh"
92 #include "G4Electron.hh"
93 #include "G4Proton.hh"
94 #include "G4VMultipleScattering.hh"
95 #include "G4VEmProcess.hh"
96 #include "G4ProductionCutsTable.hh"
97 #include "G4PhysicsTableHelper.hh"
98 #include "G4EmCorrections.hh"
99 #include "G4EmSaturation.hh"
100 #include "G4EmConfigurator.hh"
101 #include "G4ElectronIonPair.hh"
102 #include "G4EmTableType.hh"
103 #include "G4LossTableBuilder.hh"
104 #include "G4VAtomDeexcitation.hh"
105 #include "G4VSubCutProducer.hh"
106 #include "G4Region.hh"
107 #include "G4PhysicalConstants.hh"
108 #include "G4Threading.hh"
109 
110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
111 
113 
115 {
116  if(!instance) {
118  instance = inst.Instance();
119  }
120  return instance;
121 }
122 
123 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
124 
126 {
127  //G4cout << "### G4LossTableManager::~G4LossTableManager()" << G4endl;
128  for (G4int i=0; i<n_loss; ++i) {
129  //G4cout << "### eloss #" << i << G4endl;
130  if( loss_vector[i] ) {
131  delete loss_vector[i];
132  }
133  }
134  size_t msc = msc_vector.size();
135  for (size_t j=0; j<msc; ++j) {
136  if( msc_vector[j] ) { delete msc_vector[j]; }
137  }
138  size_t emp = emp_vector.size();
139  for (size_t k=0; k<emp; ++k) {
140  if( emp_vector[k] ) { delete emp_vector[k]; }
141  }
142  size_t mod = mod_vector.size();
143  size_t fmod = fmod_vector.size();
144  for (size_t a=0; a<mod; ++a) {
145  if( mod_vector[a] ) {
146  for (size_t b=0; b<fmod; ++b) {
147  if((G4VEmModel*)(fmod_vector[b]) == mod_vector[a]) {
148  fmod_vector[b] = 0;
149  }
150  }
151  delete mod_vector[a];
152  }
153  }
154  for (size_t b=0; b<fmod; ++b) {
155  if( fmod_vector[b] ) { delete fmod_vector[b]; }
156  }
157  Clear();
158  delete theMessenger;
159  delete tableBuilder;
160  delete emCorrections;
161  delete emSaturation;
162  delete emConfigurator;
163  delete emElectronIonPair;
164  delete atomDeexcitation;
165  delete subcutProducer;
166 }
167 
168 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
169 
171 {
173  n_loss = 0;
174  run = -1;
175  startInitialisation = false;
176  all_tables_are_built = false;
177  currentLoss = 0;
178  currentParticle = 0;
179  firstParticle = 0;
180  subCutoffFlag = false;
181  maxRangeVariation = 0.2;
183  integral = true;
184  integralActive = false;
185  stepFunctionActive = false;
186  isMaster = true;
190  theGenericIon= 0;
193  isMaster = false;
194  }
197  emSaturation = 0;
198  emConfigurator = 0;
199  emElectronIonPair = 0;
200  atomDeexcitation = 0;
201  subcutProducer = 0;
202 }
203 
204 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
205 
207 {
208  all_tables_are_built = false;
209  currentLoss = 0;
210  currentParticle = 0;
211  if(n_loss)
212  {
213  dedx_vector.clear();
214  range_vector.clear();
215  inv_range_vector.clear();
216  loss_map.clear();
217  loss_vector.clear();
218  part_vector.clear();
219  base_part_vector.clear();
220  tables_are_built.clear();
221  isActive.clear();
222  n_loss = 0;
223  }
224 }
225 
226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
227 
229 {
230  if(!p) { return; }
231  for (G4int i=0; i<n_loss; ++i) {
232  if(loss_vector[i] == p) { return; }
233  }
234  if(verbose > 1) {
235  G4cout << "G4LossTableManager::Register G4VEnergyLossProcess : "
236  << p->GetProcessName() << " idx= " << n_loss << G4endl;
237  }
238  ++n_loss;
239  loss_vector.push_back(p);
240  part_vector.push_back(0);
241  base_part_vector.push_back(0);
242  dedx_vector.push_back(0);
243  range_vector.push_back(0);
244  inv_range_vector.push_back(0);
245  tables_are_built.push_back(false);
246  isActive.push_back(true);
247  all_tables_are_built = false;
248  if(subCutoffFlag) { p->ActivateSubCutoff(true); }
250  maxFinalStep); }
251  if(integralActive) { p->SetIntegral(integral); }
252 }
253 
254 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
255 
257 {
259  if(!isMaster) {
261  }
268  if(atomDeexcitation) {
274  }
275 }
276 
277 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
278 
280 {
281  if(!p) { return; }
282  for (G4int i=0; i<n_loss; ++i) {
283  if(loss_vector[i] == p) { loss_vector[i] = 0; }
284  }
285 }
286 
287 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
288 
290 {
291  if(!p) { return; }
292  G4int n = msc_vector.size();
293  for (G4int i=0; i<n; ++i) {
294  if(msc_vector[i] == p) { return; }
295  }
296  if(verbose > 1) {
297  G4cout << "G4LossTableManager::Register G4VMultipleScattering : "
298  << p->GetProcessName() << " idx= " << msc_vector.size() << G4endl;
299  }
300  msc_vector.push_back(p);
301 }
302 
303 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
304 
306 {
307  if(!p) { return; }
308  size_t msc = msc_vector.size();
309  for (size_t i=0; i<msc; ++i) {
310  if(msc_vector[i] == p) { msc_vector[i] = 0; }
311  }
312 }
313 
314 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
315 
317 {
318  if(!p) { return; }
319  G4int n = emp_vector.size();
320  for (G4int i=0; i<n; ++i) {
321  if(emp_vector[i] == p) { return; }
322  }
323  if(verbose > 1) {
324  G4cout << "G4LossTableManager::Register G4VEmProcess : "
325  << p->GetProcessName() << " idx= " << emp_vector.size() << G4endl;
326  }
327  emp_vector.push_back(p);
328 }
329 
330 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
331 
333 {
334  if(!p) { return; }
335  size_t emp = emp_vector.size();
336  for (size_t i=0; i<emp; ++i) {
337  if(emp_vector[i] == p) { emp_vector[i] = 0; }
338  }
339 }
340 
341 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
342 
344 {
345  mod_vector.push_back(p);
346  if(verbose > 1) {
347  G4cout << "G4LossTableManager::Register G4VEmModel : "
348  << p->GetName() << G4endl;
349  }
350 }
351 
352 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
353 
355 {
356  size_t n = mod_vector.size();
357  for (size_t i=0; i<n; ++i) {
358  if(mod_vector[i] == p) { mod_vector[i] = 0; }
359  }
360 }
361 
362 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
363 
365 {
366  fmod_vector.push_back(p);
367  if(verbose > 1) {
368  G4cout << "G4LossTableManager::Register G4VEmFluctuationModel : "
369  << p->GetName() << G4endl;
370  }
371 }
372 
373 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
374 
376 {
377  size_t n = fmod_vector.size();
378  for (size_t i=0; i<n; ++i) {
379  if(fmod_vector[i] == p) { fmod_vector[i] = 0; }
380  }
381 }
382 
383 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
384 
386  const G4ParticleDefinition* part,
388 {
389  if(!p || !part) { return; }
390  for (G4int i=0; i<n_loss; ++i) {
391  if(loss_vector[i] == p) { return; }
392  }
393  if(verbose > 1) {
394  G4cout << "G4LossTableManager::RegisterExtraParticle "
395  << part->GetParticleName() << " G4VEnergyLossProcess : "
396  << p->GetProcessName() << " idx= " << n_loss << G4endl;
397  }
398  ++n_loss;
399  loss_vector.push_back(p);
400  part_vector.push_back(part);
401  base_part_vector.push_back(p->BaseParticle());
402  dedx_vector.push_back(0);
403  range_vector.push_back(0);
404  inv_range_vector.push_back(0);
405  tables_are_built.push_back(false);
406  all_tables_are_built = false;
407 }
408 
409 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
410 
411 void
414  G4bool theMaster)
415 {
416  if (1 < verbose) {
417  G4cout << "G4LossTableManager::PreparePhysicsTable for "
418  << particle->GetParticleName()
419  << " and " << p->GetProcessName() << " run= " << run
420  << " loss_vector " << loss_vector.size() << G4endl;
421  }
422 
423  isMaster = theMaster;
424 
425  if(!startInitialisation) {
426  ResetParameters();
427  if (1 < verbose) {
428  G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
429  << G4endl;
430  }
431  }
432 
433  // start initialisation for the first run
434  if( -1 == run ) {
435  if(emConfigurator) { emConfigurator->PrepareModels(particle, p); }
436 
437  // initialise particles for given process
438  for (G4int j=0; j<n_loss; ++j) {
439  if (p == loss_vector[j] && !part_vector[j]) {
440  part_vector[j] = particle;
441  if(particle->GetParticleName() == "GenericIon") {
442  theGenericIon = particle;
443  }
444  }
445  }
446  }
447  startInitialisation = true;
448 }
449 
450 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
451 
452 void
454  G4VEmProcess* p, G4bool theMaster)
455 {
456  if (1 < verbose) {
457  G4cout << "G4LossTableManager::PreparePhysicsTable for "
458  << particle->GetParticleName()
459  << " and " << p->GetProcessName() << G4endl;
460  }
461  isMaster = theMaster;
462 
463  if(!startInitialisation) {
464  ResetParameters();
465  if (1 < verbose) {
466  G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
467  << G4endl;
468  }
469  }
470 
471  // start initialisation for the first run
472  if( -1 == run ) {
473  if(emConfigurator) { emConfigurator->PrepareModels(particle, p); }
474  }
475  startInitialisation = true;
476 }
477 
478 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
479 
480 void
483  G4bool theMaster)
484 {
485  if (1 < verbose) {
486  G4cout << "G4LossTableManager::PreparePhysicsTable for "
487  << particle->GetParticleName()
488  << " and " << p->GetProcessName() << G4endl;
489  }
490 
491  isMaster = theMaster;
492 
493  if(!startInitialisation) {
494  ResetParameters();
495  if (1 < verbose) {
496  G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
497  << G4endl;
498  }
499  }
500 
501  // start initialisation for the first run
502  if( -1 == run ) {
503  if(emConfigurator) { emConfigurator->PrepareModels(particle, p); }
504  }
505  startInitialisation = true;
506 }
507 
508 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
509 
510 void
512 {
513  if(-1 == run && startInitialisation) {
515  }
516 }
517 
518 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
519 
521  const G4ParticleDefinition* aParticle,
523 {
524  if(1 < verbose) {
525  G4cout << "### G4LossTableManager::LocalPhysicsTable() for "
526  << aParticle->GetParticleName()
527  << " and process " << p->GetProcessName()
528  << G4endl;
529  }
530 
531  if(-1 == run && startInitialisation) {
533  firstParticle = aParticle;
534  }
535 
536  if(startInitialisation) {
537  ++run;
539  if(1 < verbose) {
540  G4cout << "===== G4LossTableManager::LocalPhysicsTable() for run "
541  << run << " =====" << G4endl;
542  }
544  if(atomDeexcitation) {
546  }
547  currentParticle = 0;
548  startInitialisation = false;
549  for (G4int i=0; i<n_loss; ++i) {
550  if(loss_vector[i]) {
551  tables_are_built[i] = false;
552  } else {
553  tables_are_built[i] = true;
554  part_vector[i] = 0;
555  }
556  }
557  }
558 
559  all_tables_are_built= true;
560  for (G4int i=0; i<n_loss; ++i) {
561  if(p == loss_vector[i]) {
562  tables_are_built[i] = true;
563  isActive[i] = true;
564  /*
565  const G4ProcessManager* pm = p->GetProcessManager();
566  isActive[i] = false;
567  if(pm) { isActive[i] = pm->GetProcessActivation(p); }
568  */
569  part_vector[i] = p->Particle();
570  base_part_vector[i] = p->BaseParticle();
571  dedx_vector[i] = p->DEDXTable();
574  if(0 == run && p->IsIonisationProcess()) {
576  }
577 
578  if(1 < verbose) {
579  G4cout << i <<". "<< p->GetProcessName();
580  if(part_vector[i]) {
581  G4cout << " for " << part_vector[i]->GetParticleName();
582  }
583  G4cout << " active= " << isActive[i]
584  << " table= " << tables_are_built[i]
585  << " isIonisation= " << p->IsIonisationProcess()
586  << G4endl;
587  }
588  break;
589  } else if(!tables_are_built[i]) {
590  all_tables_are_built = false;
591  }
592  }
593 
594  // Set run time parameters
595  SetParameters(aParticle, p);
596 
597  if(1 < verbose) {
598  G4cout << "### G4LossTableManager::LocalPhysicsTable end"
599  << G4endl;
600  }
601  if(all_tables_are_built) {
602  if(1 < verbose) {
603  G4cout << "%%%%% All dEdx and Range tables for worker are ready for run "
604  << run << " %%%%%" << G4endl;
605  }
606  }
607 }
608 
609 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
610 
612  const G4ParticleDefinition* aParticle,
614 {
615  if(1 < verbose) {
616  G4cout << "### G4LossTableManager::BuildPhysicsTable() for "
617  << aParticle->GetParticleName()
618  << " and process " << p->GetProcessName() << G4endl;
619  }
620  // clear configurator
621  if(-1 == run && startInitialisation) {
623  firstParticle = aParticle;
624  }
625  if(startInitialisation) {
626  ++run;
627  if(1 < verbose) {
628  G4cout << "===== G4LossTableManager::BuildPhysicsTable() for run "
629  << run << " ===== " << atomDeexcitation << G4endl;
630  }
631  if(atomDeexcitation) {
633  }
634  currentParticle = 0;
635  all_tables_are_built= true;
636  }
637 
638  // initialisation before any table is built
639  if ( startInitialisation && aParticle == firstParticle ) {
640 
641  startInitialisation = false;
642  if(1 < verbose) {
643  G4cout << "### G4LossTableManager start initilisation for first particle "
645  << G4endl;
646  }
647  for (G4int i=0; i<n_loss; ++i) {
649 
650  if(el) {
651  isActive[i] = true;
652  /*
653  const G4ProcessManager* pm = el->GetProcessManager();
654  isActive[i] = false;
655  if(pm) { isActive[i] = pm->GetProcessActivation(el); }
656  */
657  base_part_vector[i] = el->BaseParticle();
658  tables_are_built[i] = false;
659  all_tables_are_built= false;
660  if(!isActive[i]) {
661  el->SetIonisation(false);
662  tables_are_built[i] = true;
663  }
664 
665  if(1 < verbose) {
666  G4cout << i <<". "<< el->GetProcessName();
667  if(el->Particle()) {
668  G4cout << " for " << el->Particle()->GetParticleName();
669  }
670  G4cout << " active= " << isActive[i]
671  << " table= " << tables_are_built[i]
672  << " isIonisation= " << el->IsIonisationProcess();
673  if(base_part_vector[i]) {
674  G4cout << " base particle "
675  << base_part_vector[i]->GetParticleName();
676  }
677  G4cout << G4endl;
678  }
679  } else {
680  tables_are_built[i] = true;
681  part_vector[i] = 0;
682  isActive[i] = false;
683  }
684  }
685  }
686 
687  // Set run time parameters
688  SetParameters(aParticle, p);
689 
690  if (all_tables_are_built) { return; }
691 
692  // Build tables for given particle
693  all_tables_are_built = true;
694 
695  for(G4int i=0; i<n_loss; ++i) {
696  if(p == loss_vector[i] && !tables_are_built[i] && !base_part_vector[i]) {
697  const G4ParticleDefinition* curr_part = part_vector[i];
698  if(1 < verbose) {
699  G4cout << "### Build Table for " << p->GetProcessName()
700  << " and " << curr_part->GetParticleName() << G4endl;
701  }
702  G4VEnergyLossProcess* curr_proc = BuildTables(curr_part);
703  if(curr_proc) { CopyTables(curr_part, curr_proc); }
704  }
705  if ( !tables_are_built[i] ) { all_tables_are_built = false; }
706  }
707  if(0 == run && p->IsIonisationProcess()) { loss_map[aParticle] = p; }
708 
709  if(1 < verbose) {
710  G4cout << "### G4LossTableManager::BuildPhysicsTable end: "
711  << "all_tables_are_built= " << all_tables_are_built
712  << G4endl;
713  }
714  if(all_tables_are_built) {
715  if(1 < verbose) {
716  G4cout << "%%%%% All dEdx and Range tables are built for master run= "
717  << run << " %%%%%" << G4endl;
718  }
719  }
720 }
721 
722 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
723 
725  G4VEnergyLossProcess* base_proc)
726 {
727  for (G4int j=0; j<n_loss; ++j) {
728 
730 
731  if (!tables_are_built[j] && part == base_part_vector[j]) {
732  tables_are_built[j] = true;
733  proc->SetDEDXTable(base_proc->IonisationTable(),fRestricted);
734  proc->SetDEDXTable(base_proc->DEDXTableForSubsec(),fSubRestricted);
735  proc->SetDEDXTable(base_proc->DEDXunRestrictedTable(),fTotal);
736  proc->SetCSDARangeTable(base_proc->CSDARangeTable());
737  proc->SetRangeTableForLoss(base_proc->RangeTableForLoss());
738  proc->SetInverseRangeTable(base_proc->InverseRangeTable());
739  proc->SetLambdaTable(base_proc->LambdaTable());
740  proc->SetSubLambdaTable(base_proc->SubLambdaTable());
741  proc->SetIonisation(base_proc->IsIonisationProcess());
742  if(proc->IsIonisationProcess()) {
743  range_vector[j] = base_proc->RangeTableForLoss();
744  inv_range_vector[j] = base_proc->InverseRangeTable();
745  loss_map[part_vector[j]] = proc;
746  }
747  if (1 < verbose) {
748  G4cout << "For " << proc->GetProcessName()
749  << " for " << part_vector[j]->GetParticleName()
750  << " base_part= " << part->GetParticleName()
751  << " tables are assigned "
752  << G4endl;
753  }
754  }
755 
756  if (theElectron == part && theElectron == proc->SecondaryParticle() ) {
757  proc->SetSecondaryRangeTable(base_proc->RangeTableForLoss());
758  }
759  }
760 }
761 
762 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
763 
765  const G4ParticleDefinition* aParticle)
766 {
767  if(1 < verbose) {
768  G4cout << "G4LossTableManager::BuildTables() for "
769  << aParticle->GetParticleName() << G4endl;
770  }
771 
772  std::vector<G4PhysicsTable*> t_list;
773  std::vector<G4VEnergyLossProcess*> loss_list;
774  std::vector<G4bool> build_flags;
775  G4VEnergyLossProcess* em = 0;
776  G4VEnergyLossProcess* p = 0;
777  G4int iem = 0;
778  G4PhysicsTable* dedx = 0;
779  G4int i;
780 
781  G4ProcessVector* pvec =
782  aParticle->GetProcessManager()->GetProcessList();
783  G4int nvec = pvec->size();
784 
785  for (i=0; i<n_loss; ++i) {
786  p = loss_vector[i];
787  if (p) {
788  G4bool yes = (aParticle == part_vector[i]);
789 
790  // possible case of process sharing between particle/anti-particle
791  if(!yes) {
792  G4VProcess* ptr = static_cast<G4VProcess*>(p);
793  for(G4int j=0; j<nvec; ++j) {
794  //G4cout << "j= " << j << " " << (*pvec)[j] << " " << ptr << G4endl;
795  if(ptr == (*pvec)[j]) {
796  yes = true;
797  break;
798  }
799  }
800  }
801  // process belong to this particle
802  if(yes && isActive[i]) {
803  if (p->IsIonisationProcess() || !em) {
804  em = p;
805  iem= i;
806  }
807  // tables may be shared between particle/anti-particle
808  G4bool val = false;
809  if (!tables_are_built[i]) {
810  val = true;
811  dedx = p->BuildDEDXTable(fRestricted);
812  //G4cout << "Build DEDX table for " << p->GetProcessName()
813  // << " idx= " << i << dedx << " " << dedx->length() << G4endl;
814  p->SetDEDXTable(dedx,fRestricted);
815  tables_are_built[i] = true;
816  } else {
817  dedx = p->DEDXTable();
818  }
819  t_list.push_back(dedx);
820  loss_list.push_back(p);
821  build_flags.push_back(val);
822  }
823  }
824  }
825 
826  G4int n_dedx = t_list.size();
827  if (0 == n_dedx || !em) {
828  G4cout << "G4LossTableManager WARNING: no DEDX processes for "
829  << aParticle->GetParticleName() << G4endl;
830  return 0;
831  }
832  G4int nSubRegions = em->NumberOfSubCutoffRegions();
833 
834  if (1 < verbose) {
835  G4cout << "G4LossTableManager::BuildTables() start to build range tables"
836  << " and the sum of " << n_dedx << " processes"
837  << " iem= " << iem << " em= " << em->GetProcessName()
838  << " buildCSDARange= " << theParameters->BuildCSDARange()
839  << " nSubRegions= " << nSubRegions;
840  if(subcutProducer) {
841  G4cout << " SubCutProducer " << subcutProducer->GetName();
842  }
843  G4cout << G4endl;
844  }
845  // do not build tables if producer class is defined
846  if(subcutProducer) { nSubRegions = 0; }
847 
848  dedx = em->DEDXTable();
849  em->SetIonisation(true);
850  em->SetDEDXTable(dedx, fIsIonisation);
851 
852  if (1 < n_dedx) {
853  dedx = 0;
855  tableBuilder->BuildDEDXTable(dedx, t_list);
856  em->SetDEDXTable(dedx, fRestricted);
857  }
858 
859  /*
860  if(2==run && "e-" == aParticle->GetParticleName()) {
861  G4cout << "G4LossTableManager::BuildTables for e- " << dedx << G4endl;
862  G4cout << (*dedx) << G4endl;
863  G4cout << "%%%%% Instance ID= " << (*dedx)[0]->GetInstanceID() << G4endl;
864  G4cout << "%%%%% LastValue= " << (*dedx)[0]->GetLastValue() << G4endl;
865  G4cout << "%%%%% 1.2 " << (*(dedx))[0]->Value(1.2) << G4endl;
866  }
867  */
868  dedx_vector[iem] = dedx;
869 
870  G4PhysicsTable* range = em->RangeTableForLoss();
871  if(!range) range = G4PhysicsTableHelper::PreparePhysicsTable(range);
872  range_vector[iem] = range;
873 
874  G4PhysicsTable* invrange = em->InverseRangeTable();
875  if(!invrange) invrange = G4PhysicsTableHelper::PreparePhysicsTable(invrange);
876  inv_range_vector[iem] = invrange;
877 
878  tableBuilder->BuildRangeTable(dedx, range, true);
879  tableBuilder->BuildInverseRangeTable(range, invrange, true);
880 
881  // if(1<verbose) G4cout << *dedx << G4endl;
882 
883  em->SetRangeTableForLoss(range);
884  em->SetInverseRangeTable(invrange);
885 
886  // if(1<verbose) G4cout << *range << G4endl;
887 
888  std::vector<G4PhysicsTable*> listSub;
889  std::vector<G4PhysicsTable*> listCSDA;
890 
891  for (i=0; i<n_dedx; ++i) {
892  p = loss_list[i];
893  if(p != em) { p->SetIonisation(false); }
894  if(build_flags[i]) {
896  }
897  if (0 < nSubRegions) {
898  dedx = p->BuildDEDXTable(fSubRestricted);
899  p->SetDEDXTable(dedx,fSubRestricted);
900  listSub.push_back(dedx);
901  if(build_flags[i]) {
903  if(p != em) { em->AddCollaborativeProcess(p); }
904  }
905  }
906  if(theParameters->BuildCSDARange()) {
907  dedx = p->BuildDEDXTable(fTotal);
908  p->SetDEDXTable(dedx,fTotal);
909  listCSDA.push_back(dedx);
910  }
911  }
912 
913  if (0 < nSubRegions) {
914  G4PhysicsTable* dedxSub = em->IonisationTableForSubsec();
915  if (1 < listSub.size()) {
916  em->SetDEDXTable(dedxSub, fIsSubIonisation);
917  dedxSub = 0;
919  tableBuilder->BuildDEDXTable(dedxSub, listSub);
920  em->SetDEDXTable(dedxSub, fSubRestricted);
921  }
922  }
924  G4PhysicsTable* dedxCSDA = em->DEDXunRestrictedTable();
925  if (1 < n_dedx) {
926  dedxCSDA = 0;
927  dedxCSDA = G4PhysicsTableHelper::PreparePhysicsTable(dedxCSDA);
928  tableBuilder->BuildDEDXTable(dedxCSDA, listCSDA);
929  em->SetDEDXTable(dedxCSDA,fTotal);
930  }
931  G4PhysicsTable* rCSDA = em->CSDARangeTable();
932  if(!rCSDA) { rCSDA = G4PhysicsTableHelper::PreparePhysicsTable(rCSDA); }
933  tableBuilder->BuildRangeTable(dedxCSDA, rCSDA, true);
934  em->SetCSDARangeTable(rCSDA);
935  }
936 
937  if (1 < verbose) {
938  G4cout << "G4LossTableManager::BuildTables: Tables are built for "
939  << aParticle->GetParticleName()
940  << "; ionisation process: " << em->GetProcessName()
941  << " " << em
942  << G4endl;
943  }
944  return em;
945 }
946 
947 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
948 
950 {
951  return theMessenger;
952 }
953 
954 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
955 
957  const G4ParticleDefinition* aParticle)
958 {
960  ed << "Energy loss process not found for " << aParticle->GetParticleName()
961  << " !";
962  G4Exception("G4LossTableManager::ParticleHaveNoLoss", "em0001",
963  FatalException, ed);
964 }
965 
966 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
967 
969 {
970  return theParameters->BuildCSDARange();
971 }
972 
973 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
974 
976 {
977 }
978 
979 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
980 
982 {
983  subCutoffFlag = val;
984  for(G4int i=0; i<n_loss; ++i) {
985  if(loss_vector[i]) { loss_vector[i]->ActivateSubCutoff(val, r); }
986  }
987 }
988 
989 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
990 
992 {
993  integral = val;
994  integralActive = true;
995  for(G4int i=0; i<n_loss; ++i) {
996  if(loss_vector[i]) { loss_vector[i]->SetIntegral(val); }
997  }
998  size_t emp = emp_vector.size();
999  for (size_t k=0; k<emp; ++k) {
1000  if(emp_vector[k]) { emp_vector[k]->SetIntegral(val); }
1001  }
1002 }
1003 
1004 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1005 
1007 {
1008 }
1009 
1010 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1011 
1013 {
1014 }
1015 
1016 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1017 
1019 {
1021 }
1022 
1023 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1024 
1026 {
1028 }
1029 
1030 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1031 
1033 {
1035 }
1036 
1037 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1038 
1040 {
1042 }
1043 
1044 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1045 
1047 {
1049 }
1050 
1051 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1052 
1054 {
1055 }
1056 
1057 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1058 
1060 {
1062 }
1063 
1064 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1065 
1067 {
1069 }
1070 
1071 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1072 
1074 {
1075 }
1076 
1077 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1078 
1080 {
1081  if(0.0 < v1 && 0.0 < v2 && v2 < 1.e+50) {
1082  stepFunctionActive = true;
1083  maxRangeVariation = v1;
1084  maxFinalStep = v2;
1085  for(G4int i=0; i<n_loss; ++i) {
1086  if(loss_vector[i]) { loss_vector[i]->SetStepFunction(v1, v2); }
1087  }
1088  } else if(v1 <= 0.0) {
1089  PrintEWarning("SetStepFunction", v1);
1090  } else {
1091  PrintEWarning("SetStepFunction", v2);
1092  }
1093 }
1094 
1095 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1096 
1098 {
1099 }
1100 
1101 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1102 
1104 {
1105 }
1106 
1107 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1108 
1109 void
1112 {
1113  if(stepFunctionActive) {
1115  }
1116  if(integralActive) { p->SetIntegral(integral); }
1117 }
1118 
1119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1120 
1121 const std::vector<G4VEnergyLossProcess*>&
1123 {
1124  return loss_vector;
1125 }
1126 
1127 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1128 
1129 const std::vector<G4VEmProcess*>& G4LossTableManager::GetEmProcessVector()
1130 {
1131  return emp_vector;
1132 }
1133 
1134 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1135 
1136 const std::vector<G4VMultipleScattering*>&
1138 {
1139  return msc_vector;
1140 }
1141 
1142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1143 
1145 {
1146  return isMaster;
1147 }
1148 
1149 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1150 
1152 {
1153  return theParameters->MinKinEnergy();
1154 }
1155 
1156 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1157 
1159 {
1160  return theParameters->MaxKinEnergy();
1161 }
1162 
1163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1164 
1166 {
1167  return emCorrections;
1168 }
1169 
1170 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1171 
1173 {
1175  return emSaturation;
1176 }
1177 
1178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1179 
1181 {
1183  return emConfigurator;
1184 }
1185 
1186 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1187 
1189 {
1190  if(!emElectronIonPair) {
1192  }
1193  return emElectronIonPair;
1194 }
1195 
1196 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1197 
1199 {
1200  return atomDeexcitation;
1201 }
1202 
1203 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1204 
1206 {
1207  return tableBuilder;
1208 }
1209 
1210 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1211 
1213 {
1214  if(atomDeexcitation != p) {
1215  delete atomDeexcitation;
1216  atomDeexcitation = p;
1217  }
1218 }
1219 
1220 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1221 
1223 {
1224  if(subcutProducer != p) {
1225  delete subcutProducer;
1226  subcutProducer = p;
1227  }
1228 }
1229 
1230 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1231 
1233 {
1234  return subcutProducer;
1235 }
1236 
1237 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1238 
1241 {
1242  //G4cout << aParticle << " " << currentParticle
1243  //<< " " << currentLoss << G4endl;
1244  if(aParticle != currentParticle) {
1245  currentParticle = aParticle;
1246  std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
1247  if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
1248  currentLoss = (*pos).second;
1249  } else {
1250  currentLoss = 0;
1251  if(aParticle->GetParticleType() == "nucleus") {
1252  if ((pos = loss_map.find(theGenericIon)) != loss_map.end()) {
1253  currentLoss = (*pos).second;
1254  }
1255  }
1256  }
1257  }
1258  return currentLoss;
1259 }
1260 
1261 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1262 
1264  G4double kineticEnergy,
1265  const G4MaterialCutsCouple *couple)
1266 {
1267  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1268  G4double x = 0.0;
1269  if(currentLoss) { x = currentLoss->GetDEDX(kineticEnergy, couple); }
1270  return x;
1271 }
1272 
1273 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1274 
1276  G4double kineticEnergy,
1277  const G4MaterialCutsCouple *couple)
1278 {
1279  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1280  G4double x = 0.0;
1281  if(currentLoss) { x = currentLoss->GetDEDXForSubsec(kineticEnergy, couple); }
1282  return x;
1283 }
1284 
1285 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1286 
1288  G4double kineticEnergy,
1289  const G4MaterialCutsCouple *couple)
1290 {
1291  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1292  G4double x = DBL_MAX;
1293  if(currentLoss) { x = currentLoss->GetCSDARange(kineticEnergy, couple); }
1294  return x;
1295 }
1296 
1297 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1298 
1300  const G4ParticleDefinition *aParticle,
1301  G4double kineticEnergy,
1302  const G4MaterialCutsCouple *couple)
1303 {
1304  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1305  G4double x = DBL_MAX;
1306  if(currentLoss) { x = currentLoss->GetRangeForLoss(kineticEnergy, couple); }
1307  return x;
1308 }
1309 
1310 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1311 
1313  G4double kineticEnergy,
1314  const G4MaterialCutsCouple *couple)
1315 {
1316  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1317  G4double x = DBL_MAX;
1318  if(currentLoss) {
1319  x = currentLoss->GetRange(kineticEnergy, couple);
1320  }
1321  return x;
1322 }
1323 
1324 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1325 
1327  G4double range,
1328  const G4MaterialCutsCouple *couple)
1329 {
1330  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1331  G4double x = 0;
1332  if(currentLoss) { x = currentLoss->GetKineticEnergy(range, couple); }
1333  return x;
1334 }
1335 
1336 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1337 
1339  const G4MaterialCutsCouple *couple,
1340  const G4DynamicParticle* dp,
1341  G4double& length)
1342 {
1343  const G4ParticleDefinition* aParticle = dp->GetParticleDefinition();
1344  if(aParticle != currentParticle) { GetEnergyLossProcess(aParticle); }
1345  G4double x = 0.0;
1346  if(currentLoss) { currentLoss->GetDEDXDispersion(couple, dp, length); }
1347  return x;
1348 }
1349 
1350 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1351 
1353 {
1354  G4String ss = "G4LossTableManager::" + tit;
1356  /*
1357  ed << "Parameter is out of range: " << val
1358  << " it will have no effect!\n" << " ## "
1359  << " nbins= " << nbinsLambda
1360  << " nbinsPerDecade= " << nbinsPerDecade
1361  << " Emin(keV)= " << minKinEnergy/keV
1362  << " Emax(GeV)= " << maxKinEnergy/GeV;
1363  */
1364  G4Exception(ss, "em0044", JustWarning, ed);
1365 }
1366 
1367 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4double GetKineticEnergy(G4double &range, const G4MaterialCutsCouple *)
G4VSubCutProducer * SubCutProducer()
G4EmConfigurator * EmConfigurator()
G4int NumberOfBinsPerDecade() const
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 Spline() const
G4int WorkerVerbose() const
G4double MaxKinEnergy() const
static G4LossTableManager * Instance()
G4VAtomDeexcitation * atomDeexcitation
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void SetVerbose(G4int verb)
G4double GetDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
void SetIonisation(G4bool val)
G4VSubCutProducer * subcutProducer
G4double GetDEDX(G4double &kineticEnergy, const G4MaterialCutsCouple *)
std::vector< G4VEmModel * > mod_vector
G4PhysicsTable * SubLambdaTable() const
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double &length)
G4ElectronIonPair * emElectronIonPair
void DeRegister(G4VEnergyLossProcess *p)
G4PhysicsTable * RangeTableForLoss() const
void SetLambdaBinning(G4int val)
G4EmSaturation * emSaturation
G4double GetSubDEDX(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4bool BuildCSDARange() const
void push_back(G4PhysicsVector *)
static G4ThreadLocal G4LossTableManager * instance
G4double GetRangeForLoss(G4double &kineticEnergy, const G4MaterialCutsCouple *)
void AddCollaborativeProcess(G4VEnergyLossProcess *)
void SetStepFunction(G4double v1, G4double v2)
G4PhysicsTable * CSDARangeTable() const
const G4String & GetName() const
G4double GetCSDARange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
void SetNumberOfBins(G4int val)
G4EnergyLossMessenger * GetMessenger()
void CopyTables(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *)
void SetMaxEnergyForCSDARange(G4double val)
std::vector< G4PhysicsTable * > range_vector
G4PhysicsTable * IonisationTableForSubsec() const
G4double a
Definition: TRTMaterials.hh:39
G4double GetDEDXForSubsec(G4double &kineticEnergy, const G4MaterialCutsCouple *)
#define G4ThreadLocal
Definition: tls.hh:84
G4bool Fluo() const
G4PhysicsTable * IonisationTable() const
G4ProcessManager * GetProcessManager() const
const std::vector< G4VEmProcess * > & GetEmProcessVector()
int G4int
Definition: G4Types.hh:78
G4PhysicsTable * BuildLambdaTable(G4EmTableType tType=fRestricted)
void SetMaxEnergy(G4double val)
void SetVerbose(G4int)
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
const G4String & GetParticleName() const
G4bool DeexcitationIgnoreCut() const
G4double GetCSDARange(const G4ParticleDefinition *aParticle, G4double kineticEnergy, const G4MaterialCutsCouple *couple)
G4LossTableBuilder * GetTableBuilder()
std::map< PD, G4VEnergyLossProcess *, std::less< PD > > loss_map
G4EmConfigurator * emConfigurator
void SetBuildCSDARange(G4bool val)
void SetInitialisationFlag(G4bool flag)
G4VEnergyLossProcess * BuildTables(const G4ParticleDefinition *aParticle)
G4bool BuildCSDARange() const
G4PhysicsTable * LambdaTable() const
void ParticleHaveNoLoss(const G4ParticleDefinition *aParticle)
void SetInverseRangeTable(G4PhysicsTable *p)
const G4ParticleDefinition * SecondaryParticle() const
void SetSubCutProducer(G4VSubCutProducer *)
void SetParameters(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *)
G4GLOB_DLL std::ostream G4cout
G4int Verbose() const
G4PhysicsTable * DEDXTable() const
G4double MinKinEnergy() const
std::vector< G4VEmFluctuationModel * > fmod_vector
bool G4bool
Definition: G4Types.hh:79
G4EmCorrections * EmCorrections()
G4VEnergyLossProcess * currentLoss
G4ElectronIonPair * ElectronIonPair()
G4EmSaturation * EmSaturation()
void SetMaxEnergyForCSDARange(G4double val)
G4bool IsMaster() const
G4int NumberOfSubCutoffRegions() const
const G4String & GetParticleType() const
void Register(G4VEnergyLossProcess *p)
const G4ParticleDefinition * BaseParticle() const
G4double MinKinEnergy() 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 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
std::vector< G4PhysicsTable * > inv_range_vector
G4bool IsWorkerThread()
Definition: G4Threading.cc:128
void BuildDEDXTable(G4PhysicsTable *dedxTable, const std::vector< G4PhysicsTable * > &)
void SetVerbose(G4int value)
void SetMinSubRange(G4double val)
G4EmParameters * theParameters
G4int size() const
G4EmCorrections * emCorrections
void SetMinEnergy(G4double val)
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
void SetLambdaTable(G4PhysicsTable *p)
void SetMaxEnergy(G4double val)
const G4String & GetName() const
void SetStepFunction(G4double v1, G4double v2)
G4bool Auger() const
G4PhysicsTable * InverseRangeTable() const
void DumpG4BirksCoefficients()
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
G4LossTableBuilder * tableBuilder
const G4ParticleDefinition * Particle() const
void SetLinearLossLimit(G4double val)
G4int GetNumberOfBinsPerDecade() const
void SetMinEnergy(G4double val)
std::vector< G4PhysicsTable * > dedx_vector
void PrepareModels(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
static G4EmParameters * Instance()
std::vector< G4VEnergyLossProcess * > loss_vector
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
void SetCSDARangeTable(G4PhysicsTable *pRange)
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
void PrintEWarning(G4String, G4double)
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
std::vector< PD > part_vector
#define G4endl
Definition: G4ios.hh:61
void SetIntegral(G4bool val)
void SetLossFluctuations(G4bool val)
const G4String & GetName() const
Definition: G4VEmModel.hh:776
std::vector< G4bool > isActive
void SetSecondaryRangeTable(G4PhysicsTable *p)
G4bool Pixe() const
G4VAtomDeexcitation * AtomDeexcitation()
std::vector< G4VEmProcess * > emp_vector
G4PhysicsTable * DEDXunRestrictedTable() const
double G4double
Definition: G4Types.hh:76
void SetSubLambdaTable(G4PhysicsTable *p)
std::vector< PD > base_part_vector
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)
std::vector< G4VMultipleScattering * > msc_vector
static const double mm
Definition: G4SIunits.hh:102
void SetAtomDeexcitation(G4VAtomDeexcitation *)
G4EnergyLossMessenger * theMessenger
static const G4double pos
G4PhysicsTable * BuildDEDXTable(G4EmTableType tType=fRestricted)
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()
G4ProcessVector * GetProcessList() const
G4double GetRange(G4double &kineticEnergy, const G4MaterialCutsCouple *)
std::vector< G4bool > tables_are_built
void SetSplineFlag(G4bool flag)
G4bool IsIonisationProcess() const