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

Last change on this file since 953 was 819, checked in by garnier, 17 years ago

import all except CVS

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