Geant4-11
Public Member Functions | Static Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
G4ChannelingOptrChangeCrossSection Class Reference

#include <G4ChannelingOptrChangeCrossSection.hh>

Inheritance diagram for G4ChannelingOptrChangeCrossSection:
G4VBiasingOperator

Public Member Functions

void AttachTo (const G4LogicalVolume *)
 
virtual void Configure ()
 
virtual void ConfigureForWorker ()
 
virtual void EndTracking ()
 
void ExitingBiasing (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
 G4ChannelingOptrChangeCrossSection (G4String particleToBias, G4String name="ChannelingChangeXS")
 
const G4String GetName () const
 
G4BiasingAppliedCase GetPreviousBiasingAppliedCase () const
 
const G4VBiasingOperationGetPreviousNonPhysicsAppliedOperation ()
 
G4VBiasingOperationGetProposedFinalStateBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
G4VBiasingOperationGetProposedNonPhysicsBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
G4VBiasingOperationGetProposedOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 
void ReportOperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 
void ReportOperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void StartRun ()
 
virtual void StartTracking (const G4Track *)
 
virtual ~G4ChannelingOptrChangeCrossSection ()
 

Static Public Member Functions

static G4VBiasingOperatorGetBiasingOperator (const G4LogicalVolume *)
 
static const std::vector< G4VBiasingOperator * > & GetBiasingOperators ()
 

Protected Member Functions

virtual void ExitBiasing (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 

Private Member Functions

virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual void OperationApplied (const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
 
virtual G4VBiasingOperationProposeFinalStateBiasingOperation (const G4Track *, const G4BiasingProcessInterface *)
 
virtual G4VBiasingOperationProposeNonPhysicsBiasingOperation (const G4Track *, const G4BiasingProcessInterface *)
 
virtual G4VBiasingOperationProposeOccurenceBiasingOperation (const G4Track *track, const G4BiasingProcessInterface *callingProcess)
 

Private Attributes

std::map< const G4BiasingProcessInterface *, G4BOptnChangeCrossSection * > fChangeCrossSectionOperations
 
G4int fChannelingID
 
std::map< const G4LogicalVolume *, G4intfDepthInTree
 
G4VBiasingOperationfFinalStateBiasingOperation
 
const G4String fName
 
G4VBiasingOperationfNonPhysicsBiasingOperation
 
G4VBiasingOperationfOccurenceBiasingOperation
 
const G4ParticleDefinitionfParticleToBias
 
const G4VBiasingOperationfPreviousAppliedFinalStateBiasingOperation
 
const G4VBiasingOperationfPreviousAppliedNonPhysicsBiasingOperation
 
const G4VBiasingOperationfPreviousAppliedOccurenceBiasingOperation
 
G4BiasingAppliedCase fPreviousBiasingAppliedCase
 
const G4VBiasingOperationfPreviousProposedFinalStateBiasingOperation
 
const G4VBiasingOperationfPreviousProposedNonPhysicsBiasingOperation
 
const G4VBiasingOperationfPreviousProposedOccurenceBiasingOperation
 
std::unordered_map< std::string, G4ChannelingDensityRatiofProcessToDensity
 
std::vector< const G4LogicalVolume * > fRootVolumes
 
G4bool fSetup
 

Static Private Attributes

static G4MapCache< const G4LogicalVolume *, G4VBiasingOperator * > fLogicalToSetupMap
 
static G4VectorCache< G4VBiasingOperator * > fOperators
 
static G4Cache< G4BiasingOperatorStateNotifier * > fStateNotifier
 

Detailed Description

Definition at line 60 of file G4ChannelingOptrChangeCrossSection.hh.

Constructor & Destructor Documentation

◆ G4ChannelingOptrChangeCrossSection()

G4ChannelingOptrChangeCrossSection::G4ChannelingOptrChangeCrossSection ( G4String  particleToBias,
G4String  name = "ChannelingChangeXS" 
)

Definition at line 44 of file G4ChannelingOptrChangeCrossSection.cc.

48fSetup(true){
50
51 if ( fParticleToBias == 0 )
52 {
54 ed << "Particle `" << particleName << "' not found !" << G4endl;
55 G4Exception("G4ChannelingOptrChangeCrossSection(...)",
56 "G4Channeling",
58 ed);
59 }
60
62}
@ 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
#define G4endl
Definition: G4ios.hh:57
std::unordered_map< std::string, G4ChannelingDensityRatio > fProcessToDensity
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
static G4int GetModelID(const G4int modelIndex)
G4VBiasingOperator(G4String name)
const char * name(G4int ptype)

References fDensityRatioNone, G4ParticleTable::FindParticle(), fParticleToBias, fProcessToDensity, G4endl, G4Exception(), G4ParticleTable::GetParticleTable(), and JustWarning.

◆ ~G4ChannelingOptrChangeCrossSection()

G4ChannelingOptrChangeCrossSection::~G4ChannelingOptrChangeCrossSection ( )
virtual

Definition at line 66 of file G4ChannelingOptrChangeCrossSection.cc.

66 {
67 for ( std::map< const G4BiasingProcessInterface*, G4BOptnChangeCrossSection* >::iterator
70 it++ ) delete (*it).second;
71}
std::map< const G4BiasingProcessInterface *, G4BOptnChangeCrossSection * > fChangeCrossSectionOperations

References fChangeCrossSectionOperations.

Member Function Documentation

◆ AttachTo()

void G4VBiasingOperator::AttachTo ( const G4LogicalVolume logical)
inherited

Definition at line 59 of file G4VBiasingOperator.cc.

60{
62 it = fLogicalToSetupMap.Find(logical);
63 if ( it == fLogicalToSetupMap.End() ) fLogicalToSetupMap[logical] = this;
64 else if ( (*it).second != this )
65 {
67 ed << "Biasing operator `" << GetName()
68 << "' can not be attached to Logical volume `"
69 << logical->GetName() << "' which is already used by another operator !" << G4endl;
70 G4Exception("G4VBiasingOperator::AttachTo(...)",
71 "BIAS.MNG.01",
73 ed);
74 }
75}
const G4String & GetName() const
iterator Find(const key_type &k)
Definition: G4Cache.hh:459
iterator End()
Definition: G4Cache.hh:453
static G4MapCache< const G4LogicalVolume *, G4VBiasingOperator * > fLogicalToSetupMap
const G4String GetName() const

References G4MapCache< KEYTYPE, VALTYPE >::End(), G4MapCache< KEYTYPE, VALTYPE >::Find(), G4VBiasingOperator::fLogicalToSetupMap, G4endl, G4Exception(), G4LogicalVolume::GetName(), G4VBiasingOperator::GetName(), and JustWarning.

◆ Configure()

virtual void G4VBiasingOperator::Configure ( )
inlinevirtualinherited

Reimplemented in G4BOptrForceCollision.

Definition at line 271 of file G4VBiasingOperator.hh.

271{}

◆ ConfigureForWorker()

virtual void G4VBiasingOperator::ConfigureForWorker ( )
inlinevirtualinherited

Reimplemented in G4BOptrForceCollision.

Definition at line 274 of file G4VBiasingOperator.hh.

274{}

◆ EndTracking()

virtual void G4VBiasingOperator::EndTracking ( )
inlinevirtualinherited

Reimplemented in G4BOptrForceCollision.

Definition at line 279 of file G4VBiasingOperator.hh.

279{}

◆ ExitBiasing()

void G4VBiasingOperator::ExitBiasing ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
protectedvirtualinherited

Reimplemented in G4BOptrForceCollision.

Definition at line 173 of file G4VBiasingOperator.cc.

174{}

Referenced by G4VBiasingOperator::ExitingBiasing().

◆ ExitingBiasing()

void G4VBiasingOperator::ExitingBiasing ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
inherited

Definition at line 152 of file G4VBiasingOperator.cc.

153{
154 ExitBiasing( track, callingProcess );
155
156 // -- reset all data members:
167}
virtual void ExitBiasing(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
const G4VBiasingOperation * fPreviousProposedNonPhysicsBiasingOperation
G4VBiasingOperation * fNonPhysicsBiasingOperation
G4BiasingAppliedCase fPreviousBiasingAppliedCase
G4VBiasingOperation * fOccurenceBiasingOperation
G4VBiasingOperation * fFinalStateBiasingOperation
const G4VBiasingOperation * fPreviousProposedOccurenceBiasingOperation
const G4VBiasingOperation * fPreviousAppliedOccurenceBiasingOperation
const G4VBiasingOperation * fPreviousAppliedNonPhysicsBiasingOperation
const G4VBiasingOperation * fPreviousProposedFinalStateBiasingOperation
const G4VBiasingOperation * fPreviousAppliedFinalStateBiasingOperation

References BAC_None, G4VBiasingOperator::ExitBiasing(), G4VBiasingOperator::fFinalStateBiasingOperation, G4VBiasingOperator::fNonPhysicsBiasingOperation, G4VBiasingOperator::fOccurenceBiasingOperation, G4VBiasingOperator::fPreviousAppliedFinalStateBiasingOperation, G4VBiasingOperator::fPreviousAppliedNonPhysicsBiasingOperation, G4VBiasingOperator::fPreviousAppliedOccurenceBiasingOperation, G4VBiasingOperator::fPreviousBiasingAppliedCase, G4VBiasingOperator::fPreviousProposedFinalStateBiasingOperation, G4VBiasingOperator::fPreviousProposedNonPhysicsBiasingOperation, and G4VBiasingOperator::fPreviousProposedOccurenceBiasingOperation.

◆ GetBiasingOperator()

G4VBiasingOperator * G4VBiasingOperator::GetBiasingOperator ( const G4LogicalVolume logical)
staticinherited

◆ GetBiasingOperators()

static const std::vector< G4VBiasingOperator * > & G4VBiasingOperator::GetBiasingOperators ( )
inlinestaticinherited

◆ GetName()

const G4String G4VBiasingOperator::GetName ( ) const
inlineinherited

◆ GetPreviousBiasingAppliedCase()

G4BiasingAppliedCase G4VBiasingOperator::GetPreviousBiasingAppliedCase ( ) const
inlineinherited

◆ GetPreviousNonPhysicsAppliedOperation()

const G4VBiasingOperation * G4VBiasingOperator::GetPreviousNonPhysicsAppliedOperation ( )
inlineinherited

◆ GetProposedFinalStateBiasingOperation()

G4VBiasingOperation * G4VBiasingOperator::GetProposedFinalStateBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
inherited

Definition at line 92 of file G4VBiasingOperator.cc.

93{
96}
virtual G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

References G4VBiasingOperator::fFinalStateBiasingOperation, and G4VBiasingOperator::ProposeFinalStateBiasingOperation().

◆ GetProposedNonPhysicsBiasingOperation()

G4VBiasingOperation * G4VBiasingOperator::GetProposedNonPhysicsBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
inherited

Definition at line 98 of file G4VBiasingOperator.cc.

99{
102}
virtual G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

References G4VBiasingOperator::fNonPhysicsBiasingOperation, and G4VBiasingOperator::ProposeNonPhysicsBiasingOperation().

◆ GetProposedOccurenceBiasingOperation()

G4VBiasingOperation * G4VBiasingOperator::GetProposedOccurenceBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
inherited

Definition at line 86 of file G4VBiasingOperator.cc.

87{
90}
virtual G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)=0

References G4VBiasingOperator::fOccurenceBiasingOperation, and G4VBiasingOperator::ProposeOccurenceBiasingOperation().

Referenced by G4ChannelingOptrMultiParticleChangeCrossSection::ProposeOccurenceBiasingOperation().

◆ OperationApplied() [1/3]

void G4VBiasingOperator::OperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation occurenceOperationApplied,
G4double  weightForOccurenceInteraction,
G4VBiasingOperation finalStateOperationApplied,
const G4VParticleChange particleChangeProduced 
)
privatevirtual

Reimplemented from G4VBiasingOperator.

Definition at line 238 of file G4VBiasingOperator.cc.

182{
183}

◆ OperationApplied() [2/3]

void G4ChannelingOptrChangeCrossSection::OperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation occurenceOperationApplied,
G4double  weightForOccurenceInteraction,
G4VBiasingOperation finalStateOperationApplied,
const G4VParticleChange particleChangeProduced 
)
privatevirtual

Reimplemented from G4VBiasingOperator.

Definition at line 251 of file G4ChannelingOptrChangeCrossSection.cc.

258{
259 G4BOptnChangeCrossSection* operation = fChangeCrossSectionOperations[callingProcess];
260 if ( operation == occurenceOperationApplied ) operation->SetInteractionOccured();
261}

References fChangeCrossSectionOperations, and G4BOptnChangeCrossSection::SetInteractionOccured().

◆ OperationApplied() [3/3]

void G4VBiasingOperator::OperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation operationApplied,
const G4VParticleChange particleChangeProduced 
)
privatevirtual

