Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RPGInelastic.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// Author: D. H. Wright
28// Date: 26 May 2007
29//
30
31#ifndef G4RPGInelastic_h
32#define G4RPGInelastic_h 1
33
34// Class Description:
35//
36// Base class for re-parameterized Gheisha-style models.
37
38#include <vector>
39
40#include "globals.hh"
41#include "G4FastVector.hh"
43#include "G4ReactionProduct.hh"
44#include "Randomize.hh"
45#include "G4RPGFragmentation.hh"
46#include "G4RPGTwoCluster.hh"
47#include "G4RPGTwoBody.hh"
50
51enum{ GHADLISTSIZE=256};
52
53
55 {
56 public: // with description
57
58 G4RPGInelastic(const G4String& modelName = "RPGInelastic");
59
61 { }
62
63 protected: // with description
64
65 G4double Pmltpc(G4int np, G4int nm, G4int nz, G4int n,
66 G4double b, G4double c);
67
69
71 const G4ReactionProduct& targetParticle,
72 G4ReactionProduct& leadParticle);
73
74 void SetUpPions(const G4int np, const G4int nm, const G4int nz,
76 G4int &vecLen);
77
78 // void Rotate(G4FastVector<G4ReactionProduct,256> &vec, G4int &vecLen);
79
80 void GetNormalizationConstant(const G4double availableEnergy,
81 G4double &n,
82 G4double &anpn);
83
85 G4int &vecLen,
86 const G4HadProjectile *originalIncident,
87 const G4DynamicParticle *originalTarget,
88 G4ReactionProduct &modifiedOriginal,
89 G4Nucleus &targetNucleus,
90 G4ReactionProduct &currentParticle,
91 G4ReactionProduct &targetParticle,
92 G4bool &incidentHasChanged,
93 G4bool &targetHasChanged,
94 G4bool quasiElastic);
95
97 G4int &vecLen,
98 G4ReactionProduct &currentParticle,
99 G4ReactionProduct &targetParticle,
100 G4bool &incidentHasChanged);
101
103
105
107
109
111
112 std::pair<G4int, G4double> interpolateEnergy(G4double ke) const;
113
114 G4int sampleFlat(std::vector<G4double> sigma) const;
115
117 G4int &vecLen,
118 G4ReactionProduct &currentParticle,
119 G4ReactionProduct &targetParticle,
121
122 enum {pi0, pip, pim, kp, km, k0, k0b, pro, neu,
123 lam, sp, s0, sm, xi0, xim, om, ap, an};
124
125 protected:
126
128
129 private:
130
131 G4double cache;
132 G4ThreeVector what;
133
134 static const G4double energyScale[30];
135
136 };
137
138#endif
double S(double temp)
@ GHADLISTSIZE
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4int sampleFlat(std::vector< G4double > sigma) const
std::pair< G4int, G4double > interpolateEnergy(G4double ke) const
G4RPGPionSuppression pionSuppression
void SetUpPions(const G4int np, const G4int nm, const G4int nz, G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen)
void CheckQnums(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4double Q, G4double B, G4double S)
void GetNormalizationConstant(const G4double availableEnergy, G4double &n, G4double &anpn)
virtual ~G4RPGInelastic()
void CalculateMomenta(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, const G4HadProjectile *originalIncident, const G4DynamicParticle *originalTarget, G4ReactionProduct &modifiedOriginal, G4Nucleus &targetNucleus, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged, G4bool &targetHasChanged, G4bool quasiElastic)
G4RPGTwoBody twoBody
G4RPGFragmentation fragmentation
G4RPGTwoCluster twoCluster
void SetUpChange(G4FastVector< G4ReactionProduct, 256 > &vec, G4int &vecLen, G4ReactionProduct &currentParticle, G4ReactionProduct &targetParticle, G4bool &incidentHasChanged)
G4bool MarkLeadingStrangeParticle(const G4ReactionProduct &currentParticle, const G4ReactionProduct &targetParticle, G4ReactionProduct &leadParticle)
G4ParticleDefinition * particleDef[18]
G4int Factorial(G4int n)
G4RPGStrangeProduction strangeProduction
G4double Pmltpc(G4int np, G4int nm, G4int nz, G4int n, G4double b, G4double c)