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

Last change on this file since 1198 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

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