Reimplemented from G4VBiasingOperator.

Definition at line 232 of file G4VBiasingOperator.cc.

177{
178}

◆ ProposeFinalStateBiasingOperation()

virtual G4VBiasingOperation * G4ChannelingOptrChangeCrossSection::ProposeFinalStateBiasingOperation ( const G4Track ,
const G4BiasingProcessInterface  
)
inlineprivatevirtual

Implements G4VBiasingOperator.

Definition at line 77 of file G4ChannelingOptrChangeCrossSection.hh.

78 {return 0;}

◆ ProposeNonPhysicsBiasingOperation()

virtual G4VBiasingOperation * G4ChannelingOptrChangeCrossSection::ProposeNonPhysicsBiasingOperation ( const G4Track ,
const G4BiasingProcessInterface  
)
inlineprivatevirtual

Implements G4VBiasingOperator.

Definition at line 80 of file G4ChannelingOptrChangeCrossSection.hh.

81 {return 0;}

◆ ProposeOccurenceBiasingOperation()

G4VBiasingOperation * G4ChannelingOptrChangeCrossSection::ProposeOccurenceBiasingOperation ( const G4Track track,
const G4BiasingProcessInterface callingProcess 
)
privatevirtual

Implements G4VBiasingOperator.

Definition at line 173 of file G4ChannelingOptrChangeCrossSection.cc.

