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

Last change on this file since 1129 was 1055, checked in by garnier, 15 years ago

maj sur la beta de geant 4.9.3

File size: 28.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.96 2009/04/09 16:10:57 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-beta-cand-01 $
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 Add G4EmSaturation (V.Ivanchenko)
73//
74// Class Description:
75//
76// -------------------------------------------------------------------
77//
78//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80
81#include "G4LossTableManager.hh"
82#include "G4EnergyLossMessenger.hh"
83#include "G4PhysicsTable.hh"
84#include "G4ParticleDefinition.hh"
85#include "G4MaterialCutsCouple.hh"
86#include "G4ProcessManager.hh"
87#include "G4Electron.hh"
88#include "G4Proton.hh"
89#include "G4VMultipleScattering.hh"
90#include "G4VEmProcess.hh"
91#include "G4ProductionCutsTable.hh"
92#include "G4PhysicsTableHelper.hh"
93#include "G4EmCorrections.hh"
94#include "G4EmSaturation.hh"
95#include "G4EmConfigurator.hh"
96#include "G4EmTableType.hh"
97#include "G4LossTableBuilder.hh"
98
99G4LossTableManager* G4LossTableManager::theInstance = 0;
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
102
103G4LossTableManager* G4LossTableManager::Instance()
104{
105  if(0 == theInstance) {
106    static G4LossTableManager manager;
107    theInstance = &manager;
108  }
109  return theInstance;
110}
111
112//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
113
114G4LossTableManager::~G4LossTableManager()
115{
116  for (G4int i=0; i<n_loss; i++) {
117    if( loss_vector[i] ) delete loss_vector[i];
118  }
119  size_t msc = msc_vector.size();
120  for (size_t j=0; j<msc; j++) {
121    if( msc_vector[j] ) delete msc_vector[j];
122  }
123  size_t emp = emp_vector.size();
124  for (size_t k=0; k<emp; k++) {
125    if( emp_vector[k] ) delete emp_vector[k];
126  }
127  size_t mod = mod_vector.size();
128  for (size_t a=0; a<mod; a++) {
129    if( mod_vector[a] ) delete mod_vector[a];
130  }
131  size_t fmod = fmod_vector.size();
132  for (size_t b=0; b<fmod; b++) {
133    if( fmod_vector[b] ) delete fmod_vector[b];
134  }
135  Clear();
136  delete theMessenger;
137  delete tableBuilder;
138  delete emCorrections;
139  delete emSaturation;
140}
141
142//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
143
144G4LossTableManager::G4LossTableManager()
145{
146  n_loss = 0;
147  run = 0;
148  //  first_entry = true;
149  all_tables_are_built = false;
150  all_tables_are_stored = false;
151  currentLoss = 0;
152  currentParticle = 0;
153  lossFluctuationFlag = true;
154  subCutoffFlag = false;
155  rndmStepFlag = false;
156  minSubRange = 0.0;
157  maxRangeVariation = 1.0;
158  maxFinalStep = 0.0;
159  minKinEnergy = 0.1*keV;
160  maxKinEnergy = 100.0*TeV;
161  maxKinEnergyForMuons = 100.*TeV;
162  theMessenger = new G4EnergyLossMessenger();
163  theElectron  = G4Electron::Electron();
164  tableBuilder = new G4LossTableBuilder();
165  emCorrections= new G4EmCorrections();
166  emSaturation = new G4EmSaturation();
167  emConfigurator = new G4EmConfigurator();
168  integral = true;
169  integralActive = false;
170  buildCSDARange = false;
171  minEnergyActive = false;
172  maxEnergyActive = false;
173  maxEnergyForMuonsActive = false;
174  stepFunctionActive = false;
175  flagLPM = true;
176  splineFlag = true;
177  bremsTh = DBL_MAX;
178  verbose = 1;
179  tableBuilder->SetSplineFlag(splineFlag);
180}
181
182//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
183
184void G4LossTableManager::Clear()
185{
186  all_tables_are_built = false;
187  currentLoss = 0;
188  currentParticle = 0;
189  if(n_loss)
190    {
191      dedx_vector.clear();
192      range_vector.clear();
193      inv_range_vector.clear();
194      loss_map.clear();
195      loss_vector.clear();
196      part_vector.clear();
197      base_part_vector.clear();
198      tables_are_built.clear();
199      isActive.clear();
200      n_loss = 0;
201    }
202}
203
204//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
205
206void G4LossTableManager::Register(G4VEnergyLossProcess* p)
207{
208  n_loss++;
209  loss_vector.push_back(p);
210  part_vector.push_back(0);
211  base_part_vector.push_back(0);
212  dedx_vector.push_back(0);
213  range_vector.push_back(0);
214  inv_range_vector.push_back(0);
215  tables_are_built.push_back(false);
216  isActive.push_back(true);
217  all_tables_are_built = false;
218  if(!lossFluctuationFlag) p->SetLossFluctuations(false);
219  if(subCutoffFlag)        p->ActivateSubCutoff(true);
220  if(rndmStepFlag)         p->SetRandomStep(true);
221  if(stepFunctionActive)   p->SetStepFunction(maxRangeVariation, maxFinalStep);
222  if(integralActive)       p->SetIntegral(integral);
223  if(minEnergyActive)      p->SetMinKinEnergy(minKinEnergy);
224  if(maxEnergyActive)      p->SetMaxKinEnergy(maxKinEnergy);
225  if(verbose > 1) 
226    G4cout << "G4LossTableManager::Register G4VEnergyLossProcess : " 
227           << p->GetProcessName() << G4endl;
228}
229
230//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
231
232void G4LossTableManager::DeRegister(G4VEnergyLossProcess* p)
233{
234  for (G4int i=0; i<n_loss; i++) {
235    if(loss_vector[i] == p) loss_vector[i] = 0;
236  }
237}
238
239//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
240
241void G4LossTableManager::Register(G4VMultipleScattering* p)
242{
243  msc_vector.push_back(p);
244  if(verbose > 1) {
245    G4cout << "G4LossTableManager::Register G4VMultipleScattering : " 
246           << p->GetProcessName() << G4endl;
247  }
248}
249
250//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
251
252void G4LossTableManager::DeRegister(G4VMultipleScattering* p)
253{
254  size_t msc = msc_vector.size();
255  for (size_t i=0; i<msc; i++) {
256    if(msc_vector[i] == p) msc_vector[i] = 0;
257  }
258}
259
260//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
261
262void G4LossTableManager::Register(G4VEmProcess* p)
263{
264  emp_vector.push_back(p);
265  if(verbose > 1) {
266    G4cout << "G4LossTableManager::Register G4VEmProcess : " 
267           << p->GetProcessName() << G4endl;
268  }
269}
270
271//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
272
273void G4LossTableManager::DeRegister(G4VEmProcess* p)
274{
275  size_t emp = emp_vector.size();
276  for (size_t i=0; i<emp; i++) {
277    if(emp_vector[i] == p) emp_vector[i] = 0;
278  }
279}
280
281//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
282
283void G4LossTableManager::Register(G4VEmModel* p)
284{
285  mod_vector.push_back(p);
286  if(verbose > 1) {
287    G4cout << "G4LossTableManager::Register G4VEmModel : " 
288           << p->GetName() << G4endl;
289  }
290}
291
292//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
293
294void G4LossTableManager::DeRegister(G4VEmModel* p)
295{
296  size_t n = mod_vector.size();
297  for (size_t i=0; i<n; i++) {
298    if(mod_vector[i] == p) mod_vector[i] = 0;
299  }
300}
301
302//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
303
304void G4LossTableManager::Register(G4VEmFluctuationModel* p)
305{
306  fmod_vector.push_back(p);
307  if(verbose > 1) {
308    G4cout << "G4LossTableManager::Register G4VEmFluctuationModel : " 
309           << p->GetName() << G4endl;
310  }
311}
312
313//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
314
315void G4LossTableManager::DeRegister(G4VEmFluctuationModel* p)
316{
317  size_t n = fmod_vector.size();
318  for (size_t i=0; i<n; i++) {
319    if(fmod_vector[i] == p) fmod_vector[i] = 0;
320  }
321}
322
323//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
324
325void G4LossTableManager::RegisterIon(const G4ParticleDefinition* ion, 
326                                     G4VEnergyLossProcess* p)
327{
328  loss_map[ion] = p;
329}
330
331//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
332
333void G4LossTableManager::RegisterExtraParticle(
334     const G4ParticleDefinition* part,
335     G4VEnergyLossProcess* p)
336{
337  n_loss++;
338  loss_vector.push_back(p);
339  part_vector.push_back(part);
340  base_part_vector.push_back(p->BaseParticle());
341  dedx_vector.push_back(0);
342  range_vector.push_back(0);
343  inv_range_vector.push_back(0);
344  tables_are_built.push_back(false);
345  all_tables_are_built = false;
346}
347
348//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
349
350void G4LossTableManager::EnergyLossProcessIsInitialised(
351                   const G4ParticleDefinition* particle, 
352                   G4VEnergyLossProcess* p)
353{
354  if (run == 0 || (particle == firstParticle && all_tables_are_built) ) {
355    all_tables_are_built = true;
356
357    if(1 < verbose) 
358      G4cout << "### G4LossTableManager start initilisation of tables" 
359             << G4endl;
360    for (G4int i=0; i<n_loss; i++) {
361      G4VEnergyLossProcess* el = loss_vector[i];
362
363      if(el) {
364        const G4ProcessManager* pm = el->GetProcessManager();
365        isActive[i] = pm->GetProcessActivation(el);
366        tables_are_built[i] = false;
367        all_tables_are_built = false;
368        if(!isActive[i]) el->SetIonisation(false);
369 
370        if(1 < verbose) { 
371          G4cout << i <<".   "<< el->GetProcessName() 
372                 << "  for "  << pm->GetParticleType()->GetParticleName()
373                 << "  active= " << pm->GetProcessActivation(el)
374                 << "  table= " << tables_are_built[i]
375                 << "  isIonisation= " << el->IsIonisationProcess()
376                 << G4endl;
377        }
378      } else {
379        tables_are_built[i] = true;
380        part_vector[i] = 0;
381      }
382    }
383    if (0 == run) firstParticle = particle;
384    run++;
385  }
386
387  currentParticle = 0;
388
389  SetParameters(p);
390  for (G4int j=0; j<n_loss; j++) {
391    if (p == loss_vector[j]) {
392
393      if (!part_vector[j]) {
394        part_vector[j] = particle;
395        base_part_vector[j] = p->BaseParticle();
396      }
397      if(maxEnergyForMuonsActive) {
398        G4double dm = std::abs(particle->GetPDGMass() - 105.7*MeV);
399        if(dm < 5.*MeV) p->SetMaxKinEnergy(maxKinEnergyForMuons);
400      }
401
402      if(1 < verbose) {
403        G4cout << "For " << p->GetProcessName()
404               << " for " << part_vector[j]->GetParticleName()
405               << " tables_are_built= " << tables_are_built[j]
406               << " procFlag= " << loss_vector[j]->TablesAreBuilt()
407               << " all_tables_are_built= "     << all_tables_are_built
408               << G4endl;
409      }
410    }
411  }
412}
413
414//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
415
416G4EnergyLossMessenger* G4LossTableManager::GetMessenger()
417{
418  return theMessenger;
419}
420
421//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
422
423void G4LossTableManager::ParticleHaveNoLoss(
424     const G4ParticleDefinition* aParticle)
425{
426  G4String s = " dE/dx table not found for "
427             + aParticle->GetParticleName() + " !";
428  G4Exception("G4LossTableManager::ParticleHaveNoLoss", "EM01",
429              FatalException, s);
430
431}
432
433//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
434
435G4bool G4LossTableManager::BuildCSDARange() const
436{
437  return buildCSDARange;
438}
439
440//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
441
442void G4LossTableManager::BuildPhysicsTable(
443     const G4ParticleDefinition* aParticle,
444     G4VEnergyLossProcess* p)
445{
446  if(1 < verbose) {
447    G4cout << "### G4LossTableManager::BuildDEDXTable() is requested for "
448           << aParticle->GetParticleName()
449           << " and process " << p->GetProcessName()
450           << G4endl;
451  }
452  if (all_tables_are_built) return;
453  all_tables_are_built = true;
454
455  for(G4int i=0; i<n_loss; i++) {
456    if(!tables_are_built[i] && !base_part_vector[i]) {
457      const G4ParticleDefinition* curr_part = part_vector[i];
458      G4VEnergyLossProcess* curr_proc = BuildTables(curr_part);
459      if(curr_proc) CopyTables(curr_part, curr_proc);
460    }
461  }
462
463  for (G4int ii=0; ii<n_loss; ii++) {
464    if ( !tables_are_built[ii] ) {
465      all_tables_are_built = false;
466      break;
467    }
468  }
469
470  if(1 < verbose) {
471    G4cout << "### G4LossTableManager::BuildDEDXTable end: "
472           << "all_tables_are_built= " << all_tables_are_built
473           << G4endl;
474
475    if(all_tables_are_built)
476      G4cout << "### All dEdx and Range tables are built #####" << G4endl;
477  }
478}
479
480//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
481
482void G4LossTableManager::CopyTables(const G4ParticleDefinition* part,
483                                          G4VEnergyLossProcess* base_proc)
484{
485  for (G4int j=0; j<n_loss; j++) {
486
487    G4VEnergyLossProcess* proc = loss_vector[j];
488    if(proc == base_proc || proc->Particle() == part) 
489      tables_are_built[j] = true;
490
491    if (!tables_are_built[j] && part == base_part_vector[j]) {
492      tables_are_built[j] = true;
493      proc->SetDEDXTable(base_proc->DEDXTable(),fRestricted);
494      proc->SetDEDXTable(base_proc->DEDXTableForSubsec(),fSubRestricted);
495      proc->SetDEDXTable(base_proc->DEDXunRestrictedTable(),fTotal);
496      proc->SetCSDARangeTable(base_proc->CSDARangeTable());
497      proc->SetRangeTableForLoss(base_proc->RangeTableForLoss());
498      proc->SetInverseRangeTable(base_proc->InverseRangeTable());
499      proc->SetLambdaTable(base_proc->LambdaTable());
500      proc->SetSubLambdaTable(base_proc->SubLambdaTable());
501      proc->SetIonisation(base_proc->IsIonisationProcess());
502      loss_map[part_vector[j]] = proc;
503      if (1 < verbose) {
504         G4cout << "For " << proc->GetProcessName()
505                << " for " << part_vector[j]->GetParticleName()
506                << " base_part= " << part->GetParticleName()
507                << " tables are assigned "
508                << G4endl;
509      }
510    }
511
512    if (theElectron == part && theElectron == proc->SecondaryParticle() )
513      proc->SetSecondaryRangeTable(base_proc->RangeTableForLoss());
514  }
515}
516
517//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
518
519G4VEnergyLossProcess* G4LossTableManager::BuildTables(
520                      const G4ParticleDefinition* aParticle)
521{
522  if(1 < verbose) {
523    G4cout << "G4LossTableManager::BuildTables() for "
524           << aParticle->GetParticleName() << G4endl;
525  }
526
527  std::vector<G4PhysicsTable*> t_list; 
528  std::vector<G4VEnergyLossProcess*> loss_list;
529  loss_list.clear();
530  G4VEnergyLossProcess* em = 0;
531  G4VEnergyLossProcess* p = 0;
532  G4int iem = 0;
533  G4PhysicsTable* dedx = 0;
534  G4int i;
535
536  for (i=0; i<n_loss; i++) {
537    p = loss_vector[i];
538    if (p && aParticle == part_vector[i] && !tables_are_built[i]) {
539      if ((p->IsIonisationProcess() && isActive[i]) || 
540          !em || (em && !isActive[iem]) ) {
541        em = p;
542        iem= i;
543      }
544      dedx = p->BuildDEDXTable(fRestricted);
545      //      G4cout << "Build DEDX table for " << aParticle->GetParticleName()
546      //             << "  " << dedx << " " << dedx->length() << G4endl;
547      p->SetDEDXTable(dedx,fRestricted); 
548      t_list.push_back(dedx);
549      loss_list.push_back(p);
550      tables_are_built[i] = true;
551    }
552  }
553
554  G4int n_dedx = t_list.size();
555  if (!n_dedx) {
556    G4cout << "G4LossTableManager WARNING: no DEDX processes for " 
557           << aParticle->GetParticleName() << G4endl;
558    return 0;
559  }
560  G4int nSubRegions = em->NumberOfSubCutoffRegions();
561
562  if (1 < verbose) {
563    G4cout << "G4LossTableManager::BuildTables() start to build range tables"
564           << " and the sum of " << n_dedx << " processes"
565           << " iem= " << iem << " em= " << em->GetProcessName()
566           << " buildCSDARange= " << buildCSDARange
567           << " nSubRegions= " << nSubRegions
568           << G4endl;
569  }
570
571  dedx = em->IonisationTable();
572  if (1 < n_dedx) {
573    em->SetDEDXTable(dedx, fIsIonisation);
574    dedx = 0;
575    dedx  = G4PhysicsTableHelper::PreparePhysicsTable(dedx);
576    tableBuilder->BuildDEDXTable(dedx, t_list);
577    em->SetDEDXTable(dedx, fRestricted);
578  }
579  dedx_vector[iem] = dedx;
580
581  G4PhysicsTable* range = em->RangeTableForLoss();
582  if(!range) range  = G4PhysicsTableHelper::PreparePhysicsTable(range);
583  range_vector[iem] = range;
584
585  G4PhysicsTable* invrange = em->InverseRangeTable();
586  if(!invrange) invrange = G4PhysicsTableHelper::PreparePhysicsTable(invrange);
587  inv_range_vector[iem]  = invrange;
588
589  G4bool flag = em->IsIonisationProcess();
590  tableBuilder->BuildRangeTable(dedx, range, flag);
591  tableBuilder->BuildInverseRangeTable(range, invrange, flag);
592
593  //  if(1<verbose) G4cout << *dedx << G4endl;
594
595  em->SetRangeTableForLoss(range);
596  em->SetInverseRangeTable(invrange);
597
598  //  if(1<verbose) G4cout << *range << G4endl;
599
600  std::vector<G4PhysicsTable*> listSub;
601  std::vector<G4PhysicsTable*> listCSDA;
602
603  for (i=0; i<n_dedx; i++) {
604    p = loss_list[i];
605    p->SetIonisation(false);
606    p->SetLambdaTable(p->BuildLambdaTable(fRestricted));
607    if (0 < nSubRegions) {
608      dedx = p->BuildDEDXTable(fSubRestricted);
609      p->SetDEDXTable(dedx,fSubRestricted);
610      listSub.push_back(dedx);
611      p->SetSubLambdaTable(p->BuildLambdaTable(fSubRestricted));
612      if(p != em) em->AddCollaborativeProcess(p);
613    }
614    if(buildCSDARange) { 
615      dedx = p->BuildDEDXTable(fTotal);
616      p->SetDEDXTable(dedx,fTotal);
617      listCSDA.push_back(dedx); 
618    }     
619  }
620
621  if (0 < nSubRegions) {
622    G4PhysicsTable* dedxSub = em->IonisationTableForSubsec();
623    if (1 < listSub.size()) {
624      em->SetDEDXTable(dedxSub, fIsSubIonisation);
625      dedxSub = 0;
626      dedxSub = G4PhysicsTableHelper::PreparePhysicsTable(dedxSub);
627      tableBuilder->BuildDEDXTable(dedxSub, listSub);
628      em->SetDEDXTable(dedxSub, fSubRestricted);
629    }
630  }
631  if(buildCSDARange) {
632    G4PhysicsTable* dedxCSDA = em->DEDXunRestrictedTable();
633    if (1 < n_dedx) {
634      dedxCSDA = 0;
635      dedxCSDA  = G4PhysicsTableHelper::PreparePhysicsTable(dedxCSDA);
636      tableBuilder->BuildDEDXTable(dedxCSDA, listCSDA);
637      em->SetDEDXTable(dedxCSDA,fTotal);
638    }
639    G4PhysicsTable* rCSDA = em->CSDARangeTable();
640    if(!rCSDA) rCSDA = G4PhysicsTableHelper::PreparePhysicsTable(rCSDA);
641    tableBuilder->BuildRangeTable(dedxCSDA, rCSDA, flag);
642    em->SetCSDARangeTable(rCSDA);
643  }
644
645  em->SetIonisation(true);
646  loss_map[aParticle] = em;
647
648  if (1 < verbose) {
649    G4cout << "G4LossTableManager::BuildTables: Tables are built for "
650           << aParticle->GetParticleName()
651           << "; ionisation process: " << em->GetProcessName()
652           << G4endl;
653  }
654  return em;
655}
656
657//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
658
659void G4LossTableManager::SetLossFluctuations(G4bool val)
660{
661  lossFluctuationFlag = val;
662  for(G4int i=0; i<n_loss; i++) {
663    if(loss_vector[i]) loss_vector[i]->SetLossFluctuations(val);
664  }
665}
666
667//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
668
669void G4LossTableManager::SetSubCutoff(G4bool val)
670{
671  subCutoffFlag = val;
672  for(G4int i=0; i<n_loss; i++) {
673    if(loss_vector[i]) loss_vector[i]->ActivateSubCutoff(val);
674  }
675}
676
677//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
678
679void G4LossTableManager::SetIntegral(G4bool val)
680{
681  integral = val;
682  integralActive = true;
683  for(G4int i=0; i<n_loss; i++) {
684    if(loss_vector[i]) loss_vector[i]->SetIntegral(val);
685  }
686  size_t emp = emp_vector.size();
687  for (size_t k=0; k<emp; k++) {
688    if(emp_vector[k]) emp_vector[k]->SetIntegral(val);
689  }
690}
691
692//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
693
694void G4LossTableManager::SetMinSubRange(G4double val)
695{
696  minSubRange = val;
697  for(G4int i=0; i<n_loss; i++) {
698    if(loss_vector[i]) loss_vector[i]->SetMinSubRange(val);
699  }
700}
701
702//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
703
704void G4LossTableManager::SetRandomStep(G4bool val)
705{
706  rndmStepFlag = val;
707  for(G4int i=0; i<n_loss; i++) {
708    if(loss_vector[i]) loss_vector[i]->SetRandomStep(val);
709  }
710}
711
712//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
713
714void G4LossTableManager::SetMinEnergy(G4double val)
715{
716  minEnergyActive = true;
717  minKinEnergy = val;
718  for(G4int i=0; i<n_loss; i++) {
719    if(loss_vector[i]) loss_vector[i]->SetMinKinEnergy(val);
720  }
721  size_t msc = msc_vector.size();
722  for (size_t j=0; j<msc; j++) {
723    if(msc_vector[j]) msc_vector[j]->SetMinKinEnergy(val);
724  }
725  size_t emp = emp_vector.size();
726  for (size_t k=0; k<emp; k++) {
727    if(emp_vector[k]) emp_vector[k]->SetMinKinEnergy(val);
728  }
729}
730
731//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
732
733void G4LossTableManager::SetMaxEnergy(G4double val)
734{
735  maxEnergyActive = true;
736  maxKinEnergy = val;
737  for(G4int i=0; i<n_loss; i++) {
738    if(loss_vector[i]) loss_vector[i]->SetMaxKinEnergy(val);
739  }
740  size_t msc = msc_vector.size();
741  for (size_t j=0; j<msc; j++) {
742    if(msc_vector[j]) msc_vector[j]->SetMaxKinEnergy(val);
743  }
744  size_t emp = emp_vector.size();
745  for (size_t k=0; k<emp; k++) {
746    if(emp_vector[k]) emp_vector[k]->SetMaxKinEnergy(val);
747  }
748}
749
750//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
751
752void G4LossTableManager::SetMaxEnergyForCSDARange(G4double val)
753{
754  for(G4int i=0; i<n_loss; i++) {
755    if(loss_vector[i]) loss_vector[i]->SetMaxKinEnergyForCSDARange(val);
756  }
757}
758
759//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
760
761void G4LossTableManager::SetMaxEnergyForMuons(G4double val)
762{
763  maxEnergyForMuonsActive = true;
764  maxKinEnergyForMuons = val;
765}
766
767//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
768
769void G4LossTableManager::SetDEDXBinning(G4int val)
770{
771  for(G4int i=0; i<n_loss; i++) {
772    if(loss_vector[i]) loss_vector[i]->SetDEDXBinning(val);
773  }
774}
775
776//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
777
778void G4LossTableManager::SetDEDXBinningForCSDARange(G4int val)
779{
780  for(G4int i=0; i<n_loss; i++) {
781    if(loss_vector[i]) loss_vector[i]->SetDEDXBinningForCSDARange(val);
782  }
783}
784
785//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
786
787void G4LossTableManager::SetLambdaBinning(G4int val)
788{
789  size_t msc = msc_vector.size();
790  for (size_t j=0; j<msc; j++) {
791    if(msc_vector[j]) msc_vector[j]->SetBinning(val);
792  }
793  size_t emp = emp_vector.size();
794  for (size_t k=0; k<emp; k++) {
795    if(emp_vector[k]) emp_vector[k]->SetLambdaBinning(val);
796  }
797}
798
799//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
800
801void G4LossTableManager::SetVerbose(G4int val)
802{
803  verbose = val;
804  for(G4int i=0; i<n_loss; i++) {
805    if(loss_vector[i]) loss_vector[i]->SetVerboseLevel(val);
806  }
807  size_t msc = msc_vector.size();
808  for (size_t j=0; j<msc; j++) {
809    if(msc_vector[j]) msc_vector[j]->SetVerboseLevel(val);
810  }
811  size_t emp = emp_vector.size();
812  for (size_t k=0; k<emp; k++) {
813    if(emp_vector[k]) emp_vector[k]->SetVerboseLevel(val);
814  }
815}
816
817//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
818
819void G4LossTableManager::SetStepFunction(G4double v1, G4double v2)
820{
821  stepFunctionActive = true;
822  maxRangeVariation = v1;
823  maxFinalStep = v2;
824  for(G4int i=0; i<n_loss; i++) {
825    if(loss_vector[i]) loss_vector[i]->SetStepFunction(v1, v2);
826  }
827}
828
829//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
830
831void G4LossTableManager::SetLinearLossLimit(G4double val)
832{
833  for(G4int i=0; i<n_loss; i++) {
834    if(loss_vector[i]) loss_vector[i]->SetLinearLossLimit(val);
835  }
836}
837
838//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
839
840void G4LossTableManager::SetBuildCSDARange(G4bool val)
841{
842  buildCSDARange = val;
843}
844
845//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
846
847void G4LossTableManager::SetParameters(G4VEnergyLossProcess* p)
848{
849  if(stepFunctionActive) p->SetStepFunction(maxRangeVariation, maxFinalStep);
850  if(integralActive)     p->SetIntegral(integral);
851  if(minEnergyActive)    p->SetMinKinEnergy(minKinEnergy);
852  if(maxEnergyActive)    p->SetMaxKinEnergy(maxKinEnergy);
853  p->SetVerboseLevel(verbose);
854}
855
856//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
857
858const std::vector<G4VEnergyLossProcess*>& 
859      G4LossTableManager::GetEnergyLossProcessVector()
860{
861  return loss_vector;
862}
863
864//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
865
866const std::vector<G4VEmProcess*>& G4LossTableManager::GetEmProcessVector()
867{
868  return emp_vector;
869}
870
871//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
872
873const std::vector<G4VMultipleScattering*>& 
874      G4LossTableManager::GetMultipleScatteringVector()
875{
876  return msc_vector;
877}
878
879//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
880
881void G4LossTableManager::SetLPMFlag(G4bool val)
882{
883  flagLPM = val;
884}
885
886//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
887
888G4bool G4LossTableManager::LPMFlag() const
889{
890  return flagLPM;
891}
892
893//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
894
895void G4LossTableManager::SetSplineFlag(G4bool val)
896{
897  splineFlag = val;
898  tableBuilder->SetSplineFlag(splineFlag);
899}
900
901//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
902
903G4bool G4LossTableManager::SplineFlag() const
904{
905  return splineFlag;
906}
907
908//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
909
910void G4LossTableManager::SetBremsstrahlungTh(G4double val) 
911{
912  bremsTh = val;
913}
914
915//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
916
917G4double G4LossTableManager::BremsstrahlungTh() const
918{
919  return bremsTh;
920}
921
922//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
923
924G4EmCorrections* G4LossTableManager::EmCorrections() 
925{
926  return emCorrections;
927}
928
929//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
930
931G4EmSaturation* G4LossTableManager::EmSaturation()
932{
933  return emSaturation;
934}
935
936//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
937
938G4EmConfigurator* G4LossTableManager::EmConfigurator()
939{
940  return emConfigurator;
941}
942
943//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
Note: See TracBrowser for help on using the repository browser.