Geant4-11
G4VEnergyLossProcess.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: G4VEnergyLossProcess
32//
33// Author: Vladimir Ivanchenko
34//
35// Creation date: 03.01.2002
36//
37// Modifications: Vladimir Ivanchenko
38//
39//
40// Class Description:
41//
42// It is the unified energy loss process it calculates the continuous
43// energy loss for charged particles using a set of Energy Loss
44// models valid for different energy regions. There are a possibility
45// to create and access to dE/dx and range tables, or to calculate
46// that information on fly.
47// -------------------------------------------------------------------
48//
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
51
54#include "G4SystemOfUnits.hh"
55#include "G4ProcessManager.hh"
56#include "G4LossTableManager.hh"
57#include "G4LossTableBuilder.hh"
58#include "G4Step.hh"
60#include "G4ParticleTable.hh"
61#include "G4EmParameters.hh"
62#include "G4VEmModel.hh"
64#include "G4DataVector.hh"
65#include "G4PhysicsLogVector.hh"
66#include "G4VParticleChange.hh"
68#include "G4Gamma.hh"
69#include "G4Electron.hh"
70#include "G4Positron.hh"
71#include "G4ProcessManager.hh"
72#include "G4UnitsTable.hh"
74#include "G4Region.hh"
75#include "G4RegionStore.hh"
77#include "G4SafetyHelper.hh"
80#include "G4VSubCutProducer.hh"
81#include "G4EmBiasingManager.hh"
82#include "G4Log.hh"
83#include <iostream>
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86
88 G4ProcessType type):
90{
93
94 // low energy limit
96
97 // Size of tables
101 nBins = 84;
102 nBinsCSDA = 35;
103
104 // default linear loss limit
106
107 // particle types
111
112 // run time objects
117 ->GetSafetyHelper();
119
120 // initialise model
122 lManager->Register(this);
123
127
128 scTracks.reserve(10);
129 secParticles.reserve(12);
130}
131
132//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
133
135{
136 /*
137 G4cout << "** G4VEnergyLossProcess::~G4VEnergyLossProcess() for "
138 << GetProcessName() << " isMaster: " << isMaster
139 << " basePart: " << baseParticle
140 << G4endl;
141 G4cout << " isIonisation " << isIonisation << " "
142 << theDEDXTable << " " << theIonisationTable << G4endl;
143 */
144
145 if (isMaster && nullptr == baseParticle) {
146 if(nullptr != theDEDXTable) {
147
148 //G4cout << " theIonisationTable " << theIonisationTable << G4endl;
150 //G4cout << " delete theDEDXTable " << theDEDXTable << G4endl;
152 delete theDEDXTable;
153 theDEDXTable = nullptr;
154 }
155 //G4cout << " theIonisationTable " << theIonisationTable << G4endl;
156 if(nullptr != theIonisationTable) {
157 //G4cout << " delete theIonisationTable " << theIonisationTable << G4endl;
159 delete theIonisationTable;
160 theIonisationTable = nullptr;
161 }
162 if(nullptr != theDEDXunRestrictedTable && isIonisation) {
165 theDEDXunRestrictedTable = nullptr;
166 }
167 if(nullptr != theCSDARangeTable && isIonisation) {
169 delete theCSDARangeTable;
170 theCSDARangeTable = nullptr;
171 }
172 //G4cout << "delete RangeTable: " << theRangeTableForLoss << G4endl;
173 if(nullptr != theRangeTableForLoss && isIonisation) {
176 theRangeTableForLoss = nullptr;
177 }
178 //G4cout << "delete InvRangeTable: " << theInverseRangeTable << G4endl;
179 if(nullptr != theInverseRangeTable && isIonisation /*&& !isIon*/) {
182 theInverseRangeTable = nullptr;
183 }
184 //G4cout << "delete LambdaTable: " << theLambdaTable << G4endl;
185 if(nullptr != theLambdaTable) {
187 delete theLambdaTable;
188 theLambdaTable = nullptr;
189 }
190 if(nullptr != fXSpeaks) {
191 for(auto const & v : *fXSpeaks) { delete v; }
192 delete fXSpeaks;
193 fXSpeaks = nullptr;
194 }
195 }
196 secParticles.clear();
197 delete modelManager;
198 delete biasManager;
199 delete scoffRegions;
200 delete emModels;
201 lManager->DeRegister(this);
202 //G4cout << "** all removed" << G4endl;
203}
204
205//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
206
208 const G4Material*,
209 G4double cut)
210{
211 return cut;
212}
213
214//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
215
218 const G4Region* region)
219{
220 if(nullptr == ptr) { return; }
221 G4VEmFluctuationModel* afluc = (nullptr == fluc) ? fluctModel : fluc;
222 modelManager->AddEmModel(order, ptr, afluc, region);
224}
225
226//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
227
229{
230 if(nullptr == ptr) { return; }
231 if(nullptr == emModels) { emModels = new std::vector<G4VEmModel*>; }
232 if(!emModels->empty()) {
233 for(auto & em : *emModels) { if(em == ptr) { return; } }
234 }
235 emModels->push_back(ptr);
236}
237
238//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
239
241 G4double charge2ratio)
242{
243 massRatio = massratio;
245 fFactor = charge2ratio*biasFactor;
246 if(baseMat) { fFactor *= (*theDensityFactor)[currentCoupleIndex]; }
247 chargeSqRatio = charge2ratio;
249}
250
251//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
252
253void
255{
256 if(1 < verboseLevel) {
257 G4cout << "G4VEnergyLossProcess::PreparePhysicsTable for "
258 << GetProcessName() << " for " << part.GetParticleName()
259 << " " << this << G4endl;
260 }
262
263 // Are particle defined?
264 if(nullptr == particle) { particle = &part; }
265
266 if(part.GetParticleType() == "nucleus") {
267
269 if(pname != "deuteron" && pname != "triton" &&
270 pname != "alpha+" && pname != "alpha") {
271
272 if(nullptr == theGenericIon) {
275 }
276 isIon = true;
277 if(particle != theGenericIon) {
280 size_t n = v->size();
281 for(size_t j=0; j<n; ++j) {
282 if((*v)[j] == this) {
284 break;
285 }
286 }
287 }
288 }
289 }
290
291 if( particle != &part ) {
292 if(!isIon) {
293 lManager->RegisterExtraParticle(&part, this);
294 }
295 if(1 < verboseLevel) {
296 G4cout << "### G4VEnergyLossProcess::PreparePhysicsTable()"
297 << " interrupted for "
298 << part.GetParticleName() << " isIon=" << isIon
299 << " baseMat=" << baseMat
300 << " particle " << particle << " GenericIon " << theGenericIon
301 << G4endl;
302 }
303 return;
304 }
305
306 tablesAreBuilt = false;
307
310
311 // Base particle and set of models can be defined here
313
314 // parameters of the process
328
330
332
333 G4double initialCharge = particle->GetPDGCharge();
334 G4double initialMass = particle->GetPDGMass();
335
337
338 // integral option may be disabled
340
341 // parameters for scaling from the base particle
342 if (nullptr != baseParticle) {
343 massRatio = (baseParticle->GetPDGMass())/initialMass;
345 G4double q = initialCharge/baseParticle->GetPDGCharge();
346 chargeSqRatio = q*q;
348 }
349 lowestKinEnergy = (initialMass < CLHEP::MeV)
352
353 // Tables preparation
354 if (isMaster && nullptr == baseParticle) {
355
356 if(nullptr != theDEDXTable && isIonisation) {
359 delete theDEDXTable;
361 }
362 }
363
366
372 }
373
375
376 if(isIonisation) {
381 }
382
383 if(fXSType == fEmTwoPeaks) {
384 const G4ProductionCutsTable* theCoupleTable=
386 size_t n = theCoupleTable->GetTableSize();
387 if(nullptr == fXSpeaks) {
388 fXSpeaks = new std::vector<G4TwoPeaksXS*>;
389 }
390 fXSpeaks->resize(n, nullptr);
391 }
392 }
393 /*
394 G4cout << "** G4VEnergyLossProcess::PreparePhysicsTable() for "
395 << GetProcessName() << " and " << particle->GetParticleName()
396 << " isMaster: " << isMaster << " isIonisation: "
397 << isIonisation << G4endl;
398 G4cout << " theDEDX: " << theDEDXTable
399 << " theRange: " << theRangeTableForLoss
400 << " theInverse: " << theInverseRangeTable
401 << " theLambda: " << theLambdaTable << G4endl;
402 */
403 // forced biasing
404 if(nullptr != biasManager) {
406 biasFlag = false;
407 }
408
409 // defined ID of secondary particles
410 G4int stype = GetProcessSubType();
411 if(stype == fBremsstrahlung) {
414 } else if(stype == fPairProdByCharged) {
416 mainSecondaries = 2;
417 }
419
420 // initialisation of models
422 for(G4int i=0; i<numberOfModels; ++i) {
424 if(0 == i) { currentModel = mod; }
428 if(mod->HighEnergyLimit() > maxKinEnergy) {
430 }
432 SetEmModel(mod);
433 }
435 1.0, verboseLevel);
436
437 // subcut processor
438 if(isIonisation) {
440 }
441 if(1 == nSCoffRegions) {
442 if((*scoffRegions)[0]->GetName() == "DefaultRegionForTheWorld") {
443 delete scoffRegions;
444 scoffRegions = nullptr;
445 nSCoffRegions = 0;
446 }
447 }
448
449 if(1 < verboseLevel) {
450 G4cout << "G4VEnergyLossProcess::PrepearPhysicsTable() is done "
451 << " for local " << particle->GetParticleName()
452 << " isIon= " << isIon;
453 if(baseParticle) {
454 G4cout << "; base: " << baseParticle->GetParticleName();
455 }
456 G4cout << " chargeSqRatio= " << chargeSqRatio
457 << " massRatio= " << massRatio
458 << " reduceFactor= " << reduceFactor << G4endl;
459 if (nSCoffRegions > 0) {
460 G4cout << " SubCut secondary production is ON for regions: " << G4endl;
461 for (G4int i=0; i<nSCoffRegions; ++i) {
462 const G4Region* r = (*scoffRegions)[i];
463 G4cout << " " << r->GetName() << G4endl;
464 }
465 } else if(nullptr != subcutProducer) {
466 G4cout << " SubCut secondary production is ON for all regions" << G4endl;
467 }
468 }
469}
470
471//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
472
474{
475 if(1 < verboseLevel) {
476 G4cout << "### G4VEnergyLossProcess::BuildPhysicsTable() for "
477 << GetProcessName()
478 << " and particle " << part.GetParticleName()
479 << "; local: " << particle->GetParticleName();
480 if(baseParticle) {
481 G4cout << "; base: " << baseParticle->GetParticleName();
482 }
483 G4cout << " TablesAreBuilt= " << tablesAreBuilt
484 << " isIon= " << isIon << " " << this << G4endl;
485 }
486
487 if(&part == particle) {
488
489 if(isMaster) {
491
492 } else {
493
494 const G4VEnergyLossProcess* masterProcess =
495 static_cast<const G4VEnergyLossProcess*>(GetMasterProcess());
496
497 // copy table pointers from master thread
498 SetDEDXTable(masterProcess->DEDXTable(),fRestricted);
499 SetDEDXTable(masterProcess->DEDXunRestrictedTable(),fTotal);
501 SetRangeTableForLoss(masterProcess->RangeTableForLoss());
502 SetCSDARangeTable(masterProcess->CSDARangeTable());
504 SetInverseRangeTable(masterProcess->InverseRangeTable());
505 SetLambdaTable(masterProcess->LambdaTable());
506 SetTwoPeaksXS(masterProcess->TwoPeaksXS());
507 isIonisation = masterProcess->IsIonisationProcess();
508 baseMat = masterProcess->UseBaseMaterial();
509
510 tablesAreBuilt = true;
511 // local initialisation of models
512 G4bool printing = true;
513 for(G4int i=0; i<numberOfModels; ++i) {
514 G4VEmModel* mod = GetModelByIndex(i, printing);
515 G4VEmModel* mod0= masterProcess->GetModelByIndex(i, printing);
517 mod->InitialiseLocal(particle, mod0);
518 }
520 }
521
522 // needs to be done only once
524 }
525 // Added tracking cut to avoid tracking artifacts
526 // and identified deexcitation flag
527 if(isIonisation) {
529 if(nullptr != atomDeexcitation) {
531 }
532 }
533
534 // protection against double printout
535 if(theParameters->IsPrintLocked()) { return; }
536
537 // explicitly defined printout by particle name
538 G4String num = part.GetParticleName();
539 if(1 < verboseLevel ||
540 (0 < verboseLevel && (num == "e-" ||
541 num == "e+" || num == "mu+" ||
542 num == "mu-" || num == "proton"||
543 num == "pi+" || num == "pi-" ||
544 num == "kaon+" || num == "kaon-" ||
545 num == "alpha" || num == "anti_proton" ||
546 num == "GenericIon"|| num == "alpha+" )))
547 {
548 StreamInfo(G4cout, part);
549 }
550
551 /*
552 G4cout << "** G4VEnergyLossProcess::BuildPhysicsTable() for "
553 << GetProcessName() << " and " << particle->GetParticleName()
554 << " isMaster: " << isMaster << " isIonisation: "
555 << isIonisation << G4endl;
556 G4cout << " theDEDX: " << theDEDXTable
557 << " theRange: " << theRangeTableForLoss
558 << " theInverse: " << theInverseRangeTable
559 << " theLambda: " << theLambdaTable << G4endl;
560 */
561 //if(1 < verboseLevel || verb) {
562 if(1 < verboseLevel) {
563 G4cout << "### G4VEnergyLossProcess::BuildPhysicsTable() done for "
564 << GetProcessName()
565 << " and particle " << part.GetParticleName();
566 if(isIonisation) { G4cout << " isIonisation flag=1"; }
567 G4cout << " baseMat=" << baseMat << G4endl;
568 }
569}
570
571//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
572
574{
575 if(1 < verboseLevel ) {
576 G4cout << "G4VEnergyLossProcess::BuildDEDXTable() of type " << tType
577 << " for " << GetProcessName()
578 << " and particle " << particle->GetParticleName()
579 << G4endl;
580 }
581 G4PhysicsTable* table = nullptr;
583 G4int bin = nBins;
584
585 if(fTotal == tType) {
587 bin = nBinsCSDA;
589 } else if(fRestricted == tType) {
590 table = theDEDXTable;
591 } else {
592 G4cout << "G4VEnergyLossProcess::BuildDEDXTable WARNING: wrong type "
593 << tType << G4endl;
594 }
595
596 // Access to materials
597 const G4ProductionCutsTable* theCoupleTable=
599 size_t numOfCouples = theCoupleTable->GetTableSize();
600
601 if(1 < verboseLevel) {
602 G4cout << numOfCouples << " materials"
603 << " minKinEnergy= " << minKinEnergy
604 << " maxKinEnergy= " << emax
605 << " nbin= " << bin
606 << " EmTableType= " << tType
607 << " table= " << table << " " << this
608 << G4endl;
609 }
610 if(nullptr == table) { return table; }
611
613 G4PhysicsLogVector* aVector = nullptr;
614 G4PhysicsLogVector* bVector = nullptr;
615
616 for(size_t i=0; i<numOfCouples; ++i) {
617
618 if(1 < verboseLevel) {
619 G4cout << "G4VEnergyLossProcess::BuildDEDXVector Idx= " << i
620 << " flagTable= " << table->GetFlag(i)
621 << " Flag= " << bld->GetFlag(i) << G4endl;
622 }
623 if(bld->GetFlag(i)) {
624
625 // create physics vector and fill it
626 const G4MaterialCutsCouple* couple =
627 theCoupleTable->GetMaterialCutsCouple(i);
628 if(nullptr != (*table)[i]) { delete (*table)[i]; }
629 if(nullptr != bVector) {
630 aVector = new G4PhysicsLogVector(*bVector);
631 } else {
632 bVector = new G4PhysicsLogVector(minKinEnergy, emax, bin, spline);
633 aVector = bVector;
634 }
635
636 modelManager->FillDEDXVector(aVector, couple, tType);
637 if(spline) { aVector->FillSecondDerivatives(); }
638
639 // Insert vector for this material into the table
641 }
642 }
643
644 if(1 < verboseLevel) {
645 G4cout << "G4VEnergyLossProcess::BuildDEDXTable(): table is built for "
647 << " and process " << GetProcessName()
648 << G4endl;
649 if(2 < verboseLevel) G4cout << (*table) << G4endl;
650 }
651
652 return table;
653}
654
655//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
656
658{
659 G4PhysicsTable* table = nullptr;
660
661 if(fRestricted == tType) {
662 table = theLambdaTable;
663 } else {
664 G4cout << "G4VEnergyLossProcess::BuildLambdaTable WARNING: wrong type "
665 << tType << G4endl;
666 }
667
668 if(1 < verboseLevel) {
669 G4cout << "G4VEnergyLossProcess::BuildLambdaTable() of type "
670 << tType << " for process "
671 << GetProcessName() << " and particle "
673 << " EmTableType= " << tType
674 << " table= " << table
675 << G4endl;
676 }
677 if(nullptr == table) { return table; }
678
679 // Access to materials
680 const G4ProductionCutsTable* theCoupleTable=
682 size_t numOfCouples = theCoupleTable->GetTableSize();
684
685 G4PhysicsLogVector* aVector = nullptr;
687
688 for(size_t i=0; i<numOfCouples; ++i) {
689
690 if (bld->GetFlag(i)) {
691
692 // create physics vector and fill it
693 const G4MaterialCutsCouple* couple =
694 theCoupleTable->GetMaterialCutsCouple(i);
695 delete (*table)[i];
696
697 G4bool startNull = true;
698 G4double emin =
700 if(minKinEnergy > emin) {
701 emin = minKinEnergy;
702 startNull = false;
703 }
704
706 if(emax <= emin) { emax = 2*emin; }
707 G4int bin = G4lrint(nBins*G4Log(emax/emin)/scale);
708 bin = std::max(bin, 3);
709 aVector = new G4PhysicsLogVector(emin, emax, bin, spline);
710
711 modelManager->FillLambdaVector(aVector, couple, startNull, tType);
712 if(spline) { aVector->FillSecondDerivatives(); }
713
714 // Insert vector for this material into the table
716 }
717 }
718
719 if(1 < verboseLevel) {
720 G4cout << "Lambda table is built for "
722 << G4endl;
723 }
724
725 return table;
726}
727
728//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
729
730void G4VEnergyLossProcess::StreamInfo(std::ostream& out,
731 const G4ParticleDefinition& part, G4bool rst) const
732{
733 G4String indent = (rst ? " " : "");
734 out << std::setprecision(6);
735 out << G4endl << indent << GetProcessName() << ": ";
736 if (!rst) out << " for " << part.GetParticleName();
737 out << " XStype:" << fXSType
738 << " SubType=" << GetProcessSubType() << G4endl
739 << " dE/dx and range tables from "
740 << G4BestUnit(minKinEnergy,"Energy")
741 << " to " << G4BestUnit(maxKinEnergy,"Energy")
742 << " in " << nBins << " bins" << G4endl
743 << " Lambda tables from threshold to "
744 << G4BestUnit(maxKinEnergy,"Energy")
746 << " bins/decade, spline: " << spline
747 << G4endl;
748 if(nullptr != theRangeTableForLoss && isIonisation) {
749 out << " StepFunction=(" << dRoverRange << ", "
750 << finalRange/mm << " mm)"
751 << ", integ: " << fXSType
752 << ", fluct: " << lossFluctuationFlag
753 << ", linLossLim= " << linLossLimit
754 << G4endl;
755 }
758 if(nullptr != theCSDARangeTable && isIonisation) {
759 out << " CSDA range table up"
760 << " to " << G4BestUnit(maxKinEnergyCSDA,"Energy")
761 << " in " << nBinsCSDA << " bins" << G4endl;
762 }
763 if(nSCoffRegions>0 && isIonisation) {
764 out << " Subcutoff sampling in " << nSCoffRegions
765 << " regions" << G4endl;
766 }
767 if(2 < verboseLevel) {
768 out << " DEDXTable address= " << theDEDXTable << G4endl;
769 if(nullptr != theDEDXTable && isIonisation)
770 out << (*theDEDXTable) << G4endl;
771 out << "non restricted DEDXTable address= "
773 if(nullptr != theDEDXunRestrictedTable && isIonisation) {
774 out << (*theDEDXunRestrictedTable) << G4endl;
775 }
776 out << " CSDARangeTable address= " << theCSDARangeTable << G4endl;
777 if(nullptr != theCSDARangeTable && isIonisation) {
778 out << (*theCSDARangeTable) << G4endl;
779 }
780 out << " RangeTableForLoss address= " << theRangeTableForLoss
781 << G4endl;
782 if(nullptr != theRangeTableForLoss && isIonisation) {
783 out << (*theRangeTableForLoss) << G4endl;
784 }
785 out << " InverseRangeTable address= " << theInverseRangeTable
786 << G4endl;
787 if(nullptr != theInverseRangeTable && isIonisation) {
788 out << (*theInverseRangeTable) << G4endl;
789 }
790 out << " LambdaTable address= " << theLambdaTable << G4endl;
791 if(nullptr != theLambdaTable) {
792 out << (*theLambdaTable) << G4endl;
793 }
794 }
795}
796
797//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
798
800{
801 if(nullptr == scoffRegions) {
802 scoffRegions = new std::vector<const G4Region*>;
803 }
804 // the region is in the list
805 if(!scoffRegions->empty()) {
806 for (auto & reg : *scoffRegions) {
807 if (reg == r) { return; }
808 }
809 }
810 // new region
811 scoffRegions->push_back(r);
813}
814
815//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
816
818{
819 if(0 == nSCoffRegions) { return true; }
820 const G4Region* r = aTrack.GetVolume()->GetLogicalVolume()->GetRegion();
821 for(auto & reg : *scoffRegions) {
822 if(r == reg) { return true; }
823 }
824 return false;
825}
826
827//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
828
830{
831 /*
832 G4cout << "G4VEnergyLossProcess::StartTracking: "
833 << track->GetDefinition()->GetParticleName()
834 << " e(MeV)= " << track->GetKineticEnergy();
835 if(particle) G4cout << " " << particle->GetParticleName();
836 if(baseParticle) G4cout << " basePart: " << baseParticle->GetParticleName();
837 G4cout << " " << GetProcessName();
838 if(isIon) G4cout << " isIon: Q=" << track->GetDefinition()->GetPDGCharge()
839 << " Qdyn=" << track->GetDynamicParticle()->GetCharge();
840 G4cout << G4endl;
841 */
842 // reset parameters for the new track
845 currentCouple = nullptr;
846
847 // reset ion
848 if(isIon) {
849 const G4double newmass = track->GetDefinition()->GetPDGMass();
850 if(nullptr != baseParticle) {
851 massRatio = baseParticle->GetPDGMass()/newmass;
853 } else if(nullptr != theGenericIon) {
856 } else {
857 massRatio = 1.0;
858 logMassRatio = 0.0;
859 }
860 }
861 // forced biasing only for primary particles
862 if(biasManager) {
863 if(0 == track->GetParentID()) {
864 biasFlag = true;
866 }
867 }
868}
869
870//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
871
874 G4GPILSelection* selection)
875{
876 G4double x = DBL_MAX;
877 *selection = aGPILSelection;
880 const G4double finR = (rndmStepFlag) ? std::min(finalRange,
882 x = (fRange > finR) ?
883 fRange*dRoverRange + finR*(1.0-dRoverRange)*(2.0-finR/fRange) : fRange;
884 // if(particle->GetPDGMass() > 0.9*GeV)
885 /*
886 G4cout << GetProcessName() << ": e= " << preStepKinEnergy
887 <<" range= "<<fRange << " idx= " << basedCoupleIndex
888 << " finR= " << finR << " limit= " << x <<
889 << "\n" << "massRatio= " << massRatio << " Q^2= " << chargeSqRatio
890 << " dRoverRange= " << dRoverRange
891 << " finalRange= " << finalRange << G4endl;
892 */
893 }
894 //G4cout<<"AlongStepGPIL: " << GetProcessName()<<": e= "<<preStepKinEnergy
895 //<<" stepLimit= "<<x<<G4endl;
896 return x;
897}
898
899//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
900
902 const G4Track& track,
903 G4double previousStepSize,
905{
906 // condition is set to "Not Forced"
908 G4double x = DBL_MAX;
909
910 // initialisation of material, mass, charge, model
911 // at the beginning of the step
918
922 return x;
923 }
924
925 // change effective charge of a charged particle on fly
926 if(isIon) {
927 const G4double q2 = currentModel->ChargeSquareRatio(track);
928 if(q2 != chargeSqRatio) {
931 chargeSqRatio = q2;
932 // G4cout << "PostStepGPIL: Q^2=" << chargeSqRatio << " reducedFactor=" << reduceFactor << G4endl;
933 }
934 }
935
936 // forced biasing only for primary particles
937 if(biasManager) {
938 if(0 == track.GetParentID() && biasFlag &&
940 return biasManager->GetStepLimit(currentCoupleIndex, previousStepSize);
941 }
942 }
943
944 // compute mean free path
946
947 // zero cross section
948 if(preStepLambda <= 0.0) {
951 } else {
952
953 // non-zero cross section
955
956 // beggining of tracking (or just after DoIt of this process)
959
960 } else if(currentInteractionLength < DBL_MAX) {
961
962 // subtract NumberOfInteractionLengthLeft using previous step
964 previousStepSize/currentInteractionLength;
965
968 }
969
970 // new mean free path and step limit
973 }
974#ifdef G4VERBOSE
975 if (verboseLevel>2){
976 // if(particle->GetPDGMass() > 0.9*GeV){
977 G4cout << "G4VEnergyLossProcess::PostStepGetPhysicalInteractionLength ";
978 G4cout << "[ " << GetProcessName() << "]" << G4endl;
979 G4cout << " for " << track.GetDefinition()->GetParticleName()
980 << " in Material " << currentMaterial->GetName()
981 << " Ekin(MeV)= " << preStepKinEnergy/MeV
982 << " " << track.GetMaterial()->GetName()
983 <<G4endl;
984 G4cout << "MeanFreePath = " << currentInteractionLength/cm << "[cm]"
985 << "InteractionLength= " << x/cm <<"[cm] " <<G4endl;
986 }
987#endif
988 return x;
989}
990
991//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
992
993void
995{
996 // cross section increased with energy
997 if(fXSType == fEmIncreasing) {
998 if(e/lambdaFactor < mfpKinEnergy) {
999 mfpKinEnergy = e;
1001 }
1002
1003 // cross section has two peaks
1004 } else if(fXSType == fEmTwoPeaks) {
1005 G4TwoPeaksXS* xs = (*fXSpeaks)[basedCoupleIndex];
1006 const G4double e1peak = xs->e1peak;
1007
1008 // below the 1st peak
1009 if(e <= e1peak) {
1010 if(e/lambdaFactor < mfpKinEnergy) {
1011 mfpKinEnergy = e;
1013 }
1014 return;
1015 }
1016 const G4double e1deep = xs->e1deep;
1017 // above the 1st peak, below the deep
1018 if(e <= e1deep) {
1019 if(mfpKinEnergy >= e1deep || e <= mfpKinEnergy) {
1020 const G4double e1 = std::max(e1peak, e*lambdaFactor);
1022 mfpKinEnergy = e1;
1023 }
1024 return;
1025 }
1026 const G4double e2peak = xs->e2peak;
1027 // above the deep, below 2nd peak
1028 if(e <= e2peak) {
1029 if(e/lambdaFactor < mfpKinEnergy) {
1030 mfpKinEnergy = e;
1032 }
1033 return;
1034 }
1035 const G4double e2deep = xs->e2deep;
1036 // above the 2nd peak, below the deep
1037 if(e <= e2deep) {
1038 if(mfpKinEnergy >= e2deep || e <= mfpKinEnergy) {
1039 const G4double e1 = std::max(e2peak, e*lambdaFactor);
1041 mfpKinEnergy = e1;
1042 }
1043 return;
1044 }
1045 // above the deep, below 3d peak
1046 if(e/lambdaFactor < mfpKinEnergy) {
1047 mfpKinEnergy = e;
1049 }
1050
1051 // integral method is not used
1052 } else {
1054 }
1055}
1056
1057//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1058
1060 const G4Step& step)
1061{
1063 // The process has range table - calculate energy loss
1065 return &fParticleChange;
1066 }
1067
1068 // Get the actual (true) Step length
1069 G4double length = step.GetStepLength();
1070 if(length <= 0.0) { return &fParticleChange; }
1071 G4double eloss = 0.0;
1072
1073 /*
1074 if(-1 < verboseLevel) {
1075 const G4ParticleDefinition* d = track.GetParticleDefinition();
1076 G4cout << "AlongStepDoIt for "
1077 << GetProcessName() << " and particle "
1078 << d->GetParticleName()
1079 << " eScaled(MeV)= " << preStepScaledEnergy/MeV
1080 << " range(mm)= " << fRange/mm
1081 << " s(mm)= " << length/mm
1082 << " rf= " << reduceFactor
1083 << " q^2= " << chargeSqRatio
1084 << " md= " << d->GetPDGMass()
1085 << " status= " << track.GetTrackStatus()
1086 << " " << track.GetMaterial()->GetName()
1087 << G4endl;
1088 }
1089 */
1090
1091 const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
1092
1093 // define new weight for primary and secondaries
1095 if(weightFlag) {
1096 weight /= biasFactor;
1098 }
1099
1100 // stopping
1101 if (length >= fRange || preStepKinEnergy <= lowestKinEnergy) {
1102 eloss = preStepKinEnergy;
1103 if (useDeexcitation) {
1105 eloss, currentCoupleIndex);
1106 if(scTracks.size() > 0) { FillSecondariesAlongStep(weight); }
1107 eloss = std::max(eloss, 0.0);
1108 }
1111 return &fParticleChange;
1112 }
1113 //G4cout << theDEDXTable << " idx= " << basedCoupleIndex
1114 // << " " << GetProcessName() << " "<< currentMaterial->GetName()<<G4endl;
1115 //if(particle->GetParticleName() == "e-")G4cout << (*theDEDXTable) <<G4endl;
1116 // Short step
1119
1120 //G4cout << "Short STEP: eloss= " << eloss << G4endl;
1121
1122 // Long step
1123 if(eloss > preStepKinEnergy*linLossLimit) {
1124
1125 G4double x = (fRange - length)/reduceFactor;
1126 //G4cout << "x= " << x << " " << theInverseRangeTable << G4endl;
1128
1129 /*
1130 if(-1 < verboseLevel)
1131 G4cout << "Long STEP: rPre(mm)= "
1132 << GetScaledRangeForScaledEnergy(preStepScaledEnergy)/mm
1133 << " rPost(mm)= " << x/mm
1134 << " ePre(MeV)= " << preStepScaledEnergy/MeV
1135 << " eloss(MeV)= " << eloss/MeV
1136 << " eloss0(MeV)= "
1137 << GetDEDXForScaledEnergy(preStepScaledEnergy)*length/MeV
1138 << " lim(MeV)= " << preStepKinEnergy*linLossLimit/MeV
1139 << G4endl;
1140 */
1141 }
1142
1143 /*
1144 G4double eloss0 = eloss;
1145 if(-1 < verboseLevel ) {
1146 G4cout << "Before fluct: eloss(MeV)= " << eloss/MeV
1147 << " e-eloss= " << preStepKinEnergy-eloss
1148 << " step(mm)= " << length/mm
1149 << " range(mm)= " << fRange/mm
1150 << " fluct= " << lossFluctuationFlag
1151 << G4endl;
1152 }
1153 */
1154
1155 const G4double cut = (*theCuts)[currentCoupleIndex];
1156 G4double esec = 0.0;
1157
1158 // Corrections, which cannot be tabulated
1159 if(isIon) {
1161 length, eloss);
1162 eloss = std::max(eloss, 0.0);
1163 }
1164
1165 // Sample fluctuations if not full energy loss
1166 if(eloss >= preStepKinEnergy) {
1167 eloss = preStepKinEnergy;
1168
1169 } else if (lossFluctuationFlag) {
1170 const G4double tmax = currentModel->MaxSecondaryKinEnergy(dynParticle);
1171 const G4double tcut = std::min(cut, tmax);
1173 eloss = fluc->SampleFluctuations(currentCouple,dynParticle,
1174 tcut, tmax, length, eloss);
1175 /*
1176 if(-1 < verboseLevel)
1177 G4cout << "After fluct: eloss(MeV)= " << eloss/MeV
1178 << " fluc= " << (eloss-eloss0)/MeV
1179 << " ChargeSqRatio= " << chargeSqRatio
1180 << " massRatio= " << massRatio
1181 << " tmax= " << tmax
1182 << G4endl;
1183 */
1184 }
1185
1186 // deexcitation
1187 if (useDeexcitation) {
1188 G4double esecfluo = preStepKinEnergy;
1189 G4double de = esecfluo;
1190 //G4double eloss0 = eloss;
1191 /*
1192 G4cout << "### 1: E(keV)= " << preStepKinEnergy/keV
1193 << " Efluomax(keV)= " << de/keV
1194 << " Eloss(keV)= " << eloss/keV << G4endl;
1195 */
1197 de, currentCoupleIndex);
1198
1199 // sum of de-excitation energies
1200 esecfluo -= de;
1201
1202 // subtracted from energy loss
1203 if(eloss >= esecfluo) {
1204 esec += esecfluo;
1205 eloss -= esecfluo;
1206 } else {
1207 esec += esecfluo;
1208 eloss = 0.0;
1209 }
1210 /*
1211 if(esecfluo > 0.0) {
1212 G4cout << "### 2: E(keV)= " << preStepKinEnergy/keV
1213 << " Esec(keV)= " << esec/keV
1214 << " Esecf(kV)= " << esecfluo/keV
1215 << " Eloss0(kV)= " << eloss0/keV
1216 << " Eloss(keV)= " << eloss/keV
1217 << G4endl;
1218 }
1219 */
1220 }
1221 if(nullptr != subcutProducer && IsRegionForCubcutProcessor(track)) {
1222 subcutProducer->SampleSecondaries(step, scTracks, eloss, cut);
1223 }
1224 // secondaries from atomic de-excitation and subcut
1225 if(!scTracks.empty()) { FillSecondariesAlongStep(weight); }
1226
1227 // Energy balance
1228 G4double finalT = preStepKinEnergy - eloss - esec;
1229 if (finalT <= lowestKinEnergy) {
1230 eloss += finalT;
1231 finalT = 0.0;
1232 } else if(isIon) {
1235 currentMaterial,finalT));
1236 }
1237 eloss = std::max(eloss, 0.0);
1238
1241 /*
1242 if(-1 < verboseLevel) {
1243 G4double del = finalT + eloss + esec - preStepKinEnergy;
1244 G4cout << "Final value eloss(MeV)= " << eloss/MeV
1245 << " preStepKinEnergy= " << preStepKinEnergy
1246 << " postStepKinEnergy= " << finalT
1247 << " de(keV)= " << del/keV
1248 << " lossFlag= " << lossFluctuationFlag
1249 << " status= " << track.GetTrackStatus()
1250 << G4endl;
1251 }
1252 */
1253 return &fParticleChange;
1254}
1255
1256//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1257
1259{
1260 const G4int n0 = scTracks.size();
1261 G4double weight = wt;
1262 // weight may be changed by biasing manager
1263 if(biasManager) {
1265 weight *=
1267 }
1268 }
1269
1270 // fill secondaries
1271 const G4int n = scTracks.size();
1273
1274 for(G4int i=0; i<n; ++i) {
1275 G4Track* t = scTracks[i];
1276 if(nullptr != t) {
1277 t->SetWeight(weight);
1279 if(i >= n0) { t->SetCreatorModelID(biasID); }
1280 //G4cout << "Secondary(along step) has weight " << t->GetWeight()
1281 //<< ", kenergy " << t->GetKineticEnergy()/MeV << " MeV" <<G4endl;
1282 }
1283 }
1284 scTracks.clear();
1285}
1286
1287//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1288
1290 const G4Step& step)
1291{
1292 // In all cases clear number of interaction lengths
1295
1297 const G4double finalT = track.GetKineticEnergy();
1298
1299 const G4double postStepScaledEnergy = finalT*massRatio;
1300 SelectModel(postStepScaledEnergy);
1301
1302 if(!currentModel->IsActive(postStepScaledEnergy)) {
1303 return &fParticleChange;
1304 }
1305 /*
1306 if(-1 < verboseLevel) {
1307 G4cout << GetProcessName()
1308 << "::PostStepDoIt: E(MeV)= " << finalT/MeV
1309 << G4endl;
1310 }
1311 */
1312
1313 // forced process - should happen only once per track
1314 if(biasFlag) {
1316 biasFlag = false;
1317 }
1318 }
1319
1320 const G4DynamicParticle* dp = track.GetDynamicParticle();
1321
1322 // Integral approach
1323 if (fXSType != fEmNoIntegral) {
1324 const G4double logFinalT = dp->GetLogKineticEnergy();
1325 G4double lx = GetLambdaForScaledEnergy(postStepScaledEnergy,
1326 logFinalT + logMassRatio);
1327 lx = std::max(lx, 0.0);
1328
1329 // cache cross section useful for the false interaction
1330 const G4double lg = preStepLambda;
1331 if(postStepScaledEnergy < mfpKinEnergy) {
1332 mfpKinEnergy = postStepScaledEnergy;
1333 preStepLambda = lx;
1334 }
1335 /*
1336 if(preStepLambda<lx && 1 < verboseLevel) {
1337 G4cout << "WARNING: for " << particle->GetParticleName()
1338 << " and " << GetProcessName()
1339 << " E(MeV)= " << finalT/MeV
1340 << " preLambda= " << preStepLambda
1341 << " < " << lx << " (postLambda) "
1342 << G4endl;
1343 }
1344 */
1345 // if both lg and lx are zero then no interaction
1346 if(lg*G4UniformRand() >= lx) {
1347 return &fParticleChange;
1348 }
1349 }
1350 // define new weight for primary and secondaries
1352 if(weightFlag) {
1353 weight /= biasFactor;
1355 }
1356
1357 const G4double tcut = (*theCuts)[currentCoupleIndex];
1358
1359 // sample secondaries
1360 secParticles.clear();
1361 //G4cout<< "@@@ Eprimary= "<<dynParticle->GetKineticEnergy()/MeV
1362 // << " cut= " << tcut/MeV << G4endl;
1364
1365 const G4int num0 = secParticles.size();
1366
1367 // bremsstrahlung splitting or Russian roulette
1368 if(biasManager) {
1370 G4double eloss = 0.0;
1373 track, currentModel,
1374 &fParticleChange, eloss,
1375 currentCoupleIndex, tcut,
1376 step.GetPostStepPoint()->GetSafety());
1377 if(eloss > 0.0) {
1380 }
1381 }
1382 }
1383
1384 // save secondaries
1385 const G4int num = secParticles.size();
1386 if(num > 0) {
1387
1389 G4double time = track.GetGlobalTime();
1390
1391 G4int n1(0), n2(0);
1392 if(num0 > mainSecondaries) {
1394 }
1395
1396 for (G4int i=0; i<num; ++i) {
1397 if(nullptr != secParticles[i]) {
1398 G4Track* t = new G4Track(secParticles[i], time, track.GetPosition());
1400 if (biasManager) {
1401 t->SetWeight(weight * biasManager->GetWeight(i));
1402 } else {
1403 t->SetWeight(weight);
1404 }
1405 if(i < num0) {
1407 } else if(i < num0 + n1) {
1409 } else {
1411 }
1412
1413 //G4cout << "Secondary(post step) has weight " << t->GetWeight()
1414 // << ", kenergy " << t->GetKineticEnergy()/MeV << " MeV"
1415 // << " time= " << time/ns << " ns " << G4endl;
1417 }
1418 }
1419 }
1420
1426 }
1427
1428 /*
1429 if(-1 < verboseLevel) {
1430 G4cout << "::PostStepDoIt: Sample secondary; Efin= "
1431 << fParticleChange.GetProposedKineticEnergy()/MeV
1432 << " MeV; model= (" << currentModel->LowEnergyLimit()
1433 << ", " << currentModel->HighEnergyLimit() << ")"
1434 << " preStepLambda= " << preStepLambda
1435 << " dir= " << track.GetMomentumDirection()
1436 << " status= " << track.GetTrackStatus()
1437 << G4endl;
1438 }
1439 */
1440 return &fParticleChange;
1441}
1442
1443//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1444
1446 const G4ParticleDefinition* part, const G4String& directory,
1447 G4bool ascii)
1448{
1449 G4bool res = true;
1450 //G4cout << "G4VEnergyLossProcess::StorePhysicsTable: " << part->GetParticleName()
1451 // << " " << directory << " " << ascii << G4endl;
1452 if (!isMaster || baseParticle || part != particle ) return res;
1453
1454 if(!StoreTable(part,theDEDXTable,ascii,directory,"DEDX"))
1455 {res = false;}
1456
1457 if(!StoreTable(part,theDEDXunRestrictedTable,ascii,directory,"DEDXnr"))
1458 {res = false;}
1459
1460 if(!StoreTable(part,theIonisationTable,ascii,directory,"Ionisation"))
1461 {res = false;}
1462
1463 if(isIonisation &&
1464 !StoreTable(part,theCSDARangeTable,ascii,directory,"CSDARange"))
1465 {res = false;}
1466
1467 if(isIonisation &&
1468 !StoreTable(part,theRangeTableForLoss,ascii,directory,"Range"))
1469 {res = false;}
1470
1471 if(isIonisation &&
1472 !StoreTable(part,theInverseRangeTable,ascii,directory,"InverseRange"))
1473 {res = false;}
1474
1475 if(!StoreTable(part,theLambdaTable,ascii,directory,"Lambda"))
1476 {res = false;}
1477
1478 return res;
1479}
1480
1481//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1482
1483G4bool
1485 const G4String& directory,
1486 G4bool ascii)
1487{
1488 G4bool res = true;
1489 if (!isMaster) return res;
1490 const G4String& particleName = part->GetParticleName();
1491
1492 if(1 < verboseLevel) {
1493 G4cout << "G4VEnergyLossProcess::RetrievePhysicsTable() for "
1494 << particleName << " and process " << GetProcessName()
1495 << "; tables_are_built= " << tablesAreBuilt
1496 << G4endl;
1497 }
1498 if(particle == part) {
1499
1500 if(nullptr == baseParticle) {
1501
1502 G4bool fpi = true;
1503 if(!RetrieveTable(part,theDEDXTable,ascii,directory,"DEDX",fpi))
1504 { fpi = false; }
1505
1506 // ionisation table keeps individual dEdx and not sum of sub-processes
1507 if(!RetrieveTable(part,theDEDXTable,ascii,directory,"Ionisation",false))
1508 { fpi = false; }
1509
1510 if(!RetrieveTable(part,theRangeTableForLoss,ascii,directory,"Range",fpi))
1511 { res = false; }
1512
1513 if(!RetrieveTable(part,theDEDXunRestrictedTable,ascii,directory,
1514 "DEDXnr",false))
1515 { res = false; }
1516
1517 if(!RetrieveTable(part,theCSDARangeTable,ascii,directory,
1518 "CSDARange",false))
1519 { res = false; }
1520
1521 if(!RetrieveTable(part,theInverseRangeTable,ascii,directory,
1522 "InverseRange",fpi))
1523 { res = false; }
1524
1525 if(!RetrieveTable(part,theLambdaTable,ascii,directory,"Lambda",true))
1526 { res = false; }
1527 }
1528 }
1529 return res;
1530}
1531
1532//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1533
1535 G4PhysicsTable* aTable, G4bool ascii,
1536 const G4String& directory,
1537 const G4String& tname)
1538{
1539 G4bool res = true;
1540 if (nullptr != aTable) {
1541 const G4String& name = GetPhysicsTableFileName(part, directory, tname, ascii);
1542 if ( aTable->StorePhysicsTable(name,ascii) ) {
1543 if (0 < verboseLevel) G4cout << "Stored: " << name << G4endl;
1544 } else {
1545 res = false;
1546 G4cout << "Fail to store: " << name << G4endl;
1547 }
1548 }
1549 return res;
1550}
1551
1552//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.....
1553
1554G4bool
1556 G4PhysicsTable* aTable,
1557 G4bool ascii,
1558 const G4String& directory,
1559 const G4String& tname,
1560 G4bool mandatory)
1561{
1562 G4bool isRetrieved = false;
1563 G4String filename = GetPhysicsTableFileName(part,directory,tname,ascii);
1564 if(nullptr != aTable) {
1565 if(aTable->ExistPhysicsTable(filename)) {
1566 if(G4PhysicsTableHelper::RetrievePhysicsTable(aTable,filename,ascii,spline)) {
1567 isRetrieved = true;
1568 if(spline) {
1569 for(auto & v : *aTable) {
1570 if(nullptr != v) { v->FillSecondDerivatives(); }
1571 }
1572 }
1573 if (0 < verboseLevel) {
1574 G4cout << tname << " table for " << part->GetParticleName()
1575 << " is Retrieved from <" << filename << ">"
1576 << G4endl;
1577 }
1578 }
1579 }
1580 }
1581 if(mandatory && !isRetrieved) {
1582 if(0 < verboseLevel) {
1583 G4cout << tname << " table for " << part->GetParticleName()
1584 << " from file <"
1585 << filename << "> is not Retrieved"
1586 << G4endl;
1587 }
1588 return false;
1589 }
1590 return true;
1591}
1592
1593//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1594
1596 const G4MaterialCutsCouple *couple,
1597 const G4DynamicParticle* dp,
1598 G4double length)
1599{
1600 DefineMaterial(couple);
1601 G4double ekin = dp->GetKineticEnergy();
1602 SelectModel(ekin*massRatio);
1605 G4double d = 0.0;
1607 if(nullptr != fm) { d = fm->Dispersion(currentMaterial,dp,tcut,tmax,length); }
1608 return d;
1609}
1610
1611//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1612
1615 const G4MaterialCutsCouple* couple,
1616 G4double logKineticEnergy)
1617{
1618 // Cross section per volume is calculated
1619 DefineMaterial(couple);
1620 G4double cross = 0.0;
1621 if (nullptr != theLambdaTable) {
1622 cross = GetLambdaForScaledEnergy(kineticEnergy * massRatio,
1623 logKineticEnergy + logMassRatio);
1624 } else {
1625 SelectModel(kineticEnergy*massRatio);
1629 }
1630 return std::max(cross, 0.0);
1631}
1632
1633//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1634
1636{
1638 const G4double kinEnergy = track.GetKineticEnergy();
1639 const G4double logKinEnergy = track.GetDynamicParticle()->GetLogKineticEnergy();
1640 const G4double cs = GetLambdaForScaledEnergy(kinEnergy * massRatio,
1641 logKinEnergy + logMassRatio);
1642 return (0.0 < cs) ? 1.0/cs : DBL_MAX;
1643}
1644
1645//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1646
1648 G4double x, G4double y,
1649 G4double& z)
1650{
1652}
1653
1654//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1655
1657 const G4Track& track,
1658 G4double,
1660
1661{
1663 return MeanFreePath(track);
1664}
1665
1666//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1667
1669 const G4Track&,
1671{
1672 return DBL_MAX;
1673}
1674
1675//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1676
1679 G4double)
1680{
1681 DefineMaterial(couple);
1682 G4PhysicsVector* v = (*theLambdaTable)[basedCoupleIndex];
1683 return new G4PhysicsVector(*v);
1684}
1685
1686//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1687
1688void
1690{
1691 if(fTotal == tType) {
1693
1694 } else if(fRestricted == tType) {
1695 /*
1696 G4cout<< "G4VEnergyLossProcess::SetDEDXTable "
1697 << particle->GetParticleName()
1698 << " oldTable " << theDEDXTable << " newTable " << p
1699 << " ion " << theIonisationTable
1700 << " IsMaster " << isMaster
1701 << " " << GetProcessName() << G4endl;
1702 G4cout << (*p) << G4endl;
1703 */
1704 theDEDXTable = p;
1705 } else if(fIsIonisation == tType) {
1706 /*
1707 G4cout<< "G4VEnergyLossProcess::SetIonisationTable "
1708 << particle->GetParticleName()
1709 << " oldTable " << theDEDXTable << " newTable " << p
1710 << " ion " << theIonisationTable
1711 << " IsMaster " << isMaster
1712 << " " << GetProcessName() << G4endl;
1713 */
1715 }
1716}
1717
1718//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1719
1721{
1722 theCSDARangeTable = p;
1723 if(1 < verboseLevel) {
1724 G4cout << "### Set CSDA Range table " << p
1725 << " for " << particle->GetParticleName()
1726 << " and process " << GetProcessName() << G4endl;
1727 }
1728}
1729
1730//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1731
1733{
1735 if(1 < verboseLevel) {
1736 G4cout << "### Set Range table " << p
1737 << " for " << particle->GetParticleName()
1738 << " and process " << GetProcessName() << G4endl;
1739 }
1740}
1741
1742//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1743
1745{
1747 if(1 < verboseLevel) {
1748 G4cout << "### Set SecondaryRange table " << p
1749 << " for " << particle->GetParticleName()
1750 << " and process " << GetProcessName() << G4endl;
1751 }
1752}
1753
1754//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1755
1757{
1759 if(1 < verboseLevel) {
1760 G4cout << "### Set InverseRange table " << p
1761 << " for " << particle->GetParticleName()
1762 << " and process " << GetProcessName() << G4endl;
1763 }
1764}
1765
1766//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1767
1769{
1770 if(1 < verboseLevel) {
1771 G4cout << "### Set Lambda table " << p
1772 << " for " << particle->GetParticleName()
1773 << " and process " << GetProcessName() << G4endl;
1774 //G4cout << *p << G4endl;
1775 }
1776 theLambdaTable = p;
1777 tablesAreBuilt = true;
1778
1782
1783 if(isMaster && nullptr == baseParticle &&
1784 nullptr != theLambdaTable && fEmTwoPeaks == fXSType) {
1785
1786 size_t n = theLambdaTable->length();
1787
1788 G4double e, ss, xs, ee, e1peak, xs1peak, e1deep, e2peak, e2deep, xs2peak;
1789
1790 // first loop on existing vectors
1791 for (size_t i=0; i<n; ++i) {
1792 const G4PhysicsVector* pv = (*theLambdaTable)[i];
1793 ee = xs = xs1peak = xs2peak = 0.0;
1794 e1peak = e1deep = e2peak = e2deep = DBL_MAX;
1795 if(nullptr != pv) {
1796 size_t nb = pv->GetVectorLength();
1797 for (size_t j=0; j<nb; ++j) {
1798 e = pv->Energy(j);
1799 ss = (*pv)(j);
1800 // find out 1st peak
1801 if(e1peak == DBL_MAX) {
1802 if(ss >= xs) {
1803 xs = ss;
1804 ee = e;
1805 continue;
1806 } else {
1807 e1peak = ee;
1808 xs1peak = xs;
1809 }
1810 }
1811 // find out the deep
1812 if(e1deep == DBL_MAX) {
1813 if(ss <= xs) {
1814 xs = ss;
1815 ee = e;
1816 continue;
1817 } else {
1818 e1deep = ee;
1819 }
1820 }
1821 // find out 2nd peak
1822 if(e2peak == DBL_MAX) {
1823 if(ss >= xs) {
1824 xs = ss;
1825 ee = e;
1826 continue;
1827 } else {
1828 e2peak = ee;
1829 xs2peak = xs;
1830 }
1831 }
1832 if(e2deep == DBL_MAX) {
1833 if(ss <= xs) {
1834 xs = ss;
1835 ee = e;
1836 continue;
1837 } else {
1838 e2deep = ee;
1839 break;
1840 }
1841 }
1842 }
1843 }
1844 G4TwoPeaksXS* x = (*fXSpeaks)[i];
1845 if(nullptr == x) {
1846 x = new G4TwoPeaksXS();
1847 (*fXSpeaks)[i] = x;
1848 }
1849 x->e1peak = e1peak;
1850 x->e1deep = e1deep;
1851 x->e2peak = e2peak;
1852 x->e2deep = e2deep;
1853
1854 if(1 < verboseLevel) {
1855 G4cout << "For " << particle->GetParticleName()
1856 << " index= " << i << " data:\n" << " E1peak=" << e1peak
1857 << " xs1= " << xs1peak << " E1deep=" << e1deep
1858 << " E2peak=" << e2peak << " xs2=" << xs2peak
1859 << " E2deep=" << e2deep << G4endl;
1860 }
1861 }
1862 // second loop using base materials
1863 for (size_t i=0; i<n; ++i) {
1864 const G4PhysicsVector* pv = (*theLambdaTable)[i];
1865 if (nullptr == pv) {
1866 G4int j = (*theDensityIdx)[i];
1867 G4TwoPeaksXS* x = (*fXSpeaks)[i];
1868 G4TwoPeaksXS* y = (*fXSpeaks)[j];
1869 if(nullptr == x) {
1870 x = new G4TwoPeaksXS();
1871 (*fXSpeaks)[i] = x;
1872 }
1873 x->e1peak = y->e1peak;
1874 x->e1deep = y->e1deep;
1875 x->e2peak = y->e2peak;
1876 x->e2deep = y->e2deep;
1877 }
1878 }
1879 }
1880}
1881
1882//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1883
1884void G4VEnergyLossProcess::SetTwoPeaksXS(std::vector<G4TwoPeaksXS*>* ptr)
1885{
1886 fXSpeaks = ptr;
1887}
1888
1889//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1890
1892{
1893 return (nullptr != currentModel) ? currentModel->GetCurrentElement() : nullptr;
1894}
1895
1896//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1897
1899 G4bool flag)
1900{
1901 if(f > 0.0) {
1902 biasFactor = f;
1903 weightFlag = flag;
1904 if(1 < verboseLevel) {
1905 G4cout << "### SetCrossSectionBiasingFactor: for "
1906 << " process " << GetProcessName()
1907 << " biasFactor= " << f << " weightFlag= " << flag
1908 << G4endl;
1909 }
1910 }
1911}
1912
1913//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1914
1916 const G4String& region,
1917 G4bool flag)
1918{
1919 if(nullptr == biasManager) { biasManager = new G4EmBiasingManager(); }
1920 if(1 < verboseLevel) {
1921 G4cout << "### ActivateForcedInteraction: for "
1922 << " process " << GetProcessName()
1923 << " length(mm)= " << length/mm
1924 << " in G4Region <" << region
1925 << "> weightFlag= " << flag
1926 << G4endl;
1927 }
1928 weightFlag = flag;
1929 biasManager->ActivateForcedInteraction(length, region);
1930}
1931
1932//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1933
1934void
1936 G4double factor,
1937 G4double energyLimit)
1938{
1939 if (0.0 <= factor) {
1940 // Range cut can be applied only for e-
1941 if(0.0 == factor && secondaryParticle != G4Electron::Electron())
1942 { return; }
1943
1944 if(nullptr == biasManager) { biasManager = new G4EmBiasingManager(); }
1945 biasManager->ActivateSecondaryBiasing(region, factor, energyLimit);
1946 if(1 < verboseLevel) {
1947 G4cout << "### ActivateSecondaryBiasing: for "
1948 << " process " << GetProcessName()
1949 << " factor= " << factor
1950 << " in G4Region <" << region
1951 << "> energyLimit(MeV)= " << energyLimit/MeV
1952 << G4endl;
1953 }
1954 }
1955}
1956
1957//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1958
1960{
1961 isIonisation = val;
1963}
1964
1965//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1966
1968{
1969 if(0.0 < val && val < 1.0) {
1970 linLossLimit = val;
1971 actLinLossLimit = true;
1972 } else { PrintWarning("SetLinearLossLimit", val); }
1973}
1974
1975//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1976
1978{
1979 if(0.0 < v1 && 0.0 < v2) {
1980 dRoverRange = std::min(1.0, v1);
1981 finalRange = std::min(v2, 1.e+50);
1982 } else {
1983 PrintWarning("SetStepFunctionV1", v1);
1984 PrintWarning("SetStepFunctionV2", v2);
1985 }
1986}
1987
1988//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1989
1991{
1992 if(1.e-18 < val && val < 1.e+50) { lowestKinEnergy = val; }
1993 else { PrintWarning("SetLowestEnergyLimit", val); }
1994}
1995
1996//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
1997
1999{
2000 if(2 < n && n < 1000000000) {
2001 nBins = n;
2002 actBinning = true;
2003 } else {
2004 G4double e = (G4double)n;
2005 PrintWarning("SetDEDXBinning", e);
2006 }
2007}
2008
2009//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2010
2012{
2013 if(1.e-18 < e && e < maxKinEnergy) {
2014 minKinEnergy = e;
2015 actMinKinEnergy = true;
2016 } else { PrintWarning("SetMinKinEnergy", e); }
2017}
2018
2019//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2020
2022{
2023 if(minKinEnergy < e && e < 1.e+50) {
2024 maxKinEnergy = e;
2025 actMaxKinEnergy = true;
2026 if(e < maxKinEnergyCSDA) { maxKinEnergyCSDA = e; }
2027 } else { PrintWarning("SetMaxKinEnergy", e); }
2028}
2029
2030//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2031
2033{
2034 G4String ss = "G4VEnergyLossProcess::" + tit;
2036 ed << "Parameter is out of range: " << val
2037 << " it will have no effect!\n" << " Process "
2038 << GetProcessName() << " nbins= " << nBins
2039 << " Emin(keV)= " << minKinEnergy/keV
2040 << " Emax(GeV)= " << maxKinEnergy/GeV;
2041 G4Exception(ss, "em0044", JustWarning, ed);
2042}
2043
2044//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
2045
2046void G4VEnergyLossProcess::ProcessDescription(std::ostream& out) const
2047{
2048 if(nullptr != particle) { StreamInfo(out, *particle, true); }
2049}
2050
2051//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
static const G4double e1[44]
static const G4double emax
static const G4double reg
@ fBremsstrahlung
@ fPairProdByCharged
@ _SplitBremsstrahlung
G4EmTableType
@ fTotal
@ fRestricted
@ fIsIonisation
@ fEmNoIntegral
@ fEmTwoPeaks
@ fEmIncreasing
G4double condition(const G4ErrorSymMatrix &m)
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4ForceCondition
@ NotForced
G4GPILSelection
@ CandidateForSelection
@ NotCandidateForSelection
G4double G4Log(G4double x)
Definition: G4Log.hh:226
G4ProcessType
static constexpr double mm
Definition: G4SIunits.hh:95
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
static constexpr double cm
Definition: G4SIunits.hh:99
#define G4BestUnit(a, b)
@ fAlive
@ fStopAndKill
@ fStopButAlive
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4bool ForcedInteractionRegion(G4int coupleIdx)
G4double ApplySecondaryBiasing(std::vector< G4DynamicParticle * > &, const G4Track &track, G4VEmModel *currentModel, G4ParticleChangeForGamma *pParticleChange, G4double &eloss, G4int coupleIdx, G4double tcut, G4double safety=0.0)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
void Initialise(const G4ParticleDefinition &part, const G4String &procName, G4int verbose)
G4double GetWeight(G4int i)
void ActivateForcedInteraction(G4double length=0.0, const G4String &r="")
G4bool SecondaryBiasingRegion(G4int coupleIdx)
G4double GetStepLimit(G4int coupleIdx, G4double previousStep)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fm, const G4Region *r)
void DumpModelList(std::ostream &out, G4int verb)
G4int NumberOfModels() const
const G4DataVector * Initialise(const G4ParticleDefinition *part, const G4ParticleDefinition *secPart, G4double, G4int verb)
void FillLambdaVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4bool startFromNull=true, G4EmTableType t=fRestricted)
void FillDEDXVector(G4PhysicsVector *, const G4MaterialCutsCouple *, G4EmTableType t=fRestricted)
G4VEmModel * GetModel(G4int idx, G4bool ver=false)
G4bool IsPrintLocked() const
void DefineRegParamForLoss(G4VEnergyLossProcess *) const
void FillStepFunction(const G4ParticleDefinition *, G4VEnergyLossProcess *) const
static G4EmParameters * Instance()
G4int NumberOfBins() const
G4double MinKinEnergy() const
G4int NumberOfBinsPerDecade() const
G4bool BuildCSDARange() const
G4bool LossFluctuation() const
G4bool UseCutAsFinalRange() const
G4int Verbose() const
G4int WorkerVerbose() const
G4double MaxKinEnergy() const
G4bool Integral() const
G4double MaxEnergyForCSDARange() const
G4bool UseAngularGeneratorForIonisation() const
G4double LinearLossLimit() const
G4double LowestMuHadEnergy() const
G4double LambdaFactor() const
G4double LowestElectronEnergy() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
G4Region * GetRegion() const
const std::vector< G4double > * GetDensityFactors() const
const std::vector< G4int > * GetCoupleIndexes() const
void InitialiseBaseMaterials(const G4PhysicsTable *table=nullptr)
G4bool GetFlag(size_t idx)
static G4LossTableManager * Instance()
void LocalPhysicsTables(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
G4LossTableBuilder * GetTableBuilder()
void PreparePhysicsTable(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p, G4bool theMaster)
void BuildPhysicsTable(const G4ParticleDefinition *aParticle)
G4VSubCutProducer * SubCutProducer()
void DeRegister(G4VEnergyLossProcess *p)
void Register(G4VEnergyLossProcess *p)
G4VAtomDeexcitation * AtomDeexcitation()
void RegisterExtraParticle(const G4ParticleDefinition *aParticle, G4VEnergyLossProcess *p)
const G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
const G4String & GetName() const
Definition: G4Material.hh:173
void InitializeForPostStep(const G4Track &)
void InitializeForAlongStep(const G4Track &)
G4double GetProposedKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void SetProposedCharge(G4double theCharge)
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleType() const
G4double GetPDGCharge() const
const G4String & GetParticleName() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
static G4bool RetrievePhysicsTable(G4PhysicsTable *physTable, const G4String &fileName, G4bool ascii, G4bool spline)
static void SetPhysicsVector(G4PhysicsTable *physTable, std::size_t idx, G4PhysicsVector *vec)
void clearAndDestroy()
G4bool GetFlag(std::size_t i) const
G4bool ExistPhysicsTable(const G4String &fileName) const
G4bool StorePhysicsTable(const G4String &filename, G4bool ascii=false)
std::size_t length() const
G4double Energy(const std::size_t index) const
std::size_t GetVectorLength() const
void FillSecondDerivatives(const G4SplineType=G4SplineType::Base, const G4double dir1=0.0, const G4double dir2=0.0)
static G4Positron * Positron()
Definition: G4Positron.cc:93
G4ProcessVector * GetAlongStepProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
G4ProcessVector * GetAtRestProcessVector(G4ProcessVectorTypeIndex typ=typeGPIL) const
std::size_t size() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4double GetProductionCut(G4int index) const
const G4String & GetName() const
void InitialiseHelper()
G4double GetSafety() const
Definition: G4Step.hh:62
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
const G4ParticleDefinition * GetParticleDefinition() const
G4VPhysicalVolume * GetVolume() const
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4Material * GetMaterial() const
G4ParticleDefinition * GetDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
G4double GetKineticEnergy() const
void SetCreatorModelID(const G4int id)
G4int GetParentID() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
static G4TransportationManager * GetTransportationManager()
G4SafetyHelper * GetSafetyHelper() const
void AlongStepDeexcitation(std::vector< G4Track * > &tracks, const G4Step &step, G4double &eLoss, G4int coupleIndex)
virtual G4double Dispersion(const G4Material *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length)=0
virtual G4double SampleFluctuations(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double tcut, const G4double tmax, const G4double length, const G4double meanLoss)=0
virtual void FillNumberOfSecondaries(G4int &numberOfTriplets, G4int &numberOfRecoil)
Definition: G4VEmModel.cc:365
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:767
void SetMasterThread(G4bool val)
Definition: G4VEmModel.hh:739
G4VEmFluctuationModel * GetModelOfFluctuations()
Definition: G4VEmModel.hh:614
virtual G4double GetParticleCharge(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:391
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:237
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:505
void SetAngularGeneratorFlag(G4bool)
Definition: G4VEmModel.hh:725
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:447
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:655
G4bool IsActive(G4double kinEnergy) const
Definition: G4VEmModel.hh:795
virtual void CorrectionsAlongStep(const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
Definition: G4VEmModel.cc:399
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
void SetUseBaseMaterials(G4bool val)
Definition: G4VEmModel.hh:753
virtual void InitialiseLocal(const G4ParticleDefinition *, G4VEmModel *masterModel)
Definition: G4VEmModel.cc:204
G4double MaxSecondaryKinEnergy(const G4DynamicParticle *dynParticle)
Definition: G4VEmModel.hh:520
virtual G4double ChargeSquareRatio(const G4Track &)
Definition: G4VEmModel.cc:374
G4double GetDEDXForScaledEnergy(G4double scaledKinE)
G4double ScaledKinEnergyForLoss(G4double range)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
G4CrossSectionType fXSType
G4PhysicsTable * RangeTableForLoss() const
const G4ParticleDefinition * theGamma
G4double GetContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety) override
void SetMaxKinEnergy(G4double e)
G4ParticleChangeForLoss fParticleChange
void PreparePhysicsTable(const G4ParticleDefinition &) override
G4PhysicsTable * theDEDXTable
G4PhysicsTable * InverseRangeTable() const
G4double MeanFreePath(const G4Track &track)
const G4DataVector * theCuts
G4PhysicsTable * CSDARangeTable() const
G4PhysicsTable * theIonisationTable
virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition *, const G4ParticleDefinition *)=0
void SelectModel(G4double kinEnergy)
void SetRangeTableForLoss(G4PhysicsTable *p)
const G4ParticleDefinition * baseParticle
G4PhysicsTable * theInverseRangeTable
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *, G4double cut)
void ProcessDescription(std::ostream &outFile) const override
const G4ParticleDefinition * thePositron
G4double GetScaledRangeForScaledEnergy(G4double scaledKinE)
G4PhysicsTable * BuildDEDXTable(G4EmTableType tType=fRestricted)
G4VEmModel * GetModelByIndex(size_t idx=0, G4bool ver=false) const
const G4MaterialCutsCouple * currentCouple
const G4ParticleDefinition * theGenericIon
virtual G4double MinPrimaryEnergy(const G4ParticleDefinition *, const G4Material *, G4double cut)
void ActivateSecondaryBiasing(const G4String &region, G4double factor, G4double energyLimit)
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
const G4Element * GetCurrentElement() const
void SetDEDXBinning(G4int nbins)
void SetStepFunction(G4double v1, G4double v2)
G4VSubCutProducer * subcutProducer
void SetEmModel(G4VEmModel *, G4int index=0)
void SetCrossSectionBiasingFactor(G4double f, G4bool flag=true)
const G4Material * currentMaterial
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
void FillSecondariesAlongStep(G4double weight)
G4EmModelManager * modelManager
G4double ContinuousStepLimit(const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
std::vector< G4TwoPeaksXS * > * TwoPeaksXS() const
G4VEmFluctuationModel * fluctModel
G4EmParameters * theParameters
G4PhysicsTable * SecondaryRangeTable() const
const std::vector< G4double > * theDensityFactor
void StreamInfo(std::ostream &out, const G4ParticleDefinition &part, G4bool rst=false) const
G4PhysicsTable * theRangeTableForLoss
void SetInverseRangeTable(G4PhysicsTable *p)
G4bool RetrievePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
G4bool StoreTable(const G4ParticleDefinition *p, G4PhysicsTable *, G4bool ascii, const G4String &directory, const G4String &tname)
void ActivateForcedInteraction(G4double length, const G4String &region, G4bool flag=true)
void SetDynamicMassCharge(G4double massratio, G4double charge2ratio)
G4GPILSelection aGPILSelection
std::vector< const G4Region * > * scoffRegions
std::vector< G4TwoPeaksXS * > * fXSpeaks
G4bool IsIonisationProcess() const
G4PhysicsTable * theLambdaTable
G4double CrossSectionPerVolume(G4double kineticEnergy, const G4MaterialCutsCouple *couple)
const std::vector< G4int > * theDensityIdx
void SetDEDXTable(G4PhysicsTable *p, G4EmTableType tType)
void SetTwoPeaksXS(std::vector< G4TwoPeaksXS * > *p)
void SetSecondaryRangeTable(G4PhysicsTable *p)
G4PhysicsTable * BuildLambdaTable(G4EmTableType tType=fRestricted)
G4bool StorePhysicsTable(const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
void SetIonisation(G4bool val)
G4bool IsRegionForCubcutProcessor(const G4Track &aTrack)
void DefineMaterial(const G4MaterialCutsCouple *couple)
std::vector< G4Track * > scTracks
G4double AlongStepGetPhysicalInteractionLength(const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection) override
void SetLinearLossLimit(G4double val)
void SetLowestEnergyLimit(G4double)
G4bool RetrieveTable(const G4ParticleDefinition *p, G4PhysicsTable *, G4bool ascii, const G4String &directory, const G4String &tname, G4bool mandatory)
G4EmBiasingManager * biasManager
void SetLambdaTable(G4PhysicsTable *p)
const G4ParticleDefinition * secondaryParticle
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4VEnergyLossProcess(const G4String &name="EnergyLoss", G4ProcessType type=fElectromagnetic)
G4PhysicsTable * theSecondaryRangeTable
G4VParticleChange * AlongStepDoIt(const G4Track &, const G4Step &) override
G4VAtomDeexcitation * atomDeexcitation
G4PhysicsTable * theDEDXunRestrictedTable
std::vector< G4VEmModel * > * emModels
G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple, const G4DynamicParticle *dp, G4double length)
G4PhysicsTable * IonisationTable() const
void ActivateSubCutoff(const G4Region *region)
void PrintWarning(const G4String &, G4double val) const
G4PhysicsTable * LambdaTable() const
void SetCSDARangeTable(G4PhysicsTable *pRange)
const G4ParticleDefinition * particle
void ComputeLambdaForScaledEnergy(G4double scaledKinE, G4double logScaledKinE)
G4PhysicsTable * DEDXunRestrictedTable() const
const G4ParticleDefinition * theElectron
G4double GetLambdaForScaledEnergy(G4double scaledKinE)
G4PhysicsTable * theCSDARangeTable
G4LossTableManager * lManager
virtual void StreamProcessInfo(std::ostream &) const
G4SafetyHelper * safetyHelper
G4PhysicsTable * DEDXTable() const
G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &) override
void StartTracking(G4Track *) override
void SetMinKinEnergy(G4double e)
std::vector< G4DynamicParticle * > secParticles
G4double GetParentWeight() const
void ProposeTrackStatus(G4TrackStatus status)
void SetSecondaryWeightByProcess(G4bool)
void ProposeWeight(G4double finalWeight)
G4double GetLocalEnergyDeposit() const
void AddSecondary(G4Track *aSecondary)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4TrackStatus GetTrackStatus() const
G4LogicalVolume * GetLogicalVolume() const
G4double currentInteractionLength
Definition: G4VProcess.hh:335
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:338
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:412
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:518
G4int verboseLevel
Definition: G4VProcess.hh:356
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:331
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321
G4int GetProcessSubType() const
Definition: G4VProcess.hh:400
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:182
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
virtual void SampleSecondaries(const G4Step &step, std::vector< G4Track * > &tracks, G4double &eloss, G4double cut) const =0
static constexpr double mm
Definition: SystemOfUnits.h:96
static constexpr double proton_mass_c2
static constexpr double TeV
static constexpr double GeV
static constexpr double keV
static constexpr double MeV
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
const char * name(G4int ptype)
string pname
Definition: eplot.py:33
G4double e1peak
G4double e2deep
G4double e1deep
G4double e2peak
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62