176{
177 if ( track->GetDefinition() != fParticleToBias ) return 0;
178
179 G4double analogInteractionLength =
181 if ( analogInteractionLength > DBL_MAX/10. ) return 0;
182
183 G4double analogXS = 1./analogInteractionLength;
184
185 G4ChannelingTrackData* trackdata =
187 if(trackdata==nullptr) return 0;
188
189 G4double XStransformation = 1.;
190 auto search = fProcessToDensity.find(callingProcess->GetWrappedProcess()->GetProcessName());
191 if(search != fProcessToDensity.end()) {
192 switch (search->second) {
194 XStransformation = trackdata->GetDensity();
195 break;
196 case fDensityRatioNuD:
197 XStransformation = trackdata->GetNuD();
198 break;
199 case fDensityRatioElD:
200 XStransformation = trackdata->GetElD();
201 break;
203 return 0;
204 break;
206 return 0;
207 break;
208 default:
209 return 0;
210 break;
211 }
212 }
213 else{
214 XStransformation = trackdata->GetDensity();
215 }
216
217 G4BOptnChangeCrossSection* operation = fChangeCrossSectionOperations[callingProcess];
218 G4VBiasingOperation* previousOperation = callingProcess->GetPreviousOccurenceBiasingOperation();
219
220 if ( previousOperation == 0 ){
221 operation->SetBiasedCrossSection( XStransformation * analogXS );
222 operation->Sample();
223 }
224 else{
225 if ( previousOperation != operation ){
227 ed << " Logic problem in operation handling !" << G4endl;
228 G4Exception("G4ChannelingOptrChangeCrossSection::ProposeOccurenceBiasingOperation(...)",
229 "G4Channeling",
231 ed);
232 return 0;
233 }
234 if ( operation->GetInteractionOccured() ){
235 operation->SetBiasedCrossSection( XStransformation * analogXS );
236 operation->Sample();
237 }
238 else{
239 operation->UpdateForStep( callingProcess->GetPreviousStepSize() );
240 operation->SetBiasedCrossSection( XStransformation * analogXS );
241 operation->UpdateForStep( 0.0 );
242 }
243 }
244
245 return operation;
246
247}
double G4double
Definition: G4Types.hh:83
void UpdateForStep(G4double stepLength)
void SetBiasedCrossSection(G4double xst, bool updateInteractionLength=false)
G4VBiasingOperation * GetPreviousOccurenceBiasingOperation() const
G4VProcess * GetWrappedProcess() const
G4VAuxiliaryTrackInformation * GetAuxiliaryTrackInformation(G4int id) const
Definition: G4Track.cc:225
G4ParticleDefinition * GetDefinition() const
G4double GetCurrentInteractionLength() const
Definition: G4VProcess.hh:443
const G4String & GetProcessName() const
Definition: G4VProcess.hh:382
#define DBL_MAX
Definition: templates.hh:62

