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