Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QDiscProcessMixer.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// $Id$
27//
28// ---------------- G4QDiscProcessMixer header ----------------
29// by Mikhail Kossov, Aug 2007.
30// Header of G4QDiscProcessMixer class of the CHIPS Simulation Branch in GEANT4
31// -----------------------------------------------------------------------------
32// Short description: universal mixer of processes (NOT models as in GHAD!)
33// depending on the application energy region (defined by users).
34// ------------------------------------------------------------------------
35
36#ifndef G4QDiscProcessMixer_hh
37#define G4QDiscProcessMixer_hh
38
39// GEANT4 Headers
40#include "globals.hh"
41#include "G4ios.hh"
42#include "Randomize.hh"
43#include "G4VDiscreteProcess.hh"
44#include "G4Track.hh"
45#include "G4Step.hh"
47#include "G4Gamma.hh"
48#include "G4DynamicParticle.hh"
50
51#include <vector>
52
54{
55public:
56
57 // Constructor
58 G4QDiscProcessMixer(const G4String& processName = "Mixed Discrete Process",
60 G4ProcessType pType = fHadronic );
61 // Destructor
63
65
67 G4double previousStepSize,
69
70 G4double GetMeanFreePath(const G4Track& aTrack, G4double previousStepSize,
72
73 G4VParticleChange* PostStepDoIt(const G4Track& aTrack, const G4Step& aStep);
74
75 // DP is the name of the Discrete Process, MaxE is the maximum energy for this process
76 // The processes must be defined, starting from high energy processes (MaxE decreasing)
77 // The MaxE for the first (highest energy process) is always infinity (automatic)
78 // The MinE = MaxE of the lower energy process (MinE=0 for the lowest energy process)
79
81
82 //G4LorentzVector GetEnegryMomentumConservation();
83
84 //G4int GetNumberOfNeutronsInTarget();
85
86private:
87
88 // Hide assignment operator as private
89 G4QDiscProcessMixer& operator=(const G4QDiscProcessMixer &right);
90
91 // Copy constructor
93
94 // BODY
95 const G4ParticleDefinition* DPParticle; // Particle for DiscreteProc mixture
96 G4QDiscreteProcessVector theDPVector; // Vector of Discrete Processes
97 //G4LorentzVector EnMomConservation; // Residual of Energy/Momentum Cons.
98 //G4int nOfNeutrons; // #of neutrons in the target nucleus
99};
100#endif
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4ProcessType
@ fHadronic
std::vector< std::pair< G4VDiscreteProcess *, G4double > * > G4QDiscreteProcessVector
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4bool IsApplicable(const G4ParticleDefinition &particle)
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
void AddDiscreteProcess(G4VDiscreteProcess *DP, G4double MaxE)
Definition: G4Step.hh:78