References DBL_MAX, fChangeCrossSectionOperations, fChannelingID, fDensityRatioElD, fDensityRatioNone, fDensityRatioNotDefined, fDensityRatioNuD, fDensityRatioNuDElD, fParticleToBias, fProcessToDensity, G4endl, G4Exception(), G4Track::GetAuxiliaryTrackInformation(), G4VProcess::GetCurrentInteractionLength(), G4Track::GetDefinition(), G4ChannelingTrackData::GetDensity(), G4ChannelingTrackData::GetElD(), G4BOptnChangeCrossSection::GetInteractionOccured(), G4ChannelingTrackData::GetNuD(), G4BiasingProcessInterface::GetPreviousOccurenceBiasingOperation(), G4BiasingProcessInterface::GetPreviousStepSize(), G4VProcess::GetProcessName(), G4BiasingProcessInterface::GetWrappedProcess(), JustWarning, G4BOptnChangeCrossSection::Sample(), G4BOptnChangeCrossSection::SetBiasedCrossSection(), and G4BOptnChangeCrossSection::UpdateForStep().

◆ ReportOperationApplied() [1/2]

void G4VBiasingOperator::ReportOperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation occurenceOperationApplied,
G4double  weightForOccurenceInteraction,
G4VBiasingOperation finalStateOperationApplied,
const G4VParticleChange particleChangeProduced 
)
inherited

