source: trunk/source/processes/electromagnetic/utils/src/G4LossTableManager.cc @ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 32.4 KB
Line 
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,v 1.105 2010/11/04 12:55:09 vnivanch Exp $
27// GEANT4 tag $Name: emutils-V09-03-23 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name:     G4LossTableManager
35//
36// Author:        Vladimir Ivanchenko
37//
38// Creation date: 03.01.2002
39//
40// Modifications:
41//
42// 20-01-03 Migrade to cut per region (V.Ivanchenko)
43// 15-02-03 Lambda table can be scaled (V.Ivanchenko)
44// 17-02-03 Fix problem of store/restore tables (V.Ivanchenko)
45// 10-03-03 Add Ion registration (V.Ivanchenko)
46// 25-03-03 Add deregistration (V.Ivanchenko)
47// 02-04-03 Change messenger (V.Ivanchenko)
48// 26-04-03 Fix retrieve tables (V.Ivanchenko)
49// 13-05-03 Add calculation of precise range (V.Ivanchenko)
50// 23-07-03 Add exchange with G4EnergyLossTables (V.Ivanchenko)
51// 05-10-03 Add G4VEmProcesses registration and Verbose command (V.Ivanchenko)
52// 17-10-03 Add SetParameters method (V.Ivanchenko)
53// 23-10-03 Add control on inactive processes (V.Ivanchenko)
54// 04-11-03 Add checks in RetrievePhysicsTable (V.Ivanchenko)
55// 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
56// 14-01-04 Activate precise range calculation (V.Ivanchenko)
57// 10-03-04 Fix a problem of Precise Range table (V.Ivanchenko)
58// 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivanchenko)
59// 13-01-04 Fix problem which takes place for inactivate eIoni (V.Ivanchenko)
60// 25-01-04 Fix initialisation problem for ions (V.Ivanchenko)
61// 11-03-05 Shift verbose level by 1 (V.Ivantchenko)
62// 10-01-06 PreciseRange -> CSDARange (V.Ivantchenko)
63// 20-01-06 Introduce G4EmTableType to remove repeating code (VI)
64// 23-03-06 Set flag isIonisation (VI)
65// 10-05-06 Add methods  SetMscStepLimitation, FacRange and MscFlag (VI)
66// 22-05-06 Add methods  Set/Get bremsTh (VI)
67// 05-06-06 Do not clear loss_table map between runs (VI)
68// 16-01-07 Create new energy loss table for e+,e-,mu+,mu- and
69//          left ionisation table for further usage (VI)
70// 12-02-07 Add SetSkin, SetLinearLossLimit (V.Ivanchenko)
71// 18-06-07 Move definition of msc parameters to G4EmProcessOptions (V.Ivanchenko)
72// 21-02-08 Added G4EmSaturation (V.Ivanchenko)
73// 12-04-10 Added PreparePhsyicsTables and BuildPhysicsTables entries (V.Ivanchenko)
74//
75// Class Description:
76//
77// -------------------------------------------------------------------
78//
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81
82#include "G4LossTableManager.hh"
83#include "G4EnergyLossMessenger.hh"
84#include "G4PhysicsTable.hh"
85#include "G4ParticleDefinition.hh"
86#include "G4MaterialCutsCouple.hh"
87#include "G4ProcessManager.hh"
88#include "G4Electron.hh"
89#include "G4Proton.hh"
90#include "G4VMultipleScattering.hh"
91#include "G4VEmProcess.hh"
92#include "G4ProductionCutsTable.hh"
93#include "G4PhysicsTableHelper.hh"
94#include "G4EmCorrections.hh"
95#include "G4EmSaturation.hh"
96#include "G4EmConfigurator.hh"
97#include "G4ElectronIonPair.hh"
98#include "G4EmTableType.hh"
99#include "G4LossTableBuilder.hh"
100#include "G4VAtomDeexcitation.hh"
101#include "G4Region.hh"
102
103G4LossTableManager* G4LossTableManager::theInstance = 0;
104
105//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
106
107G4LossTableManager* G4LossTableManager::Instance()
108{
109  if(0 == theInstance) {
110    static G4LossTableManager manager;
111    theInstance = &manager;
112  }
113  return theInstance;
114}
115
116//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
117
118G4LossTableManager::~G4LossTableManager()
119{
120  for (G4int i=0; i<n_loss; ++i) {
121    if( loss_vector[i] ) { delete loss_vector[i]; }
122  }
123  size_t msc = msc_vector.size();
124  for (size_t j=0; j<msc; ++j) {
125    if( msc_vector[j] ) { delete msc_vector[j]; }
126  }
127  size_t emp = emp_vector.size();
128  for (size_t k=0; k<emp; ++k) {
129    if( emp_vector[k] ) { delete emp_vector[k]; }
130  }
131  size_t mod = mod_vector.size();
132  for (size_t a=0; a<mod; ++a) {
133    if( mod_vector[a] ) { delete mod_vector[a]; }
134  }
135  size_t fmod = fmod_vector.size();
136  for (size_t b=0; b<fmod; ++b) {
137    if( fmod_vector[b] ) { delete fmod_vector[b]; }
138  }
139  Clear();
140  delete theMessenger;
141  delete tableBuilder;
142  delete emCorrections;
143  delete emSaturation;
144  delete emConfigurator;
145  delete emElectronIonPair;
146  delete atomDeexcitation;
147}
148
149//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
150
151G4LossTableManager::G4LossTableManager()
152{
153  n_loss = 0;
154  run = 0;
155  startInitialisation = false;
156  all_tables_are_built = false;
157  currentLoss = 0;
158  currentParticle = 0;
159  firstParticle = 0;
160  lossFluctuationFlag = true;
161  subCutoffFlag = false;
162  rndmStepFlag = false;
163  minSubRange = 0.0;
164  maxRangeVariation = 1.0;
165  maxFinalStep = 0.0;
166  minKinEnergy = 0.1*keV;
167  maxKinEnergy = 10.0*TeV;
168  nbinsLambda  = 77;
169  nbinsPerDecade = 7;
170  maxKinEnergyForMuons = 10.*TeV;
171  integral = true;
172  integralActive = false;
173  buildCSDARange = false;
174  minEnergyActive = false;
175  maxEnergyActive = false;
176  maxEnergyForMuonsActive = false;
177  stepFunctionActive = false;
178  flagLPM = true;
179  splineFlag = true;
180  bremsTh = DBL_MAX;
181  factorForAngleLimit = 1.0;
182  verbose = 1;
183  theMessenger = new G4EnergyLossMessenger();
184  theElectron  = G4Electron::Electron();
185  tableBuilder = new G4LossTableBuilder();
186  emCorrections= new G4EmCorrections();
187  emSaturation = new G4EmSaturation();
188  emConfigurator = new G4EmConfigurator(verbose);
189  emElectronIonPair = new G4ElectronIonPair();
190  tableBuilder->SetSplineFlag(splineFlag);
191  atomDeexcitation = 0;
192}
193
194//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
195
196void G4LossTableManager::Clear()
197{
198  all_tables_are_built = false;
199  currentLoss = 0;
200  currentParticle = 0;
201  if(n_loss)
202    {
203      dedx_vector.clear();
204      range_vector.clear();
205      inv_range_vector.clear();
206      loss_map.clear();
207      loss_vector.clear();
208      part_vector.clear();
209      base_part_vector.clear();
210      tables_are_built.clear();
211      isActive.clear();
212      n_loss = 0;
213    }
214}
215
216//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
217
218void G4LossTableManager::Register(G4VEnergyLossProcess* p)
219{
220  ++n_loss;
221  loss_vector.push_back(p);
222  part_vector.push_back(0);
223  base_part_vector.push_back(0);
224  dedx_vector.push_back(0);
225  range_vector.push_back(0);
226  inv_range_vector.push_back(0);
227  tables_are_built.push_back(false);
228  isActive.push_back(true);
229  all_tables_are_built = false;
230  if(!lossFluctuationFlag) { p->SetLossFluctuations(false); }
231  if(subCutoffFlag)        { p->ActivateSubCutoff(true); }
232  if(rndmStepFlag)         { p->SetRandomStep(true); }
233  if(stepFunctionActive)   { p->SetStepFunction(maxRangeVariation, maxFinalStep); }
234  if(integralActive)       { p->SetIntegral(integral); }
235  if(minEnergyActive)      { p->SetMinKinEnergy(minKinEnergy); }
236  if(maxEnergyActive)      { p->SetMaxKinEnergy(maxKinEnergy); }
237  if(verbose > 1) 
238    G4cout << "G4LossTableManager::Register G4VEnergyLossProcess : " 
239           << p->GetProcessName() << G4endl;
240}
241
242//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
243
244void G4LossTableManager::DeRegister(G4VEnergyLossProcess* p)
245{
246  for (G4int i=0; i<n_loss; ++i) {
247    if(loss_vector[i] == p) { loss_vector[i] = 0; }
248  }
249}
250
251//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
252
253void G4LossTableManager::Register(G4VMultipleScattering* p)
254{
255  msc_vector.push_back(p);
256  if(verbose > 1) {
257    G4cout << "G4LossTableManager::Register G4VMultipleScattering : " 
258           << p->GetProcessName() << G4endl;
259  }
260}
261
262//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
263
264void G4LossTableManager::DeRegister(G4VMultipleScattering* p)
265{
266  size_t msc = msc_vector.size();
267  for (size_t i=0; i<msc; ++i) {
268    if(msc_vector[i] == p) { msc_vector[i] = 0; }
269  }
270}
271
272//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
273
274void G4LossTableManager::Register(G4VEmProcess* p)
275{
276  emp_vector.push_back(p);
277  if(verbose > 1) {
278    G4cout << "G4LossTableManager::Register G4VEmProcess : " 
279           << p->GetProcessName() << G4endl;
280  }
281}
282
283//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
284
285void G4LossTableManager::DeRegister(G4VEmProcess* p)
286{
287  size_t emp = emp_vector.size();
288  for (size_t i=0; i<emp; ++i) {
289    if(emp_vector[i] == p) { emp_vector[i] = 0; }
290  }
291}
292
293//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
294
295void G4LossTableManager::Register(G4VEmModel* p)
296{
297  mod_vector.push_back(p);
298  if(verbose > 1) {
299    G4cout << "G4LossTableManager::Register G4VEmModel : " 
300           << p->GetName() << G4endl;
301  }
302}
303
304//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
305
306void G4LossTableManager::DeRegister(G4VEmModel* p)
307{
308  size_t n = mod_vector.size();
309  for (size_t i=0; i<n; ++i) {
310    if(mod_vector[i] == p) { mod_vector[i] = 0; }
311  }
312}
313
314//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
315
316void G4LossTableManager::Register(G4VEmFluctuationModel* p)
317{
318  fmod_vector.push_back(p);
319  if(verbose > 1) {
320    G4cout << "G4LossTableManager::Register G4VEmFluctuationModel : " 
321           << p->GetName() << G4endl;
322  }
323}
324
325//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
326
327void G4LossTableManager::DeRegister(G4VEmFluctuationModel* p)
328{
329  size_t n = fmod_vector.size();
330  for (size_t i=0; i<n; ++i) {
331    if(fmod_vector[i] == p) { fmod_vector[i] = 0; }
332  }
333}
334
335//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
336
337void G4LossTableManager::RegisterIon(const G4ParticleDefinition* ion, 
338                                     G4VEnergyLossProcess* p)
339{
340  loss_map[ion] = p;
341}
342
343//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
344
345void G4LossTableManager::RegisterExtraParticle(
346     const G4ParticleDefinition* part,
347     G4VEnergyLossProcess* p)
348{
349  ++n_loss;
350  loss_vector.push_back(p);
351  part_vector.push_back(part);
352  base_part_vector.push_back(p->BaseParticle());
353  dedx_vector.push_back(0);
354  range_vector.push_back(0);
355  inv_range_vector.push_back(0);
356  tables_are_built.push_back(false);
357  all_tables_are_built = false;
358}
359
360//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
361
362void 
363G4LossTableManager::PreparePhysicsTable(const G4ParticleDefinition* particle,
364                                        G4VEnergyLossProcess* p)
365{
366  if (1 < verbose) {
367    G4cout << "G4LossTableManager::PreparePhysicsTable for " 
368           << particle->GetParticleName() 
369           << " and " << p->GetProcessName() << " run= " << run << G4endl;
370  }
371  // start initialisation for the first run
372  startInitialisation = true;
373
374  if( 0 == run ) {
375    emConfigurator->PrepareModels(particle, p);
376
377    // initialise particles for given process
378    for (G4int j=0; j<n_loss; ++j) {
379      if (p == loss_vector[j]) {
380        if (!part_vector[j]) { part_vector[j] = particle; }
381      }
382    }
383  }
384}
385
386//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
387
388void 
389G4LossTableManager::PreparePhysicsTable(const G4ParticleDefinition* particle,
390                                        G4VEmProcess* p)
391{
392  if (1 < verbose) {
393    G4cout << "G4LossTableManager::PreparePhysicsTable for " 
394           << particle->GetParticleName() 
395           << " and " << p->GetProcessName() << G4endl;
396  }
397  // start initialisation for the first run
398  if( 0 == run ) {
399    emConfigurator->PrepareModels(particle, p);
400  }
401  startInitialisation = true;
402}
403
404//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
405
406void 
407G4LossTableManager::PreparePhysicsTable(const G4ParticleDefinition* particle,
408                                        G4VMultipleScattering* p)
409{
410  if (1 < verbose) {
411    G4cout << "G4LossTableManager::PreparePhysicsTable for " 
412           << particle->GetParticleName() 
413           << " and " << p->GetProcessName() << G4endl;
414  }
415  // start initialisation for the first run
416  if( 0 == run ) {
417    emConfigurator->PrepareModels(particle, p);
418  }
419}
420
421//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
422
423void 
424G4LossTableManager::BuildPhysicsTable(const G4ParticleDefinition*)
425{
426  if(0 == run && startInitialisation) {
427    emConfigurator->Clear();
428  }
429}
430
431//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
432
433void G4LossTableManager::BuildPhysicsTable(
434     const G4ParticleDefinition* aParticle,
435     G4VEnergyLossProcess* p)
436{
437  if(1 < verbose) {
438    G4cout << "### G4LossTableManager::BuildDEDXTable() is requested for "
439           << aParticle->GetParticleName()
440           << " and process " << p->GetProcessName()
441           << G4endl;
442  }
443  // clear configurator
444  if(0 == run && startInitialisation) {
445    emConfigurator->Clear();
446    firstParticle = aParticle; 
447  }
448  if(startInitialisation && atomDeexcitation) {
449    atomDeexcitation->InitialiseAtomicDeexcitation();
450  }
451  startInitialisation = false;
452
453  // initialisation before any table is built
454  if ( aParticle == firstParticle ) {
455    all_tables_are_built = true;
456
457    if(1 < verbose) {
458      G4cout << "### G4LossTableManager start initilisation for first particle "
459             << firstParticle->GetParticleName() 
460             << G4endl;
461    }
462    for (G4int i=0; i<n_loss; ++i) {
463      G4VEnergyLossProcess* el = loss_vector[i];
464
465      if(el) {
466        const G4ProcessManager* pm = el->GetProcessManager();
467        isActive[i] = pm->GetProcessActivation(el);
468        if(0 == run) { base_part_vector[i] = el->BaseParticle(); }
469        tables_are_built[i] = false;
470        all_tables_are_built= false;
471        if(!isActive[i]) { el->SetIonisation(false); }
472 
473        if(1 < verbose) { 
474          G4cout << i <<".   "<< el->GetProcessName() 
475                 << "  for "  << pm->GetParticleType()->GetParticleName()
476                 << "  active= " << pm->GetProcessActivation(el)
477                 << "  table= " << tables_are_built[i]
478                 << "  isIonisation= " << el->IsIonisationProcess();
479          if(base_part_vector[i]) { 
480            G4cout << "  base particle " << base_part_vector[i]->GetParticleName();
481          }
482          G4cout << G4endl;
483        }
484      } else {
485        tables_are_built[i] = true;
486        part_vector[i] = 0;
487      }
488    }
489    ++run;
490    currentParticle = 0;
491  }
492
493  // Set run time parameters
494  SetParameters(aParticle, p);
495
496  if (all_tables_are_built) { return; }
497
498  // Build tables for given particle
499  all_tables_are_built = true;
500
501  for(G4int i=0; i<n_loss; ++i) {
502    if(p == loss_vector[i] && !tables_are_built[i] && !base_part_vector[i]) {
503      const G4ParticleDefinition* curr_part = part_vector[i];
504      if(1 < verbose) {
505        G4cout << "### BuildPhysicsTable for " << p->GetProcessName()
506               << " and " << curr_part->GetParticleName()
507               << " start BuildTable " << G4endl;
508      }
509      G4VEnergyLossProcess* curr_proc = BuildTables(curr_part);
510      if(curr_proc) { CopyTables(curr_part, curr_proc); }
511    }
512    if ( !tables_are_built[i] ) { all_tables_are_built = false; }
513  }
514
515  if(1 < verbose) {
516    G4cout << "### G4LossTableManager::BuildDEDXTable end: "
517           << "all_tables_are_built= " << all_tables_are_built
518           << G4endl;
519   
520    if(all_tables_are_built) {
521      G4cout << "### All dEdx and Range tables are built #####" << G4endl;
522    }
523  }
524}
525
526//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
527
528void G4LossTableManager::CopyTables(const G4ParticleDefinition* part,
529                                          G4VEnergyLossProcess* base_proc)
530{
531  for (G4int j=0; j<n_loss; ++j) {
532
533    G4VEnergyLossProcess* proc = loss_vector[j];
534    //if(proc == base_proc || proc->Particle() == part)
535    //  tables_are_built[j] = true;
536
537    if (!tables_are_built[j] && part == base_part_vector[j]) {
538      tables_are_built[j] = true;
539      proc->SetDEDXTable(base_proc->DEDXTable(),fRestricted);
540      proc->SetDEDXTable(base_proc->DEDXTableForSubsec(),fSubRestricted);
541      proc->SetDEDXTable(base_proc->DEDXunRestrictedTable(),fTotal);
542      proc->SetCSDARangeTable(base_proc->CSDARangeTable());
543      proc->SetRangeTableForLoss(base_proc->RangeTableForLoss());
544      proc->SetInverseRangeTable(base_proc->InverseRangeTable());
545      proc->SetLambdaTable(base_proc->LambdaTable());
546      proc->SetSubLambdaTable(base_proc->SubLambdaTable());
547      proc->SetIonisation(base_proc->IsIonisationProcess());
548      loss_map[part_vector[j]] = proc;
549      if (1 < verbose) {
550         G4cout << "For " << proc->GetProcessName()
551                << " for " << part_vector[j]->GetParticleName()
552                << " base_part= " << part->GetParticleName()
553                << " tables are assigned "
554                << G4endl;
555      }
556    }
557
558    if (theElectron == part && theElectron == proc->SecondaryParticle() ) {
559      proc->SetSecondaryRangeTable(base_proc->RangeTableForLoss());
560    }
561  }
562}
563
564//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
565
566G4VEnergyLossProcess* G4LossTableManager::BuildTables(
567                      const G4ParticleDefinition* aParticle)
568{
569  if(1 < verbose) {
570    G4cout << "G4LossTableManager::BuildTables() for "
571           << aParticle->GetParticleName() << G4endl;
572  }
573
574  std::vector<G4PhysicsTable*> t_list; 
575  std::vector<G4VEnergyLossProcess*> loss_list;
576  loss_list.clear();
577  G4VEnergyLossProcess* em = 0;
578  G4VEnergyLossProcess* p = 0;
579  G4int iem = 0;
580  G4PhysicsTable* dedx = 0;
581  G4int i;
582
583  for (i=0; i<n_loss; ++i) {
584    p = loss_vector[i];
585    if (p && aParticle == part_vector[i] && !tables_are_built[i]) {
586      if ((p->IsIonisationProcess() && isActive[i]) || 
587          !em || (em && !isActive[iem]) ) {
588        em = p;
589        iem= i;
590      }
591      dedx = p->BuildDEDXTable(fRestricted);
592      //      G4cout << "Build DEDX table for " << aParticle->GetParticleName()
593      //             << "  " << dedx << " " << dedx->length() << G4endl;
594      p->SetDEDXTable(dedx,fRestricted); 
595      t_list.push_back(dedx);
596      loss_list.push_back(p);
597      tables_are_built[i] = true;
598    }
599  }
600
601  G4int n_dedx = t_list.size();
602  if (0 == n_dedx || !em) {
603    G4cout << "G4LossTableManager WARNING: no DEDX processes for " 
604           << aParticle->GetParticleName() << G4endl;
605    return 0;
606  }
607  G4int nSubRegions = em->NumberOfSubCutoffRegions();
608
609  if (1 < verbose) {
610    G4cout << "G4LossTableManager::BuildTables() start to build range tables"
611           << " and the sum of " << n_dedx << " processes"
612           << " iem= " << iem << " em= " << em->GetProcessName()
613           << " buildCSDARange= " << buildCSDARange
614           << " nSubRegions= " << nSubRegions
615           << G4endl;
616  }
617
618  dedx = em->IonisationTable();
619  if (1 < n_dedx) {
620    em->SetDEDXTable(dedx, fIsIonisation);
621    dedx = 0;
622    dedx  = G4PhysicsTableHelper::PreparePhysicsTable(dedx);
623    tableBuilder->BuildDEDXTable(dedx, t_list);
624    em->SetDEDXTable(dedx, fRestricted);
625  }
626  dedx_vector[iem] = dedx;
627
628  G4PhysicsTable* range = em->RangeTableForLoss();
629  if(!range) range  = G4PhysicsTableHelper::PreparePhysicsTable(range);
630  range_vector[iem] = range;
631
632  G4PhysicsTable* invrange = em->InverseRangeTable();
633  if(!invrange) invrange = G4PhysicsTableHelper::PreparePhysicsTable(invrange);
634  inv_range_vector[iem]  = invrange;
635
636  G4bool flag = em->IsIonisationProcess();
637  tableBuilder->BuildRangeTable(dedx, range, flag);
638  tableBuilder->BuildInverseRangeTable(range, invrange, flag);
639
640  //  if(1<verbose) G4cout << *dedx << G4endl;
641
642  em->SetRangeTableForLoss(range);
643  em->SetInverseRangeTable(invrange);
644
645  //  if(1<verbose) G4cout << *range << G4endl;
646
647  std::vector<G4PhysicsTable*> listSub;
648  std::vector<G4PhysicsTable*> listCSDA;
649
650  for (i=0; i<n_dedx; ++i) {
651    p = loss_list[i];
652    p->SetIonisation(false);
653    p->SetLambdaTable(p->BuildLambdaTable(fRestricted));
654    if (0 < nSubRegions) {
655      dedx = p->BuildDEDXTable(fSubRestricted);
656      p->SetDEDXTable(dedx,fSubRestricted);
657      listSub.push_back(dedx);
658      p->SetSubLambdaTable(p->BuildLambdaTable(fSubRestricted));
659      if(p != em) em->AddCollaborativeProcess(p);
660    }
661    if(buildCSDARange) { 
662      dedx = p->BuildDEDXTable(fTotal);
663      p->SetDEDXTable(dedx,fTotal);
664      listCSDA.push_back(dedx); 
665    }     
666  }
667
668  if (0 < nSubRegions) {
669    G4PhysicsTable* dedxSub = em->IonisationTableForSubsec();
670    if (1 < listSub.size()) {
671      em->SetDEDXTable(dedxSub, fIsSubIonisation);
672      dedxSub = 0;
673      dedxSub = G4PhysicsTableHelper::PreparePhysicsTable(dedxSub);
674      tableBuilder->BuildDEDXTable(dedxSub, listSub);
675      em->SetDEDXTable(dedxSub, fSubRestricted);
676    }
677  }
678  if(buildCSDARange) {
679    G4PhysicsTable* dedxCSDA = em->DEDXunRestrictedTable();
680    if (1 < n_dedx) {
681      dedxCSDA = 0;
682      dedxCSDA  = G4PhysicsTableHelper::PreparePhysicsTable(dedxCSDA);
683      tableBuilder->BuildDEDXTable(dedxCSDA, listCSDA);
684      em->SetDEDXTable(dedxCSDA,fTotal);
685    }
686    G4PhysicsTable* rCSDA = em->CSDARangeTable();
687    if(!rCSDA) rCSDA = G4PhysicsTableHelper::PreparePhysicsTable(rCSDA);
688    tableBuilder->BuildRangeTable(dedxCSDA, rCSDA, flag);
689    em->SetCSDARangeTable(rCSDA);
690  }
691
692  em->SetIonisation(true);
693  loss_map[aParticle] = em;
694
695  if (1 < verbose) {
696    G4cout << "G4LossTableManager::BuildTables: Tables are built for "
697           << aParticle->GetParticleName()
698           << "; ionisation process: " << em->GetProcessName()
699           << G4endl;
700  }
701  return em;
702}
703
704//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
705
706G4EnergyLossMessenger* G4LossTableManager::GetMessenger()
707{
708  return theMessenger;
709}
710
711//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
712
713void G4LossTableManager::ParticleHaveNoLoss(
714     const G4ParticleDefinition* aParticle)
715{
716  G4String s = " dE/dx table not found for "
717             + aParticle->GetParticleName() + " !";
718  G4Exception("G4LossTableManager::ParticleHaveNoLoss", "EM01",
719              FatalException, s);
720
721}
722
723//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
724
725G4bool G4LossTableManager::BuildCSDARange() const
726{
727  return buildCSDARange;
728}
729
730//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
731
732void G4LossTableManager::SetLossFluctuations(G4bool val)
733{
734  lossFluctuationFlag = val;
735  for(G4int i=0; i<n_loss; ++i) {
736    if(loss_vector[i]) { loss_vector[i]->SetLossFluctuations(val); }
737  }
738}
739
740//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
741
742void G4LossTableManager::SetSubCutoff(G4bool val, const G4Region* r)
743{
744  subCutoffFlag = val;
745  for(G4int i=0; i<n_loss; ++i) {
746    if(loss_vector[i]) { loss_vector[i]->ActivateSubCutoff(val, r); }
747  }
748}
749
750//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
751
752void G4LossTableManager::SetIntegral(G4bool val)
753{
754  integral = val;
755  integralActive = true;
756  for(G4int i=0; i<n_loss; ++i) {
757    if(loss_vector[i]) { loss_vector[i]->SetIntegral(val); }
758  }
759  size_t emp = emp_vector.size();
760  for (size_t k=0; k<emp; ++k) {
761    if(emp_vector[k]) { emp_vector[k]->SetIntegral(val); }
762  }
763}
764
765//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
766
767void G4LossTableManager::SetMinSubRange(G4double val)
768{
769  minSubRange = val;
770  for(G4int i=0; i<n_loss; ++i) {
771    if(loss_vector[i]) { loss_vector[i]->SetMinSubRange(val); }
772  }
773}
774
775//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
776
777void G4LossTableManager::SetRandomStep(G4bool val)
778{
779  rndmStepFlag = val;
780  for(G4int i=0; i<n_loss; ++i) {
781    if(loss_vector[i]) { loss_vector[i]->SetRandomStep(val); }
782  }
783}
784
785//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
786
787void G4LossTableManager::SetMinEnergy(G4double val)
788{
789  minEnergyActive = true;
790  minKinEnergy = val;
791  for(G4int i=0; i<n_loss; ++i) {
792    if(loss_vector[i]) { loss_vector[i]->SetMinKinEnergy(val); }
793  }
794  size_t msc = msc_vector.size();
795  for (size_t j=0; j<msc; ++j) {
796    if(msc_vector[j]) { msc_vector[j]->SetMinKinEnergy(val); }
797  }
798  size_t emp = emp_vector.size();
799  for (size_t k=0; k<emp; ++k) {
800    if(emp_vector[k]) { emp_vector[k]->SetMinKinEnergy(val); }
801  }
802}
803
804//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
805
806void G4LossTableManager::SetMaxEnergy(G4double val)
807{
808  maxEnergyActive = true;
809  maxKinEnergy = val;
810  for(G4int i=0; i<n_loss; ++i) {
811    if(loss_vector[i]) { loss_vector[i]->SetMaxKinEnergy(val); }
812  }
813  size_t msc = msc_vector.size();
814  for (size_t j=0; j<msc; ++j) {
815    if(msc_vector[j]) { msc_vector[j]->SetMaxKinEnergy(val); }
816  }
817  size_t emp = emp_vector.size();
818  for (size_t k=0; k<emp; ++k) {
819    if(emp_vector[k]) { emp_vector[k]->SetMaxKinEnergy(val); }
820  }
821}
822
823//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
824
825void G4LossTableManager::SetMaxEnergyForCSDARange(G4double val)
826{
827  for(G4int i=0; i<n_loss; ++i) {
828    if(loss_vector[i]) { loss_vector[i]->SetMaxKinEnergyForCSDARange(val); }
829  }
830}
831
832//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
833
834void G4LossTableManager::SetMaxEnergyForMuons(G4double val)
835{
836  maxEnergyForMuonsActive = true;
837  maxKinEnergyForMuons = val;
838}
839
840//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
841
842void G4LossTableManager::SetDEDXBinning(G4int val)
843{
844  for(G4int i=0; i<n_loss; ++i) {
845    if(loss_vector[i]) { loss_vector[i]->SetDEDXBinning(val); }
846  }
847}
848
849//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
850
851void G4LossTableManager::SetDEDXBinningForCSDARange(G4int val)
852{
853  for(G4int i=0; i<n_loss; ++i) {
854    if(loss_vector[i]) { loss_vector[i]->SetDEDXBinningForCSDARange(val); }
855  }
856}
857
858//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
859
860void G4LossTableManager::SetLambdaBinning(G4int val)
861{
862  G4int n = val/G4int(std::log10(maxKinEnergy/minKinEnergy));
863  if(n < 5) {
864    G4cout << "G4LossTableManager::SetLambdaBinning WARNING "
865           << "too small number of bins " << val << "  ignored" 
866           << G4endl;
867    return;
868  } 
869  nbinsLambda = val;
870  nbinsPerDecade = n;
871  size_t msc = msc_vector.size();
872  for (size_t j=0; j<msc; ++j) {
873    if(msc_vector[j]) { msc_vector[j]->SetBinning(val); }
874  }
875  size_t emp = emp_vector.size();
876  for (size_t k=0; k<emp; ++k) {
877    if(emp_vector[k]) { emp_vector[k]->SetLambdaBinning(val); }
878  }
879}
880
881//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
882
883G4int G4LossTableManager::GetNumberOfBinsPerDecade() const
884{
885  return nbinsPerDecade;
886}
887
888//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
889
890void G4LossTableManager::SetVerbose(G4int val)
891{
892  verbose = val;
893  for(G4int i=0; i<n_loss; ++i) {
894    if(loss_vector[i]) { loss_vector[i]->SetVerboseLevel(val); }
895  }
896  size_t msc = msc_vector.size();
897  for (size_t j=0; j<msc; ++j) {
898    if(msc_vector[j]) { msc_vector[j]->SetVerboseLevel(val); }
899  }
900  size_t emp = emp_vector.size();
901  for (size_t k=0; k<emp; ++k) {
902    if(emp_vector[k]) { emp_vector[k]->SetVerboseLevel(val); }
903  }
904  emConfigurator->SetVerbose(val);
905  //tableBuilder->SetVerbose(val);
906  //emCorrections->SetVerbose(val);
907  emSaturation->SetVerbose(val);
908  emElectronIonPair->SetVerbose(val);
909  if(atomDeexcitation) { atomDeexcitation->SetVerboseLevel(val); }
910}
911
912//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
913
914void G4LossTableManager::SetStepFunction(G4double v1, G4double v2)
915{
916  stepFunctionActive = true;
917  maxRangeVariation = v1;
918  maxFinalStep = v2;
919  for(G4int i=0; i<n_loss; ++i) {
920    if(loss_vector[i]) { loss_vector[i]->SetStepFunction(v1, v2); }
921  }
922}
923
924//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
925
926void G4LossTableManager::SetLinearLossLimit(G4double val)
927{
928  for(G4int i=0; i<n_loss; ++i) {
929    if(loss_vector[i]) { loss_vector[i]->SetLinearLossLimit(val); }
930  }
931}
932
933//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
934
935void G4LossTableManager::SetBuildCSDARange(G4bool val)
936{
937  buildCSDARange = val;
938}
939
940//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
941
942void 
943G4LossTableManager::SetParameters(const G4ParticleDefinition* aParticle,
944                                  G4VEnergyLossProcess* p)
945{
946  if(stepFunctionActive) { p->SetStepFunction(maxRangeVariation, maxFinalStep); }
947  if(integralActive)     { p->SetIntegral(integral); }
948  if(minEnergyActive)    { p->SetMinKinEnergy(minKinEnergy); }
949  if(maxEnergyActive)    { p->SetMaxKinEnergy(maxKinEnergy); }
950  p->SetVerboseLevel(verbose);
951  if(maxEnergyForMuonsActive) {
952    G4double dm = std::abs(aParticle->GetPDGMass() - 105.7*MeV);
953    if(dm < 5.*MeV) { p->SetMaxKinEnergy(maxKinEnergyForMuons); }
954  }
955}
956
957//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
958
959const std::vector<G4VEnergyLossProcess*>& 
960G4LossTableManager::GetEnergyLossProcessVector()
961{
962  return loss_vector;
963}
964
965//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
966
967const std::vector<G4VEmProcess*>& G4LossTableManager::GetEmProcessVector()
968{
969  return emp_vector;
970}
971
972//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
973
974const std::vector<G4VMultipleScattering*>& 
975G4LossTableManager::GetMultipleScatteringVector()
976{
977  return msc_vector;
978}
979
980//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
981
982void G4LossTableManager::SetLPMFlag(G4bool val)
983{
984  flagLPM = val;
985}
986
987//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
988
989G4bool G4LossTableManager::LPMFlag() const
990{
991  return flagLPM;
992}
993
994//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
995
996void G4LossTableManager::SetSplineFlag(G4bool val)
997{
998  splineFlag = val;
999  tableBuilder->SetSplineFlag(splineFlag);
1000}
1001
1002//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1003
1004G4bool G4LossTableManager::SplineFlag() const
1005{
1006  return splineFlag;
1007}
1008
1009//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1010
1011void G4LossTableManager::SetBremsstrahlungTh(G4double val) 
1012{
1013  bremsTh = val;
1014}
1015
1016//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1017
1018G4double G4LossTableManager::BremsstrahlungTh() const
1019{
1020  return bremsTh;
1021}
1022
1023//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1024
1025void G4LossTableManager::SetFactorForAngleLimit(G4double val) 
1026{
1027  if(val > 0.0) { factorForAngleLimit = val; }
1028}
1029
1030//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1031
1032G4double G4LossTableManager::FactorForAngleLimit() const
1033{
1034  return factorForAngleLimit;
1035}
1036
1037//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1038
1039G4EmCorrections* G4LossTableManager::EmCorrections() 
1040{
1041  return emCorrections;
1042}
1043
1044//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1045
1046G4EmSaturation* G4LossTableManager::EmSaturation()
1047{
1048  return emSaturation;
1049}
1050
1051//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1052
1053G4EmConfigurator* G4LossTableManager::EmConfigurator()
1054{
1055  return emConfigurator;
1056}
1057
1058//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1059
1060G4ElectronIonPair* G4LossTableManager::ElectronIonPair()
1061{
1062  return emElectronIonPair;
1063}
1064
1065//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1066
1067G4VAtomDeexcitation* G4LossTableManager::AtomDeexcitation()
1068{
1069  return atomDeexcitation;
1070}
1071
1072//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1073 
1074void G4LossTableManager::SetAtomDeexcitation(G4VAtomDeexcitation* p)
1075{
1076  atomDeexcitation = p;
1077}
1078//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Note: See TracBrowser for help on using the repository browser.