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

Last change on this file since 968 was 961, checked in by garnier, 17 years ago

update processes

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