Definition at line 138 of file G4VBiasingOperator.cc.

144{
145 fPreviousBiasingAppliedCase = biasingCase;
146 fPreviousAppliedOccurenceBiasingOperation = occurenceOperationApplied;
147 fPreviousAppliedFinalStateBiasingOperation = finalStateOperationApplied;
148 OperationApplied( callingProcess, biasingCase, occurenceOperationApplied, weightForOccurenceInteraction, finalStateOperationApplied, particleChangeProduced );
149}
virtual void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)

References G4VBiasingOperator::fPreviousAppliedFinalStateBiasingOperation, G4VBiasingOperator::fPreviousAppliedOccurenceBiasingOperation, G4VBiasingOperator::fPreviousBiasingAppliedCase, and G4VBiasingOperator::OperationApplied().

◆ ReportOperationApplied() [2/2]

void G4VBiasingOperator::ReportOperationApplied ( const G4BiasingProcessInterface callingProcess,
G4BiasingAppliedCase  biasingCase,
G4VBiasingOperation operationApplied,
const G4VParticleChange particleChangeProduced 
)
inherited

Definition at line 104 of file G4VBiasingOperator.cc.

108{
109 fPreviousBiasingAppliedCase = biasingCase;
113 switch ( biasingCase )
114 {
115 case BAC_None:
116 break;
117 case BAC_NonPhysics:
119 break;
120 case BAC_FinalState:
122 break;
123 case BAC_Occurence:
124 G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
125 "BIAS.MNG.02",
127 "Internal logic error, please report !");
128 break;
129 default:
130 G4Exception("G4VBiasingOperator::ReportOperationApplied(...)",
131 "BIAS.MNG.03",
133 "Internal logic error, please report !");
134 }
135 OperationApplied( callingProcess, biasingCase, operationApplied, particleChangeProduced );
136}
@ BAC_Occurence
@ BAC_FinalState
@ BAC_NonPhysics

References BAC_FinalState, BAC_None, BAC_NonPhysics, BAC_Occurence, G4VBiasingOperator::fPreviousAppliedFinalStateBiasingOperation, G4VBiasingOperator::fPreviousAppliedNonPhysicsBiasingOperation, G4VBiasingOperator::fPreviousAppliedOccurenceBiasingOperation, G4VBiasingOperator::fPreviousBiasingAppliedCase, G4Exception(), JustWarning, and G4VBiasingOperator::OperationApplied().

Referenced by G4ChannelingOptrMultiParticleChangeCrossSection::OperationApplied().

◆ StartRun()

void G4ChannelingOptrChangeCrossSection::StartRun ( )
virtual

Reimplemented from G4VBiasingOperator.

Definition at line 75 of file G4ChannelingOptrChangeCrossSection.cc.

