Geant4  10.02.p02
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 95476 2016-02-12 09:39:50Z 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 
88 #include "G4VMultipleScattering.hh"
89 #include "G4VEmProcess.hh"
90 
91 #include "G4EmSaturation.hh"
92 #include "G4EmConfigurator.hh"
93 #include "G4ElectronIonPair.hh"
94 
95 #include "G4PhysicsTable.hh"
96 #include "G4ParticleDefinition.hh"
97 #include "G4MaterialCutsCouple.hh"
98 #include "G4ProcessManager.hh"
99 #include "G4Electron.hh"
100 #include "G4Proton.hh"
101 #include "G4ProductionCutsTable.hh"
102 #include "G4PhysicsTableHelper.hh"
103 #include "G4EmTableType.hh"
104 #include "G4Region.hh"
105 #include "G4PhysicalConstants.hh"
106 #include "G4Threading.hh"
107 
108 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
109 
111 
113 {
114  if(!instance) {
116  instance = inst.Instance();
117  }
118  return instance;
119 }
120 
121 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
122 
124 {
125  //G4cout << "### G4LossTableManager::~G4LossTableManager()" << G4endl;
126  for (G4int i=0; i<n_loss; ++i) {
127  //G4cout << "### eloss #" << i << G4endl;
128  if( loss_vector[i] ) {
129  delete loss_vector[i];
130  }
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  delete subcutProducer;
164 }
165 
166 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
167 
169 {
171  n_loss = 0;
172  run = -1;
173  startInitialisation = false;
174  all_tables_are_built = false;
175  currentLoss = nullptr;
176  currentParticle = nullptr;
177  firstParticle = nullptr;
178  subCutoffFlag = false;
179  maxRangeVariation = 0.2;
181  integral = true;
182  integralActive = false;
183  stepFunctionActive = false;
184  isMaster = true;
188  theGenericIon= nullptr;
191  isMaster = false;
192  }
195  emSaturation = nullptr;
196  emConfigurator = nullptr;
197  emElectronIonPair = nullptr;
198  atomDeexcitation = nullptr;
199  subcutProducer = nullptr;
200 }
201 
202 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
203 
205 {
206  all_tables_are_built = false;
207  currentLoss = nullptr;
208  currentParticle = nullptr;
209  if(n_loss)
210  {
211  dedx_vector.clear();
212  range_vector.clear();
213  inv_range_vector.clear();
214  loss_map.clear();
215  loss_vector.clear();
216  part_vector.clear();
217  base_part_vector.clear();
218  tables_are_built.clear();
219  isActive.clear();
220  n_loss = 0;
221  }
222 }
223 
224 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
225 
227 {
228  if(!p) { return; }
229  for (G4int i=0; i<n_loss; ++i) {
230  if(loss_vector[i] == p) { return; }
231  }
232  if(verbose > 1) {
233  G4cout << "G4LossTableManager::Register G4VEnergyLossProcess : "
234  << p->GetProcessName() << " idx= " << n_loss << G4endl;
235  }
236  ++n_loss;
237  loss_vector.push_back(p);
238  part_vector.push_back(nullptr);
239  base_part_vector.push_back(nullptr);
240  dedx_vector.push_back(nullptr);
241  range_vector.push_back(nullptr);
242  inv_range_vector.push_back(nullptr);
243  tables_are_built.push_back(false);
244  isActive.push_back(true);
245  all_tables_are_built = false;
246  if(subCutoffFlag) { p->ActivateSubCutoff(true); }
248  maxFinalStep); }
249  if(integralActive) { p->SetIntegral(integral); }
250 }
251 
252 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
253 
255 {
257  if(!isMaster) {
259  }
266  if(atomDeexcitation) {
269  }
270 }
271 
272 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
273 
275 {
276  if(!p) { return; }
277  for (G4int i=0; i<n_loss; ++i) {
278  if(loss_vector[i] == p) { loss_vector[i] = nullptr; }
279  }
280 }
281 
282 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
283 
285 {
286  if(!p) { return; }
287  G4int n = msc_vector.size();
288  for (G4int i=0; i<n; ++i) {
289  if(msc_vector[i] == p) { return; }
290  }
291  if(verbose > 1) {
292  G4cout << "G4LossTableManager::Register G4VMultipleScattering : "
293  << p->GetProcessName() << " idx= " << msc_vector.size() << G4endl;
294  }
295  msc_vector.push_back(p);
296 }
297 
298 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
299 
301 {
302  if(!p) { return; }
303  size_t msc = msc_vector.size();
304  for (size_t i=0; i<msc; ++i) {
305  if(msc_vector[i] == p) { msc_vector[i] = nullptr; }
306  }
307 }
308 
309 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
310 
312 {
313  if(!p) { return; }
314  G4int n = emp_vector.size();
315  for (G4int i=0; i<n; ++i) {
316  if(emp_vector[i] == p) { return; }
317  }
318  if(verbose > 1) {
319  G4cout << "G4LossTableManager::Register G4VEmProcess : "
320  << p->GetProcessName() << " idx= " << emp_vector.size() << G4endl;
321  }
322  emp_vector.push_back(p);
323 }
324 
325 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
326 
328 {
329  if(!p) { return; }
330  size_t emp = emp_vector.size();
331  for (size_t i=0; i<emp; ++i) {
332  if(emp_vector[i] == p) { emp_vector[i] = nullptr; }
333  }
334 }
335 
336 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
337 
339 {
340  mod_vector.push_back(p);
341  if(verbose > 1) {
342  G4cout << "G4LossTableManager::Register G4VEmModel : "
343  << p->GetName() << G4endl;
344  }
345 }
346 
347 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
348 
350 {
351  size_t n = mod_vector.size();
352  for (size_t i=0; i<n; ++i) {
353  if(mod_vector[i] == p) { mod_vector[i] = nullptr; }
354  }
355 }
356 
357 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
358 
360 {
361  fmod_vector.push_back(p);
362  if(verbose > 1) {
363  G4cout << "G4LossTableManager::Register G4VEmFluctuationModel : "
364  << p->GetName() << G4endl;
365  }
366 }
367 
368 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
369 
371 {
372  size_t n = fmod_vector.size();
373  for (size_t i=0; i<n; ++i) {
374  if(fmod_vector[i] == p) { fmod_vector[i] = nullptr; }
375  }
376 }
377 
378 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
379 
381  const G4ParticleDefinition* part,
383 {
384  if(!p || !part) { return; }
385  for (G4int i=0; i<n_loss; ++i) {
386  if(loss_vector[i] == p) { return; }
387  }
388  if(verbose > 1) {
389  G4cout << "G4LossTableManager::RegisterExtraParticle "
390  << part->GetParticleName() << " G4VEnergyLossProcess : "
391  << p->GetProcessName() << " idx= " << n_loss << G4endl;
392  }
393  ++n_loss;
394  loss_vector.push_back(p);
395  part_vector.push_back(part);
396  base_part_vector.push_back(p->BaseParticle());
397  dedx_vector.push_back(nullptr);
398  range_vector.push_back(nullptr);
399  inv_range_vector.push_back(nullptr);
400  tables_are_built.push_back(false);
401  all_tables_are_built = false;
402 }
403 
404 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
405 
406 void
409  G4bool theMaster)
410 {
411  if (1 < verbose) {
412  G4cout << "G4LossTableManager::PreparePhysicsTable for "
413  << particle->GetParticleName()
414  << " and " << p->GetProcessName() << " run= " << run
415  << " loss_vector " << loss_vector.size() << G4endl;
416  }
417 
418  isMaster = theMaster;
419 
420  if(!startInitialisation) {
421  ResetParameters();
422  if (1 < verbose) {
423  G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
424  << G4endl;
425  }
426  }
427 
428  // start initialisation for the first run
429  if( -1 == run ) {
430  if(emConfigurator) { emConfigurator->PrepareModels(particle, p); }
431 
432  // initialise particles for given process
433  for (G4int j=0; j<n_loss; ++j) {
434  if (p == loss_vector[j] && !part_vector[j]) {
435  part_vector[j] = particle;
436  if(particle->GetParticleName() == "GenericIon") {
437  theGenericIon = particle;
438  }
439  }
440  }
441  }
442  startInitialisation = true;
443 }
444 
445 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
446 
447 void
449  G4VEmProcess* p, G4bool theMaster)
450 {
451  if (1 < verbose) {
452  G4cout << "G4LossTableManager::PreparePhysicsTable for "
453  << particle->GetParticleName()
454  << " and " << p->GetProcessName() << G4endl;
455  }
456  isMaster = theMaster;
457 
458  if(!startInitialisation) {
459  ResetParameters();
460  if (1 < verbose) {
461  G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
462  << G4endl;
463  }
464  }
465 
466  // start initialisation for the first run
467  if( -1 == run ) {
468  if(emConfigurator) { emConfigurator->PrepareModels(particle, p); }
469  }
470  startInitialisation = true;
471 }
472 
473 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
474 
475 void
478  G4bool theMaster)
479 {
480  if (1 < verbose) {
481  G4cout << "G4LossTableManager::PreparePhysicsTable for "
482  << particle->GetParticleName()
483  << " and " << p->GetProcessName() << G4endl;
484  }
485 
486  isMaster = theMaster;
487 
488  if(!startInitialisation) {
489  ResetParameters();
490  if (1 < verbose) {
491  G4cout << "====== G4LossTableManager::PreparePhysicsTable start ====="
492  << G4endl;
493  }
494  }
495 
496  // start initialisation for the first run
497  if( -1 == run ) {
498  if(emConfigurator) { emConfigurator->PrepareModels(particle, p); }
499  }
500  startInitialisation = true;
501 }
502 
503 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
504 
505 void
507 {
508  if(-1 == run && startInitialisation) {
510  }
511 }
512 
513 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
514 
516  const G4ParticleDefinition* aParticle,
518 {
519  if(1 < verbose) {
520  G4cout << "### G4LossTableManager::LocalPhysicsTable() for "
521  << aParticle->GetParticleName()
522  << " and process " << p->GetProcessName()
523  << G4endl;
524  }
525 
526  if(-1 == run && startInitialisation) {
528  firstParticle = aParticle;
529  }
530 
531  if(startInitialisation) {
532  ++run;
534  if(1 < verbose) {
535  G4cout << "===== G4LossTableManager::LocalPhysicsTable() for run "
536  << run << " =====" << G4endl;
537  }
539  currentParticle = nullptr;
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] = nullptr;
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  part_vector[i] = p->Particle();
557  base_part_vector[i] = p->BaseParticle();
558  dedx_vector[i] = p->DEDXTable();
561  if(0 == run && p->IsIonisationProcess()) {
562  loss_map[part_vector[i]] = p;
563  //G4cout << "G4LossTableManager::LocalPhysicsTable " << part_vector[i]->GetParticleName()
564  // << " added to map " << p << G4endl;
565  }
566 
567  if(1 < verbose) {
568  G4cout << i <<". "<< p->GetProcessName();
569  if(part_vector[i]) {
570  G4cout << " for " << part_vector[i]->GetParticleName();
571  }
572  G4cout << " active= " << isActive[i]
573  << " table= " << tables_are_built[i]
574  << " isIonisation= " << p->IsIonisationProcess()
575  << G4endl;
576  }
577  break;
578  } else if(!tables_are_built[i]) {
579  all_tables_are_built = false;
580  }
581  }
582 
583  // Set run time parameters
584  SetParameters(aParticle, p);
585 
586  if(1 < verbose) {
587  G4cout << "### G4LossTableManager::LocalPhysicsTable end"
588  << G4endl;
589  }
590  if(all_tables_are_built) {
591  if(1 < verbose) {
592  G4cout << "%%%%% All dEdx and Range tables for worker are ready for run "
593  << run << " %%%%%" << G4endl;
594  }
595  }
596 }
597 
598 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
599 
601  const G4ParticleDefinition* aParticle,
603 {
604  if(1 < verbose) {
605  G4cout << "### G4LossTableManager::BuildPhysicsTable() for "
606  << aParticle->GetParticleName()
607  << " and process " << p->GetProcessName() << G4endl;
608  }
609  // clear configurator
610  if(-1 == run && startInitialisation) {
612  firstParticle = aParticle;
613  }
614  if(startInitialisation) {
615  ++run;
616  if(1 < verbose) {
617  G4cout << "===== G4LossTableManager::BuildPhysicsTable() for run "
618  << run << " ===== " << atomDeexcitation << G4endl;
619  }
620  currentParticle = nullptr;
621  all_tables_are_built= true;
622  }
623 
624  // initialisation before any table is built
625  if ( startInitialisation && aParticle == firstParticle ) {
626 
627  startInitialisation = false;
628  if(1 < verbose) {
629  G4cout << "### G4LossTableManager start initilisation for first particle "
631  << G4endl;
632  }
633  for (G4int i=0; i<n_loss; ++i) {
635 
636  if(el) {
637  isActive[i] = true;
638  base_part_vector[i] = el->BaseParticle();
639  tables_are_built[i] = false;
640  all_tables_are_built= false;
641  if(!isActive[i]) {
642  el->SetIonisation(false);
643  tables_are_built[i] = true;
644  }
645 
646  if(1 < verbose) {
647  G4cout << i <<". "<< el->GetProcessName();
648  if(el->Particle()) {
649  G4cout << " for " << el->Particle()->GetParticleName();
650  }
651  G4cout << " active= " << isActive[i]
652  << " table= " << tables_are_built[i]
653  << " isIonisation= " << el->IsIonisationProcess();
654  if(base_part_vector[i]) {
655  G4cout << " base particle "
656  << base_part_vector[i]->GetParticleName();
657  }
658  G4cout << G4endl;
659  }
660  } else {
661  tables_are_built[i] = true;
662  part_vector[i] = nullptr;
663  isActive[i] = false;
664  }
665  }
666  }
667 
668  // Set run time parameters
669  SetParameters(aParticle, p);
670 
671  if (all_tables_are_built) { return; }
672 
673  // Build tables for given particle
674  all_tables_are_built = true;
675 
676  for(G4int i=0; i<n_loss; ++i) {
677  if(p == loss_vector[i] && !tables_are_built[i] && !base_part_vector[i]) {
678  const G4ParticleDefinition* curr_part = part_vector[i];
679  if(1 < verbose) {
680  G4cout << "### Build Table for " << p->GetProcessName()
681  << " and " << curr_part->GetParticleName()
682  << " " << tables_are_built[i] << " " << base_part_vector[i] << G4endl;
683  }
684  G4VEnergyLossProcess* curr_proc = BuildTables(curr_part);
685  if(curr_proc) {
686  CopyTables(curr_part, curr_proc);
687  if(p == curr_proc && 0 == run && p->IsIonisationProcess()) {
688  loss_map[aParticle] = p;
689  //G4cout << "G4LossTableManager::BuildPhysicsTable: " << aParticle->GetParticleName()
690  // << " added to map " << p << G4endl;
691  }
692  }
693  }
694  if ( !tables_are_built[i] ) { all_tables_are_built = false; }
695  }
696  if(1 < verbose) {
697  G4cout << "### G4LossTableManager::BuildPhysicsTable end: "
698  << "all_tables_are_built= " << all_tables_are_built << " "
699  << aParticle->GetParticleName() << " proc: " << p << G4endl;
700  }
701  if(all_tables_are_built) {
702  if(1 < verbose) {
703  G4cout << "%%%%% All dEdx and Range tables are built for master run= "
704  << run << " %%%%%" << G4endl;
705  }
706  }
707 }
708 
709 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
710 
712  G4VEnergyLossProcess* base_proc)
713 {
714  for (G4int j=0; j<n_loss; ++j) {
715 
717 
718  if (!tables_are_built[j] && part == base_part_vector[j]) {
719  tables_are_built[j] = true;
720  proc->SetDEDXTable(base_proc->IonisationTable(),fRestricted);
721  proc->SetDEDXTable(base_proc->DEDXTableForSubsec(),fSubRestricted);
722  proc->SetDEDXTable(base_proc->DEDXunRestrictedTable(),fTotal);
723  proc->SetCSDARangeTable(base_proc->CSDARangeTable());
724  proc->SetRangeTableForLoss(base_proc->RangeTableForLoss());
725  proc->SetInverseRangeTable(base_proc->InverseRangeTable());
726  proc->SetLambdaTable(base_proc->LambdaTable());
727  proc->SetSubLambdaTable(base_proc->SubLambdaTable());
728  proc->SetIonisation(base_proc->IsIonisationProcess());
729  if(proc->IsIonisationProcess()) {
730  range_vector[j] = base_proc->RangeTableForLoss();
731  inv_range_vector[j] = base_proc->InverseRangeTable();
732  loss_map[part_vector[j]] = proc;
733  //G4cout << "G4LossTableManager::CopyTable " << part_vector[j]->GetParticleName()
734  // << " added to map " << proc << G4endl;
735  }
736  if (1 < verbose) {
737  G4cout << "For " << proc->GetProcessName()
738  << " for " << part_vector[j]->GetParticleName()
739  << " base_part= " << part->GetParticleName()
740  << " tables are assigned"
741  << G4endl;
742  }
743  }
744 
745  if (theElectron == part && theElectron == proc->SecondaryParticle() ) {
746  proc->SetSecondaryRangeTable(base_proc->RangeTableForLoss());
747  }
748  }
749 }
750 
751 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
752 
754  const G4ParticleDefinition* aParticle)
755 {
756  if(1 < verbose) {
757  G4cout << "G4LossTableManager::BuildTables() for "
758  << aParticle->GetParticleName() << G4endl;
759  }
760 
761  std::vector<G4PhysicsTable*> t_list;
762  std::vector<G4VEnergyLossProcess*> loss_list;
763  std::vector<G4bool> build_flags;
764  G4VEnergyLossProcess* em = nullptr;
765  G4VEnergyLossProcess* p = nullptr;
766  G4int iem = 0;
767  G4PhysicsTable* dedx = 0;
768  G4int i;
769 
770  G4ProcessVector* pvec =
771  aParticle->GetProcessManager()->GetProcessList();
772  G4int nvec = pvec->size();
773 
774  for (i=0; i<n_loss; ++i) {
775  p = loss_vector[i];
776  if (p) {
777  G4bool yes = (aParticle == part_vector[i]);
778 
779  // possible case of process sharing between particle/anti-particle
780  if(!yes) {
781  G4VProcess* ptr = static_cast<G4VProcess*>(p);
782  for(G4int j=0; j<nvec; ++j) {
783  //G4cout << "j= " << j << " " << (*pvec)[j] << " " << ptr << G4endl;
784  if(ptr == (*pvec)[j]) {
785  yes = true;
786  break;
787  }
788  }
789  }
790  // process belong to this particle
791  if(yes && isActive[i]) {
792  if (p->IsIonisationProcess() || !em) {
793  em = p;
794  iem= i;
795  }
796  // tables may be shared between particle/anti-particle
797  G4bool val = false;
798  if (!tables_are_built[i]) {
799  val = true;
800  dedx = p->BuildDEDXTable(fRestricted);
801  //G4cout << "Build DEDX table for " << p->GetProcessName()
802  // << " idx= " << i << dedx << " " << dedx->length() << G4endl;
803  p->SetDEDXTable(dedx,fRestricted);
804  tables_are_built[i] = true;
805  } else {
806  dedx = p->DEDXTable();
807  }
808  t_list.push_back(dedx);
809  loss_list.push_back(p);
810  build_flags.push_back(val);
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= " << theParameters->BuildCSDARange()
828  << " nSubRegions= " << nSubRegions;
829  if(subcutProducer) {
830  G4cout << " SubCutProducer " << subcutProducer->GetName();
831  }
832  G4cout << G4endl;
833  }
834  // do not build tables if producer class is defined
835  if(subcutProducer) { nSubRegions = 0; }
836 
837  dedx = em->DEDXTable();
838  em->SetIonisation(true);
839  em->SetDEDXTable(dedx, fIsIonisation);
840 
841  if (1 < n_dedx) {
842  dedx = 0;
844  tableBuilder->BuildDEDXTable(dedx, t_list);
845  em->SetDEDXTable(dedx, fRestricted);
846  }
847 
848  /*
849  if(2==run && "e-" == aParticle->GetParticleName()) {
850  G4cout << "G4LossTableManager::BuildTables for e- " << dedx << G4endl;
851  G4cout << (*dedx) << G4endl;
852  G4cout << "%%%%% Instance ID= " << (*dedx)[0]->GetInstanceID() << G4endl;
853  G4cout << "%%%%% LastValue= " << (*dedx)[0]->GetLastValue() << G4endl;
854  G4cout << "%%%%% 1.2 " << (*(dedx))[0]->Value(1.2) << G4endl;
855  }
856  */
857  dedx_vector[iem] = dedx;
858 
859  G4PhysicsTable* range = em->RangeTableForLoss();
860  if(!range) range = G4PhysicsTableHelper::PreparePhysicsTable(range);
861  range_vector[iem] = range;
862 
863  G4PhysicsTable* invrange = em->InverseRangeTable();
864  if(!invrange) invrange = G4PhysicsTableHelper::PreparePhysicsTable(invrange);
865  inv_range_vector[iem] = invrange;
866 
867  tableBuilder->BuildRangeTable(dedx, range, true);
868  tableBuilder->BuildInverseRangeTable(range, invrange, true);
869 
870  // if(1<verbose) G4cout << *dedx << G4endl;
871 
872  em->SetRangeTableForLoss(range);
873  em->SetInverseRangeTable(invrange);
874 
875  // if(1<verbose) G4cout << *range << G4endl;
876 
877  std::vector<G4PhysicsTable*> listSub;
878  std::vector<G4PhysicsTable*> listCSDA;
879 
880  for (i=0; i<n_dedx; ++i) {
881  p = loss_list[i];
882  if(p != em) { p->SetIonisation(false); }
883  if(build_flags[i]) {
885  }
886  if (0 < nSubRegions) {
887  dedx = p->BuildDEDXTable(fSubRestricted);
888  p->SetDEDXTable(dedx,fSubRestricted);
889  listSub.push_back(dedx);
890  if(build_flags[i]) {
892  if(p != em) { em->AddCollaborativeProcess(p); }
893  }
894  }
895  if(theParameters->BuildCSDARange()) {
896  dedx = p->BuildDEDXTable(fTotal);
897  p->SetDEDXTable(dedx,fTotal);
898  listCSDA.push_back(dedx);
899  }
900  }
901 
902  if (0 < nSubRegions) {
903  G4PhysicsTable* dedxSub = em->IonisationTableForSubsec();
904  if (1 < listSub.size()) {
905  em->SetDEDXTable(dedxSub, fIsSubIonisation);
906  dedxSub = 0;
908  tableBuilder->BuildDEDXTable(dedxSub, listSub);
909  em->SetDEDXTable(dedxSub, fSubRestricted);
910  }
911  }
913  G4PhysicsTable* dedxCSDA = em->DEDXunRestrictedTable();
914  if (1 < n_dedx) {
915  dedxCSDA = nullptr;
916  dedxCSDA = G4PhysicsTableHelper::PreparePhysicsTable(dedxCSDA);
917  tableBuilder->BuildDEDXTable(dedxCSDA, listCSDA);
918  em->SetDEDXTable(dedxCSDA,fTotal);
919  }
920  G4PhysicsTable* rCSDA = em->CSDARangeTable();
921  if(!rCSDA) { rCSDA = G4PhysicsTableHelper::PreparePhysicsTable(rCSDA); }
922  tableBuilder->BuildRangeTable(dedxCSDA, rCSDA, true);
923  em->SetCSDARangeTable(rCSDA);
924  }
925 
926  if (1 < verbose) {
927  G4cout << "G4LossTableManager::BuildTables: Tables are built for "
928  << aParticle->GetParticleName()
929  << "; ionisation process: " << em->GetProcessName()
930  << " " << em
931  << G4endl;
932  }
933  return em;
934 }
935 
936 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
937 
939 {
940  return theMessenger;
941 }
942 
943 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
944 
946  const G4ParticleDefinition* aParticle)
947 {
949  ed << "Energy loss process not found for " << aParticle->GetParticleName()
950  << " !";
951  G4Exception("G4LossTableManager::ParticleHaveNoLoss", "em0001",
952  FatalException, ed);
953 }
954 
955 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
956 
958 {
959 }
960 
961 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
962 
964 {
965  subCutoffFlag = val;
966  for(G4int i=0; i<n_loss; ++i) {
967  if(loss_vector[i]) { loss_vector[i]->ActivateSubCutoff(val, r); }
968  }
969 }
970 
971 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
972 
974 {
975  integral = val;
976  integralActive = true;
977  for(G4int i=0; i<n_loss; ++i) {
978  if(loss_vector[i]) { loss_vector[i]->SetIntegral(val); }
979  }
980  size_t emp = emp_vector.size();
981  for (size_t k=0; k<emp; ++k) {
982  if(emp_vector[k]) { emp_vector[k]->SetIntegral(val); }
983  }
984 }
985 
986 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
987 
989 {
990 }
991 
992 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
993 
995 {
996 }
997 
998 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
999 
1001 {
1003 }
1004 
1005 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1006 
1008 {
1010 }
1011 
1012 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1013 
1015 {
1017 }
1018 
1019 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1020 
1022 {
1024 }
1025 
1026 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1027 
1029 {
1031 }
1032 
1033 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1034 
1036 {
1037 }
1038 
1039 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1040 
1042 {
1044 }
1045 
1046 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1047 
1049 {
1050 }
1051 
1052 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1053 
1055 {
1056  if(0.0 < v1 && 0.0 < v2 && v2 < 1.e+50) {
1057  stepFunctionActive = true;
1058  maxRangeVariation = v1;
1059  maxFinalStep = v2;
1060  for(G4int i=0; i<n_loss; ++i) {
1061  if(loss_vector[i]) { loss_vector[i]->SetStepFunction(v1, v2); }
1062  }
1063  } else if(v1 <= 0.0) {
1064  PrintEWarning("SetStepFunction", v1);
1065  } else {
1066  PrintEWarning("SetStepFunction", v2);
1067  }
1068 }
1069 
1070 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1071 
1073 {
1074 }
1075 
1076 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1077 
1079 {
1080 }
1081 
1082 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1083 
1084 void
1087 {
1088  if(stepFunctionActive) {
1090  }
1091  if(integralActive) { p->SetIntegral(integral); }
1092 }
1093 
1094 
1095 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1096 
1097 const std::vector<G4VEnergyLossProcess*>&
1099 {
1100  return loss_vector;
1101 }
1102 
1103 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1104 
1105 const std::vector<G4VEmProcess*>& G4LossTableManager::GetEmProcessVector()
1106 {
1107  return emp_vector;
1108 }
1109 
1110 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1111 
1112 const std::vector<G4VMultipleScattering*>&
1114 {
1115  return msc_vector;
1116 }
1117 
1118 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1119 
1121 {
1123  return emSaturation;
1124 }
1125 
1126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1127 
1129 {
1131  return emConfigurator;
1132 }
1133 
1134 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1135 
1137 {
1138  if(!emElectronIonPair) {
1140  }
1141  return emElectronIonPair;
1142 }
1143 
1144 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1145 
1147 {
1148  if(atomDeexcitation != p) {
1149  delete atomDeexcitation;
1150  atomDeexcitation = p;
1151  }
1152 }
1153 
1154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1155 
1157 {
1158  if(subcutProducer != p) {
1159  delete subcutProducer;
1160  subcutProducer = p;
1161  }
1162 }
1163 
1164 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1165 
1167 {
1168  G4String ss = "G4LossTableManager::" + tit;
1170  /*
1171  ed << "Parameter is out of range: " << val
1172  << " it will have no effect!\n" << " ## "
1173  << " nbins= " << nbinsLambda
1174  << " nbinsPerDecade= " << nbinsPerDecade
1175  << " Emin(keV)= " << minKinEnergy/keV
1176  << " Emax(GeV)= " << maxKinEnergy/GeV;
1177  */
1178  G4Exception(ss, "em0044", JustWarning, ed);
1179 }
1180 
1181 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
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 Spline() const
G4int WorkerVerbose() const
static G4LossTableManager * Instance()
G4VAtomDeexcitation * atomDeexcitation
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void SetVerbose(G4int verb)
void SetIonisation(G4bool val)
G4VSubCutProducer * subcutProducer
std::vector< G4VEmModel * > mod_vector
G4PhysicsTable * SubLambdaTable() const
G4ElectronIonPair * emElectronIonPair
void DeRegister(G4VEnergyLossProcess *p)
G4PhysicsTable * RangeTableForLoss() const
void SetLambdaBinning(G4int val)
G4EmSaturation * emSaturation
void push_back(G4PhysicsVector *)
static G4ThreadLocal G4LossTableManager * instance
void AddCollaborativeProcess(G4VEnergyLossProcess *)
void SetStepFunction(G4double v1, G4double v2)
G4PhysicsTable * CSDARangeTable() const
const G4String & GetName() const
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
#define G4ThreadLocal
Definition: tls.hh:89
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
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
std::vector< G4VEmFluctuationModel * > fmod_vector
bool G4bool
Definition: G4Types.hh:79
G4VEnergyLossProcess * currentLoss
G4ElectronIonPair * ElectronIonPair()
G4EmSaturation * EmSaturation()
void SetMaxEnergyForCSDARange(G4double val)
G4int NumberOfSubCutoffRegions() const
void Register(G4VEnergyLossProcess *p)
const G4ParticleDefinition * BaseParticle() const
const G4int n
void SetDEDXBinningForCSDARange(G4int val)
void BuildInverseRangeTable(const G4PhysicsTable *rangeTable, G4PhysicsTable *invRangeTable, G4bool isIonisation=false)
G4PhysicsTable * DEDXTableForSubsec() const
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
void SetSubCutoff(G4bool val, const G4Region *r=0)
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:129
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)
void SetLambdaTable(G4PhysicsTable *p)
void SetMaxEnergy(G4double val)
const G4String & GetName() const
void SetStepFunction(G4double v1, G4double v2)
G4PhysicsTable * InverseRangeTable() const
void DumpG4BirksCoefficients()
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
G4LossTableBuilder * tableBuilder
const G4ParticleDefinition * Particle() const
void SetLinearLossLimit(G4double val)
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)
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:795
std::vector< G4bool > isActive
void SetSecondaryRangeTable(G4PhysicsTable *p)
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)
std::vector< G4VMultipleScattering * > msc_vector
static const double mm
Definition: G4SIunits.hh:114
void SetAtomDeexcitation(G4VAtomDeexcitation *)
G4EnergyLossMessenger * theMessenger
G4PhysicsTable * BuildDEDXTable(G4EmTableType tType=fRestricted)
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()
G4ProcessVector * GetProcessList() const
std::vector< G4bool > tables_are_built
void SetSplineFlag(G4bool flag)
G4bool IsIonisationProcess() const