Geant4-11
G4ChannelingOptrChangeCrossSection.hh
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// G4ChannelingOptrChangeCrossSection
29// Based On GB01BOptrChangeCrossSection
30// Class Description:
31// A G4VBiasingOperator concrete implementation example to
32// illustrate how to bias physics processes cross-section for
33// one particle type.
34// The G4VBiasingOperation G4BOptnChangeCrossSection is
35// selected by this operator, and is sent to each process
36// calling the operator.
37// A simple constant bias to the cross-section is applied,
38// but more sophisticated changes can be applied.
39//
40//---------------------------------------------------------------
41//
42
43#ifndef G4ChannelingOptrChangeCrossSection_hh
44#define G4ChannelingOptrChangeCrossSection_hh 1
45
46#include "G4VBiasingOperator.hh"
49#include <map>
50#include <unordered_map>
51
58};
59
61public:
62 // ------------------------------------------------------------
63 // -- Constructor: takes the name of the particle type to bias:
64 // ------------------------------------------------------------
65 G4ChannelingOptrChangeCrossSection(G4String particleToBias, G4String name = "ChannelingChangeXS");
67 virtual void StartRun();
68
69private:
71
72private:
73 virtual G4VBiasingOperation*
75 const G4BiasingProcessInterface* callingProcess);
76 virtual G4VBiasingOperation*
78 {return 0;}
79 virtual G4VBiasingOperation*
81 {return 0;}
82
83private:
85 virtual void OperationApplied( const G4BiasingProcessInterface* callingProcess,
86 G4BiasingAppliedCase biasingCase,
87 G4VBiasingOperation* occurenceOperationApplied,
88 G4double weightForOccurenceInteraction,
89 G4VBiasingOperation* finalStateOperationApplied,
90 const G4VParticleChange* particleChangeProduced );
91
92private:
97
98private:
99 std::unordered_map<std::string,G4ChannelingDensityRatio> fProcessToDensity;
100};
101
102#endif
G4BiasingAppliedCase
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *occurenceOperationApplied, G4double weightForOccurenceInteraction, G4VBiasingOperation *finalStateOperationApplied, const G4VParticleChange *particleChangeProduced)
virtual G4VBiasingOperation * ProposeFinalStateBiasingOperation(const G4Track *, const G4BiasingProcessInterface *)
virtual G4VBiasingOperation * ProposeOccurenceBiasingOperation(const G4Track *track, const G4BiasingProcessInterface *callingProcess)
virtual G4VBiasingOperation * ProposeNonPhysicsBiasingOperation(const G4Track *, const G4BiasingProcessInterface *)
std::unordered_map< std::string, G4ChannelingDensityRatio > fProcessToDensity
G4ChannelingOptrChangeCrossSection(G4String particleToBias, G4String name="ChannelingChangeXS")
std::map< const G4BiasingProcessInterface *, G4BOptnChangeCrossSection * > fChangeCrossSectionOperations
virtual void OperationApplied(const G4BiasingProcessInterface *callingProcess, G4BiasingAppliedCase biasingCase, G4VBiasingOperation *operationApplied, const G4VParticleChange *particleChangeProduced)
const char * name(G4int ptype)