75 {
76 if ( fSetup ){
77 const G4ProcessManager* processManager = fParticleToBias->GetProcessManager();
78 const G4BiasingProcessSharedData* sharedData =
80 if ( sharedData ){
81 for ( size_t i = 0 ; i < (sharedData->GetPhysicsBiasingProcessInterfaces()).size(); i++ ){
82 const G4BiasingProcessInterface* wrapperProcess =
83 (sharedData->GetPhysicsBiasingProcessInterfaces())[i];
84 G4String processName = wrapperProcess->GetWrappedProcess()->GetProcessName();
85 G4String operationName = "channelingChangeXS-" + processName;
86 fChangeCrossSectionOperations[wrapperProcess] =
87 new G4BOptnChangeCrossSection(operationName);
88
89 G4ProcessType type = wrapperProcess->GetWrappedProcess()->GetProcessType();
90 G4int subType = wrapperProcess->GetWrappedProcess()->GetProcessSubType();
91
92 switch (type) {
93 case fNotDefined:
95 break;
96 case fTransportation:
98 break;
100 if(subType == fCoulombScattering ||
101 subType == fMultipleScattering){
102 fProcessToDensity[processName] = fDensityRatioNuD;
103 }
104 if(subType == fIonisation ||
105 subType == fPairProdByCharged ||
106 subType == fAnnihilation ||
107 subType == fAnnihilationToMuMu ||
108 subType == fAnnihilationToHadrons){
109 fProcessToDensity[processName] = fDensityRatioElD;
110 }
111 if(subType == fBremsstrahlung ||
112 subType == fNuclearStopping){
114 }
115
116 if(subType == fCerenkov ||
117 subType == fScintillation ||
118 subType == fSynchrotronRadiation ||
119 subType == fTransitionRadiation){
121 }
122 if(subType == fRayleigh ||
123 subType == fPhotoElectricEffect ||
124 subType == fComptonScattering ||
125 subType == fGammaConversion ||
126 subType == fGammaConversionToMuMu){
128 }
129 break;
130 case fOptical:
132 break;
133 case fHadronic:
134 fProcessToDensity[processName] = fDensityRatioNuD;
135 break;
137 fProcessToDensity[processName] = fDensityRatioNuD;
138 break;
139 case fGeneral:
141 break;
142 case fDecay:
144 break;
147 break;
148 case fUserDefined:
150 break;
151 case fParallel:
153 break;
154 case fPhonon:
156 break;
157 case fUCN:
159 break;
160 default:
162 break;
163 }
164 }
165 }
166 fSetup = false;
167 }
168}
@ fGammaConversionToMuMu
@ fAnnihilationToHadrons
@ fBremsstrahlung
@ fCoulombScattering
@ fGammaConversion
@ fRayleigh
@ fIonisation
@ fPairProdByCharged
@ fSynchrotronRadiation
@ fCerenkov
@ fAnnihilationToMuMu
@ fScintillation
@ fNuclearStopping
@ fComptonScattering
@ fTransitionRadiation
@ fAnnihilation
@ fMultipleScattering
@ fPhotoElectricEffect
G4ProcessType
@ fOptical
@ fPhonon
@ fParameterisation
@ fParallel
@ fUCN
@ fGeneral
@ fDecay
@ fElectromagnetic
@ fHadronic
@ fUserDefined
@ fTransportation
@ fPhotolepton_hadron
@ fNotDefined
int G4int
Definition: G4Types.hh:85
const G4BiasingProcessSharedData * GetSharedData() const
const std::vector< const G4BiasingProcessInterface * > & GetPhysicsBiasingProcessInterfaces() const
G4ProcessManager * GetProcessManager() const
G4ProcessType GetProcessType() const
Definition: G4VProcess.hh:388
G4int GetProcessSubType() const
Definition: G4VProcess.hh:400

