Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4QDiscProcessMixer.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// $Id$
27//
28// ---------------- G4QDiscProcessMixer class -------------------
29// by Mikhail Kossov, Aug 2007.
30// 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//#define debug
37
40#include "G4SystemOfUnits.hh"
42
43
44// Constructor
46 const G4ParticleDefinition* particle,
47 G4ProcessType pType):
48 G4VDiscreteProcess(name, pType), DPParticle(particle)
49{
50 G4HadronicDeprecate("G4QDiscProcessMixer");
51
52#ifdef debug
53 G4cout<<"G4QDiscProcessMixer::Constructor is called processName="<<name<<G4endl;
54#endif
55 if (verboseLevel>0) G4cout<<GetProcessName()<<" process is created "<<G4endl;
56}
57
58// Destructor
60{
61 // Now the responsibility of deleting is deligated to the user, who created them
62 //for_each(theDPVector.begin(), theDPVector.end(), DeleteDiscreteProcess());
63}
64
66{
67 static const G4double maxEn = 1.E8*megaelectronvolt; // Conditional INF
68 if(!theDPVector.size()) // The first process in the DiscreteProcessVector (MaxE=INF)
69 {
70 std::pair<G4VDiscreteProcess*, G4double>* QDiscProc =
71 new std::pair<G4VDiscreteProcess*, G4double>(DP, maxEn);
72 theDPVector.push_back( QDiscProc );
73 }
74 else
75 {
76 if(ME < theDPVector[theDPVector.size()-1]->second)
77 {
78 std::pair<G4VDiscreteProcess*, G4double>* QDiscProc =
79 new std::pair<G4VDiscreteProcess*, G4double>(DP, ME);
80 theDPVector.push_back( QDiscProc );
81 }
82 else // Wrong Max Energy Order for the new process in the sequence of processes
83 {
84 // G4cerr<<"G4QDiscProcessMixer::AddDiscreteProcess:LastMaxE("<<theDPVector.size()-1
85 // <<")="<<theDPVector[theDPVector.size()-1]->second<<" <= MaxE="<<ME<<G4endl;
86 // G4Exception("G4QDiscProcessMixer::AddDiscreteProcess: Wrong Max Energy Order");
88 ed << " LastMaxE(" << theDPVector.size()-1 << ")="
89 << theDPVector[theDPVector.size()-1]->second << " <= MaxE=" << ME << G4endl;
90 G4Exception("G4QDiscProcessMixer::AddDiscreteProcess()", "HAD_CHPS_0000",
91 FatalException, ed);
92 }
93 }
94}
95
97{
98#ifdef debug
99 G4cout<<"G4QDiscProcessMixer::IsApplicable: part="<<particle.GetParticleName()<<" = "
100 <<DPParticle->GetParticleName()<<G4endl;
101#endif
102 if(particle == *DPParticle) return true;
103 return false;
104}
105
107 G4double PrevStSize,
109{
110 G4double kEn=Track.GetDynamicParticle()->GetKineticEnergy(); // Projectile kinetic energy
111 G4int maxDP=theDPVector.size();
112 if(maxDP) for(G4int i=maxDP-1; i>-1; i--) if(kEn < theDPVector[i]->second)
113 {
114#ifdef debug
115 G4cout<<"G4QDPMix::PSGetPIL:"<<i<<",kEn="<<kEn<<" < "<<theDPVector[i]->second<<", PIL="
116 <<theDPVector[i]->first->PostStepGetPhysicalInteractionLength(Track,PrevStSize,F)
117 <<G4endl;
118#endif
119 return theDPVector[i]->first->PostStepGetPhysicalInteractionLength(Track,PrevStSize,F);
120 }
121 return DBL_MAX;
122}
123
124// compilation fake class (length is calculated in PostStepGetPhysicalInteractionLength)
126{
127 return DBL_MAX;
128}
129
131 const G4Step& Step)
132{
133 G4double kEn=Track.GetDynamicParticle()->GetKineticEnergy(); // Projectile kinetic energy
134 G4int maxDP=theDPVector.size();
135 if(maxDP) for(G4int i=maxDP-1; i>-1; i--) if(kEn < theDPVector[i]->second)
136 {
137#ifdef debug
138 G4cout<<"G4QDPMix::PSDoIt: i="<<i<<",kEn="<<kEn<<" < "<<theDPVector[i]->second<<G4endl;
139#endif
140 //EnMomConservation= theDPVector[i]->first->GetEnegryMomentumConservation();
141 //nOfNeutrons = theDPVector[i]->first->GetNumberOfNeutronsInTarget();
142 return theDPVector[i]->first->PostStepDoIt(Track, Step);
143 }
144 return G4VDiscreteProcess::PostStepDoIt(Track, Step);
145}
146
147//G4LorentzVector G4QDiscProcessMixer::GetEnegryMomentumConservation()
148// {return EnMomConservation;}
149
150//G4int G4QDiscProcessMixer::GetNumberOfNeutronsInTarget()
151// {return nOfNeutrons;}
@ FatalException
G4ForceCondition
#define G4HadronicDeprecate(name)
G4ProcessType
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4bool IsApplicable(const G4ParticleDefinition &particle)
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4QDiscProcessMixer(const G4String &processName="Mixed Discrete Process", const G4ParticleDefinition *proj=G4Gamma::Gamma(), G4ProcessType pType=fHadronic)
G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
void AddDiscreteProcess(G4VDiscreteProcess *DP, G4double MaxE)
Definition: G4Step.hh:78
const G4DynamicParticle * GetDynamicParticle() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
#define DBL_MAX
Definition: templates.hh:83