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

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

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

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