References fAnnihilation, fAnnihilationToHadrons, fAnnihilationToMuMu, fBremsstrahlung, fCerenkov, fChangeCrossSectionOperations, fComptonScattering, fCoulombScattering, fDecay, fDensityRatioElD, fDensityRatioNone, fDensityRatioNotDefined, fDensityRatioNuD, fDensityRatioNuDElD, fElectromagnetic, fGammaConversion, fGammaConversionToMuMu, fGeneral, fHadronic, fIonisation, fMultipleScattering, fNotDefined, fNuclearStopping, fOptical, fPairProdByCharged, fParallel, fParameterisation, fParticleToBias, fPhonon, fPhotoElectricEffect, fPhotolepton_hadron, fProcessToDensity, fRayleigh, fScintillation, fSetup, fSynchrotronRadiation, fTransitionRadiation, fTransportation, fUCN, fUserDefined, G4BiasingProcessSharedData::GetPhysicsBiasingProcessInterfaces(), G4ParticleDefinition::GetProcessManager(), G4VProcess::GetProcessName(), G4VProcess::GetProcessSubType(), G4VProcess::GetProcessType(), G4BiasingProcessInterface::GetSharedData(), and G4BiasingProcessInterface::GetWrappedProcess().

◆ StartTracking()

virtual void G4VBiasingOperator::StartTracking ( const G4Track )
inlinevirtualinherited

Reimplemented in G4ChannelingOptrMultiParticleChangeCrossSection, and G4BOptrForceCollision.

Definition at line 278 of file G4VBiasingOperator.hh.

278{}

Field Documentation

◆ fChangeCrossSectionOperations

std::map< const G4BiasingProcessInterface*, G4BOptnChangeCrossSection* > G4ChannelingOptrChangeCrossSection::fChangeCrossSectionOperations
private

◆ fChannelingID

G4int G4ChannelingOptrChangeCrossSection::fChannelingID
private

◆ fDepthInTree

std::map< const G4LogicalVolume*, G4int > G4VBiasingOperator::fDepthInTree
privateinherited

Definition at line 333 of file G4VBiasingOperator.hh.

◆ fFinalStateBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fFinalStateBiasingOperation
privateinherited

◆ fLogicalToSetupMap

G4MapCache< const G4LogicalVolume *, G4VBiasingOperator * > G4VBiasingOperator::fLogicalToSetupMap
staticprivateinherited

◆ fName

const G4String G4VBiasingOperator::fName
privateinherited

Definition at line 319 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::GetName().

◆ fNonPhysicsBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fNonPhysicsBiasingOperation
privateinherited

◆ fOccurenceBiasingOperation

G4VBiasingOperation* G4VBiasingOperator::fOccurenceBiasingOperation
privateinherited

◆ fOperators

G4VectorCache< G4VBiasingOperator * > G4VBiasingOperator::fOperators
staticprivateinherited

◆ fParticleToBias

const G4ParticleDefinition* G4ChannelingOptrChangeCrossSection::fParticleToBias
private

◆ fPreviousAppliedFinalStateBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedFinalStateBiasingOperation
privateinherited

◆ fPreviousAppliedNonPhysicsBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedNonPhysicsBiasingOperation
privateinherited

◆ fPreviousAppliedOccurenceBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousAppliedOccurenceBiasingOperation
privateinherited

◆ fPreviousBiasingAppliedCase

G4BiasingAppliedCase G4VBiasingOperator::fPreviousBiasingAppliedCase
privateinherited

◆ fPreviousProposedFinalStateBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedFinalStateBiasingOperation
privateinherited

Definition at line 342 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::ExitingBiasing().

◆ fPreviousProposedNonPhysicsBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedNonPhysicsBiasingOperation
privateinherited

Definition at line 343 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::ExitingBiasing().

◆ fPreviousProposedOccurenceBiasingOperation

const G4VBiasingOperation* G4VBiasingOperator::fPreviousProposedOccurenceBiasingOperation
privateinherited

Definition at line 341 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::ExitingBiasing().

◆ fProcessToDensity

std::unordered_map<std::string,G4ChannelingDensityRatio> G4ChannelingOptrChangeCrossSection::fProcessToDensity
private

◆ fRootVolumes

std::vector< const G4LogicalVolume* > G4VBiasingOperator::fRootVolumes
privateinherited

Definition at line 332 of file G4VBiasingOperator.hh.

◆ fSetup

G4bool G4ChannelingOptrChangeCrossSection::fSetup
private

Definition at line 95 of file G4ChannelingOptrChangeCrossSection.hh.

Referenced by StartRun().

◆ fStateNotifier

G4Cache< G4BiasingOperatorStateNotifier * > G4VBiasingOperator::fStateNotifier
staticprivateinherited

Definition at line 328 of file G4VBiasingOperator.hh.

Referenced by G4VBiasingOperator::G4VBiasingOperator().


The documentation for this class was generated from the following files: