Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LossTableManager.cc
Go to the documentation of this file.
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// -------------------------------------------------------------------
27//
28// GEANT4 Class file
29//
30//
31// File name: G4LossTableManager
32//
33// Author: Vladimir Ivanchenko
34//
35// Creation date: 03.01.2002
36//
37// Modifications: by V.Ivanchenko
38//
39//
40// Class Description:
41//
42// -------------------------------------------------------------------
43//
44//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46
47#include "G4LossTableManager.hh"
48#include "G4SystemOfUnits.hh"
49
51#include "G4VEmProcess.hh"
52#include "G4VXRayModel.hh"
53
54#include "G4EmParameters.hh"
55#include "G4EmSaturation.hh"
56#include "G4EmConfigurator.hh"
57#include "G4ElectronIonPair.hh"
58#include "G4NIELCalculator.hh"
59#include "G4EmCorrections.hh"
60#include "G4LossTableBuilder.hh"
62#include "G4VSubCutProducer.hh"
63#include "G4VXRayModel.hh"
64
65#include "G4PhysicsTable.hh"
68#include "G4ProcessManager.hh"
69#include "G4Electron.hh"
70#include "G4Proton.hh"
73#include "G4EmTableType.hh"
74#include "G4Region.hh"
76
77#include "G4Gamma.hh"
78#include "G4Positron.hh"
79#include "G4OpticalPhoton.hh"
80#include "G4Neutron.hh"
81#include "G4MuonPlus.hh"
82#include "G4MuonMinus.hh"
83#include "G4GenericIon.hh"
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
86
87static std::once_flag applyOnce;
88G4ThreadLocal G4LossTableManager* G4LossTableManager::instance = nullptr;
89
91{
92 if(nullptr == instance) {
94 instance = inst.Instance();
95 }
96 return instance;
97}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
100
102{
103 for (auto const & p : loss_vector) { delete p; }
104 for (auto const & p : msc_vector) { delete p; }
105 for (auto const & p : emp_vector) { delete p; }
106 for (auto const & p : p_vector) { delete p; }
107 for (auto const & p : xray_vector) { delete p; }
108
109 std::size_t mod = mod_vector.size();
110 std::size_t fmod = fmod_vector.size();
111 for (std::size_t a=0; a<mod; ++a) {
112 if( nullptr != mod_vector[a] ) {
113 for (std::size_t b=0; b<fmod; ++b) {
114 if((G4VEmModel*)(fmod_vector[b]) == mod_vector[a]) {
115 fmod_vector[b] = nullptr;
116 }
117 }
118 delete mod_vector[a];
119 mod_vector[a] = nullptr;
120 }
121 }
122 for (auto const & p : fmod_vector) { delete p; }
123
124 Clear();
125 delete tableBuilder;
126 delete emCorrections;
127 delete emConfigurator;
128 delete emElectronIonPair;
129 delete nielCalculator;
130 delete atomDeexcitation;
131 delete subcutProducer;
132}
133
134//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
135
137{
138 theParameters = G4EmParameters::Instance();
139 theElectron = G4Electron::Electron();
140
141 // only one thread is the master
142 std::call_once(applyOnce, [this]() { isMaster = true; });
143 verbose = isMaster ? theParameters->Verbose() : theParameters->WorkerVerbose();
144
145 tableBuilder = new G4LossTableBuilder(isMaster);
146 emCorrections = new G4EmCorrections(verbose);
147
148 std::size_t n = 70;
149 loss_vector.reserve(n);
150 part_vector.reserve(n);
151 base_part_vector.reserve(n);
152 dedx_vector.reserve(n);
153 range_vector.reserve(n);
154 inv_range_vector.reserve(n);
155 tables_are_built.reserve(n);
156 isActive.reserve(n);
157 msc_vector.reserve(10);
158 emp_vector.reserve(16);
159 mod_vector.reserve(150);
160 fmod_vector.reserve(60);
161}
162
163//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
164
165void G4LossTableManager::Clear()
166{
167 all_tables_are_built = false;
168 currentLoss = nullptr;
169 currentParticle = nullptr;
170 if(n_loss) {
171 dedx_vector.clear();
172 range_vector.clear();
173 inv_range_vector.clear();
174 loss_map.clear();
175 loss_vector.clear();
176 part_vector.clear();
177 base_part_vector.clear();
178 tables_are_built.clear();
179 isActive.clear();
180 n_loss = 0;
181 }
182}
183
184//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
185
187{
188 if (nullptr == p) { return; }
189 for (G4int i=0; i<n_loss; ++i) {
190 if(loss_vector[i] == p) { return; }
191 }
192 if(verbose > 1) {
193 G4cout << "G4LossTableManager::Register G4VEnergyLossProcess : "
194 << p->GetProcessName() << " idx= " << n_loss << G4endl;
195 }
196 ++n_loss;
197 loss_vector.push_back(p);
198 part_vector.push_back(nullptr);
199 base_part_vector.push_back(nullptr);
200 dedx_vector.push_back(nullptr);
201 range_vector.push_back(nullptr);
202 inv_range_vector.push_back(nullptr);
203 tables_are_built.push_back(false);
204 isActive.push_back(true);
205 all_tables_are_built = false;
206}
207
208//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
209
211{
212 // initialisation once per run
213 if (!resetParam) { return; }
214 resetParam = false;
215 startInitialisation = true;
216 verbose = theParameters->Verbose();
217 if(!isMaster) {
218 verbose = theParameters->WorkerVerbose();
219 } else {
220 if(verbose > 0) { theParameters->Dump(); }
221 }
222
223 tableBuilder->InitialiseBaseMaterials();
224 if (nullptr != nielCalculator) { nielCalculator->Initialise(); }
225
226 emCorrections->SetVerbose(verbose);
227 if(nullptr != emConfigurator) { emConfigurator->SetVerbose(verbose); };
228 if(nullptr != emElectronIonPair) { emElectronIonPair->SetVerbose(verbose); };
229 if(nullptr != atomDeexcitation) {
230 atomDeexcitation->SetVerboseLevel(verbose);
231 atomDeexcitation->InitialiseAtomicDeexcitation();
232 }
233 if (1 < verbose) {
234 G4cout << "====== G4LossTableManager::ResetParameters "
235 << " Nloss=" << loss_vector.size()
236 << " run=" << run << " master=" << isMaster
237 << G4endl;
238 }
239}
240
241//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
242
244{
245 if (nullptr == p) { return; }
246 for (G4int i=0; i<n_loss; ++i) {
247 if(loss_vector[i] == p) {
248 loss_vector[i] = nullptr;
249 break;
250 }
251 }
252}
253
254//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
255
257{
258 if (nullptr == p) { return; }
259 std::size_t n = msc_vector.size();
260 for (std::size_t i=0; i<n; ++i) {
261 if(msc_vector[i] == p) { return; }
262 }
263 if(verbose > 1) {
264 G4cout << "G4LossTableManager::Register G4VMultipleScattering : "
265 << p->GetProcessName() << " idx= " << msc_vector.size() << G4endl;
266 }
267 msc_vector.push_back(p);
268}
269
270//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
271
273{
274 if (nullptr == p) { return; }
275 std::size_t msc = msc_vector.size();
276 for (std::size_t i=0; i<msc; ++i) {
277 if(msc_vector[i] == p) {
278 msc_vector[i] = nullptr;
279 break;
280 }
281 }
282}
283
284//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
285
287{
288 if (nullptr == p) { return; }
289 std::size_t n = emp_vector.size();
290 for (std::size_t i=0; i<n; ++i) {
291 if(emp_vector[i] == p) { return; }
292 }
293 if(verbose > 1) {
294 G4cout << "G4LossTableManager::Register G4VEmProcess : "
295 << p->GetProcessName() << " idx= " << emp_vector.size() << G4endl;
296 }
297 emp_vector.push_back(p);
298}
299
300//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
301
303{
304 if (nullptr == p) { return; }
305 std::size_t emp = emp_vector.size();
306 for (std::size_t i=0; i<emp; ++i) {
307 if(emp_vector[i] == p) {
308 emp_vector[i] = nullptr;
309 break;
310 }
311 }
312}
313
314//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
315
317{
318 if (nullptr == p) { return; }
319 std::size_t n = p_vector.size();
320 for (std::size_t i=0; i<n; ++i) {
321 if(p_vector[i] == p) { return; }
322 }
323 if(verbose > 1) {
324 G4cout << "G4LossTableManager::Register G4VProcess : "
325 << p->GetProcessName() << " idx= " << p_vector.size() << G4endl;
326 }
327 p_vector.push_back(p);
328}
329
330//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
331
333{
334 if (nullptr == p) { return; }
335 std::size_t emp = p_vector.size();
336 for (std::size_t i=0; i<emp; ++i) {
337 if(p_vector[i] == p) {
338 p_vector[i] = nullptr;
339 break;
340 }
341 }
342}
343
344//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
345
347{
348 if (nullptr == p) { return; }
349 std::size_t n = mod_vector.size();
350 for (std::size_t i=0; i<n; ++i) { if (mod_vector[i] == p) { return; } }
351 mod_vector.push_back(p);
352 if(verbose > 1) {
353 G4cout << "G4LossTableManager::Register G4VEmModel : "
354 << p->GetName() << " " << p << " " << mod_vector.size() << G4endl;
355 }
356}
357
358//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
359
361{
362 std::size_t n = mod_vector.size();
363 for (std::size_t i=0; i<n; ++i) {
364 if(mod_vector[i] == p) {
365 mod_vector[i] = nullptr;
366 break;
367 }
368 }
369}
370
371//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
372
374{
375 if (nullptr == p) { return; }
376 std::size_t n = fmod_vector.size();
377 for (std::size_t i=0; i<n; ++i) { if (fmod_vector[i] == p) { return; } }
378 fmod_vector.push_back(p);
379 if(verbose > 1) {
380 G4cout << "G4LossTableManager::Register G4VEmFluctuationModel : "
381 << p->GetName() << " " << fmod_vector.size() << G4endl;
382 }
383}
384
385//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
386
388{
389 std::size_t n = fmod_vector.size();
390 for (std::size_t i=0; i<n; ++i) {
391 if(fmod_vector[i] == p) {
392 fmod_vector[i] = nullptr;
393 break;
394 }
395 }
396}
397
398//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
399
401{
402 if (nullptr == p) { return; }
403 std::size_t n = xray_vector.size();
404 for (std::size_t i=0; i<n; ++i) { if (xray_vector[i] == p) { return; } }
405 xray_vector.push_back(p);
406 if (verbose > 1) {
407 G4cout << "G4LossTableManager::Register G4VXRayModel : "
408 << p->GetName() << " " << xray_vector.size() << G4endl;
409 }
410}
411
412//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
413
415{
416 std::size_t n = xray_vector.size();
417 for (std::size_t i=0; i<n; ++i) {
418 if (xray_vector[i] == p) {
419 xray_vector[i] = nullptr;
420 break;
421 }
422 }
423}
424
425//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
426
428 const G4ParticleDefinition* part,
430{
431 if (nullptr == p || nullptr == part) { return; }
432 for (G4int i=0; i<n_loss; ++i) {
433 if(loss_vector[i] == p) { return; }
434 }
435 if(verbose > 1) {
436 G4cout << "G4LossTableManager::RegisterExtraParticle "
437 << part->GetParticleName() << " G4VEnergyLossProcess : "
438 << p->GetProcessName() << " idx= " << n_loss << G4endl;
439 }
440 ++n_loss;
441 loss_vector.push_back(p);
442 part_vector.push_back(part);
443 base_part_vector.push_back(p->BaseParticle());
444 dedx_vector.push_back(nullptr);
445 range_vector.push_back(nullptr);
446 inv_range_vector.push_back(nullptr);
447 tables_are_built.push_back(false);
448 all_tables_are_built = false;
449}
450
451//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
452
455{
456 if(aParticle != currentParticle) {
457 currentParticle = aParticle;
458 std::map<PD,G4VEnergyLossProcess*,std::less<PD> >::const_iterator pos;
459 if ((pos = loss_map.find(aParticle)) != loss_map.end()) {
460 currentLoss = (*pos).second;
461 } else {
462 currentLoss = nullptr;
463 if(0.0 != aParticle->GetPDGCharge() &&
464 (pos = loss_map.find(theGenericIon)) != loss_map.end()) {
465 currentLoss = (*pos).second;
466 }
467 }
468 }
469 return currentLoss;
470}
471
472//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
473
474void
477{
478 if (1 < verbose) {
479 G4cout << "G4LossTableManager::PreparePhysicsTable for "
480 << particle->GetParticleName()
481 << " and " << p->GetProcessName() << " run= " << run
482 << " loss_vector " << loss_vector.size()
483 << " run=" << run << " master=" << isMaster
484 << G4endl;
485 }
486
487 // start initialisation for the first run
488 if ( -1 == run ) {
489 if (nullptr != emConfigurator) {
490 emConfigurator->PrepareModels(particle, p);
491 }
492
493 // initialise particles for given process
494 for (G4int j=0; j<n_loss; ++j) {
495 if (p == loss_vector[j] && nullptr == part_vector[j]) {
496 part_vector[j] = particle;
497 if (particle->GetParticleName() == "GenericIon") {
498 theGenericIon = particle;
499 }
500 }
501 }
502 }
504}
505
506//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
507
508void
510 G4VEmProcess* p)
511{
512 if (1 < verbose) {
513 G4cout << "G4LossTableManager::PreparePhysicsTable for "
514 << particle->GetParticleName()
515 << " and " << p->GetProcessName()
516 << " run=" << run << " master=" << isMaster
517 << G4endl;
518 }
519
520 // start initialisation for the first run
521 if( -1 == run ) {
522 if (nullptr != emConfigurator) { emConfigurator->PrepareModels(particle, p); }
523 }
524
526}
527
528//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
529
530void
533{
534 if (1 < verbose) {
535 G4cout << "G4LossTableManager::PreparePhysicsTable for "
536 << particle->GetParticleName()
537 << " and " << p->GetProcessName()
538 << " run=" << run << " master=" << isMaster
539 << G4endl;
540 }
541
542 // start initialisation for the first run
543 if ( -1 == run ) {
544 if (nullptr != emConfigurator) {
545 emConfigurator->PrepareModels(particle, p);
546 }
547 }
548
550}
551
552//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
553
554void
556{
557 if(-1 == run && startInitialisation) {
558 if (nullptr != emConfigurator) { emConfigurator->Clear(); }
559 }
560 if (startInitialisation) { resetParam = true; }
561}
562
563//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
564
566 const G4ParticleDefinition* aParticle,
568{
569 if (1 < verbose) {
570 G4cout << "### G4LossTableManager::LocalPhysicsTable() for "
571 << aParticle->GetParticleName()
572 << " and process " << p->GetProcessName()
573 << G4endl;
574 }
575
576 if(-1 == run && startInitialisation) {
577 if (nullptr != emConfigurator) { emConfigurator->Clear(); }
578 firstParticle = aParticle;
579 }
580
581 if (startInitialisation) {
582 ++run;
583 if (1 < verbose) {
584 G4cout << "===== G4LossTableManager::LocalPhysicsTable() for run "
585 << run << " =====" << G4endl;
586 }
587 currentParticle = nullptr;
588 startInitialisation = false;
589 resetParam = true;
590 for (G4int i=0; i<n_loss; ++i) {
591 if (nullptr != loss_vector[i]) {
592 tables_are_built[i] = false;
593 } else {
594 tables_are_built[i] = true;
595 part_vector[i] = nullptr;
596 }
597 }
598 }
599
600 all_tables_are_built= true;
601 for (G4int i=0; i<n_loss; ++i) {
602 if(p == loss_vector[i]) {
603 tables_are_built[i] = true;
604 isActive[i] = true;
605 part_vector[i] = p->Particle();
606 base_part_vector[i] = p->BaseParticle();
607 dedx_vector[i] = p->DEDXTable();
608 range_vector[i] = p->RangeTableForLoss();
609 inv_range_vector[i] = p->InverseRangeTable();
610 if(0 == run && p->IsIonisationProcess()) {
611 loss_map[part_vector[i]] = p;
612 }
613
614 if(1 < verbose) {
615 G4cout << i <<". "<< p->GetProcessName();
616 if(part_vector[i]) {
617 G4cout << " for " << part_vector[i]->GetParticleName();
618 }
619 G4cout << " active= " << isActive[i]
620 << " table= " << tables_are_built[i]
621 << " isIonisation= " << p->IsIonisationProcess()
622 << G4endl;
623 }
624 break;
625 } else if(!tables_are_built[i]) {
626 all_tables_are_built = false;
627 }
628 }
629
630 if(1 < verbose) {
631 G4cout << "### G4LossTableManager::LocalPhysicsTable end"
632 << G4endl;
633 }
634 if(all_tables_are_built) {
635 if(1 < verbose) {
636 G4cout << "%%%%% All dEdx and Range tables for worker are ready for run "
637 << run << " %%%%%" << G4endl;
638 }
639 }
640}
641
642//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
643
645 const G4ParticleDefinition* aParticle,
647{
648 if(1 < verbose) {
649 G4cout << "### G4LossTableManager::BuildPhysicsTable() for "
650 << aParticle->GetParticleName()
651 << " and process " << p->GetProcessName() << G4endl;
652 }
653 // clear configurator
654 if(-1 == run && startInitialisation) {
655 if( nullptr != emConfigurator) { emConfigurator->Clear(); }
656 firstParticle = aParticle;
657 }
658 if(startInitialisation) {
659 ++run;
660 resetParam = true;
661 startInitialisation = false;
662 if(1 < verbose) {
663 G4cout << "===== G4LossTableManager::BuildPhysicsTable() for run "
664 << run << " ===== " << atomDeexcitation << G4endl;
665 }
666 currentParticle = nullptr;
667 all_tables_are_built = false;
668
669 for (G4int i=0; i<n_loss; ++i) {
670 G4VEnergyLossProcess* el = loss_vector[i];
671
672 if(nullptr != el) {
673 isActive[i] = true;
674 part_vector[i] = el->Particle();
675 base_part_vector[i] = el->BaseParticle();
676 tables_are_built[i] = false;
677 if(1 < verbose) {
678 G4cout << i <<". "<< el->GetProcessName();
679 if(el->Particle()) {
680 G4cout << " for " << el->Particle()->GetParticleName();
681 }
682 G4cout << " active= " << isActive[i]
683 << " table= " << tables_are_built[i]
684 << " isIonisation= " << el->IsIonisationProcess();
685 if(base_part_vector[i]) {
686 G4cout << " base particle "
687 << base_part_vector[i]->GetParticleName();
688 }
689 G4cout << G4endl;
690 }
691 } else {
692 tables_are_built[i] = true;
693 part_vector[i] = nullptr;
694 isActive[i] = false;
695 }
696 }
697 }
698
699 if (all_tables_are_built) {
700 theParameters->SetIsPrintedFlag(true);
701 return;
702 }
703
704 // Build tables for given particle
705 all_tables_are_built = true;
706
707 for(G4int i=0; i<n_loss; ++i) {
708 if(p == loss_vector[i] && !tables_are_built[i] && nullptr == base_part_vector[i]) {
709 const G4ParticleDefinition* curr_part = part_vector[i];
710 if(1 < verbose) {
711 G4cout << "### Build Table for " << p->GetProcessName()
712 << " and " << curr_part->GetParticleName()
713 << " " << tables_are_built[i] << " " << base_part_vector[i]
714 << G4endl;
715 }
716 G4VEnergyLossProcess* curr_proc = BuildTables(curr_part);
717 if(curr_proc) {
718 CopyTables(curr_part, curr_proc);
719 if(p == curr_proc && 0 == run && p->IsIonisationProcess()) {
720 loss_map[aParticle] = p;
721 //G4cout << "G4LossTableManager::BuildPhysicsTable: "
722 // << aParticle->GetParticleName()
723 // << " added to map " << p << G4endl;
724 }
725 }
726 }
727 if ( !tables_are_built[i] ) { all_tables_are_built = false; }
728 }
729 if(1 < verbose) {
730 G4cout << "### G4LossTableManager::BuildPhysicsTable end: "
731 << "all_tables_are_built= " << all_tables_are_built << " "
732 << aParticle->GetParticleName() << " proc: " << p << G4endl;
733 }
734}
735
736//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
737
738void G4LossTableManager::CopyTables(const G4ParticleDefinition* part,
739 G4VEnergyLossProcess* base_proc)
740{
741 for (G4int j=0; j<n_loss; ++j) {
742
743 G4VEnergyLossProcess* proc = loss_vector[j];
744
745 if (!tables_are_built[j] && part == base_part_vector[j]) {
746 tables_are_built[j] = true;
747 // for base particle approach only ionisation table should be used
748 proc->SetDEDXTable(base_proc->IonisationTable(),fRestricted);
749 proc->SetDEDXTable(base_proc->DEDXunRestrictedTable(),fTotal);
750 proc->SetCSDARangeTable(base_proc->CSDARangeTable());
751 proc->SetRangeTableForLoss(base_proc->RangeTableForLoss());
752 proc->SetInverseRangeTable(base_proc->InverseRangeTable());
753 proc->SetLambdaTable(base_proc->LambdaTable());
754 if(proc->IsIonisationProcess()) {
755 range_vector[j] = base_proc->RangeTableForLoss();
756 inv_range_vector[j] = base_proc->InverseRangeTable();
757 loss_map[part_vector[j]] = proc;
758 //G4cout << "G4LossTableManager::CopyTable "
759 // << part_vector[j]->GetParticleName()
760 // << " added to map " << proc << G4endl;
761 }
762 if (1 < verbose) {
763 G4cout << " CopyTables for " << proc->GetProcessName()
764 << " for " << part_vector[j]->GetParticleName()
765 << " base_part= " << part->GetParticleName()
766 << " tables are assigned"
767 << G4endl;
768 }
769 }
770 }
771}
772
773//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
774
775G4VEnergyLossProcess* G4LossTableManager::BuildTables(
776 const G4ParticleDefinition* aParticle)
777{
778 if(1 < verbose) {
779 G4cout << " G4LossTableManager::BuildTables(part) for "
780 << aParticle->GetParticleName() << G4endl;
781 }
782
783 std::vector<G4PhysicsTable*> t_list;
784 std::vector<G4VEnergyLossProcess*> loss_list;
785 std::vector<G4bool> build_flags;
786 G4VEnergyLossProcess* em = nullptr;
787 G4VEnergyLossProcess* p = nullptr;
788 G4int iem = 0;
789 G4PhysicsTable* dedx = nullptr;
790 G4int i;
791
792 G4ProcessVector* pvec =
793 aParticle->GetProcessManager()->GetProcessList();
794 G4int nvec = (G4int)pvec->size();
795
796 for (i=0; i<n_loss; ++i) {
797 p = loss_vector[i];
798 if (nullptr != p) {
799 G4bool yes = (aParticle == part_vector[i]);
800
801 // possible case of process sharing between particle/anti-particle
802 if(!yes) {
803 auto ptr = static_cast<G4VProcess*>(p);
804 for(G4int j=0; j<nvec; ++j) {
805 //G4cout << "j= " << j << " " << (*pvec)[j] << " " << ptr << G4endl;
806 if(ptr == (*pvec)[j]) {
807 yes = true;
808 break;
809 }
810 }
811 }
812 // process belong to this particle
813 if(yes && isActive[i]) {
814 if (p->IsIonisationProcess() || !em) {
815 em = p;
816 iem= i;
817 }
818 // tables may be shared between particle/anti-particle
819 G4bool val = false;
820 if (!tables_are_built[i]) {
821 val = true;
822 dedx = p->BuildDEDXTable(fRestricted);
823 //G4cout << "===Build DEDX table for " << p->GetProcessName()
824 // << " idx= " << i << " dedx:" << dedx << " " << dedx->length() << G4endl;
825 p->SetDEDXTable(dedx,fRestricted);
826 tables_are_built[i] = true;
827 } else {
828 dedx = p->DEDXTable();
829 }
830 t_list.push_back(dedx);
831 loss_list.push_back(p);
832 build_flags.push_back(val);
833 }
834 }
835 }
836
837 G4int n_dedx = (G4int)t_list.size();
838 if (0 == n_dedx || !em) {
839 G4cout << "G4LossTableManager WARNING: no DEDX processes for "
840 << aParticle->GetParticleName() << G4endl;
841 return nullptr;
842 }
843 G4int nSubRegions = em->NumberOfSubCutoffRegions();
844
845 if (1 < verbose) {
846 G4cout << " Start to build the sum of " << n_dedx << " processes"
847 << " iem= " << iem << " em= " << em->GetProcessName()
848 << " buildCSDARange= " << theParameters->BuildCSDARange()
849 << " nSubRegions= " << nSubRegions;
850 if(subcutProducer) {
851 G4cout << " SubCutProducer " << subcutProducer->GetName();
852 }
853 G4cout << G4endl;
854 }
855 // do not build tables if producer class is defined
856 if(subcutProducer) { nSubRegions = 0; }
857
858 dedx = em->DEDXTable();
859 em->SetDEDXTable(dedx, fIsIonisation);
860
861 if (1 < n_dedx) {
862 dedx = nullptr;
864 tableBuilder->BuildDEDXTable(dedx, t_list);
865 em->SetDEDXTable(dedx, fRestricted);
866 }
867
868 dedx_vector[iem] = dedx;
869
870 G4PhysicsTable* range = em->RangeTableForLoss();
871 if(!range) range = G4PhysicsTableHelper::PreparePhysicsTable(range);
872 range_vector[iem] = range;
873
874 G4PhysicsTable* invrange = em->InverseRangeTable();
875 if(!invrange) invrange = G4PhysicsTableHelper::PreparePhysicsTable(invrange);
876 inv_range_vector[iem] = invrange;
877
878 tableBuilder->BuildRangeTable(dedx, range);
879 tableBuilder->BuildInverseRangeTable(range, invrange);
880
881 em->SetRangeTableForLoss(range);
882 em->SetInverseRangeTable(invrange);
883
884 std::vector<G4PhysicsTable*> listCSDA;
885
886 for (i=0; i<n_dedx; ++i) {
887 p = loss_list[i];
888 if(build_flags[i]) {
890 }
891 if(theParameters->BuildCSDARange()) {
892 dedx = p->BuildDEDXTable(fTotal);
893 p->SetDEDXTable(dedx,fTotal);
894 listCSDA.push_back(dedx);
895 }
896 }
897
898 if(theParameters->BuildCSDARange()) {
899 G4PhysicsTable* dedxCSDA = em->DEDXunRestrictedTable();
900 if (1 < n_dedx) {
902 tableBuilder->BuildDEDXTable(dedxCSDA, listCSDA);
903 em->SetDEDXTable(dedxCSDA,fTotal);
904 }
905 G4PhysicsTable* rCSDA = em->CSDARangeTable();
906 if(!rCSDA) { rCSDA = G4PhysicsTableHelper::PreparePhysicsTable(rCSDA); }
907 tableBuilder->BuildRangeTable(dedxCSDA, rCSDA);
908 em->SetCSDARangeTable(rCSDA);
909 }
910
911 if (1 < verbose) {
912 G4cout << "G4LossTableManager::BuildTables: Tables are built for "
913 << aParticle->GetParticleName()
914 << "; ionisation process: " << em->GetProcessName()
915 << " " << em
916 << G4endl;
917 }
918 return em;
919}
920
921//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
922
923void G4LossTableManager::ParticleHaveNoLoss(
924 const G4ParticleDefinition* aParticle)
925{
927 ed << "Energy loss process not found for " << aParticle->GetParticleName()
928 << " !";
929 G4Exception("G4LossTableManager::ParticleHaveNoLoss", "em0001",
930 FatalException, ed);
931}
932
933//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
934
936{
937 verbose = val;
938}
939
940//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
941
942const std::vector<G4VEnergyLossProcess*>&
944{
945 return loss_vector;
946}
947
948//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
949
950const std::vector<G4VEmProcess*>& G4LossTableManager::GetEmProcessVector()
951{
952 return emp_vector;
953}
954
955//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
956
957const std::vector<G4VMultipleScattering*>&
959{
960 return msc_vector;
961}
962
963//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
964
966{
967 return theParameters->GetEmSaturation();
968}
969
970//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
971
973{
974 if(!emConfigurator) {
975 emConfigurator = new G4EmConfigurator(verbose);
976 }
977 return emConfigurator;
978}
979
980//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
981
983{
984 if(!emElectronIonPair) {
985 emElectronIonPair = new G4ElectronIonPair(verbose);
986 }
987 return emElectronIonPair;
988}
989
990//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
991
993{
994 if(nullptr != ptr && ptr != nielCalculator) {
995 delete nielCalculator;
996 nielCalculator = ptr;
997 }
998}
999
1000//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
1001
1003{
1004 if(!nielCalculator) {
1005 nielCalculator = new G4NIELCalculator(nullptr, verbose);
1006 }
1007 return nielCalculator;
1008}
1009
1010//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1011
1013{
1014 if(atomDeexcitation != p) {
1015 delete atomDeexcitation;
1016 atomDeexcitation = p;
1017 }
1018}
1019
1020//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1021
1023{
1024 if(subcutProducer != p) {
1025 delete subcutProducer;
1026 subcutProducer = p;
1027 }
1028}
1029
1030//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1031
1032void G4LossTableManager::PrintEWarning(G4String tit, G4double /*val*/)
1033{
1034 G4String ss = "G4LossTableManager::" + tit;
1036 /*
1037 ed << "Parameter is out of range: " << val
1038 << " it will have no effect!\n" << " ## "
1039 << " nbins= " << nbinsLambda
1040 << " nbinsPerDecade= " << nbinsPerDecade
1041 << " Emin(keV)= " << minKinEnergy/keV
1042 << " Emax(GeV)= " << maxKinEnergy/GeV;
1043 */
1044 G4Exception(ss, "em0044", JustWarning, ed);
1045}
1046
1047//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1048
1050{
1051 // Automatic generation of html documentation page for physics lists
1052 // List processes and models for the most important
1053 // particles in descending order of importance
1054 // NB. for model names with length > 18 characters the .rst file needs
1055 // to be edited by hand. Or modify G4EmModelManager::DumpModelList
1056
1057 char* dirName = std::getenv("G4PhysListDocDir");
1058 char* physList = std::getenv("G4PhysListName");
1059 if (dirName && physList) {
1060 G4String physListName = G4String(physList);
1061 G4String pathName = G4String(dirName) + "/" + physListName + ".rst";
1062
1063 std::ofstream outFile;
1064 outFile.open(pathName);
1065
1066 outFile << physListName << G4endl;
1067 outFile << std::string(physListName.length(), '=') << G4endl;
1068
1069 std::vector<G4ParticleDefinition*> particles {
1076 };
1077
1078 std::vector<G4VEmProcess*> emproc_vector = GetEmProcessVector();
1079 std::vector<G4VEnergyLossProcess*> enloss_vector =
1081 std::vector<G4VMultipleScattering*> mscat_vector =
1083
1084 for (auto theParticle : particles) {
1085 outFile << G4endl << "**" << theParticle->GetParticleName()
1086 << "**" << G4endl << G4endl << " .. code-block:: none" << G4endl;
1087
1088 G4ProcessManager* pm = theParticle->GetProcessManager();
1089 G4ProcessVector* pv = pm->GetProcessList();
1090 G4int plen = pm->GetProcessListLength();
1091
1092 for (auto emproc : emproc_vector) {
1093 for (G4int i = 0; i < plen; ++i) {
1094 G4VProcess* proc = (*pv)[i];
1095 if (proc == emproc) {
1096 outFile << G4endl;
1097 proc->ProcessDescription(outFile);
1098 break;
1099 }
1100 }
1101 }
1102
1103 for (auto mscproc : mscat_vector) {
1104 for (G4int i = 0; i < plen; ++i) {
1105 G4VProcess* proc = (*pv)[i];
1106 if (proc == mscproc) {
1107 outFile << G4endl;
1108 proc->ProcessDescription(outFile);
1109 break;
1110 }
1111 }
1112 }
1113
1114 for (auto enlossproc : enloss_vector) {
1115 for (G4int i = 0; i < plen; ++i) {
1116 G4VProcess* proc = (*pv)[i];
1117 if (proc == enlossproc) {
1118 outFile << G4endl;
1119 proc->ProcessDescription(outFile);
1120 break;
1121 }
1122 }
1123 }
1124 }
1125 outFile.close();
1126 }
1127}
1128
1129//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1130
@ fTotal
@ fRestricted
@ fIsIonisation
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4EmParameters * Instance()
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
void SetAtomDeexcitation(G4VAtomDeexcitation *)
static G4LossTableManager * Instance()
const std::vector< G4VEmProcess * > & GetEmProcessVector()
G4VEnergyLossProcess * GetEnergyLossProcess(const G4ParticleDefinition *)
void LocalPhysicsTables(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
const std::vector< G4VMultipleScattering * > & GetMultipleScatteringVector()
void SetVerbose(G4int val)
void DeRegister(G4VEnergyLossProcess *p)
G4NIELCalculator * NIELCalculator()
void SetNIELCalculator(G4NIELCalculator *)
G4EmConfigurator * EmConfigurator()
void Register(G4VEnergyLossProcess *p)
G4ElectronIonPair * ElectronIonPair()
G4EmSaturation * EmSaturation()
const std::vector< G4VEnergyLossProcess * > & GetEnergyLossProcessVector()
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
void SetSubCutProducer(G4VSubCutProducer *)
void RegisterExtraParticle(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
G4LossTableManager(G4LossTableManager &)=delete
static G4MuonMinus * MuonMinusDefinition()
static G4MuonPlus * MuonPlusDefinition()
Definition G4MuonPlus.cc:93
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
static G4Positron * Positron()
Definition G4Positron.cc:90
G4int GetProcessListLength() const
G4ProcessVector * GetProcessList() const
std::size_t size() const
static G4Proton * ProtonDefinition()
Definition G4Proton.cc:85
const G4String & GetName() const
const G4String & GetName() const
const G4ParticleDefinition * BaseParticle() const
G4PhysicsTable * RangeTableForLoss() const
G4PhysicsTable * InverseRangeTable() const
G4PhysicsTable * CSDARangeTable() const
void SetRangeTableForLoss(G4PhysicsTable *p)
G4int NumberOfSubCutoffRegions() const
G4PhysicsTable * BuildDEDXTable(G4EmTableType tType=fRestricted)
const G4ParticleDefinition * Particle() const
void SetInverseRangeTable(G4PhysicsTable *p)
void SetDEDXTable(G4PhysicsTable *p, G4EmTableType tType)
G4PhysicsTable * BuildLambdaTable(G4EmTableType tType=fRestricted)
void SetLambdaTable(G4PhysicsTable *p)
G4PhysicsTable * IonisationTable() const
G4PhysicsTable * LambdaTable() const
void SetCSDARangeTable(G4PhysicsTable *pRange)
G4PhysicsTable * DEDXunRestrictedTable() const
G4PhysicsTable * DEDXTable() const
virtual void ProcessDescription(std::ostream &outfile) const
const G4String & GetProcessName() const
const G4String & GetName() const
#define G4ThreadLocal
Definition tls.hh:77